MOM6
ideal_age_example Module Reference

Detailed Description

A tracer package of ideal age tracers.

Originally by Robert Hallberg, 2002

This file contains an example of the code that is needed to set up and use a set (in this case two) of dynamically passive tracers for diagnostic purposes. The tracers here are an ideal age tracer that ages at a rate of 1/year once it is isolated from the surface, and a vintage tracer, whose surface concentration grows exponen- with time with a 30-year timescale (similar to CFCs).

Data Types

type  ideal_age_tracer_cs
 The control structure for the ideal_age_tracer package. More...
 

Functions/Subroutines

logical function, public register_ideal_age_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Register the ideal age tracer fields to be used with MOM. More...
 
subroutine, public initialize_ideal_age_tracer (restart, day, G, GV, US, h, diag, OBC, CS, sponge_CSp)
 Sets the ideal age traces to their initial values and sets up the tracer output. More...
 
subroutine, public ideal_age_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
 Applies diapycnal diffusion, aging and regeneration at the surface to the ideal age tracers. More...
 
integer function, public ideal_age_stock (h, stocks, G, GV, CS, names, units, stock_index)
 Calculates the mass-weighted integral of all tracer stocks, returning the number of stocks it has calculated. If stock_index is present, only the stock corresponding to that coded index is found. More...
 
subroutine, public ideal_age_tracer_surface_state (state, h, G, CS)
 This subroutine extracts the surface fields from this tracer package that are to be shared with the atmosphere in coupled configurations. This particular tracer package does not report anything back to the coupler. More...
 
subroutine, public ideal_age_example_end (CS)
 Deallocate any memory associated with this tracer package. More...
 

Variables

integer, parameter ntr_max = 3
 the maximum number of tracers in this module.
 

Function/Subroutine Documentation

◆ ideal_age_example_end()

subroutine, public ideal_age_example::ideal_age_example_end ( type(ideal_age_tracer_cs), pointer  CS)

Deallocate any memory associated with this tracer package.

Parameters
csThe control structure returned by a previous call to register_ideal_age_tracer.

Definition at line 453 of file ideal_age_example.F90.

453  type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
454  !! call to register_ideal_age_tracer.
455 
456  integer :: m
457 
458  if (associated(cs)) then
459  if (associated(cs%tr)) deallocate(cs%tr)
460  deallocate(cs)
461  endif

◆ ideal_age_stock()

integer function, public ideal_age_example::ideal_age_stock ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
real, dimension(:), intent(out)  stocks,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(ideal_age_tracer_cs), pointer  CS,
character(len=*), dimension(:), intent(out)  names,
character(len=*), dimension(:), intent(out)  units,
integer, intent(in), optional  stock_index 
)

Calculates the mass-weighted integral of all tracer stocks, returning the number of stocks it has calculated. If stock_index is present, only the stock corresponding to that coded index is found.

Parameters
[in]gThe ocean's grid structure
[in]hLayer thicknesses [H ~> m or kg m-2]
[out]stocksthe mass-weighted integrated amount of each tracer, in kg times concentration units [kg conc].
[in]gvThe ocean's vertical grid structure
csThe control structure returned by a previous call to register_ideal_age_tracer.
[out]namesthe names of the stocks calculated.
[out]unitsthe units of the stocks calculated.
[in]stock_indexthe coded index of a specific stock being sought.
Returns
The number of stocks calculated here.

Definition at line 373 of file ideal_age_example.F90.

373  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
374  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
375  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
376  real, dimension(:), intent(out) :: stocks !< the mass-weighted integrated amount of each
377  !! tracer, in kg times concentration units [kg conc].
378  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
379  type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
380  !! call to register_ideal_age_tracer.
381  character(len=*), dimension(:), intent(out) :: names !< the names of the stocks calculated.
382  character(len=*), dimension(:), intent(out) :: units !< the units of the stocks calculated.
383  integer, optional, intent(in) :: stock_index !< the coded index of a specific stock
384  !! being sought.
385  integer :: ideal_age_stock !< The number of stocks calculated here.
386 ! This function calculates the mass-weighted integral of all tracer stocks,
387 ! returning the number of stocks it has calculated. If the stock_index
388 ! is present, only the stock corresponding to that coded index is returned.
389 
390  integer :: i, j, k, is, ie, js, je, nz, m
391  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
392 
393  ideal_age_stock = 0
394  if (.not.associated(cs)) return
395  if (cs%ntr < 1) return
396 
397  if (present(stock_index)) then ; if (stock_index > 0) then
398  ! Check whether this stock is available from this routine.
399 
400  ! No stocks from this routine are being checked yet. Return 0.
401  return
402  endif ; endif
403 
404  do m=1,cs%ntr
405  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="ideal_age_stock")
406  units(m) = trim(units(m))//" kg"
407  stocks(m) = 0.0
408  do k=1,nz ; do j=js,je ; do i=is,ie
409  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
410  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
411  enddo ; enddo ; enddo
412  stocks(m) = gv%H_to_kg_m2 * stocks(m)
413  enddo
414  ideal_age_stock = cs%ntr
415 

