MOM6
mom_continuity Module Reference

Detailed Description

Solve the layer continuity equation.

Data Types

type  continuity_cs
 Control structure for mom_continuity. More...
 

Functions/Subroutines

subroutine, public continuity (u, v, hin, h, uh, vh, dt, G, GV, US, CS, uhbt, vhbt, OBC, visc_rem_u, visc_rem_v, u_cor, v_cor, BT_cont)
 Time steps the layer thicknesses, using a monotonically limited, directionally split PPM scheme, based on Lin (1994). More...
 
subroutine, public continuity_init (Time, G, GV, US, param_file, diag, CS)
 Initializes continuity_cs. More...
 
integer function, public continuity_stencil (CS)
 continuity_stencil returns the continuity solver stencil size More...
 
subroutine, public continuity_end (CS)
 Destructor for continuity_cs. More...
 

Variables

integer, parameter ppm_scheme = 1
 Enumerated constant to select PPM. More...
 
character(len=20), parameter ppm_string = "PPM"
 String to select PPM. More...
 

Function/Subroutine Documentation

◆ continuity()

subroutine, public mom_continuity::continuity ( real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(in)  v,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  hin,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(inout)  h,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(out)  uh,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(out)  vh,
real, intent(in)  dt,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(continuity_cs), pointer  CS,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed), intent(in), optional  uhbt,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb), intent(in), optional  vhbt,
type(ocean_obc_type), optional, pointer  OBC,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(in), optional  visc_rem_u,
real, dimension(szi_(g),szjb_(g),szk_(g)), intent(in), optional  visc_rem_v,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(out), optional  u_cor,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(out), optional  v_cor,
type(bt_cont_type), optional, pointer  BT_cont 
)

Time steps the layer thicknesses, using a monotonically limited, directionally split PPM scheme, based on Lin (1994).

Parameters
[in,out]gOcean grid structure.
[in]gvVertical grid structure.
[in]uZonal velocity [L T-1 ~> m s-1].
[in]vMeridional velocity [L T-1 ~> m s-1].
[in]hinInitial layer thickness [H ~> m or kg m-2].
[in,out]hFinal layer thickness [H ~> m or kg m-2].
[out]uhVolume flux through zonal faces =
[out]vhVolume flux through meridional faces =
[in]dtTime increment [T ~> s].
[in]usA dimensional unit scaling type
csControl structure for mom_continuity.
[in]uhbtThe vertically summed volume
[in]vhbtThe vertically summed volume
obcOpen boundaries control structure.
[in]visc_rem_uBoth the fraction of
[in]visc_rem_vBoth the fraction of
[out]u_corThe zonal velocities that
[out]v_corThe meridional velocities that
bt_contA structure with elements

Definition at line 44 of file MOM_continuity.F90.

