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, 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.
 
character(len=20), parameter ppm_string = "PPM"
 String to select PPM.
 

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,
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 [m s-1].
[in]vMeridional velocity [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 [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
[in]uhbt_auxA second summed zonal
[in]vhbt_auxA second summed meridional
[in,out]u_cor_auxThe zonal velocities
[in,out]v_cor_auxThe meridional velocities
bt_contA structure with elements

Definition at line 45 of file MOM_continuity.F90.

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

◆ 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 192 of file MOM_continuity.F90.

192  type(continuity_CS), pointer :: CS !< Control structure for mom_continuity.
193 
194  if (cs%continuity_scheme == ppm_scheme) then
195  call continuity_ppm_end(cs%PPM_CSp)
196  endif
197 
198  deallocate(cs)
199 

◆ 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 133 of file MOM_continuity.F90.

133  type(time_type), target, intent(in) :: Time !< Current model time.
134  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
135  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
136  type(param_file_type), intent(in) :: param_file !< Parameter file handles.
137  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
138  type(continuity_CS), pointer :: CS !< Control structure for mom_continuity.
139 ! This include declares and sets the variable "version".
140 #include "version_variable.h"
141  character(len=40) :: mdl = "MOM_continuity" ! This module's name.
142  character(len=20) :: tmpstr
143 
144  if (associated(cs)) then
145  call mom_error(warning, "continuity_init called with associated control structure.")
146  return
147  endif
148  allocate(cs)
149 
150  ! Read all relevant parameters and write them to the model log.
151  call log_version(param_file, mdl, version, "")
152  call get_param(param_file, mdl, "CONTINUITY_SCHEME", tmpstr, &
153  "CONTINUITY_SCHEME selects the discretization for the "//&
154  "continuity solver. The only valid value currently is: \n"//&
155  "\t PPM - use a positive-definite (or monotonic) \n"//&
156  "\t piecewise parabolic reconstruction solver.", &
157  default=ppm_string)
158 
159  tmpstr = uppercase(tmpstr) ; cs%continuity_scheme = 0
160  select case (trim(tmpstr))
161  case (ppm_string) ; cs%continuity_scheme = ppm_scheme
162  case default
163  call mom_mesg('continuity_init: CONTINUITY_SCHEME ="'//trim(tmpstr)//'"', 0)
164  call mom_mesg("continuity_init: The only valid value is currently "// &
165  trim(ppm_string), 0)
166  call mom_error(fatal, "continuity_init: Unrecognized setting "// &
167  "#define CONTINUITY_SCHEME "//trim(tmpstr)//" found in input file.")
168  end select
169 
170  if (cs%continuity_scheme == ppm_scheme) then
171  call continuity_ppm_init(time, g, gv, param_file, diag, cs%PPM_CSp)
172  endif
173 

◆ 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 179 of file MOM_continuity.F90.

179  type(continuity_CS), pointer :: CS !< Module's control structure.
180  integer :: stencil !< The continuity solver stencil size with the current settings.
181 
182  stencil = 1
183 
184  if (cs%continuity_scheme == ppm_scheme) then
185  stencil = continuity_ppm_stencil(cs%PPM_CSp)
186  endif
187