MOM6
dome_tracer Module Reference

Detailed Description

A tracer package that is used as a diagnostic in the DOME experiments.

By Robert Hallberg, 2002

This file contains an example of the code that is needed to set up and use a set (in this case eleven) of dynamically passive tracers. These tracers dye the inflowing water or water initially within a range of latitudes or water initially in a range of depths.

A single subroutine is called from within each file to register each of the tracers for reinitialization and advection and to register the subroutine that initializes the tracers and set up their output and the subroutine that does any tracer physics or chemistry along with diapycnal mixing (included here because some tracers may float or swim vertically or dye diapycnal processes).

Data Types

type  dome_tracer_cs
 The DOME_tracer control structure. More...
 

Functions/Subroutines

logical function, public register_dome_tracer (HI, GV, param_file, CS, tr_Reg, restart_CS)
 Register tracer fields and subroutines to be used with MOM. More...
 
subroutine, public initialize_dome_tracer (restart, day, G, GV, US, h, diag, OBC, CS, sponge_CSp, param_file)
 Initializes the NTR tracer fields in tr(:,:,:,:) and sets up the tracer output. More...
 
subroutine, public dome_tracer_column_physics (h_old, h_new, ea, eb, fluxes, dt, G, GV, US, CS, evap_CFL_limit, minimum_forcing_depth)
 This subroutine applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers from this file. This is a simple example of a set of advected passive tracers. More...
 
subroutine, public dome_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 dome_tracer_end (CS)
 Clean up memory allocations, if any. More...
 

Variables

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

Function/Subroutine Documentation

◆ dome_tracer_column_physics()

subroutine, public dome_tracer::dome_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(unit_scale_type), intent(in)  US,
type(dome_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 tracer physics or chemistry to the tracers from this file. This is a simple example of a set of advected passive tracers.

The arguments to this subroutine are redundant in that h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)

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 DOME_register_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 [H ~> m or kg m-2]

Definition at line 287 of file DOME_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)), &
290  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
291  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
292  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
293  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
294  intent(in) :: ea !< an array to which the amount of fluid entrained
295  !! from the layer above during this call will be
296  !! added [H ~> m or kg m-2].
297  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
298  intent(in) :: eb !< an array to which the amount of fluid entrained
299  !! from the layer below during this call will be
300  !! added [H ~> m or kg m-2].
301  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
302  !! and tracer forcing fields. Unused fields have NULL ptrs.
303  real, intent(in) :: dt !< The amount of time covered by this call [T ~> s]
304  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
305  type(DOME_tracer_CS), pointer :: CS !< The control structure returned by a previous
306  !! call to DOME_register_tracer.
307  real, optional, intent(in) :: evap_CFL_limit !< Limit on the fraction of the water that can
308  !! be fluxed out of the top layer in a timestep [nondim]
309  real, optional, intent(in) :: minimum_forcing_depth !< The smallest depth over which
310  !! fluxes can be applied [H ~> m or kg m-2]
311 
312 ! Local variables
313  real :: b1(SZI_(G)) ! b1 and c1 are variables used by the
314  real :: c1(SZI_(G),SZK_(G)) ! tridiagonal solver.
315  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h_work ! Used so that h can be modified
316  integer :: i, j, k, is, ie, js, je, nz, m
317  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
318 
319  if (.not.associated(cs)) return
320 
321  if (present(evap_cfl_limit) .and. present(minimum_forcing_depth)) then
322  do m=1,ntr
323  do k=1,nz ;do j=js,je ; do i=is,ie
324  h_work(i,j,k) = h_old(i,j,k)
325  enddo ; enddo ; enddo
326  call applytracerboundaryfluxesinout(g, gv, cs%tr(:,:,:,m), dt, fluxes, h_work, &
327  evap_cfl_limit, minimum_forcing_depth)
328  call tracer_vertdiff(h_work, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
329  enddo
330  else
331  do m=1,ntr
332  call tracer_vertdiff(h_old, ea, eb, dt, cs%tr(:,:,:,m), g, gv)
333  enddo
334  endif
335 

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

Here is the call graph for this function:

◆ dome_tracer_end()

subroutine, public dome_tracer::dome_tracer_end ( type(dome_tracer_cs), pointer  CS)

Clean up memory allocations, if any.

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

Definition at line 373 of file DOME_tracer.F90.

373  type(DOME_tracer_CS), pointer :: CS !< The control structure returned by a previous
374  !! call to DOME_register_tracer.
375  integer :: m
376 
377  if (associated(cs)) then
378  if (associated(cs%tr)) deallocate(cs%tr)
379  deallocate(cs)
380  endif

◆ dome_tracer_surface_state()

subroutine, public dome_tracer::dome_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(dome_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 DOME_register_tracer.

Definition at line 342 of file DOME_tracer.F90.

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

References ntr.

◆ initialize_dome_tracer()

subroutine, public dome_tracer::initialize_dome_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,
type(unit_scale_type), intent(in)  US,
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(dome_tracer_cs), pointer  CS,
type(sponge_cs), pointer  sponge_CSp,
type(param_file_type), intent(in)  param_file 
)

Initializes the NTR tracer fields in tr(:,:,:,:) and sets up the tracer output.

Parameters
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in]restart.true. if the fields have already been read from a restart file.
[in]dayTime of the start of the run.
[in]hLayer thicknesses [H ~> m or kg m-2]
[in]diagStructure used to regulate diagnostic output.
obcStructure specifying open boundary options.
csThe control structure returned by a previous call to DOME_register_tracer.
sponge_cspA pointer to the control structure for the sponges, if they are in use.
[in]param_fileA structure to parse for run-time parameters

