5 !----------------------------------------------------------------------------------------------
7 ! Purpose: Contains data and functions for manipulating advected and non-advected constituents.
10 ! B.A. Boville Original version
11 ! June 2003 P. Rasch Add wet/dry m.r. specifier
12 ! 2004-08-28 B. Eaton Add query function to allow turning off the default CAM output of
13 ! constituents so that chemistry module can make the outfld calls.
14 ! Allow cnst_get_ind to return without aborting when constituent not
16 ! 2006-10-31 B. Eaton Remove 'non-advected' constituent functionality.
17 !----------------------------------------------------------------------------------------------
18 use shr_kind_mod, only: r8 => shr_kind_r8
19 use physconst, only: r_universal
22 use spmd_utils, only: masterproc
23 use abortutils, only: endrun
24 use cam_logfile, only: iulog
26 use module_cam_support, only: masterproc,endrun,iulog,pcnst
35 public cnst_add ! add a constituent to the list of advected constituents
36 public cnst_num_avail ! returns the number of available slots in the constituent array
37 public cnst_get_ind ! get the index of a constituent
38 public cnst_get_type_byind ! get the type of a constituent
39 public cnst_get_type_byname ! get the type of a constituent
40 public cnst_read_iv ! query whether constituent initial values are read from initial file
41 public cnst_chk_dim ! check that number of constituents added equals dimensions (pcnst)
42 public cnst_cam_outfld ! Returns true if default CAM output was specified in the cnst_add calls.
46 integer, parameter, public :: pcnst = PCNST ! number of advected constituents (including water vapor)
49 character(len=16), public :: cnst_name(pcnst) ! constituent names
50 character(len=128),public :: cnst_longname(pcnst) ! long name of constituents
53 logical, public :: readtrace = .true. ! true => obtain initial tracer data from IC file
56 ! Constants for each tracer
57 real(r8), public :: cnst_cp (pcnst) ! specific heat at constant pressure (J/kg/K)
58 real(r8), public :: cnst_cv (pcnst) ! specific heat at constant volume (J/kg/K)
59 real(r8), public :: cnst_mw (pcnst) ! molecular weight (kg/kmole)
60 character*3, public :: cnst_type(pcnst) ! wet or dry mixing ratio
61 real(r8), public :: cnst_rgas(pcnst) ! gas constant ()
62 real(r8), public :: qmin (pcnst) ! minimum permitted constituent concentration (kg/kg)
63 real(r8), public :: qmincg (pcnst) ! for backward compatibility only
64 logical, public :: cnst_fixed_ubc(pcnst) = .false. ! upper bndy condition = fixed ?
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 :: dmetendnam(pcnst) ! names of dme adjusted tracers (FV)
77 character(len=16), public :: sflxnam (pcnst) ! names of surface fluxes of species
78 character(len=16), public :: tottnam (pcnst) ! names for horz + vert + fixer tendencies
82 integer :: padv = 0 ! index pointer to last advected tracer
83 logical :: read_init_vals(pcnst) ! true => read initial values from initial file
84 logical :: cam_outfld_(pcnst) ! true => default CAM output of constituents in kg/kg
85 ! false => chemistry is responsible for making outfld
86 ! calls for constituents
88 !==============================================================================================
90 !==============================================================================================
92 subroutine cnst_add (name, mwc, cpc, qminc, &
93 ind, longname, readiv, mixtype, cam_outfld, fixed_ubc)
94 !-----------------------------------------------------------------------
96 ! Purpose: Register a constituent to be advected by the large scale winds and transported by
97 ! subgrid scale processes.
99 !---------------------------------------------------------------------------------
101 character(len=*), intent(in) :: &
102 name ! constituent name used as variable name in history file output (8 char max)
103 real(r8),intent(in) :: mwc ! constituent molecular weight (kg/kmol)
104 real(r8),intent(in) :: cpc ! constituent specific heat at constant pressure (J/kg/K)
105 real(r8),intent(in) :: qminc ! minimum value of mass mixing ratio (kg/kg)
106 ! normally 0., except water 1.E-12, for radiation.
107 integer, intent(out) :: ind ! global constituent index (in q array)
109 character(len=*), intent(in), optional :: &
110 longname ! value for long_name attribute in netcdf output (128 char max, defaults to name)
111 logical, intent(in), optional :: &
112 readiv ! true => read initial values from initial file (default: true)
113 character(len=*), intent(in), optional :: &
114 mixtype ! mixing ratio type (dry, wet)
115 logical, intent(in), optional :: &
116 cam_outfld ! true => default CAM output of constituent in kg/kg
117 logical, intent(in), optional :: &
118 fixed_ubc ! true => const has a fixed upper bndy condition
120 !-----------------------------------------------------------------------
122 ! set tracer index and check validity, advected tracer
125 if (padv > pcnst) then
126 write(iulog,*) 'CNST_ADD: advected tracer index greater than pcnst = ', pcnst
130 ! set tracer name and constants
131 cnst_name(ind) = name
132 if ( present(longname) )then
133 cnst_longname(ind) = longname
135 cnst_longname(ind) = name
138 ! set whether to read initial values from initial file
139 if ( present(readiv) ) then
140 read_init_vals(ind) = readiv
142 read_init_vals(ind) = readtrace
145 ! set constituent mixing ratio type
146 if ( present(mixtype) )then
147 cnst_type(ind) = mixtype
149 cnst_type(ind) = 'wet'
153 ! (false: the module declaring the constituent is responsible for outfld calls)
154 if ( present(cam_outfld) ) then
155 cam_outfld_(ind) = cam_outfld
157 cam_outfld_(ind) = .true.
160 ! set upper boundary condition type
161 if ( present(fixed_ubc) ) then
162 cnst_fixed_ubc(ind) = fixed_ubc
164 cnst_fixed_ubc(ind) = .false.
171 if (ind == 1) qmincg = 0._r8 ! This crap is replicate what was there before ****
173 cnst_rgas(ind) = r_universal * mwc
174 cnst_cv (ind) = cpc - cnst_rgas(ind)
177 end subroutine cnst_add
179 !==============================================================================
181 function cnst_num_avail()
183 ! return number of available slots in the constituent array
185 integer cnst_num_avail
187 cnst_num_avail = pcnst - padv
189 end function cnst_num_avail
191 !==============================================================================
193 subroutine cnst_get_ind (name, ind, abort)
194 !-----------------------------------------------------------------------
196 ! Purpose: Get the index of a constituent
198 ! Author: B.A. Boville
200 !-----------------------------Arguments---------------------------------
202 character(len=*), intent(in) :: name ! constituent name
203 integer, intent(out) :: ind ! global constituent index (in q array)
204 logical, optional, intent(in) :: abort ! optional flag controlling abort
206 !---------------------------Local workspace-----------------------------
207 integer :: m ! tracer index
208 logical :: abort_on_error
209 !-----------------------------------------------------------------------
211 ! Find tracer name in list
213 if (name == cnst_name(m)) then
220 abort_on_error = .true.
221 if ( present(abort) ) abort_on_error = abort
223 if ( abort_on_error ) then
224 write(iulog,*) 'CNST_GET_IND, name:', name, ' not found in list:', cnst_name(:)
226 call wrf_message(iulog)
228 call endrun('CNST_GET_IND: name not found')
234 end subroutine cnst_get_ind
236 !==============================================================================================
238 character*3 function cnst_get_type_byind (ind)
239 !-----------------------------------------------------------------------
241 ! Purpose: Get the type of a constituent
244 ! <Describe the algorithm(s) used in the routine.>
245 ! <Also include any applicable external references.>
247 ! Author: P. J. Rasch
249 !-----------------------------Arguments---------------------------------
251 integer, intent(in) :: ind ! global constituent index (in q array)
253 !---------------------------Local workspace-----------------------------
254 integer :: m ! tracer index
256 !-----------------------------------------------------------------------
258 if (ind.le.pcnst) then
259 cnst_get_type_byind = cnst_type(ind)
262 write(iulog,*) 'CNST_GET_TYPE_BYIND, ind:', ind
267 end function cnst_get_type_byind
269 !==============================================================================================
271 character*3 function cnst_get_type_byname (name)
272 !-----------------------------------------------------------------------
274 ! Purpose: Get the type of a constituent
277 ! <Describe the algorithm(s) used in the routine.>
278 ! <Also include any applicable external references.>
280 ! Author: P. J. Rasch
282 !-----------------------------Arguments---------------------------------
284 character(len=*), intent(in) :: name ! constituent name
286 !---------------------------Local workspace-----------------------------
287 integer :: m ! tracer index
289 !-----------------------------------------------------------------------
292 if (name == cnst_name(m)) then
293 cnst_get_type_byname = cnst_type(m)
299 write(iulog,*) 'CNST_GET_TYPE_BYNAME, name:', name, ' not found in list:', cnst_name(:)
302 end function cnst_get_type_byname
304 !==============================================================================
305 function cnst_read_iv(m)
306 !-----------------------------------------------------------------------
308 ! Purpose: Query whether constituent initial values are read from initial file.
312 !-----------------------------Arguments---------------------------------
314 integer, intent(in) :: m ! constituent index
316 logical :: cnst_read_iv ! true => read initial values from inital file
317 !-----------------------------------------------------------------------
319 cnst_read_iv = read_init_vals(m)
320 end function cnst_read_iv
322 !==============================================================================
323 subroutine cnst_chk_dim
324 !-----------------------------------------------------------------------
326 ! Purpose: Check that the number of registered constituents of each type is the
327 ! same as the dimension
330 ! <Describe the algorithm(s) used in the routine.>
331 ! <Also include any applicable external references.>
333 ! Author: B.A. Boville
336 !-----------------------------------------------------------------------
338 if (padv /= pcnst) then
339 write(iulog,*)'CNST_CHK_DIM: number of advected tracer ',padv, ' not equal to pcnst = ',pcnst
344 write(iulog,*) 'Advected constituent list:'
346 write(iulog,'(i4,2x,a8,2x,a128,2x,a3)') i, cnst_name(i), cnst_longname(i), cnst_type(i)
350 ! Set names of advected tracer diagnostics
352 apcnst (m) = trim(cnst_name(m))//'AP'
353 bpcnst (m) = trim(cnst_name(m))//'BP'
354 hadvnam (m) = 'HA'//cnst_name(m)
355 vadvnam (m) = 'VA'//cnst_name(m)
356 fixcnam (m) = 'DF'//cnst_name(m)
357 tendnam (m) = 'TE'//cnst_name(m)
358 ptendnam (m) = 'PTE'//cnst_name(m)
359 dmetendnam(m) = 'DME'//cnst_name(m)
360 tottnam (m) = 'TA'//cnst_name(m)
361 sflxnam(m) = 'SF'//cnst_name(m)
365 end subroutine cnst_chk_dim
367 !==============================================================================
369 function cnst_cam_outfld(m)
370 !-----------------------------------------------------------------------
373 ! Query whether default CAM outfld calls should be made.
375 !-----------------------------------------------------------------------
376 integer, intent(in) :: m ! constituent index
377 logical :: cnst_cam_outfld ! true => use default CAM outfld calls
378 !-----------------------------------------------------------------------
380 cnst_cam_outfld = cam_outfld_(m)
382 end function cnst_cam_outfld
384 !==============================================================================
386 end module constituents