MOM6
isomip_tracer Module Reference

Detailed Description

This module contains the routines used to set up and use a set of (one for now) dynamically passive tracers. For now, just one passive tracer is injected in the sponge layer. Set up and use passive tracers requires the following: (1) register_ISOMIP_tracer (2)

Data Types

type  isomip_tracer_cs
 tracer control structure More...
 
type  p3d
 

Functions/Subroutines

logical function, public register_isomip_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 This subroutine is used to register tracer fields. More...
 
subroutine, public initialize_isomip_tracer (restart, day, G, GV, h, diag, OBC, CS, ALE_sponge_CSp, diag_to_Z_CSp)
 Initializes the NTR tracer fields in tr(:,:,:,:) More...
 
subroutine, public isomip_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, evap_CFL_limit, minimum_forcing_depth)
 This subroutine applies diapycnal diffusion and any other column. More...
 
subroutine, public isomip_tracer_surface_state (state, h, G, CS)
 This particular tracer package does not report anything back to the coupler. More...
 
subroutine, public isomip_tracer_end (CS)
 

Variables

integer, parameter ntr = 1
 

Function/Subroutine Documentation

◆ initialize_isomip_tracer()

subroutine, public isomip_tracer::initialize_isomip_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( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(diag_ctrl), intent(in), target  diag,
type(ocean_obc_type), pointer  OBC,
type(isomip_tracer_cs), pointer  CS,
type(ale_sponge_cs), pointer  ALE_sponge_CSp,
type(diag_to_z_cs), pointer  diag_to_Z_CSp 
)

Initializes the NTR tracer fields in tr(:,:,:,:)

Parameters
[in]gGrid structure.
[in]gvThe ocean's vertical grid structure.
[in]restart.true. if the fields have already been read from a restart file.
[in]dayTime of the start of the run.
[in]hLayer thickness, in m or kg m-2.
obcThis open boundary condition type specifies whether, where, and what open boundary conditions are used. This is not being used for now.
csThe control structure returned by a previous call to ISOMIP_register_tracer.
ale_sponge_cspA pointer to the control structure for the sponges, if they are in use. Otherwise this may be unassociated.
diag_to_z_cspA pointer to the control structure for diagnostics in depth space.

Definition at line 172 of file ISOMIP_tracer.F90.

References mom_tracer_registry::add_tracer_diagnostics(), mom_error_handler::mom_error(), ntr, mom_io::query_vardesc(), and mom_diag_to_z::register_z_tracer().

