MOM6
mom_variables Module Reference

Data Types

type  accel_diag_ptrs
 The accel_diag_ptrs structure contains pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances. More...
 
type  bt_cont_type
 The BT_cont_type structure contains information about the summed layer transports and how they will vary as the barotropic velocity is changed. More...
 
type  cont_diag_ptrs
 The cont_diag_ptrs structure contains pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances. More...
 
type  ocean_internal_state
 The ocean_internal_state structure contains pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90. It is useful for sending these variables for diagnostics, and in preparation for ensembles later on. All variables have the same names as the local (public) variables they refer to in MOM.F90. More...
 
type  p2d
 
type  p3d
 
type  surface
 The following structure contains pointers to various fields which may be used describe the surface state of MOM, and which will be returned to a the calling program. More...
 
type  thermo_var_ptrs
 The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be available, including potential temperature, salinity, heat capacity, and the equation of state control structure. More...
 
type  vertvisc_type
 The vertvisc_type structure contains vertical viscosities, drag coefficients, and related fields. More...
 

Functions/Subroutines

subroutine, public alloc_bt_cont_type (BT_cont, G, alloc_faces)
 alloc_BT_cont_type allocates the arrays contained within a BT_cont_type and initializes them to 0. More...
 
subroutine, public dealloc_bt_cont_type (BT_cont)
 dealloc_BT_cont_type deallocates the arrays contained within a BT_cont_type. More...
 
subroutine, public mom_thermovar_chksum (mesg, tv, G)
 MOM_thermovar_chksum does diagnostic checksums on various elements of a thermo_var_ptrs type for debugging. More...
 

Function/Subroutine Documentation

◆ alloc_bt_cont_type()

subroutine, public mom_variables::alloc_bt_cont_type ( type(bt_cont_type), pointer  BT_cont,
type(ocean_grid_type), intent(in)  G,
logical, intent(in), optional  alloc_faces 
)

alloc_BT_cont_type allocates the arrays contained within a BT_cont_type and initializes them to 0.

Parameters
[in]gThe ocean's grid structure

Definition at line 286 of file MOM_variables.F90.

References mom_error_handler::mom_error().

Referenced by mom_barotropic::barotropic_init().

