MOM6
mom_cvmix_conv Module Reference

Detailed Description

Interface to CVMix convection scheme.

Data Types

type  cvmix_conv_cs
 Control structure including parameters for CVMix convection. More...
 
character(len=40) mdl = "MOM_CVMix_conv"
 This module's name. More...
 
logical function, public cvmix_conv_init (Time, G, GV, US, param_file, diag, CS)
 Initialized the CVMix convection mixing routine. More...
 
subroutine, public calculate_cvmix_conv (h, tv, G, GV, US, CS, hbl)
 Subroutine for calculating enhanced diffusivity/viscosity due to convection via CVMix. More...
 
logical function, public cvmix_conv_is_used (param_file)
 Reads the parameter "USE_CVMix_CONVECTION" and returns state. This function allows other modules to know whether this parameterization will be used without needing to duplicate the log entry. More...
 
subroutine, public cvmix_conv_end (CS)
 Clear pointers and dealocate memory. More...
 

Function/Subroutine Documentation

◆ calculate_cvmix_conv()

subroutine, public mom_cvmix_conv::calculate_cvmix_conv ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(cvmix_conv_cs), pointer  CS,
real, dimension(:,:), optional, pointer  hbl 
)

Subroutine for calculating enhanced diffusivity/viscosity due to convection via CVMix.

Parameters
[in]gGrid structure.
[in]gvVertical grid structure.
[in]usA dimensional unit scaling type
[in]hLayer thickness [H ~> m or kg m-2].
[in]tvThermodynamics structure.
csThe control structure returned by a previous call to CVMix_conv_init.
hblDepth of ocean boundary layer [m]

Definition at line 151 of file MOM_CVMix_conv.F90.

151 
152  type(ocean_grid_type), intent(in) :: G !< Grid structure.
153  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
154  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
155  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
156  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure.
157  type(CVMix_conv_cs), pointer :: CS !< The control structure returned
158  !! by a previous call to CVMix_conv_init.
159  real, dimension(:,:), optional, pointer :: hbl!< Depth of ocean boundary layer [m]
160  ! local variables
161  real, dimension(SZK_(G)) :: rho_lwr !< Adiabatic Water Density, this is a dummy
162  !! variable since here convection is always
163  !! computed based on Brunt Vaisala.
164  real, dimension(SZK_(G)) :: rho_1d !< water density in a column, this is also
165  !! a dummy variable, same reason as above.
166  real, dimension(SZK_(G)+1) :: kv_col !< Viscosities at interfaces in the column [m2 s-1]
167  real, dimension(SZK_(G)+1) :: kd_col !< Diffusivities at interfaces in the column [m2 s-1]
168  real, dimension(SZK_(G)+1) :: iFaceHeight !< Height of interfaces [m]
169  real, dimension(SZK_(G)) :: cellHeight !< Height of cell centers [m]
170  integer :: kOBL !< level of OBL extent
171  real :: g_o_rho0 ! Gravitational acceleration divided by density in MKS units [m4 s-2]
172  real :: pref, rhok, rhokm1, dz, dh, hcorr
173  integer :: i, j, k
174 
175  g_o_rho0 = gv%mks_g_Earth / (us%R_to_kg_m3*gv%Rho0)
176 
177  ! initialize dummy variables
178  rho_lwr(:) = 0.0; rho_1d(:) = 0.0
179 
180  if (.not. associated(hbl)) then
181  allocate(hbl(szi_(g), szj_(g)))
182  hbl(:,:) = 0.0
183  endif
184 
185  do j = g%jsc, g%jec
186  do i = g%isc, g%iec
187 
188  ! set N2 to zero at the top- and bottom-most interfaces
189  cs%N2(i,j,1) = 0.
190  cs%N2(i,j,g%ke+1) = 0.
191 
192  ! skip calling at land points
193  !if (G%mask2dT(i,j) == 0.) cycle
194 
195  pref = 0.
196  ! Compute Brunt-Vaisala frequency (static stability) on interfaces
197  do k=2,g%ke
198 
199  ! pRef is pressure at interface between k and km1.
200  pref = pref + gv%H_to_Pa * h(i,j,k)
201  call calculate_density(tv%t(i,j,k), tv%s(i,j,k), pref, rhok, tv%eqn_of_state)
202  call calculate_density(tv%t(i,j,k-1), tv%s(i,j,k-1), pref, rhokm1, tv%eqn_of_state)
203 
204  dz = ((0.5*(h(i,j,k-1) + h(i,j,k))+gv%H_subroundoff)*gv%H_to_m)
205  cs%N2(i,j,k) = g_o_rho0 * (rhok - rhokm1) / dz ! Can be negative
206 
207  enddo
208 
209  ifaceheight(1) = 0.0 ! BBL is all relative to the surface
210  hcorr = 0.0
211  ! compute heights at cell center and interfaces
212  do k=1,g%ke
213  dh = h(i,j,k) * gv%H_to_m ! Nominal thickness to use for increment
214  dh = dh + hcorr ! Take away the accumulated error (could temporarily make dh<0)
215  hcorr = min( dh - cs%min_thickness, 0. ) ! If inflating then hcorr<0
216  dh = max( dh, cs%min_thickness ) ! Limit increment dh>=min_thickness
217  cellheight(k) = ifaceheight(k) - 0.5 * dh
218  ifaceheight(k+1) = ifaceheight(k) - dh
219  enddo
220 
221  ! gets index of the level and interface above hbl
222  kobl = cvmix_kpp_compute_kobl_depth(ifaceheight, cellheight,hbl(i,j))
223 
224  kv_col(:) = 0.0 ; kd_col(:) = 0.0
225  call cvmix_coeffs_conv(mdiff_out=kv_col(:), &
226  tdiff_out=kd_col(:), &
227  nsqr=cs%N2(i,j,:), &
228  dens=rho_1d(:), &
229  dens_lwr=rho_lwr(:), &
230  nlev=g%ke, &
231  max_nlev=g%ke, &
232  obl_ind=kobl)
233 
234  do k=1,g%ke+1
235  cs%kv_conv(i,j,k) = us%m2_s_to_Z2_T * kv_col(k)
236  cs%Kd_conv(i,j,k) = us%m2_s_to_Z2_T * kd_col(k)
237  enddo
238  ! Do not apply mixing due to convection within the boundary layer
239  do k=1,kobl
240  cs%kv_conv(i,j,k) = 0.0
241  cs%kd_conv(i,j,k) = 0.0
242  enddo
243 
244  enddo
245  enddo
246 
247  if (cs%debug) then
248  call hchksum(cs%N2, "MOM_CVMix_conv: N2",g%HI,haloshift=0)
249  call hchksum(cs%kd_conv, "MOM_CVMix_conv: kd_conv",g%HI,haloshift=0,scale=us%Z2_T_to_m2_s)
250  call hchksum(cs%kv_conv, "MOM_CVMix_conv: kv_conv",g%HI,haloshift=0,scale=us%m2_s_to_Z2_T)
251  endif
252 
253  ! send diagnostics to post_data
254  if (cs%id_N2 > 0) call post_data(cs%id_N2, cs%N2, cs%diag)
255  if (cs%id_kd_conv > 0) call post_data(cs%id_kd_conv, cs%kd_conv, cs%diag)
256  if (cs%id_kv_conv > 0) call post_data(cs%id_kv_conv, cs%kv_conv, cs%diag)
257 

