MOM6
mom_cvmix_shear Module Reference

Detailed Description

Interface to CVMix interior shear schemes.

Data Types

type  cvmix_shear_cs
 Control structure including parameters for CVMix interior shear schemes. More...
 
character(len=40) mdl = "MOM_CVMix_shear"
 This module's name.
 
subroutine, public calculate_cvmix_shear (u_H, v_H, h, tv, kd, kv, G, GV, US, CS)
 Subroutine for calculating (internal) vertical diffusivities/viscosities. More...
 
logical function, public cvmix_shear_init (Time, G, GV, US, param_file, diag, CS)
 Initialized the CVMix internal shear mixing routine. More...
 
logical function, public cvmix_shear_is_used (param_file)
 Reads the parameters "LMD94" and "PP81" 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_shear_end (CS)
 Clear pointers and dealocate memory. More...
 

Function/Subroutine Documentation

◆ calculate_cvmix_shear()

subroutine, public mom_cvmix_shear::calculate_cvmix_shear ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  u_H,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  v_H,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
type(thermo_var_ptrs), intent(in)  tv,
real, dimension(szi_(g),szj_(g),szk_(g)+1), intent(out)  kd,
real, dimension(szi_(g),szj_(g),szk_(g)+1), intent(out)  kv,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(cvmix_shear_cs), pointer  CS 
)

Subroutine for calculating (internal) vertical diffusivities/viscosities.

Parameters
[in]gGrid structure.
[in]gvVertical grid structure.
[in]usA dimensional unit scaling type
[in]u_hInitial zonal velocity on T points [m s-1].
[in]v_hInitial meridional velocity on T points [m s-1].
[in]hLayer thickness [H ~> m or kg m-2].
[in]tvThermodynamics structure.
[out]kdThe vertical diffusivity at each interface (not layer!) [Z2 T-1 ~> m2 s-1].
[out]kvThe vertical viscosity at each interface (not layer!) [Z2 T-1 ~> m2 s-1].
csThe control structure returned by a previous call to CVMix_shear_init.

Definition at line 60 of file MOM_CVMix_shear.F90.

