Line data Source code
1 : !===============================================================================
2 : module calc_pressure
3 :
4 : implicit none
5 :
6 : public :: init_pressure, & ! Procedure(s)
7 : calculate_thvm
8 :
9 : private ! default scope
10 :
11 : contains
12 :
13 : !=============================================================================
14 0 : subroutine init_pressure( gr, thvm, p_sfc, &
15 0 : p_in_Pa, exner, p_in_Pa_zm, exner_zm )
16 :
17 : ! Description:
18 : ! Calculates the initial pressure according to the hydrostatic
19 : ! approximation. Combining the moist equation of state and the hydrostatic
20 : ! approximation, the change of pressure with respect to height can be
21 : ! calculated based on theta_v, such that:
22 : !
23 : ! dp/dz = - p * grav / ( Rd * theta_v * exner );
24 : !
25 : ! where exner = ( p / p0 )^(Rd/Cp);
26 : !
27 : ! and where p0 is a reference pressure of 100000 Pa.
28 : !
29 : ! The integral equation is set up to integrate over p on the left-hand side
30 : ! and integrate over z on the right-hand side. The equation is:
31 : !
32 : ! INT(p1:p2) p^(Rd/Cp-1) dp
33 : ! = - p0^(Rd/Cp) * ( grav / Rd ) * INT(z1:z2) (1/thvm) dz.
34 : !
35 : ! The value of mean theta_v (thvm) is calculated at each thermodynamic grid
36 : ! level, and linear interpolation is used in the integral equation for all
37 : ! altitude in-between successive thermodynamic levels, such that:
38 : !
39 : ! thvm(z) = ( ( thvm2 - thvm1 ) / ( z2 - z1 ) ) * ( z - z1 ) + thvm1.
40 : !
41 : ! The integrals are solved, and the results for pressure can be rewritten
42 : ! in terms of exner, such that:
43 : !
44 : ! exner2 - exner1
45 : ! | - ( grav / Cp )
46 : ! | * ( ( z2 - z1 ) / ( thvm2 - thvm1 ) ) * ln( thvm2 / thvm1 );
47 : ! = | where thvm2 /= thvm1;
48 : ! |
49 : ! | - ( grav / Cp ) * ( z2 - z1 ) / thvm; where thvm2 = thvm1 (= thvm).
50 : !
51 : ! The value of pressure (exner) can be calculated using the above equation
52 : ! at all vertical levels once the value of pressure is known at one level.
53 : ! Since the surface pressure is known at the initial time, that allows
54 : ! pressure to be calculated for the rest of the vertical profile.
55 :
56 : ! References:
57 : !-----------------------------------------------------------------------
58 :
59 : use grid_class, only: &
60 : grid, & ! Type
61 : zt2zm ! Procedure(s)
62 :
63 : use constants_clubb, only: &
64 : one, & ! 1
65 : Cp, & ! Specific heat of dry air [J/(kg K)]
66 : kappa, & ! Rd/Cp [-]
67 : p0, & ! Reference pressure of 100000 Pa [Pa]
68 : grav ! Acceleration of gravity (9.81 m/s^2) [m/s^2]
69 :
70 : use clubb_precision, only: &
71 : core_rknd ! Variable(s)
72 :
73 : implicit none
74 :
75 : type (grid), target, intent(in) :: gr
76 :
77 : ! Input Variables
78 : real( kind = core_rknd ), dimension(gr%nz), intent(in) :: &
79 : thvm ! Mean theta_v (thermodynamic levels) [K]
80 :
81 : real( kind = core_rknd ), intent(in) :: &
82 : p_sfc ! Surface pressure [Pa]
83 :
84 : ! Output Variables
85 : real( kind = core_rknd ), dimension(gr%nz), intent(out) :: &
86 : p_in_Pa, & ! Pressure (thermodynamic levels) [Pa]
87 : exner, & ! Exner function (thermodynamic levels) [-]
88 : p_in_Pa_zm, & ! Pressure on momentum levels [Pa]
89 : exner_zm ! Exner function on momentum levels [-]
90 :
91 : ! Local Variables
92 : real( kind = core_rknd ), dimension(gr%nz) :: &
93 0 : thvm_zm ! Mean theta_v interpolated to momentum grid levels [K]
94 :
95 : real( kind = core_rknd ), parameter :: &
96 : g_ov_Cp = grav / Cp ! g / Cp [K/m]
97 :
98 : integer :: k ! Vertical level index
99 :
100 :
101 : ! The pressure (and exner) at the lowest momentum level is the surface
102 : ! pressure (and exner based on the surface pressure).
103 0 : p_in_Pa_zm(1) = p_sfc
104 0 : exner_zm(1) = ( p_sfc / p0 )**kappa
105 :
106 : ! Set the pressure (and exner) at the lowest thermodynamic level, which is
107 : ! below the model lower boundary, to their values at the model lower
108 : ! boundary or surface.
109 0 : p_in_Pa(1) = p_in_Pa_zm(1)
110 0 : exner(1) = exner_zm(1)
111 :
112 : ! Interpolate theta_v to momentum levels.
113 0 : thvm_zm = zt2zm( gr, thvm )
114 :
115 : ! Calculate exner at all other thermodynamic and momentum grid levels.
116 : ! exner2
117 : ! = exner1
118 : ! | ( grav / Cp )
119 : ! | * ( ( z2 - z1 ) / ( thvm2 - thvm1 ) ) * ln( thvm2 / thvm1 );
120 : ! - | where thvm2 /= thvm1;
121 : ! |
122 : ! | ( grav / Cp ) * ( z2 - z1 ) / thvm; where thvm2 = thvm1 (= thvm).
123 :
124 : ! Calculate exner at thermodynamic level 2 (first thermodynamic grid level
125 : ! that is above the lower boundary).
126 0 : if ( abs( thvm(2) - thvm_zm(1) ) > epsilon( thvm ) * thvm(2) ) then
127 :
128 : exner(2) &
129 : = exner_zm(1) &
130 0 : - g_ov_Cp * ( gr%zt(1,2) - gr%zm(1,1) ) / ( thvm(2) - thvm_zm(1) ) &
131 0 : * log( thvm(2) / thvm_zm(1) )
132 :
133 : else ! thvm(2) = thvm_zm(1)
134 :
135 0 : exner(2) = exner_zm(1) - g_ov_Cp * ( gr%zt(1,2) - gr%zm(1,1) ) / thvm(2)
136 :
137 : endif
138 :
139 : ! Calculate pressure on thermodynamic levels.
140 0 : p_in_Pa(2) = p0 * exner(2)**(one/kappa)
141 :
142 : ! Loop over all other thermodynamic vertical grid levels.
143 0 : do k = 3, gr%nz, 1
144 :
145 : ! Calculate exner on thermodynamic levels.
146 0 : if ( abs( thvm(k) - thvm(k-1) ) > epsilon( thvm ) * thvm(k) ) then
147 :
148 : exner(k) &
149 : = exner(k-1) &
150 0 : - g_ov_Cp * ( gr%zt(1,k) - gr%zt(1,k-1) ) / ( thvm(k) - thvm(k-1) ) &
151 0 : * log( thvm(k) / thvm(k-1) )
152 :
153 : else ! thvm(k+1) = thvm(k)
154 :
155 0 : exner(k) = exner(k-1) - g_ov_Cp * ( gr%zt(1,k) - gr%zt(1,k-1) ) / thvm(k)
156 :
157 : endif
158 :
159 : ! Calculate pressure on thermodynamic levels.
160 0 : p_in_Pa(k) = p0 * exner(k)**(one/kappa)
161 :
162 : enddo ! k = 2, gr%nz, 1
163 :
164 : ! Loop over all momentum grid levels.
165 0 : do k = 2, gr%nz, 1
166 :
167 : ! Calculate exner on momentum levels.
168 0 : if ( abs( thvm(k) - thvm_zm(k) ) > epsilon( thvm ) * thvm(k) ) then
169 :
170 : exner_zm(k) &
171 : = exner(k) &
172 0 : - g_ov_Cp * ( gr%zm(1,k) - gr%zt(1,k) ) / ( thvm_zm(k) - thvm(k) ) &
173 0 : * log( thvm_zm(k) / thvm(k) )
174 :
175 : else ! thvm(k) = thvm_zm(k)
176 :
177 : exner_zm(k) &
178 0 : = exner(k) - g_ov_Cp * ( gr%zm(1,k) - gr%zt(1,k) ) / thvm_zm(k)
179 :
180 : endif
181 :
182 : ! Calculate pressure on momentum levels.
183 0 : p_in_Pa_zm(k) = p0 * exner_zm(k)**(one/kappa)
184 :
185 : enddo ! k = 2, gr%nz, 1
186 :
187 :
188 0 : return
189 :
190 : end subroutine init_pressure
191 :
192 : !=============================================================================
193 352944 : subroutine calculate_thvm( nz, ngrdcol, &
194 352944 : thlm, rtm, rcm, exner, thv_ds_zt, &
195 352944 : thvm)
196 : ! Description:
197 : ! Calculates mean theta_v based on a linearized approximation to the theta_v
198 : ! equation. This equation also includes liquid water loading.
199 :
200 : ! References:
201 : !-----------------------------------------------------------------------
202 :
203 : use constants_clubb, only: &
204 : Lv, & ! Latent Heat of Vaporizaion [J/kg]
205 : Cp, & ! Specific Heat of Dry Air [J/(kg K)]
206 : ep1, & ! Rv/Rd - 1 [-]
207 : ep2 ! Rv/Rd [-]
208 :
209 : use clubb_precision, only: &
210 : core_rknd
211 :
212 : implicit none
213 :
214 : ! ---------------------------- Input Variables ----------------------------
215 : integer, intent(in) :: &
216 : nz, &
217 : ngrdcol
218 :
219 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
220 : thlm, & ! Mean theta_l (thermodynamic levels) [K]
221 : rtm, & ! Mean total water (thermodynamic levels) [kg/kg]
222 : rcm, & ! Mean cloud water (thermodynamic levels) [kg/kg]
223 : exner, & ! Exner function (thermodynamic levels) [-]
224 : thv_ds_zt ! Reference theta_v on thermodynamic levels [K]
225 :
226 : ! ---------------------------- Return Variable ----------------------------
227 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
228 : thvm ! Mean theta_v (thermodynamic levels) [K]
229 :
230 : ! ---------------------------- Local Variables ----------------------------
231 : integer :: i, k
232 :
233 : ! ---------------------------- Begin Code ----------------------------
234 :
235 : !$acc data copyin( thlm, rtm, rcm, exner, thv_ds_zt ) &
236 : !$acc copyout( thvm )
237 :
238 : ! Calculate mean theta_v
239 : !$acc parallel loop gang vector collapse(2) default(present)
240 30353184 : do k = 1, nz
241 501287184 : do i = 1, ngrdcol
242 941868000 : thvm(i,k) = thlm(i,k) + ep1 * thv_ds_zt(i,k) * rtm(i,k) &
243 1442802240 : + ( Lv / ( Cp * exner(i,k) ) - ep2 * thv_ds_zt(i,k) ) * rcm(i,k)
244 : end do
245 : end do
246 : !$acc end parallel loop
247 :
248 : !$acc end data
249 :
250 352944 : return
251 :
252 : end subroutine calculate_thvm
253 : !=============================================================================
254 :
255 : end module calc_pressure
|