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