44  type(ocean_grid_type), intent(inout) :: G !< Ocean grid structure.
45  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
46  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
47  intent(in) :: u !< Zonal velocity [L T-1 ~> m s-1].
48  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
49  intent(in) :: v !< Meridional velocity [L T-1 ~> m s-1].
50  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
51  intent(in) :: hin !< Initial layer thickness [H ~> m or kg m-2].
52  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
53  intent(inout) :: h !< Final layer thickness [H ~> m or kg m-2].
54  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
55  intent(out) :: uh !< Volume flux through zonal faces =
56  !! u*h*dy [H L2 T-1 ~> m3 s-1 or kg s-1].
57  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
58  intent(out) :: vh !< Volume flux through meridional faces =
59  !! v*h*dx [H L2 T-1 ~> m3 s-1 or kg s-1].
60  real, intent(in) :: dt !< Time increment [T ~> s].
61  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
62  type(continuity_CS), pointer :: CS !< Control structure for mom_continuity.
63  real, dimension(SZIB_(G),SZJ_(G)), &
64  optional, intent(in) :: uhbt !< The vertically summed volume
65  !! flux through zonal faces [H L2 T-1 ~> m3 s-1 or kg s-1].
66  real, dimension(SZI_(G),SZJB_(G)), &
67  optional, intent(in) :: vhbt !< The vertically summed volume
68  !! flux through meridional faces [H L2 T-1 ~> m3 s-1 or kg s-1].
69  type(ocean_OBC_type), &
70  optional, pointer :: OBC !< Open boundaries control structure.
71  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
72  optional, intent(in) :: visc_rem_u !< Both the fraction of
73  !! zonal momentum that remains after a time-step of viscosity, and the fraction of a time-step's
74  !! worth of a barotropic acceleration that a layer experiences after viscosity is applied.
75  !! Non-dimensional between 0 (at the bottom) and 1 (far above the bottom).
76  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
77  optional, intent(in) :: visc_rem_v !< Both the fraction of
78  !! meridional momentum that remains after a time-step of viscosity, and the fraction of a time-step's
79  !! worth of a barotropic acceleration that a layer experiences after viscosity is applied.
80  !! Non-dimensional between 0 (at the bottom) and 1 (far above the bottom).
81  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
82  optional, intent(out) :: u_cor !< The zonal velocities that
83  !! give uhbt as the depth-integrated transport [L T-1 ~> m s-1].
84  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
85  optional, intent(out) :: v_cor !< The meridional velocities that
86  !! give vhbt as the depth-integrated transport [L T-1 ~> m s-1].
87  type(BT_cont_type), &
88  optional, pointer :: BT_cont !< A structure with elements
89  !! that describe the effective open face areas as a function of barotropic flow.
90 
91  if (present(visc_rem_u) .neqv. present(visc_rem_v)) call mom_error(fatal, &
92  "MOM_continuity: Either both visc_rem_u and visc_rem_v or neither"// &
93  " one must be present in call to continuity.")
94  if (present(u_cor) .neqv. present(v_cor)) call mom_error(fatal, &
95  "MOM_continuity: Either both u_cor and v_cor or neither"// &
96  " one must be present in call to continuity.")
97 
98  if (cs%continuity_scheme == ppm_scheme) then
99  call continuity_ppm(u, v, hin, h, uh, vh, dt, g, gv, us, cs%PPM_CSp, uhbt, vhbt, obc, &
100  visc_rem_u, visc_rem_v, u_cor, v_cor, bt_cont=bt_cont)
101  else
102  call mom_error(fatal, "continuity: Unrecognized value of continuity_scheme")
103  endif
104 

References mom_error_handler::mom_error(), and ppm_scheme.

Referenced by mom_dynamics_unsplit::step_mom_dyn_unsplit(), and mom_dynamics_unsplit_rk2::step_mom_dyn_unsplit_rk2().

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

◆ continuity_end()

subroutine, public mom_continuity::continuity_end ( type(continuity_cs), pointer  CS)

Destructor for continuity_cs.

Parameters
csControl structure for mom_continuity.

Definition at line 169 of file MOM_continuity.F90.

169  type(continuity_CS), pointer :: CS !< Control structure for mom_continuity.
170 
171  if (cs%continuity_scheme == ppm_scheme) then
172  call continuity_ppm_end(cs%PPM_CSp)
173  endif
174 
175  deallocate(cs)
176 

References mom_continuity_ppm::continuity_ppm_end(), and ppm_scheme.

Here is the call graph for this function:

◆ continuity_init()

subroutine, public mom_continuity::continuity_init ( type(time_type), intent(in), target  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(diag_ctrl), intent(inout), target  diag,
type(continuity_cs), pointer  CS 
)

Initializes continuity_cs.

Parameters
[in]timeCurrent model time.
[in]gOcean grid structure.
[in]gvVertical grid structure.
[in]usA dimensional unit scaling type
[in]param_fileParameter file handles.
[in,out]diagDiagnostics control structure.
csControl structure for mom_continuity.

