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

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

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

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

Definition at line 303 of file MOM_variables.F90.

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

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

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

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

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

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

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