MOM6
mom_variables Module Reference

Detailed Description

Provides transparent structures with groups of MOM6 variables and supporting routines.

Data Types

type  accel_diag_ptrs
 Pointers to arrays with accelerations, which can later be used for derived diagnostics, like energy balances. More...
 
type  bt_cont_type
 Container for information about the summed layer transports and how they will vary as the barotropic velocity is changed. More...
 
type  cont_diag_ptrs
 Pointers to arrays with transports, which can later be used for derived diagnostics, like energy balances. More...
 
type  ocean_internal_state
 Pointers to all of the prognostic variables allocated in MOM_variables.F90 and MOM.F90. More...
 
type  p2d
 A structure for creating arrays of pointers to 2D arrays. More...
 
type  p3d
 A structure for creating arrays of pointers to 3D arrays. More...
 
type  surface
 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
 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
 Vertical viscosities, drag coefficients, and related fields. More...
 

Functions/Subroutines

subroutine, public allocate_surface_state (sfc_state, G, use_temperature, do_integrals, gas_fields_ocn, use_meltpot, use_iceshelves)
 Allocates the fields for the surface (return) properties of the ocean model. Unused fields are unallocated. More...
 
subroutine, public deallocate_surface_state (sfc_state)
 Deallocates the elements of a surface state type. More...
 
subroutine, public alloc_bt_cont_type (BT_cont, G, alloc_faces)
 Allocates the arrays contained within a BT_cont_type and initializes them to 0. More...
 
subroutine, public dealloc_bt_cont_type (BT_cont)
 Deallocates the arrays contained within a BT_cont_type. More...
 
subroutine, public mom_thermovar_chksum (mesg, tv, G)
 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 
)

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

Parameters
bt_contThe BT_cont_type whose elements will be allocated
[in]gThe ocean's grid structure
[in]alloc_facesIf present and true, allocate memory for effective face thicknesses.

Definition at line 395 of file MOM_variables.F90.

