Line data Source code
1 :
2 : module constituents
3 :
4 : ! Metadata manager for the advected constituents.
5 :
6 : use shr_kind_mod, only: r8 => shr_kind_r8
7 : use shr_const_mod, only: shr_const_rgas
8 : use spmd_utils, only: masterproc
9 : use cam_abortutils, only: endrun
10 : use cam_logfile, only: iulog
11 :
12 : implicit none
13 : private
14 : save
15 :
16 : ! Public interfaces
17 : public :: &
18 : cnst_readnl, &! read namelist
19 : cnst_add, &! add a constituent to the list of advected constituents
20 : cnst_num_avail, &! returns the number of available slots in the constituent array
21 : cnst_get_ind, &! get the index of a constituent
22 : cnst_get_type_byind, &! return mixing ratio type of a constituent
23 : cnst_get_molec_byind,&! return molecular diffusion type of a constituent
24 : cnst_read_iv, &! query whether constituent initial values are read from initial file
25 : cnst_chk_dim, &! check that number of constituents added equals dimensions (pcnst)
26 : cnst_cam_outfld, &! Returns true if default CAM output was specified in the cnst_add calls.
27 : cnst_set_spec_class, &! Sets the type of species class
28 : cnst_is_a_water_species,&! Returns true for constituents identified as water species
29 : cnst_set_convtran2 ! Override for convtran2 values set by the cnst_add routine
30 :
31 : ! Public data
32 :
33 : integer, parameter, public :: pcnst = PCNST ! number of advected constituents (including water vapor)
34 :
35 : character(len=16), public, protected :: cnst_name(pcnst) ! constituent names
36 : character(len=128),public, protected :: cnst_longname(pcnst) ! long name of constituents
37 :
38 : ! Namelist variables
39 : logical, public, protected :: readtrace = .true. ! true => obtain initial tracer data from IC file
40 :
41 : integer, public, parameter :: cnst_spec_class_undefined = 0
42 : integer, public, parameter :: cnst_spec_class_cldphysics = 1
43 : integer, public, parameter :: cnst_spec_class_aerosol = 2
44 : integer, public, parameter :: cnst_spec_class_gas = 3
45 : integer, public, parameter :: cnst_spec_class_other = 4
46 :
47 : !
48 : ! Constants for each tracer
49 :
50 : integer, public, protected :: cnst_species_class(pcnst) = cnst_spec_class_undefined ! indicates species class &
51 : ! (cldphysics, aerosol, gas )
52 : real(r8), public :: cnst_cp (pcnst) ! specific heat at constant pressure (J/kg/K)
53 : real(r8), public :: cnst_cv (pcnst) ! specific heat at constant volume (J/kg/K)
54 : real(r8), public :: cnst_mw (pcnst) ! molecular weight (kg/kmole)
55 : character*3, public, protected :: cnst_type(pcnst)! wet or dry mixing ratio
56 : character*5, public :: cnst_molec(pcnst) ! major or minor species molecular diffusion
57 : logical, public, protected :: cnst_ndropmixed(pcnst) = .false. ! vertically mixed by ndrop activation process
58 : real(r8), public :: cnst_rgas(pcnst) ! gas constant ()
59 : real(r8), public :: qmin (pcnst) ! minimum permitted constituent concentration (kg/kg)
60 : real(r8), public :: qmincg (pcnst) ! for backward compatibility only
61 : logical, public, protected :: cnst_fixed_ubc(pcnst) = .false. ! upper boundary condition (concentration)
62 : logical, public, protected :: cnst_fixed_ubflx(pcnst) = .false. ! upper boundary non-zero fixed constituent flux
63 : logical, public, protected :: cnst_is_convtran1(pcnst) = .false. ! do convective transport in phase 1
64 : logical, public, protected :: cnst_is_convtran2(pcnst) = .false. ! do convective transport in phase 2
65 :
66 : !++bee - temporary... These names should be declared in the module that makes the addfld and outfld calls.
67 : ! Lists of tracer names and diagnostics
68 : character(len=16), public :: apcnst (pcnst) ! constituents after physics (FV core only)
69 : character(len=16), public :: bpcnst (pcnst) ! constituents before physics (FV core only)
70 : character(len=16), public :: hadvnam (pcnst) ! names of horizontal advection tendencies
71 : character(len=16), public :: vadvnam (pcnst) ! names of vertical advection tendencies
72 : character(len=16), public :: dcconnam (pcnst) ! names of convection tendencies
73 : character(len=16), public :: fixcnam (pcnst) ! names of species slt fixer tendencies
74 : character(len=16), public :: tendnam (pcnst) ! names of total tendencies of species
75 : character(len=16), public :: ptendnam (pcnst) ! names of total physics tendencies of species
76 : character(len=16), public :: sflxnam (pcnst) ! names of surface fluxes of species
77 : character(len=16), public :: tottnam (pcnst) ! names for horz + vert + fixer tendencies
78 :
79 : ! Private data
80 :
81 : integer :: padv = 0 ! index pointer to last advected tracer
82 : logical :: read_init_vals(pcnst) ! true => read initial values from initial file
83 : logical :: cam_outfld_(pcnst) ! true => default CAM output of constituents in kg/kg
84 : ! false => chemistry is responsible for making outfld
85 : ! calls for constituents
86 :
87 : !==============================================================================================
88 : CONTAINS
89 : !==============================================================================================
90 :
91 2304 : subroutine cnst_readnl(nlfile)
92 :
93 : use namelist_utils, only: find_group_name
94 : use units, only: getunit, freeunit
95 : use spmd_utils, only: mpicom, mstrid=>masterprocid, mpi_logical
96 :
97 :
98 : character(len=*), intent(in) :: nlfile ! filepath for file containing namelist input
99 :
100 : ! Local variables
101 : integer :: unitn, ierr
102 : character(len=*), parameter :: sub = 'cnst_readnl'
103 :
104 : namelist /constituents_nl/ readtrace
105 : !-----------------------------------------------------------------------------
106 :
107 2304 : if (masterproc) then
108 3 : unitn = getunit()
109 3 : open( unitn, file=trim(nlfile), status='old' )
110 3 : call find_group_name(unitn, 'constituents_nl', status=ierr)
111 3 : if (ierr == 0) then
112 0 : read(unitn, constituents_nl, iostat=ierr)
113 0 : if (ierr /= 0) then
114 0 : call endrun(sub//': FATAL: reading namelist')
115 : end if
116 : end if
117 3 : close(unitn)
118 3 : call freeunit(unitn)
119 : end if
120 :
121 2304 : call mpi_bcast(readtrace, 1, mpi_logical, mstrid, mpicom, ierr)
122 2304 : if (ierr /= 0) call endrun(sub//": FATAL: mpi_bcast: readtrace")
123 :
124 2304 : if (masterproc) then
125 3 : write(iulog,*)'Summary of constituent module options:'
126 3 : write(iulog,*)' Read constituent initial values from initial file by default: ', readtrace
127 : end if
128 :
129 2304 : end subroutine cnst_readnl
130 :
131 : !=========================================================================================
132 :
133 :
134 94464 : subroutine cnst_add (name, mwc, cpc, qminc, &
135 : ind, longname, readiv, mixtype, molectype, ndropmixed, cam_outfld, &
136 : fixed_ubc, fixed_ubflx, is_convtran1, is_convtran2, cnst_spec_class)
137 :
138 : ! Register a constituent.
139 :
140 : character(len=*), intent(in) :: &
141 : name ! constituent name used as variable name in history file output (8 char max)
142 : real(r8),intent(in) :: mwc ! constituent molecular weight (kg/kmol)
143 : real(r8),intent(in) :: cpc ! constituent specific heat at constant pressure (J/kg/K)
144 : real(r8),intent(in) :: qminc ! minimum value of mass mixing ratio (kg/kg)
145 : ! normally 0., except water 1.E-12, for radiation.
146 : integer, intent(out) :: ind ! global constituent index (in q array)
147 :
148 : character(len=*), intent(in), optional :: &
149 : longname ! value for long_name attribute in netcdf output (128 char max, defaults to name)
150 : logical, intent(in), optional :: &
151 : readiv ! true => read initial values from initial file (default: true)
152 : character(len=*), intent(in), optional :: &
153 : mixtype ! mixing ratio type (dry, wet)
154 : character(len=*), intent(in), optional :: &
155 : molectype ! molecular diffusion type (minor, major)
156 : logical, intent(in), optional :: &
157 : ndropmixed ! vertically mixed by ndrop activation process
158 : logical, intent(in), optional :: &
159 : cam_outfld ! true => default CAM output of constituent in kg/kg
160 : logical, intent(in), optional :: &
161 : fixed_ubc ! true => const has a fixed upper bndy condition
162 : logical, intent(in), optional :: &
163 : fixed_ubflx ! true => const has a non-zero fixed upper bndy flux value
164 : logical, intent(in), optional :: &
165 : is_convtran1 ! true => convective transport in convtran1
166 : logical, intent(in), optional :: &
167 : is_convtran2 ! true => convective transport in convtran2
168 : integer, intent(in), optional :: &
169 : cnst_spec_class ! type of species class
170 :
171 : character(len=*), parameter :: sub='cnst_add'
172 : character(len=128) :: errmsg
173 : !-----------------------------------------------------------------------
174 :
175 : ! set tracer index and check validity
176 94464 : padv = padv+1
177 94464 : ind = padv
178 94464 : if (padv > pcnst) then
179 0 : write(errmsg, *) sub//': FATAL: advected tracer (', trim(name), ') index is greater than number of constituents'
180 0 : call endrun(errmsg)
181 : end if
182 :
183 : ! set tracer name and constants
184 94464 : cnst_name(ind) = name
185 94464 : if (present(longname)) then
186 94464 : cnst_longname(ind) = longname
187 : else
188 0 : cnst_longname(ind) = name
189 : end if
190 :
191 : ! set whether to read initial values from initial file
192 94464 : if (present(readiv)) then
193 71424 : read_init_vals(ind) = readiv
194 : else
195 23040 : read_init_vals(ind) = readtrace
196 : end if
197 :
198 : ! set constituent mixing ratio type
199 94464 : if (present(mixtype)) then
200 69120 : cnst_type(ind) = mixtype
201 : else
202 25344 : cnst_type(ind) = 'wet'
203 : end if
204 :
205 : ! set constituent molecular diffusion type
206 94464 : if (present(molectype)) then
207 69120 : cnst_molec(ind) = molectype
208 : else
209 25344 : cnst_molec(ind) = 'minor'
210 : end if
211 :
212 : ! vertically mixed by ndrop activation process
213 94464 : if (present(ndropmixed)) then
214 71424 : cnst_ndropmixed(ind) = ndropmixed
215 : else
216 23040 : cnst_ndropmixed(ind) = .false.
217 : end if
218 :
219 : ! set outfld type
220 : ! (false: the module declaring the constituent is responsible for outfld calls)
221 94464 : if (present(cam_outfld)) then
222 69120 : cam_outfld_(ind) = cam_outfld
223 : else
224 25344 : cam_outfld_(ind) = .true.
225 : end if
226 :
227 : ! set upper boundary condition type
228 94464 : if (present(fixed_ubc)) then
229 69120 : cnst_fixed_ubc(ind) = fixed_ubc
230 : else
231 25344 : cnst_fixed_ubc(ind) = .false.
232 : end if
233 :
234 : ! set upper boundary flux type
235 94464 : if (present(fixed_ubflx)) then
236 69120 : cnst_fixed_ubflx(ind) = fixed_ubflx
237 : else
238 25344 : cnst_fixed_ubflx(ind) = .false.
239 : end if
240 :
241 : ! Set flag for convective transport by first call to convtran (phase 1).
242 94464 : if (present(is_convtran1)) then
243 25344 : cnst_is_convtran1(ind) = is_convtran1
244 : else
245 69120 : cnst_is_convtran1(ind) = .false.
246 : end if
247 : ! Set flag for convective transport after wetdep (phase 2).
248 94464 : if (present(is_convtran2)) then
249 0 : cnst_is_convtran2(ind) = is_convtran2
250 : else
251 : ! The default is that all constituents except water vapor are transported in phase 2
252 : ! unless they were transported in phase 1 (typically the condensed water species)
253 94464 : if (ind > 1) cnst_is_convtran2(ind) = .not. cnst_is_convtran1(ind)
254 : end if
255 : ! consistency check -- It is OK to completely turn off the deep scheme transport by setting
256 : ! both cnst_is_convtran1 and cnst_is_convtran2 to FALSE. But it is an error to
257 : ! have both set TRUE.
258 94464 : if (cnst_is_convtran1(ind) .and. cnst_is_convtran2(ind)) then
259 0 : call endrun(sub//': FATAL: cannot set both cnst_is_convtran1 and cnst_is_convtran2 to TRUE')
260 : end if
261 :
262 : ! Set type for species class
263 94464 : if ( present(cnst_spec_class) ) then
264 0 : cnst_species_class(ind) = cnst_spec_class
265 : else
266 94464 : cnst_species_class(ind) = cnst_spec_class_undefined
267 : end if
268 :
269 94464 : cnst_cp (ind) = cpc
270 94464 : cnst_mw (ind) = mwc
271 94464 : qmin (ind) = qminc
272 94464 : if (ind == 1) then
273 : ! qmincg for water vapor set to zero
274 2304 : qmincg(ind) = 0._r8
275 : else
276 92160 : qmincg(ind) = qminc
277 : end if
278 :
279 94464 : cnst_rgas(ind) = shr_const_rgas * mwc
280 94464 : cnst_cv (ind) = cpc - cnst_rgas(ind)
281 :
282 94464 : end subroutine cnst_add
283 :
284 : !----------------------------------------------------------------------------------------------
285 :
286 43776 : subroutine cnst_set_convtran2(ind, is_convtran2)
287 :
288 : ! Allow user to override the value of cnst_is_convtran2 set by a previous cnst_add call.
289 :
290 : integer, intent(in) :: ind ! global constituent index (in q array)
291 : logical, intent(in) :: is_convtran2 ! true => convect in convtran2
292 :
293 : character(len=*), parameter :: sub = 'cnst_set_convtran2'
294 : character(len=128) :: errmsg
295 : !-----------------------------------------------------------------------
296 :
297 : ! check index
298 43776 : if (ind <= 0 .or. ind > padv) then
299 0 : write(errmsg,*) sub//': FATAL: bad tracer index: padv, ind = ', padv, ind
300 0 : call endrun(errmsg)
301 : end if
302 :
303 : ! Set flag for convective transport after wetdep (phase 2).
304 43776 : cnst_is_convtran2(ind) = is_convtran2
305 :
306 : ! consistency check -- It is OK to completely turn off the tracer convection by setting
307 : ! both cnst_is_convtran1 and cnst_is_convtran2 to FALSE. But it is an error to
308 : ! have both set TRUE.
309 43776 : if (cnst_is_convtran1(ind) .and. cnst_is_convtran2(ind)) then
310 0 : call endrun(sub//': FATAL: cannot set both cnst_is_convtran1 and cnst_is_convtran2 to TRUE')
311 : end if
312 :
313 43776 : end subroutine cnst_set_convtran2
314 :
315 : !----------------------------------------------------------------------------------------------
316 :
317 69120 : subroutine cnst_set_spec_class(ind, cnst_spec_class_in)
318 :
319 : ! Allow user to override the value of cnst_spec_class set by a previous cnst_add call.
320 :
321 : integer, intent(in) :: ind ! global constituent index (in q array)
322 : integer, intent(in) :: cnst_spec_class_in ! species class designator
323 :
324 : character(len=*), parameter :: subname = 'cnst_set_spec_class'
325 : !-----------------------------------------------------------------------
326 :
327 : ! check index
328 69120 : if (ind <= 0 .or. ind > padv) then
329 0 : write(iulog,*) subname//': illegal tracer index: padv, ind = ', padv, ind
330 0 : call endrun(subname//': illegal tracer index')
331 : end if
332 :
333 : ! Check designator
334 : if (cnst_spec_class_in /= cnst_spec_class_undefined .and. &
335 : cnst_spec_class_in /= cnst_spec_class_cldphysics .and. &
336 : cnst_spec_class_in /= cnst_spec_class_aerosol .and. &
337 69120 : cnst_spec_class_in /= cnst_spec_class_gas .and. &
338 : cnst_spec_class_in /= cnst_spec_class_other ) then
339 0 : write(iulog,*) subname//': trying to use invalid cnst_spec_class designator', cnst_spec_class_in
340 0 : call endrun(subname//': invalid cnst_spec_class designator')
341 : end if
342 :
343 : ! Set flag for convective transport after wetdep (phase 2).
344 69120 : cnst_species_class(ind) = cnst_spec_class_in
345 :
346 69120 : end subroutine cnst_set_spec_class
347 :
348 : !==============================================================================
349 :
350 0 : function cnst_num_avail()
351 :
352 : ! return number of available slots in the constituent array
353 :
354 : integer cnst_num_avail
355 :
356 0 : cnst_num_avail = pcnst - padv
357 :
358 0 : end function cnst_num_avail
359 :
360 : !==============================================================================
361 :
362 91871736 : subroutine cnst_get_ind (name, ind, abort)
363 :
364 : ! Get the index of a constituent. Optional abort argument allows returning
365 : ! control to caller when constituent name is not found. Default behavior is
366 : ! to call endrun when name is not found.
367 :
368 : !-----------------------------Arguments---------------------------------
369 : character(len=*), intent(in) :: name ! constituent name
370 : integer, intent(out) :: ind ! global constituent index (in q array)
371 : logical, optional, intent(in) :: abort ! optional flag controlling abort
372 :
373 : !---------------------------Local workspace-----------------------------
374 : integer :: m ! tracer index
375 : logical :: abort_on_error
376 : character(len=*), parameter :: sub='cnst_get_ind'
377 : !-----------------------------------------------------------------------
378 :
379 : ! Find tracer name in list
380 1976943408 : do m = 1, pcnst
381 1976943408 : if (name == cnst_name(m)) then
382 84218760 : ind = m
383 84218760 : return
384 : end if
385 : end do
386 :
387 : ! Unrecognized name
388 7652976 : abort_on_error = .true.
389 7652976 : if (present(abort)) abort_on_error = abort
390 :
391 7652976 : if (abort_on_error) then
392 0 : write(iulog, *) sub//': FATAL: name:', name, ' not found in constituent list: ', cnst_name(:)
393 0 : call endrun(sub//': FATAL: name not found')
394 : end if
395 :
396 : ! error return
397 7652976 : ind = -1
398 :
399 : end subroutine cnst_get_ind
400 :
401 : !==============================================================================================
402 :
403 346368 : character*3 function cnst_get_type_byind(ind)
404 :
405 : ! Return the mixing ratio type of a constituent
406 :
407 : !-----------------------------Arguments---------------------------------
408 : integer, intent(in) :: ind ! global constituent index (in q array)
409 :
410 : !---------------------------Local workspace-----------------------------
411 : character(len=*), parameter :: sub='cnst_get_type_byind'
412 : character(len=128) :: errmsg
413 : !-----------------------------------------------------------------------
414 :
415 346368 : if (ind > 0 .and. ind <= pcnst) then
416 346368 : cnst_get_type_byind = cnst_type(ind)
417 : else
418 : ! index out of range
419 0 : write(errmsg,*) sub//': FATAL: bad value for constituent index=', ind
420 0 : call endrun(errmsg)
421 : end if
422 :
423 346368 : end function cnst_get_type_byind
424 :
425 : !==============================================================================================
426 :
427 48384 : character*5 function cnst_get_molec_byind (ind)
428 :
429 : ! Return the molecular diffusion type of a constituent
430 :
431 : !-----------------------------Arguments---------------------------------
432 : integer, intent(in) :: ind ! global constituent index (in q array)
433 :
434 : !---------------------------Local workspace-----------------------------
435 : character(len=*), parameter :: sub='cnst_get_molec_byind'
436 : character(len=128) :: errmsg
437 : !-----------------------------------------------------------------------
438 :
439 48384 : if (ind > 0 .and. ind <= pcnst) then
440 48384 : cnst_get_molec_byind = cnst_molec(ind)
441 : else
442 : ! index out of range
443 0 : write(errmsg,*) sub//': FATAL: bad value for constituent index=', ind
444 0 : call endrun(errmsg)
445 : end if
446 :
447 48384 : end function cnst_get_molec_byind
448 :
449 : !==============================================================================
450 :
451 69120 : function cnst_read_iv(m)
452 :
453 : ! Query whether to attempt to read constituent initial values from initial file.
454 :
455 : !-----------------------------Arguments---------------------------------
456 : integer, intent(in) :: m ! constituent index
457 :
458 : logical :: cnst_read_iv ! true => read initial values from inital file
459 :
460 : character(len=*), parameter :: sub='cnst_read_iv'
461 : character(len=128) :: errmsg
462 : !-----------------------------------------------------------------------
463 :
464 69120 : if (m > 0 .and. m <= pcnst) then
465 69120 : cnst_read_iv = read_init_vals(m)
466 : else
467 : ! index out of range
468 0 : write(errmsg,*) sub//': FATAL: bad value for constiuent index=', m
469 0 : call endrun(errmsg)
470 : end if
471 :
472 69120 : end function cnst_read_iv
473 :
474 : !==============================================================================
475 :
476 2304 : subroutine cnst_chk_dim
477 :
478 : ! Check that the number of registered constituents is pcnst
479 : ! Write constituent list to log file.
480 :
481 : integer :: i, m
482 : character(len=*), parameter :: sub='cnst_chk_dim'
483 : character(len=128) :: errmsg
484 : !-----------------------------------------------------------------------
485 :
486 2304 : if (padv /= pcnst) then
487 0 : write(errmsg, *) sub//': FATAL: number of advected tracer (',padv, &
488 0 : ') not equal to pcnst (', pcnst, ')'
489 0 : call endrun (errmsg)
490 : endif
491 :
492 2304 : if (masterproc) then
493 3 : write(iulog,*) 'Advected constituent list:'
494 126 : do i = 1, pcnst
495 246 : write(iulog,'(2x,i4,2x,a8,2x,a128,2x,a3)') i, cnst_name(i), cnst_longname(i), &
496 372 : cnst_type(i)
497 : end do
498 : end if
499 :
500 : ! Set names of advected tracer diagnostics
501 96768 : do m=1,pcnst
502 94464 : apcnst (m) = trim(cnst_name(m))//'AP'
503 94464 : bpcnst (m) = trim(cnst_name(m))//'BP'
504 94464 : hadvnam (m) = 'HA'//cnst_name(m)
505 94464 : vadvnam (m) = 'VA'//cnst_name(m)
506 94464 : fixcnam (m) = 'DF'//cnst_name(m)
507 94464 : tendnam (m) = 'TE'//cnst_name(m)
508 94464 : ptendnam (m) = 'PTE'//cnst_name(m)
509 94464 : tottnam (m) = 'TA'//cnst_name(m)
510 96768 : sflxnam(m) = 'SF'//cnst_name(m)
511 : end do
512 :
513 2304 : end subroutine cnst_chk_dim
514 :
515 : !==============================================================================
516 :
517 4606848 : function cnst_cam_outfld(m)
518 :
519 : ! Query whether default CAM outfld calls should be made.
520 :
521 : !-----------------------------------------------------------------------
522 : integer, intent(in) :: m ! constituent index
523 :
524 : logical :: cnst_cam_outfld ! true => use default CAM outfld calls
525 :
526 : character(len=*), parameter :: sub='cnst_cam_outfld'
527 : character(len=128) :: errmsg
528 : !-----------------------------------------------------------------------
529 :
530 4606848 : if (m > 0 .and. m <= pcnst) then
531 4606848 : cnst_cam_outfld = cam_outfld_(m)
532 : else
533 : ! index out of range
534 0 : write(errmsg,*) sub//': FATAL: bad value for constiuent index=', m
535 0 : call endrun(errmsg)
536 : end if
537 :
538 4606848 : end function cnst_cam_outfld
539 :
540 : !==============================================================================
541 :
542 69120 : pure logical function cnst_is_a_water_species(name)
543 :
544 : ! test whether the input name matches the name of a water species
545 :
546 : character(len=*), intent(in) :: name
547 : !-------------------------------------------------------------------------
548 :
549 69120 : cnst_is_a_water_species = .false.
550 :
551 : if (name == 'Q' .or. &
552 : name == 'CLDLIQ' .or. &
553 : name == 'CLDICE' .or. &
554 : name == 'RAINQM' .or. &
555 69120 : name == 'SNOWQM' .or. &
556 13824 : name == 'GRAUQM' ) cnst_is_a_water_species = .true.
557 :
558 69120 : end function cnst_is_a_water_species
559 :
560 : !==============================================================================
561 :
562 : end module constituents
|