Referenced by mom_diabatic_driver::diabatic_ale(), mom_diabatic_driver::diabatic_ale_legacy(), and mom_diabatic_driver::layered_diabatic().

Here is the caller graph for this function:

◆ cvmix_conv_end()

subroutine, public mom_cvmix_conv::cvmix_conv_end ( type(cvmix_conv_cs), pointer  CS)

Clear pointers and dealocate memory.

Parameters
csControl structure for this module that will be deallocated in this subroutine

Definition at line 272 of file MOM_CVMix_conv.F90.

272  type(CVMix_conv_cs), pointer :: CS !< Control structure for this module that
273  !! will be deallocated in this subroutine
274 
275  if (.not. associated(cs)) return
276 
277  deallocate(cs%N2)
278  deallocate(cs%kd_conv)
279  deallocate(cs%kv_conv)
280  deallocate(cs)
281 

Referenced by mom_diabatic_driver::diabatic_driver_end().

Here is the caller graph for this function:

◆ cvmix_conv_init()

logical function, public mom_cvmix_conv::cvmix_conv_init ( type(time_type), intent(in)  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(cvmix_conv_cs), pointer  CS 
)

Initialized the CVMix convection mixing routine.

Parameters
[in]timeThe current time.
[in]gGrid structure.
[in]gvVertical grid structure.
[in]usA dimensional unit scaling type
[in]param_fileRun-time parameter file handle
[in,out]diagDiagnostics control structure.
csThis module's control structure.

Definition at line 56 of file MOM_CVMix_conv.F90.

