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, CS, uhbt, vhbt, OBC, visc_rem_u, visc_rem_v, u_cor, v_cor, uhbt_aux, vhbt_aux, u_cor_aux, v_cor_aux, 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, 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(inout)  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(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,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed), intent(in), optional  uhbt_aux,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb), intent(in), optional  vhbt_aux,
real, dimension( g %isdb: g %iedb, g %jsd: g %jed, g %ke), intent(inout), optional  u_cor_aux,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, g %ke), intent(inout), optional  v_cor_aux,
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, in m/s.
[in]vMeridional velocity, in m/s.
[in,out]hinInitial layer thickness, in m or kg/m2.
[in,out]hFinal layer thickness, in m or kg/m2.
[out]uhVolume flux through zonal faces = u*h*dy, in m3/s.
[out]vhVolume flux through meridional faces = v*h*dx, in m3/s.
[in]dtTime increment, in s.
csControl structure for mom_continuity.
[in]uhbtThe vertically summed volume flux through zonal faces, in m3/s.
[in]vhbtThe vertically summed volume flux through meridional faces, in m3/s.
obcOpen boundaries control structure.
[in]visc_rem_uBoth the fraction of zonal momentum that remains after a time-step of viscosity, and the fraction of a time-step's worth of a barotropic acceleration that a layer experiences after viscosity is applied. Non-dimensional between 0 (at the bottom) and 1 (far above the bottom).
[in]visc_rem_vBoth the fraction of meridional momentum that remains after a time-step of viscosity, and the fraction of a time-step's worth of a barotropic acceleration that a layer experiences after viscosity is applied. Non-dimensional between 0 (at the bottom) and 1 (far above the bottom).
[out]u_corThe zonal velocities that give uhbt as the depth-integrated transport, in m/s.
[out]v_corThe meridional velocities that give vhbt as the depth-integrated transport, in m/s.
[in]uhbt_auxA second summed zonal volume flux in m3/s.
[in]vhbt_auxA second summed meridional volume flux in m3/s.
[in,out]u_cor_auxThe zonal velocities that give uhbt_aux as the depth-integrated transport, in m/s.
[in,out]v_cor_auxThe meridional velocities that give vhbt_aux as the depth-integrated transport, in m/s.
bt_contA structure with elements that describe the effective open face areas as a function of barotropic flow.

Definition at line 44 of file MOM_continuity.F90.