60  type(ocean_grid_type), intent(in) :: G !< Grid structure.
61  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
62  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
63  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: u_H !< Initial zonal velocity on T points [m s-1].
64  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: v_H !< Initial meridional velocity on T points [m s-1].
65  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
66  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure.
67  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(out) :: kd !< The vertical diffusivity at each interface
68  !! (not layer!) [Z2 T-1 ~> m2 s-1].
69  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(out) :: kv !< The vertical viscosity at each interface
70  !! (not layer!) [Z2 T-1 ~> m2 s-1].
71  type(CVMix_shear_cs), pointer :: CS !< The control structure returned by a previous call to
72  !! CVMix_shear_init.
73  ! Local variables
74  integer :: i, j, k, kk, km1
75  real :: GoRho ! Gravitational acceleration divided by density in MKS units [m4 s-2]
76  real :: pref, DU, DV, dRho, DZ, N2, S2, dummy
77  real, dimension(2*(G%ke)) :: pres_1d, temp_1d, salt_1d, rho_1d
78  real, dimension(G%ke+1) :: Ri_Grad !< Gradient Richardson number
79  real, dimension(G%ke+1) :: Kvisc !< Vertical viscosity at interfaces [m2 s-1]
80  real, dimension(G%ke+1) :: Kdiff !< Diapycnal diffusivity at interfaces [m2 s-1]
81  real, parameter :: epsln = 1.e-10 !< Threshold to identify vanished layers
82 
83  ! some constants
84  gorho = gv%mks_g_Earth / gv%Rho0
85 
86  do j = g%jsc, g%jec
87  do i = g%isc, g%iec
88 
89  ! skip calling for land points
90  if (g%mask2dT(i,j)==0.) cycle
91 
92  ! Richardson number computed for each cell in a column.
93  pref = 0.
94  ri_grad(:)=1.e8 !Initialize w/ large Richardson value
95  do k=1,g%ke
96  ! pressure, temp, and saln for EOS
97  ! kk+1 = k fields
98  ! kk+2 = km1 fields
99  km1 = max(1, k-1)
100  kk = 2*(k-1)
101  pres_1d(kk+1) = pref
102  pres_1d(kk+2) = pref
103  temp_1d(kk+1) = tv%T(i,j,k)
104  temp_1d(kk+2) = tv%T(i,j,km1)
105  salt_1d(kk+1) = tv%S(i,j,k)
106  salt_1d(kk+2) = tv%S(i,j,km1)
107 
108  ! pRef is pressure at interface between k and km1.
109  ! iterate pRef for next pass through k-loop.
110  pref = pref + gv%H_to_Pa * h(i,j,k)
111 
112  enddo ! k-loop finishes
113 
114  ! compute in-situ density
115  call calculate_density(temp_1d, salt_1d, pres_1d, rho_1d, 1, 2*g%ke, tv%EQN_OF_STATE)
116 
117  ! N2 (can be negative) on interface
118  do k = 1, g%ke
119  km1 = max(1, k-1)
120  kk = 2*(k-1)
121  du = (u_h(i,j,k))-(u_h(i,j,km1))
122  dv = (v_h(i,j,k))-(v_h(i,j,km1))
123  drho = (gorho * (rho_1d(kk+1) - rho_1d(kk+2)) )
124  dz = ((0.5*(h(i,j,km1) + h(i,j,k))+gv%H_subroundoff)*gv%H_to_m)
125  n2 = drho/dz
126  s2 = (du*du+dv*dv)/(dz*dz)
127  ri_grad(k) = max(0.,n2)/max(s2,1.e-10)
128 
129  ! fill 3d arrays, if user asks for diagsnostics
130  if (cs%id_N2 > 0) cs%N2(i,j,k) = n2
131  if (cs%id_S2 > 0) cs%S2(i,j,k) = s2
132 
133  enddo
134 
135  ri_grad(g%ke+1) = ri_grad(g%ke)
136 
137  if (cs%id_ri_grad > 0) cs%ri_grad(i,j,:) = ri_grad(:)
138 
139  if (cs%smooth_ri) then
140  ! 1) fill Ri_grad in vanished layers with adjacent value
141  do k = 2, g%ke
142  if (h(i,j,k) .le. epsln) ri_grad(k) = ri_grad(k-1)
143  enddo
144 
145  ri_grad(g%ke+1) = ri_grad(g%ke)
146 
147  ! 2) vertically smooth Ri with 1-2-1 filter
148  dummy = 0.25 * ri_grad(2)
149  ri_grad(g%ke+1) = ri_grad(g%ke)
150  do k = 3, g%ke
151  ri_grad(k) = dummy + 0.5 * ri_grad(k) + 0.25 * ri_grad(k+1)
152  dummy = 0.25 * ri_grad(k)
153  enddo
154 
155  if (cs%id_ri_grad_smooth > 0) cs%ri_grad_smooth(i,j,:) = ri_grad(:)
156  endif
157 
158  do k=1,g%ke+1
159  kvisc(k) = us%Z2_T_to_m2_s * kv(i,j,k)
160  kdiff(k) = us%Z2_T_to_m2_s * kd(i,j,k)
161  enddo
162 
163  ! Call to CVMix wrapper for computing interior mixing coefficients.
164  call cvmix_coeffs_shear(mdiff_out=kvisc(:), &
165  tdiff_out=kdiff(:), &
166  rich=ri_grad(:), &
167  nlev=g%ke, &
168  max_nlev=g%ke)
169  do k=1,g%ke+1
170  kv(i,j,k) = us%m2_s_to_Z2_T * kvisc(k)
171  kd(i,j,k) = us%m2_s_to_Z2_T * kdiff(k)
172  enddo
173  enddo
174  enddo
175 
176  ! write diagnostics
177  if (cs%id_kd > 0) call post_data(cs%id_kd, kd, cs%diag)
178  if (cs%id_kv > 0) call post_data(cs%id_kv, kv, cs%diag)
179  if (cs%id_N2 > 0) call post_data(cs%id_N2, cs%N2, cs%diag)
180  if (cs%id_S2 > 0) call post_data(cs%id_S2, cs%S2, cs%diag)
181  if (cs%id_ri_grad > 0) call post_data(cs%id_ri_grad, cs%ri_grad, cs%diag)
182  if (cs%id_ri_grad_smooth > 0) call post_data(cs%id_ri_grad_smooth ,cs%ri_grad_smooth, cs%diag)
183 

