Line data Source code
1 : ! $Id$
2 : !===============================================================================
3 : module adg1_adg2_3d_luhar_pdf
4 :
5 : ! Description:
6 : ! Contains code to calculate the PDF parameters using the Analytic Double
7 : ! Gaussian 1 (ADG1) PDF closure, the Analytic Double Gaussian 2 (ADG2) PDF
8 : ! closure, or the 3D Luhar PDF closure. The ADG1 PDF is the original PDF used
9 : ! by CLUBB in interactive runs.
10 :
11 : ! References:
12 : !-------------------------------------------------------------------------
13 :
14 : implicit none
15 :
16 : public :: ADG1_pdf_driver, & ! Procedure(s)
17 : ADG2_pdf_driver, &
18 : Luhar_3D_pdf_driver, &
19 : ADG1_w_closure, &
20 : calc_Luhar_params, &
21 : close_Luhar_pdf
22 :
23 : private :: ADG1_ADG2_responder_params, & ! Procedure(s)
24 : backsolve_Luhar_params, &
25 : max_cubic_root
26 :
27 : private
28 :
29 : contains
30 :
31 : !=============================================================================
32 705888 : subroutine ADG1_pdf_driver( nz, ngrdcol, & ! In
33 705888 : wm, rtm, thlm, um, vm, & ! In
34 705888 : wp2, rtp2, thlp2, up2, vp2, & ! In
35 705888 : Skw, wprtp, wpthlp, upwp, vpwp, sqrt_wp2, & ! In
36 705888 : sigma_sqd_w, beta, mixt_frac_max_mag, & ! In
37 705888 : sclrm, sclrp2, wpsclrp, l_scalar_calc, & ! In
38 705888 : w_1, w_2, & ! Out
39 705888 : rt_1, rt_2, & ! Out
40 705888 : thl_1, thl_2, & ! Out
41 705888 : u_1, u_2, v_1, v_2, & ! Out
42 705888 : varnce_w_1, varnce_w_2, & ! Out
43 705888 : varnce_rt_1, varnce_rt_2, & ! Out
44 705888 : varnce_thl_1, varnce_thl_2, & ! Out
45 705888 : varnce_u_1, varnce_u_2, & ! Out
46 705888 : varnce_v_1, varnce_v_2, & ! Out
47 705888 : mixt_frac, & ! Out
48 705888 : alpha_rt, alpha_thl, & ! Out
49 705888 : alpha_u, alpha_v, & ! Out
50 705888 : sclr_1, sclr_2, varnce_sclr_1, & ! Out
51 705888 : varnce_sclr_2, alpha_sclr ) ! Out
52 :
53 : ! Calculates the PDF parameters using the Analytic Double Gaussian 1 (ADG1)
54 : ! PDF closure.
55 :
56 : ! References:
57 : !-----------------------------------------------------------------------
58 :
59 : use constants_clubb, only: &
60 : rt_tol, & ! Constant(s)
61 : thl_tol
62 :
63 : use parameters_model, only: &
64 : sclr_dim, & ! Variable(s)
65 : sclr_tol
66 :
67 : use clubb_precision, only: &
68 : core_rknd ! Variable(s)
69 :
70 : implicit none
71 :
72 : integer, intent(in) :: &
73 : ngrdcol, & ! Number of grid columns
74 : nz ! Number of vertical level
75 :
76 : ! Input Variables
77 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
78 : wm, & ! Mean of w-wind comp. (vert. vel.) [m/s]
79 : rtm, & ! Mean of total water mixing ratio [kg/kg]
80 : thlm, & ! Mean of liquid water potential temp. [K]
81 : um, & ! Mean of eastward wind [m/s]
82 : vm, & ! Mean of northward wind [m/s]
83 : wp2, & ! Variance of w (overall) [m^2/s^2]
84 : rtp2, & ! Variance of r_t (overall) [(kg/kg)^2]
85 : thlp2, & ! Variance of th_l (overall) [K^2]
86 : up2, & ! Variance of eastward wind (overall) [m^2/s^2]
87 : vp2, & ! Variance of northward wind (overall) [m^2/s^2]
88 : Skw, & ! Skewness of w [-]
89 : wprtp, & ! Covariance of w and r_t [(kg/kg)(m/s)]
90 : wpthlp, & ! Covariance of w and th_l [K(m/s)]
91 : upwp, & ! Covariance of u and w [m^2/s^2]
92 : vpwp, & ! Covariance of v and w [m^2/s^2]
93 : sqrt_wp2, & ! Square root of variance of w [m/s]
94 : sigma_sqd_w ! Width of individual w plumes [-]
95 :
96 : real( kind = core_rknd ), intent(in) :: &
97 : beta, & ! CLUBB tunable parameter beta [-]
98 : mixt_frac_max_mag ! Maximum allowable mag. of mixt_frac [-]
99 :
100 : real( kind = core_rknd ), dimension(ngrdcol, nz, sclr_dim), intent(in) :: &
101 : sclrm, & ! Mean of passive scalar (overall) [units vary]
102 : sclrp2, & ! Variance of passive scalar (overall) [(units vary)^2]
103 : wpsclrp ! Covariance of w and passive scalar [m/s (units vary)]
104 :
105 : logical, intent(in) :: &
106 : l_scalar_calc ! Flag to perform calculations for passive scalars
107 :
108 : ! Output Variables
109 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
110 : w_1, & ! Mean of w (1st PDF component) [m/s]
111 : w_2, & ! Mean of w (2nd PDF component) [m/s]
112 : rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
113 : rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
114 : thl_1, & ! Mean of th_l (1st PDF component) [K]
115 : thl_2, & ! Mean of th_l (2nd PDF component) [K]
116 : u_1, & ! Mean of u (eastward wind) (1st PDF component) [m/s]
117 : u_2, & ! Mean of u (eastward wind) (2nd PDF component) [m/s]
118 : v_1, & ! Mean of v (northward wind) (1st PDF component) [m/s]
119 : v_2, & ! Mean of v (northward wind) (2nd PDF component) [m/s]
120 : varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
121 : varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
122 : varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
123 : varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
124 : varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
125 : varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
126 : varnce_u_1, & ! Variance of u wind (1st PDF component) [m^2/s^2]
127 : varnce_u_2, & ! Variance of u wind (2nd PDF component) [m^2/s^2]
128 : varnce_v_1, & ! Variance of v wind (1st PDF component) [m^2/s^2]
129 : varnce_v_2, & ! Variance of v wind (2nd PDF component) [m^2/s^2]
130 : mixt_frac, & ! Mixture fraction (weight of 1st PDF component) [-]
131 : alpha_thl, & ! Factor relating to normalized variance for th_l [-]
132 : alpha_rt, & ! Factor relating to normalized variance for r_t [-]
133 : alpha_u, & ! Factor relating to normalized variance for u wind [-]
134 : alpha_v ! Factor relating to normalized variance for v wind [-]
135 :
136 : real( kind = core_rknd ), dimension(ngrdcol, nz, sclr_dim), intent(out) :: &
137 : sclr_1, & ! Mean of passive scalar (1st PDF component) [units vary]
138 : sclr_2, & ! Mean of passive scalar (2nd PDF component) [units vary]
139 : varnce_sclr_1, & ! Variance of pass. sclr (1st PDF comp.) [(units vary)^2]
140 : varnce_sclr_2, & ! Variance of pass. sclr (2nd PDF comp.) [(units vary)^2]
141 : alpha_sclr ! Factor relating to normalized variance for sclr [-]
142 :
143 : ! Local Variables
144 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
145 1411776 : w_1_n, & ! Normalized mean of w (1st PDF component) [-]
146 1411776 : w_2_n ! Normalized mean of w (2nd PDF component) [-]
147 :
148 : integer :: j ! Loop index
149 :
150 : !$acc enter data create( w_1_n, w_2_n )
151 :
152 : ! Calculate the mixture fraction and the PDF component means and variances
153 : ! of w.
154 : call ADG1_w_closure( nz, ngrdcol, wm, wp2, Skw, sigma_sqd_w, & ! In
155 : sqrt_wp2, mixt_frac_max_mag, & ! In
156 : w_1, w_2, w_1_n, w_2_n, & ! Out
157 705888 : varnce_w_1, varnce_w_2, mixt_frac ) ! Out
158 :
159 : ! Calculate the PDF component means and variances of rt.
160 : call ADG1_ADG2_responder_params( nz, ngrdcol, &
161 : rtm, rtp2, wp2, sqrt_wp2, & ! In
162 : wprtp, w_1_n, w_2_n, mixt_frac, & ! In
163 : sigma_sqd_w, beta, rt_tol, & ! In
164 : rt_1, rt_2, varnce_rt_1, & ! Out
165 705888 : varnce_rt_2, alpha_rt ) ! Out
166 :
167 : ! Calculate the PDF component means and variances of thl.
168 : call ADG1_ADG2_responder_params( nz, ngrdcol, &
169 : thlm, thlp2, wp2, sqrt_wp2, & ! In
170 : wpthlp, w_1_n, w_2_n, mixt_frac, & ! In
171 : sigma_sqd_w, beta, thl_tol, & ! In
172 : thl_1, thl_2, varnce_thl_1, & ! Out
173 705888 : varnce_thl_2, alpha_thl ) ! Out
174 :
175 : ! Calculate the PDF component means and variances of u wind.
176 : call ADG1_ADG2_responder_params( nz, ngrdcol, &
177 : um, up2, wp2, sqrt_wp2, & ! In
178 : upwp, w_1_n, w_2_n, mixt_frac, & ! In
179 : sigma_sqd_w, beta, thl_tol, & ! In
180 : u_1, u_2, varnce_u_1, & ! Out
181 705888 : varnce_u_2, alpha_u ) ! Out
182 :
183 : ! Calculate the PDF component means and variances of v wind.
184 : call ADG1_ADG2_responder_params( nz, ngrdcol, &
185 : vm, vp2, wp2, sqrt_wp2, & ! In
186 : vpwp, w_1_n, w_2_n, mixt_frac, & ! In
187 : sigma_sqd_w, beta, thl_tol, & ! In
188 : v_1, v_2, varnce_v_1, & ! Out
189 705888 : varnce_v_2, alpha_v ) ! Out
190 :
191 : ! Calculate the PDF component means and variances of passive scalars.
192 705888 : if ( l_scalar_calc ) then
193 0 : do j = 1, sclr_dim
194 : call ADG1_ADG2_responder_params( nz, ngrdcol, &
195 : sclrm(:,:,j), sclrp2(:,:,j), & ! In
196 : wp2, sqrt_wp2, wpsclrp(:,:,j), & ! In
197 : w_1_n, w_2_n, mixt_frac, & ! In
198 : sigma_sqd_w, beta, & ! In
199 0 : sclr_tol(j), & ! In
200 : sclr_1(:,:,j), sclr_2(:,:,j), & ! Out
201 : varnce_sclr_1(:,:,j), & ! Out
202 : varnce_sclr_2(:,:,j), & ! Out
203 0 : alpha_sclr(:,:,j) ) ! Out
204 : enddo ! i=1, sclr_dim
205 : endif ! l_scalar_calc
206 :
207 : !$acc exit data delete( w_1_n, w_2_n )
208 :
209 705888 : return
210 :
211 : end subroutine ADG1_pdf_driver
212 :
213 : !=============================================================================
214 0 : subroutine ADG2_pdf_driver( nz, ngrdcol, & ! In
215 0 : wm, rtm, thlm, wp2, rtp2, thlp2, & ! In
216 0 : Skw, wprtp, wpthlp, sqrt_wp2, beta, & ! In
217 0 : sclrm, sclrp2, wpsclrp, l_scalar_calc, & ! In
218 0 : w_1, w_2, & ! Out
219 0 : rt_1, rt_2, & ! Out
220 0 : thl_1, thl_2, & ! Out
221 0 : varnce_w_1, varnce_w_2, & ! Out
222 0 : varnce_rt_1, varnce_rt_2, & ! Out
223 0 : varnce_thl_1, varnce_thl_2, & ! Out
224 0 : mixt_frac, & ! Out
225 0 : alpha_rt, alpha_thl, & ! Out
226 0 : sigma_sqd_w, sclr_1, sclr_2, & ! Out
227 0 : varnce_sclr_1, varnce_sclr_2, alpha_sclr ) ! Out
228 :
229 : ! Calculates the PDF parameters using the Analytic Double Gaussian 2 (ADG2)
230 : ! PDF closure.
231 :
232 : ! References:
233 : !-----------------------------------------------------------------------
234 :
235 : use grid_class, only: &
236 : grid ! Type
237 :
238 : use constants_clubb, only: &
239 : one, & ! Constant(s)
240 : w_tol_sqd, &
241 : rt_tol, &
242 : thl_tol
243 :
244 : use parameters_model, only: &
245 : sclr_dim, & ! Variable(s)
246 : sclr_tol
247 :
248 : use clubb_precision, only: &
249 : core_rknd ! Variable(s)
250 :
251 : implicit none
252 :
253 : integer, intent(in) :: &
254 : ngrdcol, & ! Number of grid columns
255 : nz ! Number of vertical level
256 :
257 : ! Input Variables
258 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
259 : wm, & ! Mean of w-wind component (vertical velocity) [m/s]
260 : rtm, & ! Mean of total water mixing ratio [kg/kg]
261 : thlm, & ! Mean of liquid water potential temperature [K]
262 : wp2, & ! Variance of w (overall) [m^2/s^2]
263 : rtp2, & ! Variance of r_t (overall) [(kg/kg)^2]
264 : thlp2, & ! Variance of th_l (overall) [K^2]
265 : Skw, & ! Skewness of w [-]
266 : wprtp, & ! Covariance of w and r_t [(kg/kg)(m/s)]
267 : wpthlp, & ! Covariance of w and th_l [K(m/s)]
268 : sqrt_wp2 ! Square root of variance of w [m/s]
269 :
270 : real ( kind = core_rknd ), intent(in) :: &
271 : beta ! CLUBB tunable parameter beta [-]
272 :
273 : real( kind = core_rknd ), dimension(ngrdcol,nz, sclr_dim), intent(in) :: &
274 : sclrm, & ! Mean of passive scalar (overall) [units vary]
275 : sclrp2, & ! Variance of passive scalar (overall) [(units vary)^2]
276 : wpsclrp ! Covariance of w and passive scalar [m/s (units vary)]
277 :
278 : logical, intent(in) :: &
279 : l_scalar_calc ! Flag to perform calculations for passive scalars
280 :
281 : ! Output Variables
282 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
283 : w_1, & ! Mean of w (1st PDF component) [m/s]
284 : w_2, & ! Mean of w (2nd PDF component) [m/s]
285 : rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
286 : rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
287 : thl_1, & ! Mean of th_l (1st PDF component) [K]
288 : thl_2, & ! Mean of th_l (2nd PDF component) [K]
289 : varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
290 : varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
291 : varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
292 : varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
293 : varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
294 : varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
295 : mixt_frac, & ! Mixture fraction (weight of 1st PDF component) [-]
296 : alpha_thl, & ! Factor relating to normalized variance for th_l [-]
297 : alpha_rt, & ! Factor relating to normalized variance for r_t [-]
298 : sigma_sqd_w ! Width of individual w plumes (ADG1 parameter) [-]
299 :
300 : real( kind = core_rknd ), dimension(ngrdcol,nz, sclr_dim), intent(out) :: &
301 : sclr_1, & ! Mean of passive scalar (1st PDF component) [units vary]
302 : sclr_2, & ! Mean of passive scalar (2nd PDF component) [units vary]
303 : varnce_sclr_1, & ! Variance of pass. sclr (1st PDF comp.) [(units vary)^2]
304 : varnce_sclr_2, & ! Variance of pass. sclr (2nd PDF comp.) [(units vary)^2]
305 : alpha_sclr ! Factor relating to normalized variance for sclr [-]
306 :
307 : ! Local Variables
308 : real( kind = core_rknd ), dimension(ngrdcol,nz) :: &
309 0 : w_1_n, & ! Normalized mean of w (1st PDF component) [-]
310 0 : w_2_n, & ! Normalized mean of w (2nd PDF component) [-]
311 0 : small_m_w, & ! Luhar's m (tunable parameter) [-]
312 0 : big_m_w, & ! Luhar's M (a function of Luhar's m) [-]
313 0 : sigma_sqd_w_1, & ! Normalized width parameter of w (1st PDF component) [-]
314 0 : sigma_sqd_w_2 ! Normalized width parameter of w (2nd PDF component) [-]
315 :
316 : integer :: j, i ! Loop index
317 :
318 :
319 : ! Calculate the mixture fraction and the PDF component means and variances
320 : ! of w.
321 :
322 : ! Reproduce ADG2_w_closure using separate functions
323 0 : do i = 1, ngrdcol
324 :
325 0 : call calc_Luhar_params( nz, Skw(i,:), wp2(i,:), & ! intent(in)
326 0 : wp2(i,:), w_tol_sqd, & ! intent(in)
327 0 : mixt_frac(i,:), big_m_w(i,:), small_m_w(i,:) ) ! intent(out)
328 :
329 0 : call close_Luhar_pdf( nz, wm(i,:), wp2(i,:), mixt_frac(i,:), & ! intent(in)
330 0 : small_m_w(i,:), wp2(i,:), w_tol_sqd, & ! intent(in)
331 0 : sigma_sqd_w_1(i,:), sigma_sqd_w_2(i,:), & ! intent(out)
332 0 : varnce_w_1(i,:), varnce_w_2(i,:), & ! intent(out)
333 0 : w_1_n(i,:), w_2_n(i,:), w_1(i,:), w_2(i,:) ) ! intent(out)
334 :
335 : ! Overwrite sigma_sqd_w for consistency with ADG1
336 0 : sigma_sqd_w(i,:) = min( one / ( one + small_m_w(i,:)**2 ), 0.99_core_rknd )
337 :
338 : end do
339 :
340 : ! Calculate the PDF component means and variances of rt.
341 : call ADG1_ADG2_responder_params( nz, ngrdcol, & ! In
342 : rtm, rtp2, wp2, sqrt_wp2, & ! In
343 : wprtp, w_1_n, w_2_n, mixt_frac, & ! In
344 : sigma_sqd_w, beta, rt_tol, & ! In
345 : rt_1, rt_2, varnce_rt_1, & ! Out
346 0 : varnce_rt_2, alpha_rt ) ! Out
347 :
348 : ! Calculate the PDF component means and variances of thl.
349 : call ADG1_ADG2_responder_params( nz, ngrdcol, & ! In
350 : thlm, thlp2, wp2, sqrt_wp2, & ! In
351 : wpthlp, w_1_n, w_2_n, mixt_frac, & ! In
352 : sigma_sqd_w, beta, thl_tol, & ! In
353 : thl_1, thl_2, varnce_thl_1, & ! Out
354 0 : varnce_thl_2, alpha_thl ) ! Out
355 :
356 : ! Calculate the PDF component means and variances of passive scalars.
357 0 : if ( l_scalar_calc ) then
358 0 : do j = 1, sclr_dim
359 : call ADG1_ADG2_responder_params( nz, ngrdcol, & ! In
360 : sclrm(:,:,j), sclrp2(:,:,j), & ! In
361 : wp2, sqrt_wp2, wpsclrp(:,:,j), & ! In
362 : w_1_n, w_2_n, mixt_frac, & ! In
363 : sigma_sqd_w, beta, & ! In
364 0 : sclr_tol(j), & ! In
365 : sclr_1(:,:,j), sclr_2(:,:,j), & ! Out
366 : varnce_sclr_1(:,:,j), & ! Out
367 : varnce_sclr_2(:,:,j), & ! Out
368 0 : alpha_sclr(:,:,j) ) ! Out
369 : enddo ! i=1, sclr_dim
370 : endif ! l_scalar_calc
371 :
372 0 : return
373 :
374 : end subroutine ADG2_pdf_driver
375 :
376 : !=============================================================================
377 0 : subroutine Luhar_3D_pdf_driver( nz, wm, rtm, thlm, wp2, rtp2, thlp2, & ! In
378 0 : Skw, Skrt, Skthl, wprtp, wpthlp, & ! In
379 0 : w_1, w_2, & ! Out
380 0 : rt_1, rt_2, & ! Out
381 0 : thl_1, thl_2, & ! Out
382 0 : varnce_w_1, varnce_w_2, & ! Out
383 0 : varnce_rt_1, varnce_rt_2, & ! Out
384 0 : varnce_thl_1, varnce_thl_2, & ! Out
385 0 : mixt_frac ) ! Out
386 :
387 : ! Calculates the PDF parameters using the 3D Luhar PDF closure.
388 :
389 : ! References:
390 : !-----------------------------------------------------------------------
391 :
392 : use constants_clubb, only: &
393 : w_tol_sqd, & ! Constant(s)
394 : rt_tol, &
395 : thl_tol
396 :
397 : use clubb_precision, only: &
398 : core_rknd ! Variable(s)
399 :
400 : implicit none
401 :
402 : integer, intent(in) :: &
403 : nz
404 :
405 : ! Input Variables
406 : real( kind = core_rknd ), dimension(nz), intent(in) :: &
407 : wm, & ! Mean of w-wind component (vertical velocity) [m/s]
408 : rtm, & ! Mean of total water mixing ratio [kg/kg]
409 : thlm, & ! Mean of liquid water potential temperature [K]
410 : wp2, & ! Variance of w (overall) [m^2/s^2]
411 : rtp2, & ! Variance of r_t (overall) [(kg/kg)^2]
412 : thlp2, & ! Variance of th_l (overall) [K^2]
413 : Skw, & ! Skewness of w [-]
414 : Skrt, & ! Skewness of rt [-]
415 : Skthl, & ! Skewness of th_l [-]
416 : wprtp, & ! Covariance of w and r_t [(kg/kg)(m/s)]
417 : wpthlp ! Covariance of w and th_l [K(m/s)]
418 :
419 : ! Output Variables
420 : real( kind = core_rknd ), dimension(nz), intent(out) :: &
421 : w_1, & ! Mean of w (1st PDF component) [m/s]
422 : w_2, & ! Mean of w (2nd PDF component) [m/s]
423 : rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
424 : rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
425 : thl_1, & ! Mean of th_l (1st PDF component) [K]
426 : thl_2, & ! Mean of th_l (2nd PDF component) [K]
427 : varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
428 : varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
429 : varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
430 : varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
431 : varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
432 : varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
433 : mixt_frac ! Mixture fraction (weight of 1st PDF component) [-]
434 :
435 : real( kind = core_rknd ), dimension(nz) :: &
436 0 : w_1_n, & ! Normalized mean of w (1st PDF component) [-]
437 0 : w_2_n, & ! Normalized mean of w (2nd PDF component) [-]
438 0 : rt_1_n, & ! Normalized mean of rt (1st PDF component) [-]
439 0 : rt_2_n, & ! Normalized mean of rt (2nd PDF component) [-]
440 0 : thl_1_n, & ! Normalized mean of thl (1st PDF component) [-]
441 0 : thl_2_n, & ! Normalized mean of thl (2nd PDF component) [-]
442 0 : small_m_w, & ! Luhar's m (tunable parameter) for w [-]
443 0 : big_m_w, & ! Luhar's M (a function of Luhar's m) for w [-]
444 0 : small_m_rt, & ! Luhar's m (tunable parameter) for rt [-]
445 0 : big_m_rt, & ! Luhar's M (a function of Luhar's m) for rt [-]
446 0 : small_m_thl, & ! Luhar's m (tunable parameter) for thl [-]
447 0 : big_m_thl, & ! Luhar's M (a function of Luhar's m) for thl [-]
448 0 : sigma_sqd_w_1, & ! Normalized width parameter of w (1st PDF comp.) [-]
449 0 : sigma_sqd_w_2, & ! Normalized width parameter of w (2nd PDF comp.) [-]
450 0 : sigma_sqd_rt_1, & ! Normalized width parameter of rt (1st PDF comp.) [-]
451 0 : sigma_sqd_rt_2, & ! Normalized width parameter of rt (2nd PDF comp.) [-]
452 0 : sigma_sqd_thl_1, & ! Normalized width parameter of thl (1st PDF comp.) [-]
453 0 : sigma_sqd_thl_2 ! Normalized width parameter of thl (2nd PDF comp.) [-]
454 :
455 : integer :: k ! Vertical loop index
456 :
457 :
458 0 : do k = 1, nz, 1
459 :
460 0 : if ( ( abs( Skw(k) ) >= abs( Skthl(k) ) ) &
461 0 : .and. ( abs( Skw(k) ) >= abs( Skrt(k) ) ) ) then
462 :
463 : ! w has the greatest magnitude of skewness.
464 :
465 : ! Solve for the w PDF
466 : call calc_Luhar_params( 1, Skw(k), wp2(k), & ! In
467 : wp2(k), w_tol_sqd, & ! In
468 0 : mixt_frac(k), big_m_w(k), small_m_w(k) ) ! Out
469 :
470 : call close_Luhar_pdf( 1, wm(k), wp2(k), mixt_frac(k), & ! In
471 : small_m_w(k), wp2(k), w_tol_sqd, & ! In
472 : sigma_sqd_w_1(k), sigma_sqd_w_2(k), & ! Out
473 : varnce_w_1(k), varnce_w_2(k), & ! Out
474 0 : w_1_n(k), w_2_n(k), w_1(k), w_2(k) ) ! Out
475 :
476 : ! Solve for the thl PDF
477 : call backsolve_Luhar_params( Skw(k), Skthl(k), & ! In
478 : big_m_w(k), mixt_frac(k), & ! In
479 : big_m_thl(k), small_m_thl(k) ) ! Out
480 :
481 : call close_Luhar_pdf( 1, thlm(k), thlp2(k), mixt_frac(k), & ! In
482 : small_m_thl(k), wpthlp(k), thl_tol**2, & ! In
483 : sigma_sqd_thl_1(k), sigma_sqd_thl_2(k), & ! Out
484 : varnce_thl_1(k), varnce_thl_2(k), & ! Out
485 : thl_1_n(k), thl_2_n(k), & ! Out
486 0 : thl_1(k), thl_2(k) ) ! Out
487 :
488 : ! Solve for the rt PDF
489 : call backsolve_Luhar_params( Skw(k), Skrt(k), & ! In
490 : big_m_w(k), mixt_frac(k), & ! In
491 : big_m_rt(k), small_m_rt(k) ) ! Out
492 :
493 : call close_Luhar_pdf( 1, rtm(k), rtp2(k), mixt_frac(k), & ! In
494 : small_m_rt(k), wprtp(k), rt_tol**2, & ! In
495 : sigma_sqd_rt_1(k), sigma_sqd_rt_2(k), & ! Out
496 : varnce_rt_1(k), varnce_rt_2(k), & ! Out
497 0 : rt_1_n(k), rt_2_n(k), rt_1(k), rt_2(k) ) ! Out
498 :
499 : elseif ( ( abs( Skthl(k) ) > abs( Skw(k) ) ) &
500 0 : .and. ( abs( Skthl(k) ) >= abs( Skrt(k) ) ) ) then
501 :
502 : ! theta-l has the greatest magnitude of skewness.
503 :
504 : ! Solve for the thl PDF
505 : call calc_Luhar_params( 1, Skthl(k), wpthlp(k), & ! In
506 : thlp2(k), thl_tol**2, & ! In
507 : mixt_frac(k), big_m_thl(k), & ! Out
508 0 : small_m_thl(k) ) ! Out
509 :
510 : ! Solve for the thl PDF
511 : call close_Luhar_pdf( 1, thlm(k), thlp2(k), mixt_frac(k), & ! In
512 : small_m_thl(k), wpthlp(k), thl_tol**2, & ! In
513 : sigma_sqd_thl_1(k), sigma_sqd_thl_2(k), & ! Out
514 : varnce_thl_1(k), varnce_thl_2(k), & ! Out
515 : thl_1_n(k), thl_2_n(k), & ! Out
516 0 : thl_1(k), thl_2(k) ) ! Out
517 :
518 : ! Solve for the w PDF
519 : call backsolve_Luhar_params( Skthl(k), Skw(k), & ! In
520 : big_m_thl(k), mixt_frac(k), & ! In
521 : big_m_w(k), small_m_w(k) ) ! Out
522 :
523 : call close_Luhar_pdf( 1, wm(k), wp2(k), mixt_frac(k), & ! In
524 : small_m_w(k), wp2(k), w_tol_sqd, & ! In
525 : sigma_sqd_w_1(k), sigma_sqd_w_2(k), & ! Out
526 : varnce_w_1(k), varnce_w_2(k), & ! Out
527 0 : w_1_n(k), w_2_n(k), w_1(k), w_2(k) ) ! Out
528 :
529 : ! Solve for the rt PDF
530 : call backsolve_Luhar_params( Skthl(k), Skrt(k), & ! In
531 : big_m_thl(k), mixt_frac(k), & ! In
532 : big_m_rt(k), small_m_rt(k) ) ! Out
533 :
534 : call close_Luhar_pdf( 1, rtm(k), rtp2(k), mixt_frac(k), & ! In
535 : small_m_rt(k), wprtp(k), rt_tol**2, & ! In
536 : sigma_sqd_rt_1(k), sigma_sqd_rt_2(k), & ! Out
537 : varnce_rt_1(k), varnce_rt_2(k), & ! Out
538 0 : rt_1_n(k), rt_2_n(k), rt_1(k), rt_2(k) ) ! Out
539 :
540 : else
541 :
542 : ! rt has the greatest magnitude of skewness.
543 :
544 : ! Solve for the rt PDF
545 : call calc_Luhar_params( 1, Skrt(k), wprtp(k), & ! In
546 : rtp2(k), rt_tol**2, & ! In
547 : mixt_frac(k), big_m_rt(k), & ! Out
548 0 : small_m_rt(k) ) ! Out
549 :
550 : ! Solve for the rt PDF
551 : call close_Luhar_pdf( 1, rtm(k), rtp2(k), mixt_frac(k), & ! In
552 : small_m_rt(k), wprtp(k), rt_tol**2, & ! In
553 : sigma_sqd_rt_1(k), sigma_sqd_rt_2(k), & ! Out
554 : varnce_rt_1(k), varnce_rt_2(k), & ! Out
555 0 : rt_1_n(k), rt_2_n(k), rt_1(k), rt_2(k) ) ! Out
556 :
557 : ! Solve for the w PDF
558 : call backsolve_Luhar_params( Skrt(k), Skw(k), & ! In
559 : big_m_rt(k), mixt_frac(k), & ! In
560 : big_m_w(k), small_m_w(k) ) ! Out
561 :
562 : call close_Luhar_pdf( 1, wm(k), wp2(k), mixt_frac(k), & ! In
563 : small_m_w(k), wp2(k), w_tol_sqd, & ! In
564 : sigma_sqd_w_1(k), sigma_sqd_w_2(k), & ! Out
565 : varnce_w_1(k), varnce_w_2(k), & ! Out
566 0 : w_1_n(k), w_2_n(k), w_1(k), w_2(k) ) ! OUt
567 :
568 : ! Solve for the thl PDF
569 : call backsolve_Luhar_params( Skrt(k), Skthl(k), & ! In
570 : big_m_rt(k), mixt_frac(k), & ! In
571 : big_m_thl(k), small_m_thl(k) ) ! Out
572 :
573 : call close_Luhar_pdf( 1, thlm(k), thlp2(k), mixt_frac(k), & ! In
574 : small_m_thl(k), wpthlp(k), thl_tol**2, & ! In
575 : sigma_sqd_thl_1(k), sigma_sqd_thl_2(k), & ! Out
576 : varnce_thl_1(k), varnce_thl_2(k), & ! Out
577 : thl_1_n(k), thl_2_n(k), & ! Out
578 0 : thl_1(k), thl_2(k) ) ! Out
579 :
580 : endif
581 :
582 : enddo ! k = 1, nz, 1
583 :
584 :
585 0 : return
586 :
587 : end subroutine Luhar_3D_pdf_driver
588 :
589 : !======================================================================
590 705888 : subroutine ADG1_w_closure( nz, ngrdcol, wm, wp2, Skw, sigma_sqd_w, & !In
591 705888 : sqrt_wp2, mixt_frac_max_mag, & !In
592 705888 : w_1, w_2, w_1_n, w_2_n, & !Out
593 705888 : varnce_w_1, varnce_w_2, mixt_frac ) !Out
594 :
595 :
596 : ! Description:
597 : ! Calculates the mixture fraction, the PDF component means of w, and the PDF
598 : ! component variances of w for the Analytic Double Gaussian 1 (ADG1)
599 : ! closure. It sets the widths of both w Gaussians to be the same, and
600 : ! furthermore, relates them to the overall variance of w (<w'^2>) by a
601 : ! parameter, sigma_sqd_w. The equation is:
602 : !
603 : ! sigma_w_1^2 = sigma_w_2^2 = sigma_sqd_w * <w'^2>;
604 : !
605 : ! where sigma_w_1^2 is the variance of w in the 1st PDF component,
606 : ! sigma_w_2^2 is the variance of w in the 2nd PDF component, and parameter
607 : ! sigma_sqd_w must have a value within the range 0 <= sigma_sqd_w <= 1.
608 : !
609 : ! References:
610 : ! Golaz, J-C., V. E. Larson, and W. R. Cotton, 2002a: A PDF-based model for
611 : ! boundary layer clouds. Part I: Method and model description. J. Atmos.
612 : ! Sci., 59, 3540–3551.
613 : !
614 : ! Vincent E. Larson and Jean-Christophe Golaz, 2005: Using Probability
615 : ! Density Functions to Derive Consistent Closure Relationships among
616 : ! Higher-Order Moments. Mon. Wea. Rev., 133, 1023–1042.
617 : !-----------------------------------------------------------------------
618 :
619 : use constants_clubb, only: &
620 : four, & ! Constant(s)
621 : one, &
622 : one_half, &
623 : zero, &
624 : w_tol_sqd
625 :
626 : use clubb_precision, only: &
627 : core_rknd ! Precision
628 :
629 : implicit none
630 :
631 : integer, intent(in) :: &
632 : ngrdcol, & ! Number of grid columns
633 : nz ! Number of vertical level
634 :
635 : ! Input Variables
636 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
637 : wm, & ! Mean of w (overall) [m/s]
638 : wp2, & ! Variance of w (overall) [m^2/s^2]
639 : Skw, & ! Skewness of w [-]
640 : sigma_sqd_w, & ! Widths of each w Gaussian [-]
641 : sqrt_wp2 ! Square root of the variance of w [m/s]
642 :
643 : real( kind = core_rknd ), intent(in) :: &
644 : mixt_frac_max_mag ! Maximum allowable mag. of mixt_frac [-]
645 :
646 : ! Output Variables
647 : real( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
648 : w_1, & ! Mean of w (1st PDF component) [m/s]
649 : w_2, & ! Mean of w (2nd PDF component) [m/s]
650 : w_1_n, & ! Normalized mean of w (1st PDF component) [-]
651 : w_2_n, & ! Normalized mean of w (2nd PDF component) [-]
652 : varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
653 : varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
654 : mixt_frac ! Mixture fraction [-]
655 :
656 : ! Local Variable
657 : integer :: k, i
658 :
659 : !----- Begin Code -----
660 :
661 : !$acc parallel loop gang vector collapse(2) default(present)
662 60706368 : do k = 1, nz
663 1002574368 : do i = 1, ngrdcol
664 :
665 1001868480 : if ( wp2(i,k) > w_tol_sqd ) then
666 :
667 : ! Width (standard deviation) parameters are non-zero
668 :
669 : ! The variable "mixt_frac" is the weight of the 1st PDF component. The
670 : ! weight of the 2nd PDF component is "1-mixt_frac". If there isn't any
671 : ! skewness of w (Sk_w = 0 because w'^3 = 0), mixt_frac = 0.5, and both
672 : ! PDF components are equally weighted. If there is positive skewness of
673 : ! w (Sk_w > 0 because w'^3 > 0), 0 < mixt_frac < 0.5, and the 2nd PDF
674 : ! component has greater weight than does the 1st PDF component. If there
675 : ! is negative skewness of w (Sk_w < 0 because w'^3 < 0),
676 : ! 0.5 < mixt_frac < 1, and the 1st PDF component has greater weight than
677 : ! does the 2nd PDF component.
678 229926100 : if ( abs( Skw(i,k) ) <= 1.0e-5_core_rknd ) then
679 6512716 : mixt_frac(i,k) = one_half
680 : else
681 : mixt_frac(i,k) &
682 : = one_half &
683 223413384 : * ( one - Skw(i,k) / sqrt( four * ( one - sigma_sqd_w(i,k) )**3 + Skw(i,k)**2 ) )
684 : endif
685 :
686 : ! Clip mixt_frac, and 1 - mixt_frac, to avoid dividing by a small number.
687 : ! Formula for mixt_frac_max_mag =
688 : ! 1 - ( 1/2 * ( 1 - Skw_max
689 : ! / sqrt( 4*( 1 - sigma_sqd_w )^3 + Skw_max^2 ) ) ),
690 : ! where sigma_sqd_w is fixed at 0.4 for this calculation.
691 : mixt_frac(i,k) = min( max( mixt_frac(i,k), one - mixt_frac_max_mag ), &
692 229926100 : mixt_frac_max_mag )
693 :
694 : ! The normalized mean of w for Gaussian "plume" 1 is w_1_n. It's value
695 : ! will always be greater than 0. As an example, a value of 1.0 would
696 : ! indicate that the actual mean of w for Gaussian "plume" 1 is found 1.0
697 : ! standard deviation above the overall mean for w.
698 : w_1_n(i,k) &
699 229926100 : = sqrt( ( ( one - mixt_frac(i,k) ) / mixt_frac(i,k) ) * ( one - sigma_sqd_w(i,k) ) )
700 : ! The normalized mean of w for Gaussian "plume" 2 is w_2_n. It's value
701 : ! will always be less than 0. As an example, a value of -0.5 would
702 : ! indicate that the actual mean of w for Gaussian "plume" 2 is found 0.5
703 : ! standard deviations below the overall mean for w.
704 : w_2_n(i,k) &
705 229926100 : = -sqrt( ( mixt_frac(i,k) / ( one - mixt_frac(i,k) ) ) * ( one - sigma_sqd_w(i,k) ) )
706 : ! The mean of w for Gaussian "plume" 1 is w_1.
707 229926100 : w_1(i,k) = wm(i,k) + sqrt_wp2(i,k) * w_1_n(i,k)
708 : ! The mean of w for Gaussian "plume" 2 is w_2.
709 229926100 : w_2(i,k) = wm(i,k) + sqrt_wp2(i,k) * w_2_n(i,k)
710 :
711 : ! The variance of w for Gaussian "plume" 1 for varnce_w_1.
712 229926100 : varnce_w_1(i,k) = sigma_sqd_w(i,k) * wp2(i,k)
713 : ! The variance of w for Gaussian "plume" 2 for varnce_w_2.
714 : ! The variance in both Gaussian "plumes" is defined to be the same.
715 229926100 : varnce_w_2(i,k) = sigma_sqd_w(i,k) * wp2(i,k)
716 :
717 : else
718 :
719 : ! Vertical velocity doesn't vary.
720 711941900 : mixt_frac(i,k) = one_half
721 711941900 : w_1_n(i,k) = sqrt( one - sigma_sqd_w(i,k) )
722 711941900 : w_2_n(i,k) = -sqrt( one - sigma_sqd_w(i,k) )
723 711941900 : w_1(i,k) = wm(i,k)
724 711941900 : w_2(i,k) = wm(i,k)
725 711941900 : varnce_w_1(i,k) = zero
726 711941900 : varnce_w_2(i,k) = zero
727 :
728 : endif ! Widths non-zero
729 : end do
730 : end do
731 : !$acc end parallel loop
732 705888 : return
733 :
734 :
735 :
736 : end subroutine ADG1_w_closure
737 :
738 :
739 : !=============================================================================
740 0 : subroutine calc_Luhar_params( nz, Skx, wpxp, & ! In
741 0 : xp2, x_tol_sqd, & ! In
742 0 : mixt_frac, big_m, small_m ) ! Out
743 :
744 : ! Description:
745 : ! For the Luhar closure, this subroutine takes Skx (and Skw) as input and
746 : ! outputs the mixture fraction, big_m, and small_m. This code was written
747 : ! using the equations and nomenclature of Larson et al. (2002) Appendix
748 : ! section e.
749 : !
750 : ! The relationship between skewness of x (Skx), mixture fraction (a), and
751 : ! Luhar's small m (m) is given by:
752 : !
753 : ! Skx^2 = ( m^2 * ( m^2 + 3 )^2 / ( m^2 + 1 )^3 )
754 : ! * ( 1 - 2*a )^2 / ( a * ( 1 - a ) ).
755 : !
756 : ! Luhar's large M (M) is used to more easily express the factor involving
757 : ! the m's:
758 : !
759 : ! M = ( m^2 + 1 )^3 / ( m^2 * ( m^2 + 3 )^2 ).
760 : !
761 : ! The equation involving skewness of x becomes:
762 : !
763 : ! Skx^2 = ( 1 / M ) * ( 1 - 2*a )^2 / ( a * ( 1 - a ) );
764 : !
765 : ! or:
766 : !
767 : ! M * Skx^2 = ( 1 - 2*a )^2 / ( a * ( 1 - a ) ).
768 : !
769 : ! This equation can be rewritten as:
770 : !
771 : ! ( a * ( 1 - a ) ) * M * Skx^2 = ( 1 - 2*a )^2;
772 : !
773 : ! as well as:
774 : !
775 : ! ( a - a^2 ) * M * Skx^2 = 1 - 4*a + 4*a^2;
776 : !
777 : ! and eventually as:
778 : !
779 : ! ( 4 + M * Skx^2 ) * a^2 - ( 4 + M * Skx^2 ) * a + 1 = 0.
780 : !
781 : ! Solving the quadratic equation for a:
782 : !
783 : ! a = (1/2) * ( 1 +- Skx * sqrt( 1 / ( 4/M + Skx^2 ) ) ).
784 : !
785 : ! Since by definition, mu_w_1 >= mu_w_2, a < 0.5 when Skw > 0, the equation
786 : ! for mixture fraction is:
787 : !
788 : ! a = (1/2) * ( 1 - Skx * sqrt( 1 / ( 4/M + Skx^2 ) ) ).
789 : !
790 : ! For 3-D Luhar, the variable (w, rt, or theta-l) with the greatest
791 : ! magnitude of skewness is used to calculate mixture fraction. Since it is
792 : ! desirable to still have a < 0.5 when Skw > 0 and a > 0.5 when Skw < 0, the
793 : ! sign function is used. The value of Skx is replaced by:
794 : !
795 : ! Skx|_adj = sgn( <w'x'> ) * Skx;
796 : !
797 : ! where
798 : !
799 : ! sgn( <w'x'> ) = | 1 when <w'x'> >= 0
800 : ! | -1 when <w'x'> < 0.
801 : !
802 : ! Since Skx|_adj^2 = ( sgn( <w'x'> ) * Skx )^2 = ( sgn( <w'x'> ) )^2 * Skx^2
803 : ! = Skx^2, the equation for mixture fraction is:
804 : !
805 : ! a = (1/2) * ( 1 - sgn( <w'x'> ) * Skx * sqrt( 1 / ( 4/M + Skx^2 ) ) ).
806 : !
807 : ! When using the ADG2 closure or when using the 3-D Luhar closure when the
808 : ! variable with the greatest magnitude of skewness is w, Skw = Skx and
809 : ! sgn( <w'^2> ) is always equal to 1, reducing the equation to its previous
810 : ! form.
811 :
812 : ! References:
813 : ! Vincent E. Larson, Jean-Christophe Golaz, and William R. Cotton, 2002:
814 : ! Small-Scale and Mesoscale Variability in Cloudy Boundary Layers: Joint
815 : ! Probability Density Functions. J. Atmos. Sci., 59, 3519–3539.
816 : !-----------------------------------------------------------------------
817 :
818 : use constants_clubb, only: &
819 : four, & ! Constant(s)
820 : three, &
821 : one, &
822 : two_thirds, &
823 : one_half, &
824 : one_third, &
825 : zero
826 :
827 : use clubb_precision, only: &
828 : core_rknd ! Precision
829 :
830 : implicit none
831 :
832 : ! Input Variables
833 : integer, intent(in) :: &
834 : nz ! Number of vertical levels
835 :
836 : real( kind = core_rknd ), dimension(nz), intent(in) :: &
837 : Skx, & ! Skewness of x ( <x'^3> / <x'2>^(3/2) ) [-]
838 : wpxp, & ! Covariance of w and x [m/s (x units)]
839 : xp2 ! Variance of x [(x units)^2]
840 :
841 : real( kind = core_rknd ), intent(in) :: &
842 : x_tol_sqd ! Tolerance value of x squared [(x units)^2]
843 :
844 : ! Output Variables
845 : real( kind = core_rknd ), dimension(nz), intent(out) :: &
846 : mixt_frac, & ! Mixture fraction [-]
847 : big_m, & ! Luhar's M [-]
848 : small_m ! Luhar's m [-]
849 :
850 : ! Local Variables
851 : real( kind = core_rknd ), dimension(nz) :: &
852 0 : small_m_sqd, & ! Luhar's m^2 [-]
853 0 : sgn_wpxp ! Sgn(<w'x'>); 1 when <w'x'> >= 0 or -1 when <w'x'> < 0 [-]
854 :
855 :
856 0 : where ( xp2 > x_tol_sqd )
857 :
858 : ! Width (standard deviation) parameters are non-zero
859 :
860 : ! Calculate Luhar's m (small m).
861 : ! If Skx is very small, then small_m will tend to zero which risks
862 : ! divide-by-zero. To ameliorate this problem, we enforce abs( x_1_n )
863 : ! and abs( x_2_n ) > 0.05.
864 : ! Note: Luhar's small_m (m) is the only tunable parameter in the Luhar
865 : ! closure, so this equation can be changed. However, the value
866 : ! of m should go toward 0 as Skx goes toward 0 so that the double
867 : ! Gaussian reduces to a single Gaussian when the distribution is
868 : ! unskewed.
869 : small_m = max( two_thirds * abs( Skx )**one_third, 0.05_core_rknd )
870 :
871 : ! Calculate m^2.
872 : small_m_sqd = small_m**2
873 :
874 : ! Calculate Luhar's M (big M).
875 : big_m = ( one + small_m_sqd )**3 &
876 : / ( ( three + small_m_sqd )**2 * small_m_sqd )
877 :
878 : ! Calculate sgn( <w'x'> ).
879 : where ( wpxp >= zero )
880 : sgn_wpxp = one
881 : elsewhere ! <w'x'> < 0
882 : sgn_wpxp = -one
883 : endwhere ! <w'x'> >= 0
884 :
885 : ! Calculate mixture fraction.
886 : mixt_frac = one_half &
887 : * ( one - sgn_wpxp * Skx &
888 : * sqrt( one / ( ( four / big_m ) + Skx**2 ) ) )
889 :
890 : elsewhere
891 :
892 : ! Variable x doesn't vary.
893 : mixt_frac = one_half
894 :
895 : ! For output purposes, set small_m and big_m to 0.
896 : small_m = zero
897 : big_m = zero
898 :
899 : endwhere ! Widths non-zero
900 :
901 :
902 0 : return
903 :
904 : end subroutine calc_Luhar_params
905 :
906 : !=============================================================================
907 0 : subroutine close_Luhar_pdf( nz, xm, xp2, mixt_frac, & ! In
908 0 : small_m, wpxp, x_tol_sqd, & ! In
909 0 : sigma_sqd_x_1, sigma_sqd_x_2, & ! Out
910 0 : varnce_x_1, varnce_x_2, & ! Out
911 0 : x_1_n, x_2_n, x_1, x_2 ) ! Out
912 :
913 : ! Description:
914 : ! For the Luhar closure, this subroutine takes Skx, xm, xp2, and mixt_frac,
915 : ! big_m, and small_m (calculated in calc_Luhar_params) as input and outputs
916 : ! the PDF component means and variances of a variable x in the joint-PDF
917 : ! according to Luhar et al. (1996). This code was written using the
918 : ! equations and nomenclature of Larson et al. (2002) Appendix section e.
919 :
920 : ! References:
921 : ! Vincent E. Larson, Jean-Christophe Golaz, and William R. Cotton, 2002:
922 : ! Small-Scale and Mesoscale Variability in Cloudy Boundary Layers: Joint
923 : ! Probability Density Functions. J. Atmos. Sci., 59, 3519–3539.
924 : !-----------------------------------------------------------------------
925 :
926 : use constants_clubb, only: &
927 : one, & ! Constant(s)
928 : zero
929 :
930 : use clubb_precision, only: &
931 : core_rknd ! Precision
932 :
933 : implicit none
934 :
935 : ! Input Variables
936 : integer, intent(in) :: &
937 : nz ! Number of vertical levels
938 :
939 : real( kind = core_rknd ), dimension(nz), intent(in) :: &
940 : xm, & ! Mean (overall) of x, <x> [(x units)]
941 : xp2, & ! Variance (overall) of x, <x'^2> [(x units)^2]
942 : mixt_frac, & ! Mixture fraction [-]
943 : small_m, & ! Luhar's small m [-]
944 : wpxp ! Covariance of w and x, <w'x'> [m/s (x units)]
945 :
946 : real( kind = core_rknd ), intent(in) :: &
947 : x_tol_sqd ! Tolerance value of x squared [(x units)^2]
948 :
949 : ! Output Variables
950 : real( kind = core_rknd ), dimension(nz), intent(out) :: &
951 : sigma_sqd_x_1, & ! Normalized width parameter of x (1st PDF component) [-]
952 : sigma_sqd_x_2, & ! Normalized width parameter of x (1st PDF component) [-]
953 : varnce_x_1, & ! Variance of x (1st PDF component) [(x units)^2]
954 : varnce_x_2, & ! Variance of x (2nd PDF component) [(x units)^2]
955 : x_1_n, & ! Normalized mean of x (1st PDF component) [-]
956 : x_2_n, & ! Normalized mean of x (2nd PDF component) [-]
957 : x_1, & ! Mean of x (1st PDF component) [(x units)]
958 : x_2 ! Mean of x (2nd PDF component) [(x units)]
959 :
960 : ! Local Variables
961 : real( kind = core_rknd), dimension(nz) :: &
962 0 : sqrt_xp2, & ! Square root of the variance of x [(x units)]
963 0 : sgn_wpxp ! Sgn( <w'x'> ); 1 when <w'x'> >= 0 or -1 when <w'x'> < 0 [-]
964 :
965 :
966 0 : where ( xp2 > x_tol_sqd )
967 :
968 : ! Width (standard deviation) parameters are non-zero
969 :
970 : ! Calculate sgn( <w'x'> ).
971 : where ( wpxp >= zero )
972 : sgn_wpxp = one
973 : elsewhere ! <w'x'> < 0
974 : sgn_wpxp = -one
975 : endwhere ! <w'x'> >= 0
976 :
977 : ! Calculate the square root of the overall variance of x.
978 : sqrt_xp2 = sqrt( xp2 )
979 :
980 : ! Normalized width parameter of x in the 1st PDF component.
981 : sigma_sqd_x_1 &
982 : = ( one - mixt_frac ) / ( mixt_frac * ( one + small_m**2 ) )
983 :
984 : ! The variance of x in the 1st PDF component.
985 : varnce_x_1 = sigma_sqd_x_1 * xp2
986 :
987 : ! Normalized width parameter of x in the 2nd PDF component.
988 : sigma_sqd_x_2 &
989 : = mixt_frac / ( ( one - mixt_frac ) * ( one + small_m**2 ) )
990 :
991 : ! The variance of x in the 2nd PDF component.
992 : varnce_x_2 = sigma_sqd_x_2 * xp2
993 :
994 : ! Normalized mean of x in the 1st PDF component.
995 : x_1_n = sgn_wpxp * small_m * sqrt( sigma_sqd_x_1 )
996 :
997 : ! Normalized mean of x in the 2nd PDF component.
998 : x_2_n = -sgn_wpxp * small_m * sqrt( sigma_sqd_x_2 )
999 :
1000 : ! The mean of x in the 1st PDF component.
1001 : x_1 = xm + sqrt_xp2 * x_1_n
1002 :
1003 : ! The mean of x in the 2nd PDF component.
1004 : x_2 = xm + sqrt_xp2 * x_2_n
1005 :
1006 : elsewhere
1007 :
1008 : ! Variable x doesn't vary.
1009 : sigma_sqd_x_1 = ( one / ( one + small_m**2 ) )
1010 : sigma_sqd_x_2 = ( one / ( one + small_m**2 ) )
1011 : varnce_x_1 = zero
1012 : varnce_x_2 = zero
1013 : x_1_n = sgn_wpxp * small_m * sqrt( sigma_sqd_x_1 )
1014 : x_2_n = -sgn_wpxp * small_m * sqrt( sigma_sqd_x_2 )
1015 : x_1 = xm
1016 : x_2 = xm
1017 :
1018 : endwhere ! Widths non-zero
1019 :
1020 :
1021 0 : return
1022 :
1023 : end subroutine close_Luhar_pdf
1024 :
1025 : !=============================================================================
1026 2823552 : subroutine ADG1_ADG2_responder_params( nz, ngrdcol, & ! In
1027 2823552 : xm, xp2, wp2, sqrt_wp2, & ! In
1028 2823552 : wpxp, w_1_n, w_2_n, mixt_frac, & ! In
1029 2823552 : sigma_sqd_w, beta, x_tol, & ! In
1030 2823552 : x_1, x_2, varnce_x_1, & ! Out
1031 2823552 : varnce_x_2, alpha_x ) ! Out
1032 :
1033 : ! Description:
1034 : ! Calculates the PDF component means and PDF component variances for thl,
1035 : ! rt, or sclr when either the ADG1 PDF or ADG2 PDF are used.
1036 : !
1037 : ! The normalized variance for thl, rt, and sclr for "plume" 1 is:
1038 : !
1039 : ! { 1 - [1/(1-sigma_sqd_w)]*[ <w'x'>^2 / (<w'^2> * <x'^2>) ] / mixt_frac }
1040 : ! * { (1/3)*beta + mixt_frac*( 1 - (2/3)*beta ) };
1041 : !
1042 : ! where "x" stands for thl, rt, or sclr; "mixt_frac" is the weight of
1043 : ! Gaussian "plume" 1, and 0 <= beta <= 3.
1044 : !
1045 : ! The factor { (1/3)*beta + mixt_frac*( 1 - (2/3)*beta ) } does not depend
1046 : ! on which varable "x" stands for. The factor is multiplied by 2 and
1047 : ! defined as width_factor_1.
1048 : !
1049 : ! The factor
1050 : ! { 1 - [1/(1-sigma_sqd_w)]*[ (w'x')^2 / (w'^2 * x'^2) ] / mixt_frac }
1051 : ! depends on which variable "x" stands for. It is multiplied by 1/2 and
1052 : ! defined as alpha_x, where "x" stands for thl, rt, or sclr.
1053 :
1054 : ! Vince Larson added a dimensionless factor so that the
1055 : ! width of plumes in theta_l, rt can vary.
1056 : ! beta is a constant defined in CLUBB's tunable parameters.
1057 : ! Set 0<beta<3.
1058 : ! beta=1.5_core_rknd recovers Chris Golaz' simplified formula.
1059 : ! 3 Nov 2003
1060 :
1061 : ! References:
1062 :
1063 : use constants_clubb, only: &
1064 : two, & ! Constant(s)
1065 : one, &
1066 : two_thirds, &
1067 : one_half, &
1068 : zero_threshold, &
1069 : w_tol_sqd
1070 :
1071 : use clubb_precision, only: &
1072 : core_rknd ! Variable(s)
1073 :
1074 : implicit none
1075 :
1076 : integer, intent(in) :: &
1077 : ngrdcol, & ! Number of grid columns
1078 : nz ! Number of vertical level
1079 :
1080 : ! Input Variables
1081 : real ( kind = core_rknd ), dimension(ngrdcol,nz), intent(in) :: &
1082 : xm, & ! Mean of x (overall) [units vary]
1083 : xp2, & ! Variance of x (overall) [(units vary)^2]
1084 : wp2, & ! Variance of w (overall) [m^2/s^2]
1085 : sqrt_wp2, & ! Square root of the variance of w [m/s]
1086 : wpxp, & ! Covariance of w and x [m/s (units vary)]
1087 : w_1_n, & ! Normalized mean of w (1st PDF comp.) [-]
1088 : w_2_n, & ! Normalized mean of w (2nd PDF comp.) [-]
1089 : mixt_frac, & ! Mixture fraction [-]
1090 : sigma_sqd_w ! Width of individual w plumes [-]
1091 :
1092 : real ( kind = core_rknd ), intent(in) :: &
1093 : beta, & ! CLUBB tunable parameter beta [-]
1094 : x_tol ! Tolerance value for x [units vary]
1095 :
1096 : ! Output Variables
1097 : real ( kind = core_rknd ), dimension(ngrdcol,nz), intent(out) :: &
1098 : x_1, & ! Mean of x (1st PDF component) [units vary]
1099 : x_2, & ! Mean of x (2nd PDF component) [units vary]
1100 : varnce_x_1, & ! Variance of x (1st PDF component) [(units vary)^2]
1101 : varnce_x_2, & ! Variance of x (2nd PDF component) [(units vary)^2]
1102 : alpha_x ! Factor relating to normalized variance for x [-]
1103 :
1104 : ! Local Variables
1105 : ! variables for a generalization of Chris Golaz' closure
1106 : ! varies width of plumes in theta_l, rt
1107 : real ( kind = core_rknd ) :: &
1108 : width_factor_1 ! Width factor relating to PDF component 1 [-]
1109 :
1110 : integer :: &
1111 : k, i ! Vertical loop index
1112 :
1113 : !----- Begin Code -----
1114 : !$acc parallel loop gang vector collapse(2) default(present)
1115 242825472 : do k = 1, nz
1116 4010297472 : do i = 1, ngrdcol
1117 :
1118 4007473920 : if ( wp2(i,k) <= w_tol_sqd .or. xp2(i,k) <= x_tol**2 ) then
1119 :
1120 : ! Vertical velocity doesn't vary. Variable x is treated as a single
1121 : ! Gaussian in this situation.
1122 2865231606 : x_1(i,k) = xm(i,k)
1123 2865231606 : x_2(i,k) = xm(i,k)
1124 2865231606 : varnce_x_1(i,k) = xp2(i,k)
1125 2865231606 : varnce_x_2(i,k) = xp2(i,k)
1126 2865231606 : alpha_x(i,k) = one_half
1127 :
1128 : else
1129 :
1130 902240394 : x_1(i,k) = xm(i,k) - wpxp(i,k) / ( sqrt_wp2(i,k) * w_2_n(i,k) )
1131 902240394 : x_2(i,k) = xm(i,k) - wpxp(i,k) / ( sqrt_wp2(i,k) * w_1_n(i,k) )
1132 :
1133 : alpha_x(i,k) = one_half &
1134 : * ( one - wpxp(i,k) * wpxp(i,k) &
1135 902240394 : / ( ( one - sigma_sqd_w(i,k) ) * wp2(i,k) * xp2(i,k) ) )
1136 :
1137 902240394 : alpha_x(i,k) = max( min( alpha_x(i,k), one ), zero_threshold )
1138 :
1139 : width_factor_1 = two_thirds * beta &
1140 902240394 : + two * mixt_frac(i,k) * ( one - two_thirds * beta )
1141 :
1142 : ! Vince Larson multiplied original expressions by width_factor_1,2
1143 : ! to generalize scalar skewnesses. 05 Nov 03
1144 902240394 : varnce_x_1(i,k) = width_factor_1 * xp2(i,k) * alpha_x(i,k) / mixt_frac(i,k)
1145 : varnce_x_2(i,k) = ( two - width_factor_1 ) * xp2(i,k) &
1146 902240394 : * alpha_x(i,k) / ( one - mixt_frac(i,k) )
1147 :
1148 : end if
1149 :
1150 : end do
1151 : end do
1152 : !$acc end parallel loop
1153 :
1154 2823552 : return
1155 :
1156 : end subroutine ADG1_ADG2_responder_params
1157 :
1158 : !=============================================================================
1159 0 : subroutine backsolve_Luhar_params( Sk_max, Skx, & ! In
1160 : big_m_max, mixt_frac, & ! In
1161 : big_m_x, small_m_x ) ! Out
1162 :
1163 : ! Description:
1164 : ! This subroutine calculates Luhar's big_m and small_m for the variate 'x'
1165 : ! consistent with the mixture fraction of the variate with the largest
1166 : ! skewness.
1167 : !
1168 : ! The relationship between skewness of x (Skx), mixture fraction (a), and
1169 : ! Luhar's small m (m) is given by:
1170 : !
1171 : ! Skx^2 = ( m^2 * ( m^2 + 3 )^2 / ( m^2 + 1 )^3 )
1172 : ! * ( 1 - 2*a )^2 / ( a * ( 1 - a ) ).
1173 : !
1174 : ! Moving the factor involving mixture fraction to the right-hand side:
1175 : !
1176 : ! ( ( a * ( 1 - a ) ) / ( 1 - 2*a )^2 ) * Skx^2
1177 : ! = m^2 * ( m^2 + 3 )^2 / ( m^2 + 1 )^3.
1178 : !
1179 : ! This can be rewritten as:
1180 : !
1181 : ! ( ( a * ( 1 - a ) ) / ( 1 - 2*a )^2 ) * Skx^2
1182 : ! = ( m^6 + 6*m^4 + 9*m^2 ) / ( m^6 + 3*m^4 + 3*m^2 + 1 ).
1183 : !
1184 : ! Setting alpha = ( ( a * ( 1 - a ) ) / ( 1 - 2*a )^2 ) * Skx^2, the
1185 : ! equation can be rewritten as:
1186 : !
1187 : ! ( m^6 + 3*m^4 + 3*m^2 + 1 ) * alpha = m^6 + 6*m^4 + 9*m^2.
1188 : !
1189 : ! This can be rearranged and rewritten as:
1190 : !
1191 : ! ( alpha - 1 ) * m^6 + ( 3 * alpha - 6 ) * m^4
1192 : ! + ( 3 * alpha - 9 ) * m^2 + alpha = 0.
1193 : !
1194 : ! This can be rewritten again as:
1195 : !
1196 : ! ( alpha - 1 ) * (m^2)^3 + ( 3 * alpha - 6 ) * (m^2)^2
1197 : ! + ( 3 * alpha - 9 ) * (m^2) + alpha = 0.
1198 : !
1199 : ! The goal is to solve for m^2, and then take the square root of m^2 to
1200 : ! solve for m. This can be accomplished by using the cubic formula (with
1201 : ! the l_use_cubic_backsolve option), or else by a quadratic approximation.
1202 :
1203 : ! References:
1204 : !-----------------------------------------------------------------------
1205 :
1206 : use constants_clubb, only: &
1207 : nine, & ! Constant(s)
1208 : six, &
1209 : five, &
1210 : four, &
1211 : three, &
1212 : two, &
1213 : one, &
1214 : one_half, &
1215 : zero, &
1216 : eps
1217 :
1218 : use clubb_precision, only: &
1219 : core_rknd ! Precision
1220 :
1221 : implicit none
1222 :
1223 : ! Input Variables
1224 : real( kind = core_rknd ), intent(in) :: &
1225 : Sk_max, & ! Maximum skewness [-]
1226 : Skx, & ! Skewness of the variate solving small_m and big_m for [-]
1227 : big_m_max, & ! Luhar's big_m of the variate with maximum skewness [-]
1228 : mixt_frac ! Mixture fraction [-]
1229 :
1230 : ! Output Variables
1231 : real( kind = core_rknd ), intent(out) :: &
1232 : big_m_x, & ! Luhar's big_m for the variate being solved for [-]
1233 : small_m_x ! Luhar's small_m for the variate being solved for [-]
1234 :
1235 : ! Local Variables
1236 : real( kind = core_rknd ) :: &
1237 : alpha, & ! 1 / big_m_x
1238 : a, & ! For readability, quadratic equation
1239 : b, &
1240 : c, &
1241 : alpha_upr, &
1242 : alpha_low, &
1243 : discrim
1244 :
1245 : real( kind = core_rknd ), dimension(1) :: &
1246 : a_coef_in, & ! Coefficient a (of x^3) in a*x^3 + b*x^2 + c^x + d = 0 [-]
1247 : b_coef_in, & ! Coefficient b (of x^2) in a*x^3 + b*x^2 + c^x + d = 0 [-]
1248 : c_coef_in, & ! Coefficient c (of x) in a*x^3 + b*x^2 + c^x + d = 0 [-]
1249 : d_coef_in ! Coefficient d in a*x^3 + b*x^2 + c^x + d = 0 [-]
1250 :
1251 : ! Flag to backsolve for m^2 using cubic formula
1252 : logical, parameter :: &
1253 : l_use_cubic_backsolve = .true.
1254 :
1255 :
1256 : if ( l_use_cubic_backsolve ) then
1257 :
1258 0 : if ( abs( mixt_frac - one_half ) < 0.001_core_rknd ) then
1259 :
1260 : ! When mixture fraction = 0.5 (based on the variable with the largest
1261 : ! magnitude of skewness), all variables must have a skewness of 0.
1262 : ! Set m to the minimum threshold of 0.05.
1263 0 : small_m_x = 0.05_core_rknd
1264 :
1265 : ! Calculate the corresponding value of big_m_x.
1266 : big_m_x = ( one + small_m_x**2 )**3 &
1267 0 : / ( ( three + small_m_x**2 )**2 * small_m_x**2 )
1268 :
1269 0 : elseif ( abs(Skx) < eps ) then
1270 :
1271 : ! Mixture fraction /= 0.5 because the variable with the largest
1272 : ! magnitude of skewness has a skewness /= 0. However, variable x has
1273 : ! a skewness of 0. In order to reproduce the correct skewness for
1274 : ! variable x, set m to 0 (regardless of minimum thresholds used in
1275 : ! other parts of the code).
1276 0 : small_m_x = zero
1277 :
1278 : ! The value of big_m_x should be inf. Set it to huge. This is not
1279 : ! used in any calculation, anyway.
1280 0 : big_m_x = huge( big_m_x )
1281 :
1282 : else ! mixt_frac /= 0.5 and Skx /= 0
1283 :
1284 : ! Backsolve for m, given mixt_frac and Skx.
1285 :
1286 : ! alpha = 1/M is given by:
1287 : ! [ mixt_frac * ( 1 - mixt_frac ) / ( 1 - 2 * mixt_frac )^2 ] * Skx^2.
1288 : alpha = ( mixt_frac * ( one - mixt_frac ) &
1289 0 : / ( one - two * mixt_frac )**2 ) * Skx**2
1290 :
1291 : ! Calculate big_m_x.
1292 0 : big_m_x = one / alpha
1293 :
1294 : ! Solve the cubic equation for m^2:
1295 : ! ( alpha - 1 ) * (m^2)^3 + ( 3 * alpha - 6 ) * (m^2)^2
1296 : ! + ( 3 * alpha - 9 ) * (m^2) + alpha = 0.
1297 : ! The largest root is preferred.
1298 0 : a_coef_in(1) = alpha - one
1299 0 : b_coef_in(1) = three * alpha - six
1300 0 : c_coef_in(1) = three * alpha - nine
1301 0 : d_coef_in(1) = alpha
1302 : small_m_x &
1303 : = sqrt( max( max_cubic_root( a_coef_in, b_coef_in, &
1304 : c_coef_in, d_coef_in ), &
1305 0 : 0.05_core_rknd**2 ) )
1306 :
1307 : endif
1308 :
1309 : else ! original formualation
1310 :
1311 : alpha = ( Skx**2 / (max(Sk_max**2, eps) * big_m_max) ) ! 1 / big_m_x
1312 :
1313 : ! This limit keeps the discriminant >= 0
1314 : alpha_upr = two * sqrt( 13.0_core_rknd ) - five
1315 :
1316 : alpha_low = eps
1317 :
1318 : ! For this approximation, alpha must be less than 2*sqrt(13) - 5 to get a
1319 : ! real ans.
1320 : alpha = min( alpha, alpha_upr )
1321 :
1322 : ! For testing, eliminate possibility of divide by zero
1323 : alpha = max( alpha, alpha_low )
1324 :
1325 : ! Use a piece-wise approximation
1326 : if ( alpha < one ) then
1327 :
1328 : a = max( three * alpha - six, eps) ! Prevent divide by zero
1329 : b = three * alpha - nine
1330 : c = alpha
1331 :
1332 : discrim = b**2 - four * a * c
1333 : small_m_x = sqrt( ( - b - sqrt( discrim ) ) / ( two * a ) )
1334 :
1335 : else ! alpha >= 1
1336 :
1337 : ! For this approximation, alpha must be less than 2*sqrt(13) - 5
1338 : ! to get a real ans.
1339 : alpha = min( alpha, two )
1340 :
1341 : a = max( six * alpha - nine, eps) ! Prevent divide by zero
1342 : b = -six
1343 : c = two * alpha - one
1344 :
1345 : discrim = b**2 - four * a * c
1346 : small_m_x = sqrt( ( - b - sqrt( discrim ) ) / ( two * a ) )
1347 :
1348 : endif ! alpha < 1
1349 :
1350 : ! Clip consistently with subroutine calc_Luhar_params
1351 : small_m_x = max( 5e-2_core_rknd, small_m_x)
1352 :
1353 : big_m_x = one / alpha
1354 :
1355 : endif ! l_use_cubic_backsolve
1356 :
1357 :
1358 0 : return
1359 :
1360 : end subroutine backsolve_Luhar_params
1361 :
1362 : !=============================================================================
1363 0 : function max_cubic_root( a_coef, b_coef, c_coef, d_coef ) &
1364 : result( max_root )
1365 :
1366 : ! Description:
1367 : ! Calculates the largest root that results from solving a cubic equation of
1368 : ! the form a*x^3 + b*x^2 + c*x + d = 0.
1369 : !
1370 : ! This is done to backsolve for m^2 for the 3-D Luhar closure, given the
1371 : ! values of mixt_frac and Skx.
1372 :
1373 : ! References:
1374 : !-----------------------------------------------------------------------
1375 :
1376 : use constants_clubb, only: &
1377 : eps ! Constant(s)
1378 :
1379 : use calc_roots, only: &
1380 : cubic_solve, & ! Procedure(s)
1381 : quadratic_solve
1382 :
1383 : use clubb_precision, only: &
1384 : core_rknd ! Variable(s)
1385 :
1386 : implicit none
1387 :
1388 : ! Input Variables
1389 : real( kind = core_rknd ), dimension(1), intent(in) :: &
1390 : a_coef, & ! Coefficient a (of x^3) in a*x^3 + b*x^2 + c^x + d = 0 [-]
1391 : b_coef, & ! Coefficient b (of x^2) in a*x^3 + b*x^2 + c^x + d = 0 [-]
1392 : c_coef, & ! Coefficient c (of x) in a*x^3 + b*x^2 + c^x + d = 0 [-]
1393 : d_coef ! Coefficient d in a*x^3 + b*x^2 + c^x + d = 0 [-]
1394 :
1395 : ! Return Variable
1396 : real( kind = core_rknd ) :: &
1397 : max_root ! Maximum root that solves the cubic equation [-]
1398 :
1399 : ! Local Variables
1400 : complex( kind = core_rknd ), dimension(1,3) :: &
1401 : cubic_roots ! Roots of x that satisfy a*x^3 + b*x^2 + c*x + d = 0 [-]
1402 :
1403 : complex( kind = core_rknd ), dimension(1,2) :: &
1404 : quadratic_roots ! Roots of x that satisfy b*x^2 + c*x + d = 0 [-]
1405 :
1406 : real( kind = core_rknd ) :: &
1407 : a_coef_thresh, & ! Minimum threshold of |a| to use cubic solver [-]
1408 : b_coef_thresh ! Minimum threshold of |b| to use quadratic solver [-]
1409 :
1410 :
1411 : ! Calculate a minimum threshold for |a| to call this a cubic equation.
1412 : a_coef_thresh = 0.001_core_rknd &
1413 0 : * max( abs(b_coef(1)), abs(c_coef(1)), abs(d_coef(1)) )
1414 :
1415 : ! Calculate a minimum threshold for |b| to call this a quadratic equation.
1416 : ! This only matters when |a| <= a_coef_thresh.
1417 0 : b_coef_thresh = 0.001_core_rknd * max( abs(c_coef(1)), abs(d_coef(1)) )
1418 :
1419 0 : if ( abs( a_coef(1) ) > a_coef_thresh ) then
1420 :
1421 : ! The equation is a cubic equation.
1422 0 : cubic_roots = cubic_solve( 1, a_coef, b_coef, c_coef, d_coef )
1423 :
1424 0 : if ( abs(aimag( cubic_roots(1,2) )) < eps .and. &
1425 : abs(aimag( cubic_roots(1,3) )) < eps ) then
1426 :
1427 : ! Find the maximum root of the three roots.
1428 : max_root = max( real( cubic_roots(1,1), kind = core_rknd ), &
1429 : real( cubic_roots(1,2), kind = core_rknd ), &
1430 0 : real( cubic_roots(1,3), kind = core_rknd ) )
1431 :
1432 : else ! cubic_roots(2) and cubic_roots(3) are complex.
1433 :
1434 0 : max_root = real( cubic_roots(1,1), kind = core_rknd )
1435 :
1436 : endif
1437 :
1438 0 : elseif ( abs( b_coef(1) ) > b_coef_thresh ) then
1439 :
1440 : ! The equation is a quadratic equation, since a = 0, but b /= 0.
1441 : ! This should very rarely occur for 3-D Luhar. When it does, the result
1442 : ! will always be two real-valued roots.
1443 0 : quadratic_roots = quadratic_solve( 1, b_coef, c_coef, d_coef )
1444 :
1445 : ! Find the maximum root of the two roots.
1446 : max_root = max( real( quadratic_roots(1,1), kind = core_rknd ), &
1447 0 : real( quadratic_roots(1,2), kind = core_rknd ) )
1448 :
1449 : else ! |a| = 0 and |b| = 0
1450 :
1451 : ! The equation is a linear equation.
1452 : ! This won't happen for 3-D Luhar.
1453 0 : max_root = - d_coef(1) / c_coef(1)
1454 :
1455 : endif ! |a| > 0
1456 :
1457 :
1458 : return
1459 :
1460 : end function max_cubic_root
1461 :
1462 : !=============================================================================
1463 :
1464 : end module adg1_adg2_3d_luhar_pdf
|