ideal_age_example.F90
1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> A tracer package of ideal age tracers
6module ideal_age_example
7
9use mom_coupler_types, only : set_coupler_type_data, atmos_ocn_coupler_flux
10use mom_diag_mediator, only : diag_ctrl
11use mom_error_handler, only : mom_error, fatal, warning, note
12use mom_file_parser, only : get_param, log_param, log_version, param_file_type
13use mom_forcing_type, only : forcing
14use mom_grid, only : ocean_grid_type
15use mom_hor_index, only : hor_index_type
16use mom_io, only : file_exists, mom_read_data, slasher, vardesc, var_desc, query_vardesc
17use mom_open_boundary, only : ocean_obc_type
18use mom_restart, only : query_initialized, set_initialized, mom_restart_cs
19use mom_spatial_means, only : global_mass_int_efp
20use mom_sponge, only : set_up_sponge_field, sponge_cs
21use mom_time_manager, only : time_type, time_to_real
22use mom_tracer_registry, only : register_tracer, tracer_registry_type
23use mom_tracer_diabatic, only : tracer_vertdiff, applytracerboundaryfluxesinout
24use mom_tracer_z_init, only : tracer_z_init
25use mom_unit_scaling, only : unit_scale_type
26use mom_variables, only : surface
27use mom_verticalgrid, only : verticalgrid_type
28
29implicit none ; private
30
31#include <MOM_memory.h>
32
35public ideal_age_stock, ideal_age_example_end
36public count_bl_layers
37
38integer, parameter :: ntr_max = 4 !< the maximum number of tracers in this module.
39
40!> The control structure for the ideal_age_tracer package
41type, public :: ideal_age_tracer_cs ; private
42 integer :: ntr !< The number of tracers that are actually used.
43 logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
44 integer :: nkbl !< The number of layers in the boundary layer. The ideal
45 !1 age tracers are reset in the top nkbl layers.
46 character(len=200) :: ic_file !< The file in which the age-tracer initial values
47 !! can be found, or an empty string for internal initialization.
48 logical :: z_ic_file !< If true, the IC_file is in Z-space. The default is false.
49 type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
50 type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the tracer registry
51 real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this package [years] or other units
52 real, dimension(NTR_MAX) :: ic_val = 0.0 !< The (uniform) initial condition value [years] or other units
53 real, dimension(NTR_MAX) :: young_val = 0.0 !< The value assigned to tr at the surface [years] or other units
54 real, dimension(NTR_MAX) :: land_val = -1.0 !< The value of tr used where land is masked out [years] or other units
55 real, dimension(NTR_MAX) :: growth_rate !< The exponential growth rate for the young value [year-1]
56 real, dimension(NTR_MAX) :: tracer_start_year !< The year in which tracers start aging, or at which the
57 !! surface value equals young_val [years].
58 logical :: use_real_bl_depth !< If true, uses the BL scheme to determine the number of
59 !! layers above the BL depth instead of the fixed nkbl value.
60 integer :: bl_residence_num !< The tracer number assigned to the BL residence tracer in this module
61 logical :: tracers_may_reinit !< If true, these tracers be set up via the initialization code if
62 !! they are not found in the restart files.
63 logical :: tracer_ages(ntr_max) !< Indicates whether each tracer ages.
64
65 integer, dimension(NTR_MAX) :: ind_tr !< Indices returned by atmos_ocn_coupler_flux if it is used and the
66 !! surface tracer concentrations are to be provided to the coupler.
67
68 type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
69 !! regulate the timing of diagnostic output.
70 type(mom_restart_cs), pointer :: restart_csp => null() !< A pointer to the restart controls structure
71
72 type(vardesc) :: tr_desc(ntr_max) !< Descriptions and metadata for the tracers
73
74end type ideal_age_tracer_cs
75
76contains
77
78!> Register the ideal age tracer fields to be used with MOM.
79function register_ideal_age_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
80 type(hor_index_type), intent(in) :: hi !< A horizontal index type structure
81 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
82 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
83 type(ideal_age_tracer_cs), pointer :: cs !< The control structure returned by a previous
84 !! call to register_ideal_age_tracer.
85 type(tracer_registry_type), pointer :: tr_reg !< A pointer that is set to point to the control
86 !! structure for the tracer advection and
87 !! diffusion module
88 type(mom_restart_cs), target, intent(inout) :: restart_cs !< MOM restart control struct
89
90 ! This include declares and sets the variable "version".
91# include "version_variable.h"
92 character(len=40) :: mdl = "ideal_age_example" ! This module's name.
93 character(len=200) :: inputdir ! The directory where the input files are.
94 character(len=48) :: var_name ! The variable's name.
95 real, pointer :: tr_ptr(:,:,:) => null() ! The tracer concentration [years]
96 logical :: register_ideal_age_tracer
97 logical :: do_ideal_age, do_vintage, do_ideal_age_dated, do_bl_residence
98 integer :: isd, ied, jsd, jed, nz, m
99 isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
100
101 if (associated(cs)) then
102 call mom_error(fatal, "register_ideal_age_tracer called with an "// &
103 "associated control structure.")
104 endif
105 allocate(cs)
106
107 ! Read all relevant parameters and write them to the model log.
108 call log_version(param_file, mdl, version, "")
109 call get_param(param_file, mdl, "DO_IDEAL_AGE", do_ideal_age, &
110 "If true, use an ideal age tracer that is set to 0 age "//&
111 "in the boundary layer and ages at unit rate in the interior.", &
112 default=.true.)
113 call get_param(param_file, mdl, "DO_IDEAL_VINTAGE", do_vintage, &
114 "If true, use an ideal vintage tracer that is set to an "//&
115 "exponentially increasing value in the boundary layer and "//&
116 "is conserved thereafter.", default=.false.)
117 call get_param(param_file, mdl, "DO_IDEAL_AGE_DATED", do_ideal_age_dated, &
118 "If true, use an ideal age tracer that is everywhere 0 "//&
119 "before IDEAL_AGE_DATED_START_YEAR, but the behaves like "//&
120 "the standard ideal age tracer - i.e. is set to 0 age in "//&
121 "the boundary layer and ages at unit rate in the interior.", &
122 default=.false.)
123 call get_param(param_file, mdl, "DO_BL_RESIDENCE", do_bl_residence, &
124 "If true, use a residence tracer that is set to 0 age "//&
125 "in the interior and ages at unit rate in the boundary layer.", &
126 default=.false.)
127 call get_param(param_file, mdl, "USE_REAL_BL_DEPTH", cs%use_real_BL_depth, &
128 "If true, the ideal age tracers will use the boundary layer "//&
129 "depth diagnosed from the BL or bulkmixedlayer scheme.", &
130 default=.false.)
131 call get_param(param_file, mdl, "AGE_IC_FILE", cs%IC_file, &
132 "The file in which the age-tracer initial values can be "//&
133 "found, or an empty string for internal initialization.", &
134 default=" ")
135 if ((len_trim(cs%IC_file) > 0) .and. (scan(cs%IC_file,'/') == 0)) then
136 ! Add the directory if CS%IC_file is not already a complete path.
137 call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
138 cs%IC_file = trim(slasher(inputdir))//trim(cs%IC_file)
139 call log_param(param_file, mdl, "INPUTDIR/AGE_IC_FILE", cs%IC_file)
140 endif
141 call get_param(param_file, mdl, "AGE_IC_FILE_IS_Z", cs%Z_IC_file, &
142 "If true, AGE_IC_FILE is in depth space, not layer space", &
143 default=.false.)
144 call get_param(param_file, mdl, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
145 "If true, tracers may go through the initialization code "//&
146 "if they are not found in the restart files. Otherwise "//&
147 "it is a fatal error if the tracers are not found in the "//&
148 "restart files of a restarted run.", default=.false.)
149
150 cs%ntr = 0
151 if (do_ideal_age) then
152 cs%ntr = cs%ntr + 1 ; m = cs%ntr
153 cs%tr_desc(m) = var_desc("age", "yr", "Ideal Age Tracer", cmor_field_name="agessc", caller=mdl)
154 cs%tracer_ages(m) = .true. ; cs%growth_rate(m) = 0.0
155 cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
156 endif
157
158 if (do_vintage) then
159 cs%ntr = cs%ntr + 1 ; m = cs%ntr
160 cs%tr_desc(m) = var_desc("vintage", "yr", "Exponential Vintage Tracer", &
161 caller=mdl)
162 cs%tracer_ages(m) = .false. ; cs%growth_rate(m) = 1.0/30.0
163 cs%IC_val(m) = 0.0 ; cs%young_val(m) = 1e-20 ; cs%tracer_start_year(m) = 0.0
164 call get_param(param_file, mdl, "IDEAL_VINTAGE_START_YEAR", cs%tracer_start_year(m), &
165 "The date at which the ideal vintage tracer starts.", &
166 units="years", default=0.0)
167 endif
168
169 if (do_ideal_age_dated) then
170 cs%ntr = cs%ntr + 1 ; m = cs%ntr
171 cs%tr_desc(m) = var_desc("age_dated","yr","Ideal Age Tracer with a Start Date",&
172 caller=mdl)
173 cs%tracer_ages(m) = .true. ; cs%growth_rate(m) = 0.0
174 cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
175 call get_param(param_file, mdl, "IDEAL_AGE_DATED_START_YEAR", cs%tracer_start_year(m), &
176 "The date at which the dated ideal age tracer starts.", &
177 units="years", default=0.0)
178 endif
179
180 cs%BL_residence_num = 0
181 if (do_bl_residence) then
182 cs%ntr = cs%ntr + 1 ; m = cs%ntr ; cs%BL_residence_num = cs%ntr
183 cs%tr_desc(m) = var_desc("BL_age", "yr", "BL Residence Time Tracer", caller=mdl)
184 cs%tracer_ages(m) = .true. ; cs%growth_rate(m) = 0.0
185 cs%IC_val(m) = 0.0 ; cs%young_val(m) = 0.0 ; cs%tracer_start_year(m) = 0.0
186 endif
187
188 allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr), source=0.0)
189
190 do m=1,cs%ntr
191 ! This is needed to force the compiler not to do a copy in the registration
192 ! calls. Curses on the designers and implementers of Fortran90.
193 tr_ptr => cs%tr(:,:,:,m)
194 call query_vardesc(cs%tr_desc(m), name=var_name, &
195 caller="register_ideal_age_tracer")
196 ! Register the tracer for horizontal advection, diffusion, and restarts.
197 call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, tr_desc=cs%tr_desc(m), &
198 registry_diags=.true., restart_cs=restart_cs, &
199 mandatory=.not.cs%tracers_may_reinit, &
200 flux_scale=gv%H_to_m)
201
202 ! Set coupled_tracers to be true (hard-coded above) to provide the surface
203 ! values to the coupler (if any). This is meta-code and its arguments will
204 ! currently (deliberately) give fatal errors if it is used.
205 if (cs%coupled_tracers) &
206 cs%ind_tr(m) = atmos_ocn_coupler_flux(trim(var_name)//'_flux', &
207 flux_type=' ', implementation=' ', caller="register_ideal_age_tracer")
208 enddo
209
210 cs%tr_Reg => tr_reg
211 cs%restart_CSp => restart_cs
212 register_ideal_age_tracer = .true.
213end function register_ideal_age_tracer
214
215!> Sets the ideal age traces to their initial values and sets up the tracer output
216subroutine initialize_ideal_age_tracer(restart, day, G, GV, US, h, diag, OBC, CS, &
217 sponge_CSp)
218 logical, intent(in) :: restart !< .true. if the fields have already
219 !! been read from a restart file.
220 type(time_type), target, intent(in) :: day !< Time of the start of the run.
221 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
222 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
223 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
224 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
225 intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
226 type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
227 !! diagnostic output.
228 type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
229 !! whether, where, and what open boundary
230 !! conditions are used.
231 type(ideal_age_tracer_cs), pointer :: cs !< The control structure returned by a previous
232 !! call to register_ideal_age_tracer.
233 type(sponge_cs), pointer :: sponge_csp !< Pointer to the control structure for the sponges.
234
235! This subroutine initializes the CS%ntr tracer fields in tr(:,:,:,:)
236! and it sets up the tracer output.
237
238 ! Local variables
239 character(len=24) :: name ! A variable's name in a NetCDF file.
240 logical :: ok
241 integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
242 integer :: isdb, iedb, jsdb, jedb
243 logical :: use_real_bl_depth
244
245 if (.not.associated(cs)) return
246 if (cs%ntr < 1) return
247 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
248 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
249 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
250
251 cs%Time => day
252 cs%diag => diag
253 cs%nkbl = max(gv%nkml,1)
254
255 do m=1,cs%ntr
256 call query_vardesc(cs%tr_desc(m), name=name, &
257 caller="initialize_ideal_age_tracer")
258 if ((.not.restart) .or. (cs%tracers_may_reinit .and. .not. &
259 query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
260
261 if (len_trim(cs%IC_file) > 0) then
262 ! Read the tracer concentrations from a netcdf file.
263 if (.not.file_exists(cs%IC_file, g%Domain)) &
264 call mom_error(fatal, "initialize_ideal_age_tracer: "// &
265 "Unable to open "//cs%IC_file)
266
267 if (cs%Z_IC_file) then
268 ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, name,&
269 g, gv, us, -1e34, 0.0) ! CS%land_val(m))
270 if (.not.ok) then
271 ok = tracer_z_init(cs%tr(:,:,:,m), h, cs%IC_file, &
272 trim(name), g, gv, us, -1e34, 0.0) ! CS%land_val(m))
273 if (.not.ok) call mom_error(fatal,"initialize_ideal_age_tracer: "//&
274 "Unable to read "//trim(name)//" from "//&
275 trim(cs%IC_file)//".")
276 endif
277 else
278 call mom_read_data(cs%IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
279 endif
280 else
281 do k=1,nz ; do j=js,je ; do i=is,ie
282 if (g%mask2dT(i,j) < 0.5) then
283 cs%tr(i,j,k,m) = cs%land_val(m)
284 else
285 cs%tr(i,j,k,m) = cs%IC_val(m)
286 endif
287 enddo ; enddo ; enddo
288 endif
289
290 call set_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp)
291 endif ! restart
292 enddo ! Tracer loop
293
294 if (associated(obc)) then
295 ! Steal from updated DOME in the fullness of time.
296 endif
297
298end subroutine initialize_ideal_age_tracer
299
300!> Applies diapycnal diffusion, aging and regeneration at the surface to the ideal age tracers
301subroutine ideal_age_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, &
302 evap_CFL_limit, minimum_forcing_depth, Hbl)
303 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
304 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
305 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
306 intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
307 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
308 intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
309 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
310 intent(in) :: ea !< an array to which the amount of fluid entrained
311 !! from the layer above during this call will be
312 !! added [H ~> m or kg m-2].
313 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
314 intent(in) :: eb !< an array to which the amount of fluid entrained
315 !! from the layer below during this call will be
316 !! added [H ~> m or kg m-2].
317 type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
318 !! and tracer forcing fields. Unused fields have NULL ptrs.
319 real, intent(in) :: dt !< The amount of time covered by this call [T ~> s]
320 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
321 type(ideal_age_tracer_cs), pointer :: cs !< The control structure returned by a previous
322 !! call to register_ideal_age_tracer.
323 real, optional, intent(in) :: evap_cfl_limit !< Limit on the fraction of the water that can
324 !! be fluxed out of the top layer in a timestep [nondim]
325 real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
326 !! fluxes can be applied [H ~> m or kg m-2]
327 real, dimension(SZI_(G),SZJ_(G)), optional, intent(in) :: hbl !< Boundary layer thickness [H ~> m or kg m-2]
328
329! This subroutine applies diapycnal diffusion and any other column
330! tracer physics or chemistry to the tracers from this file.
331! This is a simple example of a set of advected passive tracers.
332
333! The arguments to this subroutine are redundant in that
334! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
335 ! Local variables
336 real, dimension(SZI_(G),SZJ_(G)) :: bl_layers ! Stores number of layers in boundary layer [nondim]
337 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h_work ! Used so that h can be modified [H ~> m or kg m-2]
338 real :: young_val ! The "young" value for the tracers [years] or other units
339 real :: isecs_per_year ! The inverse of the amount of time in a year [T-1 ~> s-1]
340 real :: year ! The time in years [years]
341 real :: layer_frac ! The fraction of the current layer that is within the mixed layer [nondim]
342 integer :: i, j, k, is, ie, js, je, nz, m, nk
343 character(len=255) :: msg
344 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
345
346 if (cs%use_real_BL_depth .and. .not. present(hbl)) then
347 call mom_error(fatal, "Attempting to use real boundary layer depth for ideal age tracers, " &
348 // "but no valid boundary layer scheme was found")
349 endif
350
351 if (cs%use_real_BL_depth .and. present(hbl)) then
352 call count_bl_layers(g, gv, h_old, hbl, bl_layers)
353 endif
354
355 if (.not.associated(cs)) return
356 if (cs%ntr < 1) return
357
358 if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
359 do m=1,cs%ntr
360 do k=1,nz ;do j=js,je ; do i=is,ie
361 h_work(i,j,k) = h_old(i,j,k)
362 enddo ; enddo ; enddo
363 call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m), dt, fluxes, h_work, &
364 evap_cfl_limit, minimum_forcing_depth)
365 call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
366 enddo
367 else
368 do m=1,cs%ntr
369 call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
370 enddo
371 endif
372
373 isecs_per_year = 1.0 / (365.0*86400.0*us%s_to_T)
374 ! Set the surface value of tracer 1 to increase exponentially
375 ! with a 30 year time scale.
376 year = time_to_real(cs%Time, scale=us%s_to_T) * isecs_per_year
377
378 do m=1,cs%ntr
379
380 if (cs%growth_rate(m) == 0.0) then
381 young_val = cs%young_val(m)
382 else
383 young_val = cs%young_val(m) * &
384 exp((year-cs%tracer_start_year(m)) * cs%growth_rate(m))
385 endif
386
387 if (m == cs%BL_residence_num) then
388
389 if (cs%use_real_BL_depth) then
390 do j=js,je ; do i=is,ie
391 nk = floor(bl_layers(i,j))
392
393 do k=1,nk
394 if (g%mask2dT(i,j) > 0.0) then
395 cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt*isecs_per_year
396 else
397 cs%tr(i,j,k,m) = cs%land_val(m)
398 endif
399 enddo
400
401 k = min(nk+1,nz)
402
403 if (g%mask2dT(i,j) > 0.0) then
404 layer_frac = bl_layers(i,j)-nk
405 cs%tr(i,j,k,m) = layer_frac * (cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt &
406 *isecs_per_year) + (1.-layer_frac) * young_val
407 else
408 cs%tr(i,j,k,m) = cs%land_val(m)
409 endif
410
411
412 do k=nk+2,nz
413 if (g%mask2dT(i,j) > 0.0) then
414 cs%tr(i,j,k,m) = young_val
415 else
416 cs%tr(i,j,k,m) = cs%land_val(m)
417 endif
418 enddo
419 enddo ; enddo
420
421 else ! use real BL depth
422 do j=js,je ; do i=is,ie
423 do k=1,cs%nkbl
424 if (g%mask2dT(i,j) > 0.0) then
425 cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt*isecs_per_year
426 else
427 cs%tr(i,j,k,m) = cs%land_val(m)
428 endif
429 enddo
430
431 do k=cs%nkbl+1,nz
432 if (g%mask2dT(i,j) > 0.0) then
433 cs%tr(i,j,k,m) = young_val
434 else
435 cs%tr(i,j,k,m) = cs%land_val(m)
436 endif
437 enddo
438 enddo ; enddo
439
440 endif ! use real BL depth
441
442 else ! if BL residence tracer
443
444 if (cs%use_real_BL_depth) then
445 do j=js,je ; do i=is,ie
446 nk = floor(bl_layers(i,j))
447 do k=1,nk
448 if (g%mask2dT(i,j) > 0.0) then
449 cs%tr(i,j,k,m) = young_val
450 else
451 cs%tr(i,j,k,m) = cs%land_val(m)
452 endif
453 enddo
454
455 k = min(nk+1,nz)
456 if (g%mask2dT(i,j) > 0.0) then
457 layer_frac = bl_layers(i,j)-nk
458 cs%tr(i,j,k,m) = (1.-layer_frac) * (cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt &
459 *isecs_per_year) + layer_frac * young_val
460 else
461 cs%tr(i,j,k,m) = cs%land_val(m)
462 endif
463
464 do k=nk+2,nz
465 if (g%mask2dT(i,j) > 0.0) then
466 cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt*isecs_per_year
467 else
468 cs%tr(i,j,k,m) = cs%land_val(m)
469 endif
470 enddo
471 enddo ; enddo
472
473 else ! use real BL depth
474 do k=1,cs%nkbl ; do j=js,je ; do i=is,ie
475 if (g%mask2dT(i,j) > 0.0) then
476 cs%tr(i,j,k,m) = young_val
477 else
478 cs%tr(i,j,k,m) = cs%land_val(m)
479 endif
480 enddo ; enddo ; enddo
481
482 if (cs%tracer_ages(m) .and. (year>=cs%tracer_start_year(m))) then
483 !$OMP parallel do default(none) shared(is,ie,js,je,CS,nz,G,dt,Isecs_per_year,m)
484 do k=cs%nkbl+1,nz ; do j=js,je ; do i=is,ie
485 cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + g%mask2dT(i,j)*dt*isecs_per_year
486 enddo ; enddo ; enddo
487 endif
488
489
490 endif ! if use real BL depth
491 endif ! if BL residence tracer
492
493 enddo ! loop over all tracers
494
495end subroutine ideal_age_tracer_column_physics
496
497!> Calculates the mass-weighted integral of all tracer stocks, returning the number of stocks it
498!! has calculated. If stock_index is present, only the stock corresponding to that coded index is found.
499function ideal_age_stock(h, stocks, G, GV, CS, names, units, stock_index)
500 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
501 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
502 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
503 intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
504 type(efp_type), dimension(:), intent(out) :: stocks !< the mass-weighted integrated amount of each
505 !! tracer, in kg times concentration units [kg conc].
506 type(ideal_age_tracer_cs), pointer :: cs !< The control structure returned by a previous
507 !! call to register_ideal_age_tracer.
508 character(len=*), dimension(:), intent(out) :: names !< the names of the stocks calculated.
509 character(len=*), dimension(:), intent(out) :: units !< the units of the stocks calculated.
510 integer, optional, intent(in) :: stock_index !< the coded index of a specific stock
511 !! being sought.
512 integer :: ideal_age_stock !< The number of stocks calculated here.
513
514 ! Local variables
515 integer :: m
516
517 ideal_age_stock = 0
518 if (.not.associated(cs)) return
519 if (cs%ntr < 1) return
520
521 if (present(stock_index)) then ; if (stock_index > 0) then
522 ! Check whether this stock is available from this routine.
523
524 ! No stocks from this routine are being checked yet. Return 0.
525 return
526 endif ; endif
527
528 do m=1,cs%ntr
529 call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="ideal_age_stock")
530 units(m) = trim(units(m))//" kg"
531 stocks(m) = global_mass_int_efp(h, g, gv, cs%tr(:,:,:,m), on_pe_only=.true.)
532 enddo
533 ideal_age_stock = cs%ntr
534
535end function ideal_age_stock
536
537!> This subroutine extracts the surface fields from this tracer package that
538!! are to be shared with the atmosphere in coupled configurations.
539!! This particular tracer package does not report anything back to the coupler.
540subroutine ideal_age_tracer_surface_state(sfc_state, h, G, GV, CS)
541 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
542 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
543 type(surface), intent(inout) :: sfc_state !< A structure containing fields that
544 !! describe the surface state of the ocean.
545 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
546 intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
547 type(ideal_age_tracer_cs), pointer :: cs !< The control structure returned by a previous
548 !! call to register_ideal_age_tracer.
549
550 ! This particular tracer package does not report anything back to the coupler.
551 ! The code that is here is just a rough guide for packages that would.
552
553 integer :: m, is, ie, js, je, isd, ied, jsd, jed
554 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
555 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
556
557 if (.not.associated(cs)) return
558
559 if (cs%coupled_tracers) then
560 do m=1,cs%ntr
561 ! This call loads the surface values into the appropriate array in the
562 ! coupler-type structure.
563 call set_coupler_type_data(cs%tr(:,:,1,m), cs%ind_tr(m), sfc_state%tr_fields, &
564 idim=(/isd, is, ie, ied/), jdim=(/jsd, js, je, jed/), turns=g%HI%turns)
565 enddo
566 endif
567
568end subroutine ideal_age_tracer_surface_state
569
570!> Deallocate any memory associated with this tracer package
571subroutine ideal_age_example_end(CS)
572 type(ideal_age_tracer_cs), pointer :: cs !< The control structure returned by a previous
573 !! call to register_ideal_age_tracer.
574
575 if (associated(cs)) then
576 if (associated(cs%tr)) deallocate(cs%tr)
577 deallocate(cs)
578 endif
579end subroutine ideal_age_example_end
580
581subroutine count_bl_layers(G, GV, h, Hbl, BL_layers)
582 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
583 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
584 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
585 intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
586 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: hbl !< Boundary layer thickness [H ~> m or kg m-2]
587 real, dimension(SZI_(G),SZJ_(G)), intent(out) :: bl_layers !< Number of model layers in the boundary layer [nondim]
588
589 real :: current_depth ! Distance from the free surface [H ~> m or kg m-2]
590 integer :: i, j, k, is, ie, js, je, nz, m, nk
591 character(len=255) :: msg
592 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
593
594 bl_layers(:,:) = 0.
595 do j=js,je
596 do i=is,ie
597 current_depth = 0.
598 do k=1,nz
599 current_depth = current_depth + h(i,j,k)
600 if (hbl(i,j) <= current_depth) then
601 bl_layers(i,j) = bl_layers(i,j) + (1.0 - (current_depth - hbl(i,j)) / h(i,j,k))
602 exit
603 else
604 bl_layers(i,j) = bl_layers(i,j) + 1.0
605 endif
606 enddo
607 enddo
608 enddo
609
610end subroutine count_bl_layers
611
612!> \namespace ideal_age_example
613!!
614!! Originally by Robert Hallberg, 2002
615!!
616!! This file contains an example of the code that is needed to set
617!! up and use a set (in this case two) of dynamically passive tracers
618!! for diagnostic purposes. The tracers here are an ideal age tracer
619!! that ages at a rate of 1/year once it is isolated from the surface,
620!! and a vintage tracer, whose surface concentration grows exponen-
621!! with time with a 30-year timescale (similar to CFCs).
622
623end module ideal_age_example