MOM6
MOM_tracer_registry.F90
1 !> This module contains the tracer_registry_type and the subroutines
2 !! that handle registration of tracers and related subroutines.
3 !! The primary subroutine, register_tracer, is called to indicate the
4 !! tracers advected and diffused.
6 
7 ! This file is part of MOM6. See LICENSE.md for the license.
8 
9 ! use MOM_diag_mediator, only : diag_ctrl
10 use mom_coms, only : reproducing_sum
11 use mom_debugging, only : hchksum
12 use mom_diag_mediator, only : diag_ctrl, register_diag_field, post_data, safe_alloc_ptr
14 use mom_diag_mediator, only : diag_copy_storage_to_diag, diag_save_grids, diag_restore_grids
15 use mom_error_handler, only : mom_error, fatal, warning, mom_mesg, is_root_pe
17 use mom_hor_index, only : hor_index_type
18 use mom_grid, only : ocean_grid_type
19 use mom_io, only : vardesc, query_vardesc, cmor_long_std
21 use mom_string_functions, only : lowercase
22 use mom_time_manager, only : time_type
24 
25 implicit none ; private
26 
27 #include <MOM_memory.h>
28 
29 public register_tracer
30 public mom_tracer_chksum, mom_tracer_chkinv
31 public register_tracer_diagnostics, post_tracer_diagnostics, post_tracer_transport_diagnostics
32 public preale_tracer_diagnostics, postale_tracer_diagnostics
33 public tracer_registry_init, lock_tracer_registry, tracer_registry_end
34 public tracer_name_lookup
35 
36 !> The tracer type
37 type, public :: tracer_type
38 
39  real, dimension(:,:,:), pointer :: t => null() !< tracer concentration array [conc]
40 ! real :: OBC_inflow_conc= 0.0 !< tracer concentration for generic inflows [conc]
41 ! real, dimension(:,:,:), pointer :: OBC_in_u => NULL() !< structured values for flow into the domain
42 ! !! specified in OBCs through u-face of cell
43 ! real, dimension(:,:,:), pointer :: OBC_in_v => NULL() !< structured values for flow into the domain
44 ! !! specified in OBCs through v-face of cell
45 
46  real, dimension(:,:,:), pointer :: ad_x => null() !< diagnostic array for x-advective tracer flux
47  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
48  real, dimension(:,:,:), pointer :: ad_y => null() !< diagnostic array for y-advective tracer flux
49  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
50  real, dimension(:,:), pointer :: ad2d_x => null() !< diagnostic vertical sum x-advective tracer flux
51  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
52  real, dimension(:,:), pointer :: ad2d_y => null() !< diagnostic vertical sum y-advective tracer flux
53  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
54 
55  real, dimension(:,:,:), pointer :: df_x => null() !< diagnostic array for x-diffusive tracer flux
56  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
57  real, dimension(:,:,:), pointer :: df_y => null() !< diagnostic array for y-diffusive tracer flux
58  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
59  real, dimension(:,:), pointer :: df2d_x => null() !< diagnostic vertical sum x-diffusive flux
60  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
61  real, dimension(:,:), pointer :: df2d_y => null() !< diagnostic vertical sum y-diffusive flux
62  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
63  real, dimension(:,:), pointer :: df2d_conc_x => null() !< diagnostic vertical sum x-diffusive content flux
64  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
65  real, dimension(:,:), pointer :: df2d_conc_y => null() !< diagnostic vertical sum y-diffusive content flux
66  !! [conc H m2 s-1 ~> conc m3 s-1 or conc kg s-1]
67 
68  real, dimension(:,:,:), pointer :: advection_xy => null() !< convergence of lateral advective tracer fluxes
69  !! [conc H s-1 ~> conc m s-1 or conc kg m-2 s-1]
70  real, dimension(:,:,:), pointer :: diff_cont_xy => null() !< convergence of lateral diffusive tracer fluxes
71  !! [conc H s-1 ~> conc m s-1 or conc kg m-2 s-1]
72  real, dimension(:,:,:), pointer :: diff_conc_xy => null() !< convergence of lateral diffusive tracer fluxes
73  !! expressed as a change in concentration [conc s-1]
74  real, dimension(:,:,:), pointer :: t_prev => null() !< tracer concentration array at a previous
75  !! timestep used for diagnostics [conc]
76  real, dimension(:,:,:), pointer :: trxh_prev => null() !< layer integrated tracer concentration array
77  !! at a previous timestep used for diagnostics
78 
79  character(len=32) :: name !< tracer name used for diagnostics and error messages
80  character(len=64) :: units !< Physical dimensions of the tracer concentration
81  character(len=240) :: longname !< Long name of the variable
82 ! type(vardesc), pointer :: vd => NULL() !< metadata describing the tracer
83  logical :: registry_diags = .false. !< If true, use the registry to set up the
84  !! diagnostics associated with this tracer.
85  character(len=64) :: cmor_name !< CMOR name of this tracer
86  character(len=64) :: cmor_units !< CMOR physical dimensions of the tracer
87  character(len=240) :: cmor_longname !< CMOR long name of the tracer
88  character(len=32) :: flux_nameroot = "" !< Short tracer name snippet used construct the
89  !! names of flux diagnostics.
90  character(len=64) :: flux_longname = "" !< A word or phrase used construct the long
91  !! names of flux diagnostics.
92  real :: flux_scale= 1.0 !< A scaling factor used to convert the fluxes
93  !! of this tracer to its desired units.
94  character(len=48) :: flux_units = "" !< The units for fluxes of this variable.
95  character(len=48) :: conv_units = "" !< The units for the flux convergence of this tracer.
96  real :: conv_scale = 1.0 !< A scaling factor used to convert the flux
97  !! convergence of this tracer to its desired units.
98  character(len=48) :: cmor_tendprefix = "" !< The CMOR variable prefix for tendencies of this
99  !! tracer, required because CMOR does not follow any
100  !! discernable pattern for these names.
101  integer :: ind_tr_squared = -1 !< The tracer registry index for the square of this tracer
102 
103  !### THESE CAPABILITIES HAVE NOT YET BEEN IMPLEMENTED.
104  logical :: advect_tr = .true. !< If true, this tracer should be advected
105  logical :: hordiff_tr = .true. !< If true, this tracer should experience epineutral diffusion
106  logical :: remap_tr = .true. !< If true, this tracer should be vertically remapped
107 
108  integer :: diag_form = 1 !< An integer indicating which template is to be used to label diagnostics.
109  !>@{ Diagnostic IDs
110  integer :: id_tr = -1
111  integer :: id_adx = -1, id_ady = -1, id_dfx = -1, id_dfy = -1
112  integer :: id_adx_2d = -1, id_ady_2d = -1, id_dfx_2d = -1, id_dfy_2d = -1
113  integer :: id_adv_xy = -1, id_adv_xy_2d = -1
114  integer :: id_dfxy_cont = -1, id_dfxy_cont_2d = -1, id_dfxy_conc = -1
115  integer :: id_remap_conc = -1, id_remap_cont = -1, id_remap_cont_2d = -1
116  integer :: id_tendency = -1, id_trxh_tendency = -1, id_trxh_tendency_2d = -1
117  integer :: id_tr_vardec = -1
118  !!@}
119 end type tracer_type
120 
121 !> Type to carry basic tracer information
122 type, public :: tracer_registry_type
123  integer :: ntr = 0 !< number of registered tracers
124  type(tracer_type) :: tr(max_fields_) !< array of registered tracers
125 ! type(diag_ctrl), pointer :: diag !< structure to regulate timing of diagnostics
126  logical :: locked = .false. !< New tracers may be registered if locked=.false.
127  !! When locked=.true., no more tracers can be registered,
128  !! at which point common diagnostics can be set up
129  !! for the registered tracers
130 end type tracer_registry_type
131 
132 contains
133 
134 !> This subroutine registers a tracer to be advected and laterally diffused.
135 subroutine register_tracer(tr_ptr, Reg, param_file, HI, GV, name, longname, units, &
136  cmor_name, cmor_units, cmor_longname, tr_desc, &
137  OBC_inflow, OBC_in_u, OBC_in_v, ad_x, ad_y, df_x, df_y, &
138  ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy, registry_diags, &
139  flux_nameroot, flux_longname, flux_units, flux_scale, &
140  convergence_units, convergence_scale, cmor_tendprefix, diag_form, &
141  restart_CS, mandatory)
142  type(hor_index_type), intent(in) :: hi !< horizontal index type
143  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
144  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
145  real, dimension(SZI_(HI),SZJ_(HI),SZK_(GV)), &
146  target :: tr_ptr !< target or pointer to the tracer array
147  type(param_file_type), intent(in) :: param_file !< file to parse for model parameter values
148  character(len=*), optional, intent(in) :: name !< Short tracer name
149  character(len=*), optional, intent(in) :: longname !< The long tracer name
150  character(len=*), optional, intent(in) :: units !< The units of this tracer
151  character(len=*), optional, intent(in) :: cmor_name !< CMOR name
152  character(len=*), optional, intent(in) :: cmor_units !< CMOR physical dimensions of variable
153  character(len=*), optional, intent(in) :: cmor_longname !< CMOR long name
154  type(vardesc), optional, intent(in) :: tr_desc !< A structure with metadata about the tracer
155 
156  real, optional, intent(in) :: obc_inflow !< the tracer for all inflows via OBC for which OBC_in_u
157  !! or OBC_in_v are not specified (units of tracer CONC)
158  real, dimension(:,:,:), optional, pointer :: obc_in_u !< tracer at inflows through u-faces of
159  !! tracer cells (units of tracer CONC)
160  real, dimension(:,:,:), optional, pointer :: obc_in_v !< tracer at inflows through v-faces of
161  !! tracer cells (units of tracer CONC)
162 
163  ! The following are probably not necessary if registry_diags is present and true.
164  real, dimension(:,:,:), optional, pointer :: ad_x !< diagnostic x-advective flux (CONC m3/s or CONC*kg/s)
165  real, dimension(:,:,:), optional, pointer :: ad_y !< diagnostic y-advective flux (CONC m3/s or CONC*kg/s)
166  real, dimension(:,:,:), optional, pointer :: df_x !< diagnostic x-diffusive flux (CONC m3/s or CONC*kg/s)
167  real, dimension(:,:,:), optional, pointer :: df_y !< diagnostic y-diffusive flux (CONC m3/s or CONC*kg/s)
168  real, dimension(:,:), optional, pointer :: ad_2d_x !< vert sum of diagnostic x-advect flux
169  !! (CONC m3/s or CONC*kg/s)
170  real, dimension(:,:), optional, pointer :: ad_2d_y !< vert sum of diagnostic y-advect flux
171  !! (CONC m3/s or CONC*kg/s)
172  real, dimension(:,:), optional, pointer :: df_2d_x !< vert sum of diagnostic x-diffuse flux
173  !! (CONC m3/s or CONC*kg/s)
174  real, dimension(:,:), optional, pointer :: df_2d_y !< vert sum of diagnostic y-diffuse flux
175  !! (CONC m3/s or CONC*kg/s)
176 
177  real, dimension(:,:,:), optional, pointer :: advection_xy !< convergence of lateral advective tracer fluxes
178  logical, optional, intent(in) :: registry_diags !< If present and true, use the registry for
179  !! the diagnostics of this tracer.
180  character(len=*), optional, intent(in) :: flux_nameroot !< Short tracer name snippet used construct the
181  !! names of flux diagnostics.
182  character(len=*), optional, intent(in) :: flux_longname !< A word or phrase used construct the long
183  !! names of flux diagnostics.
184  character(len=*), optional, intent(in) :: flux_units !< The units for the fluxes of this tracer.
185  real, optional, intent(in) :: flux_scale !< A scaling factor used to convert the fluxes
186  !! of this tracer to its desired units.
187  character(len=*), optional, intent(in) :: convergence_units !< The units for the flux convergence of
188  !! this tracer.
189  real, optional, intent(in) :: convergence_scale !< A scaling factor used to convert the flux
190  !! convergence of this tracer to its desired units.
191  character(len=*), optional, intent(in) :: cmor_tendprefix !< The CMOR name for the layer-integrated
192  !! tendencies of this tracer.
193  integer, optional, intent(in) :: diag_form !< An integer (1 or 2, 1 by default) indicating the
194  !! character string template to use in
195  !! labeling diagnostics
196  type(mom_restart_cs), optional, pointer :: restart_cs !< A pointer to the restart control structure
197  !! this tracer will be registered for
198  !! restarts if this argument is present
199  logical, optional, intent(in) :: mandatory !< If true, this tracer must be read
200  !! from a restart file.
201 
202  logical :: mand
203  type(tracer_type), pointer :: tr=>null()
204  character(len=256) :: mesg ! Message for error messages.
205 
206  if (.not. associated(reg)) call tracer_registry_init(param_file, reg)
207 
208  if (reg%ntr>=max_fields_) then
209  write(mesg,'("Increase MAX_FIELDS_ in MOM_memory.h to at least ",I3," to allow for &
210  &all the tracers being registered via register_tracer.")') reg%ntr+1
211  call mom_error(fatal,"MOM register_tracer: "//mesg)
212  endif
213  reg%ntr = reg%ntr + 1
214 
215  tr => reg%Tr(reg%ntr)
216 
217  if (present(name)) then
218  tr%name = name
219  tr%longname = name ; if (present(longname)) tr%longname = longname
220  tr%units = "Conc" ; if (present(units)) tr%units = units
221 
222  tr%cmor_name = ""
223  if (present(cmor_name)) tr%cmor_name = cmor_name
224 
225  tr%cmor_units = tr%units
226  if (present(cmor_units)) tr%cmor_units = cmor_units
227 
228  tr%cmor_longname = ""
229  if (present(cmor_longname)) tr%cmor_longname = cmor_longname
230 
231  if (present(tr_desc)) call mom_error(warning, "MOM register_tracer: "//&
232  "It is a bad idea to use both name and tr_desc when registring "//trim(name))
233  elseif (present(tr_desc)) then
234  call query_vardesc(tr_desc, name=tr%name, units=tr%units, &
235  longname=tr%longname, cmor_field_name=tr%cmor_name, &
236  cmor_longname=tr%cmor_longname, caller="register_tracer")
237  tr%cmor_units = tr%units
238  else
239  call mom_error(fatal,"MOM register_tracer: Either name or "//&
240  "tr_desc must be present when registering a tracer.")
241  endif
242 
243  if (reg%locked) call mom_error(fatal, &
244  "MOM register_tracer was called for variable "//trim(tr%name)//&
245  " with a locked tracer registry.")
246 
247  tr%flux_nameroot = tr%name
248  if (present(flux_nameroot)) then
249  if (len_trim(flux_nameroot) > 0) tr%flux_nameroot = flux_nameroot
250  endif
251 
252  tr%flux_longname = tr%longname
253  if (present(flux_longname)) then
254  if (len_trim(flux_longname) > 0) tr%flux_longname = flux_longname
255  endif
256 
257  tr%flux_units = ""
258  if (present(flux_units)) tr%flux_units = flux_units
259 
260  tr%flux_scale = 1.0
261  if (present(flux_scale)) tr%flux_scale = flux_scale
262 
263  tr%conv_units = ""
264  if (present(convergence_units)) tr%conv_units = convergence_units
265 
266  tr%cmor_tendprefix = ""
267  if (present(cmor_tendprefix)) tr%cmor_tendprefix = cmor_tendprefix
268 
269  tr%conv_scale = 1.0
270  if (present(convergence_scale)) then
271  tr%conv_scale = convergence_scale
272  elseif (present(flux_scale)) then
273  tr%conv_scale = flux_scale
274  endif
275 
276  tr%diag_form = 1
277  if (present(diag_form)) tr%diag_form = diag_form
278 
279  tr%t => tr_ptr
280 
281  if (present(registry_diags)) tr%registry_diags = registry_diags
282 
283  if (present(ad_x)) then ; if (associated(ad_x)) tr%ad_x => ad_x ; endif
284  if (present(ad_y)) then ; if (associated(ad_y)) tr%ad_y => ad_y ; endif
285  if (present(df_x)) then ; if (associated(df_x)) tr%df_x => df_x ; endif
286  if (present(df_y)) then ; if (associated(df_y)) tr%df_y => df_y ; endif
287 ! if (present(OBC_inflow)) Tr%OBC_inflow_conc = OBC_inflow
288 ! if (present(OBC_in_u)) then ; if (associated(OBC_in_u)) &
289 ! Tr%OBC_in_u => OBC_in_u ; endif
290 ! if (present(OBC_in_v)) then ; if (associated(OBC_in_v)) &
291 ! Tr%OBC_in_v => OBC_in_v ; endif
292  if (present(ad_2d_x)) then ; if (associated(ad_2d_x)) tr%ad2d_x => ad_2d_x ; endif
293  if (present(ad_2d_y)) then ; if (associated(ad_2d_y)) tr%ad2d_y => ad_2d_y ; endif
294  if (present(df_2d_x)) then ; if (associated(df_2d_x)) tr%df2d_x => df_2d_x ; endif
295 
296  if (present(advection_xy)) then ; if (associated(advection_xy)) tr%advection_xy => advection_xy ; endif
297 
298  if (present(restart_cs)) then ; if (associated(restart_cs)) then
299  ! Register this tracer to be read from and written to restart files.
300  mand = .true. ; if (present(mandatory)) mand = mandatory
301 
302  call register_restart_field(tr_ptr, tr%name, mand, restart_cs, &
303  longname=tr%longname, units=tr%units)
304  endif ; endif
305 
306 end subroutine register_tracer
307 
308 
309 !> This subroutine locks the tracer registry to prevent the addition of more
310 !! tracers. After locked=.true., can then register common diagnostics.
311 subroutine lock_tracer_registry(Reg)
312  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
313 
314  if (.not. associated(reg)) call mom_error(warning, &
315  "lock_tracer_registry called with an unassociated registry.")
316 
317  reg%locked = .true.
318 
319 end subroutine lock_tracer_registry
320 
321 !> register_tracer_diagnostics does a set of register_diag_field calls for any previously
322 !! registered in a tracer registry with a value of registry_diags set to .true.
323 subroutine register_tracer_diagnostics(Reg, h, Time, diag, G, GV, use_ALE)
324  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
325  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
326  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
327  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
328  intent(in) :: h !< Layer thicknesses
329  type(time_type), intent(in) :: time !< current model time
330  type(diag_ctrl), intent(in) :: diag !< structure to regulate diagnostic output
331  logical, intent(in) :: use_ale !< If true active diagnostics that only
332  !! apply to ALE configurations
333 
334  character(len=24) :: name ! A variable's name in a NetCDF file.
335  character(len=24) :: shortnm ! A shortened version of a variable's name for
336  ! creating additional diagnostics.
337  character(len=72) :: longname ! The long name of that tracer variable.
338  character(len=72) :: flux_longname ! The tracer name in the long names of fluxes.
339  character(len=48) :: units ! The dimensions of the tracer.
340  character(len=48) :: flux_units ! The units for fluxes, either
341  ! [units] m3 s-1 or [units] kg s-1.
342  character(len=48) :: conv_units ! The units for flux convergences, either
343  ! [units] m2 s-1 or [units] kg s-1.
344  character(len=48) :: unit2 ! The dimensions of the tracer squared
345  character(len=72) :: cmorname ! The CMOR name of this tracer.
346  character(len=120) :: cmor_longname ! The CMOR long name of that variable.
347  character(len=120) :: var_lname ! A temporary longname for a diagnostic.
348  character(len=120) :: cmor_var_lname ! The temporary CMOR long name for a diagnostic
349  character(len=72) :: cmor_varname ! The temporary CMOR name for a diagnostic
350  type(tracer_type), pointer :: tr=>null()
351  integer :: i, j, k, is, ie, js, je, nz, m, m2, ntr_in
352  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
353  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
354  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
355  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
356 
357  if (.not. associated(reg)) call mom_error(fatal, "register_tracer_diagnostics: "//&
358  "register_tracer must be called before register_tracer_diagnostics")
359 
360  ntr_in = reg%ntr
361 
362  do m=1,ntr_in ; if (reg%Tr(m)%registry_diags) then
363  tr => reg%Tr(m)
364 ! call query_vardesc(Tr%vd, name, units=units, longname=longname, &
365 ! cmor_field_name=cmorname, cmor_longname=cmor_longname, &
366 ! caller="register_tracer_diagnostics")
367  name = tr%name ; units=adjustl(tr%units) ; longname = tr%longname
368  cmorname = tr%cmor_name ; cmor_longname = tr%cmor_longname
369  shortnm = tr%flux_nameroot
370  flux_longname = tr%flux_longname
371  if (len_trim(cmor_longname) == 0) cmor_longname = longname
372 
373  if (len_trim(tr%flux_units) > 0) then ; flux_units = tr%flux_units
374  elseif (gv%Boussinesq) then ; flux_units = trim(units)//" m3 s-1"
375  else ; flux_units = trim(units)//" kg s-1" ; endif
376 
377  if (len_trim(tr%conv_units) > 0) then ; conv_units = tr%conv_units
378  elseif (gv%Boussinesq) then ; conv_units = trim(units)//" m s-1"
379  else ; conv_units = trim(units)//" kg m-2 s-1" ; endif
380 
381  if (len_trim(cmorname) == 0) then
382  tr%id_tr = register_diag_field("ocean_model", trim(name), diag%axesTL, &
383  time, trim(longname), trim(units))
384  else
385  tr%id_tr = register_diag_field("ocean_model", trim(name), diag%axesTL, &
386  time, trim(longname), trim(units), cmor_field_name=cmorname, &
387  cmor_long_name=cmor_longname, cmor_units=tr%cmor_units, &
388  cmor_standard_name=cmor_long_std(cmor_longname))
389  endif
390  if (tr%diag_form == 1) then
391  tr%id_adx = register_diag_field("ocean_model", trim(shortnm)//"_adx", &
392  diag%axesCuL, time, trim(flux_longname)//" advective zonal flux" , &
393  trim(flux_units), v_extensive = .true., y_cell_method = 'sum')
394  tr%id_ady = register_diag_field("ocean_model", trim(shortnm)//"_ady", &
395  diag%axesCvL, time, trim(flux_longname)//" advective meridional flux" , &
396  trim(flux_units), v_extensive = .true., x_cell_method = 'sum')
397  tr%id_dfx = register_diag_field("ocean_model", trim(shortnm)//"_dfx", &
398  diag%axesCuL, time, trim(flux_longname)//" diffusive zonal flux" , &
399  trim(flux_units), v_extensive = .true., y_cell_method = 'sum')
400  tr%id_dfy = register_diag_field("ocean_model", trim(shortnm)//"_dfy", &
401  diag%axesCvL, time, trim(flux_longname)//" diffusive zonal flux" , &
402  trim(flux_units), v_extensive = .true., x_cell_method = 'sum')
403  else
404  tr%id_adx = register_diag_field("ocean_model", trim(shortnm)//"_adx", &
405  diag%axesCuL, time, "Advective (by residual mean) Zonal Flux of "//trim(flux_longname), &
406  flux_units, v_extensive=.true., conversion=tr%flux_scale, y_cell_method = 'sum')
407  tr%id_ady = register_diag_field("ocean_model", trim(shortnm)//"_ady", &
408  diag%axesCvL, time, "Advective (by residual mean) Meridional Flux of "//trim(flux_longname), &
409  flux_units, v_extensive=.true., conversion=tr%flux_scale, x_cell_method = 'sum')
410  tr%id_dfx = register_diag_field("ocean_model", trim(shortnm)//"_diffx", &
411  diag%axesCuL, time, "Diffusive Zonal Flux of "//trim(flux_longname), &
412  flux_units, v_extensive=.true., conversion=tr%flux_scale, y_cell_method = 'sum')
413  tr%id_dfy = register_diag_field("ocean_model", trim(shortnm)//"_diffy", &
414  diag%axesCvL, time, "Diffusive Meridional Flux of "//trim(flux_longname), &
415  flux_units, v_extensive=.true., conversion=tr%flux_scale, x_cell_method = 'sum')
416  endif
417  if (tr%id_adx > 0) call safe_alloc_ptr(tr%ad_x,isdb,iedb,jsd,jed,nz)
418  if (tr%id_ady > 0) call safe_alloc_ptr(tr%ad_y,isd,ied,jsdb,jedb,nz)
419  if (tr%id_dfx > 0) call safe_alloc_ptr(tr%df_x,isdb,iedb,jsd,jed,nz)
420  if (tr%id_dfy > 0) call safe_alloc_ptr(tr%df_y,isd,ied,jsdb,jedb,nz)
421 
422  tr%id_adx_2d = register_diag_field("ocean_model", trim(shortnm)//"_adx_2d", &
423  diag%axesCu1, time, &
424  "Vertically Integrated Advective Zonal Flux of "//trim(flux_longname), &
425  flux_units, conversion=tr%flux_scale, y_cell_method = 'sum')
426  tr%id_ady_2d = register_diag_field("ocean_model", trim(shortnm)//"_ady_2d", &
427  diag%axesCv1, time, &
428  "Vertically Integrated Advective Meridional Flux of "//trim(flux_longname), &
429  flux_units, conversion=tr%flux_scale, x_cell_method = 'sum')
430  tr%id_dfx_2d = register_diag_field("ocean_model", trim(shortnm)//"_diffx_2d", &
431  diag%axesCu1, time, &
432  "Vertically Integrated Diffusive Zonal Flux of "//trim(flux_longname), &
433  flux_units, conversion=tr%flux_scale, y_cell_method = 'sum')
434  tr%id_dfy_2d = register_diag_field("ocean_model", trim(shortnm)//"_diffy_2d", &
435  diag%axesCv1, time, &
436  "Vertically Integrated Diffusive Meridional Flux of "//trim(flux_longname), &
437  flux_units, conversion=tr%flux_scale, x_cell_method = 'sum')
438 
439  if (tr%id_adx_2d > 0) call safe_alloc_ptr(tr%ad2d_x,isdb,iedb,jsd,jed)
440  if (tr%id_ady_2d > 0) call safe_alloc_ptr(tr%ad2d_y,isd,ied,jsdb,jedb)
441  if (tr%id_dfx_2d > 0) call safe_alloc_ptr(tr%df2d_x,isdb,iedb,jsd,jed)
442  if (tr%id_dfy_2d > 0) call safe_alloc_ptr(tr%df2d_y,isd,ied,jsdb,jedb)
443 
444  tr%id_adv_xy = register_diag_field('ocean_model', trim(shortnm)//"_advection_xy", &
445  diag%axesTL, time, &
446  'Horizontal convergence of residual mean advective fluxes of '//&
447  trim(lowercase(flux_longname)), conv_units, v_extensive=.true., &
448  conversion=tr%conv_scale)
449  tr%id_adv_xy_2d = register_diag_field('ocean_model', trim(shortnm)//"_advection_xy_2d", &
450  diag%axesT1, time, &
451  'Vertical sum of horizontal convergence of residual mean advective fluxes of '//&
452  trim(lowercase(flux_longname)), conv_units, conversion=tr%conv_scale)
453  if ((tr%id_adv_xy > 0) .or. (tr%id_adv_xy_2d > 0)) &
454  call safe_alloc_ptr(tr%advection_xy,isd,ied,jsd,jed,nz)
455 
456  tr%id_tendency = register_diag_field('ocean_model', trim(shortnm)//'_tendency', &
457  diag%axesTL, time, &
458  'Net time tendency for '//trim(lowercase(longname)), trim(units)//' s-1')
459 
460  if (tr%id_tendency > 0) then
461  call safe_alloc_ptr(tr%t_prev,isd,ied,jsd,jed,nz)
462  do k=1,nz ; do j=js,je ; do i=is,ie
463  tr%t_prev(i,j,k) = tr%t(i,j,k)
464  enddo ; enddo ; enddo
465  endif
466 
467  ! Lateral diffusion convergence tendencies
468  if (tr%diag_form == 1) then
469  tr%id_dfxy_cont = register_diag_field("ocean_model", trim(shortnm)//'_dfxy_cont_tendency', &
470  diag%axesTL, time, "Lateral or neutral diffusion tracer content tendency for "//trim(shortnm), &
471  conv_units, conversion=tr%conv_scale, x_cell_method='sum', y_cell_method='sum', v_extensive=.true.)
472 
473  tr%id_dfxy_cont_2d = register_diag_field("ocean_model", trim(shortnm)//'_dfxy_cont_tendency_2d', &
474  diag%axesT1, time, "Depth integrated lateral or neutral diffusion tracer concentration "//&
475  "tendency for "//trim(shortnm), conv_units, conversion = tr%conv_scale, &
476  x_cell_method = 'sum', y_cell_method = 'sum')
477  else
478  cmor_var_lname = 'Tendency of '//trim(lowercase(cmor_longname))//&
479  ' expressed as '//trim(lowercase(flux_longname))//&
480  ' content due to parameterized mesoscale diffusion'
481  tr%id_dfxy_cont = register_diag_field("ocean_model", trim(shortnm)//'_dfxy_cont_tendency', &
482  diag%axesTL, time, "Lateral or neutral diffusion tracer concentration tendency for "//trim(shortnm), &
483  conv_units, conversion = tr%conv_scale, cmor_field_name = trim(tr%cmor_tendprefix)//'pmdiff', &
484  cmor_long_name = trim(cmor_var_lname), cmor_standard_name = trim(cmor_long_std(cmor_var_lname)), &
485  x_cell_method = 'sum', y_cell_method = 'sum', v_extensive = .true.)
486 
487  cmor_var_lname = 'Tendency of '//trim(lowercase(cmor_longname))//' expressed as '//&
488  trim(lowercase(flux_longname))//' content due to parameterized mesoscale diffusion'
489  tr%id_dfxy_cont_2d = register_diag_field("ocean_model", trim(shortnm)//'_dfxy_cont_tendency_2d', &
490  diag%axesT1, time, "Depth integrated lateral or neutral diffusion tracer "//&
491  "concentration tendency for "//trim(shortnm), conv_units, &
492  conversion=tr%conv_scale, cmor_field_name=trim(tr%cmor_tendprefix)//'pmdiff_2d', &
493  cmor_long_name=trim(cmor_var_lname), cmor_standard_name=trim(cmor_long_std(cmor_var_lname)), &
494  x_cell_method='sum', y_cell_method='sum')
495  endif
496  tr%id_dfxy_conc = register_diag_field("ocean_model", trim(shortnm)//'_dfxy_conc_tendency', &
497  diag%axesTL, time, "Lateral (neutral) tracer concentration tendency for "//trim(shortnm), &
498  trim(units)//' s-1')
499 
500  var_lname = "Net time tendency for "//lowercase(flux_longname)
501  if (len_trim(tr%cmor_tendprefix) == 0) then
502  tr%id_trxh_tendency = register_diag_field('ocean_model', trim(shortnm)//'h_tendency', &
503  diag%axesTL, time, var_lname, conv_units, v_extensive=.true.)
504  tr%id_trxh_tendency_2d = register_diag_field('ocean_model', trim(shortnm)//'h_tendency_2d', &
505  diag%axesT1, time, "Vertical sum of "//trim(lowercase(var_lname)), conv_units)
506  else
507  cmor_var_lname = "Tendency of "//trim(cmor_longname)//" Expressed as "//&
508  trim(flux_longname)//" Content"
509  tr%id_trxh_tendency = register_diag_field('ocean_model', trim(shortnm)//'h_tendency', &
510  diag%axesTL, time, var_lname, conv_units, &
511  cmor_field_name=trim(tr%cmor_tendprefix)//"tend", &
512  cmor_standard_name=cmor_long_std(cmor_var_lname), cmor_long_name=cmor_var_lname, &
513  v_extensive=.true., conversion=tr%conv_scale)
514  cmor_var_lname = trim(cmor_var_lname)//" Vertical Sum"
515  tr%id_trxh_tendency_2d = register_diag_field('ocean_model', trim(shortnm)//'h_tendency_2d', &
516  diag%axesT1, time, "Vertical sum of "//trim(lowercase(var_lname)), conv_units, &
517  cmor_field_name=trim(tr%cmor_tendprefix)//"tend_2d", &
518  cmor_standard_name=cmor_long_std(cmor_var_lname), cmor_long_name=cmor_var_lname, &
519  conversion=tr%conv_scale)
520  endif
521  if ((tr%id_trxh_tendency > 0) .or. (tr%id_trxh_tendency_2d > 0)) then
522  call safe_alloc_ptr(tr%Trxh_prev,isd,ied,jsd,jed,nz)
523  do k=1,nz ; do j=js,je ; do i=is,ie
524  tr%Trxh_prev(i,j,k) = tr%t(i,j,k) * h(i,j,k)
525  enddo ; enddo ; enddo
526  endif
527 
528  ! Vertical regridding/remapping tendencies
529  if (use_ale .and. tr%remap_tr) then
530  var_lname = "Vertical remapping tracer concentration tendency for "//trim(reg%Tr(m)%name)
531  tr%id_remap_conc= register_diag_field('ocean_model', &
532  trim(tr%flux_nameroot)//'_tendency_vert_remap', diag%axesTL, time, var_lname, &
533  trim(units)//' s-1')
534 
535  var_lname = "Vertical remapping tracer content tendency for "//trim(reg%Tr(m)%flux_longname)
536  tr%id_remap_cont = register_diag_field('ocean_model', &
537  trim(tr%flux_nameroot)//'h_tendency_vert_remap', &
538  diag%axesTL, time, var_lname, flux_units, v_extensive=.true., conversion = tr%conv_scale)
539 
540  var_lname = "Vertical sum of vertical remapping tracer content tendency for "//&
541  trim(reg%Tr(m)%flux_longname)
542  tr%id_remap_cont_2d = register_diag_field('ocean_model', &
543  trim(tr%flux_nameroot)//'h_tendency_vert_remap_2d', &
544  diag%axesT1, time, var_lname, flux_units, conversion = tr%conv_scale)
545 
546  endif
547 
548  if (use_ale .and. (reg%ntr<max_fields_) .and. tr%remap_tr) then
549  unit2 = trim(units)//"2"
550  if (index(units(1:len_trim(units))," ") > 0) unit2 = "("//trim(units)//")2"
551  tr%id_tr_vardec = register_diag_field('ocean_model', trim(shortnm)//"_vardec", diag%axesTL, &
552  time, "ALE variance decay for "//lowercase(longname), trim(unit2)//" s-1")
553  if (tr%id_tr_vardec > 0) then
554  ! Set up a new tracer for this tracer squared
555  m2 = reg%ntr+1
556  tr%ind_tr_squared = m2
557  call safe_alloc_ptr(reg%Tr(m2)%t,isd,ied,jsd,jed,nz) ; reg%Tr(m2)%t(:,:,:) = 0.0
558  reg%Tr(m2)%name = trim(shortnm)//"2"
559  reg%Tr(m2)%longname = "Squared "//trim(longname)
560  reg%Tr(m2)%units = unit2
561  reg%Tr(m2)%registry_diags = .false.
562  reg%Tr(m2)%ind_tr_squared = -1
563  ! Augment the total number of tracers, including the squared tracers.
564  reg%ntr = reg%ntr + 1
565  endif
566  endif
567 
568  endif ; enddo
569 
570 end subroutine register_tracer_diagnostics
571 
572 subroutine preale_tracer_diagnostics(Reg, G, GV)
573  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
574  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
575  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
576 
577  integer :: i, j, k, is, ie, js, je, nz, m, m2
578  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
579 
580  do m=1,reg%ntr ; if (reg%Tr(m)%ind_tr_squared > 0) then
581  m2 = reg%Tr(m)%ind_tr_squared
582  ! Update squared quantities
583  do k=1,nz ; do j=js,je ; do i=is,ie
584  reg%Tr(m2)%T(i,j,k) = reg%Tr(m)%T(i,j,k)**2
585  enddo ; enddo ; enddo
586  endif ; enddo
587 
588 end subroutine preale_tracer_diagnostics
589 
590 subroutine postale_tracer_diagnostics(Reg, G, GV, diag, dt)
591  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
592  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
593  type(verticalgrid_type), intent(in) :: gv !< ocean vertical grid structure
594  type(diag_ctrl), intent(in) :: diag !< regulates diagnostic output
595  real, intent(in) :: dt !< total time interval for these diagnostics
596 
597  real :: work(szi_(g),szj_(g),szk_(g))
598  real :: idt
599  integer :: i, j, k, is, ie, js, je, nz, m, m2
600  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
601 
602  ! The "if" is to avoid NaNs if the diagnostic is called for a zero length interval
603  idt = 0.0 ; if (dt /= 0.0) idt = 1.0 / dt
604 
605  do m=1,reg%ntr ; if (reg%Tr(m)%id_tr_vardec > 0) then
606  m2 = reg%Tr(m)%ind_tr_squared
607  if (m2 < 1) call mom_error(fatal, "Bad value of Tr%ind_tr_squared for "//trim(reg%Tr(m)%name))
608  ! Update squared quantities
609  do k=1,nz ; do j=js,je ; do i=is,ie
610  work(i,j,k) = (reg%Tr(m2)%T(i,j,k) - reg%Tr(m)%T(i,j,k)**2) * idt
611  enddo ; enddo ; enddo
612  call post_data(reg%Tr(m)%id_tr_vardec, work, diag)
613  endif ; enddo
614 
615 end subroutine postale_tracer_diagnostics
616 
617 !> post_tracer_diagnostics does post_data calls for any diagnostics that are
618 !! being handled via the tracer registry.
619 subroutine post_tracer_diagnostics(Reg, h, diag_prev, diag, G, GV, dt)
620  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
621  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
622  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
623  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
624  intent(in) :: h !< Layer thicknesses
625  type(diag_grid_storage), intent(in) :: diag_prev !< Contains diagnostic grids from previous timestep
626  type(diag_ctrl), intent(inout) :: diag !< structure to regulate diagnostic output
627  real, intent(in) :: dt !< total time step for tracer updates
628 
629  real :: work3d(szi_(g),szj_(g),szk_(gv))
630  real :: work2d(szi_(g),szj_(g))
631  real :: idt
632  type(tracer_type), pointer :: tr=>null()
633  integer :: i, j, k, is, ie, js, je, nz, m
634  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
635 
636  idt = 0.; if (dt/=0.) idt = 1.0 / dt ! The "if" is in case the diagnostic is called for a zero length interval
637 
638  ! Tendency diagnostics need to be posted on the grid from the last call to this routine
639  call diag_save_grids(diag)
640  call diag_copy_storage_to_diag(diag, diag_prev)
641  do m=1,reg%ntr ; if (reg%Tr(m)%registry_diags) then
642  tr => reg%Tr(m)
643  if (tr%id_tendency > 0) then
644  work3d(:,:,:) = 0.0
645  do k=1,nz ; do j=js,je ; do i=is,ie
646  work3d(i,j,k) = (tr%t(i,j,k) - tr%t_prev(i,j,k))*idt
647  tr%t_prev(i,j,k) = tr%t(i,j,k)
648  enddo ; enddo ; enddo
649  call post_data(tr%id_tendency, work3d, diag, alt_h = diag_prev%h_state)
650  endif
651  if ((tr%id_trxh_tendency > 0) .or. (tr%id_trxh_tendency_2d > 0)) then
652  do k=1,nz ; do j=js,je ; do i=is,ie
653  work3d(i,j,k) = (tr%t(i,j,k)*h(i,j,k) - tr%Trxh_prev(i,j,k)) * idt
654  tr%Trxh_prev(i,j,k) = tr%t(i,j,k) * h(i,j,k)
655  enddo ; enddo ; enddo
656  if (tr%id_trxh_tendency > 0) call post_data(tr%id_trxh_tendency, work3d, diag, alt_h = diag_prev%h_state)
657  if (tr%id_trxh_tendency_2d > 0) then
658  work2d(:,:) = 0.0
659  do k=1,nz ; do j=js,je ; do i=is,ie
660  work2d(i,j) = work2d(i,j) + work3d(i,j,k)
661  enddo ; enddo ; enddo
662  call post_data(tr%id_trxh_tendency_2d, work2d, diag)
663  endif
664  endif
665  endif ; enddo
666  call diag_restore_grids(diag)
667 
668 end subroutine post_tracer_diagnostics
669 
670 !> Post the advective and diffusive tendencies
671 subroutine post_tracer_transport_diagnostics(G, GV, Reg, h_diag, diag)
672  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
673  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
674  type(tracer_registry_type), pointer :: reg !< pointer to the tracer registry
675  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
676  intent(in) :: h_diag !< Layer thicknesses on which to post fields
677  type(diag_ctrl), intent(in) :: diag !< structure to regulate diagnostic output
678 
679  integer :: i, j, k, is, ie, js, je, nz, m
680  real :: work2d(szi_(g),szj_(g))
681  type(tracer_type), pointer :: tr=>null()
682 
683  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
684 
685  do m=1,reg%ntr ; if (reg%Tr(m)%registry_diags) then
686  tr => reg%Tr(m)
687  if (tr%id_tr > 0) call post_data(tr%id_tr, tr%t, diag)
688  if (tr%id_adx > 0) call post_data(tr%id_adx, tr%ad_x, diag, alt_h = h_diag)
689  if (tr%id_ady > 0) call post_data(tr%id_ady, tr%ad_y, diag, alt_h = h_diag)
690  if (tr%id_dfx > 0) call post_data(tr%id_dfx, tr%df_x, diag, alt_h = h_diag)
691  if (tr%id_dfy > 0) call post_data(tr%id_dfy, tr%df_y, diag, alt_h = h_diag)
692  if (tr%id_adx_2d > 0) call post_data(tr%id_adx_2d, tr%ad2d_x, diag)
693  if (tr%id_ady_2d > 0) call post_data(tr%id_ady_2d, tr%ad2d_y, diag)
694  if (tr%id_dfx_2d > 0) call post_data(tr%id_dfx_2d, tr%df2d_x, diag)
695  if (tr%id_dfy_2d > 0) call post_data(tr%id_dfy_2d, tr%df2d_y, diag)
696  if (tr%id_adv_xy > 0) call post_data(tr%id_adv_xy, tr%advection_xy, diag, alt_h = h_diag)
697  if (tr%id_adv_xy_2d > 0) then
698  work2d(:,:) = 0.0
699  do k=1,nz ; do j=js,je ; do i=is,ie
700  work2d(i,j) = work2d(i,j) + tr%advection_xy(i,j,k)
701  enddo ; enddo ; enddo
702  call post_data(tr%id_adv_xy_2d, work2d, diag)
703  endif
704  endif ; enddo
705 
706 end subroutine post_tracer_transport_diagnostics
707 
708 !> This subroutine writes out chksums for tracers.
709 subroutine mom_tracer_chksum(mesg, Tr, ntr, G)
710  character(len=*), intent(in) :: mesg !< message that appears on the chksum lines
711  type(tracer_type), intent(in) :: tr(:) !< array of all of registered tracers
712  integer, intent(in) :: ntr !< number of registered tracers
713  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
714 
715  integer :: is, ie, js, je, nz
716  integer :: i, j, k, m
717 
718  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
719  do m=1,ntr
720  call hchksum(tr(m)%t, mesg//trim(tr(m)%name), g%HI)
721  enddo
722 
723 end subroutine mom_tracer_chksum
724 
725 !> Calculates and prints the global inventory of all tracers in the registry.
726 subroutine mom_tracer_chkinv(mesg, G, h, Tr, ntr)
727  character(len=*), intent(in) :: mesg !< message that appears on the chksum lines
728  type(ocean_grid_type), intent(in) :: g !< ocean grid structure
729  type(tracer_type), dimension(:), intent(in) :: tr !< array of all of registered tracers
730  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses
731  integer, intent(in) :: ntr !< number of registered tracers
732 
733  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: tr_inv !< Tracer inventory
734  real :: total_inv
735  integer :: is, ie, js, je, nz
736  integer :: i, j, k, m
737 
738  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
739  do m=1,ntr
740  do k=1,nz ; do j=js,je ; do i=is,ie
741  tr_inv(i,j,k) = tr(m)%t(i,j,k)*h(i,j,k)*g%areaT(i,j)*g%mask2dT(i,j)
742  enddo ; enddo ; enddo
743  total_inv = reproducing_sum(tr_inv, is+(1-g%isd), ie+(1-g%isd), js+(1-g%jsd), je+(1-g%jsd))
744  if (is_root_pe()) write(0,'(A,1X,A5,1X,ES25.16,1X,A)') "h-point: inventory", tr(m)%name, total_inv, mesg
745  enddo
746 
747 end subroutine mom_tracer_chkinv
748 
749 !> Find a tracer in the tracer registry by name.
750 subroutine tracer_name_lookup(Reg, tr_ptr, name)
751  type(tracer_registry_type), pointer :: reg !< pointer to tracer registry
752  type(tracer_type), pointer :: tr_ptr !< target or pointer to the tracer array
753  character(len=32), intent(in) :: name !< tracer name
754 
755  integer n
756  do n=1,reg%ntr
757  if (lowercase(reg%Tr(n)%name) == lowercase(name)) tr_ptr => reg%Tr(n)
758  enddo
759 
760 end subroutine tracer_name_lookup
761 
762 !> Initialize the tracer registry.
763 subroutine tracer_registry_init(param_file, Reg)
764  type(param_file_type), intent(in) :: param_file !< open file to parse for model parameters
765  type(tracer_registry_type), pointer :: reg !< pointer to tracer registry
766 
767  integer, save :: init_calls = 0
768 
769 ! This include declares and sets the variable "version".
770 #include "version_variable.h"
771  character(len=40) :: mdl = "MOM_tracer_registry" ! This module's name.
772  character(len=256) :: mesg ! Message for error messages.
773 
774  if (.not.associated(reg)) then ; allocate(reg)
775  else ; return ; endif
776 
777  ! Read all relevant parameters and write them to the model log.
778  call log_version(param_file, mdl, version, "")
779 
780  init_calls = init_calls + 1
781  if (init_calls > 1) then
782  write(mesg,'("tracer_registry_init called ",I3, &
783  &" times with different registry pointers.")') init_calls
784  if (is_root_pe()) call mom_error(warning,"MOM_tracer"//mesg)
785  endif
786 
787 end subroutine tracer_registry_init
788 
789 
790 !> This routine closes the tracer registry module.
791 subroutine tracer_registry_end(Reg)
792  type(tracer_registry_type), pointer :: reg !< The tracer registry that will be deallocated
793  if (associated(reg)) deallocate(reg)
794 end subroutine tracer_registry_end
795 
796 end module mom_tracer_registry
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_tracer_registry::tracer_type
The tracer type.
Definition: MOM_tracer_registry.F90:37
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_tracer_registry
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Definition: MOM_tracer_registry.F90:5
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_hor_index
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Definition: MOM_hor_index.F90:2
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_restart::mom_restart_cs
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:72
mom_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_diag_mediator::diag_grid_storage
Stores all the remapping grids and the model's native space thicknesses.
Definition: MOM_diag_mediator.F90:146
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_restart
The MOM6 facility for reading and writing restart files, and querying what has been read.
Definition: MOM_restart.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_tracer_registry::tracer_registry_type
Type to carry basic tracer information.
Definition: MOM_tracer_registry.F90:122
mom_hor_index::hor_index_type
Container for horizontal index ranges for data, computational and global domains.
Definition: MOM_hor_index.F90:15
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_restart::register_restart_field
Register fields for restarts.
Definition: MOM_restart.F90:107
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_io::vardesc
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
mom_coms::reproducing_sum
Find an accurate and order-invariant sum of distributed 2d or 3d fields.
Definition: MOM_coms.F90:54
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239