286  type(bt_cont_type), pointer :: bt_cont
287  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
288  logical, optional, intent(in) :: alloc_faces
289 
290  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
291  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
292  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
293 
294  if (associated(bt_cont)) call mom_error(fatal, &
295  "alloc_BT_cont_type called with an associated BT_cont_type pointer.")
296 
297  allocate(bt_cont)
298  allocate(bt_cont%FA_u_WW(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_WW(:,:) = 0.0
299  allocate(bt_cont%FA_u_W0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_W0(:,:) = 0.0
300  allocate(bt_cont%FA_u_E0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_E0(:,:) = 0.0
301  allocate(bt_cont%FA_u_EE(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_EE(:,:) = 0.0
302  allocate(bt_cont%uBT_WW(isdb:iedb,jsd:jed)) ; bt_cont%uBT_WW(:,:) = 0.0
303  allocate(bt_cont%uBT_EE(isdb:iedb,jsd:jed)) ; bt_cont%uBT_EE(:,:) = 0.0
304 
305  allocate(bt_cont%FA_v_SS(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_SS(:,:) = 0.0
306  allocate(bt_cont%FA_v_S0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_S0(:,:) = 0.0
307  allocate(bt_cont%FA_v_N0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_N0(:,:) = 0.0
308  allocate(bt_cont%FA_v_NN(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_NN(:,:) = 0.0
309  allocate(bt_cont%vBT_SS(isd:ied,jsdb:jedb)) ; bt_cont%vBT_SS(:,:) = 0.0
310  allocate(bt_cont%vBT_NN(isd:ied,jsdb:jedb)) ; bt_cont%vBT_NN(:,:) = 0.0
311 
312  if (present(alloc_faces)) then ; if (alloc_faces) then
313  allocate(bt_cont%h_u(isdb:iedb,jsd:jed,1:g%ke)) ; bt_cont%h_u(:,:,:) = 0.0
314  allocate(bt_cont%h_v(isd:ied,jsdb:jedb,1:g%ke)) ; bt_cont%h_v(:,:,:) = 0.0
315  endif ; endif
316 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the call graph for this function:
Here is the caller graph for this function:

◆ dealloc_bt_cont_type()

subroutine, public mom_variables::dealloc_bt_cont_type ( type(bt_cont_type), pointer  BT_cont)

dealloc_BT_cont_type deallocates the arrays contained within a BT_cont_type.

Definition at line 321 of file MOM_variables.F90.

321  type(bt_cont_type), pointer :: bt_cont
322 
323  if (.not.associated(bt_cont)) return
324 
325  deallocate(bt_cont%FA_u_WW) ; deallocate(bt_cont%FA_u_W0)
326  deallocate(bt_cont%FA_u_E0) ; deallocate(bt_cont%FA_u_EE)
327  deallocate(bt_cont%uBT_WW) ; deallocate(bt_cont%uBT_EE)
328 
329  deallocate(bt_cont%FA_v_SS) ; deallocate(bt_cont%FA_v_S0)
330  deallocate(bt_cont%FA_v_N0) ; deallocate(bt_cont%FA_v_NN)
331  deallocate(bt_cont%vBT_SS) ; deallocate(bt_cont%vBT_NN)
332 
333  if (associated(bt_cont%h_u)) deallocate(bt_cont%h_u)
334  if (associated(bt_cont%h_v)) deallocate(bt_cont%h_v)
335 
336  deallocate(bt_cont)
337 

◆ mom_thermovar_chksum()

subroutine, public mom_variables::mom_thermovar_chksum ( character(len=*), intent(in)  mesg,
type(thermo_var_ptrs), intent(in)  tv,
type(ocean_grid_type), intent(in)  G 
)

MOM_thermovar_chksum does diagnostic checksums on various elements of a thermo_var_ptrs type for debugging.

Parameters
[in]tvA structure pointing to various thermodynamic variables
[in]gThe ocean's grid structure

Definition at line 343 of file MOM_variables.F90.

Referenced by mom_diabatic_driver::diabatic().

343  character(len=*), intent(in) :: mesg
344  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
345  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
346 ! This subroutine writes out chksums for the model's basic state variables.
347 ! Arguments: mesg - A message that appears on the chksum lines.
348 ! (in) u - Zonal velocity, in m s-1.
349 ! (in) v - Meridional velocity, in m s-1.
350 ! (in) h - Layer thickness, in m.
351 ! (in) uh - Volume flux through zonal faces = u*h*dy, m3 s-1.
352 ! (in) vh - Volume flux through meridional faces = v*h*dx, in m3 s-1.
353 ! (in) G - The ocean's grid structure.
354  integer :: is, ie, js, je, nz
355  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
356 
357  ! Note that for the chksum calls to be useful for reproducing across PE
358  ! counts, there must be no redundant points, so all variables use is..ie
359  ! and js...je as their extent.
360  if (associated(tv%T)) &
361  call hchksum(tv%T, mesg//" tv%T",g%HI)
362  if (associated(tv%S)) &
363  call hchksum(tv%S, mesg//" tv%S",g%HI)
364  if (associated(tv%frazil)) &
365  call hchksum(tv%frazil, mesg//" tv%frazil",g%HI)
366  if (associated(tv%salt_deficit)) &
367  call hchksum(tv%salt_deficit, mesg//" tv%salt_deficit",g%HI)
368  if (associated(tv%TempxPmE)) &
369  call hchksum(tv%TempxPmE, mesg//" tv%TempxPmE",g%HI)
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Here is the caller graph for this function: