MOM6
mom_cvmix_ddiff Module Reference

Detailed Description

Interface to CVMix double diffusion scheme.

Data Types

type  cvmix_ddiff_cs
 Control structure including parameters for CVMix double diffusion. More...
 

Functions/Subroutines

logical function, public cvmix_ddiff_init (Time, G, GV, US, param_file, diag, CS)
 Initialized the CVMix double diffusion module. More...
 
subroutine, public compute_ddiff_coeffs (h, tv, G, GV, US, j, Kd_T, Kd_S, CS, R_rho)
 Subroutine for computing vertical diffusion coefficients for the double diffusion mixing parameterization. More...
 
logical function, public cvmix_ddiff_is_used (param_file)
 Reads the parameter "USE_CVMIX_DDIFF" 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_ddiff_end (CS)
 Clear pointers and dealocate memory. More...
 

Variables

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

Function/Subroutine Documentation

◆ compute_ddiff_coeffs()

subroutine, public mom_cvmix_ddiff::compute_ddiff_coeffs ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %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,
integer, intent(in)  j,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke+1), intent(inout)  Kd_T,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke+1), intent(inout)  Kd_S,
type(cvmix_ddiff_cs), pointer  CS,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke+1), intent(inout), optional  R_rho 
)

Subroutine for computing vertical diffusion coefficients for the double diffusion mixing parameterization.

Parameters
[in]gGrid structure.
[in]gvVertical grid structure.
[in]hLayer thickness [H ~> m or kg m-2].
[in]tvThermodynamics structure.
[in]usA dimensional unit scaling type
[in]jMeridional grid index to work on.
[in,out]kd_tInterface double diffusion diapycnal diffusivity for temp [Z2 T-1 ~> m2 s-1].
[in,out]kd_sInterface double diffusion diapycnal diffusivity for salt [Z2 T-1 ~> m2 s-1].
csThe control structure returned by a previous call to CVMix_ddiff_init.
[in,out]r_rhoThe density ratios at interfaces [nondim].

Definition at line 143 of file MOM_CVMix_ddiff.F90.

143 
144  type(ocean_grid_type), intent(in) :: G !< Grid structure.
145  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
146  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2].
147  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure.
148  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
149  integer, intent(in) :: j !< Meridional grid index to work on.
150  ! Kd_T and Kd_S are intent inout because only one j-row is set here, but they are essentially outputs.
151  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(inout) :: Kd_T !< Interface double diffusion diapycnal
152  !! diffusivity for temp [Z2 T-1 ~> m2 s-1].
153  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(inout) :: Kd_S !< Interface double diffusion diapycnal
154  !! diffusivity for salt [Z2 T-1 ~> m2 s-1].
155  type(CVMix_ddiff_cs), pointer :: CS !< The control structure returned
156  !! by a previous call to CVMix_ddiff_init.
157  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), &
158  optional, intent(inout) :: R_rho !< The density ratios at interfaces [nondim].
159 
160  ! Local variables
161  real, dimension(SZK_(GV)) :: &
162  cellHeight, & !< Height of cell centers [m]
163  dRho_dT, & !< partial derivatives of density wrt temp [R degC-1 ~> kg m-3 degC-1]
164  dRho_dS, & !< partial derivatives of density wrt saln [R ppt-1 ~> kg m-3 ppt-1]
165  pres_int, & !< pressure at each interface [R L2 T-2 ~> Pa]
166  temp_int, & !< temp and at interfaces [degC]
167  salt_int, & !< salt at at interfaces [ppt]
168  alpha_dT, & !< alpha*dT across interfaces [kg m-3]
169  beta_dS, & !< beta*dS across interfaces [kg m-3]
170  dT, & !< temp. difference between adjacent layers [degC]
171  dS !< salt difference between adjacent layers [ppt]
172  real, dimension(SZK_(GV)+1) :: &
173  Kd1_T, & !< Diapycanal diffusivity of temperature [m2 s-1].
174  Kd1_S !< Diapycanal diffusivity of salinity [m2 s-1].
175 
176  real, dimension(SZK_(GV)+1) :: iFaceHeight !< Height of interfaces [m]
177  integer :: kOBL !< level of OBL extent
178  real :: dh, hcorr
179  integer :: i, k
180 
181  ! initialize dummy variables
182  pres_int(:) = 0.0; temp_int(:) = 0.0; salt_int(:) = 0.0
183  alpha_dt(:) = 0.0; beta_ds(:) = 0.0; drho_dt(:) = 0.0
184  drho_ds(:) = 0.0; dt(:) = 0.0; ds(:) = 0.0
185 
186 
187  ! GMM, I am leaving some code commented below. We need to pass BLD to
188  ! this soubroutine to avoid adding diffusivity above that. This needs
189  ! to be done once we re-structure the order of the calls.
190  !if (.not. associated(hbl)) then
191  ! allocate(hbl(SZI_(G), SZJ_(G)));
192  ! hbl(:,:) = 0.0
193  !endif
194 
195  do i = g%isc, g%iec
196 
197  ! skip calling at land points
198  if (g%mask2dT(i,j) == 0.) cycle
199 
200  pres_int(1) = 0. ; if (associated(tv%p_surf)) pres_int(1) = tv%p_surf(i,j)
201  ! we don't have SST and SSS, so let's use values at top-most layer
202  temp_int(1) = tv%T(i,j,1); salt_int(1) = tv%S(i,j,1)
203  do k=2,g%ke
204  ! pressure at interface
205  pres_int(k) = pres_int(k-1) + (gv%g_Earth * gv%H_to_RZ) * h(i,j,k-1)
206  ! temp and salt at interface
207  ! for temp: (t1*h1 + t2*h2)/(h1+h2)
208  temp_int(k) = (tv%T(i,j,k-1)*h(i,j,k-1) + tv%T(i,j,k)*h(i,j,k)) / (h(i,j,k-1)+h(i,j,k))
209  salt_int(k) = (tv%S(i,j,k-1)*h(i,j,k-1) + tv%S(i,j,k)*h(i,j,k)) / (h(i,j,k-1)+h(i,j,k))
210  ! dT and dS
211  dt(k) = (tv%T(i,j,k-1)-tv%T(i,j,k))
212  ds(k) = (tv%S(i,j,k-1)-tv%S(i,j,k))
213  enddo ! k-loop finishes
214 
215  call calculate_density_derivs(temp_int, salt_int, pres_int, drho_dt, drho_ds, tv%eqn_of_state)
216 
217  ! The "-1.0" below is needed so that the following criteria is satisfied:
218  ! if ((alpha_dT > beta_dS) .and. (beta_dS > 0.0)) then "salt finger"
219  ! if ((alpha_dT < 0.) .and. (beta_dS < 0.) .and. (alpha_dT > beta_dS)) then "diffusive convection"
220  do k=1,g%ke
221  alpha_dt(k) = -1.0*us%R_to_kg_m3*drho_dt(k) * dt(k)
222  beta_ds(k) = us%R_to_kg_m3*drho_ds(k) * ds(k)
223  enddo
224 
225  if (present(r_rho)) then
226  do k=1,g%ke
227  ! Set R_rho using Adcroft's rule of reciprocals.
228  r_rho(i,j,k) = 0.0 ; if (abs(beta_ds(k)) > 0.0) r_rho(i,j,k) = alpha_dt(k) / beta_ds(k)
229  ! avoid NaN's again for safety, perhaps unnecessarily.
230  if (r_rho(i,j,k) /= r_rho(i,j,k)) r_rho(i,j,k) = 0.0
231  enddo
232  endif
233 
234  ifaceheight(1) = 0.0 ! BBL is all relative to the surface
235  hcorr = 0.0
236  ! compute heights at cell center and interfaces
237  do k=1,g%ke
238  dh = h(i,j,k) * gv%H_to_m ! Nominal thickness to use for increment
239  dh = dh + hcorr ! Take away the accumulated error (could temporarily make dh<0)
240  hcorr = min( dh - cs%min_thickness, 0. ) ! If inflating then hcorr<0
241  dh = max( dh, cs%min_thickness ) ! Limit increment dh>=min_thickness
242  cellheight(k) = ifaceheight(k) - 0.5 * dh
243  ifaceheight(k+1) = ifaceheight(k) - dh
244  enddo
245 
246  ! gets index of the level and interface above hbl
247  !kOBL = CVmix_kpp_compute_kOBL_depth(iFaceHeight, cellHeight,hbl(i,j))
248 
249  kd1_t(:) = 0.0 ; kd1_s(:) = 0.0
250  call cvmix_coeffs_ddiff(tdiff_out=kd1_t(:), &
251  sdiff_out=kd1_s(:), &
252  strat_param_num=alpha_dt(:), &
253  strat_param_denom=beta_ds(:), &
254  nlev=g%ke, &
255  max_nlev=g%ke)
256  do k=1,g%ke+1
257  kd_t(i,j,k) = us%m2_s_to_Z2_T * kd1_t(k)
258  kd_s(i,j,k) = us%m2_s_to_Z2_T * kd1_s(k)
259  enddo
260 
261  ! Do not apply mixing due to convection within the boundary layer
262  !do k=1,kOBL
263  ! Kd_T(i,j,k) = 0.0
264  ! Kd_S(i,j,k) = 0.0
265  !enddo
266 
267  enddo ! i-loop
268 

◆ cvmix_ddiff_end()

subroutine, public mom_cvmix_ddiff::cvmix_ddiff_end ( type(cvmix_ddiff_cs), pointer  CS)

Clear pointers and dealocate memory.

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

Definition at line 283 of file MOM_CVMix_ddiff.F90.

283  type(CVMix_ddiff_cs), pointer :: CS !< Control structure for this module that
284  !! will be deallocated in this subroutine
285 
286  deallocate(cs)
287 

◆ cvmix_ddiff_init()

logical function, public mom_cvmix_ddiff::cvmix_ddiff_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_ddiff_cs), pointer  CS 
)

