1 | #define WRF_PORT |
---|
2 | |
---|
3 | module trb_mtn_stress |
---|
4 | |
---|
5 | #ifndef WRF_PORT |
---|
6 | use cam_logfile, only : iulog |
---|
7 | #else |
---|
8 | use module_cam_support, only: iulog |
---|
9 | #endif |
---|
10 | use shr_kind_mod, only : r8 => shr_kind_r8 |
---|
11 | |
---|
12 | |
---|
13 | implicit none |
---|
14 | private |
---|
15 | save |
---|
16 | |
---|
17 | public init_tms ! Initialization |
---|
18 | public compute_tms ! Full routine |
---|
19 | |
---|
20 | ! ------------ ! |
---|
21 | ! Private data ! |
---|
22 | ! ------------ ! |
---|
23 | |
---|
24 | #if( defined WACCM_PHYS) |
---|
25 | real(r8), parameter :: z0fac = 0.1_r8 ! Factor determining z_0 from orographic standard deviation [ no unit ] |
---|
26 | real(r8), parameter :: horomin= 1._r8 ! Minimum value of subgrid orographic height for mountain stress [ m ] |
---|
27 | #else |
---|
28 | real(r8), parameter :: z0fac = 0.075_r8 ! Factor determining z_0 from orographic standard deviation [ no unit ] |
---|
29 | real(r8), parameter :: horomin= 1._r8 ! Minimum value of subgrid orographic height for mountain stress [ m ] |
---|
30 | #endif |
---|
31 | |
---|
32 | real(r8), parameter :: z0max = 100._r8 ! Maximum value of z_0 for orography [ m ] |
---|
33 | real(r8), parameter :: dv2min = 0.01_r8 ! Minimum shear squared [ m2/s2 ] |
---|
34 | real(r8) :: oroconst ! Converts from standard deviation to height [ no unit ] |
---|
35 | real(r8) :: karman ! von Karman constant |
---|
36 | real(r8) :: gravit ! Acceleration due to gravity |
---|
37 | real(r8) :: rair ! Gas constant for dry air |
---|
38 | |
---|
39 | contains |
---|
40 | |
---|
41 | !============================================================================ ! |
---|
42 | ! ! |
---|
43 | !============================================================================ ! |
---|
44 | |
---|
45 | subroutine init_tms( kind, oro_in, karman_in, gravit_in, rair_in ) |
---|
46 | |
---|
47 | integer, intent(in) :: kind |
---|
48 | real(r8), intent(in) :: oro_in, karman_in, gravit_in, rair_in |
---|
49 | |
---|
50 | if( kind .ne. r8 ) then |
---|
51 | write(iulog,*) 'KIND of reals passed to init_tms -- exiting.' |
---|
52 | stop 'compute_tms' |
---|
53 | endif |
---|
54 | |
---|
55 | oroconst = oro_in |
---|
56 | karman = karman_in |
---|
57 | gravit = gravit_in |
---|
58 | rair = rair_in |
---|
59 | |
---|
60 | return |
---|
61 | end subroutine init_tms |
---|
62 | |
---|
63 | !============================================================================ ! |
---|
64 | ! ! |
---|
65 | !============================================================================ ! |
---|
66 | |
---|
67 | subroutine compute_tms( pcols , pver , ncol , & |
---|
68 | u , v , t , pmid , exner , & |
---|
69 | zm , sgh , ksrf , taux , tauy , & |
---|
70 | landfrac ) |
---|
71 | |
---|
72 | !------------------------------------------------------------------------------ ! |
---|
73 | ! Turbulent mountain stress parameterization ! |
---|
74 | ! ! |
---|
75 | ! Returns surface drag coefficient and stress associated with subgrid mountains ! |
---|
76 | ! For points where the orographic variance is small ( including ocean ), ! |
---|
77 | ! the returned surface drag coefficient and stress is zero. ! |
---|
78 | ! ! |
---|
79 | ! Lastly arranged : Sungsu Park. Jan. 2010. ! |
---|
80 | !------------------------------------------------------------------------------ ! |
---|
81 | |
---|
82 | ! ---------------------- ! |
---|
83 | ! Input-Output Arguments ! |
---|
84 | ! ---------------------- ! |
---|
85 | |
---|
86 | integer, intent(in) :: pcols ! Number of columns dimensioned |
---|
87 | integer, intent(in) :: pver ! Number of model layers |
---|
88 | integer, intent(in) :: ncol ! Number of columns actually used |
---|
89 | |
---|
90 | real(r8), intent(in) :: u(pcols,pver) ! Layer mid-point zonal wind [ m/s ] |
---|
91 | real(r8), intent(in) :: v(pcols,pver) ! Layer mid-point meridional wind [ m/s ] |
---|
92 | real(r8), intent(in) :: t(pcols,pver) ! Layer mid-point temperature [ K ] |
---|
93 | real(r8), intent(in) :: pmid(pcols,pver) ! Layer mid-point pressure [ Pa ] |
---|
94 | real(r8), intent(in) :: exner(pcols,pver) ! Layer mid-point exner function [ no unit ] |
---|
95 | real(r8), intent(in) :: zm(pcols,pver) ! Layer mid-point height [ m ] |
---|
96 | real(r8), intent(in) :: sgh(pcols) ! Standard deviation of orography [ m ] |
---|
97 | real(r8), intent(in) :: landfrac(pcols) ! Land fraction [ fraction ] |
---|
98 | |
---|
99 | real(r8), intent(out) :: ksrf(pcols) ! Surface drag coefficient [ kg/s/m2 ] |
---|
100 | real(r8), intent(out) :: taux(pcols) ! Surface zonal wind stress [ N/m2 ] |
---|
101 | real(r8), intent(out) :: tauy(pcols) ! Surface meridional wind stress [ N/m2 ] |
---|
102 | |
---|
103 | ! --------------- ! |
---|
104 | ! Local Variables ! |
---|
105 | ! --------------- ! |
---|
106 | |
---|
107 | integer :: i ! Loop index |
---|
108 | integer :: kb, kt ! Bottom and top of source region |
---|
109 | |
---|
110 | real(r8) :: horo ! Orographic height [ m ] |
---|
111 | real(r8) :: z0oro ! Orographic z0 for momentum [ m ] |
---|
112 | real(r8) :: dv2 ! (delta v)**2 [ m2/s2 ] |
---|
113 | real(r8) :: ri ! Richardson number [ no unit ] |
---|
114 | real(r8) :: stabfri ! Instability function of Richardson number [ no unit ] |
---|
115 | real(r8) :: rho ! Density [ kg/m3 ] |
---|
116 | real(r8) :: cd ! Drag coefficient [ no unit ] |
---|
117 | real(r8) :: vmag ! Velocity magnitude [ m /s ] |
---|
118 | |
---|
119 | ! ----------------------- ! |
---|
120 | ! Main Computation Begins ! |
---|
121 | ! ----------------------- ! |
---|
122 | |
---|
123 | do i = 1, ncol |
---|
124 | |
---|
125 | ! determine subgrid orgraphic height ( mean to peak ) |
---|
126 | |
---|
127 | horo = oroconst * sgh(i) |
---|
128 | |
---|
129 | ! No mountain stress if horo is too small |
---|
130 | |
---|
131 | if( horo < horomin ) then |
---|
132 | |
---|
133 | ksrf(i) = 0._r8 |
---|
134 | taux(i) = 0._r8 |
---|
135 | tauy(i) = 0._r8 |
---|
136 | |
---|
137 | else |
---|
138 | |
---|
139 | ! Determine z0m for orography |
---|
140 | |
---|
141 | z0oro = min( z0fac * horo, z0max ) |
---|
142 | |
---|
143 | ! Calculate neutral drag coefficient |
---|
144 | |
---|
145 | cd = ( karman / log( ( zm(i,pver) + z0oro ) / z0oro) )**2 |
---|
146 | |
---|
147 | ! Calculate the Richardson number over the lowest 2 layers |
---|
148 | |
---|
149 | kt = pver - 1 |
---|
150 | kb = pver |
---|
151 | dv2 = max( ( u(i,kt) - u(i,kb) )**2 + ( v(i,kt) - v(i,kb) )**2, dv2min ) |
---|
152 | |
---|
153 | ! Modification : Below computation of Ri is wrong. Note that 'Exner' function here is |
---|
154 | ! inverse exner function. Here, exner function is not multiplied in |
---|
155 | ! the denominator. Also, we should use moist Ri not dry Ri. |
---|
156 | ! Also, this approach using the two lowest model layers can be potentially |
---|
157 | ! sensitive to the vertical resolution. |
---|
158 | ! OK. I only modified the part associated with exner function. |
---|
159 | |
---|
160 | ri = 2._r8 * gravit * ( t(i,kt) * exner(i,kt) - t(i,kb) * exner(i,kb) ) * ( zm(i,kt) - zm(i,kb) ) & |
---|
161 | / ( ( t(i,kt) * exner(i,kt) + t(i,kb) * exner(i,kb) ) * dv2 ) |
---|
162 | |
---|
163 | ! ri = 2._r8 * gravit * ( t(i,kt) * exner(i,kt) - t(i,kb) * exner(i,kb) ) * ( zm(i,kt) - zm(i,kb) ) & |
---|
164 | ! / ( ( t(i,kt) + t(i,kb) ) * dv2 ) |
---|
165 | |
---|
166 | ! Calculate the instability function and modify the neutral drag cofficient. |
---|
167 | ! We should probably follow more elegant approach like Louis et al (1982) or Bretherton and Park (2009) |
---|
168 | ! but for now we use very crude approach : just 1 for ri < 0, 0 for ri > 1, and linear ramping. |
---|
169 | |
---|
170 | stabfri = max( 0._r8, min( 1._r8, 1._r8 - ri ) ) |
---|
171 | cd = cd * stabfri |
---|
172 | |
---|
173 | ! Compute density, velocity magnitude and stress using bottom level properties |
---|
174 | |
---|
175 | rho = pmid(i,pver) / ( rair * t(i,pver) ) |
---|
176 | vmag = sqrt( u(i,pver)**2 + v(i,pver)**2 ) |
---|
177 | ksrf(i) = rho * cd * vmag * landfrac(i) |
---|
178 | taux(i) = -ksrf(i) * u(i,pver) |
---|
179 | tauy(i) = -ksrf(i) * v(i,pver) |
---|
180 | |
---|
181 | end if |
---|
182 | |
---|
183 | end do |
---|
184 | |
---|
185 | return |
---|
186 | end subroutine compute_tms |
---|
187 | |
---|
188 | end module trb_mtn_stress |
---|