MOM6
advection_test_tracer Module Reference

Detailed Description

This tracer package is used to test advection schemes.

Data Types

type  advection_test_tracer_cs
 The control structure for the advect_test_tracer module. More...
 

Functions/Subroutines

logical function, public register_advection_test_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Register tracer fields and subroutines to be used with MOM. More...
 
subroutine, public initialize_advection_test_tracer (restart, day, G, GV, h, diag, OBC, CS, sponge_CSp)
 Initializes the NTR tracer fields in tr(:,:,:,:) and it sets up the tracer output. More...
 
subroutine, public advection_test_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
 Applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers from this package. This is a simple example of a set of advected passive tracers. More...
 
subroutine, public advection_test_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...
 
integer function, public advection_test_stock (h, stocks, G, GV, CS, names, units, stock_index)
 Calculate the mass-weighted integral of all tracer stocks, returning the number of stocks it has calculated. If the stock_index is present, only the stock corresponding to that coded index is returned. More...
 
subroutine, public advection_test_tracer_end (CS)
 Deallocate memory associated with this module. More...
 

Variables

integer, parameter ntr = 11
 The number of tracers in this module.
 

Function/Subroutine Documentation

◆ advection_test_stock()

integer function, public advection_test_tracer::advection_test_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(advection_test_tracer_cs), pointer  CS,
character(len=*), dimension(:), intent(out)  names,
character(len=*), dimension(:), intent(out)  units,
integer, intent(in), optional  stock_index 
)

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

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_advection_test_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 353 of file advection_test_tracer.F90.

353  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
354  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
355  real, dimension(:), intent(out) :: stocks !< the mass-weighted integrated amount of each
356  !! tracer, in kg times concentration units [kg conc].
357  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
358  type(advection_test_tracer_CS), pointer :: CS !< The control structure returned by a previous
359  !! call to register_advection_test_tracer.
360  character(len=*), dimension(:), intent(out) :: names !< the names of the stocks calculated.
361  character(len=*), dimension(:), intent(out) :: units !< the units of the stocks calculated.
362  integer, optional, intent(in) :: stock_index !< the coded index of a specific stock being sought.
363  integer :: advection_test_stock !< the number of stocks calculated here.
364 
365  integer :: i, j, k, is, ie, js, je, nz, m
366  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
367 
368  advection_test_stock = 0
369  if (.not.associated(cs)) return
370  if (cs%ntr < 1) return
371 
372  if (present(stock_index)) then ; if (stock_index > 0) then
373  ! Check whether this stock is available from this routine.
374 
375  ! No stocks from this routine are being checked yet. Return 0.
376  return
377  endif ; endif
378 
379  do m=1,cs%ntr
380  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="advection_test_stock")
381  stocks(m) = 0.0
382  do k=1,nz ; do j=js,je ; do i=is,ie
383  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
384  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
385  enddo ; enddo ; enddo
386  stocks(m) = gv%H_to_kg_m2 * stocks(m)
387  enddo
388  advection_test_stock = cs%ntr
389 

◆ advection_test_tracer_column_physics()

subroutine, public advection_test_tracer::advection_test_tracer_column_physics ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h_old,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h_new,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  ea,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), 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(advection_test_tracer_cs), pointer  CS,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)

Applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers from this package. This is a simple example of a set of advected passive 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_advection_test_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 262 of file advection_test_tracer.F90.

262  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
263  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
264  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
265  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
266  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
267  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
268  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
269  intent(in) :: ea !< an array to which the amount of fluid entrained
270  !! from the layer above during this call will be
271  !! added [H ~> m or kg m-2].
272  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
273  intent(in) :: eb !< an array to which the amount of fluid entrained
274  !! from the layer below during this call will be
275  !! added [H ~> m or kg m-2].
276  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
277  !! and tracer forcing fields. Unused fields have NULL ptrs.
278  real, intent(in) :: dt !< The amount of time covered by this call [s]
279  type(advection_test_tracer_CS), pointer :: CS !< The control structure returned by a previous
280  !! call to register_advection_test_tracer.
281  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
282  !! be fluxed out of the top layer in a timestep [nondim]
283  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
284  !! fluxes can be applied [m]
285 ! This subroutine applies diapycnal diffusion and any other column
286 ! tracer physics or chemistry to the tracers from this file.
287 ! This is a simple example of a set of advected passive tracers.
288 
289 ! The arguments to this subroutine are redundant in that
290 ! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
291  ! Local variables
292  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
293  real :: b1(SZI_(G)) ! b1 and c1 are variables used by the
294  real :: c1(SZI_(G),SZK_(G)) ! tridiagonal solver.
295  integer :: i, j, k, is, ie, js, je, nz, m
296  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
297 
298  if (.not.associated(cs)) return
299 
300  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
301  do m=1,cs%ntr
302  do k=1,nz ;do j=js,je ; do i=is,ie
303  h_work(i,j,k) = h_old(i,j,k)
304  enddo ; enddo ; enddo
305  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
306  evap_cfl_limit, minimum_forcing_depth)
307  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
308  enddo
309  else
310  do m=1,ntr
311  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
312  enddo
313  endif
314 

◆ advection_test_tracer_end()

subroutine, public advection_test_tracer::advection_test_tracer_end ( type(advection_test_tracer_cs), pointer  CS)

Deallocate memory associated with this module.

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

Definition at line 394 of file advection_test_tracer.F90.

394  type(advection_test_tracer_CS), pointer :: CS !< The control structure returned by a previous
395  !! call to register_advection_test_tracer.
396  integer :: m
397 
398  if (associated(cs)) then
399  if (associated(cs%tr)) deallocate(cs%tr)
400  deallocate(cs)
401  endif

◆ advection_test_tracer_surface_state()

subroutine, public advection_test_tracer::advection_test_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(advection_test_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_advection_test_tracer.

Definition at line 321 of file advection_test_tracer.F90.

321  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
322  type(surface), intent(inout) :: state !< A structure containing fields that
323  !! describe the surface state of the ocean.
324  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
325  intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
326  type(advection_test_tracer_CS), pointer :: CS !< The control structure returned by a previous
327  !! call to register_advection_test_tracer.
328 
329  ! This particular tracer package does not report anything back to the coupler.
330  ! The code that is here is just a rough guide for packages that would.
331 
332  integer :: m, is, ie, js, je, isd, ied, jsd, jed
333  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
334  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
335 
336  if (.not.associated(cs)) return
337 
338  if (cs%coupled_tracers) then
339  do m=1,cs%ntr
340  ! This call loads the surface values into the appropriate array in the
341  ! coupler-type structure.
342  call coupler_type_set_data(cs%tr(:,:,1,m), cs%ind_tr(m), ind_csurf, &
343  state%tr_fields, idim=(/isd, is, ie, ied/), &
344  jdim=(/jsd, js, je, jed/) )
345  enddo
346  endif
347 

◆ initialize_advection_test_tracer()

subroutine, public advection_test_tracer::initialize_advection_test_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(advection_test_tracer_cs), pointer  CS,
type(sponge_cs), pointer  sponge_CSp 
)

Initializes the NTR tracer fields in tr(:,:,:,:) and it 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]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_advection_test_tracer.
sponge_cspPointer to the control structure for the sponges.

Definition at line 167 of file advection_test_tracer.F90.