◆ ideal_age_tracer_column_physics()

subroutine, public ideal_age_example::ideal_age_tracer_column_physics ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h_old,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h_new,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  ea,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  eb,
type(forcing), intent(in)  fluxes,
real, intent(in)  dt,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(ideal_age_tracer_cs), pointer  CS,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)

Applies diapycnal diffusion, aging and regeneration at the surface to the ideal age tracers.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]h_oldLayer thickness before entrainment [H ~> m or kg m-2].
[in]h_newLayer thickness after entrainment [H ~> m or kg m-2].
[in]eaan array to which the amount of fluid entrained
[in]eban array to which the amount of fluid entrained
[in]fluxesA structure containing pointers to thermodynamic and tracer forcing fields. Unused fields have NULL ptrs.
[in]dtThe amount of time covered by this call [s]
csThe control structure returned by a previous call to register_ideal_age_tracer.
[in]evap_cfl_limitLimit on the fraction of the water that can be fluxed out of the top layer in a timestep [nondim]
[in]minimum_forcing_depthThe smallest depth over which fluxes can be applied [m]

Definition at line 285 of file ideal_age_example.F90.

285  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
286  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
287  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
288  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
289  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
290  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
291  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
292  intent(in) :: ea !< an array to which the amount of fluid entrained
293  !! from the layer above during this call will be
294  !! added [H ~> m or kg m-2].
295  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
296  intent(in) :: eb !< an array to which the amount of fluid entrained
297  !! from the layer below during this call will be
298  !! added [H ~> m or kg m-2].
299  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
300  !! and tracer forcing fields. Unused fields have NULL ptrs.
301  real, intent(in) :: dt !< The amount of time covered by this call [s]
302  type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
303  !! call to register_ideal_age_tracer.
304  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
305  !! be fluxed out of the top layer in a timestep [nondim]
306  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
307  !! fluxes can be applied [m]
308 ! This subroutine applies diapycnal diffusion and any other column
309 ! tracer physics or chemistry to the tracers from this file.
310 ! This is a simple example of a set of advected passive tracers.
311 
312 ! The arguments to this subroutine are redundant in that
313 ! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
314  ! Local variables
315  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
316  real :: sfc_val ! The surface value for the tracers.
317  real :: Isecs_per_year ! The number of seconds in a year.
318  real :: year ! The time in years.
319  integer :: i, j, k, is, ie, js, je, nz, m
320  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
321 
322  if (.not.associated(cs)) return
323  if (cs%ntr < 1) return
324 
325  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
326  do m=1,cs%ntr
327  do k=1,nz ;do j=js,je ; do i=is,ie
328  h_work(i,j,k) = h_old(i,j,k)
329  enddo ; enddo ; enddo
330  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
331  evap_cfl_limit, minimum_forcing_depth)
332  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
333  enddo
334  else
335  do m=1,cs%ntr
336  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
337  enddo
338  endif
339 
340  isecs_per_year = 1.0 / (365.0*86400.0)
341  ! Set the surface value of tracer 1 to increase exponentially
342  ! with a 30 year time scale.
343  year = time_type_to_real(cs%Time) * isecs_per_year
344 
345  do m=1,cs%ntr
346  if (cs%sfc_growth_rate(m) == 0.0) then
347  sfc_val = cs%young_val(m)
348  else
349  sfc_val = cs%young_val(m) * &
350  exp((year-cs%tracer_start_year(m)) * cs%sfc_growth_rate(m))
351  endif
352  do k=1,cs%nkml ; do j=js,je ; do i=is,ie
353  if (g%mask2dT(i,j) > 0.5) then
354  cs%tr(i,j,k,m) = sfc_val
355  else
356  cs%tr(i,j,k,m) = cs%land_val(m)
357  endif
358  enddo ; enddo ; enddo
359  enddo
360  do m=1,cs%ntr ; if (cs%tracer_ages(m) .and. &
361  (year>=cs%tracer_start_year(m))) then
362 !$OMP parallel do default(none) shared(is,ie,js,je,CS,nz,G,dt,Isecs_per_year,m)
363  do k=cs%nkml+1,nz ; do j=js,je ; do i=is,ie
364  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt*isecs_per_year
365  enddo ; enddo ; enddo
366  endif ; enddo
367 

