MOM6
DOME_tracer.F90
1 !> A tracer package that is used as a diagnostic in the DOME experiments
2 module dome_tracer
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_hor_index, only : hor_index_type
11 use mom_grid, only : ocean_grid_type
12 use mom_io, only : file_exists, mom_read_data, slasher, vardesc, var_desc, query_vardesc
14 use mom_open_boundary, only : obc_segment_type, register_segment_tracer
15 use mom_restart, only : mom_restart_cs
16 use mom_sponge, only : set_up_sponge_field, sponge_cs
17 use mom_time_manager, only : time_type
18 use mom_tracer_registry, only : register_tracer, tracer_registry_type
19 use mom_tracer_diabatic, only : tracer_vertdiff, applytracerboundaryfluxesinout
21 use mom_variables, only : surface
23 
24 use coupler_types_mod, only : coupler_type_set_data, ind_csurf
25 use atmos_ocean_fluxes_mod, only : aof_set_coupler_flux
26 
27 implicit none ; private
28 
29 #include <MOM_memory.h>
30 
31 public register_dome_tracer, initialize_dome_tracer
32 public dome_tracer_column_physics, dome_tracer_surface_state, dome_tracer_end
33 
34 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
35 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
36 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
37 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
38 
39 integer, parameter :: ntr = 11 !< The number of tracers in this module.
40 
41 !> The DOME_tracer control structure
42 type, public :: dome_tracer_cs ; private
43  logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
44  character(len=200) :: tracer_ic_file !< The full path to the IC file, or " " to initialize internally.
45  type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
46  type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the tracer registry
47  real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this package, in g m-3?
48  real :: land_val(ntr) = -1.0 !< The value of tr used where land is masked out.
49  logical :: use_sponge !< If true, sponges may be applied somewhere in the domain.
50 
51  integer, dimension(NTR) :: ind_tr !< Indices returned by aof_set_coupler_flux if it is used and the
52  !! surface tracer concentrations are to be provided to the coupler.
53 
54  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
55  !! regulate the timing of diagnostic output.
56 
57  type(vardesc) :: tr_desc(ntr) !< Descriptions and metadata for the tracers
58 end type dome_tracer_cs
59 
60 contains
61 
62 !> Register tracer fields and subroutines to be used with MOM.
63 function register_dome_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
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., flux_units=flux_units, &
127  restart_cs=restart_cs)
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.
139 end function register_dome_tracer
140 
141 !> Initializes the NTR tracer fields in tr(:,:,:,:) and sets up the tracer output.
142 subroutine initialize_dome_tracer(restart, day, G, GV, US, h, diag, OBC, CS, &
143  sponge_CSp, param_file)
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 :: dist2 ! The distance squared from a line [m2].
175  real :: h_neglect ! A thickness that is so small it is usually lost
176  ! in roundoff and can be neglected [H ~> m or kg m-2].
177  real :: e(szk_(g)+1), e_top, e_bot ! Heights [Z ~> m].
178  real :: d_tr ! A change in tracer concentraions, in tracer units.
179  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz, m
180  integer :: isdb, iedb, jsdb, jedb
181 
182  if (.not.associated(cs)) return
183  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
184  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
185  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
186  h_neglect = gv%H_subroundoff
187 
188  cs%Time => day
189  cs%diag => diag
190 
191  if (.not.restart) then
192  if (len_trim(cs%tracer_IC_file) >= 1) then
193  ! Read the tracer concentrations from a netcdf file.
194  if (.not.file_exists(cs%tracer_IC_file, g%Domain)) &
195  call mom_error(fatal, "DOME_initialize_tracer: Unable to open "// &
196  cs%tracer_IC_file)
197  do m=1,ntr
198  call query_vardesc(cs%tr_desc(m), name, caller="initialize_DOME_tracer")
199  call mom_read_data(cs%tracer_IC_file, trim(name), cs%tr(:,:,:,m), g%Domain)
200  enddo
201  else
202  do m=1,ntr
203  do k=1,nz ; do j=js,je ; do i=is,ie
204  cs%tr(i,j,k,m) = 1.0e-20 ! This could just as well be 0.
205  enddo ; enddo ; enddo
206  enddo
207 
208 ! This sets a stripe of tracer across the basin.
209  do m=2,ntr ; do j=js,je ; do i=is,ie
210  tr_y = 0.0
211  if ((m <= 6) .and. (g%geoLatT(i,j) > (300.0+50.0*real(m-1))) .and. &
212  (g%geoLatT(i,j) < (350.0+50.0*real(m-1)))) tr_y = 1.0
213  do k=1,nz
214 ! This adds the stripes of tracer to every layer.
215  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + tr_y
216  enddo
217  enddo ; enddo ; enddo
218 
219  if (ntr > 7) then
220  do j=js,je ; do i=is,ie
221  e(nz+1) = -g%bathyT(i,j)
222  do k=nz,1,-1
223  e(k) = e(k+1) + h(i,j,k)*gv%H_to_Z
224  do m=7,ntr
225  e_top = (-600.0*real(m-1) + 3000.0) * us%m_to_Z
226  e_bot = (-600.0*real(m-1) + 2700.0) * us%m_to_Z
227  if (e_top < e(k)) then
228  if (e_top < e(k+1)) then ; d_tr = 0.0
229  elseif (e_bot < e(k+1)) then
230  d_tr = 1.0 * (e_top-e(k+1)) / ((h(i,j,k)+h_neglect)*gv%H_to_Z)
231  else ; d_tr = 1.0 * (e_top-e_bot) / ((h(i,j,k)+h_neglect)*gv%H_to_Z)
232  endif
233  elseif (e_bot < e(k)) then
234  if (e_bot < e(k+1)) then ; d_tr = 1.0
235  else ; d_tr = 1.0 * (e(k)-e_bot) / ((h(i,j,k)+h_neglect)*gv%H_to_Z)
236  endif
237  else
238  d_tr = 0.0
239  endif
240  if (h(i,j,k) < 2.0*gv%Angstrom_H) d_tr=0.0
241  cs%tr(i,j,k,m) = cs%tr(i,j,k,m) + d_tr
242  enddo
243  enddo
244  enddo ; enddo
245  endif
246 
247  endif
248  endif ! restart
249 
250  if ( cs%use_sponge ) then
251 ! If sponges are used, this example damps tracers in sponges in the
252 ! northern half of the domain to 1 and tracers in the southern half
253 ! to 0. For any tracers that are not damped in the sponge, the call
254 ! to set_up_sponge_field can simply be omitted.
255  if (.not.associated(sponge_csp)) &
256  call mom_error(fatal, "DOME_initialize_tracer: "// &
257  "The pointer to sponge_CSp must be associated if SPONGE is defined.")
258 
259  allocate(temp(g%isd:g%ied,g%jsd:g%jed,nz))
260  do k=1,nz ; do j=js,je ; do i=is,ie
261  if (g%geoLatT(i,j) > 700.0 .and. (k > nz/2)) then
262  temp(i,j,k) = 1.0
263  else
264  temp(i,j,k) = 0.0
265  endif
266  enddo ; enddo ; enddo
267 
268 ! do m=1,NTR
269  do m=1,1
270  ! This is needed to force the compiler not to do a copy in the sponge
271  ! calls. Curses on the designers and implementers of Fortran90.
272  tr_ptr => cs%tr(:,:,:,m)
273  call set_up_sponge_field(temp, tr_ptr, g, nz, sponge_csp)
274  enddo
275  deallocate(temp)
276  endif
277 
278 end subroutine initialize_dome_tracer
279 
280 !> This subroutine applies diapycnal diffusion and any other column
281 !! tracer physics or chemistry to the tracers from this file.
282 !! This is a simple example of a set of advected passive tracers.
283 !!
284 !! The arguments to this subroutine are redundant in that
285 !! h_new(k) = h_old(k) + ea(k) - eb(k-1) + eb(k) - ea(k+1)
286 subroutine dome_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, &
287  evap_CFL_limit, minimum_forcing_depth)
288  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
289  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
290  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
291  intent(in) :: h_old !< Layer thickness before entrainment [H ~> m or kg m-2].
292  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
293  intent(in) :: h_new !< Layer thickness after entrainment [H ~> m or kg m-2].
294  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
295  intent(in) :: ea !< an array to which the amount of fluid entrained
296  !! from the layer above during this call will be
297  !! added [H ~> m or kg m-2].
298  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
299  intent(in) :: eb !< an array to which the amount of fluid entrained
300  !! from the layer below during this call will be
301  !! added [H ~> m or kg m-2].
302  type(forcing), intent(in) :: fluxes !< A structure containing pointers to thermodynamic
303  !! and tracer forcing fields. Unused fields have NULL ptrs.
304  real, intent(in) :: dt !< The amount of time covered by this call [s]
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 [m]
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 
336 end subroutine dome_tracer_column_physics
337 
338 !> This subroutine extracts the surface fields from this tracer package that
339 !! are to be shared with the atmosphere in coupled configurations.
340 !! This particular tracer package does not report anything back to the coupler.
341 subroutine dome_tracer_surface_state(state, h, G, CS)
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 
369 end subroutine dome_tracer_surface_state
370 
371 !> Clean up memory allocations, if any.
372 subroutine dome_tracer_end(CS)
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
381 end subroutine dome_tracer_end
382 
383 !> \namespace dome_tracer
384 !!
385 !! By Robert Hallberg, 2002
386 !!
387 !! This file contains an example of the code that is needed to set
388 !! up and use a set (in this case eleven) of dynamically passive
389 !! tracers. These tracers dye the inflowing water or water initially
390 !! within a range of latitudes or water initially in a range of
391 !! depths.
392 !!
393 !! A single subroutine is called from within each file to register
394 !! each of the tracers for reinitialization and advection and to
395 !! register the subroutine that initializes the tracers and set up
396 !! their output and the subroutine that does any tracer physics or
397 !! chemistry along with diapycnal mixing (included here because some
398 !! tracers may float or swim vertically or dye diapycnal processes).
399 
400 end module dome_tracer
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_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_tracer_diabatic
This module contains routines that implement physical fluxes of tracers (e.g. due to surface fluxes o...
Definition: MOM_tracer_diabatic.F90:4
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_open_boundary::obc_segment_tracer_type
Tracer segment data structure, for putting into an array of objects, not all the same shape.
Definition: MOM_open_boundary.F90:88
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
dome_tracer
A tracer package that is used as a diagnostic in the DOME experiments.
Definition: DOME_tracer.F90:2
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
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_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
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
dome_tracer::dome_tracer_cs
The DOME_tracer control structure.
Definition: DOME_tracer.F90:42
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_open_boundary::obc_segment_type
Open boundary segment data structure.
Definition: MOM_open_boundary.F90:107
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