167  logical, intent(in) :: restart !< .true. if the fields have already
168  !! been read from a restart file.
169  type(time_type), target, intent(in) :: day !< Time of the start of the run.
170  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
171  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
172  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
173  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
174  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
175  !! diagnostic output.
176  type(ocean_OBC_type), pointer :: OBC !< This open boundary condition type specifies
177  !! whether, where, and what open boundary
178  !! conditions are used.
179  type(advection_test_tracer_CS), pointer :: CS !< The control structure returned by a previous
180  !! call to register_advection_test_tracer.
181  type(sponge_CS), pointer :: sponge_CSp !< Pointer to the control structure for the sponges.
182 
183  ! Local variables
184  real, allocatable :: temp(:,:,:)
185  real, pointer, dimension(:,:,:) :: &
186  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
187  obc_tr1_v => null() ! specify the values of tracer 1 that should come
188  ! in through u- and v- points through the open
189  ! boundary conditions, in the same units as tr.
190  character(len=16) :: name ! A variable's name in a NetCDF file.
191  character(len=72) :: longname ! The long name of that variable.
192  character(len=48) :: units ! The dimensions of the variable.
193  character(len=48) :: flux_units ! The units for tracer fluxes, usually
194  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
195  real, pointer :: tr_ptr(:,:,:) => null()
196  real :: PI ! 3.1415926... calculated as 4*atan(1)
197  real :: tr_y ! Initial zonally uniform tracer concentrations.
198  real :: dist2 ! The distance squared from a line [m2].
199  real :: h_neglect ! A thickness that is so small it is usually lost
200  ! in roundoff and can be neglected [H ~> m or kg m-2].
201  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
202  integer :: IsdB, IedB, JsdB, JedB
203  real :: tmpx, tmpy, locx, locy
204 
205  if (.not.associated(cs)) return
206  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
207  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
208  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
209  h_neglect = gv%H_subroundoff
210 
211  cs%diag => diag
212  cs%ntr = ntr
213  do m=1,ntr
214  call query_vardesc(cs%tr_desc(m), name=name, &
215  caller="initialize_advection_test_tracer")
216  if ((.not.restart) .or. (cs%tracers_may_reinit .and. .not. &
217  query_initialized(cs%tr(:,:,:,m), name, cs%restart_CSp))) then
218  do k=1,nz ; do j=js,je ; do i=is,ie
219  cs%tr(i,j,k,m) = 0.0
220  enddo ; enddo ; enddo
221  k=1 ! Square wave
222  do j=js,je ; do i=is,ie
223  if (abs(g%geoLonT(i,j)-cs%x_origin)<0.5*cs%x_width .and. &
224  abs(g%geoLatT(i,j)-cs%y_origin)<0.5*cs%y_width) cs%tr(i,j,k,m) = 1.0
225  enddo ; enddo
226  k=2 ! Triangle wave
227  do j=js,je ; do i=is,ie
228  locx=abs(g%geoLonT(i,j)-cs%x_origin)/cs%x_width
229  locy=abs(g%geoLatT(i,j)-cs%y_origin)/cs%y_width
230  cs%tr(i,j,k,m) = max(0.0, 1.0-locx)*max(0.0, 1.0-locy)
231  enddo ; enddo
232  k=3 ! Cosine bell
233  do j=js,je ; do i=is,ie
234  locx=min(1.0, abs(g%geoLonT(i,j)-cs%x_origin)/cs%x_width)*(acos(0.0)*2.)
235  locy=min(1.0, abs(g%geoLatT(i,j)-cs%y_origin)/cs%y_width)*(acos(0.0)*2.)
236  cs%tr(i,j,k,m) = (1.0+cos(locx))*(1.0+cos(locy))*0.25
237  enddo ; enddo
238  k=4 ! Cylinder
239  do j=js,je ; do i=is,ie
240  locx=abs(g%geoLonT(i,j)-cs%x_origin)/cs%x_width
241  locy=abs(g%geoLatT(i,j)-cs%y_origin)/cs%y_width
242  if (locx**2+locy**2<=1.0) cs%tr(i,j,k,m) = 1.0
243  enddo ; enddo
244  k=5 ! Cut cylinder
245  do j=js,je ; do i=is,ie
246  locx=(g%geoLonT(i,j)-cs%x_origin)/cs%x_width
247  locy=(g%geoLatT(i,j)-cs%y_origin)/cs%y_width
248  if (locx**2+locy**2<=1.0) cs%tr(i,j,k,m) = 1.0
249  if (locx>0.0.and.abs(locy)<0.2) cs%tr(i,j,k,m) = 0.0
250  enddo ; enddo
251  endif ! restart
252  enddo
253 
254 

◆ register_advection_test_tracer()