◆ ideal_age_tracer_surface_state()

subroutine, public ideal_age_example::ideal_age_tracer_surface_state ( type(surface), intent(inout)  state,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(ideal_age_tracer_cs), pointer  CS 
)

This subroutine extracts the surface fields from this tracer package that are to be shared with the atmosphere in coupled configurations. This particular tracer package does not report anything back to the coupler.

Parameters
[in]gThe ocean's grid structure.
[in,out]stateA structure containing fields that describe the surface state of the ocean.
[in]hLayer thickness [H ~> m or kg m-2].
csThe control structure returned by a previous call to register_ideal_age_tracer.

Definition at line 422 of file ideal_age_example.F90.

422  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
423  type(surface), intent(inout) :: state !< A structure containing fields that
424  !! describe the surface state of the ocean.
425  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
426  intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
427  type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
428  !! call to register_ideal_age_tracer.
429 
430  ! This particular tracer package does not report anything back to the coupler.
431  ! The code that is here is just a rough guide for packages that would.
432 
433  integer :: m, is, ie, js, je, isd, ied, jsd, jed
434  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
435  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
436 
437  if (.not.associated(cs)) return
438 
439  if (cs%coupled_tracers) then
440  do m=1,cs%ntr
441  ! This call loads the surface values into the appropriate array in the
442  ! coupler-type structure.
443  call coupler_type_set_data(cs%tr(:,:,1,m), cs%ind_tr(m), ind_csurf, &
444  state%tr_fields, idim=(/isd, is, ie, ied/), &
445  jdim=(/jsd, js, je, jed/) )
446  enddo
447  endif
448 

◆ initialize_ideal_age_tracer()

subroutine, public ideal_age_example::initialize_ideal_age_tracer ( logical, intent(in)  restart,
type(time_type), intent(in), target  day,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(ideal_age_tracer_cs), pointer  CS,
type(sponge_cs), pointer  sponge_CSp 
)

Sets the ideal age traces to their initial values and sets up the tracer output.

Parameters
[in]restart.true. if the fields have already been read from a restart file.
[in]dayTime of the start of the run.
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]diagA structure that is used to regulate diagnostic output.
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used.
csThe control structure returned by a previous call to register_ideal_age_tracer.
sponge_cspPointer to the control structure for the sponges.

Definition at line 197 of file ideal_age_example.F90.

197  logical, intent(in) :: restart !< .true. if the fields have already
198  !! been read from a restart file.
199  type(time_type), target, intent(in) :: day !< Time of the start of the run.
200  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
201  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
202  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
203  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
204  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
205  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
206  !! diagnostic output.
207  type(ocean_OBC_type), pointer :: OBC !< This open boundary condition type specifies
208  !! whether, where, and what open boundary
209  !! conditions are used.
210  type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
211  !! call to register_ideal_age_tracer.
212  type(sponge_CS), pointer :: sponge_CSp !< Pointer to the control structure for the sponges.
213 
214 ! This subroutine initializes the CS%ntr tracer fields in tr(:,:,:,:)
215 ! and it sets up the tracer output.
216 
217  ! Local variables
218  character(len=24) :: name ! A variable's name in a NetCDF file.
219  character(len=72) :: longname ! The long name of that variable.
220  character(len=48) :: units ! The dimensions of the variable.
221  character(len=48) :: flux_units ! The units for age tracer fluxes, either
222  ! years m3 s-1 or years kg s-1.
223  character(len=72) :: cmorname ! The CMOR name of that variable.
224  logical :: OK
225  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
226  integer :: IsdB, IedB, JsdB, JedB
227 
228  if (.not.associated(cs)) return
229  if (cs%ntr < 1) return
230  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
231  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
232  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
233 
234  cs%Time => day
235  cs%diag => diag
236  cs%nkml = max(gv%nkml,1)
237 
238  do m=1,cs%ntr
239  call query_vardesc(cs%tr_desc(m), name=name, &
240  caller="initialize_ideal_age_tracer")
241  if ((.not.restart) .or. (cs%tracers_may_reinit .and. .not. &
242  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
243 
244  if (len_trim(cs%IC_file) > 0) then
245  ! Read the tracer concentrations from a netcdf file.
246  if (.not.file_exists(cs%IC_file, g%Domain)) &
247  call mom_error(fatal, "initialize_ideal_age_tracer: "// &
248  "Unable to open "//cs%IC_file)
249 
250  if (cs%Z_IC_file) then
251  ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, name,&
252  g, us, -1e34, 0.0) ! CS%land_val(m))
253  if (.not.ok) then
254  ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, &
255  trim(name), g, us, -1e34, 0.0) ! CS%land_val(m))
256  if (.not.ok) call mom_error(fatal,"initialize_ideal_age_tracer: "//&
257  "Unable to read "//trim(name)//" from "//&
258  trim(cs%IC_file)//".")
259  endif
260  else
261  call mom_read_data(cs%IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
262  endif
263  else
264  do k=1,nz ; do j=js,je ; do i=is,ie
265  if (g%mask2dT(i,j) < 0.5) then
266  cs%tr(i,j,k,m) = cs%land_val(m)
267  else
268  cs%tr(i,j,k,m) = cs%IC_val(m)
269  endif
270  enddo ; enddo ; enddo
271  endif
272 
273  endif ! restart
274  enddo ! Tracer loop
275 
276  if (associated(obc)) then
277  ! Steal from updated DOME in the fullness of time.
278  endif
279 

◆ register_ideal_age_tracer()

logical function, public ideal_age_example::register_ideal_age_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(ideal_age_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

Register the ideal age tracer fields to be used with MOM.

Parameters
[in]hiA horizontal index type structure
[in]gvThe ocean's vertical grid structure
[in]param_fileA structure to parse for run-time parameters
csThe control structure returned by a previous call to register_ideal_age_tracer.
tr_regA pointer that is set to point to the control structure for the tracer advection and diffusion module
restart_csA pointer to the restart control structure

Definition at line 73 of file ideal_age_example.F90.

73  type(hor_index_type), intent(in) :: HI !< A horizontal index type structure
74  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
75  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
76  type(ideal_age_tracer_CS), pointer :: CS !< The control structure returned by a previous
77  !! call to register_ideal_age_tracer.
78  type(tracer_registry_type), pointer :: tr_Reg !< A pointer that is set to point to the control
79  !! structure for the tracer advection and
80  !! diffusion module
81  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure
82 
83 ! This include declares and sets the variable "version".
84 #include "version_variable.h"
85  character(len=40) :: mdl = "ideal_age_example" ! This module's name.
86  character(len=200) :: inputdir ! The directory where the input files are.
87  character(len=48) :: var_name ! The variable's name.
88  real, pointer :: tr_ptr(:,:,:) => null()
89  logical :: register_ideal_age_tracer
90  logical :: do_ideal_age, do_vintage, do_ideal_age_dated
91  integer :: isd, ied, jsd, jed, nz, m
92  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
93 
94  if (associated(cs)) then
95  call mom_error(warning, "register_ideal_age_tracer called with an "// &
96  "associated control structure.")
97  return
98  endif
99  allocate(cs)
100 
101  ! Read all relevant parameters and write them to the model log.
102  call log_version(param_file, mdl, version, "")
103  call get_param(param_file, mdl, "DO_IDEAL_AGE", do_ideal_age, &
104  "If true, use an ideal age tracer that is set to 0 age "//&
105  "in the mixed layer and ages at unit rate in the interior.", &
106  default=.true.)
107  call get_param(param_file, mdl, "DO_IDEAL_VINTAGE", do_vintage, &
108  "If true, use an ideal vintage tracer that is set to an "//&
109  "exponentially increasing value in the mixed layer and "//&
110  "is conserved thereafter.", default=.false.)
111  call get_param(param_file, mdl, "DO_IDEAL_AGE_DATED", do_ideal_age_dated, &
112  "If true, use an ideal age tracer that is everywhere 0 "//&
113  "before IDEAL_AGE_DATED_START_YEAR, but the behaves like "//&
114  "the standard ideal age tracer - i.e. is set to 0 age in "//&
115  "the mixed layer and ages at unit rate in the interior.", &
116  default=.false.)
117 
118 
119  call get_param(param_file, mdl, "AGE_IC_FILE", cs%IC_file, &
120  "The file in which the age-tracer initial values can be "//&
121  "found, or an empty string for internal initialization.", &
122  default=" ")
123  if ((len_trim(cs%IC_file) > 0) .and. (scan(cs%IC_file,'/') == 0)) then
124  ! Add the directory if CS%IC_file is not already a complete path.
125  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
126  cs%IC_file = trim(slasher(inputdir))//trim(cs%IC_file)
127  call log_param(param_file, mdl, "INPUTDIR/AGE_IC_FILE", cs%IC_file)
128  endif
129  call get_param(param_file, mdl, "AGE_IC_FILE_IS_Z", cs%Z_IC_file, &
130  "If true, AGE_IC_FILE is in depth space, not layer space", &
131  default=.false.)
132  call get_param(param_file, mdl, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
133  "If true, tracers may go through the initialization code "//&
134  "if they are not found in the restart files. Otherwise "//&
135  "it is a fatal error if the tracers are not found in the "//&
136  "restart files of a restarted run.", default=.false.)
137 
138  cs%ntr = 0
139  if (do_ideal_age) then
140  cs%ntr = cs%ntr + 1 ; m = cs%ntr
141  cs%tr_desc(m) = var_desc("age", "yr", "Ideal Age Tracer", cmor_field_name="agessc", caller=mdl)
142  cs%tracer_ages(m) = .true. ; cs%sfc_growth_rate(m) = 0.0
143  cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
144  endif
145 
146  if (do_vintage) then
147  cs%ntr = cs%ntr + 1 ; m = cs%ntr
148  cs%tr_desc(m) = var_desc("vintage", "yr", "Exponential Vintage Tracer", &
149  caller=mdl)
150  cs%tracer_ages(m) = .false. ; cs%sfc_growth_rate(m) = 1.0/30.0
151  cs%IC_val(m) = 0.0 ; cs%young_val(m) = 1e-20 ; cs%tracer_start_year(m) = 0.0
152  call get_param(param_file, mdl, "IDEAL_VINTAGE_START_YEAR", cs%tracer_start_year(m), &
153  "The date at which the ideal vintage tracer starts.", &
154  units="years", default=0.0)
155  endif
156 
157  if (do_ideal_age_dated) then
158  cs%ntr = cs%ntr + 1 ; m = cs%ntr
159  cs%tr_desc(m) = var_desc("age_dated","yr","Ideal Age Tracer with a Start Date",&
160  caller=mdl)
161  cs%tracer_ages(m) = .true. ; cs%sfc_growth_rate(m) = 0.0
162  cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
163  call get_param(param_file, mdl, "IDEAL_AGE_DATED_START_YEAR", cs%tracer_start_year(m), &
164  "The date at which the dated ideal age tracer starts.", &
165  units="years", default=0.0)
166  endif
167 
168  allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
169 
170  do m=1,cs%ntr
171  ! This is needed to force the compiler not to do a copy in the registration
172  ! calls. Curses on the designers and implementers of Fortran90.
173  tr_ptr => cs%tr(:,:,:,m)
174  call query_vardesc(cs%tr_desc(m), name=var_name, &
175  caller="register_ideal_age_tracer")
176  ! Register the tracer for horizontal advection, diffusion, and restarts.
177  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, tr_desc=cs%tr_desc(m), &
178  registry_diags=.true., restart_cs=restart_cs, &
179  mandatory=.not.cs%tracers_may_reinit)
180 
181  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
182  ! values to the coupler (if any). This is meta-code and its arguments will
183  ! currently (deliberately) give fatal errors if it is used.
184  if (cs%coupled_tracers) &
185  cs%ind_tr(m) = aof_set_coupler_flux(trim(var_name)//'_flux', &
186  flux_type=' ', implementation=' ', caller="register_ideal_age_tracer")
187  enddo
188 
189  cs%tr_Reg => tr_reg
190  cs%restart_CSp => restart_cs
191  register_ideal_age_tracer = .true.