Line data Source code
1 : !-----------------------------------------------------------------------
2 : !$Id$
3 : !-------------------------------------------------------------------------------
4 : module corr_varnce_module
5 :
6 : use clubb_precision, only: &
7 : core_rknd
8 :
9 : use array_index, only: &
10 : iiPDF_chi, &
11 : iiPDF_eta, &
12 : iiPDF_w, &
13 : iiPDF_rr, &
14 : iiPDF_rs, &
15 : iiPDF_ri, &
16 : iiPDF_rg, &
17 : iiPDF_Nr, &
18 : iiPDF_Ns, &
19 : iiPDF_Ni, &
20 : iiPDF_Ng, &
21 : iiPDF_Ncn
22 :
23 : use error_code, only: &
24 : clubb_at_least_debug_level, & ! Procedure
25 : clubb_fatal_error, & ! Constant
26 : err_code ! Error indicator
27 :
28 : implicit none
29 :
30 : type hmp2_ip_on_hmm2_ip_ratios_type
31 :
32 : ! Prescribed parameters for hydrometeor values of <hm|_ip'^2> / <hm|_ip>^2,
33 : ! where <hm|_ip'^2> is the in-precip. variance of the hydrometeor and
34 : ! <hm|_ip> is the in-precip. mean of the hydrometeor
35 : !
36 : ! These values are dependent on the horizontal grid spacing of the run, and are calculated
37 : ! using a slope and intercept corresponding to each hydrometer.
38 : real( kind = core_rknd ) :: &
39 : rr = 1.0_core_rknd, & ! For rain water mixing ratio [-]
40 : Nr = 1.0_core_rknd, & ! For rain drop concentration [-]
41 : ri = 1.0_core_rknd, & ! For ice mixing ratio [-]
42 : Ni = 1.0_core_rknd, & ! For ice concentration [-]
43 : rs = 1.0_core_rknd, & ! For snow mixing ratio [-]
44 : Ns = 1.0_core_rknd, & ! For snow concentration [-]
45 : rg = 1.0_core_rknd, & ! For graupel mixing ratio [-]
46 : Ng = 1.0_core_rknd ! For graupel concentration [-]
47 :
48 : end type hmp2_ip_on_hmm2_ip_ratios_type
49 :
50 : ! These slopes and intercepts below are used to calculate the hmp2_ip_on_hmm2_ip_ratios_type
51 : ! values that are defined above. This functionality is described by equations 8, 10, & 11
52 : ! in ``Parameterization of the Spatial Variability of Rain for Large-Scale Models and
53 : ! Remote Sensing`` Lebo, et al, October 2015
54 : ! https://journals.ametsoc.org/doi/pdf/10.1175/JAMC-D-15-0066.1
55 : ! see clubb:ticket:830 for more detail
56 : !
57 : ! hmp2_ip_on_hmm2_ip(iirr) = hmp2_ip_on_hmm2_ip_intrcpt%rr + &
58 : ! hmp2_ip_on_hmm2_ip_slope%rr * max( host_dx, host_dy )
59 : !
60 : ! In Lebo et al. the suggested values were
61 : ! slope = 2.12e-5 [1/m]
62 : ! intercept = 0.54 [-]
63 : !
64 : ! In CLUBB standalone, these parameters can be set based on the value for a
65 : ! given case in the CASE_model.in file.
66 : type hmp2_ip_on_hmm2_ip_slope_type
67 :
68 : real( kind = core_rknd ) :: &
69 : rr = 2.12e-5_core_rknd, & ! For rain water mixing ratio [1/m]
70 : Nr = 2.12e-5_core_rknd, & ! For rain drop concentration [1/m]
71 : ri = 2.12e-5_core_rknd, & ! For ice mixing ratio [1/m]
72 : Ni = 2.12e-5_core_rknd, & ! For ice concentration [1/m]
73 : rs = 2.12e-5_core_rknd, & ! For snow mixing ratio [1/m]
74 : Ns = 2.12e-5_core_rknd, & ! For snow concentration [1/m]
75 : rg = 2.12e-5_core_rknd, & ! For graupel mixing ratio [1/m]
76 : Ng = 2.12e-5_core_rknd ! For graupel concentration [1/m]
77 :
78 : end type hmp2_ip_on_hmm2_ip_slope_type
79 :
80 : type hmp2_ip_on_hmm2_ip_intrcpt_type
81 :
82 : real( kind = core_rknd ) :: &
83 : rr = 0.54_core_rknd, & ! For rain water mixing ratio [-]
84 : Nr = 0.54_core_rknd, & ! For rain drop concentration [-]
85 : ri = 0.54_core_rknd, & ! For ice mixing ratio [-]
86 : Ni = 0.54_core_rknd, & ! For ice concentration [-]
87 : rs = 0.54_core_rknd, & ! For snow mixing ratio [-]
88 : Ns = 0.54_core_rknd, & ! For snow concentration [-]
89 : rg = 0.54_core_rknd, & ! For graupel mixing ratio [-]
90 : Ng = 0.54_core_rknd ! For graupel concentration [-]
91 :
92 : end type hmp2_ip_on_hmm2_ip_intrcpt_type
93 :
94 :
95 : ! Prescribed parameter for <N_cn'^2> / <N_cn>^2.
96 : ! NOTE: In the case that l_const_Nc_in_cloud is true, Ncn is constant
97 : ! throughout the entire grid box, so the parameter below should be
98 : ! ignored.
99 : real( kind = core_rknd ), public :: &
100 : Ncnp2_on_Ncnm2 = 1.0_core_rknd ! Prescribed ratio <N_cn'^2> / <N_cn>^2 [-]
101 :
102 : !$omp threadprivate(Ncnp2_on_Ncnm2)
103 :
104 : integer, parameter, public :: &
105 : d_var_total = 12 ! Size of the default correlation arrays
106 :
107 : integer, public :: &
108 : pdf_dim
109 : !$omp threadprivate(pdf_dim)
110 :
111 : real( kind = core_rknd ), dimension(:), allocatable, public :: &
112 : hmp2_ip_on_hmm2_ip
113 :
114 : !$omp threadprivate(hmp2_ip_on_hmm2_ip)
115 :
116 : real( kind = core_rknd ), public, dimension(:,:), allocatable :: &
117 : corr_array_n_cloud, &
118 : corr_array_n_below
119 : !$omp threadprivate(corr_array_n_cloud, corr_array_n_below)
120 :
121 : real( kind = core_rknd ), public, dimension(:,:), allocatable :: &
122 : corr_array_n_cloud_def, &
123 : corr_array_n_below_def
124 : !$omp threadprivate( corr_array_n_cloud_def, corr_array_n_below_def )
125 :
126 : private
127 :
128 : public :: hmp2_ip_on_hmm2_ip_ratios_type, &
129 : hmp2_ip_on_hmm2_ip_slope_type, &
130 : hmp2_ip_on_hmm2_ip_intrcpt_type, &
131 : read_correlation_matrix, init_pdf_indices, &
132 : setup_corr_varnce_array, cleanup_corr_matrix_arrays, &
133 : assert_corr_symmetric, print_corr_matrix, init_hydromet_arrays
134 :
135 : private :: get_corr_var_index, def_corr_idx
136 :
137 :
138 : contains
139 :
140 : !-----------------------------------------------------------------------------
141 0 : subroutine init_default_corr_arrays( )
142 :
143 : ! Description:
144 : ! Initializes the default correlation arrays.
145 : !---------------------------------------------------------------------------
146 :
147 : use constants_clubb, only: &
148 : one, & ! Constant(s)
149 : zero
150 :
151 : implicit none
152 :
153 : integer:: indx
154 :
155 : ! This "renaming" is used to shorten the matrix declarations below.
156 : integer, parameter :: c = core_rknd
157 :
158 : ! ---- Begin Code ----
159 :
160 : ! Allocate Arrays.
161 0 : allocate( corr_array_n_cloud_def(d_var_total,d_var_total) )
162 0 : allocate( corr_array_n_below_def(d_var_total,d_var_total) )
163 :
164 : ! Initialize all values to 0.
165 0 : corr_array_n_cloud_def = zero
166 0 : corr_array_n_below_def = zero
167 :
168 : ! Set the correlation of any variable with itself to 1.
169 0 : do indx = 1, d_var_total, 1
170 0 : corr_array_n_cloud_def(indx,indx) = one
171 0 : corr_array_n_below_def(indx,indx) = one
172 : enddo
173 :
174 : ! Set up default normal space correlation arrays.
175 : ! The default normal space correlation arrays used here are the normal space
176 : ! correlation arrays used for the ARM 97 case. Any changes should be made
177 : ! concurrently here and in
178 : ! ../../input/case_setups/arm_97_corr_array_cloud.in (for "in-cloud") and
179 : ! in ../../input/case_setups/arm_97_corr_array_cloud.in (for "below-cloud").
180 : corr_array_n_cloud_def = reshape( &
181 :
182 : (/1._c,-.6_c, .09_c , .09_c , .788_c, .675_c, .240_c, .222_c, .240_c, .222_c, .240_c, .222_c, &! chi
183 : 0._c, 1._c, .027_c, .027_c, .114_c, .115_c,-.029_c, .093_c, .022_c, .013_c, 0._c , 0._c , &! eta
184 : 0._c, 0._c, 1._c , .34_c , .315_c, .270_c, .120_c, .167_c, 0._c , 0._c , 0._c , 0._c , &! w
185 : 0._c, 0._c, 0._c , 1._c , 0._c , 0._c , .464_c, .320_c, .168_c, .232_c, 0._c , 0._c , &! Ncn
186 : 0._c, 0._c, 0._c , 0._c , 1._c , .821_c, 0._c , 0._c , .173_c, .164_c, .319_c, .308_c, &! rr
187 : 0._c, 0._c, 0._c , 0._c , 0._c , 1._c , .152_c, .143_c, 0._c , 0._c , .285_c, .273_c, &! Nr
188 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 1._c , .758_c, .585_c, .571_c, .379_c, .363_c, &! ri
189 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .571_c, .550_c, .363_c, .345_c, &! Ni
190 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .758_c, .485_c, .470_c, &! rs
191 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .470_c, .450_c, &! Ns
192 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .758_c, &! rg
193 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c/), &! Ng
194 :
195 0 : shape(corr_array_n_cloud_def) )
196 : ! chi eta w Ncn rr Nr ri Ni rs Ns rg Ng
197 :
198 0 : corr_array_n_cloud_def = transpose( corr_array_n_cloud_def )
199 :
200 :
201 : corr_array_n_below_def = reshape( &
202 :
203 : (/1._c, .3_c, .09_c , .09_c , .788_c, .675_c, .240_c, .222_c, .240_c, .222_c, .240_c, .222_c, &! chi
204 : 0._c, 1._c, .027_c, .027_c, .114_c, .115_c,-.029_c, .093_c, .022_c, .013_c, 0._c , 0._c , &! eta
205 : 0._c, 0._c, 1._c , .34_c , .315_c, .270_c, .120_c, .167_c, 0._c , 0._c , 0._c , 0._c , &! w
206 : 0._c, 0._c, 0._c , 1._c , 0._c , 0._c , .464_c, .320_c, .168_c, .232_c, 0._c , 0._c , &! Ncn
207 : 0._c, 0._c, 0._c , 0._c , 1._c , .821_c, 0._c , 0._c , .173_c, .164_c, .319_c, .308_c, &! rr
208 : 0._c, 0._c, 0._c , 0._c , 0._c , 1._c , .152_c, .143_c, 0._c , 0._c , .285_c, .273_c, &! Nr
209 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 1._c , .758_c, .585_c, .571_c, .379_c, .363_c, &! ri
210 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .571_c, .550_c, .363_c, .345_c, &! Ni
211 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .758_c, .485_c, .470_c, &! rs
212 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .470_c, .450_c, &! Ns
213 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c , .758_c, &! rg
214 : 0._c, 0._c, 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 0._c , 1._c/), &! Ng
215 :
216 0 : shape(corr_array_n_below_def) )
217 : ! chi eta w Ncn rr Nr ri Ni rs Ns rg Ng
218 :
219 0 : corr_array_n_below_def = transpose( corr_array_n_below_def )
220 :
221 :
222 0 : return
223 :
224 : end subroutine init_default_corr_arrays
225 :
226 : !-----------------------------------------------------------------------------
227 0 : function def_corr_idx( iiPDF_x ) result(ii_def_corr)
228 :
229 : ! Description:
230 : ! Map from a iiPDF index to the corresponding index in the default
231 : ! correlation arrays.
232 : !-----------------------------------------------------------------------------
233 :
234 : implicit none
235 :
236 : ! Constant Parameters
237 :
238 : ! Indices that represent the order in the default corr arrays
239 : ! (chi (old s), eta (old t), w, Ncn, rr, Nr, ri, Ni, rs, Ns, rg, Ng)
240 : integer, parameter :: &
241 : ii_chi = 1, &
242 : ii_eta = 2, &
243 : ii_w = 3, &
244 : ii_Ncn = 4, &
245 : ii_rr = 5, &
246 : ii_Nr = 6, &
247 : ii_ri = 7, &
248 : ii_Ni = 8, &
249 : ii_rs = 9, &
250 : ii_Ns = 10, &
251 : ii_rg = 11, &
252 : ii_Ng = 12
253 :
254 : ! Input Variables
255 :
256 : integer, intent(in) :: iiPDF_x
257 :
258 : ! Return Variable
259 :
260 : integer :: ii_def_corr
261 :
262 : ! ---- Begin Code ----
263 :
264 0 : ii_def_corr = -1
265 :
266 0 : if (iiPDF_x == iiPDF_chi) then
267 : ii_def_corr = ii_chi
268 :
269 0 : elseif (iiPDF_x == iiPDF_eta) then
270 : ii_def_corr = ii_eta
271 :
272 0 : elseif (iiPDF_x == iiPDF_w) then
273 : ii_def_corr = ii_w
274 :
275 0 : elseif (iiPDF_x == iiPDF_Ncn) then
276 : ii_def_corr = ii_Ncn
277 :
278 0 : elseif (iiPDF_x == iiPDF_rr) then
279 : ii_def_corr = ii_rr
280 :
281 0 : elseif (iiPDF_x == iiPDF_Nr) then
282 : ii_def_corr = ii_Nr
283 :
284 0 : elseif (iiPDF_x == iiPDF_ri) then
285 : ii_def_corr = ii_ri
286 :
287 0 : elseif (iiPDF_x == iiPDF_Ni) then
288 : ii_def_corr = ii_Ni
289 :
290 0 : elseif (iiPDF_x == iiPDF_rs) then
291 : ii_def_corr = ii_rs
292 :
293 0 : elseif (iiPDF_x == iiPDF_Ns) then
294 : ii_def_corr = ii_Ns
295 :
296 0 : elseif (iiPDF_x == iiPDF_rg) then
297 : ii_def_corr = ii_rg
298 :
299 0 : elseif (iiPDF_x == iiPDF_Ng) then
300 0 : ii_def_corr = ii_Ng
301 :
302 : endif
303 0 : end function def_corr_idx
304 :
305 : !-----------------------------------------------------------------------------
306 0 : subroutine set_corr_arrays_to_default( )
307 :
308 : ! Description:
309 : ! If there are no corr_array.in files for the current case, default
310 : ! correlations are used.
311 : !-----------------------------------------------------------------------------
312 :
313 : use constants_clubb, only: &
314 : zero, &
315 : one
316 :
317 : implicit none
318 :
319 : ! Local Variables
320 : integer :: i, j ! Loop iterators
321 :
322 :
323 : ! ---- Begin Code ----
324 :
325 0 : corr_array_n_cloud = zero
326 0 : corr_array_n_below = zero
327 :
328 0 : do i = 1, pdf_dim
329 0 : corr_array_n_cloud(i,i) = one
330 0 : corr_array_n_below(i,i) = one
331 : enddo
332 :
333 0 : do i = 1, pdf_dim-1
334 0 : do j = i+1, pdf_dim
335 0 : if ( def_corr_idx(i) > def_corr_idx(j) ) then
336 0 : corr_array_n_cloud(j, i) = corr_array_n_cloud_def(def_corr_idx(j), def_corr_idx(i))
337 0 : corr_array_n_below(j, i) = corr_array_n_below_def(def_corr_idx(j), def_corr_idx(i))
338 : else
339 0 : corr_array_n_cloud(j, i) = corr_array_n_cloud_def(def_corr_idx(i), def_corr_idx(j))
340 0 : corr_array_n_below(j, i) = corr_array_n_below_def(def_corr_idx(i), def_corr_idx(j))
341 : endif
342 : enddo
343 : enddo
344 :
345 0 : end subroutine set_corr_arrays_to_default
346 :
347 :
348 : !-----------------------------------------------------------------------------
349 0 : subroutine read_correlation_matrix( iunit, input_file, pdf_dim, &
350 0 : corr_array_n )
351 :
352 : ! Description:
353 : ! Reads a correlation variance array from a file and stores it in an array.
354 : !-----------------------------------------------------------------------------
355 :
356 : use input_reader, only: &
357 : one_dim_read_var, & ! Variable(s)
358 : read_one_dim_file, deallocate_one_dim_vars, count_columns ! Procedure(s)
359 :
360 : use matrix_operations, only: set_lower_triangular_matrix ! Procedure(s)
361 :
362 : use constants_clubb, only: fstderr ! Variable(s)
363 :
364 : use clubb_precision, only: &
365 : core_rknd ! Variable(s)
366 :
367 : implicit none
368 :
369 : ! Input Variable(s)
370 : integer, intent(in) :: &
371 : iunit, & ! File I/O unit
372 : pdf_dim! number of variables in the array
373 :
374 : character(len=*), intent(in) :: input_file ! Path to the file
375 :
376 : ! Input/Output Variable(s)
377 : real( kind = core_rknd ), dimension(pdf_dim,pdf_dim), intent(inout) :: &
378 : corr_array_n ! Normal space correlation array
379 :
380 : ! Local Variable(s)
381 :
382 : type(one_dim_read_var), allocatable, dimension(:) :: &
383 0 : retVars ! stores the variables read in from the corr_varnce.in file
384 :
385 : integer :: &
386 : var_index1, & ! variable index
387 : var_index2, & ! variable index
388 : nCols, & ! the number of columns in the file
389 : i, j ! Loop index
390 :
391 :
392 : !--------------------------- BEGIN CODE -------------------------
393 :
394 0 : nCols = count_columns( iunit, input_file )
395 :
396 : ! Allocate all arrays based on pdf_dim
397 0 : allocate( retVars(1:nCols) )
398 :
399 : ! Initializing to zero means that correlations we don't have are assumed to be 0.
400 0 : corr_array_n(:,:) = 0.0_core_rknd
401 :
402 : ! Set main diagonal to 1
403 0 : do i=1, pdf_dim
404 0 : corr_array_n(i,i) = 1.0_core_rknd
405 : end do
406 :
407 : ! Read the values from the specified file
408 : call read_one_dim_file( iunit, nCols, input_file, & ! intent(in)
409 0 : retVars ) ! intent(out)
410 :
411 0 : if( size( retVars(1)%values ) /= nCols ) then
412 0 : write(fstderr, *) "Correlation matrix must have an equal number of rows and cols in file ", &
413 0 : input_file
414 0 : error stop "Bad data in correlation file."
415 : end if
416 :
417 : ! Start at 2 because the first index is always just 1.0 in the first row
418 : ! and the rest of the rows are ignored
419 0 : do i=2, nCols
420 0 : var_index1 = get_corr_var_index( retVars(i)%name )
421 0 : if( var_index1 > -1 ) then
422 0 : do j=1, (i-1)
423 0 : var_index2 = get_corr_var_index( retVars(j)%name )
424 0 : if( var_index2 > -1 ) then
425 : call set_lower_triangular_matrix &
426 0 : ( pdf_dim, var_index1, var_index2, retVars(i)%values(j), & ! intent(in)
427 0 : corr_array_n ) ! intent(inout)
428 : end if
429 : end do
430 : end if
431 : end do
432 :
433 : call deallocate_one_dim_vars( nCols, & ! intent(in)
434 0 : retVars ) ! intent(inout)
435 :
436 0 : return
437 0 : end subroutine read_correlation_matrix
438 :
439 : !--------------------------------------------------------------------------
440 0 : function get_corr_var_index( var_name ) result( i )
441 :
442 : ! Definition:
443 : ! Returns the index for a variable based on its name.
444 : !--------------------------------------------------------------------------
445 :
446 : implicit none
447 :
448 : character(len=*), intent(in) :: var_name ! The name of the variable
449 :
450 : ! Output variable
451 : integer :: i
452 :
453 : !------------------ BEGIN CODE -----------------------------
454 0 : i = -1
455 :
456 0 : select case( trim(var_name) )
457 :
458 : case( "chi" )
459 0 : i = iiPDF_chi
460 :
461 : case( "eta" )
462 0 : i = iiPDF_eta
463 :
464 : case( "w" )
465 0 : i = iiPDF_w
466 :
467 : case( "Ncn" )
468 0 : i = iiPDF_Ncn
469 :
470 : case( "rr" )
471 0 : i = iiPDF_rr
472 :
473 : case( "Nr" )
474 0 : i = iiPDF_Nr
475 :
476 : case( "ri" )
477 0 : i = iiPDF_ri
478 :
479 : case( "Ni" )
480 0 : i = iiPDF_Ni
481 :
482 : case( "rs" )
483 0 : i = iiPDF_rs
484 :
485 : case( "Ns" )
486 0 : i = iiPDF_Ns
487 :
488 : case( "rg" )
489 0 : i = iiPDF_rg
490 :
491 : case( "Ng" )
492 0 : i = iiPDF_Ng
493 :
494 : end select
495 :
496 : return
497 :
498 : end function get_corr_var_index
499 :
500 : !===============================================================================
501 0 : subroutine init_pdf_indices( hydromet_dim, iirr, iiNr, & ! intent(in)
502 : iiri, iiNi, iirs, iiNs, & ! intent(in)
503 : iirg, iiNg ) ! intent(in)
504 :
505 : ! Description:
506 : !
507 : ! Setup for the iiPDF indices. These indices are used to address
508 : ! chi(s), eta(t), w and the hydrometeors in the mean/stdev/corr arrays
509 : !
510 : ! References:
511 : !-------------------------------------------------------------------------------
512 :
513 :
514 :
515 : implicit none
516 :
517 : ! Input Variables
518 : integer, intent(in) :: &
519 : hydromet_dim, & ! Total number of hydrometeor species.
520 : iirr, & ! Index of rain water mixing ratio
521 : iiNr, & ! Index of rain drop concentration
522 : iiri, & ! Index of ice mixing ratio
523 : iiNi, & ! Index of ice crystal concentration
524 : iirs, & ! Index of snow mixing ratio
525 : iiNs, & ! Index of snow flake concentration
526 : iirg, & ! Index of graupel mixing ratio
527 : iiNg ! Index of graupel concentration
528 :
529 : ! Local Variables
530 : integer :: &
531 : pdf_count, & ! Count number of PDF variables
532 : i ! Hydrometeor loop index
533 :
534 : !--------------------- Begin Code --------------------------------
535 :
536 0 : iiPDF_chi = 1 ! Extended liquid water mixing ratio, chi
537 0 : iiPDF_eta = 2 ! 'eta' orthogonal to 'chi'
538 0 : iiPDF_w = 3 ! vertical velocity
539 0 : iiPDF_Ncn = 4 ! Simplified cloud nuclei concentration or extended Nc.
540 :
541 0 : pdf_count = iiPDF_Ncn
542 :
543 : ! Loop over hydrometeors.
544 : ! Hydrometeor indices in the PDF arrays should be in the same order as
545 : ! found in the hydrometeor arrays.
546 0 : if ( hydromet_dim > 0 ) then
547 :
548 0 : do i = 1, hydromet_dim, 1
549 :
550 0 : if ( i == iirr ) then
551 0 : pdf_count = pdf_count + 1
552 0 : iiPDF_rr = pdf_count
553 : endif
554 :
555 0 : if ( i == iiNr ) then
556 0 : pdf_count = pdf_count + 1
557 0 : iiPDF_Nr = pdf_count
558 : endif
559 :
560 0 : if ( i == iiri ) then
561 0 : pdf_count = pdf_count + 1
562 0 : iiPDF_ri = pdf_count
563 : endif
564 :
565 0 : if ( i == iiNi ) then
566 0 : pdf_count = pdf_count + 1
567 0 : iiPDF_Ni = pdf_count
568 : endif
569 :
570 0 : if ( i == iirs ) then
571 0 : pdf_count = pdf_count + 1
572 0 : iiPDF_rs = pdf_count
573 : endif
574 :
575 0 : if ( i == iiNs ) then
576 0 : pdf_count = pdf_count + 1
577 0 : iiPDF_Ns = pdf_count
578 : endif
579 :
580 0 : if ( i == iirg ) then
581 0 : pdf_count = pdf_count + 1
582 0 : iiPDF_rg = pdf_count
583 : endif
584 :
585 0 : if ( i == iiNg ) then
586 0 : pdf_count = pdf_count + 1
587 0 : iiPDF_Ng = pdf_count
588 : endif
589 :
590 : enddo ! i = 1, hydromet_dim, 1
591 :
592 : endif ! hydromet_dim > 0
593 :
594 0 : pdf_dim = pdf_count
595 :
596 :
597 0 : return
598 :
599 : end subroutine init_pdf_indices
600 :
601 : !===============================================================================
602 0 : subroutine init_hydromet_arrays( hydromet_dim, iirr, iiNr, & ! intent(in)
603 : iiri, iiNi, iirs, iiNs, & ! intent(in)
604 : iirg, iiNg ) ! intent(in)
605 : ! Description:
606 : !
607 : ! Initialization for the hydromet arrays. How the arrays are initialized is
608 : ! determined by how the hydromet indicies (iirr, iiNr, etc.) have been set,
609 : ! which are determined by the microphysics scheme.
610 : !-------------------------------------------------------------------------------
611 :
612 : use constants_clubb, only: &
613 : rr_tol, & ! Constants, tolerances for each hydrometeor
614 : ri_tol, &
615 : rs_tol, &
616 : rg_tol, &
617 : Nr_tol, &
618 : Ni_tol, &
619 : Ns_tol, &
620 : Ng_tol
621 :
622 : use array_index, only: &
623 : hydromet_list, & ! Names of the hydrometeor species
624 : hydromet_tol, & ! List of tolerances for each enabled hydrometeor
625 : l_frozen_hm, & ! True means hydrometeor is frozen
626 : l_mix_rat_hm ! True means hydrometeor is a mixing ratio
627 :
628 : implicit none
629 :
630 : ! Input Variables
631 : integer, intent(in) :: &
632 : hydromet_dim, & ! Total number of hydrometeor species.
633 : iirr, & ! Index of rain water mixing ratio
634 : iiNr, & ! Index of rain drop concentration
635 : iiri, & ! Index of ice mixing ratio
636 : iiNi, & ! Index of ice crystal concentration
637 : iirs, & ! Index of snow mixing ratio
638 : iiNs, & ! Index of snow flake concentration
639 : iirg, & ! Index of graupel mixing ratio
640 : iiNg ! Index of graupel concentration
641 :
642 : !--------------------- Begin Code --------------------------------
643 :
644 :
645 : ! Set up predictive precipitating hydrometeor arrays.
646 0 : allocate( hydromet_list(hydromet_dim) )
647 0 : allocate( hydromet_tol(hydromet_dim) )
648 0 : allocate( l_mix_rat_hm(hydromet_dim) )
649 0 : allocate( l_frozen_hm(hydromet_dim) )
650 :
651 0 : if ( iirr > 0 ) then
652 : ! The microphysics scheme predicts rain water mixing ratio, rr.
653 0 : hydromet_list(iirr) = "rrm"
654 0 : l_mix_rat_hm(iirr) = .true.
655 0 : l_frozen_hm(iirr) = .false.
656 0 : hydromet_tol(iirr) = rr_tol
657 : endif
658 :
659 0 : if ( iiri > 0 ) then
660 : ! The microphysics scheme predicts ice mixing ratio, ri.
661 0 : hydromet_list(iiri) = "rim"
662 0 : l_mix_rat_hm(iiri) = .true.
663 0 : l_frozen_hm(iiri) = .true.
664 0 : hydromet_tol(iiri) = ri_tol
665 : endif
666 :
667 0 : if ( iirs > 0 ) then
668 : ! The microphysics scheme predicts snow mixing ratio, rs.
669 0 : hydromet_list(iirs) = "rsm"
670 0 : l_mix_rat_hm(iirs) = .true.
671 0 : l_frozen_hm(iirs) = .true.
672 0 : hydromet_tol(iirs) = rs_tol
673 : endif
674 :
675 0 : if ( iirg > 0 ) then
676 : ! The microphysics scheme predicts graupel mixing ratio, rg.
677 0 : hydromet_list(iirg) = "rgm"
678 0 : l_mix_rat_hm(iirg) = .true.
679 0 : l_frozen_hm(iirg) = .true.
680 0 : hydromet_tol(iirg) = rg_tol
681 : endif
682 :
683 0 : if ( iiNr > 0 ) then
684 : ! The microphysics scheme predicts rain drop concentration, Nr.
685 0 : hydromet_list(iiNr) = "Nrm"
686 0 : l_frozen_hm(iiNr) = .false.
687 0 : l_mix_rat_hm(iiNr) = .false.
688 0 : hydromet_tol(iiNr) = Nr_tol
689 : endif
690 :
691 0 : if ( iiNi > 0 ) then
692 : ! The microphysics scheme predicts ice concentration, Ni.
693 0 : hydromet_list(iiNi) = "Nim"
694 0 : l_mix_rat_hm(iiNi) = .false.
695 0 : l_frozen_hm(iiNi) = .true.
696 0 : hydromet_tol(iiNi) = Ni_tol
697 : endif
698 :
699 0 : if ( iiNs > 0 ) then
700 : ! The microphysics scheme predicts snowflake concentration, Ns.
701 0 : hydromet_list(iiNs) = "Nsm"
702 0 : l_mix_rat_hm(iiNs) = .false.
703 0 : l_frozen_hm(iiNs) = .true.
704 0 : hydromet_tol(iiNs) = Ns_tol
705 : endif
706 :
707 0 : if ( iiNg > 0 ) then
708 : ! The microphysics scheme predicts graupel concentration, Ng.
709 0 : hydromet_list(iiNg) = "Ngm"
710 0 : l_mix_rat_hm(iiNg) = .false.
711 0 : l_frozen_hm(iiNg) = .true.
712 0 : hydromet_tol(iiNg) = Ng_tol
713 : endif
714 :
715 0 : return
716 :
717 : end subroutine init_hydromet_arrays
718 :
719 : !===============================================================================
720 0 : subroutine setup_corr_varnce_array( input_file_cloud, input_file_below, &
721 : iunit, &
722 : l_fix_w_chi_eta_correlations )
723 :
724 : ! Description:
725 : ! Setup an array with the x'^2/xm^2 variables on the diagonal and the other
726 : ! elements to be correlations between various variables.
727 :
728 : ! References:
729 : ! None.
730 : !-------------------------------------------------------------------------------
731 :
732 : use matrix_operations, only: mirror_lower_triangular_matrix ! Procedure
733 :
734 : use constants_clubb, only: &
735 : fstderr, & ! Constant(s)
736 : eps
737 :
738 : implicit none
739 :
740 : ! External
741 : intrinsic :: max, epsilon, trim
742 :
743 : character(len=*), intent(in) :: &
744 : input_file_cloud, & ! Path to the in cloud correlation file
745 : input_file_below ! Path to the out of cloud correlation file
746 :
747 : ! Input Variables
748 : integer, intent(in) :: &
749 : iunit ! The file unit
750 :
751 : logical, intent(in) :: &
752 : l_fix_w_chi_eta_correlations ! Use a fixed correlation for s and t Mellor(chi/eta)
753 :
754 : ! Local variables
755 : logical :: l_warning, l_corr_file_1_exist, l_corr_file_2_exist
756 : integer :: i
757 :
758 : ! ---- Begin Code ----
759 :
760 0 : allocate( corr_array_n_cloud(pdf_dim,pdf_dim) )
761 0 : allocate( corr_array_n_below(pdf_dim,pdf_dim) )
762 :
763 0 : inquire( file = input_file_cloud, exist = l_corr_file_1_exist )
764 0 : inquire( file = input_file_below, exist = l_corr_file_2_exist )
765 :
766 0 : if ( l_corr_file_1_exist .and. l_corr_file_2_exist ) then
767 :
768 : call read_correlation_matrix( iunit, trim( input_file_cloud ), pdf_dim, & ! In
769 0 : corr_array_n_cloud ) ! Out
770 :
771 : call read_correlation_matrix( iunit, trim( input_file_below ), pdf_dim, & ! In
772 0 : corr_array_n_below ) ! Out
773 :
774 : else ! Read in default correlation matrices
775 :
776 0 : if ( clubb_at_least_debug_level( 1 ) ) then
777 : write(fstderr,*) "Warning: "//trim( input_file_cloud )//" was not found! " // &
778 0 : "The default correlation arrays will be used."
779 : end if
780 :
781 0 : call init_default_corr_arrays( )
782 :
783 0 : call set_corr_arrays_to_default( )
784 :
785 : endif
786 :
787 : ! Mirror the correlation matrices
788 : call mirror_lower_triangular_matrix( pdf_dim, & ! intent(in)
789 0 : corr_array_n_cloud ) ! intent(inout)
790 : call mirror_lower_triangular_matrix( pdf_dim, & ! intent(in)
791 0 : corr_array_n_below ) ! intent(inout)
792 :
793 : ! Sanity check to avoid confusing non-convergence results.
794 0 : if ( clubb_at_least_debug_level( 2 ) ) then
795 :
796 0 : if ( .not. l_fix_w_chi_eta_correlations .and. iiPDF_Ncn > 0 ) then
797 0 : l_warning = .false.
798 0 : do i = 1, pdf_dim
799 0 : if ( ( (abs(corr_array_n_cloud(i,iiPDF_Ncn)) > eps) .or. &
800 0 : (abs(corr_array_n_below(i,iiPDF_Ncn))) > eps) .and. &
801 0 : i /= iiPDF_Ncn ) then
802 0 : l_warning = .true.
803 : end if
804 : end do ! 1..pdf_dim
805 0 : if ( l_warning ) then
806 : write(fstderr,*) "Warning: the specified correlations for chi" &
807 0 : // " (old s) and Ncn are non-zero."
808 : write(fstderr,*) "The latin hypercube code will not converge to" &
809 0 : // " the analytic solution using these settings."
810 : end if
811 : end if ! l_fix_w_chi_eta_correlations .and. iiPDF_Ncn > 0
812 :
813 : end if ! clubb_at_least_debug_level( 2 )
814 :
815 :
816 0 : return
817 :
818 : end subroutine setup_corr_varnce_array
819 :
820 : !-----------------------------------------------------------------------------
821 0 : subroutine cleanup_corr_matrix_arrays( )
822 :
823 : ! Description:
824 : ! De-allocate latin hypercube arrays
825 : ! References:
826 : ! None
827 : !---------------------------------------------------------------------------
828 : implicit none
829 :
830 : ! External
831 : intrinsic :: allocated
832 :
833 : ! ---- Begin Code ----
834 :
835 0 : if ( allocated( corr_array_n_cloud ) ) then
836 0 : deallocate( corr_array_n_cloud )
837 : end if
838 :
839 0 : if ( allocated( corr_array_n_below ) ) then
840 0 : deallocate( corr_array_n_below )
841 : end if
842 :
843 0 : if ( allocated( corr_array_n_cloud_def ) ) then
844 0 : deallocate( corr_array_n_cloud_def )
845 : end if
846 :
847 0 : if ( allocated( corr_array_n_below_def ) ) then
848 0 : deallocate( corr_array_n_below_def )
849 : end if
850 :
851 :
852 0 : return
853 :
854 : end subroutine cleanup_corr_matrix_arrays
855 :
856 : !-----------------------------------------------------------------------------
857 0 : subroutine assert_corr_symmetric( corr_array_n, & ! intent(in)
858 : pdf_dim) ! intent(in)
859 :
860 : ! Description:
861 : ! Asserts that corr_matrix(i,j) == corr_matrix(j,i) for all indeces
862 : ! in the correlation array. If this is not the case, stops the program.
863 : ! References:
864 : ! None
865 : !---------------------------------------------------------------------------
866 :
867 : use constants_clubb, only: fstderr, eps, one ! Constant(s)
868 :
869 : implicit none
870 :
871 : ! Input Variables
872 : integer, intent(in) :: &
873 : pdf_dim ! Number of variables in the correlation array
874 :
875 : real( kind = core_rknd ), dimension(pdf_dim, pdf_dim), &
876 : intent(in) :: corr_array_n ! Normal space correlation array to be checked
877 :
878 : ! Local Variables
879 :
880 : ! tolerance used for real precision testing
881 : real( kind = core_rknd ), parameter :: tol = 1.0e-6_core_rknd
882 :
883 : integer :: n_row, n_col ! Indices
884 :
885 : !----- Begin Code -----
886 :
887 : !Do the check
888 0 : do n_col = 1, pdf_dim
889 0 : do n_row = 1, pdf_dim
890 :
891 0 : if ( abs(corr_array_n(n_col, n_row) - corr_array_n(n_row, n_col)) > tol ) then
892 0 : err_code = clubb_fatal_error
893 0 : write(fstderr,*) "in subroutine assert_corr_symmetric: ", &
894 0 : "Correlation array is non symmetric."
895 0 : write(fstderr,*) corr_array_n
896 0 : return
897 : end if
898 :
899 0 : if ( n_col == n_row .and. abs(corr_array_n(n_col, n_row)-one) > eps ) then
900 0 : err_code = clubb_fatal_error
901 0 : write(fstderr,*) "in subroutine assert_corr_symmetric: ", &
902 0 : "Correlation array is formatted incorrectly."
903 0 : write(fstderr,*) corr_array_n
904 0 : return
905 : end if
906 : end do
907 : end do
908 :
909 : end subroutine assert_corr_symmetric
910 :
911 : !-----------------------------------------------------------------------------
912 0 : subroutine print_corr_matrix( pdf_dim, & ! intent(in)
913 0 : corr_array_n ) ! intent(in)
914 :
915 : ! Description:
916 : ! Prints the correlation matrix to the console.
917 : ! References:
918 : ! None
919 : !---------------------------------------------------------------------------
920 :
921 : use clubb_precision, only: core_rknd
922 :
923 : implicit none
924 :
925 : ! Input Variables
926 : integer, intent(in) :: &
927 : pdf_dim ! Number of variables in the correlation array
928 :
929 : real( kind = core_rknd ), dimension(pdf_dim, pdf_dim), &
930 : intent(in) :: corr_array_n ! Normal space correlation array to be printed
931 :
932 : ! Local Variables
933 : integer :: n, & ! Loop indeces
934 : m, &
935 : current_character_index ! keeps track of the position in the string
936 :
937 : character(LEN=72) :: current_line ! The current line to be printed
938 : character(LEN=10) :: str_array_value
939 :
940 : !----- Begin Code -----
941 :
942 0 : current_character_index = 0
943 :
944 0 : do n = 1, pdf_dim
945 0 : do m = 1, pdf_dim
946 0 : write(str_array_value,'(F5.2)') corr_array_n(m,n)
947 0 : current_line = current_line(1:current_character_index)//str_array_value
948 0 : current_character_index = current_character_index + 6
949 : end do
950 0 : write(*, *) current_line
951 0 : current_line = ""
952 0 : current_character_index = 0
953 : end do
954 :
955 0 : end subroutine print_corr_matrix
956 : !-----------------------------------------------------------------------------
957 :
958 0 : end module corr_varnce_module
|