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, omit_frazil)
 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 rotate_surface_state (sfc_state_in, G_in, sfc_state, G, turns)
 Rotate the surface state fields from the input to the model indices.
 
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 477 of file MOM_variables.F90.

477  type(bt_cont_type), pointer :: bt_cont !< The BT_cont_type whose elements will be allocated
478  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
479  logical, optional, intent(in) :: alloc_faces !< If present and true, allocate
480  !! memory for effective face thicknesses.
481 
482  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
483  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
484  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
485 
486  if (associated(bt_cont)) call mom_error(fatal, &
487  "alloc_BT_cont_type called with an associated BT_cont_type pointer.")
488 
489  allocate(bt_cont)
490  allocate(bt_cont%FA_u_WW(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_WW(:,:) = 0.0
491  allocate(bt_cont%FA_u_W0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_W0(:,:) = 0.0
492  allocate(bt_cont%FA_u_E0(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_E0(:,:) = 0.0
493  allocate(bt_cont%FA_u_EE(isdb:iedb,jsd:jed)) ; bt_cont%FA_u_EE(:,:) = 0.0
494  allocate(bt_cont%uBT_WW(isdb:iedb,jsd:jed)) ; bt_cont%uBT_WW(:,:) = 0.0
495  allocate(bt_cont%uBT_EE(isdb:iedb,jsd:jed)) ; bt_cont%uBT_EE(:,:) = 0.0
496 
497  allocate(bt_cont%FA_v_SS(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_SS(:,:) = 0.0
498  allocate(bt_cont%FA_v_S0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_S0(:,:) = 0.0
499  allocate(bt_cont%FA_v_N0(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_N0(:,:) = 0.0
500  allocate(bt_cont%FA_v_NN(isd:ied,jsdb:jedb)) ; bt_cont%FA_v_NN(:,:) = 0.0
501  allocate(bt_cont%vBT_SS(isd:ied,jsdb:jedb)) ; bt_cont%vBT_SS(:,:) = 0.0
502  allocate(bt_cont%vBT_NN(isd:ied,jsdb:jedb)) ; bt_cont%vBT_NN(:,:) = 0.0
503 
504  if (present(alloc_faces)) then ; if (alloc_faces) then
505  allocate(bt_cont%h_u(isdb:iedb,jsd:jed,1:g%ke)) ; bt_cont%h_u(:,:,:) = 0.0
506  allocate(bt_cont%h_v(isd:ied,jsdb:jedb,1:g%ke)) ; bt_cont%h_v(:,:,:) = 0.0
507  endif ; endif
508 

◆ 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,
logical, intent(in), optional  omit_frazil 
)

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.
[in]omit_frazilIf present and false, do not allocate the space to pass frazil fluxes to the coupler

Definition at line 300 of file MOM_variables.F90.

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

◆ 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 513 of file MOM_variables.F90.

513  type(bt_cont_type), pointer :: bt_cont !< The BT_cont_type whose elements will be deallocated.
514 
515  if (.not.associated(bt_cont)) return
516 
517  deallocate(bt_cont%FA_u_WW) ; deallocate(bt_cont%FA_u_W0)
518  deallocate(bt_cont%FA_u_E0) ; deallocate(bt_cont%FA_u_EE)
519  deallocate(bt_cont%uBT_WW) ; deallocate(bt_cont%uBT_EE)
520 
521  deallocate(bt_cont%FA_v_SS) ; deallocate(bt_cont%FA_v_S0)
522  deallocate(bt_cont%FA_v_N0) ; deallocate(bt_cont%FA_v_NN)
523  deallocate(bt_cont%vBT_SS) ; deallocate(bt_cont%vBT_NN)
524 
525  if (allocated(bt_cont%h_u)) deallocate(bt_cont%h_u)
526  if (allocated(bt_cont%h_v)) deallocate(bt_cont%h_v)
527 
528  deallocate(bt_cont)
529 

◆ 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 379 of file MOM_variables.F90.

379  type(surface), intent(inout) :: sfc_state !< ocean surface state type to be deallocated here.
380 
381  if (.not.sfc_state%arrays_allocated) return
382 
383  if (allocated(sfc_state%melt_potential)) deallocate(sfc_state%melt_potential)
384  if (allocated(sfc_state%SST)) deallocate(sfc_state%SST)
385  if (allocated(sfc_state%SSS)) deallocate(sfc_state%SSS)
386  if (allocated(sfc_state%sfc_density)) deallocate(sfc_state%sfc_density)
387  if (allocated(sfc_state%sea_lev)) deallocate(sfc_state%sea_lev)
388  if (allocated(sfc_state%Hml)) deallocate(sfc_state%Hml)
389  if (allocated(sfc_state%u)) deallocate(sfc_state%u)
390  if (allocated(sfc_state%v)) deallocate(sfc_state%v)
391  if (allocated(sfc_state%ocean_mass)) deallocate(sfc_state%ocean_mass)
392  if (allocated(sfc_state%ocean_heat)) deallocate(sfc_state%ocean_heat)
393  if (allocated(sfc_state%ocean_salt)) deallocate(sfc_state%ocean_salt)
394  if (allocated(sfc_state%salt_deficit)) deallocate(sfc_state%salt_deficit)
395 
396  call coupler_type_destructor(sfc_state%tr_fields)
397 
398  sfc_state%arrays_allocated = .false.
399 

◆ 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 534 of file MOM_variables.F90.

534  character(len=*), intent(in) :: mesg !< A message that appears in the checksum lines
535  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
536  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
537 
538  ! Note that for the chksum calls to be useful for reproducing across PE
539  ! counts, there must be no redundant points, so all variables use is..ie
540  ! and js...je as their extent.
541  if (associated(tv%T)) &
542  call hchksum(tv%T, mesg//" tv%T", g%HI)
543  if (associated(tv%S)) &
544  call hchksum(tv%S, mesg//" tv%S", g%HI)
545  if (associated(tv%frazil)) &
546  call hchksum(tv%frazil, mesg//" tv%frazil", g%HI, scale=g%US%Q_to_J_kg*g%US%RZ_to_kg_m2)
547  if (associated(tv%salt_deficit)) &
548  call hchksum(tv%salt_deficit, mesg//" tv%salt_deficit", g%HI, scale=g%US%RZ_to_kg_m2)
549  if (associated(tv%TempxPmE)) &
550  call hchksum(tv%TempxPmE, mesg//" tv%TempxPmE", g%HI, scale=g%US%RZ_to_kg_m2)