Definition at line 144 of file DOME_tracer.F90.

144  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
145  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
146  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
147  logical, intent(in) :: restart !< .true. if the fields have already
148  !! been read from a restart file.
149  type(time_type), target, intent(in) :: day !< Time of the start of the run.
150  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
151  type(diag_ctrl), target, intent(in) :: diag !< Structure used to regulate diagnostic output.
152  type(ocean_OBC_type), pointer :: OBC !< Structure specifying open boundary options.
153  type(DOME_tracer_CS), pointer :: CS !< The control structure returned by a previous
154  !! call to DOME_register_tracer.
155  type(sponge_CS), pointer :: sponge_CSp !< A pointer to the control structure
156  !! for the sponges, if they are in use.
157  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
158 
159 ! Local variables
160  real, allocatable :: temp(:,:,:)
161  real, pointer, dimension(:,:,:) :: &
162  OBC_tr1_u => null(), & ! These arrays should be allocated and set to
163  obc_tr1_v => null() ! specify the values of tracer 1 that should come
164  ! in through u- and v- points through the open
165  ! boundary conditions, in the same units as tr.
166  character(len=16) :: name ! A variable's name in a NetCDF file.
167  character(len=72) :: longname ! The long name of that variable.
168  character(len=48) :: units ! The dimensions of the variable.
169  character(len=48) :: flux_units ! The units for tracer fluxes, usually
170  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
171  real, pointer :: tr_ptr(:,:,:) => null()
172  real :: PI ! 3.1415926... calculated as 4*atan(1)
173  real :: tr_y ! Initial zonally uniform tracer concentrations.
174  real :: h_neglect ! A thickness that is so small it is usually lost
175  ! in roundoff and can be neglected [H ~> m or kg m-2].
176  real :: e(SZK_(G)+1), e_top, e_bot ! Heights [Z ~> m].
177  real :: d_tr ! A change in tracer concentraions, in tracer units.
178  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
179  integer :: IsdB, IedB, JsdB, JedB
180 
181  if (.not.associated(cs)) 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  h_neglect = gv%H_subroundoff
186 
187  cs%Time => day
188  cs%diag => diag
189 
190  if (.not.restart) then
191  if (len_trim(cs%tracer_IC_file) >= 1) then
192  ! Read the tracer concentrations from a netcdf file.
193  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
194  call mom_error(fatal, "DOME_initialize_tracer: Unable to open "// &
195  cs%tracer_IC_file)
196  do m=1,ntr
197  call query_vardesc(cs%tr_desc(m), name, caller="initialize_DOME_tracer")
198  call mom_read_data(cs%tracer_IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
199  enddo
200  else
201  do m=1,ntr
202  do k=1,nz ; do j=js,je ; do i=is,ie
203  cs%tr(i,j,k,m) = 1.0e-20 ! This could just as well be 0.
204  enddo ; enddo ; enddo
205  enddo
206 
207 ! This sets a stripe of tracer across the basin.
208  do m=2,ntr ; do j=js,je ; do i=is,ie
209  tr_y = 0.0
210  if ((m <= 6) .and. (g%geoLatT(i,j) > (300.0+50.0*real(m-1))) .and. &
211  (g%geoLatT(i,j) < (350.0+50.0*real(m-1)))) tr_y = 1.0
212  do k=1,nz
213 ! This adds the stripes of tracer to every layer.
214  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + tr_y
215  enddo
216  enddo ; enddo ; enddo
217 
218  if (ntr > 7) then
219  do j=js,je ; do i=is,ie
220  e(nz+1) = -g%bathyT(i,j)
221  do k=nz,1,-1
222  e(k) = e(k+1) + h(i,j,k)*gv%H_to_Z
223  do m=7,ntr
224  e_top = (-600.0*real(m-1) + 3000.0) * us%m_to_Z
225  e_bot = (-600.0*real(m-1) + 2700.0) * us%m_to_Z
226  if (e_top < e(k)) then
227  if (e_top < e(k+1)) then ; d_tr = 0.0
228  elseif (e_bot < e(k+1)) then
229  d_tr = 1.0 * (e_top-e(k+1)) / ((h(i,j,k)+h_neglect)*gv%H_to_Z)
230  else ; d_tr = 1.0 * (e_top-e_bot) / ((h(i,j,k)+h_neglect)*gv%H_to_Z)
231  endif
232  elseif (e_bot < e(k)) then
233  if (e_bot < e(k+1)) then ; d_tr = 1.0
234  else ; d_tr = 1.0 * (e(k)-e_bot) / ((h(i,j,k)+h_neglect)*gv%H_to_Z)
235  endif
236  else
237  d_tr = 0.0
238  endif
239  if (h(i,j,k) < 2.0*gv%Angstrom_H) d_tr=0.0
240  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + d_tr
241  enddo
242  enddo
243  enddo ; enddo
244  endif
245 
246  endif
247  endif ! restart
248 
249  if ( cs%use_sponge ) then
250 ! If sponges are used, this example damps tracers in sponges in the
251 ! northern half of the domain to 1 and tracers in the southern half
252 ! to 0. For any tracers that are not damped in the sponge, the call
253 ! to set_up_sponge_field can simply be omitted.
254  if (.not.associated(sponge_csp)) &
255  call mom_error(fatal, "DOME_initialize_tracer: "// &
256  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
257 
258  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
259  do k=1,nz ; do j=js,je ; do i=is,ie
260  if (g%geoLatT(i,j) > 700.0 .and. (k > nz/2)) then
261  temp(i,j,k) = 1.0
262  else
263  temp(i,j,k) = 0.0
264  endif
265  enddo ; enddo ; enddo
266 
267 ! do m=1,NTR
268  do m=1,1
269  ! This is needed to force the compiler not to do a copy in the sponge
270  ! calls. Curses on the designers and implementers of Fortran90.
271  tr_ptr => cs%tr(:,:,:,m)
272  call set_up_sponge_field(temp, tr_ptr, g, nz, sponge_csp)
273  enddo
274  deallocate(temp)
275  endif
276 

References mom_error_handler::mom_error(), ntr, mom_io::query_vardesc(), and mom_sponge::set_up_sponge_field().

Here is the call graph for this function:

◆ register_dome_tracer()

logical function, public dome_tracer::register_dome_tracer ( type(hor_index_type), intent(in)  HI,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(dome_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
csA pointer that is set to point to the control structure for this module
tr_regA pointer to the tracer registry.
restart_csA pointer to the restart control structure.

Definition at line 64 of file DOME_tracer.F90.

64  type(hor_index_type), intent(in) :: HI !< A horizontal index type structure.
65  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
66  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
67  type(DOME_tracer_CS), pointer :: CS !< A pointer that is set to point to the
68  !! control structure for this module
69  type(tracer_registry_type), pointer :: tr_Reg !< A pointer to the tracer registry.
70  type(MOM_restart_CS), pointer :: restart_CS !< A pointer to the restart control structure.
71 
72 ! Local variables
73  character(len=80) :: name, longname
74 ! This include declares and sets the variable "version".
75 #include "version_variable.h"
76  character(len=40) :: mdl = "DOME_tracer" ! This module's name.
77  character(len=48) :: flux_units ! The units for tracer fluxes, usually
78  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
79  character(len=200) :: inputdir
80  real, pointer :: tr_ptr(:,:,:) => null()
81  logical :: register_DOME_tracer
82  integer :: isd, ied, jsd, jed, nz, m
83  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
84 
85  if (associated(cs)) then
86  call mom_error(warning, "DOME_register_tracer called with an "// &
87  "associated control structure.")
88  return
89  endif
90  allocate(cs)
91 
92  ! Read all relevant parameters and write them to the model log.
93  call log_version(param_file, mdl, version, "")
94  call get_param(param_file, mdl, "DOME_TRACER_IC_FILE", cs%tracer_IC_file, &
95  "The name of a file from which to read the initial "//&
96  "conditions for the DOME tracers, or blank to initialize "//&
97  "them internally.", default=" ")
98  if (len_trim(cs%tracer_IC_file) >= 1) then
99  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
100  inputdir = slasher(inputdir)
101  cs%tracer_IC_file = trim(inputdir)//trim(cs%tracer_IC_file)
102  call log_param(param_file, mdl, "INPUTDIR/DOME_TRACER_IC_FILE", &
103  cs%tracer_IC_file)
104  endif
105  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
106  "If true, sponges may be applied anywhere in the domain. "//&
107  "The exact location and properties of those sponges are "//&
108  "specified from MOM_initialization.F90.", default=.false.)
109 
110  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
111 
112  do m=1,ntr
113  if (m < 10) then ; write(name,'("tr_D",I1.1)') m
114  else ; write(name,'("tr_D",I2.2)') m ; endif
115  write(longname,'("Concentration of DOME Tracer ",I2.2)') m
116  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
117  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
118  else ; flux_units = "kg s-1" ; endif
119 
120  ! This is needed to force the compiler not to do a copy in the registration
121  ! calls. Curses on the designers and implementers of Fortran90.
122  tr_ptr => cs%tr(:,:,:,m)
123  ! Register the tracer for horizontal advection, diffusion, and restarts.
124  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
125  name=name, longname=longname, units="kg kg-1", &
126  registry_diags=.true., restart_cs=restart_cs, &
127  flux_units=trim(flux_units), flux_scale=gv%H_to_MKS)
128 
129  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
130  ! values to the coupler (if any). This is meta-code and its arguments will
131  ! currently (deliberately) give fatal errors if it is used.
132  if (cs%coupled_tracers) &
133  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
134  flux_type=' ', implementation=' ', caller="register_DOME_tracer")
135  enddo
136 
137  cs%tr_Reg => tr_reg
138  register_dome_tracer = .true.

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

Here is the call graph for this function:

Variable Documentation

◆ ntr

integer, parameter dome_tracer::ntr = 11
private

The number of tracers in this module.

Definition at line 39 of file DOME_tracer.F90.

39 integer, parameter :: ntr = 11 !< The number of tracers in this module.

Referenced by dome_tracer_column_physics(), dome_tracer_surface_state(), initialize_dome_tracer(), and register_dome_tracer().