logical function, public advection_test_tracer::register_advection_test_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(advection_test_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

Register tracer fields and subroutines 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_advection_test_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 advection_test_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(advection_test_tracer_CS), pointer :: CS !< The control structure returned by a previous
71  !! call to register_advection_test_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=80) :: name, longname
79 ! This include declares and sets the variable "version".
80 #include "version_variable.h"
81  character(len=40) :: mdl = "advection_test_tracer" ! This module's name.
82  character(len=200) :: inputdir
83  character(len=48) :: flux_units ! The units for tracer fluxes, usually
84  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
85  real, pointer :: tr_ptr(:,:,:) => null()
86  logical :: register_advection_test_tracer
87  integer :: isd, ied, jsd, jed, nz, m
88  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
89 
90  if (associated(cs)) then
91  call mom_error(warning, "register_advection_test_tracer called with an "// &
92  "associated control structure.")
93  return
94  endif
95  allocate(cs)
96 
97  ! Read all relevant parameters and write them to the model log.
98  call log_version(param_file, mdl, version, "")
99 
100  call get_param(param_file, mdl, "ADVECTION_TEST_X_ORIGIN", cs%x_origin, &
101  "The x-coorindate of the center of the test-functions.", default=0.)
102  call get_param(param_file, mdl, "ADVECTION_TEST_Y_ORIGIN", cs%y_origin, &
103  "The y-coorindate of the center of the test-functions.", default=0.)
104  call get_param(param_file, mdl, "ADVECTION_TEST_X_WIDTH", cs%x_width, &
105  "The x-width of the test-functions.", default=0.)
106  call get_param(param_file, mdl, "ADVECTION_TEST_Y_WIDTH", cs%y_width, &
107  "The y-width of the test-functions.", default=0.)
108  call get_param(param_file, mdl, "ADVECTION_TEST_TRACER_IC_FILE", cs%tracer_IC_file, &
109  "The name of a file from which to read the initial "//&
110  "conditions for the tracers, or blank to initialize "//&
111  "them internally.", default=" ")
112 
113  if (len_trim(cs%tracer_IC_file) >= 1) then
114  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
115  cs%tracer_IC_file = trim(slasher(inputdir))//trim(cs%tracer_IC_file)
116  call log_param(param_file, mdl, "INPUTDIR/ADVECTION_TEST_TRACER_IC_FILE", &
117  cs%tracer_IC_file)
118  endif
119  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
120  "If true, sponges may be applied anywhere in the domain. "//&
121  "The exact location and properties of those sponges are "//&
122  "specified from MOM_initialization.F90.", default=.false.)
123 
124  call get_param(param_file, mdl, "TRACERS_MAY_REINIT", cs%tracers_may_reinit, &
125  "If true, tracers may go through the initialization code "//&
126  "if they are not found in the restart files. Otherwise "//&
127  "it is a fatal error if the tracers are not found in the "//&
128  "restart files of a restarted run.", default=.false.)
129 
130 
131  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
132 
133  do m=1,ntr
134  if (m < 10) then ; write(name,'("tr",I1.1)') m
135  else ; write(name,'("tr",I2.2)') m ; endif
136  write(longname,'("Concentration of Tracer ",I2.2)') m
137  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
138  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
139  else ; flux_units = "kg s-1" ; endif
140 
141 
142  ! This is needed to force the compiler not to do a copy in the registration
143  ! calls. Curses on the designers and implementers of Fortran90.
144  tr_ptr => cs%tr(:,:,:,m)
145  ! Register the tracer for horizontal advection, diffusion, and restarts.
146  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
147  name=name, longname=longname, units="kg kg-1", &
148  registry_diags=.true., flux_units=flux_units, &
149  restart_cs=restart_cs, mandatory=.not.cs%tracers_may_reinit)
150 
151  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
152  ! values to the coupler (if any). This is meta-code and its arguments will
153  ! currently (deliberately) give fatal errors if it is used.
154  if (cs%coupled_tracers) &
155  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
156  flux_type=' ', implementation=' ', caller="register_advection_test_tracer")
157  enddo
158 
159  cs%tr_Reg => tr_reg
160  cs%restart_CSp => restart_cs
161  register_advection_test_tracer = .true.