Definition at line 109 of file MOM_continuity.F90.

109  type(time_type), target, intent(in) :: Time !< Current model time.
110  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
111  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
112  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
113  type(param_file_type), intent(in) :: param_file !< Parameter file handles.
114  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
115  type(continuity_CS), pointer :: CS !< Control structure for mom_continuity.
116 ! This include declares and sets the variable "version".
117 #include "version_variable.h"
118  character(len=40) :: mdl = "MOM_continuity" ! This module's name.
119  character(len=20) :: tmpstr
120 
121  if (associated(cs)) then
122  call mom_error(warning, "continuity_init called with associated control structure.")
123  return
124  endif
125  allocate(cs)
126 
127  ! Read all relevant parameters and write them to the model log.
128  call log_version(param_file, mdl, version, "")
129  call get_param(param_file, mdl, "CONTINUITY_SCHEME", tmpstr, &
130  "CONTINUITY_SCHEME selects the discretization for the "//&
131  "continuity solver. The only valid value currently is: \n"//&
132  "\t PPM - use a positive-definite (or monotonic) \n"//&
133  "\t piecewise parabolic reconstruction solver.", &
134  default=ppm_string)
135 
136  tmpstr = uppercase(tmpstr) ; cs%continuity_scheme = 0
137  select case (trim(tmpstr))
138  case (ppm_string) ; cs%continuity_scheme = ppm_scheme
139  case default
140  call mom_mesg('continuity_init: CONTINUITY_SCHEME ="'//trim(tmpstr)//'"', 0)
141  call mom_mesg("continuity_init: The only valid value is currently "// &
142  trim(ppm_string), 0)
143  call mom_error(fatal, "continuity_init: Unrecognized setting "// &
144  "#define CONTINUITY_SCHEME "//trim(tmpstr)//" found in input file.")
145  end select
146 
147  if (cs%continuity_scheme == ppm_scheme) then
148  call continuity_ppm_init(time, g, gv, us, param_file, diag, cs%PPM_CSp)
149  endif
150 

References mom_error_handler::mom_error(), mom_error_handler::mom_mesg(), ppm_scheme, ppm_string, and mom_string_functions::uppercase().

Referenced by mom_dynamics_unsplit::initialize_dyn_unsplit(), and mom_dynamics_unsplit_rk2::initialize_dyn_unsplit_rk2().

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

◆ continuity_stencil()

integer function, public mom_continuity::continuity_stencil ( type(continuity_cs), pointer  CS)

continuity_stencil returns the continuity solver stencil size

Parameters
csModule's control structure.
Returns
The continuity solver stencil size with the current settings.

Definition at line 156 of file MOM_continuity.F90.

156  type(continuity_CS), pointer :: CS !< Module's control structure.
157  integer :: stencil !< The continuity solver stencil size with the current settings.
158 
159  stencil = 1
160 
161  if (cs%continuity_scheme == ppm_scheme) then
162  stencil = continuity_ppm_stencil(cs%PPM_CSp)
163  endif
164 

References ppm_scheme.

Referenced by mom_dynamics_split_rk2::step_mom_dyn_split_rk2().

Here is the caller graph for this function:

Variable Documentation

◆ ppm_scheme

integer, parameter mom_continuity::ppm_scheme = 1

Enumerated constant to select PPM.

Definition at line 35 of file MOM_continuity.F90.

35 integer, parameter :: PPM_SCHEME = 1 !< Enumerated constant to select PPM

Referenced by continuity(), continuity_end(), continuity_init(), and continuity_stencil().

◆ ppm_string

character(len=20), parameter mom_continuity::ppm_string = "PPM"
private

String to select PPM.

Definition at line 36 of file MOM_continuity.F90.

36 character(len=20), parameter :: PPM_STRING = "PPM" !< String to select PPM

Referenced by continuity_init().