◆ cvmix_shear_end()

subroutine, public mom_cvmix_shear::cvmix_shear_end ( type(cvmix_shear_cs), pointer  CS)

Clear pointers and dealocate memory.

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

Definition at line 314 of file MOM_CVMix_shear.F90.

314  type(CVMix_shear_cs), pointer :: CS !< Control structure for this module that
315  !! will be deallocated in this subroutine
316 
317  if (.not. associated(cs)) return
318 
319  if (cs%id_N2 > 0) deallocate(cs%N2)
320  if (cs%id_S2 > 0) deallocate(cs%S2)
321  if (cs%id_ri_grad > 0) deallocate(cs%ri_grad)
322  deallocate(cs)
323 

◆ cvmix_shear_init()

logical function, public mom_cvmix_shear::cvmix_shear_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_shear_cs), pointer  CS 
)

Initialized the CVMix internal shear mixing routine.

Note
*This is where we test to make sure multiple internal shear mixing routines (including JHL) are not enabled at the same time. (returns) CVMix_shear_init - True if module is to be used, False otherwise
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 192 of file MOM_CVMix_shear.F90.

192  type(time_type), intent(in) :: Time !< The current time.
193  type(ocean_grid_type), intent(in) :: G !< Grid structure.
194  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
195  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
196  type(param_file_type), intent(in) :: param_file !< Run-time parameter file handle
197  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
198  type(CVMix_shear_cs), pointer :: CS !< This module's control structure.
199  ! Local variables
200  integer :: NumberTrue=0
201  logical :: use_JHL
202 ! This include declares and sets the variable "version".
203 #include "version_variable.h"
204 
205  if (associated(cs)) then
206  call mom_error(warning, "CVMix_shear_init called with an associated "// &
207  "control structure.")
208  return
209  endif
210  allocate(cs)
211 
212 ! Set default, read and log parameters
213  call log_version(param_file, mdl, version, &
214  "Parameterization of shear-driven turbulence via CVMix (various options)")
215  call get_param(param_file, mdl, "USE_LMD94", cs%use_LMD94, &
216  "If true, use the Large-McWilliams-Doney (JGR 1994) "//&
217  "shear mixing parameterization.", default=.false.)
218  if (cs%use_LMD94) then
219  numbertrue=numbertrue + 1
220  cs%Mix_Scheme='KPP'
221  endif
222  call get_param(param_file, mdl, "USE_PP81", cs%use_PP81, &
223  "If true, use the Pacanowski and Philander (JPO 1981) "//&
224  "shear mixing parameterization.", default=.false.)
225  if (cs%use_PP81) then
226  numbertrue = numbertrue + 1
227  cs%Mix_Scheme='PP'
228  endif
229  use_jhl=kappa_shear_is_used(param_file)
230  if (use_jhl) numbertrue = numbertrue + 1
231  ! After testing for interior schemes, make sure only 0 or 1 are enabled.
232  ! Otherwise, warn user and kill job.
233  if ((numbertrue) > 1) then
234  call mom_error(fatal, 'MOM_CVMix_shear_init: '// &
235  'Multiple shear driven internal mixing schemes selected,'//&
236  ' please disable all but one scheme to proceed.')
237  endif
238  cvmix_shear_init=(cs%use_PP81.or.cs%use_LMD94)
239 
240 ! Forego remainder of initialization if not using this scheme
241  if (.not. cvmix_shear_init) return
242  call get_param(param_file, mdl, "NU_ZERO", cs%Nu_Zero, &
243  "Leading coefficient in KPP shear mixing.", &
244  units="nondim", default=5.e-3)
245  call get_param(param_file, mdl, "RI_ZERO", cs%Ri_Zero, &
246  "Critical Richardson for KPP shear mixing, "// &
247  "NOTE this the internal mixing and this is "// &
248  "not for setting the boundary layer depth." &
249  ,units="nondim", default=0.8)
250  call get_param(param_file, mdl, "KPP_EXP", cs%KPP_exp, &
251  "Exponent of unitless factor of diffusivities, "// &
252  "for KPP internal shear mixing scheme." &
253  ,units="nondim", default=3.0)
254  call get_param(param_file, mdl, "SMOOTH_RI", cs%smooth_ri, &
255  "If true, vertically smooth the Richardson "// &
256  "number by applying a 1-2-1 filter once.", &
257  default = .false.)
258  call cvmix_init_shear(mix_scheme=cs%Mix_Scheme, &
259  kpp_nu_zero=cs%Nu_Zero, &
260  kpp_ri_zero=cs%Ri_zero, &
261  kpp_exp=cs%KPP_exp)
262 
263  ! Register diagnostics; allocation and initialization
264  cs%diag => diag
265 
266  cs%id_N2 = register_diag_field('ocean_model', 'N2_shear', diag%axesTi, time, &
267  'Square of Brunt-Vaisala frequency used by MOM_CVMix_shear module', '1/s2')
268  if (cs%id_N2 > 0) then
269  allocate( cs%N2( szi_(g), szj_(g), szk_(g)+1 ) ) ; cs%N2(:,:,:) = 0.
270  endif
271 
272  cs%id_S2 = register_diag_field('ocean_model', 'S2_shear', diag%axesTi, time, &
273  'Square of vertical shear used by MOM_CVMix_shear module','1/s2')
274  if (cs%id_S2 > 0) then
275  allocate( cs%S2( szi_(g), szj_(g), szk_(g)+1 ) ) ; cs%S2(:,:,:) = 0.
276  endif
277 
278  cs%id_ri_grad = register_diag_field('ocean_model', 'ri_grad_shear', diag%axesTi, time, &
279  'Gradient Richarson number used by MOM_CVMix_shear module','nondim')
280  if (cs%id_ri_grad > 0) then !Initialize w/ large Richardson value
281  allocate( cs%ri_grad( szi_(g), szj_(g), szk_(g)+1 )) ; cs%ri_grad(:,:,:) = 1.e8
282  endif
283 
284  cs%id_ri_grad_smooth = register_diag_field('ocean_model', 'ri_grad_shear_smooth', &
285  diag%axesTi, time, &
286  'Smoothed gradient Richarson number used by MOM_CVMix_shear module','nondim')
287  if (cs%id_ri_grad_smooth > 0) then !Initialize w/ large Richardson value
288  allocate( cs%ri_grad_smooth( szi_(g), szj_(g), szk_(g)+1 )) ; cs%ri_grad_smooth(:,:,:) = 1.e8
289  endif
290 
291  cs%id_kd = register_diag_field('ocean_model', 'kd_shear_CVMix', diag%axesTi, time, &
292  'Vertical diffusivity added by MOM_CVMix_shear module', 'm2/s', conversion=us%Z2_T_to_m2_s)
293  cs%id_kv = register_diag_field('ocean_model', 'kv_shear_CVMix', diag%axesTi, time, &
294  'Vertical viscosity added by MOM_CVMix_shear module', 'm2/s', conversion=us%Z2_T_to_m2_s)
295 

◆ cvmix_shear_is_used()

logical function, public mom_cvmix_shear::cvmix_shear_is_used ( type(param_file_type), intent(in)  param_file)

Reads the parameters "LMD94" and "PP81" 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_fileRun-time parameter files handle.

Definition at line 302 of file MOM_CVMix_shear.F90.

302  type(param_file_type), intent(in) :: param_file !< Run-time parameter files handle.
303  ! Local variables
304  logical :: LMD94, PP81
305  call get_param(param_file, mdl, "USE_LMD94", lmd94, &
306  default=.false., do_not_log = .true.)
307  call get_param(param_file, mdl, "Use_PP81", pp81, &
308  default=.false., do_not_log = .true.)
309  cvmix_shear_is_used = (lmd94 .or. pp81)