References mom_error_handler::mom_error(), and ppm_scheme.

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

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)), intent(in) :: u !< Zonal velocity, in m/s.
47  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(in) :: v !< Meridional velocity, in m/s.
48  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: hin !< Initial layer thickness, in m or kg/m2.
49  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< Final layer thickness, in m or kg/m2.
50  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: uh !< Volume flux through zonal faces =
51  !! u*h*dy, in m3/s.
52  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: vh !< Volume flux through meridional faces =
53  !! v*h*dx, in m3/s.
54  real, intent(in) :: dt !< Time increment, in s.
55  type(continuity_cs), pointer :: cs !< Control structure for mom_continuity.
56  real, dimension(SZIB_(G),SZJ_(G)), intent(in), optional :: uhbt !< The vertically summed volume
57  !! flux through zonal faces, in m3/s.
58  real, dimension(SZI_(G),SZJB_(G)), intent(in), optional :: vhbt !< The vertically summed volume
59  !! flux through meridional faces, in m3/s.
60  type(ocean_obc_type), pointer, optional :: obc !< Open boundaries control structure.
61  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(in), optional :: visc_rem_u !< Both the fraction of
62  !! zonal momentum that remains after a time-step of viscosity, and the fraction of a time-step's
63  !! worth of a barotropic acceleration that a layer experiences after viscosity is applied.
64  !! Non-dimensional between 0 (at the bottom) and 1 (far above the bottom).
65  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(in), optional :: visc_rem_v !< Both the fraction of
66  !! meridional momentum that remains after a time-step of viscosity, and the fraction of a time-step's
67  !! worth of a barotropic acceleration that a layer experiences after viscosity is applied.
68  !! Non-dimensional between 0 (at the bottom) and 1 (far above the bottom).
69  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out), optional :: u_cor !< The zonal velocities that
70  !! give uhbt as the depth-integrated transport, in m/s.
71  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out), optional :: v_cor !< The meridional velocities that
72  !! give vhbt as the depth-integrated transport, in m/s.
73  real, dimension(SZIB_(G),SZJ_(G)), intent(in), optional :: uhbt_aux !< A second summed zonal
74  !! volume flux in m3/s.
75  real, dimension(SZI_(G),SZJB_(G)), intent(in), optional :: vhbt_aux !< A second summed meridional
76  !! volume flux in m3/s.
77  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(inout), optional :: u_cor_aux !< The zonal velocities
78  !! that give uhbt_aux as the depth-integrated transport, in m/s.
79  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(inout), optional :: v_cor_aux !< The meridional velocities
80  !! that give vhbt_aux as the depth-integrated transport, in m/s.
81  type(bt_cont_type), pointer, optional :: bt_cont !< A structure with elements
82  !! that describe the effective open face areas as a function of barotropic flow.
83 
84  if (present(visc_rem_u) .neqv. present(visc_rem_v)) call mom_error(fatal, &
85  "MOM_continuity: Either both visc_rem_u and visc_rem_v or neither"// &
86  " one must be present in call to continuity.")
87  if (present(u_cor) .neqv. present(v_cor)) call mom_error(fatal, &
88  "MOM_continuity: Either both u_cor and v_cor or neither"// &
89  " one must be present in call to continuity.")
90  if (present(uhbt_aux) .neqv. present(vhbt_aux)) call mom_error(fatal, &
91  "MOM_continuity: Either both uhbt_aux and uhbt_aux or neither"// &
92  " one must be present in call to continuity.")
93  if (present(u_cor_aux) .neqv. present(v_cor_aux)) call mom_error(fatal, &
94  "MOM_continuity: Either both u_cor_aux and v_cor_aux or neither"// &
95  " one must be present in call to continuity.")
96  if (present(u_cor_aux) .neqv. present(uhbt_aux)) call mom_error(fatal, &
97  "MOM_continuity: u_cor_aux can only be calculated if uhbt_aux is"// &
98  " provided, and uhbt_aux has no other purpose. Include both arguments"//&
99  " or neither.")
100 
101  if (cs%continuity_scheme == ppm_scheme) then
102  call continuity_ppm(u, v, hin, h, uh, vh, dt, g, gv, cs%PPM_CSp, uhbt, vhbt, obc, &
103  visc_rem_u, visc_rem_v, u_cor, v_cor, &
104  uhbt_aux, vhbt_aux, u_cor_aux, v_cor_aux, bt_cont)
105  else
106  call mom_error(fatal, "continuity: Unrecognized value of continuity_scheme")
107  endif
108 
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 172 of file MOM_continuity.F90.

References mom_continuity_ppm::continuity_ppm_end(), and ppm_scheme.

172  type(continuity_cs), pointer :: cs !< Control structure for mom_continuity.
173 
174  if (cs%continuity_scheme == ppm_scheme) then
175  call continuity_ppm_end(cs%PPM_CSp)
176  endif
177 
178  deallocate(cs)
179 
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(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]param_fileParameter file handles.
[in,out]diagDiagnostics control structure.
csControl structure for mom_continuity.

Definition at line 113 of file MOM_continuity.F90.

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

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

References ppm_scheme.

Referenced by mom_dynamics_split_rk2::step_mom_dyn_split_rk2().

159  type(continuity_cs), pointer :: cs !< Module's control structure.
160  integer :: stencil !< The continuity solver stencil size with the current settings.
161 
162  stencil = 1
163 
164  if (cs%continuity_scheme == ppm_scheme) then
165  stencil = continuity_ppm_stencil(cs%PPM_CSp)
166  endif
167 
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 34 of file MOM_continuity.F90.

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

34 integer, parameter :: ppm_scheme = 1 !< Enumerated constant to select PPM

◆ ppm_string

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

String to select PPM.

Definition at line 35 of file MOM_continuity.F90.

Referenced by continuity_init().

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