CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
cires_ugwp_triggers.F90
1
3
7 contains
8!
9
11 subroutine slat_geos5_tamp_v0(im, tau_amp, xlatdeg, tau_gw)
12!=================
13! GEOS-5 & MERRA-2 lat-dependent GW-source function tau(z=Zlaunch) =rho*<u'w'>
14!=================
15 implicit none
16 integer :: im
17 real :: tau_amp, xlatdeg(im), tau_gw(im)
18 real :: latdeg, flat_gw, tem
19 integer :: i
20
21!
22! if-lat
23!
24 do i=1, im
25 latdeg = abs(xlatdeg(i))
26 if (latdeg < 15.3) then
27 tem = (latdeg-3.0) / 8.0
28 flat_gw = 0.75 * exp(-tem * tem)
29 if (flat_gw < 1.2 .and. latdeg <= 3.0) flat_gw = 0.75
30 elseif (latdeg < 31.0 .and. latdeg >= 15.3) then
31 flat_gw = 0.10
32 elseif (latdeg < 60.0 .and. latdeg >= 31.0) then
33 tem = (latdeg-60.0) / 23.0
34 flat_gw = 0.50 * exp(- tem * tem)
35 elseif (latdeg >= 60.0) then
36 tem = (latdeg-60.0) / 70.0
37 flat_gw = 0.50 * exp(- tem * tem)
38 endif
39 tau_gw(i) = tau_amp*flat_gw
40 enddo
41!
42 end subroutine slat_geos5_tamp_v0
43
45 subroutine slat_geos5_v0(im, xlatdeg, tau_gw)
46!=================
47! GEOS-5 & MERRA-2 lat-dependent GW-source function tau(z=Zlaunch) =rho*<u'w'>
48!=================
49 implicit none
50 integer :: im
51 real :: xlatdeg(im)
52 real :: tau_gw(im)
53 real :: latdeg
54 real, parameter :: tau_amp = 100.e-3
55 real :: trop_gw, flat_gw
56 integer :: i
57!
58! if-lat
59!
60 trop_gw = 0.75
61 do i=1, im
62 latdeg = xlatdeg(i)
63 if (-15.3 < latdeg .and. latdeg < 15.3) then
64 flat_gw = trop_gw*exp(-( (abs(latdeg)-3.)/8.0)**2)
65 if (flat_gw < 1.2 .and. abs(latdeg) <= 3.) flat_gw = trop_gw
66 else if (latdeg > -31. .and. latdeg <= -15.3) then
67 flat_gw = 0.10
68 else if (latdeg < 31. .and. latdeg >= 15.3) then
69 flat_gw = 0.10
70 else if (latdeg > -60. .and. latdeg <= -31.) then
71 flat_gw = 0.50*exp(-((abs(latdeg)-60.)/23.)**2)
72 else if (latdeg < 60. .and. latdeg >= 31.) then
73 flat_gw = 0.50*exp(-((abs(latdeg)-60.)/23.)**2)
74 else if (latdeg <= -60.) then
75 flat_gw = 0.50*exp(-((abs(latdeg)-60.)/70.)**2)
76 else if (latdeg >= 60.) then
77 flat_gw = 0.50*exp(-((abs(latdeg)-60.)/70.)**2)
78 end if
79 tau_gw(i) = tau_amp*flat_gw
80 enddo
81!
82 end subroutine slat_geos5_v0
83
85 subroutine init_nazdir_v0(naz, xaz, yaz)
86 use ugwp_common_v0 , only : pi2
87 implicit none
88 integer :: naz
89 real, dimension(naz) :: xaz, yaz
90 integer :: idir
91 real :: phic, drad
92 drad = pi2/float(naz)
93 if (naz.ne.4) then
94 do idir =1, naz
95 phic = drad*(float(idir)-1.0)
96 xaz(idir) = cos(phic)
97 yaz(idir) = sin(phic)
98 enddo
99 else
100! if (naz.eq.4) then
101 xaz(1) = 1.0 !E
102 yaz(1) = 0.0
103 xaz(2) = 0.0
104 yaz(2) = 1.0 !N
105 xaz(3) =-1.0 !W
106 yaz(3) = 0.0
107 xaz(4) = 0.0
108 yaz(4) =-1.0 !S
109 endif
110 end subroutine init_nazdir_v0
111 end module cires_ugwp_triggers
This module contains routines describing the the latitudinal shape of vertical momentum flux function...
This module contains UGWP v0 initialization schemes.