395  type(BT_cont_type), pointer :: BT_cont !< The BT_cont_type whose elements will be allocated
396  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
397  logical, optional, intent(in) :: alloc_faces !< If present and true, allocate
398  !! memory for effective face thicknesses.
399 
400  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
401  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
402  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
403 
404  if (associated(bt_cont)) call mom_error(fatal, &
405  "alloc_BT_cont_type called with an associated BT_cont_type pointer.")
406 
407  allocate(bt_cont)
408  allocate(bt_cont%FA_u_WW(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_WW(:,:) = 0.0
409  allocate(bt_cont%FA_u_W0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_W0(:,:) = 0.0
410  allocate(bt_cont%FA_u_E0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_E0(:,:) = 0.0
411  allocate(bt_cont%FA_u_EE(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_EE(:,:) = 0.0
412  allocate(bt_cont%uBT_WW(isdb:iedb,jsd:jed)) ; bt_cont%uBT_WW(:,:) = 0.0
413  allocate(bt_cont%uBT_EE(isdb:iedb,jsd:jed)) ; bt_cont%uBT_EE(:,:) = 0.0
414 
415  allocate(bt_cont%FA_v_SS(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_SS(:,:) = 0.0
416  allocate(bt_cont%FA_v_S0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_S0(:,:) = 0.0
417  allocate(bt_cont%FA_v_N0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_N0(:,:) = 0.0
418  allocate(bt_cont%FA_v_NN(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_NN(:,:) = 0.0
419  allocate(bt_cont%vBT_SS(isd:ied,jsdb:jedb)) ; bt_cont%vBT_SS(:,:) = 0.0
420  allocate(bt_cont%vBT_NN(isd:ied,jsdb:jedb)) ; bt_cont%vBT_NN(:,:) = 0.0
421 
422  if (present(alloc_faces)) then ; if (alloc_faces) then
423  allocate(bt_cont%h_u(isdb:iedb,jsd:jed,1:g%ke)) ; bt_cont%h_u(:,:,:) = 0.0
424  allocate(bt_cont%h_v(isd:ied,jsdb:jedb,1:g%ke)) ; bt_cont%h_v(:,:,:) = 0.0
425  endif ; endif
426 

References mom_error_handler::mom_error().

Referenced by mom_barotropic::barotropic_init().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ allocate_surface_state()

subroutine, public mom_variables::allocate_surface_state ( type(surface), intent(inout)  sfc_state,
type(ocean_grid_type), intent(in)  G,
logical, intent(in), optional  use_temperature,
logical, intent(in), optional  do_integrals,
type(coupler_1d_bc_type), intent(in), optional  gas_fields_ocn,
logical, intent(in), optional  use_meltpot,
logical, intent(in), optional  use_iceshelves 
)

Allocates the fields for the surface (return) properties of the ocean model. Unused fields are unallocated.

Parameters
[in]gocean grid structure
[in,out]sfc_stateocean surface state type to be allocated.
[in]use_temperatureIf true, allocate the space for thermodynamic variables.
[in]do_integralsIf true, allocate the space for vertically integrated fields.
[in]gas_fields_ocnIf present, this type describes the ocean
[in]use_meltpotIf true, allocate the space for melt potential
[in]use_iceshelvesIf true, allocate the space for the stresses under ice shelves.

Definition at line 297 of file MOM_variables.F90.

297  type(ocean_grid_type), intent(in) :: G !< ocean grid structure
298  type(surface), intent(inout) :: sfc_state !< ocean surface state type to be allocated.
299  logical, optional, intent(in) :: use_temperature !< If true, allocate the space for thermodynamic variables.
300  logical, optional, intent(in) :: do_integrals !< If true, allocate the space for vertically
301  !! integrated fields.
302  type(coupler_1d_bc_type), &
303  optional, intent(in) :: gas_fields_ocn !< If present, this type describes the ocean
304  !! ocean and surface-ice fields that will participate
305  !! in the calculation of additional gas or other
306  !! tracer fluxes, and can be used to spawn related
307  !! internal variables in the ice model.
308  logical, optional, intent(in) :: use_meltpot !< If true, allocate the space for melt potential
309  logical, optional, intent(in) :: use_iceshelves !< If true, allocate the space for the stresses
310  !! under ice shelves.
311 
312  ! local variables
313  logical :: use_temp, alloc_integ, use_melt_potential, alloc_iceshelves
314  integer :: is, ie, js, je, isd, ied, jsd, jed
315  integer :: isdB, iedB, jsdB, jedB
316 
317  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
318  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
319  isdb = g%isdB ; iedb = g%iedB; jsdb = g%jsdB ; jedb = g%jedB
320 
321  use_temp = .true. ; if (present(use_temperature)) use_temp = use_temperature
322  alloc_integ = .true. ; if (present(do_integrals)) alloc_integ = do_integrals
323  use_melt_potential = .false. ; if (present(use_meltpot)) use_melt_potential = use_meltpot
324  alloc_iceshelves = .false. ; if (present(use_iceshelves)) alloc_iceshelves = use_iceshelves
325 
326  if (sfc_state%arrays_allocated) return
327 
328  if (use_temp) then
329  allocate(sfc_state%SST(isd:ied,jsd:jed)) ; sfc_state%SST(:,:) = 0.0
330  allocate(sfc_state%SSS(isd:ied,jsd:jed)) ; sfc_state%SSS(:,:) = 0.0
331  else
332  allocate(sfc_state%sfc_density(isd:ied,jsd:jed)) ; sfc_state%sfc_density(:,:) = 0.0
333  endif
334  allocate(sfc_state%sea_lev(isd:ied,jsd:jed)) ; sfc_state%sea_lev(:,:) = 0.0
335  allocate(sfc_state%Hml(isd:ied,jsd:jed)) ; sfc_state%Hml(:,:) = 0.0
336  allocate(sfc_state%u(isdb:iedb,jsd:jed)) ; sfc_state%u(:,:) = 0.0
337  allocate(sfc_state%v(isd:ied,jsdb:jedb)) ; sfc_state%v(:,:) = 0.0
338 
339  if (use_melt_potential) then
340  allocate(sfc_state%melt_potential(isd:ied,jsd:jed)) ; sfc_state%melt_potential(:,:) = 0.0
341  endif
342 
343  if (alloc_integ) then
344  ! Allocate structures for the vertically integrated ocean_mass, ocean_heat, and ocean_salt.
345  allocate(sfc_state%ocean_mass(isd:ied,jsd:jed)) ; sfc_state%ocean_mass(:,:) = 0.0
346  if (use_temp) then
347  allocate(sfc_state%ocean_heat(isd:ied,jsd:jed)) ; sfc_state%ocean_heat(:,:) = 0.0
348  allocate(sfc_state%ocean_salt(isd:ied,jsd:jed)) ; sfc_state%ocean_salt(:,:) = 0.0
349  allocate(sfc_state%TempxPmE(isd:ied,jsd:jed)) ; sfc_state%TempxPmE(:,:) = 0.0
350  allocate(sfc_state%salt_deficit(isd:ied,jsd:jed)) ; sfc_state%salt_deficit(:,:) = 0.0
351  allocate(sfc_state%internal_heat(isd:ied,jsd:jed)) ; sfc_state%internal_heat(:,:) = 0.0
352  endif
353  endif
354 
355  if (alloc_iceshelves) then
356  allocate(sfc_state%taux_shelf(isdb:iedb,jsd:jed)) ; sfc_state%taux_shelf(:,:) = 0.0
357  allocate(sfc_state%tauy_shelf(isd:ied,jsdb:jedb)) ; sfc_state%tauy_shelf(:,:) = 0.0
358  endif
359 
360  if (present(gas_fields_ocn)) &
361  call coupler_type_spawn(gas_fields_ocn, sfc_state%tr_fields, &
362  (/is,is,ie,ie/), (/js,js,je,je/), as_needed=.true.)
363 
364  sfc_state%arrays_allocated = .true.
365 

◆ dealloc_bt_cont_type()

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

Deallocates the arrays contained within a BT_cont_type.

Parameters
bt_contThe BT_cont_type whose elements will be deallocated.

Definition at line 431 of file MOM_variables.F90.

431  type(BT_cont_type), pointer :: BT_cont !< The BT_cont_type whose elements will be deallocated.
432 
433  if (.not.associated(bt_cont)) return
434 
435  deallocate(bt_cont%FA_u_WW) ; deallocate(bt_cont%FA_u_W0)
436  deallocate(bt_cont%FA_u_E0) ; deallocate(bt_cont%FA_u_EE)
437  deallocate(bt_cont%uBT_WW) ; deallocate(bt_cont%uBT_EE)
438 
439  deallocate(bt_cont%FA_v_SS) ; deallocate(bt_cont%FA_v_S0)
440  deallocate(bt_cont%FA_v_N0) ; deallocate(bt_cont%FA_v_NN)
441  deallocate(bt_cont%vBT_SS) ; deallocate(bt_cont%vBT_NN)
442 
443  if (allocated(bt_cont%h_u)) deallocate(bt_cont%h_u)
444  if (allocated(bt_cont%h_v)) deallocate(bt_cont%h_v)
445 
446  deallocate(bt_cont)
447 

◆ deallocate_surface_state()

subroutine, public mom_variables::deallocate_surface_state ( type(surface), intent(inout)  sfc_state)

Deallocates the elements of a surface state type.

Parameters
[in,out]sfc_stateocean surface state type to be deallocated here.

Definition at line 370 of file MOM_variables.F90.

370  type(surface), intent(inout) :: sfc_state !< ocean surface state type to be deallocated here.
371 
372  if (.not.sfc_state%arrays_allocated) return
373 
374  if (allocated(sfc_state%melt_potential)) deallocate(sfc_state%melt_potential)
375  if (allocated(sfc_state%SST)) deallocate(sfc_state%SST)
376  if (allocated(sfc_state%SSS)) deallocate(sfc_state%SSS)
377  if (allocated(sfc_state%sfc_density)) deallocate(sfc_state%sfc_density)
378  if (allocated(sfc_state%sea_lev)) deallocate(sfc_state%sea_lev)
379  if (allocated(sfc_state%Hml)) deallocate(sfc_state%Hml)
380  if (allocated(sfc_state%u)) deallocate(sfc_state%u)
381  if (allocated(sfc_state%v)) deallocate(sfc_state%v)
382  if (allocated(sfc_state%ocean_mass)) deallocate(sfc_state%ocean_mass)
383  if (allocated(sfc_state%ocean_heat)) deallocate(sfc_state%ocean_heat)
384  if (allocated(sfc_state%ocean_salt)) deallocate(sfc_state%ocean_salt)
385  if (allocated(sfc_state%salt_deficit)) deallocate(sfc_state%salt_deficit)
386 
387  call coupler_type_destructor(sfc_state%tr_fields)
388 
389  sfc_state%arrays_allocated = .false.
390 

◆ 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 
)

Diagnostic checksums on various elements of a thermo_var_ptrs type for debugging.

Parameters
[in]mesgA message that appears in the checksum lines
[in]tvA structure pointing to various thermodynamic variables
[in]gThe ocean's grid structure

Definition at line 452 of file MOM_variables.F90.

452  character(len=*), intent(in) :: mesg !< A message that appears in the checksum lines
453  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
454  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
455 
456  ! Note that for the chksum calls to be useful for reproducing across PE
457  ! counts, there must be no redundant points, so all variables use is..ie
458  ! and js...je as their extent.
459  if (associated(tv%T)) &
460  call hchksum(tv%T, mesg//" tv%T", g%HI)
461  if (associated(tv%S)) &
462  call hchksum(tv%S, mesg//" tv%S", g%HI)
463  if (associated(tv%frazil)) &
464  call hchksum(tv%frazil, mesg//" tv%frazil", g%HI)
465  if (associated(tv%salt_deficit)) &
466  call hchksum(tv%salt_deficit, mesg//" tv%salt_deficit", g%HI, scale=g%US%R_to_kg_m3*g%US%Z_to_m)
467  if (associated(tv%TempxPmE)) &
468  call hchksum(tv%TempxPmE, mesg//" tv%TempxPmE", g%HI, scale=g%US%R_to_kg_m3*g%US%Z_to_m)

Referenced by mom_diabatic_driver::diabatic_ale(), mom_diabatic_driver::diabatic_ale_legacy(), and mom_diabatic_driver::layered_diabatic().

Here is the caller graph for this function: