ROMS
Loading...
Searching...
No Matches
tl_lmd_swfrac.F File Reference
#include "cppdefs.h"
#include "set_bounds.h"
Include dependency graph for tl_lmd_swfrac.F:

Go to the source code of this file.

Functions/Subroutines

subroutine tl_lmd_swfrac_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, zscale, z, tl_z, tl_swdk)
 

Function/Subroutine Documentation

◆ tl_lmd_swfrac_tile()

subroutine tl_lmd_swfrac_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
real(r8), intent(in) zscale,
real(r8), dimension(imins:imaxs,jmins:jmaxs), intent(in) z,
real(r8), dimension(imins:imaxs,jmins:jmaxs), intent(in) tl_z,
real(r8), dimension(imins:imaxs,jmins:jmaxs), intent(out) tl_swdk )

Definition at line 6 of file tl_lmd_swfrac.F.

10!
11!git $Id$
12!================================================== Hernan G. Arango ===
13! Copyright (c) 2002-2025 The ROMS Group !
14! Licensed under a MIT/X style license !
15! See License_ROMS.md !
16!=======================================================================
17! !
18! This routine computes the tangent linear fraction of solar !
19! shortwave flux penetrating to specified depth (times Zscale) !
20! due to exponential decay in Jerlov water type. !
21! !
22! On Input: !
23! !
24! Zscale Scale factor to apply to depth array. !
25! Z Vertical height (meters, negative) for !
26! desired solar short-wave fraction. !
27! tl_Z Tangent linear vertical height for !
28! desired solar short-wave fraction. !
29! !
30! On Output: !
31! !
32! tl_swdk Tangent linear shortwave (radiation) fractional decay. !
33! !
34! Reference: !
35! !
36! Paulson, C.A., and J.J. Simpson, 1977: Irradiance meassurements !
37! in the upper ocean, J. Phys. Oceanogr., 7, 952-956. !
38! !
39! This routine was adapted from Bill Large 1995 code. !
40! !
41!=======================================================================
42!
43 USE mod_param
44 USE mod_mixing
45 USE mod_scalars
46!
47 implicit none
48!
49! Imported variable declarations.
50!
51 integer, intent(in) :: ng, tile
52 integer, intent(in) :: LBi, UBi, LBj, UBj
53 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
54
55 real(r8), intent(in) :: Zscale
56
57 real(r8), intent(in) :: Z(IminS:ImaxS,JminS:JmaxS)
58 real(r8), intent(in) :: tl_Z(IminS:ImaxS,JminS:JmaxS)
59
60 real(r8), intent(out) :: tl_swdk(IminS:ImaxS,JminS:JmaxS)
61!
62! Local variable declarations.
63!
64 integer :: Jindex, i, j
65
66 real(r8) :: cff1, cff2
67 real(r8) :: tl_cff1, tl_cff2
68
69 real(r8), dimension(IminS:ImaxS) :: fac1, fac2, fac3
70
71# include "set_bounds.h"
72!
73!-----------------------------------------------------------------------
74! Use Paulson and Simpson (1977) two wavelength bands solar
75! absorption model.
76!-----------------------------------------------------------------------
77!
78 DO j=jstr,jend
79 DO i=istr,iend
80 jindex=int(mixing(ng)%Jwtype(i,j))
81 fac1(i)=zscale/lmd_mu1(jindex)
82 fac2(i)=zscale/lmd_mu2(jindex)
83 fac3(i)=lmd_r1(jindex)
84 END DO
85!!DIR$ VECTOR ALWAYS
86 DO i=istr,iend
87 cff1=exp(z(i,j)*fac1(i))
88 tl_cff1=fac1(i)*tl_z(i,j)*cff1
89 cff2=exp(z(i,j)*fac2(i))
90 tl_cff2=fac2(i)*tl_z(i,j)*cff2
91!^ swdk(i,j)=cff1*fac3(i)+ &
92!^ & cff2*(1.0_r8-fac3(i))
93!^
94 tl_swdk(i,j)=tl_cff1*fac3(i)+ &
95 & tl_cff2*(1.0_r8-fac3(i))
96 END DO
97 END DO
98 RETURN
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
real(r8), dimension(9) lmd_mu1
real(r8), dimension(9) lmd_mu2
real(r8), dimension(9) lmd_r1

References mod_scalars::lmd_mu1, mod_scalars::lmd_mu2, mod_scalars::lmd_r1, and mod_mixing::mixing.

Referenced by tl_pre_step3d_mod::tl_pre_step3d_tile().

Here is the caller graph for this function: