MOM6
boundary_impulse_tracer Module Reference

Detailed Description

Implements a boundary impulse response tracer to calculate Green's functions.

Boundary Impulse Response Tracer and Transit Time Distributions

Transit time distributions (TTD) are the Green's function solution of the passive tracer equation between the oceanic surface and interior. The name derives from the idea that the 'age' (e.g. time since last contact with the atmosphere) of a water parcel is best characterized as a distribution of ages because water parcels leaving the surface arrive at a particular interior point at different times. The more commonly used ideal age tracer is the first moment of the TTD, equivalently referred to as the mean age.

A boundary impulse response (BIR) is a passive tracer whose surface boundary condition is a rectangle function with width \(\Delta t\). In the case of unsteady flow, multiple BIRs, initiated at different times in the model can be used to infer the transit time distribution or Green's function between the oceanic surface and interior. In the case of steady or cyclostationary flow, a single BIR is sufficient.

In the References section, both the theoretical discussion of TTDs and BIRs are listed along with modeling studies which have this used framework in scientific investigations

Run-time parameters

-DO_BOUNDARY_IMPULSE_TRACER: Enables the boundary impulse tracer model -IMPULSE_SOURCE_TIME: Length of time that the surface layer acts as a source of the BIR tracer

References

and BIR Theory

-Holzer, M., and T.M. Hall, 2000: Transit-time and tracer-age distributions in geophysical flows. J. Atmos. Sci., 57, 3539-3558, doi:10.1175/1520-0469(2000)057<3539:TTATAD>2.0.CO;2. -T.W.N. Haine, H. Zhang, D.W. Waugh, M. Holzer, On transit-time distributions in unsteady circulation models, Ocean Modelling, Volume 21, Issues 1–2, 2008, Pages 35-45, ISSN 1463-5003 http://dx.doi.org/10.1016/j.ocemod.2007.11.004.

Modelling applications

-Peacock, S., and M. Maltrud (2006), Transit-time distributions in a global ocean model, J. Phys. Oceanogr., 36(3), 474–495, doi:10.1175/JPO2860.1. -Maltrud, M., Bryan, F. & Peacock, Boundary impulse response functions in a century-long eddying global ocean simulation, S. Environ Fluid Mech (2010) 10: 275. doi:10.1007/s10652-009-9154-3

Data Types

type  boundary_impulse_tracer_cs
 The control structure for the boundary impulse tracer package. More...
 

Functions/Subroutines

logical function, public register_boundary_impulse_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Read in runtime options and add boundary impulse tracer to tracer registry. More...
 
subroutine, public initialize_boundary_impulse_tracer (restart, day, G, GV, h, diag, OBC, CS, sponge_CSp, tv)
 Initialize tracer from restart or set to 1 at surface to initialize. More...
 
subroutine, public boundary_impulse_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, tv, debug, evap_CFL_limit, minimum_forcing_depth)
 Apply source or sink at boundary and do vertical diffusion. More...
 
integer function, public boundary_impulse_stock (h, stocks, G, GV, CS, names, units, stock_index)
 Calculate total inventory of tracer. More...
 
subroutine, public boundary_impulse_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 boundary_impulse_tracer_end (CS)
 Performs finalization of boundary impulse tracer. More...
 

Variables

integer, parameter ntr_max = 1
 NTR_MAX is the maximum number of tracers in this module. More...
 

Function/Subroutine Documentation

◆ boundary_impulse_stock()

integer function, public boundary_impulse_tracer::boundary_impulse_stock ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
real, dimension(:), intent(out)  stocks,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(boundary_impulse_tracer_cs), pointer  CS,
character(len=*), dimension(:), intent(out)  names,
character(len=*), dimension(:), intent(out)  units,
integer, intent(in), optional  stock_index 
)

Calculate total inventory of tracer.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical 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].
csThe control structure returned by a previous call to register_boundary_impulse_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
Return value: the number of stocks calculated here.

Definition at line 287 of file boundary_impulse_tracer.F90.

287  type(ocean_grid_type), intent(in ) :: G !< The ocean's grid structure
288  type(verticalGrid_type), intent(in ) :: GV !< The ocean's vertical grid structure
289  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in ) :: h !< Layer thicknesses [H ~> m or kg m-2]
290  real, dimension(:), intent( out) :: stocks !< the mass-weighted integrated amount of each
291  !! tracer, in kg times concentration units [kg conc].
292  type(boundary_impulse_tracer_CS), pointer :: CS !< The control structure returned by a previous
293  !! call to register_boundary_impulse_tracer.
294  character(len=*), dimension(:), intent( out) :: names !< The names of the stocks calculated.
295  character(len=*), dimension(:), intent( out) :: units !< The units of the stocks calculated.
296  integer, optional, intent(in ) :: stock_index !< The coded index of a specific stock
297  !! being sought.
298  integer :: boundary_impulse_stock !< Return value: the number of stocks calculated here.
299 
300 ! This function calculates the mass-weighted integral of all tracer stocks,
301 ! returning the number of stocks it has calculated. If the stock_index
302 ! is present, only the stock corresponding to that coded index is returned.
303 
304  ! Local variables
305  integer :: i, j, k, is, ie, js, je, nz, m
306  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
307 
308  boundary_impulse_stock = 0
309  if (.not.associated(cs)) return
310  if (cs%ntr < 1) return
311 
312  if (present(stock_index)) then ; if (stock_index > 0) then
313  ! Check whether this stock is available from this routine.
314 
315  ! No stocks from this routine are being checked yet. Return 0.
316  return
317  endif ; endif
318 
319  do m=1,1
320  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="boundary_impulse_stock")
321  units(m) = trim(units(m))//" kg"
322  stocks(m) = 0.0
323  do k=1,nz ; do j=js,je ; do i=is,ie
324  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
325  (g%mask2dT(i,j) * g%US%L_to_m**2*g%areaT(i,j) * h(i,j,k))
326  enddo ; enddo ; enddo
327  stocks(m) = gv%H_to_kg_m2 * stocks(m)
328  enddo
329 
330  boundary_impulse_stock = cs%ntr
331 

References mom_io::query_vardesc().

Here is the call graph for this function:

◆ boundary_impulse_tracer_column_physics()

subroutine, public boundary_impulse_tracer::boundary_impulse_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(unit_scale_type), intent(in)  US,
type(boundary_impulse_tracer_cs), pointer  CS,
type(thermo_var_ptrs), intent(in)  tv,
logical, intent(in)  debug,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)

Apply source or sink at boundary and do vertical diffusion.

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 [T ~> s]
[in]usA dimensional unit scaling type
csThe control structure returned by a previous call to register_boundary_impulse_tracer.
[in]tvA structure pointing to various thermodynamic variables
[in]debugIf true calculate checksums
[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 [H ~> m or kg m-2]

Definition at line 209 of file boundary_impulse_tracer.F90.

209  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
210  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
211  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
212  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
213  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
214  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
215  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
216  intent(in) :: ea !< an array to which the amount of fluid entrained
217  !! from the layer above during this call will be
218  !! added [H ~> m or kg m-2].
219  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
220  intent(in) :: eb !< an array to which the amount of fluid entrained
221  !! from the layer below during this call will be
222  !! added [H ~> m or kg m-2].
223  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
224  !! and tracer forcing fields. Unused fields have NULL ptrs.
225  real, intent(in) :: dt !< The amount of time covered by this call [T ~> s]
226  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
227  type(boundary_impulse_tracer_CS), pointer :: CS !< The control structure returned by a previous
228  !! call to register_boundary_impulse_tracer.
229  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
230  !! thermodynamic variables
231  logical, intent(in) :: debug !< If true calculate checksums
232  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
233  !! be fluxed out of the top layer in a timestep [nondim]
234  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
235  !! fluxes can be applied [H ~> m or kg m-2]
236 
237 ! This subroutine applies diapycnal diffusion and any other column
238 ! tracer physics or chemistry to the tracers from this file.
239 ! This is a simple example of a set of advected passive tracers.
240 
241 ! The arguments to this subroutine are redundant in that
242 ! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
243 
244  ! Local variables
245  real :: Isecs_per_year = 1.0 / (365.0*86400.0)
246  real :: year, h_total, scale, htot, Ih_limit
247  integer :: secs, days
248  integer :: i, j, k, is, ie, js, je, nz, m, k_max
249  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
250 
251  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
252 
253  if (.not.associated(cs)) return
254  if (cs%ntr < 1) return
255 
256  ! This uses applyTracerBoundaryFluxesInOut, usually in ALE mode
257  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
258  do k=1,nz ;do j=js,je ; do i=is,ie
259  h_work(i,j,k) = h_old(i,j,k)
260  enddo ; enddo ; enddo
261  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,1), dt, fluxes, h_work, &
262  evap_cfl_limit, minimum_forcing_depth)
263  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,1), g, gv)
264  else
265  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,1), g, gv)
266  endif
267 
268  ! Set surface conditions
269  do m=1,1
270  if (cs%remaining_source_time>0.0) then
271  do k=1,cs%nkml ; do j=js,je ; do i=is,ie
272  cs%tr(i,j,k,m) = 1.0
273  enddo ; enddo ; enddo
274  cs%remaining_source_time = cs%remaining_source_time-us%T_to_s*dt
275  else
276  do k=1,cs%nkml ; do j=js,je ; do i=is,ie
277  cs%tr(i,j,k,m) = 0.0
278  enddo ; enddo ; enddo
279  endif
280 
281  enddo
282 

References mom_tracer_diabatic::applytracerboundaryfluxesinout(), and mom_tracer_diabatic::tracer_vertdiff().

Here is the call graph for this function:

◆ boundary_impulse_tracer_end()

subroutine, public boundary_impulse_tracer::boundary_impulse_tracer_end ( type(boundary_impulse_tracer_cs), pointer  CS)

Performs finalization of boundary impulse tracer.

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

Definition at line 369 of file boundary_impulse_tracer.F90.

369  type(boundary_impulse_tracer_CS), pointer :: CS !< The control structure returned by a previous
370  !! call to register_boundary_impulse_tracer.
371  integer :: m
372 
373  if (associated(cs)) then
374  if (associated(cs%tr)) deallocate(cs%tr)
375  deallocate(cs)
376  endif

◆ boundary_impulse_tracer_surface_state()

subroutine, public boundary_impulse_tracer::boundary_impulse_tracer_surface_state ( type(surface), intent(inout)  state,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(ocean_grid_type), intent(in)  G,
type(boundary_impulse_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_boundary_impulse_tracer.

Definition at line 338 of file boundary_impulse_tracer.F90.

338  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
339  type(surface), intent(inout) :: state !< A structure containing fields that
340  !! describe the surface state of the ocean.
341  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
342  intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
343  type(boundary_impulse_tracer_CS), pointer :: CS !< The control structure returned by a previous
344  !! call to register_boundary_impulse_tracer.
345 
346  ! This particular tracer package does not report anything back to the coupler.
347  ! The code that is here is just a rough guide for packages that would.
348 
349  integer :: m, is, ie, js, je, isd, ied, jsd, jed
350  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
351  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
352 
353  if (.not.associated(cs)) return
354 
355  if (cs%coupled_tracers) then
356  do m=1,cs%ntr
357  ! This call loads the surface values into the appropriate array in the
358  ! coupler-type structure.
359  call coupler_type_set_data(cs%tr(:,:,1,m), cs%ind_tr(m), ind_csurf, &
360  state%tr_fields, idim=(/isd, is, ie, ied/), &
361  jdim=(/jsd, js, je, jed/) )
362  enddo
363  endif
364 

◆ initialize_boundary_impulse_tracer()

subroutine, public boundary_impulse_tracer::initialize_boundary_impulse_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,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(boundary_impulse_tracer_cs), pointer  CS,
type(sponge_cs), pointer  sponge_CSp,
type(thermo_var_ptrs), intent(in)  tv 
)

Initialize tracer from restart or set to 1 at surface to initialize.

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]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_boundary_impulse_tracer.
sponge_cspPointer to the control structure for the sponges.
[in]tvA structure pointing to various thermodynamic variables

Definition at line 153 of file boundary_impulse_tracer.F90.

153  logical, intent(in) :: restart !< .true. if the fields have already
154  !! been read from a restart file.
155  type(time_type), target, intent(in) :: day !< Time of the start of the run.
156  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
157  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
158  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
159  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
160  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
161  !! diagnostic output.
162  type(ocean_OBC_type), pointer :: OBC !< This open boundary condition type specifies
163  !! whether, where, and what open boundary
164  !! conditions are used.
165  type(boundary_impulse_tracer_CS), pointer :: CS !< The control structure returned by a previous
166  !! call to register_boundary_impulse_tracer.
167  type(sponge_CS), pointer :: sponge_CSp !< Pointer to the control structure for the sponges.
168  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
169  !! thermodynamic variables
170  ! Local variables
171  character(len=16) :: name ! A variable's name in a NetCDF file.
172  character(len=72) :: longname ! The long name of that variable.
173  character(len=48) :: units ! The dimensions of the variable.
174  character(len=48) :: flux_units ! The units for age tracer fluxes, either
175  ! years m3 s-1 or years kg s-1.
176  logical :: OK
177  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
178  integer :: IsdB, IedB, JsdB, JedB
179 
180  if (.not.associated(cs)) return
181  if (cs%ntr < 1) return
182  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
183  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
184  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
185 
186  cs%Time => day
187  cs%diag => diag
188  name = "boundary_impulse"
189 
190  do m=1,cs%ntr
191  call query_vardesc(cs%tr_desc(m), name=name, caller="initialize_boundary_impulse_tracer")
192  if ((.not.restart) .or. (.not. &
193  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
194  do k=1,cs%nkml ; do j=jsd,jed ; do i=isd,ied
195  cs%tr(i,j,k,m) = 1.0
196  enddo ; enddo ; enddo
197  endif
198  enddo ! Tracer loop
199 
200  if (associated(obc)) then
201  ! Steal from updated DOME in the fullness of time.
202  endif
203 

References mom_io::query_vardesc().

Here is the call graph for this function:

◆ register_boundary_impulse_tracer()

logical function, public boundary_impulse_tracer::register_boundary_impulse_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(boundary_impulse_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

Read in runtime options and add boundary impulse tracer to tracer registry.

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_boundary_impulse_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 67 of file boundary_impulse_tracer.F90.

67  type(hor_index_type), intent(in ) :: HI !< A horizontal index type structure
68  type(verticalGrid_type), intent(in ) :: GV !< The ocean's vertical grid structure
69  type(param_file_type), intent(in ) :: param_file !< A structure to parse for run-time parameters
70  type(boundary_impulse_tracer_CS), pointer :: CS !< The control structure returned by a previous
71  !! call to register_boundary_impulse_tracer.
72  type(tracer_registry_type), pointer :: tr_Reg !< A pointer that is set to point to the control
73  !! structure for the tracer advection and
74  !! diffusion module
75  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure
76 
77  ! Local variables
78  character(len=40) :: mdl = "boundary_impulse_tracer" ! This module's name.
79  character(len=200) :: inputdir ! The directory where the input files are.
80  character(len=48) :: var_name ! The variable's name.
81  character(len=3) :: name_tag ! String for creating identifying boundary_impulse
82  character(len=48) :: flux_units ! The units for tracer fluxes, usually
83  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
84  ! This include declares and sets the variable "version".
85 #include "version_variable.h"
86  real, pointer :: tr_ptr(:,:,:) => null()
87  real, pointer :: rem_time_ptr => null()
88  logical :: register_boundary_impulse_tracer
89  integer :: isd, ied, jsd, jed, nz, m, i, j
90  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
91 
92  if (associated(cs)) then
93  call mom_error(warning, "register_boundary_impulse_tracer called with an "// &
94  "associated control structure.")
95  return
96  endif
97  allocate(cs)
98 
99  ! Read all relevant parameters and write them to the model log.
100  call log_version(param_file, mdl, version, "")
101  call get_param(param_file, mdl, "IMPULSE_SOURCE_TIME", cs%remaining_source_time, &
102  "Length of time for the boundary tracer to be injected "//&
103  "into the mixed layer. After this time has elapsed, the "//&
104  "surface becomes a sink for the boundary impulse tracer.", &
105  default=31536000.0)
106  call get_param(param_file, mdl, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
107  "If true, tracers may go through the initialization code "//&
108  "if they are not found in the restart files. Otherwise "//&
109  "it is a fatal error if the tracers are not found in the "//&
110  "restart files of a restarted run.", default=.false.)
111  cs%ntr = ntr_max
112  allocate(cs%tr(isd:ied,jsd:jed,nz,cs%ntr)) ; cs%tr(:,:,:,:) = 0.0
113 
114  cs%nkml = max(gv%nkml,1)
115 
116  do m=1,cs%ntr
117  ! This is needed to force the compiler not to do a copy in the registration
118  ! calls. Curses on the designers and implementers of Fortran90.
119  cs%tr_desc(m) = var_desc(trim("boundary_impulse"), "kg kg-1", &
120  "Boundary impulse tracer", caller=mdl)
121  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
122  else ; flux_units = "kg s-1" ; endif
123 
124  tr_ptr => cs%tr(:,:,:,m)
125  call query_vardesc(cs%tr_desc(m), name=var_name, caller="register_boundary_impulse_tracer")
126  ! Register the tracer for horizontal advection, diffusion, and restarts.
127  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, tr_desc=cs%tr_desc(m), &
128  registry_diags=.true., flux_units=flux_units, &
129  restart_cs=restart_cs, mandatory=.not.cs%tracers_may_reinit)
130 
131  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
132  ! values to the coupler (if any). This is meta-code and its arguments will
133  ! currently (deliberately) give fatal errors if it is used.
134  if (cs%coupled_tracers) &
135  cs%ind_tr(m) = aof_set_coupler_flux(trim(var_name)//'_flux', &
136  flux_type=' ', implementation=' ', caller="register_boundary_impulse_tracer")
137  enddo
138  ! Register remaining source time as a restart field
139  rem_time_ptr => cs%remaining_source_time
140  call register_restart_field(rem_time_ptr, "bir_remain_time", &
141  .not.cs%tracers_may_reinit, restart_cs, &
142  "Remaining time to apply BIR source", "s")
143 
144  cs%tr_Reg => tr_reg
145  cs%restart_CSp => restart_cs
146  register_boundary_impulse_tracer = .true.
147 

References atmos_ocean_fluxes_mod::aof_set_coupler_flux(), mom_error_handler::mom_error(), ntr_max, mom_io::query_vardesc(), mom_tracer_registry::register_tracer(), and mom_io::var_desc().

Here is the call graph for this function:

Variable Documentation

◆ ntr_max

integer, parameter boundary_impulse_tracer::ntr_max = 1
private

NTR_MAX is the maximum number of tracers in this module.

Definition at line 37 of file boundary_impulse_tracer.F90.

37 integer, parameter :: NTR_MAX = 1

Referenced by register_boundary_impulse_tracer().