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

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

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

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

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

522  type(bt_cont_type), pointer :: bt_cont !< The BT_cont_type whose elements will be deallocated.
523 
524  if (.not.associated(bt_cont)) return
525 
526  deallocate(bt_cont%FA_u_WW) ; deallocate(bt_cont%FA_u_W0)
527  deallocate(bt_cont%FA_u_E0) ; deallocate(bt_cont%FA_u_EE)
528  deallocate(bt_cont%uBT_WW) ; deallocate(bt_cont%uBT_EE)
529 
530  deallocate(bt_cont%FA_v_SS) ; deallocate(bt_cont%FA_v_S0)
531  deallocate(bt_cont%FA_v_N0) ; deallocate(bt_cont%FA_v_NN)
532  deallocate(bt_cont%vBT_SS) ; deallocate(bt_cont%vBT_NN)
533 
534  if (allocated(bt_cont%h_u)) deallocate(bt_cont%h_u)
535  if (allocated(bt_cont%h_v)) deallocate(bt_cont%h_v)
536 
537  deallocate(bt_cont)
538 

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

388  type(surface), intent(inout) :: sfc_state !< ocean surface state type to be deallocated here.
389 
390  if (.not.sfc_state%arrays_allocated) return
391 
392  if (allocated(sfc_state%melt_potential)) deallocate(sfc_state%melt_potential)
393  if (allocated(sfc_state%SST)) deallocate(sfc_state%SST)
394  if (allocated(sfc_state%SSS)) deallocate(sfc_state%SSS)
395  if (allocated(sfc_state%sfc_density)) deallocate(sfc_state%sfc_density)
396  if (allocated(sfc_state%sea_lev)) deallocate(sfc_state%sea_lev)
397  if (allocated(sfc_state%Hml)) deallocate(sfc_state%Hml)
398  if (allocated(sfc_state%u)) deallocate(sfc_state%u)
399  if (allocated(sfc_state%v)) deallocate(sfc_state%v)
400  if (allocated(sfc_state%ocean_mass)) deallocate(sfc_state%ocean_mass)
401  if (allocated(sfc_state%ocean_heat)) deallocate(sfc_state%ocean_heat)
402  if (allocated(sfc_state%ocean_salt)) deallocate(sfc_state%ocean_salt)
403  if (allocated(sfc_state%salt_deficit)) deallocate(sfc_state%salt_deficit)
404 
405  call coupler_type_destructor(sfc_state%tr_fields)
406 
407  sfc_state%arrays_allocated = .false.
408 

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

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