MOM6
ideal_age_example.F90
Go to the documentation of this file.
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !***********************************************************************
21 
22 !********+*********+*********+*********+*********+*********+*********+**
23 !* *
24 !* By Robert Hallberg, 2002 *
25 !* *
26 !* This file contains an example of the code that is needed to set *
27 !* up and use a set (in this case two) of dynamically passive tracers *
28 !* for diagnostic purposes. The tracers here are an ideal age tracer *
29 !* that ages at a rate of 1/year once it is isolated from the surface,*
30 !* and a vintage tracer, whose surface concentration grows exponen- *
31 !* with time with a 30-year timescale (similar to CFCs). *
32 !* *
33 !* A single subroutine is called from within each file to register *
34 !* each of the tracers for reinitialization and advection and to *
35 !* register the subroutine that initializes the tracers and set up *
36 !* their output and the subroutine that does any tracer physics or *
37 !* chemistry along with diapycnal mixing (included here because some *
38 !* tracers may float or swim vertically or dye diapycnal processes). *
39 !* *
40 !* *
41 !* Macros written all in capital letters are defined in MOM_memory.h. *
42 !* *
43 !* A small fragment of the grid is shown below: *
44 !* *
45 !* j+1 x ^ x ^ x At x: q *
46 !* j+1 > o > o > At ^: v *
47 !* j x ^ x ^ x At >: u *
48 !* j > o > o > At o: h, tr *
49 !* j-1 x ^ x ^ x *
50 !* i-1 i i+1 At x & ^: *
51 !* i i+1 At > & o: *
52 !* *
53 !* The boundaries always run through q grid points (x). *
54 !* *
55 !********+*********+*********+*********+*********+*********+*********+**
56 
57 use mom_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
58 use mom_diag_mediator, only : diag_ctrl
60 use mom_error_handler, only : mom_error, fatal, warning
62 use mom_forcing_type, only : forcing
63 use mom_grid, only : ocean_grid_type
64 use mom_hor_index, only : hor_index_type
65 use mom_io, only : file_exists, read_data, slasher, vardesc, var_desc, query_vardesc
69 use mom_time_manager, only : time_type, get_time
74 use mom_variables, only : surface
76 
77 use coupler_util, only : set_coupler_values, ind_csurf
79 
80 implicit none ; private
81 
82 #include <MOM_memory.h>
83 
87 
88 ! NTR_MAX is the maximum number of tracers in this module.
89 integer, parameter :: ntr_max = 3
90 
91 type p3d
92  real, dimension(:,:,:), pointer :: p => null()
93 end type p3d
94 
95 type, public :: ideal_age_tracer_cs ; private
96  integer :: ntr ! The number of tracers that are actually used.
97  logical :: coupled_tracers = .false. ! These tracers are not offered to the
98  ! coupler.
99  integer :: nkml ! The number of layers in the mixed layer. The ideal
100  ! age tracers are reset in the top nkml layers.
101  character(len=200) :: ic_file ! The file in which the age-tracer initial values
102  ! can be found, or an empty string for internal initialization.
103  logical :: z_ic_file ! If true, the IC_file is in Z-space. The default is false.
104  type(time_type), pointer :: time ! A pointer to the ocean model's clock.
105  type(tracer_registry_type), pointer :: tr_reg => null()
106  real, pointer :: tr(:,:,:,:) => null() ! The array of tracers used in this
107  ! subroutine, in g m-3?
108  real, pointer :: tr_aux(:,:,:,:) => null() ! The masked tracer concentration
109  ! for output, in g m-3.
110  type(p3d), dimension(NTR_MAX) :: &
111  tr_adx, &! Tracer zonal advective fluxes in g m-3 m3 s-1.
112  tr_ady, &! Tracer meridional advective fluxes in g m-3 m3 s-1.
113  tr_dfx, &! Tracer zonal diffusive fluxes in g m-3 m3 s-1.
114  tr_dfy ! Tracer meridional diffusive fluxes in g m-3 m3 s-1.
115  real, dimension(NTR_MAX) :: &
116  ic_val = 0.0, & ! The (uniform) initial condition value.
117  young_val = 0.0, & ! The value assigned to tr at the surface.
118  land_val = -1.0, & ! The value of tr used where land is masked out.
119  sfc_growth_rate, & ! The exponential growth rate for the surface value,
120  ! in units of year-1.
121  tracer_start_year ! The year in which tracers start aging, or at which the
122  ! surface value equals young_val, in years.
123  logical :: mask_tracers ! If true, tracers are masked out in massless layers.
124  logical :: tracers_may_reinit ! If true, tracers may go through the
125  ! initialization code if they are not found in the
126  ! restart files.
127  logical :: tracer_ages(ntr_max)
128 
129  integer, dimension(NTR_MAX) :: &
130  ind_tr, & ! Indices returned by aof_set_coupler_flux if it is used and the
131  ! surface tracer concentrations are to be provided to the coupler.
132  id_tracer = -1, id_tr_adx = -1, id_tr_ady = -1, &
133  id_tr_dfx = -1, id_tr_dfy = -1
134 
135  type(diag_ctrl), pointer :: diag ! A structure that is used to regulate the
136  ! timing of diagnostic output.
137  type(mom_restart_cs), pointer :: restart_csp => null()
138 
139  type(vardesc) :: tr_desc(ntr_max)
140 end type ideal_age_tracer_cs
141 
142 contains
143 
144 function register_ideal_age_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
145  type(hor_index_type), intent(in) :: HI
146  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
147  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
148  type(ideal_age_tracer_cs), pointer :: CS
149  type(tracer_registry_type), pointer :: tr_Reg
150  type(mom_restart_cs), pointer :: restart_CS
151 ! This subroutine is used to register tracer fields and subroutines
152 ! to be used with MOM.
153 ! Arguments: HI - A horizontal index type structure.
154 ! (in) param_file - A structure indicating the open file to parse for
155 ! model parameter values.
156 ! (in/out) CS - A pointer that is set to point to the control structure
157 ! for this module
158 ! (in/out) tr_Reg - A pointer that is set to point to the control structure
159 ! for the tracer advection and diffusion module.
160 ! (in) restart_CS - A pointer to the restart control structure.
161 
162 ! This include declares and sets the variable "version".
163 #include "version_variable.h"
164  character(len=40) :: mdl = "ideal_age_example" ! This module's name.
165  character(len=200) :: inputdir ! The directory where the input files are.
166  character(len=48) :: var_name ! The variable's name.
167  real, pointer :: tr_ptr(:,:,:) => null()
168  logical :: register_ideal_age_tracer
169  logical :: do_ideal_age, do_vintage, do_ideal_age_dated
170  integer :: isd, ied, jsd, jed, nz, m
171  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
172 
173  if (associated(cs)) then
174  call mom_error(warning, "register_ideal_age_tracer called with an "// &
175  "associated control structure.")
176  return
177  endif
178  allocate(cs)
179 
180  ! Read all relevant parameters and write them to the model log.
181  call log_version(param_file, mdl, version, "")
182  call get_param(param_file, mdl, "DO_IDEAL_AGE", do_ideal_age, &
183  "If true, use an ideal age tracer that is set to 0 age \n"//&
184  "in the mixed layer and ages at unit rate in the interior.", &
185  default=.true.)
186  call get_param(param_file, mdl, "DO_IDEAL_VINTAGE", do_vintage, &
187  "If true, use an ideal vintage tracer that is set to an \n"//&
188  "exponentially increasing value in the mixed layer and \n"//&
189  "is conserved thereafter.", default=.false.)
190  call get_param(param_file, mdl, "DO_IDEAL_AGE_DATED", do_ideal_age_dated, &
191  "If true, use an ideal age tracer that is everywhere 0 \n"//&
192  "before IDEAL_AGE_DATED_START_YEAR, but the behaves like \n"//&
193  "the standard ideal age tracer - i.e. is set to 0 age in \n"//&
194  "the mixed layer and ages at unit rate in the interior.", &
195  default=.false.)
196 
197 
198  call get_param(param_file, mdl, "AGE_IC_FILE", cs%IC_file, &
199  "The file in which the age-tracer initial values can be \n"//&
200  "found, or an empty string for internal initialization.", &
201  default=" ")
202  if ((len_trim(cs%IC_file) > 0) .and. (scan(cs%IC_file,'/') == 0)) then
203  ! Add the directory if CS%IC_file is not already a complete path.
204  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
205  cs%IC_file = trim(slasher(inputdir))//trim(cs%IC_file)
206  call log_param(param_file, mdl, "INPUTDIR/AGE_IC_FILE", cs%IC_file)
207  endif
208  call get_param(param_file, mdl, "AGE_IC_FILE_IS_Z", cs%Z_IC_file, &
209  "If true, AGE_IC_FILE is in depth space, not layer space", &
210  default=.false.)
211  call get_param(param_file, mdl, "MASK_MASSLESS_TRACERS", cs%mask_tracers, &
212  "If true, the tracers are masked out in massless layer. \n"//&
213  "This can be a problem with time-averages.", default=.false.)
214  call get_param(param_file, mdl, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
215  "If true, tracers may go through the initialization code \n"//&
216  "if they are not found in the restart files. Otherwise \n"//&
217  "it is a fatal error if the tracers are not found in the \n"//&
218  "restart files of a restarted run.", default=.false.)
219 
220  cs%ntr = 0
221  if (do_ideal_age) then
222  cs%ntr = cs%ntr + 1 ; m = cs%ntr
223  cs%tr_desc(m) = var_desc("age", "years", "Ideal Age Tracer", cmor_field_name="agessc", caller=mdl)
224  cs%tracer_ages(m) = .true. ; cs%sfc_growth_rate(m) = 0.0
225  cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
226  endif
227 
228  if (do_vintage) then
229  cs%ntr = cs%ntr + 1 ; m = cs%ntr
230  cs%tr_desc(m) = var_desc("vintage", "years", "Exponential Vintage Tracer", &
231  caller=mdl)
232  cs%tracer_ages(m) = .false. ; cs%sfc_growth_rate(m) = 1.0/30.0
233  cs%IC_val(m) = 0.0 ; cs%young_val(m) = 1e-20 ; cs%tracer_start_year(m) = 0.0
234  call get_param(param_file, mdl, "IDEAL_VINTAGE_START_YEAR", cs%tracer_start_year(m), &
235  "The date at which the ideal vintage tracer starts.", &
236  units="years", default=0.0)
237  endif
238 
239  if (do_ideal_age_dated) then
240  cs%ntr = cs%ntr + 1 ; m = cs%ntr
241  cs%tr_desc(m) = var_desc("age_dated","years","Ideal Age Tracer with a Start Date",&
242  caller=mdl)
243  cs%tracer_ages(m) = .true. ; cs%sfc_growth_rate(m) = 0.0
244  cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
245  call get_param(param_file, mdl, "IDEAL_AGE_DATED_START_YEAR", cs%tracer_start_year(m), &
246  "The date at which the dated ideal age tracer starts.", &
247  units="years", default=0.0)
248  endif
249 
250  allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
251  if (cs%mask_tracers) then
252  allocate(cs%tr_aux(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr_aux(:,:,:,:) = 0.0
253  endif
254 
255  do m=1,cs%ntr
256  ! This is needed to force the compiler not to do a copy in the registration
257  ! calls. Curses on the designers and implementers of Fortran90.
258  tr_ptr => cs%tr(:,:,:,m)
259  call query_vardesc(cs%tr_desc(m), name=var_name, &
260  caller="register_ideal_age_tracer")
261  ! Register the tracer for the restart file.
262  call register_restart_field(tr_ptr, cs%tr_desc(m), &
263  .not.cs%tracers_may_reinit, restart_cs)
264  ! Register the tracer for horizontal advection & diffusion.
265  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
266  tr_desc_ptr=cs%tr_desc(m))
267 
268  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
269  ! values to the coupler (if any). This is meta-code and its arguments will
270  ! currently (deliberately) give fatal errors if it is used.
271  if (cs%coupled_tracers) &
272  cs%ind_tr(m) = aof_set_coupler_flux(trim(var_name)//'_flux', &
273  flux_type=' ', implementation=' ', caller="register_ideal_age_tracer")
274  enddo
275 
276  cs%tr_Reg => tr_reg
277  cs%restart_CSp => restart_cs
278  register_ideal_age_tracer = .true.
279 end function register_ideal_age_tracer
280 
281 subroutine initialize_ideal_age_tracer(restart, day, G, GV, h, diag, OBC, CS, &
282  sponge_CSp, diag_to_Z_CSp)
283  logical, intent(in) :: restart
284  type(time_type), target, intent(in) :: day
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)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
288  type(diag_ctrl), target, intent(in) :: diag
289  type(ocean_obc_type), pointer :: OBC
290  type(ideal_age_tracer_cs), pointer :: CS
291  type(sponge_cs), pointer :: sponge_CSp
292  type(diag_to_z_cs), pointer :: diag_to_Z_CSp
293 ! This subroutine initializes the CS%ntr tracer fields in tr(:,:,:,:)
294 ! and it sets up the tracer output.
295 
296 ! Arguments: restart - .true. if the fields have already been read from
297 ! a restart file.
298 ! (in) day - Time of the start of the run.
299 ! (in) G - The ocean's grid structure.
300 ! (in) GV - The ocean's vertical grid structure.
301 ! (in) h - Layer thickness, in m or kg m-2.
302 ! (in) diag - A structure that is used to regulate diagnostic output.
303 ! (in) OBC - This open boundary condition type specifies whether, where,
304 ! and what open boundary conditions are used.
305 ! (in/out) CS - The control structure returned by a previous call to
306 ! register_ideal_age_tracer.
307 ! (in/out) sponge_CSp - A pointer to the control structure for the sponges, if
308 ! they are in use. Otherwise this may be unassociated.
309 ! (in/out) diag_to_Z_Csp - A pointer to the control structure for diagnostics
310 ! in depth space.
311  character(len=24) :: name ! A variable's name in a NetCDF file.
312  character(len=72) :: longname ! The long name of that variable.
313  character(len=48) :: units ! The dimensions of the variable.
314  character(len=48) :: flux_units ! The units for age tracer fluxes, either
315  ! years m3 s-1 or years kg s-1.
316  character(len=72) :: cmorname ! The CMOR name of that variable.
317  logical :: OK
318  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
319  integer :: IsdB, IedB, JsdB, JedB
320 
321  if (.not.associated(cs)) return
322  if (cs%ntr < 1) return
323  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
324  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
325  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
326 
327  cs%Time => day
328  cs%diag => diag
329  cs%nkml = max(gv%nkml,1)
330 
331  do m=1,cs%ntr
332  call query_vardesc(cs%tr_desc(m), name=name, &
333  caller="initialize_ideal_age_tracer")
334  if ((.not.restart) .or. (cs%tracers_may_reinit .and. .not. &
335  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
336 
337  if (len_trim(cs%IC_file) > 0) then
338  ! Read the tracer concentrations from a netcdf file.
339  if (.not.file_exists(cs%IC_file, g%Domain)) &
340  call mom_error(fatal, "initialize_ideal_age_tracer: "// &
341  "Unable to open "//cs%IC_file)
342 
343  if (cs%Z_IC_file) then
344  ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, name,&
345  g, -1e34, 0.0) ! CS%land_val(m))
346  if (.not.ok) then
347  ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, &
348  trim(name), g, -1e34, 0.0) ! CS%land_val(m))
349  if (.not.ok) call mom_error(fatal,"initialize_ideal_age_tracer: "//&
350  "Unable to read "//trim(name)//" from "//&
351  trim(cs%IC_file)//".")
352  endif
353  else
354  call read_data(cs%IC_file, trim(name), cs%tr(:,:,:,m), &
355  domain=g%Domain%mpp_domain)
356  endif
357  else
358  do k=1,nz ; do j=js,je ; do i=is,ie
359  if (g%mask2dT(i,j) < 0.5) then
360  cs%tr(i,j,k,m) = cs%land_val(m)
361  else
362  cs%tr(i,j,k,m) = cs%IC_val(m)
363  endif
364  enddo ; enddo ; enddo
365  endif
366 
367  endif ! restart
368  enddo ! Tracer loop
369 
370  if (associated(obc)) then
371  ! All tracers but the first have 0 concentration in their inflows. As this
372  ! is the default value, the following calls are unnecessary.
373  ! do m=1,CS%ntr
374  ! call add_tracer_OBC_values(trim(CS%tr_desc(m)%name), CS%tr_Reg, 0.0)
375  ! enddo
376  endif
377 
378  ! This needs to be changed if the units of tracer are changed above.
379  if (gv%Boussinesq) then ; flux_units = "years m3 s-1"
380  else ; flux_units = "years kg s-1" ; endif
381 
382  do m=1,cs%ntr
383  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
384  cmor_field_name=cmorname, caller="initialize_ideal_age_tracer")
385  if (len_trim(cmorname)==0) then
386  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
387  day, trim(longname) , trim(units))
388  else
389  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
390  day, trim(longname) , trim(units), cmor_field_name=cmorname)
391  endif
392  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
393  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
394  trim(flux_units))
395  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
396  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
397  trim(flux_units))
398  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
399  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
400  trim(flux_units))
401  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
402  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
403  trim(flux_units))
404  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
405  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
406  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
407  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
408 
409 ! Register the tracer for horizontal advection & diffusion.
410  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
411  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
412  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
413  cs%tr_ady(m)%p,cs%tr_dfx(m)%p,cs%tr_dfy(m)%p)
414 
415  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
416  day, g, diag_to_z_csp)
417  enddo
418 
419 end subroutine initialize_ideal_age_tracer
420 
421 subroutine ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, &
422  evap_CFL_limit, minimum_forcing_depth)
423  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
424  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
425  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
426  type(forcing), intent(in) :: fluxes
427  real, intent(in) :: dt !< The amount of time covered by this call, in s
428  type(ideal_age_tracer_cs), pointer :: CS
429  real, optional,intent(in) :: evap_CFL_limit
430  real, optional,intent(in) :: minimum_forcing_depth
431 ! This subroutine applies diapycnal diffusion and any other column
432 ! tracer physics or chemistry to the tracers from this file.
433 ! This is a simple example of a set of advected passive tracers.
434 
435 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
436 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
437 ! (in) ea - an array to which the amount of fluid entrained
438 ! from the layer above during this call will be
439 ! added, in m or kg m-2.
440 ! (in) eb - an array to which the amount of fluid entrained
441 ! from the layer below during this call will be
442 ! added, in m or kg m-2.
443 ! (in) fluxes - A structure containing pointers to any possible
444 ! forcing fields. Unused fields have NULL ptrs.
445 ! (in) dt - The amount of time covered by this call, in s.
446 ! (in) G - The ocean's grid structure.
447 ! (in) GV - The ocean's vertical grid structure.
448 ! (in) CS - The control structure returned by a previous call to
449 ! register_ideal_age_tracer.
450 !
451 ! The arguments to this subroutine are redundant in that
452 ! h_new[k] = h_old[k] + ea[k] - eb[k-1] + eb[k] - ea[k+1]
453  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
454  real :: sfc_val ! The surface value for the tracers.
455  real :: Isecs_per_year ! The number of seconds in a year.
456  real :: year ! The time in years.
457  integer :: secs, days ! Integer components of the time type.
458  integer :: i, j, k, is, ie, js, je, nz, m
459  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
460 
461  if (.not.associated(cs)) return
462  if (cs%ntr < 1) return
463 
464  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
465  do m=1,cs%ntr
466  do k=1,nz ;do j=js,je ; do i=is,ie
467  h_work(i,j,k) = h_old(i,j,k)
468  enddo ; enddo ; enddo;
469  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
470  evap_cfl_limit, minimum_forcing_depth)
471  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
472  enddo
473  else
474  do m=1,cs%ntr
475  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
476  enddo
477  endif
478 
479  isecs_per_year = 1.0 / (365.0*86400.0)
480  ! Set the surface value of tracer 1 to increase exponentially
481  ! with a 30 year time scale.
482  call get_time(cs%Time, secs, days)
483  year = (86400.0*days + real(secs)) * isecs_per_year
484 
485  do m=1,cs%ntr
486  if (cs%sfc_growth_rate(m) == 0.0) then
487  sfc_val = cs%young_val(m)
488  else
489  sfc_val = cs%young_val(m) * &
490  exp((year-cs%tracer_start_year(m)) * cs%sfc_growth_rate(m))
491  endif
492  do k=1,cs%nkml ; do j=js,je ; do i=is,ie
493  if (g%mask2dT(i,j) > 0.5) then
494  cs%tr(i,j,k,m) = sfc_val
495  else
496  cs%tr(i,j,k,m) = cs%land_val(m)
497  endif
498  enddo ; enddo ; enddo
499  enddo
500  do m=1,cs%ntr ; if (cs%tracer_ages(m) .and. &
501  (year>=cs%tracer_start_year(m))) then
502 !$OMP parallel do default(none) shared(is,ie,js,je,CS,nz,G,dt,Isecs_per_year,m)
503  do k=cs%nkml+1,nz ; do j=js,je ; do i=is,ie
504  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt*isecs_per_year
505  enddo ; enddo ; enddo
506  endif ; enddo
507 
508  if (cs%mask_tracers) then
509  do m=1,cs%ntr ; if (cs%id_tracer(m) > 0) then
510  do k=1,nz ; do j=js,je ; do i=is,ie
511  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
512  cs%tr_aux(i,j,k,m) = cs%land_val(m)
513  else
514  cs%tr_aux(i,j,k,m) = cs%tr(i,j,k,m)
515  endif
516  enddo ; enddo ; enddo
517  endif ; enddo
518  endif
519 
520  do m=1,cs%ntr
521  if (cs%mask_tracers) then
522  if (cs%id_tracer(m)>0) &
523  call post_data(cs%id_tracer(m),cs%tr_aux(:,:,:,m),cs%diag)
524  else
525  if (cs%id_tracer(m)>0) &
526  call post_data(cs%id_tracer(m),cs%tr(:,:,:,m),cs%diag)
527  endif
528  if (cs%id_tr_adx(m)>0) &
529  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
530  if (cs%id_tr_ady(m)>0) &
531  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
532  if (cs%id_tr_dfx(m)>0) &
533  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
534  if (cs%id_tr_dfy(m)>0) &
535  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
536  enddo
537 
538 end subroutine ideal_age_tracer_column_physics
539 
540 function ideal_age_stock(h, stocks, G, GV, CS, names, units, stock_index)
541  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
542  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
543  real, dimension(:), intent(out) :: stocks
544  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
545  type(ideal_age_tracer_cs), pointer :: CS
546  character(len=*), dimension(:), intent(out) :: names
547  character(len=*), dimension(:), intent(out) :: units
548  integer, optional, intent(in) :: stock_index
549  integer :: ideal_age_stock
550 ! This function calculates the mass-weighted integral of all tracer stocks,
551 ! returning the number of stocks it has calculated. If the stock_index
552 ! is present, only the stock corresponding to that coded index is returned.
553 
554 ! Arguments: h - Layer thickness, in m or kg m-2.
555 ! (out) stocks - the mass-weighted integrated amount of each tracer,
556 ! in kg times concentration units.
557 ! (in) G - The ocean's grid structure.
558 ! (in) GV - The ocean's vertical grid structure.
559 ! (in) CS - The control structure returned by a previous call to
560 ! register_ideal_age_tracer.
561 ! (out) names - the names of the stocks calculated.
562 ! (out) units - the units of the stocks calculated.
563 ! (in,opt) stock_index - the coded index of a specific stock being sought.
564 ! Return value: the number of stocks calculated here.
565 
566  integer :: i, j, k, is, ie, js, je, nz, m
567  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
568 
569  ideal_age_stock = 0
570  if (.not.associated(cs)) return
571  if (cs%ntr < 1) return
572 
573  if (present(stock_index)) then ; if (stock_index > 0) then
574  ! Check whether this stock is available from this routine.
575 
576  ! No stocks from this routine are being checked yet. Return 0.
577  return
578  endif ; endif
579 
580  do m=1,cs%ntr
581  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="ideal_age_stock")
582  units(m) = trim(units(m))//" kg"
583  stocks(m) = 0.0
584  do k=1,nz ; do j=js,je ; do i=is,ie
585  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
586  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
587  enddo ; enddo ; enddo
588  stocks(m) = gv%H_to_kg_m2 * stocks(m)
589  enddo
590  ideal_age_stock = cs%ntr
591 
592 end function ideal_age_stock
593 
594 subroutine ideal_age_tracer_surface_state(state, h, G, CS)
595  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
596  type(surface), intent(inout) :: state
597  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses, in H (usually m or kg m-2)
598  type(ideal_age_tracer_cs), pointer :: CS
599 ! This particular tracer package does not report anything back to the coupler.
600 ! The code that is here is just a rough guide for packages that would.
601 ! Arguments: state - A structure containing fields that describe the
602 ! surface state of the ocean.
603 ! (in) h - Layer thickness, in m or kg m-2.
604 ! (in) G - The ocean's grid structure.
605 ! (in) CS - The control structure returned by a previous call to
606 ! register_ideal_age_tracer.
607  integer :: i, j, m, is, ie, js, je
608  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
609 
610  if (.not.associated(cs)) return
611 
612  if (cs%coupled_tracers) then
613  do m=1,cs%ntr
614  ! This call loads the surface vlues into the appropriate array in the
615  ! coupler-type structure.
616  call set_coupler_values(cs%tr(:,:,1,m), state%tr_fields, cs%ind_tr(m), &
617  ind_csurf, is, ie, js, je)
618  enddo
619  endif
620 
621 end subroutine ideal_age_tracer_surface_state
622 
623 subroutine ideal_age_example_end(CS)
624  type(ideal_age_tracer_cs), pointer :: CS
625  integer :: m
626 
627  if (associated(cs)) then
628  if (associated(cs%tr)) deallocate(cs%tr)
629  if (associated(cs%tr_aux)) deallocate(cs%tr_aux)
630  do m=1,cs%ntr
631  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
632  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
633  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
634  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
635  enddo
636 
637  deallocate(cs)
638  endif
639 end subroutine ideal_age_example_end
640 
641 end module ideal_age_example
The following structure contains pointers to various fields which may be used describe the surface st...
subroutine, public ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
type(vardesc) function, public var_desc(name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
Returns a vardesc type whose elements have been filled with the provided fields. The argument name is...
Definition: MOM_io.F90:585
subroutine, public ideal_age_tracer_surface_state(state, h, G, CS)
This module implements boundary forcing for MOM6.
integer function, public aof_set_coupler_flux(name, flux_type, implementation, atm_tr_index, param, flag, ice_restart_file, ocean_restart_file, units, caller)
subroutine, public initialize_ideal_age_tracer(restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, diag_to_Z_CSp)
logical function, public tracer_z_init(tr, h, filename, tr_name, G, missing_val, land_val)
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
subroutine, public ideal_age_example_end(CS)
Provides the ocean grid type.
Definition: MOM_grid.F90:2
This module contains I/O framework code.
Definition: MOM_io.F90:2
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Container for horizontal index ranges for data, computational and global domains. ...
subroutine, public register_tracer(tr1, tr_desc, param_file, HI, GV, Reg, tr_desc_ptr, ad_x, ad_y, df_x, df_y, OBC_inflow, OBC_in_u, OBC_in_v, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine registers a tracer to be advected and laterally diffused.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
subroutine, public applytracerboundaryfluxesinout(G, GV, Tr, dt, fluxes, h, evap_CFL_limit, minimum_forcing_depth, in_flux_optional, out_flux_optional, update_h_opt)
This routine is modeled after applyBoundaryFluxesInOut in MOM_diabatic_aux.F90 NOTE: Please note that...
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
Type to carry basic tracer information.
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
subroutine, public register_z_tracer(tr_ptr, name, long_name, units, Time, G, CS, standard_name, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name)
This subroutine registers a tracer to be output in depth space.
subroutine, public add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
This subroutine adds open boundary condition concentrations for a tracer that has previously been reg...
subroutine, public add_tracer_diagnostics(name, Reg, ad_x, ad_y, df_x, df_y, ad_2d_x, ad_2d_y, df_2d_x, df_2d_y, advection_xy)
This subroutine adds diagnostic arrays for a tracer that has previously been registered by a call to ...
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
subroutine, public set_coupler_values(array_in, BC_struc, BC_index, BC_element, is, ie, js, je, conversion)
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
subroutine, public tracer_vertdiff(h_old, ea, eb, dt, tr, G, GV, sfc_flux, btm_flux, btm_reservoir, sink_rate, convert_flux_in)
This subroutine solves a tridiagonal equation for the final tracer concentrations after the dual-entr...
logical function, public register_ideal_age_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
integer, parameter ntr_max
subroutine, public query_vardesc(vd, name, units, longname, hor_grid, z_grid, t_grid, cmor_field_name, cmor_units, conversion, caller)
This routine queries vardesc.
Definition: MOM_io.F90:664
Controls where open boundary conditions are applied.
integer function, public register_diag_field(module_name, field_name, axes, init_time, long_name, units, missing_value, range, mask_variant, standard_name, verbose, do_not_log, err_msg, interp_method, tile_count, cmor_field_name, cmor_long_name, cmor_units, cmor_standard_name, cell_methods, x_cell_method, y_cell_method, v_cell_method, conversion, v_extensive)
Returns the "diag_mediator" handle for a group (native, CMOR, z-coord, ...) of diagnostics derived fr...
subroutine, public mom_error(level, message, all_print)
integer function, public ideal_age_stock(h, stocks, G, GV, CS, names, units, stock_index)