-
Notifications
You must be signed in to change notification settings - Fork 12
/
pdf_parameter_module.F90
544 lines (466 loc) · 25.4 KB
/
pdf_parameter_module.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
!-----------------------------------------------------------------------
! $Id$
!===============================================================================
module pdf_parameter_module
! Description:
! This module defines the derived type pdf_parameter.
! References:
! None
!-----------------------------------------------------------------------
use clubb_precision, only: &
core_rknd
implicit none
private ! Default scope
public :: pdf_parameter, & ! Variable Type(s)
implicit_coefs_terms, &
init_pdf_params, & ! Procedure(s)
init_pdf_implicit_coefs_terms
! CLUBB's PDF parameters.
type pdf_parameter
real( kind = core_rknd ), dimension(:), allocatable :: &
w_1, & ! Mean of w (1st PDF component) [m/s]
w_2, & ! Mean of w (2nd PDF component) [m/s]
varnce_w_1, & ! Variance of w (1st PDF component) [m^2/s^2]
varnce_w_2, & ! Variance of w (2nd PDF component) [m^2/s^2]
rt_1, & ! Mean of r_t (1st PDF component) [kg/kg]
rt_2, & ! Mean of r_t (2nd PDF component) [kg/kg]
varnce_rt_1, & ! Variance of r_t (1st PDF component) [kg^2/kg^2]
varnce_rt_2, & ! Variance of r_t (2nd PDF component) [kg^2/kg^2]
thl_1, & ! Mean of th_l (1st PDF component) [K]
thl_2, & ! Mean of th_l (2nd PDF component) [K]
varnce_thl_1, & ! Variance of th_l (1st PDF component) [K^2]
varnce_thl_2, & ! Variance of th_l (2nd PDF component) [K^2]
corr_w_rt_1, & ! Correlation of w and r_t (1st PDF component) [-]
corr_w_rt_2, & ! Correlation of w and r_t (2nd PDF component) [-]
corr_w_thl_1, & ! Correlation of w and th_l (1st PDF component) [-]
corr_w_thl_2, & ! Correlation of w and th_l (2nd PDF component) [-]
corr_rt_thl_1, & ! Correlation of r_t and th_l (1st PDF component) [-]
corr_rt_thl_2, & ! Correlation of r_t and th_l (2nd PDF component) [-]
alpha_thl, & ! Factor relating to normalized variance for th_l [-]
alpha_rt, & ! Factor relating to normalized variance for r_t [-]
crt_1, & ! r_t coef. in chi/eta eqns. (1st PDF comp.) [-]
crt_2, & ! r_t coef. in chi/eta eqns. (2nd PDF comp.) [-]
cthl_1, & ! th_l coef.: chi/eta eqns. (1st PDF comp.) [(kg/kg)/K]
cthl_2, & ! th_l coef.: chi/eta eqns. (2nd PDF comp.) [(kg/kg)/K]
chi_1, & ! Mean of chi (old s) (1st PDF component) [kg/kg]
chi_2, & ! Mean of chi (old s) (2nd PDF component) [kg/kg]
stdev_chi_1, & ! Standard deviation of chi (1st PDF component) [kg/kg]
stdev_chi_2, & ! Standard deviation of chi (2nd PDF component) [kg/kg]
stdev_eta_1, & ! Standard dev. of eta (old t) (1st PDF comp.) [kg/kg]
stdev_eta_2, & ! Standard dev. of eta (old t) (2nd PDF comp.) [kg/kg]
covar_chi_eta_1, & ! Covariance of chi and eta (1st PDF comp.) [kg^2/kg^2]
covar_chi_eta_2, & ! Covariance of chi and eta (2nd PDF comp.) [kg^2/kg^2]
corr_w_chi_1, & ! Correlation of w and chi (1st PDF component) [-]
corr_w_chi_2, & ! Correlation of w and chi (2nd PDF component) [-]
corr_w_eta_1, & ! Correlation of w and eta (1st PDF component) [-]
corr_w_eta_2, & ! Correlation of w and eta (2nd PDF component) [-]
corr_chi_eta_1, & ! Correlation of chi and eta (1st PDF component) [-]
corr_chi_eta_2, & ! Correlation of chi and eta (2nd PDF component) [-]
rsatl_1, & ! Saturation mixing ratio r_sat(mu_Tl_1,p) [kg/kg]
rsatl_2, & ! Saturation mixing ratio r_sat(mu_Tl_2,p) [kg/kg]
rc_1, & ! Mean of r_c (1st PDF component) [kg/kg]
rc_2, & ! Mean of r_c (2nd PDF component) [kg/kg]
cloud_frac_1, & ! Cloud fraction (1st PDF component) [-]
cloud_frac_2, & ! Cloud fraction (2nd PDF component) [-]
mixt_frac, & ! Weight of 1st PDF component (Sk_w dependent) [-]
ice_supersat_frac_1, & ! Ice supersaturation fraction (1st PDF comp.) [-]
ice_supersat_frac_2 ! Ice supersaturation fraction (2nd PDF comp.) [-]
end type pdf_parameter
! The implicit coefficients, semi-implicit coefficients and terms, and
! explicit terms for turbulent advection of turbulent fields are calculated
! from the PDF and the resulting PDF parameters.
type implicit_coefs_terms
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wp4_implicit ! <w'^4> = coef_wp4_implicit * <w'^2>^2 [-]
! <w'^2 rt'> = coef_wp2rtp_implicit * <w'rt'> + term_wp2rtp_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wp2rtp_implicit, & ! Coefficient that is multiplied by <w'rt'> [m/s]
term_wp2rtp_explicit ! Term that is on the RHS [m^2/s^2 kg/kg]
! <w'^2 thl'> = coef_wp2thlp_implicit * <w'thl'> + term_wp2thlp_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wp2thlp_implicit, & ! Coef. that is multiplied by <w'thl'> [m/s]
term_wp2thlp_explicit ! Term that is on the RHS [m^2/s^2 K]
! <w'^2 u'> = coef_wp2up_implicit * <u'w'> + term_wp2up_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wp2up_implicit, & ! Coefficient that is multiplied by <u'w'> [m/s]
term_wp2up_explicit ! Term that is on the RHS [m^3/s^3]
! <w'^2 v'> = coef_wp2vp_implicit * <v'w'> + term_wp2vp_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wp2vp_implicit, & ! Coefficient that is multiplied by <v'w'> [m/s]
term_wp2vp_explicit ! Term that is on the RHS [m^3/s^3]
! <w'rt'^2> = coef_wprtp2_implicit * <rt'^2> + term_wprtp2_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wprtp2_implicit, & ! Coefficient that is multiplied by <rt'^2> [m/s]
term_wprtp2_explicit ! Term that is on the RHS [m/s kg^2/kg^2]
! <w'thl'^2> = coef_wpthlp2_implicit * <thl'^2> + term_wpthlp2_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wpthlp2_implicit, & ! Coef. that is multiplied by <thl'^2> [m/s]
term_wpthlp2_explicit ! Term that is on the RHS [m/s K^2]
! <w'rt'thl'> = coef_wprtpthlp_implicit*<rt'thl'> + term_wprtpthlp_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wprtpthlp_implicit, & ! Coef. that is multiplied by <rt'thl'> [m/s]
term_wprtpthlp_explicit ! Term that is on the RHS [m/s(kg/kg)K]
! <w'u'^2> = coef_wpup2_implicit * <u'^2> + term_wpup2_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wpup2_implicit, & ! Coefficient that is multiplied by <u'^2> [m/s]
term_wpup2_explicit ! Term that is on the RHS [m^3/s^3]
! <w'v'^2> = coef_wpvp2_implicit * <v'^2> + term_wpvp2_explicit
real( kind = core_rknd ), dimension(:), allocatable :: &
coef_wpvp2_implicit, & ! Coefficient that is multiplied by <v'^2> [m/s]
term_wpvp2_explicit ! Term that is on the RHS [m^3/s^3]
! <w'^2 sclr'> = coef_wp2sclrp_implicit * <w'sclr'> + term_wp2sclrp_explicit
real( kind = core_rknd ), dimension(:,:), allocatable :: &
coef_wp2sclrp_implicit, & ! Coef. that is multiplied by <w'sclr'> [m/s]
term_wp2sclrp_explicit ! Term that is on the RHS [m^2/s^2 (un. vary)]
! <w'sclr'^2> = coef_wpsclrp2_implicit * <sclr'^2> + term_wpsclrp2_explicit
real( kind = core_rknd ), dimension(:,:), allocatable :: &
coef_wpsclrp2_implicit, & ! Coef. that is multiplied by <sclr'^2> [m/s]
term_wpsclrp2_explicit ! Term that is on the RHS [m/s(units vary)^2]
! <w'rt'sclr'> = coef_wprtpsclrp_implicit * <sclr'rt'>
! + term_wprtpsclrp_explicit
real( kind = core_rknd ), dimension(:,:), allocatable :: &
coef_wprtpsclrp_implicit, & ! Coef. that is multiplied by <sclr'rt'> [m/s]
term_wprtpsclrp_explicit ! Term that is on the RHS [m/s(kg/kg)(un. v.)]
! <w'thl'sclr'> = coef_wpthlpsclrp_implicit * <sclr'thl'>
! + term_wpthlpsclrp_explicit
real( kind = core_rknd ), dimension(:,:), allocatable :: &
coef_wpthlpsclrp_implicit, & ! Coef. that is mult. by <sclr'thl'> [m/s]
term_wpthlpsclrp_explicit ! Term that is on the RHS [(m/s)K(un. vary)]
end type implicit_coefs_terms
! The CLUBB_CAM preprocessor directives are being commented out because this
! code is now also used for WRF-CLUBB.
!#ifdef CLUBB_CAM /* Code for storing pdf_parameter structs in pbuf as array */
public :: pack_pdf_params, unpack_pdf_params
integer, public, parameter :: num_pdf_params = 47
!#endif /* CLUBB_CAM */
contains
!=============================================================================
subroutine init_pdf_params( nz, pdf_params )
! Description:
! Initializes all PDF parameters in the variable type pdf_parameter.
! References:
!--------------------------------------------------------------------
use constants_clubb, only: &
zero ! Constant(s)
implicit none
! Input Variable(s)
integer, intent(in) :: &
nz ! Number of vertical grid levels [-]
! Output Variable(s)
type(pdf_parameter), intent(out) :: &
pdf_params ! PDF parameters [units vary]
allocate( pdf_params%w_1(nz), &
pdf_params%w_2(nz), &
pdf_params%varnce_w_1(nz), &
pdf_params%varnce_w_2(nz), &
pdf_params%rt_1(nz), &
pdf_params%rt_2(nz), &
pdf_params%varnce_rt_1(nz), &
pdf_params%varnce_rt_2(nz), &
pdf_params%thl_1(nz), &
pdf_params%thl_2(nz), &
pdf_params%varnce_thl_1(nz), &
pdf_params%varnce_thl_2(nz), &
pdf_params%corr_w_rt_1(nz), &
pdf_params%corr_w_rt_2(nz), &
pdf_params%corr_w_thl_1(nz), &
pdf_params%corr_w_thl_2(nz), &
pdf_params%corr_rt_thl_1(nz), &
pdf_params%corr_rt_thl_2(nz), &
pdf_params%alpha_thl(nz), &
pdf_params%alpha_rt(nz), &
pdf_params%crt_1(nz), &
pdf_params%crt_2(nz), &
pdf_params%cthl_1(nz), &
pdf_params%cthl_2(nz), &
pdf_params%chi_1(nz), &
pdf_params%chi_2(nz), &
pdf_params%stdev_chi_1(nz), &
pdf_params%stdev_chi_2(nz), &
pdf_params%stdev_eta_1(nz), &
pdf_params%stdev_eta_2(nz), &
pdf_params%covar_chi_eta_1(nz), &
pdf_params%covar_chi_eta_2(nz), &
pdf_params%corr_w_chi_1(nz), &
pdf_params%corr_w_chi_2(nz), &
pdf_params%corr_w_eta_1(nz), &
pdf_params%corr_w_eta_2(nz), &
pdf_params%corr_chi_eta_1(nz), &
pdf_params%corr_chi_eta_2(nz), &
pdf_params%rsatl_1(nz), &
pdf_params%rsatl_2(nz), &
pdf_params%rc_1(nz), &
pdf_params%rc_2(nz), &
pdf_params%cloud_frac_1(nz), &
pdf_params%cloud_frac_2(nz), &
pdf_params%mixt_frac(nz), &
pdf_params%ice_supersat_frac_1(nz), &
pdf_params%ice_supersat_frac_2(nz) )
pdf_params%w_1 = zero
pdf_params%w_2 = zero
pdf_params%varnce_w_1 = zero
pdf_params%varnce_w_2 = zero
pdf_params%rt_1 = zero
pdf_params%rt_2 = zero
pdf_params%varnce_rt_1 = zero
pdf_params%varnce_rt_2 = zero
pdf_params%thl_1 = zero
pdf_params%thl_2 = zero
pdf_params%varnce_thl_1 = zero
pdf_params%varnce_thl_2 = zero
pdf_params%corr_w_rt_1 = zero
pdf_params%corr_w_rt_2 = zero
pdf_params%corr_w_thl_1 = zero
pdf_params%corr_w_thl_2 = zero
pdf_params%corr_rt_thl_1 = zero
pdf_params%corr_rt_thl_2 = zero
pdf_params%alpha_thl = zero
pdf_params%alpha_rt = zero
pdf_params%crt_1 = zero
pdf_params%crt_2 = zero
pdf_params%cthl_1 = zero
pdf_params%cthl_2 = zero
pdf_params%chi_1 = zero
pdf_params%chi_2 = zero
pdf_params%stdev_chi_1 = zero
pdf_params%stdev_chi_2 = zero
pdf_params%stdev_eta_1 = zero
pdf_params%stdev_eta_2 = zero
pdf_params%covar_chi_eta_1 = zero
pdf_params%covar_chi_eta_2 = zero
pdf_params%corr_w_chi_1 = zero
pdf_params%corr_w_chi_2 = zero
pdf_params%corr_w_eta_1 = zero
pdf_params%corr_w_eta_2 = zero
pdf_params%corr_chi_eta_1 = zero
pdf_params%corr_chi_eta_2 = zero
pdf_params%rsatl_1 = zero
pdf_params%rsatl_2 = zero
pdf_params%rc_1 = zero
pdf_params%rc_2 = zero
pdf_params%cloud_frac_1 = zero
pdf_params%cloud_frac_2 = zero
pdf_params%mixt_frac = zero
pdf_params%ice_supersat_frac_1 = zero
pdf_params%ice_supersat_frac_2 = zero
return
end subroutine init_pdf_params
!=============================================================================
subroutine init_pdf_implicit_coefs_terms( nz, sclr_dim, &
pdf_implicit_coefs_terms )
! Description:
! Initializes all PDF implicit coefficients and explicit terms in the
! variable type implicit_coefs_terms.
! References:
!--------------------------------------------------------------------
use constants_clubb, only: &
zero ! Constant(s)
implicit none
! Input Variables
integer, intent(in) :: &
nz, & ! Number of vertical grid levels [-]
sclr_dim ! Number of scalar variables [-]
! Output Variable
type(implicit_coefs_terms), intent(out) :: &
pdf_implicit_coefs_terms ! Implicit coefs / explicit terms [units vary]
! Allocate pdf_implicit_coefs_terms
allocate( pdf_implicit_coefs_terms%coef_wp4_implicit(1:nz), &
pdf_implicit_coefs_terms%coef_wp2rtp_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wp2rtp_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wp2thlp_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wp2thlp_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wp2up_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wp2up_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wp2vp_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wp2vp_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wprtp2_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wprtp2_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wpthlp2_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wpthlp2_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wprtpthlp_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wprtpthlp_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wpup2_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wpup2_explicit(1:nz), &
pdf_implicit_coefs_terms%coef_wpvp2_implicit(1:nz), &
pdf_implicit_coefs_terms%term_wpvp2_explicit(1:nz) )
if ( sclr_dim > 0 ) then
allocate( &
pdf_implicit_coefs_terms%coef_wp2sclrp_implicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%term_wp2sclrp_explicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%coef_wpsclrp2_implicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%term_wpsclrp2_explicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%coef_wprtpsclrp_implicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%term_wprtpsclrp_explicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%coef_wpthlpsclrp_implicit(1:nz,1:sclr_dim), &
pdf_implicit_coefs_terms%term_wpthlpsclrp_explicit(1:nz,1:sclr_dim) )
endif ! sclr_dim > 0
! Initialize pdf_implicit_coefs_terms
pdf_implicit_coefs_terms%coef_wp4_implicit = zero
pdf_implicit_coefs_terms%coef_wp2rtp_implicit = zero
pdf_implicit_coefs_terms%term_wp2rtp_explicit = zero
pdf_implicit_coefs_terms%coef_wp2thlp_implicit = zero
pdf_implicit_coefs_terms%term_wp2thlp_explicit = zero
pdf_implicit_coefs_terms%coef_wp2up_implicit = zero
pdf_implicit_coefs_terms%term_wp2up_explicit = zero
pdf_implicit_coefs_terms%coef_wp2vp_implicit = zero
pdf_implicit_coefs_terms%term_wp2vp_explicit = zero
pdf_implicit_coefs_terms%coef_wprtp2_implicit = zero
pdf_implicit_coefs_terms%term_wprtp2_explicit = zero
pdf_implicit_coefs_terms%coef_wpthlp2_implicit = zero
pdf_implicit_coefs_terms%term_wpthlp2_explicit = zero
pdf_implicit_coefs_terms%coef_wprtpthlp_implicit = zero
pdf_implicit_coefs_terms%term_wprtpthlp_explicit = zero
pdf_implicit_coefs_terms%coef_wpup2_implicit = zero
pdf_implicit_coefs_terms%term_wpup2_explicit = zero
pdf_implicit_coefs_terms%coef_wpvp2_implicit = zero
pdf_implicit_coefs_terms%term_wpvp2_explicit = zero
if ( sclr_dim > 0 ) then
pdf_implicit_coefs_terms%coef_wp2sclrp_implicit = zero
pdf_implicit_coefs_terms%term_wp2sclrp_explicit = zero
pdf_implicit_coefs_terms%coef_wpsclrp2_implicit = zero
pdf_implicit_coefs_terms%term_wpsclrp2_explicit = zero
pdf_implicit_coefs_terms%coef_wprtpsclrp_implicit = zero
pdf_implicit_coefs_terms%term_wprtpsclrp_explicit = zero
pdf_implicit_coefs_terms%coef_wpthlpsclrp_implicit = zero
pdf_implicit_coefs_terms%term_wpthlpsclrp_explicit = zero
endif ! sclr_dim > 0
return
end subroutine init_pdf_implicit_coefs_terms
!=============================================================================
! The CLUBB_CAM preprocessor directives are being commented out because this
! code is now also used for WRF-CLUBB.
!#ifdef CLUBB_CAM /* Code for storing pdf_parameter structs in pbuf as array */
subroutine pack_pdf_params(pdf_params, nz, r_param_array, &
k_start_in, k_end_in )
implicit none
integer, intent(in) :: nz ! Num Vert Model Levs
! Input a pdf_parameter array with nz instances of pdf_parameter
type (pdf_parameter), intent(in) :: pdf_params
! Output a two dimensional real array with all values
real (kind = core_rknd), dimension(nz,num_pdf_params), intent(out) :: &
r_param_array
integer, optional, intent(in) :: k_start_in, k_end_in
integer :: k_start, k_end
if( present( k_start_in ) .and. present( k_end_in ) ) then
k_start = k_start_in
k_end = k_end_in
else
k_start = 1
k_end = nz
end if
r_param_array(:,1) = pdf_params%w_1(k_start:k_end)
r_param_array(:,2) = pdf_params%w_2(k_start:k_end)
r_param_array(:,3) = pdf_params%varnce_w_1(k_start:k_end)
r_param_array(:,4) = pdf_params%varnce_w_2(k_start:k_end)
r_param_array(:,5) = pdf_params%rt_1(k_start:k_end)
r_param_array(:,6) = pdf_params%rt_2(k_start:k_end)
r_param_array(:,7) = pdf_params%varnce_rt_1(k_start:k_end)
r_param_array(:,8) = pdf_params%varnce_rt_2(k_start:k_end)
r_param_array(:,9) = pdf_params%thl_1(k_start:k_end)
r_param_array(:,10) = pdf_params%thl_2(k_start:k_end)
r_param_array(:,11) = pdf_params%varnce_thl_1(k_start:k_end)
r_param_array(:,12) = pdf_params%varnce_thl_2(k_start:k_end)
r_param_array(:,13) = pdf_params%corr_w_rt_1(k_start:k_end)
r_param_array(:,14) = pdf_params%corr_w_rt_2(k_start:k_end)
r_param_array(:,15) = pdf_params%corr_w_thl_1(k_start:k_end)
r_param_array(:,16) = pdf_params%corr_w_thl_2(k_start:k_end)
r_param_array(:,17) = pdf_params%corr_rt_thl_1(k_start:k_end)
r_param_array(:,18) = pdf_params%corr_rt_thl_2(k_start:k_end)
r_param_array(:,19) = pdf_params%alpha_thl(k_start:k_end)
r_param_array(:,20) = pdf_params%alpha_rt(k_start:k_end)
r_param_array(:,21) = pdf_params%crt_1(k_start:k_end)
r_param_array(:,22) = pdf_params%crt_2(k_start:k_end)
r_param_array(:,23) = pdf_params%cthl_1(k_start:k_end)
r_param_array(:,24) = pdf_params%cthl_2(k_start:k_end)
r_param_array(:,25) = pdf_params%chi_1(k_start:k_end)
r_param_array(:,26) = pdf_params%chi_2(k_start:k_end)
r_param_array(:,27) = pdf_params%stdev_chi_1(k_start:k_end)
r_param_array(:,28) = pdf_params%stdev_chi_2(k_start:k_end)
r_param_array(:,29) = pdf_params%stdev_eta_1(k_start:k_end)
r_param_array(:,30) = pdf_params%stdev_eta_2(k_start:k_end)
r_param_array(:,31) = pdf_params%covar_chi_eta_1(k_start:k_end)
r_param_array(:,32) = pdf_params%covar_chi_eta_2(k_start:k_end)
r_param_array(:,33) = pdf_params%corr_w_chi_1(k_start:k_end)
r_param_array(:,34) = pdf_params%corr_w_chi_2(k_start:k_end)
r_param_array(:,35) = pdf_params%corr_w_eta_1(k_start:k_end)
r_param_array(:,36) = pdf_params%corr_w_eta_2(k_start:k_end)
r_param_array(:,37) = pdf_params%corr_chi_eta_1(k_start:k_end)
r_param_array(:,38) = pdf_params%corr_chi_eta_2(k_start:k_end)
r_param_array(:,39) = pdf_params%rsatl_1(k_start:k_end)
r_param_array(:,40) = pdf_params%rsatl_2(k_start:k_end)
r_param_array(:,41) = pdf_params%rc_1(k_start:k_end)
r_param_array(:,42) = pdf_params%rc_2(k_start:k_end)
r_param_array(:,43) = pdf_params%cloud_frac_1(k_start:k_end)
r_param_array(:,44) = pdf_params%cloud_frac_2(k_start:k_end)
r_param_array(:,45) = pdf_params%mixt_frac(k_start:k_end)
r_param_array(:,46) = pdf_params%ice_supersat_frac_1(k_start:k_end)
r_param_array(:,47) = pdf_params%ice_supersat_frac_2(k_start:k_end)
end subroutine pack_pdf_params
subroutine unpack_pdf_params(r_param_array, nz, pdf_params, &
k_start_in, k_end_in )
implicit none
integer, intent(in) :: nz ! Num Vert Model Levs
! Input a two dimensional real array with pdf values
real (kind = core_rknd), dimension(nz,num_pdf_params), intent(in) :: &
r_param_array
! Output a pdf_parameter array with nz instances of pdf_parameter
type (pdf_parameter), intent(inout) :: pdf_params
integer, optional, intent(in) :: k_start_in, k_end_in
integer :: k_start, k_end
if( present( k_start_in ) .and. present( k_end_in ) ) then
k_start = k_start_in
k_end = k_end_in
else
k_start = 1
k_end = nz
end if
pdf_params%w_1(k_start:k_end) = r_param_array(:,1)
pdf_params%w_2(k_start:k_end) = r_param_array(:,2)
pdf_params%varnce_w_1(k_start:k_end) = r_param_array(:,3)
pdf_params%varnce_w_2(k_start:k_end) = r_param_array(:,4)
pdf_params%rt_1(k_start:k_end) = r_param_array(:,5)
pdf_params%rt_2(k_start:k_end) = r_param_array(:,6)
pdf_params%varnce_rt_1(k_start:k_end) = r_param_array(:,7)
pdf_params%varnce_rt_2(k_start:k_end)= r_param_array(:,8)
pdf_params%thl_1(k_start:k_end) = r_param_array(:,9)
pdf_params%thl_2(k_start:k_end) = r_param_array(:,10)
pdf_params%varnce_thl_1(k_start:k_end) = r_param_array(:,11)
pdf_params%varnce_thl_2(k_start:k_end) = r_param_array(:,12)
pdf_params%corr_w_rt_1(k_start:k_end) = r_param_array(:,13)
pdf_params%corr_w_rt_2(k_start:k_end) = r_param_array(:,14)
pdf_params%corr_w_thl_1(k_start:k_end) = r_param_array(:,15)
pdf_params%corr_w_thl_2(k_start:k_end) = r_param_array(:,16)
pdf_params%corr_rt_thl_1(k_start:k_end) = r_param_array(:,17)
pdf_params%corr_rt_thl_2(k_start:k_end) = r_param_array(:,18)
pdf_params%alpha_thl(k_start:k_end) = r_param_array(:,19)
pdf_params%alpha_rt(k_start:k_end) = r_param_array(:,20)
pdf_params%crt_1(k_start:k_end) = r_param_array(:,21)
pdf_params%crt_2(k_start:k_end) = r_param_array(:,22)
pdf_params%cthl_1(k_start:k_end) = r_param_array(:,23)
pdf_params%cthl_2(k_start:k_end) = r_param_array(:,24)
pdf_params%chi_1(k_start:k_end) = r_param_array(:,25)
pdf_params%chi_2(k_start:k_end) = r_param_array(:,26)
pdf_params%stdev_chi_1(k_start:k_end) = r_param_array(:,27)
pdf_params%stdev_chi_2(k_start:k_end) = r_param_array(:,28)
pdf_params%stdev_eta_1(k_start:k_end) = r_param_array(:,29)
pdf_params%stdev_eta_2(k_start:k_end) = r_param_array(:,30)
pdf_params%covar_chi_eta_1(k_start:k_end) = r_param_array(:,31)
pdf_params%covar_chi_eta_2(k_start:k_end) = r_param_array(:,32)
pdf_params%corr_w_chi_1(k_start:k_end) = r_param_array(:,33)
pdf_params%corr_w_chi_2(k_start:k_end) = r_param_array(:,34)
pdf_params%corr_w_eta_1(k_start:k_end) = r_param_array(:,35)
pdf_params%corr_w_eta_2(k_start:k_end) = r_param_array(:,36)
pdf_params%corr_chi_eta_1(k_start:k_end) = r_param_array(:,37)
pdf_params%corr_chi_eta_2(k_start:k_end) = r_param_array(:,38)
pdf_params%rsatl_1(k_start:k_end) = r_param_array(:,39)
pdf_params%rsatl_2(k_start:k_end) = r_param_array(:,40)
pdf_params%rc_1(k_start:k_end) = r_param_array(:,41)
pdf_params%rc_2(k_start:k_end) = r_param_array(:,42)
pdf_params%cloud_frac_1(k_start:k_end) = r_param_array(:,43)
pdf_params%cloud_frac_2(k_start:k_end) = r_param_array(:,44)
pdf_params%mixt_frac(k_start:k_end) = r_param_array(:,45)
pdf_params%ice_supersat_frac_1(k_start:k_end) = r_param_array(:,46)
pdf_params%ice_supersat_frac_2(k_start:k_end) = r_param_array(:,47)
end subroutine unpack_pdf_params
!#endif /* CLUBB_CAM */
end module pdf_parameter_module