56 
57  type(time_type), intent(in) :: Time !< The current time.
58  type(ocean_grid_type), intent(in) :: G !< Grid structure.
59  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
60  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
61  type(param_file_type), intent(in) :: param_file !< Run-time parameter file handle
62  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
63  type(CVMix_conv_cs), pointer :: CS !< This module's control structure.
64  ! Local variables
65  real :: prandtl_conv !< Turbulent Prandtl number used in convective instabilities.
66  logical :: useEPBL !< If True, use the ePBL boundary layer scheme.
67 
68 ! This include declares and sets the variable "version".
69 #include "version_variable.h"
70 
71  if (associated(cs)) then
72  call mom_error(warning, "CVMix_conv_init called with an associated "// &
73  "control structure.")
74  return
75  endif
76  allocate(cs)
77 
78  ! Read parameters
79  call log_version(param_file, mdl, version, &
80  "Parameterization of enhanced mixing due to convection via CVMix")
81  call get_param(param_file, mdl, "USE_CVMix_CONVECTION", cvmix_conv_init, &
82  "If true, turns on the enhanced mixing due to convection "//&
83  "via CVMix. This scheme increases diapycnal diffs./viscs. "//&
84  "at statically unstable interfaces. Relevant parameters are "//&
85  "contained in the CVMix_CONVECTION% parameter block.", &
86  default=.false.)
87 
88  if (.not. cvmix_conv_init) return
89 
90  call get_param(param_file, mdl, "ENERGETICS_SFC_PBL", useepbl, default=.false., &
91  do_not_log=.true.)
92 
93  ! Warn user if EPBL is being used, since in this case mixing due to convection will
94  ! be aplied in the boundary layer
95  if (useepbl) then
96  call mom_error(warning, 'MOM_CVMix_conv_init: '// &
97  'CVMix convection may not be properly applied when ENERGETICS_SFC_PBL = True'//&
98  'as convective mixing might occur in the boundary layer.')
99  endif
100 
101  call get_param(param_file, mdl, 'DEBUG', cs%debug, default=.false., do_not_log=.true.)
102 
103  call get_param(param_file, mdl, 'MIN_THICKNESS', cs%min_thickness, default=0.001, do_not_log=.true.)
104 
105  call openparameterblock(param_file,'CVMix_CONVECTION')
106 
107  call get_param(param_file, mdl, "PRANDTL_CONV", prandtl_conv, &
108  "The turbulent Prandtl number applied to convective "//&
109  "instabilities (i.e., used to convert KD_CONV into KV_CONV)", &
110  units="nondim", default=1.0)
111 
112  call get_param(param_file, mdl, 'KD_CONV', cs%kd_conv_const, &
113  "Diffusivity used in convective regime. Corresponding viscosity "//&
114  "(KV_CONV) will be set to KD_CONV * PRANDTL_TURB.", &
115  units='m2/s', default=1.00)
116 
117  call get_param(param_file, mdl, 'BV_SQR_CONV', cs%bv_sqr_conv, &
118  "Threshold for squared buoyancy frequency needed to trigger "//&
119  "Brunt-Vaisala parameterization.", &
120  units='1/s^2', default=0.0)
121 
122  call closeparameterblock(param_file)
123 
124  ! set kv_conv_const based on kd_conv_const and prandtl_conv
125  cs%kv_conv_const = cs%kd_conv_const * prandtl_conv
126 
127  ! allocate arrays and set them to zero
128  allocate(cs%N2(szi_(g), szj_(g), szk_(g)+1)); cs%N2(:,:,:) = 0.
129  allocate(cs%kd_conv(szi_(g), szj_(g), szk_(g)+1)); cs%kd_conv(:,:,:) = 0.
130  allocate(cs%kv_conv(szi_(g), szj_(g), szk_(g)+1)); cs%kv_conv(:,:,:) = 0.
131 
132  ! Register diagnostics
133  cs%diag => diag
134  cs%id_N2 = register_diag_field('ocean_model', 'N2_conv', diag%axesTi, time, &
135  'Square of Brunt-Vaisala frequency used by MOM_CVMix_conv module', '1/s2')
136  cs%id_kd_conv = register_diag_field('ocean_model', 'kd_conv', diag%axesTi, time, &
137  'Additional diffusivity added by MOM_CVMix_conv module', 'm2/s', conversion=us%Z2_T_to_m2_s)
138  cs%id_kv_conv = register_diag_field('ocean_model', 'kv_conv', diag%axesTi, time, &
139  'Additional viscosity added by MOM_CVMix_conv module', 'm2/s', conversion=us%Z2_T_to_m2_s)
140 
141  call cvmix_init_conv(convect_diff=cs%kd_conv_const, &
142  convect_visc=cs%kv_conv_const, &
143  lbruntvaisala=.true., &
144  bvsqr_convect=cs%bv_sqr_conv)
145 

References mdl, and mom_error_handler::mom_error().

Here is the call graph for this function:

◆ cvmix_conv_is_used()

logical function, public mom_cvmix_conv::cvmix_conv_is_used ( type(param_file_type), intent(in)  param_file)

Reads the parameter "USE_CVMix_CONVECTION" and returns state. This function allows other modules to know whether this parameterization will be used without needing to duplicate the log entry.

Parameters
[in]param_fileA structure to parse for run-time parameters

Definition at line 264 of file MOM_CVMix_conv.F90.

264  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
265  call get_param(param_file, mdl, "USE_CVMix_CONVECTION", cvmix_conv_is_used, &
266  default=.false., do_not_log = .true.)
267 

References mdl.

Referenced by mom_set_visc::set_visc_register_restarts().

Here is the caller graph for this function:

Variable Documentation

◆ mdl

character(len=40) mom_cvmix_conv::mdl = "MOM_CVMix_conv"

This module's name.

Definition at line 50 of file MOM_CVMix_conv.F90.

50 character(len=40) :: mdl = "MOM_CVMix_conv" !< This module's name.

Referenced by cvmix_conv_init(), and cvmix_conv_is_used().