172 
173  type(ocean_grid_type), intent(in) :: g !< Grid structure.
174  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
175  logical, intent(in) :: restart !< .true. if the fields have already been read from a restart file.
176  type(time_type), target, intent(in) :: day !< Time of the start of the run.
177  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness, in m or kg m-2.
178  type(diag_ctrl), target, intent(in) :: diag
179  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies whether, where, and what open boundary conditions are used. This is not being used for now.
180  type(isomip_tracer_cs), pointer :: cs !< The control structure returned by a previous call to ISOMIP_register_tracer.
181  type(ale_sponge_cs), pointer :: ale_sponge_csp !< A pointer to the control structure for the sponges, if they are in use. Otherwise this may be unassociated.
182  type(diag_to_z_cs), pointer :: diag_to_z_csp !< A pointer to the control structure for diagnostics in depth space.
183 
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, in m2.
199  real :: h_neglect ! A thickness that is so small it is usually lost
200  ! in roundoff and can be neglected, in m.
201  real :: e(szk_(g)+1), e_top, e_bot, d_tr
202  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
203  integer :: isdb, iedb, jsdb, jedb
204 
205  if (.not.associated(cs)) return
206  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%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%Time => day
212  cs%diag => diag
213 
214  if (.not.restart) then
215  if (len_trim(cs%tracer_IC_file) >= 1) then
216  ! Read the tracer concentrations from a netcdf file.
217  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
218  call mom_error(fatal, "ISOMIP_initialize_tracer: Unable to open "// &
219  cs%tracer_IC_file)
220  do m=1,ntr
221  call query_vardesc(cs%tr_desc(m), name, caller="initialize_ISOMIP_tracer")
222  call read_data(cs%tracer_IC_file, trim(name), &
223  cs%tr(:,:,:,m), domain=g%Domain%mpp_domain)
224  enddo
225  else
226  do m=1,ntr
227  do k=1,nz ; do j=js,je ; do i=is,ie
228  cs%tr(i,j,k,m) = 0.0
229  enddo ; enddo ; enddo
230  enddo
231  endif
232  endif ! restart
233 
234 ! the following does not work in layer mode yet
235 !! if ( CS%use_sponge ) then
236  ! If sponges are used, this example damps tracers in sponges in the
237  ! northern half of the domain to 1 and tracers in the southern half
238  ! to 0. For any tracers that are not damped in the sponge, the call
239  ! to set_up_sponge_field can simply be omitted.
240 ! if (.not.associated(ALE_sponge_CSp)) &
241 ! call MOM_error(FATAL, "ISOMIP_initialize_tracer: "// &
242 ! "The pointer to ALEsponge_CSp must be associated if SPONGE is defined.")
243 
244 ! allocate(temp(G%isd:G%ied,G%jsd:G%jed,nz))
245 
246 ! do j=js,je ; do i=is,ie
247 ! if (G%geoLonT(i,j) >= 790.0 .AND. G%geoLonT(i,j) <= 800.0) then
248 ! temp(i,j,:) = 1.0
249 ! else
250 ! temp(i,j,:) = 0.0
251 ! endif
252 ! enddo ; enddo
253 
254  ! do m=1,NTR
255 ! do m=1,1
256  ! This is needed to force the compiler not to do a copy in the sponge
257  ! calls. Curses on the designers and implementers of Fortran90.
258 ! tr_ptr => CS%tr(:,:,:,m)
259 ! call set_up_ALE_sponge_field(temp, G, tr_ptr, ALE_sponge_CSp)
260 ! enddo
261 ! deallocate(temp)
262 ! endif
263 
264  ! This needs to be changed if the units of tracer are changed above.
265  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
266  else ; flux_units = "kg s-1" ; endif
267 
268  do m=1,ntr
269  ! Register the tracer for the restart file.
270  call query_vardesc(cs%tr_desc(m), name, units=units, longname=longname, &
271  caller="initialize_ISOMIP_tracer")
272  cs%id_tracer(m) = register_diag_field("ocean_model", trim(name), cs%diag%axesTL, &
273  day, trim(longname) , trim(units))
274  cs%id_tr_adx(m) = register_diag_field("ocean_model", trim(name)//"_adx", &
275  cs%diag%axesCuL, day, trim(longname)//" advective zonal flux" , &
276  trim(flux_units))
277  cs%id_tr_ady(m) = register_diag_field("ocean_model", trim(name)//"_ady", &
278  cs%diag%axesCvL, day, trim(longname)//" advective meridional flux" , &
279  trim(flux_units))
280  cs%id_tr_dfx(m) = register_diag_field("ocean_model", trim(name)//"_dfx", &
281  cs%diag%axesCuL, day, trim(longname)//" diffusive zonal flux" , &
282  trim(flux_units))
283  cs%id_tr_dfy(m) = register_diag_field("ocean_model", trim(name)//"_dfy", &
284  cs%diag%axesCvL, day, trim(longname)//" diffusive zonal flux" , &
285  trim(flux_units))
286  if (cs%id_tr_adx(m) > 0) call safe_alloc_ptr(cs%tr_adx(m)%p,isdb,iedb,jsd,jed,nz)
287  if (cs%id_tr_ady(m) > 0) call safe_alloc_ptr(cs%tr_ady(m)%p,isd,ied,jsdb,jedb,nz)
288  if (cs%id_tr_dfx(m) > 0) call safe_alloc_ptr(cs%tr_dfx(m)%p,isdb,iedb,jsd,jed,nz)
289  if (cs%id_tr_dfy(m) > 0) call safe_alloc_ptr(cs%tr_dfy(m)%p,isd,ied,jsdb,jedb,nz)
290 
291 ! Register the tracer for horizontal advection & diffusion.
292  if ((cs%id_tr_adx(m) > 0) .or. (cs%id_tr_ady(m) > 0) .or. &
293  (cs%id_tr_dfx(m) > 0) .or. (cs%id_tr_dfy(m) > 0)) &
294  call add_tracer_diagnostics(name, cs%tr_Reg, cs%tr_adx(m)%p, &
295  cs%tr_ady(m)%p, cs%tr_dfx(m)%p, cs%tr_dfy(m)%p)
296 
297  call register_z_tracer(cs%tr(:,:,:,m), trim(name), longname, units, &
298  day, g, diag_to_z_csp)
299  enddo
300 
Here is the call graph for this function:

