MOM6
tracer_example.F90
1 !> A sample tracer package that has striped initial conditions
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
7 use mom_error_handler, only : mom_error, fatal, warning
9 use mom_forcing_type, only : forcing
10 use mom_grid, only : ocean_grid_type
11 use mom_hor_index, only : hor_index_type
12 use mom_io, only : file_exists, mom_read_data, slasher, vardesc, var_desc, query_vardesc
14 use mom_restart, only : mom_restart_cs
15 use mom_sponge, only : set_up_sponge_field, sponge_cs
16 use mom_time_manager, only : time_type
17 use mom_tracer_registry, only : register_tracer, tracer_registry_type
18 use mom_variables, only : surface
20 
21 use coupler_types_mod, only : coupler_type_set_data, ind_csurf
22 use atmos_ocean_fluxes_mod, only : aof_set_coupler_flux
23 
24 implicit none ; private
25 
26 #include <MOM_memory.h>
27 
28 public user_register_tracer_example, user_initialize_tracer, user_tracer_stock
29 public tracer_column_physics, user_tracer_surface_state, user_tracer_example_end
30 
31 integer, parameter :: ntr = 1 !< The number of tracers in this module.
32 
33 !> The control structure for the USER_tracer_example module
34 type, public :: user_tracer_example_cs ; private
35  logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
36  character(len=200) :: tracer_ic_file !< The full path to the IC file, or " "
37  !! to initialize internally.
38  type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
39  type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the tracer registry
40  real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this subroutine, in g m-3?
41  real :: land_val(ntr) = -1.0 !< The value of tr that is used where land is masked out.
42  logical :: use_sponge !< If true, sponges may be applied somewhere in the domain.
43 
44  integer, dimension(NTR) :: ind_tr !< Indices returned by aof_set_coupler_flux if it is used and the
45  !! surface tracer concentrations are to be provided to the coupler.
46 
47  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to regulate the timing of diagnostic output.
48 
49  type(vardesc) :: tr_desc(ntr) !< Descriptions of each of the tracers.
51 
52 contains
53 
54 !> This subroutine is used to register tracer fields and subroutines to be used with MOM.
55 function user_register_tracer_example(HI, GV, param_file, CS, tr_Reg, restart_CS)
56  type(hor_index_type), intent(in) :: hi !< A horizontal index type structure
57  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
58  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
59  type(user_tracer_example_cs), pointer :: cs !< A pointer that is set to point to the control
60  !! structure for this module
61  type(tracer_registry_type), pointer :: tr_reg !< A pointer that is set to point to the control
62  !! structure for the tracer advection and
63  !! diffusion module
64  type(mom_restart_cs), pointer :: restart_cs !< A pointer to the restart control structure
65 
66 ! Local variables
67  character(len=80) :: name, longname
68 ! This include declares and sets the variable "version".
69 #include "version_variable.h"
70  character(len=40) :: mdl = "tracer_example" ! This module's name.
71  character(len=200) :: inputdir
72  character(len=48) :: flux_units ! The units for tracer fluxes, usually
73  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
74  real, pointer :: tr_ptr(:,:,:) => null()
75  logical :: user_register_tracer_example
76  integer :: isd, ied, jsd, jed, nz, m
77  isd = hi%isd ; ied = hi%ied ; jsd = hi%jsd ; jed = hi%jed ; nz = gv%ke
78 
79  if (associated(cs)) then
80  call mom_error(warning, "USER_register_tracer_example called with an "// &
81  "associated control structure.")
82  return
83  endif
84  allocate(cs)
85 
86  ! Read all relevant parameters and write them to the model log.
87  call log_version(param_file, mdl, version, "")
88  call get_param(param_file, mdl, "TRACER_EXAMPLE_IC_FILE", cs%tracer_IC_file, &
89  "The name of a file from which to read the initial "//&
90  "conditions for the DOME tracers, or blank to initialize "//&
91  "them internally.", default=" ")
92  if (len_trim(cs%tracer_IC_file) >= 1) then
93  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
94  cs%tracer_IC_file = trim(slasher(inputdir))//trim(cs%tracer_IC_file)
95  call log_param(param_file, mdl, "INPUTDIR/TRACER_EXAMPLE_IC_FILE", &
96  cs%tracer_IC_file)
97  endif
98  call get_param(param_file, mdl, "SPONGE", cs%use_sponge, &
99  "If true, sponges may be applied anywhere in the domain. "//&
100  "The exact location and properties of those sponges are "//&
101  "specified from MOM_initialization.F90.", default=.false.)
102 
103  allocate(cs%tr(isd:ied,jsd:jed,nz,ntr)) ; cs%tr(:,:,:,:) = 0.0
104 
105  do m=1,ntr
106  if (m < 10) then ; write(name,'("tr",I1.1)') m
107  else ; write(name,'("tr",I2.2)') m ; endif
108  write(longname,'("Concentration of Tracer ",I2.2)') m
109  cs%tr_desc(m) = var_desc(name, units="kg kg-1", longname=longname, caller=mdl)
110 
111  ! This needs to be changed if the units of tracer are changed above.
112  if (gv%Boussinesq) then ; flux_units = "kg kg-1 m3 s-1"
113  else ; flux_units = "kg s-1" ; endif
114 
115  ! This is needed to force the compiler not to do a copy in the registration
116  ! calls. Curses on the designers and implementers of Fortran90.
117  tr_ptr => cs%tr(:,:,:,m)
118  ! Register the tracer for horizontal advection, diffusion, and restarts.
119  call register_tracer(tr_ptr, tr_reg, param_file, hi, gv, &
120  name=name, longname=longname, units="kg kg-1", &
121  registry_diags=.true., flux_units=flux_units, &
122  restart_cs=restart_cs)
123 
124  ! Set coupled_tracers to be true (hard-coded above) to provide the surface
125  ! values to the coupler (if any). This is meta-code and its arguments will
126  ! currently (deliberately) give fatal errors if it is used.
127  if (cs%coupled_tracers) &
128  cs%ind_tr(m) = aof_set_coupler_flux(trim(name)//'_flux', &
129  flux_type=' ', implementation=' ', caller="USER_register_tracer_example")
130  enddo
131 
132  cs%tr_Reg => tr_reg
133  user_register_tracer_example = .true.
134 end function user_register_tracer_example
135 
136 !> This subroutine initializes the NTR tracer fields in tr(:,:,:,:)
137 !! and it sets up the tracer output.
138 subroutine user_initialize_tracer(restart, day, G, GV, h, diag, OBC, CS, &
139  sponge_CSp)
140  logical, intent(in) :: restart !< .true. if the fields have already
141  !! been read from a restart file.
142  type(time_type), target, intent(in) :: day !< Time of the start of the run.
143  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
144  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
145  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
146  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
147  type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate
148  !! diagnostic output.
149  type(ocean_obc_type), pointer :: obc !< This open boundary condition type specifies
150  !! whether, where, and what open boundary
151  !! conditions are used.
152  type(user_tracer_example_cs), pointer :: cs !< The control structure returned by a previous
153  !! call to USER_register_tracer_example.
154  type(sponge_cs), pointer :: sponge_csp !< A pointer to the control structure
155  !! for the sponges, if they are in use.
156 
157 ! Local variables
158  real, allocatable :: temp(:,:,:)
159  character(len=32) :: name ! A variable's name in a NetCDF file.
160  character(len=72) :: longname ! The long name of that variable.
161  character(len=48) :: units ! The dimensions of the variable.
162  character(len=48) :: flux_units ! The units for tracer fluxes, usually
163  ! kg(tracer) kg(water)-1 m3 s-1 or kg(tracer) s-1.
164  real, pointer :: tr_ptr(:,:,:) => null()
165  real :: pi ! 3.1415926... calculated as 4*atan(1)
166  real :: tr_y ! Initial zonally uniform tracer concentrations.
167  real :: dist2 ! The distance squared from a line [m2].
168  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
169  integer :: isdb, iedb, jsdb, jedb, lntr
170 
171  if (.not.associated(cs)) return
172  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
173  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
174  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
175 
176  lntr = ntr ! Avoids compile-time warning when NTR<2
177  cs%Time => day
178  cs%diag => diag
179 
180  if (.not.restart) then
181  if (len_trim(cs%tracer_IC_file) >= 1) then
182 ! Read the tracer concentrations from a netcdf file.
183  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
184  call mom_error(fatal, "USER_initialize_tracer: Unable to open "// &
185  cs%tracer_IC_file)
186  do m=1,ntr
187  call query_vardesc(cs%tr_desc(m), name, caller="USER_initialize_tracer")
188  call mom_read_data(cs%tracer_IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
189  enddo
190  else
191  do m=1,ntr
192  do k=1,nz ; do j=js,je ; do i=is,ie
193  cs%tr(i,j,k,m) = 1.0e-20 ! This could just as well be 0.
194  enddo ; enddo ; enddo
195  enddo
196 
197 ! This sets a stripe of tracer across the basin.
198  pi = 4.0*atan(1.0)
199  do j=js,je
200  dist2 = (g%Rad_Earth * pi / 180.0)**2 * &
201  (g%geoLatT(i,j) - 40.0) * (g%geoLatT(i,j) - 40.0)
202  tr_y = 0.5*exp(-dist2/(1.0e5*1.0e5))
203 
204  do k=1,nz ; do i=is,ie
205 ! This adds the stripes of tracer to every layer.
206  cs%tr(i,j,k,1) = cs%tr(i,j,k,1) + tr_y
207  enddo ; enddo
208  enddo
209  endif
210  endif ! restart
211 
212  if ( cs%use_sponge ) then
213 ! If sponges are used, this example damps tracers in sponges in the
214 ! northern half of the domain to 1 and tracers in the southern half
215 ! to 0. For any tracers that are not damped in the sponge, the call
216 ! to set_up_sponge_field can simply be omitted.
217  if (.not.associated(sponge_csp)) &
218  call mom_error(fatal, "USER_initialize_tracer: "// &
219  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
220 
221  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
222  do k=1,nz ; do j=js,je ; do i=is,ie
223  if (g%geoLatT(i,j) > 700.0 .and. (k > nz/2)) then
224  temp(i,j,k) = 1.0
225  else
226  temp(i,j,k) = 0.0
227  endif
228  enddo ; enddo ; enddo
229 
230 ! do m=1,NTR
231  do m=1,1
232  ! This is needed to force the compiler not to do a copy in the sponge
233  ! calls. Curses on the designers and implementers of Fortran90.
234  tr_ptr => cs%tr(:,:,:,m)
235  call set_up_sponge_field(temp, tr_ptr, g, nz, sponge_csp)
236  enddo
237  deallocate(temp)
238  endif
239 
240  if (associated(obc)) then
241  call query_vardesc(cs%tr_desc(1), name, caller="USER_initialize_tracer")
242  if (obc%specified_v_BCs_exist_globally) then
243  ! Steal from updated DOME in the fullness of time.
244  else
245  ! Steal from updated DOME in the fullness of time.
246  endif
247  ! All tracers but the first have 0 concentration in their inflows. As this
248  ! is the default value, the following calls are unnecessary.
249  do m=2,lntr
250  call query_vardesc(cs%tr_desc(m), name, caller="USER_initialize_tracer")
251  ! Steal from updated DOME in the fullness of time.
252  enddo
253  endif
254 
255 end subroutine user_initialize_tracer
256 
257 !> This subroutine applies diapycnal diffusion and any other column
258 !! tracer physics or chemistry to the tracers from this file.
259 !! This is a simple example of a set of advected passive tracers.
260 !! The arguments to this subroutine are redundant in that
261 !! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
262 subroutine tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS)
263  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
264  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
265  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
266  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
267  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
268  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
269  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
270  intent(in) :: ea !< an array to which the amount of fluid entrained
271  !! from the layer above during this call will be
272  !! added [H ~> m or kg m-2].
273  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
274  intent(in) :: eb !< an array to which the amount of fluid entrained
275  !! from the layer below during this call will be
276  !! added [H ~> m or kg m-2].
277  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
278  !! and tracer forcing fields. Unused fields have NULL ptrs.
279  real, intent(in) :: dt !< The amount of time covered by this call [s]
280  type(user_tracer_example_cs), pointer :: cs !< The control structure returned by a previous
281  !! call to USER_register_tracer_example.
282 
283 ! Local variables
284  real :: hold0(szi_(g)) ! The original topmost layer thickness,
285  ! with surface mass fluxes added back, m.
286  real :: b1(szi_(g)) ! b1 and c1 are variables used by the
287  real :: c1(szi_(g),szk_(g)) ! tridiagonal solver.
288  real :: d1(szi_(g)) ! d1=1-c1 is used by the tridiagonal solver.
289  real :: h_neglect ! A thickness that is so small it is usually lost
290  ! in roundoff and can be neglected [H ~> m or kg m-2].
291  real :: b_denom_1 ! The first term in the denominator of b1 [H ~> m or kg m-2].
292  integer :: i, j, k, is, ie, js, je, nz, m
293 
294 ! The following array (trdc) determines the behavior of the tracer
295 ! diapycnal advection. The first element is 1 if tracers are
296 ! passively advected. The second and third are the concentrations
297 ! to which downwelling and upwelling water are set, respectively.
298 ! For most (normal) tracers, the appropriate vales are {1,0,0}.
299 
300  real :: trdc(3)
301 ! Uncomment the following line to dye both upwelling and downwelling.
302 ! data trdc / 0.0,1.0,1.0 /
303 ! Uncomment the following line to dye downwelling.
304 ! data trdc / 0.0,1.0,0.0 /
305 ! Uncomment the following line to dye upwelling.
306 ! data trdc / 0.0,0.0,1.0 /
307 ! Uncomment the following line for tracer concentrations to be set
308 ! to zero in any diapycnal motions.
309 ! data trdc / 0.0,0.0,0.0 /
310 ! Uncomment the following line for most "physical" tracers, which
311 ! are advected diapycnally in the usual manner.
312  data trdc / 1.0,0.0,0.0 /
313  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
314 
315  if (.not.associated(cs)) return
316  h_neglect = gv%H_subroundoff
317 
318  do j=js,je
319  do i=is,ie
320 ! The following line is appropriate for quantities like salinity
321 ! that are left behind by evaporation, and any surface fluxes would
322 ! be explicitly included in the flux structure.
323  hold0(i) = h_old(i,j,1)
324 ! The following line is appropriate for quantities like temperature
325 ! that can be assumed to have the same concentration in evaporation
326 ! as they had in the water. The explicit surface fluxes here would
327 ! reflect differences in concentration from the ambient water, not
328 ! the absolute fluxes.
329  ! hold0(i) = h_old(i,j,1) + ea(i,j,1)
330  b_denom_1 = h_old(i,j,1) + ea(i,j,1) + h_neglect
331  b1(i) = 1.0 / (b_denom_1 + eb(i,j,1))
332 ! d1(i) = b_denom_1 * b1(i)
333  d1(i) = trdc(1) * (b_denom_1 * b1(i)) + (1.0 - trdc(1))
334  do m=1,ntr
335  cs%tr(i,j,1,m) = b1(i)*(hold0(i)*cs%tr(i,j,1,m) + trdc(3)*eb(i,j,1))
336  ! Add any surface tracer fluxes to the preceding line.
337  enddo
338  enddo
339  do k=2,nz ; do i=is,ie
340  c1(i,k) = trdc(1) * eb(i,j,k-1) * b1(i)
341  b_denom_1 = h_old(i,j,k) + d1(i)*ea(i,j,k) + h_neglect
342  b1(i) = 1.0 / (b_denom_1 + eb(i,j,k))
343  d1(i) = trdc(1) * (b_denom_1 * b1(i)) + (1.0 - trdc(1))
344  do m=1,ntr
345  cs%tr(i,j,k,m) = b1(i) * (h_old(i,j,k)*cs%tr(i,j,k,m) + &
346  ea(i,j,k)*(trdc(1)*cs%tr(i,j,k-1,m)+trdc(2)) + &
347  eb(i,j,k)*trdc(3))
348  enddo
349  enddo ; enddo
350  do m=1,ntr ; do k=nz-1,1,-1 ; do i=is,ie
351  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + c1(i,k+1)*cs%tr(i,j,k+1,m)
352  enddo ; enddo ; enddo
353  enddo
354 
355 end subroutine tracer_column_physics
356 
357 !> This function calculates the mass-weighted integral of all tracer stocks,
358 !! returning the number of stocks it has calculated. If the stock_index
359 !! is present, only the stock corresponding to that coded index is returned.
360 function user_tracer_stock(h, stocks, G, GV, CS, names, units, stock_index)
361  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
362  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
363  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
364  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
365  real, dimension(:), intent(out) :: stocks !< the mass-weighted integrated amount of each
366  !! tracer, in kg times concentration units [kg conc].
367  type(user_tracer_example_cs), pointer :: cs !< The control structure returned by a
368  !! previous call to register_USER_tracer.
369  character(len=*), dimension(:), intent(out) :: names !< The names of the stocks calculated.
370  character(len=*), dimension(:), intent(out) :: units !< The units of the stocks calculated.
371  integer, optional, intent(in) :: stock_index !< The coded index of a specific stock
372  !! being sought.
373  integer :: user_tracer_stock !< Return value: the number of
374  !! stocks calculated here.
375 
376 ! Local variables
377  integer :: i, j, k, is, ie, js, je, nz, m
378  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
379 
380  user_tracer_stock = 0
381  if (.not.associated(cs)) return
382 
383  if (present(stock_index)) then ; if (stock_index > 0) then
384  ! Check whether this stock is available from this routine.
385 
386  ! No stocks from this routine are being checked yet. Return 0.
387  return
388  endif ; endif
389 
390  do m=1,ntr
391  call query_vardesc(cs%tr_desc(m), name=names(m), units=units(m), caller="USER_tracer_stock")
392  units(m) = trim(units(m))//" kg"
393  stocks(m) = 0.0
394  do k=1,nz ; do j=js,je ; do i=is,ie
395  stocks(m) = stocks(m) + cs%tr(i,j,k,m) * &
396  (g%mask2dT(i,j) * g%areaT(i,j) * h(i,j,k))
397  enddo ; enddo ; enddo
398  stocks(m) = gv%H_to_kg_m2 * stocks(m)
399  enddo
400  user_tracer_stock = ntr
401 
402 end function user_tracer_stock
403 
404 !> This subroutine extracts the surface fields from this tracer package that
405 !! are to be shared with the atmosphere in coupled configurations.
406 subroutine user_tracer_surface_state(state, h, G, CS)
407  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
408  type(surface), intent(inout) :: state !< A structure containing fields that
409  !! describe the surface state of the ocean.
410  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
411  intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
412  type(user_tracer_example_cs), pointer :: cs !< The control structure returned by a previous
413  !! call to register_USER_tracer.
414 
415  ! This particular tracer package does not report anything back to the coupler.
416  ! The code that is here is just a rough guide for packages that would.
417 
418  integer :: m, is, ie, js, je, isd, ied, jsd, jed
419  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
420  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
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 values into the appropriate array in the
427  ! coupler-type structure.
428  call coupler_type_set_data(cs%tr(:,:,1,m), cs%ind_tr(m), ind_csurf, &
429  state%tr_fields, idim=(/isd, is, ie, ied/), &
430  jdim=(/jsd, js, je, jed/) )
431  enddo
432  endif
433 
434 end subroutine user_tracer_surface_state
435 
436 !> Clean up allocated memory at the end.
437 subroutine user_tracer_example_end(CS)
438  type(user_tracer_example_cs), pointer :: cs !< The control structure returned by a previous
439  !! call to register_USER_tracer.
440  integer :: m
441 
442  if (associated(cs)) then
443  if (associated(cs%tr)) deallocate(cs%tr)
444  deallocate(cs)
445  endif
446 end subroutine user_tracer_example_end
447 
448 !> \namespace user_tracer_example
449 !!
450 !! Original by Robert Hallberg, 2002
451 !!
452 !! This file contains an example of the code that is needed to set
453 !! up and use a set (in this case one) of dynamically passive tracers.
454 !!
455 !! A single subroutine is called from within each file to register
456 !! each of the tracers for reinitialization and advection and to
457 !! register the subroutine that initializes the tracers and set up
458 !! their output and the subroutine that does any tracer physics or
459 !! chemistry along with diapycnal mixing (included here because some
460 !! tracers may float or swim vertically or dye diapycnal processes).
461 end module user_tracer_example
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
atmos_ocean_fluxes_mod
A dummy version of atmos_ocean_fluxes_mod module for use when the vastly larger FMS package is not ne...
Definition: atmos_ocean_fluxes.F90:3
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
mom_tracer_registry
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Definition: MOM_tracer_registry.F90:5
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_hor_index
Defines the horizontal index type (hor_index_type) used for providing index ranges.
Definition: MOM_hor_index.F90:2
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_restart::mom_restart_cs
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:72
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_restart
The MOM6 facility for reading and writing restart files, and querying what has been read.
Definition: MOM_restart.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_open_boundary
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.F90:2
user_tracer_example::user_tracer_example_cs
The control structure for the USER_tracer_example module.
Definition: tracer_example.F90:34
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_sponge
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
user_tracer_example
A sample tracer package that has striped initial conditions.
Definition: tracer_example.F90:2
mom_tracer_registry::tracer_registry_type
Type to carry basic tracer information.
Definition: MOM_tracer_registry.F90:122
mom_hor_index::hor_index_type
Container for horizontal index ranges for data, computational and global domains.
Definition: MOM_hor_index.F90:15
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_open_boundary::ocean_obc_type
Open-boundary data.
Definition: MOM_open_boundary.F90:186
mom_sponge::sponge_cs
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:40
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:49
mom_io::vardesc
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239