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...
 

Functions/Subroutines

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...
 

Variables

character(len=40) mdl = "MOM_CVMix_conv"
 This module's name.
 

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( g %isd: g %ied, g %jsd: g %jed), intent(in)  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.
[in]hblDepth of ocean boundary layer [Z ~> m]

Definition at line 153 of file MOM_CVMix_conv.F90.

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

◆ 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 275 of file MOM_CVMix_conv.F90.

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

◆ 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 get_param(param_file, mdl, "USE_CVMix_CONVECTION", cvmix_conv_init, default=.false., do_not_log=.true.)
80  call log_version(param_file, mdl, version, &
81  "Parameterization of enhanced mixing due to convection via CVMix", &
82  all_default=.not.cvmix_conv_init)
83  call get_param(param_file, mdl, "USE_CVMix_CONVECTION", cvmix_conv_init, &
84  "If true, turns on the enhanced mixing due to convection "//&
85  "via CVMix. This scheme increases diapycnal diffs./viscs. "//&
86  "at statically unstable interfaces. Relevant parameters are "//&
87  "contained in the CVMix_CONVECTION% parameter block.", &
88  default=.false.)
89 
90  if (.not. cvmix_conv_init) return
91 
92  call get_param(param_file, mdl, "ENERGETICS_SFC_PBL", useepbl, default=.false., &
93  do_not_log=.true.)
94 
95  ! Warn user if EPBL is being used, since in this case mixing due to convection will
96  ! be aplied in the boundary layer
97  if (useepbl) then
98  call mom_error(warning, 'MOM_CVMix_conv_init: '// &
99  'CVMix convection may not be properly applied when ENERGETICS_SFC_PBL = True'//&
100  'as convective mixing might occur in the boundary layer.')
101  endif
102 
103  call get_param(param_file, mdl, 'DEBUG', cs%debug, default=.false., do_not_log=.true.)
104 
105  call get_param(param_file, mdl, 'MIN_THICKNESS', cs%min_thickness, default=0.001, do_not_log=.true.)
106 
107  call openparameterblock(param_file,'CVMix_CONVECTION')
108 
109  call get_param(param_file, mdl, "PRANDTL_CONV", prandtl_conv, &
110  "The turbulent Prandtl number applied to convective "//&
111  "instabilities (i.e., used to convert KD_CONV into KV_CONV)", &
112  units="nondim", default=1.0)
113 
114  call get_param(param_file, mdl, 'KD_CONV', cs%kd_conv_const, &
115  "Diffusivity used in convective regime. Corresponding viscosity "//&
116  "(KV_CONV) will be set to KD_CONV * PRANDTL_TURB.", &
117  units='m2/s', default=1.00)
118 
119  call get_param(param_file, mdl, 'BV_SQR_CONV', cs%bv_sqr_conv, &
120  "Threshold for squared buoyancy frequency needed to trigger "//&
121  "Brunt-Vaisala parameterization.", &
122  units='1/s^2', default=0.0)
123 
124  call closeparameterblock(param_file)
125 
126  ! set kv_conv_const based on kd_conv_const and prandtl_conv
127  cs%kv_conv_const = cs%kd_conv_const * prandtl_conv
128 
129  ! allocate arrays and set them to zero
130  allocate(cs%N2(szi_(g), szj_(g), szk_(g)+1)); cs%N2(:,:,:) = 0.
131  allocate(cs%kd_conv(szi_(g), szj_(g), szk_(g)+1)); cs%kd_conv(:,:,:) = 0.
132  allocate(cs%kv_conv(szi_(g), szj_(g), szk_(g)+1)); cs%kv_conv(:,:,:) = 0.
133 
134  ! Register diagnostics
135  cs%diag => diag
136  cs%id_N2 = register_diag_field('ocean_model', 'N2_conv', diag%axesTi, time, &
137  'Square of Brunt-Vaisala frequency used by MOM_CVMix_conv module', '1/s2')
138  cs%id_kd_conv = register_diag_field('ocean_model', 'kd_conv', diag%axesTi, time, &
139  'Additional diffusivity added by MOM_CVMix_conv module', 'm2/s', conversion=us%Z2_T_to_m2_s)
140  cs%id_kv_conv = register_diag_field('ocean_model', 'kv_conv', diag%axesTi, time, &
141  'Additional viscosity added by MOM_CVMix_conv module', 'm2/s', conversion=us%Z2_T_to_m2_s)
142 
143  call cvmix_init_conv(convect_diff=cs%kd_conv_const, &
144  convect_visc=cs%kv_conv_const, &
145  lbruntvaisala=.true., &
146  bvsqr_convect=cs%bv_sqr_conv)
147 

◆ 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 267 of file MOM_CVMix_conv.F90.

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