◆ isomip_tracer_column_physics()

subroutine, public isomip_tracer::isomip_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(isomip_tracer_cs), pointer  CS,
real, intent(in), optional  evap_CFL_limit,
real, intent(in), optional  minimum_forcing_depth 
)

This subroutine applies diapycnal diffusion and any other column.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]dtThe amount of time covered by this call, in s

Definition at line 308 of file ISOMIP_tracer.F90.

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

308  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
309  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
310  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h_old, h_new, ea, eb
311  type(forcing), intent(in) :: fluxes
312  real, intent(in) :: dt !< The amount of time covered by this call, in s
313  type(isomip_tracer_cs), pointer :: cs
314  real, optional,intent(in) :: evap_cfl_limit
315  real, optional,intent(in) :: minimum_forcing_depth
316 
317 ! Arguments: h_old - Layer thickness before entrainment, in m or kg m-2.
318 ! (in) h_new - Layer thickness after entrainment, in m or kg m-2.
319 ! (in) ea - an array to which the amount of fluid entrained
320 ! from the layer above during this call will be
321 ! added, in m or kg m-2.
322 ! (in) eb - an array to which the amount of fluid entrained
323 ! from the layer below during this call will be
324 ! added, in m or kg m-2.
325 ! (in) fluxes - A structure containing pointers to any possible
326 ! forcing fields. Unused fields have NULL ptrs.
327 ! (in) dt - The amount of time covered by this call, in s.
328 ! (in) G - The ocean's grid structure.
329 ! (in) GV - The ocean's vertical grid structure.
330 ! (in) CS - The control structure returned by a previous call to
331 ! ISOMIP_register_tracer.
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 
336  real :: mmax
337  real :: b1(szi_(g)) ! b1 and c1 are variables used by the
338  real :: c1(szi_(g),szk_(g)) ! tridiagonal solver.
339  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
340  real :: melt(szi_(g),szj_(g)) ! melt water (positive for melting
341  ! negative for freezing)
342  integer :: i, j, k, is, ie, js, je, nz, m
343  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
344 
345  if (.not.associated(cs)) return
346 
347  melt(:,:) = fluxes%iceshelf_melt
348 
349  ! max. melt
350  mmax = maxval(melt(is:ie,js:je))
351  call max_across_pes(mmax)
352  !write(*,*)'max melt', mmax
353  ! dye melt water (m=1), dye = 1 if melt=max(melt)
354  do m=1,ntr
355  do j=js,je ; do i=is,ie
356  if (melt(i,j) > 0.0) then ! melting
357  !write(*,*)'i,j,melt,melt/mmax',i,j,melt(i,j),melt(i,j)/mmax
358  cs%tr(i,j,1:2,m) = melt(i,j)/mmax ! inject dye in the ML
359  else ! freezing
360  cs%tr(i,j,1:2,m) = 0.0
361  endif
362  enddo ; enddo
363  enddo
364 
365  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
366  do m=1,ntr
367  do k=1,nz ;do j=js,je ; do i=is,ie
368  h_work(i,j,k) = h_old(i,j,k)
369  enddo ; enddo ; enddo;
370  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m) , dt, fluxes, h_work, &
371  evap_cfl_limit, minimum_forcing_depth)
372  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
373  enddo
374  else
375  do m=1,ntr
376  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
377  enddo
378  endif
379 
380  if (cs%mask_tracers) then
381  do m = 1,ntr ; if (cs%id_tracer(m) > 0) then
382  do k=1,nz ; do j=js,je ; do i=is,ie
383  if (h_new(i,j,k) < 1.1*gv%Angstrom) then
384  cs%tr_aux(i,j,k,m) = cs%land_val(m)
385  else
386  cs%tr_aux(i,j,k,m) = cs%tr(i,j,k,m)
387  endif
388  enddo ; enddo ; enddo
389  endif ; enddo
390  endif
391 
392  do m=1,ntr
393  if (cs%mask_tracers) then
394  if (cs%id_tracer(m)>0) &
395  call post_data(cs%id_tracer(m),cs%tr_aux(:,:,:,m),cs%diag)
396  else
397  if (cs%id_tracer(m)>0) &
398  call post_data(cs%id_tracer(m),cs%tr(:,:,:,m),cs%diag)
399  endif
400  if (cs%id_tr_adx(m)>0) &
401  call post_data(cs%id_tr_adx(m),cs%tr_adx(m)%p(:,:,:),cs%diag)
402  if (cs%id_tr_ady(m)>0) &
403  call post_data(cs%id_tr_ady(m),cs%tr_ady(m)%p(:,:,:),cs%diag)
404  if (cs%id_tr_dfx(m)>0) &
405  call post_data(cs%id_tr_dfx(m),cs%tr_dfx(m)%p(:,:,:),cs%diag)
406  if (cs%id_tr_dfy(m)>0) &
407  call post_data(cs%id_tr_dfy(m),cs%tr_dfy(m)%p(:,:,:),cs%diag)
408  enddo
409 
Here is the call graph for this function:

◆ isomip_tracer_end()

subroutine, public isomip_tracer::isomip_tracer_end ( type(isomip_tracer_cs), pointer  CS)

Definition at line 436 of file ISOMIP_tracer.F90.

References ntr.

Referenced by mom_tracer_flow_control::tracer_flow_control_end().

436  type(isomip_tracer_cs), pointer :: cs
437  integer :: m
438 
439  if (associated(cs)) then
440  if (associated(cs%tr)) deallocate(cs%tr)
441  if (associated(cs%tr_aux)) deallocate(cs%tr_aux)
442  do m=1,ntr
443  if (associated(cs%tr_adx(m)%p)) deallocate(cs%tr_adx(m)%p)
444  if (associated(cs%tr_ady(m)%p)) deallocate(cs%tr_ady(m)%p)
445  if (associated(cs%tr_dfx(m)%p)) deallocate(cs%tr_dfx(m)%p)
446  if (associated(cs%tr_dfy(m)%p)) deallocate(cs%tr_dfy(m)%p)
447  enddo
448 
449  deallocate(cs)
450  endif
Here is the caller graph for this function:

◆ isomip_tracer_surface_state()

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

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, in m or kg m-2.
csThe control structure returned by a previous call to ISOMIP_register_tracer.

Definition at line 415 of file ISOMIP_tracer.F90.

References ntr, and coupler_util::set_coupler_values().

415  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
416  type(surface), intent(inout) :: state !< A structure containing fields that describe the surface state of the ocean.
417  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness, in m or kg m-2.
418  type(isomip_tracer_cs), pointer :: cs !< The control structure returned by a previous call to ISOMIP_register_tracer.
419  integer :: i, j, m, is, ie, js, je, nz
420  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
421 
422  if (.not.associated(cs)) return
423 
424  if (cs%coupled_tracers) then
425  do m=1,ntr
426  ! This call loads the surface vlues into the appropriate array in the
427  ! coupler-type structure.
428  call set_coupler_values(cs%tr(:,:,1,1), state%tr_fields, cs%ind_tr(m), &
429  ind_csurf, is, ie, js, je)
430  enddo
431  endif
432 
Here is the call graph for this function:

◆ register_isomip_tracer()

logical function, public isomip_tracer::register_isomip_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(isomip_tracer_cs), pointer  CS,
type(tracer_registry_type), pointer  tr_Reg,
type(mom_restart_cs), pointer  restart_CS 
)

This subroutine is used to register tracer fields.

