MOM6
advection_test_tracer.F90
1 !> This tracer package is used to test advection schemes
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, read_data, slasher, vardesc, var_desc, query_vardesc
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_tracer_diabatic, only : tracer_vertdiff, applytracerboundaryfluxesinout
19 use mom_variables, only : surface
21 
22 use coupler_types_mod, only : coupler_type_set_data, ind_csurf
23 use atmos_ocean_fluxes_mod, only : aof_set_coupler_flux
24 
25 implicit none ; private
26 
27 #include <MOM_memory.h>
28 
29 public register_advection_test_tracer, initialize_advection_test_tracer
30 public advection_test_tracer_surface_state, advection_test_tracer_end
31 public advection_test_tracer_column_physics, advection_test_stock
32 
33 integer, parameter :: ntr = 11 !< The number of tracers in this module.
34 
35 !> The control structure for the advect_test_tracer module
36 type, public :: advection_test_tracer_cs ; private
37  integer :: ntr = ntr !< Number of tracers in this module
38  logical :: coupled_tracers = .false. !< These tracers are not offered to the coupler.
39  character(len=200) :: tracer_ic_file !< The full path to the IC file, or " " to initialize internally.
40  type(time_type), pointer :: time => null() !< A pointer to the ocean model's clock.
41  type(tracer_registry_type), pointer :: tr_reg => null() !< A pointer to the MOM tracer registry
42  real, pointer :: tr(:,:,:,:) => null() !< The array of tracers used in this subroutine, in g m-3?
43  real :: land_val(ntr) = -1.0 !< The value of tr used where land is masked out.
44  logical :: use_sponge !< If true, sponges may be applied somewhere in the domain.
45  logical :: tracers_may_reinit !< If true, the tracers may be set up via the initialization code if
46  !! they are not found in the restart files. Otherwise it is a fatal error
47  !! if the tracers are not found in the restart files of a restarted run.
48  real :: x_origin !< Parameters describing the test functions
49  real :: x_width !< Parameters describing the test functions
50  real :: y_origin !< Parameters describing the test functions
51  real :: y_width !< Parameters describing the test functions
52 
53  integer, dimension(NTR) :: ind_tr !< Indices returned by aof_set_coupler_flux if it is used and
54  !! the surface tracer concentrations are to be provided to the coupler.
55 
56  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
57  !! regulate the timing of diagnostic output.
58  type(mom_restart_cs), pointer :: restart_csp => null() !< A pointer to the restart control structure.
59 
60  type(vardesc) :: tr_desc(ntr) !< Descriptions and metadata for the tracers
62 
63 contains
64 
65 !> Register tracer fields and subroutines to be used with MOM.
66 function register_advection_test_tracer(HI, GV, param_file, CS, tr_Reg, restart_CS)
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.
162 end function register_advection_test_tracer
163 
164 !> Initializes the NTR tracer fields in tr(:,:,:,:) and it sets up the tracer output.
165 subroutine initialize_advection_test_tracer(restart, day, G, GV, h,diag, OBC, CS, &
166  sponge_CSp)
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 
255 end subroutine initialize_advection_test_tracer
256 
257 
258 !> Applies diapycnal diffusion and any other column tracer physics or chemistry to the tracers
259 !! from this package. This is a simple example of a set of advected passive tracers.
260 subroutine advection_test_tracer_column_physics(h_old, h_new, ea, eb, fluxes, dt, G, GV, CS, &
261  evap_CFL_limit, minimum_forcing_depth)
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 
315 end subroutine advection_test_tracer_column_physics
316 
317 !> This subroutine extracts the surface fields from this tracer package that
318 !! are to be shared with the atmosphere in coupled configurations.
319 !! This particular tracer package does not report anything back to the coupler.
320 subroutine advection_test_tracer_surface_state(state, h, G, CS)
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 
348 end subroutine advection_test_tracer_surface_state
349 
350 !> Calculate the mass-weighted integral of all tracer stocks, returning the number of stocks it has calculated.
351 !! If the stock_index is present, only the stock corresponding to that coded index is returned.
352 function advection_test_stock(h, stocks, G, GV, CS, names, units, stock_index)
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 
390 end function advection_test_stock
391 
392 !> Deallocate memory associated with this module
393 subroutine advection_test_tracer_end(CS)
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
402 end subroutine advection_test_tracer_end
403 
404 end module advection_test_tracer
advection_test_tracer::advection_test_tracer_cs
The control structure for the advect_test_tracer module.
Definition: advection_test_tracer.F90:36
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_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
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.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
advection_test_tracer
This tracer package is used to test advection schemes.
Definition: advection_test_tracer.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_restart::query_initialized
Indicate whether a field has been read from a restart file.
Definition: MOM_restart.F90:116
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