Initialized the CVMix double diffusion module.

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 51 of file MOM_CVMix_ddiff.F90.

51 
52  type(time_type), intent(in) :: Time !< The current time.
53  type(ocean_grid_type), intent(in) :: G !< Grid structure.
54  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
55  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
56  type(param_file_type), intent(in) :: param_file !< Run-time parameter file handle
57  type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
58  type(CVMix_ddiff_cs), pointer :: CS !< This module's control structure.
59 
60 ! This include declares and sets the variable "version".
61 #include "version_variable.h"
62 
63  if (associated(cs)) then
64  call mom_error(warning, "CVMix_ddiff_init called with an associated "// &
65  "control structure.")
66  return
67  endif
68  allocate(cs)
69 
70  ! Read parameters
71  call get_param(param_file, mdl, "USE_CVMIX_DDIFF", cvmix_ddiff_init, default=.false., do_not_log=.true.)
72  call log_version(param_file, mdl, version, &
73  "Parameterization of mixing due to double diffusion processes via CVMix", &
74  all_default=.not.cvmix_ddiff_init)
75  call get_param(param_file, mdl, "USE_CVMIX_DDIFF", cvmix_ddiff_init, &
76  "If true, turns on double diffusive processes via CVMix. "//&
77  "Note that double diffusive processes on viscosity are ignored "//&
78  "in CVMix, see http://cvmix.github.io/ for justification.", &
79  default=.false.)
80 
81  if (.not. cvmix_ddiff_init) return
82 
83  call get_param(param_file, mdl, 'DEBUG', cs%debug, default=.false., do_not_log=.true.)
84 
85  call get_param(param_file, mdl, 'MIN_THICKNESS', cs%min_thickness, default=0.001, do_not_log=.true.)
86 
87  call openparameterblock(param_file,'CVMIX_DDIFF')
88 
89  call get_param(param_file, mdl, "STRAT_PARAM_MAX", cs%strat_param_max, &
90  "The maximum value for the double dissusion stratification parameter", &
91  units="nondim", default=2.55)
92 
93  call get_param(param_file, mdl, "KAPPA_DDIFF_S", cs%kappa_ddiff_s, &
94  "Leading coefficient in formula for salt-fingering regime "//&
95  "for salinity diffusion.", units="m2 s-1", default=1.0e-4)
96 
97  call get_param(param_file, mdl, "DDIFF_EXP1", cs%ddiff_exp1, &
98  "Interior exponent in salt-fingering regime formula.", &
99  units="nondim", default=1.0)
100 
101  call get_param(param_file, mdl, "DDIFF_EXP2", cs%ddiff_exp2, &
102  "Exterior exponent in salt-fingering regime formula.", &
103  units="nondim", default=3.0)
104 
105  call get_param(param_file, mdl, "KAPPA_DDIFF_PARAM1", cs%kappa_ddiff_param1, &
106  "Exterior coefficient in diffusive convection regime.", &
107  units="nondim", default=0.909)
108 
109  call get_param(param_file, mdl, "KAPPA_DDIFF_PARAM2", cs%kappa_ddiff_param2, &
110  "Middle coefficient in diffusive convection regime.", &
111  units="nondim", default=4.6)
112 
113  call get_param(param_file, mdl, "KAPPA_DDIFF_PARAM3", cs%kappa_ddiff_param3, &
114  "Interior coefficient in diffusive convection regime.", &
115  units="nondim", default=-0.54)
116 
117  call get_param(param_file, mdl, "MOL_DIFF", cs%mol_diff, &
118  "Molecular diffusivity used in CVMix double diffusion.", &
119  units="m2 s-1", default=1.5e-6)
120 
121  call get_param(param_file, mdl, "DIFF_CONV_TYPE", cs%diff_conv_type, &
122  "type of diffusive convection to use. Options are Marmorino \n" //&
123  "and Caldwell 1976 (MC76) and Kelley 1988, 1990 (K90).", &
124  default="MC76")
125 
126  call closeparameterblock(param_file)
127 
128  call cvmix_init_ddiff(strat_param_max=cs%strat_param_max, &
129  kappa_ddiff_s=cs%kappa_ddiff_s, &
130  ddiff_exp1=cs%ddiff_exp1, &
131  ddiff_exp2=cs%ddiff_exp2, &
132  mol_diff=cs%mol_diff, &
133  kappa_ddiff_param1=cs%kappa_ddiff_param1, &
134  kappa_ddiff_param2=cs%kappa_ddiff_param2, &
135  kappa_ddiff_param3=cs%kappa_ddiff_param3, &
136  diff_conv_type=cs%diff_conv_type)
137 

◆ cvmix_ddiff_is_used()

logical function, public mom_cvmix_ddiff::cvmix_ddiff_is_used ( type(param_file_type), intent(in)  param_file)

Reads the parameter "USE_CVMIX_DDIFF" 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 275 of file MOM_CVMix_ddiff.F90.

275  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
276  call get_param(param_file, mdl, "USE_CVMIX_DDIFF", cvmix_ddiff_is_used, &
277  default=.false., do_not_log = .true.)
278