Parameters
[in]hiA horizontal index type structure.
[in]gvThe ocean's vertical grid structure.
[in]param_fileA structure indicating the open file to parse for model parameter values.
csA pointer that is set to point to the control structure for this module (in/out).
tr_regA pointer to the tracer registry.
restart_csA pointer to the restart control structure.

Definition at line 94 of file ISOMIP_tracer.F90.

References atmos_ocean_fluxes_mod::aof_set_coupler_flux(), mom_error_handler::mom_error(), ntr, and mom_io::var_desc().

94  type(hor_index_type), intent(in) :: hi !<A horizontal index type structure.
95  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
96  type(param_file_type), intent(in) :: param_file !<A structure indicating the open file to parse for model parameter values.
97  type(isomip_tracer_cs), pointer :: cs !<A pointer that is set to point to the control structure for this module (in/out).
98  type(tracer_registry_type), pointer :: tr_reg !<A pointer to the tracer registry.
99  type(mom_restart_cs), pointer :: restart_cs !<A pointer to the restart control structure.
100 
101  character(len=80) :: name, longname
102 ! This include declares and sets the variable "version".
103 #include "version_variable.h"
104  character(len=40) :: mdl = "ISOMIP_tracer" ! This module's name.
105  character(len=200) :: inputdir
106  real, pointer :: tr_ptr(:,:,:) => null()
107  logical :: register_isomip_tracer
108  integer :: isd, ied, jsd, jed, nz, m
109  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
110 
111  if (associated(cs)) then
112  call mom_error(warning, "ISOMIP_register_tracer called with an "// &
113  "associated control structure.")
114  return
115  endif
116  allocate(cs)
117 
118  ! Read all relevant parameters and write them to the model log.
119  call log_version(param_file, mdl, version, "")
120  call get_param(param_file, mdl, "ISOMIP_TRACER_IC_FILE", cs%tracer_IC_file, &
121  "The name of a file from which to read the initial \n"//&
122  "conditions for the ISOMIP tracers, or blank to initialize \n"//&
123  "them internally.", default=" ")
124  if (len_trim(cs%tracer_IC_file) >= 1) then
125  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
126  inputdir = slasher(inputdir)
127  cs%tracer_IC_file = trim(inputdir)//trim(cs%tracer_IC_file)
128  call log_param(param_file, mdl, "INPUTDIR/ISOMIP_TRACER_IC_FILE", &
129  cs%tracer_IC_file)
130  endif
131  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
132  "If true, sponges may be applied anywhere in the domain. \n"//&
133  "The exact location and properties of those sponges are \n"//&
134  "specified from MOM_initialization.F90.", default=.false.)
135 
136  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
137  if (cs%mask_tracers) then
138  allocate(cs%tr_aux(isd:ied,jsd:jed,nz,ntr)) ; cs%tr_aux(:,:,:,:) = 0.0
139  endif
140 
141  do m=1,ntr
142  if (m < 10) then ; write(name,'("tr_D",I1.1)') m
143  else ; write(name,'("tr_D",I2.2)') m ; endif
144  write(longname,'("Concentration of ISOMIP Tracer ",I2.2)') m
145  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
146 
147  ! This is needed to force the compiler not to do a copy in the registration
148  ! calls. Curses on the designers and implementers of Fortran90.
149  tr_ptr => cs%tr(:,:,:,m)
150  ! Register the tracer for the restart file.
151  call register_restart_field(tr_ptr, cs%tr_desc(m), .true., restart_cs)
152  ! Register the tracer for horizontal advection & diffusion.
153  call register_tracer(tr_ptr, cs%tr_desc(m), param_file, hi, gv, tr_reg, &
154  tr_desc_ptr=cs%tr_desc(m))
155 
156  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
157  ! values to the coupler (if any). This is meta-code and its arguments will
158  ! currently (deliberately) give fatal errors if it is used.
159  if (cs%coupled_tracers) &
160  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
161  flux_type=' ', implementation=' ', caller="register_ISOMIP_tracer")
162  enddo
163 
164  cs%tr_Reg => tr_reg
165  register_isomip_tracer = .true.
Here is the call graph for this function:

Variable Documentation

◆ ntr

integer, parameter isomip_tracer::ntr = 1
private