MOM6
mom_diffconvection Module Reference

Data Types

type  diffconvection_cs
 

Functions/Subroutines

logical function, public diffconvection_init (paramFile, G, diag, Time, CS)
 
subroutine, public diffconvection_calculate (CS, G, GV, h, Temp, Salt, EOS, Kd_int)
 
subroutine, public diffconvection_end (CS)
 

Variables

logical, parameter verbose = .False.
 

Function/Subroutine Documentation

◆ diffconvection_calculate()

subroutine, public mom_diffconvection::diffconvection_calculate ( type(diffconvection_cs), pointer  CS,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  h,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  Temp,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(in)  Salt,
type(eos_type), pointer  EOS,
real, dimension(szi_(g),szj_(g),szk_(g)+1), intent(inout)  Kd_int 
)
Parameters
csControl structure
[in]gOcean grid
[in]gvOcean vertical grid
[in]hLayer/level thicknesses (units of H)
[in]tempPot. temperature (degrees C)
[in]saltSalinity (ppt)
eosEquation of state
[in,out]kd_int(in) Vertical diffusivity on interfaces (m2/s) (out) Modified vertical diffusivity (m2/s)

Definition at line 99 of file MOM_diffConvection.F90.

Referenced by mom_diabatic_driver::diabatic().

99 !< Calculates diffusivity and non-local transport for KPP parameterization
100 
101 ! Arguments
102  type(diffconvection_cs), pointer :: cs !< Control structure
103  type(ocean_grid_type), intent(in) :: g !< Ocean grid
104  type(verticalgrid_type), intent(in) :: gv !< Ocean vertical grid
105  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer/level thicknesses (units of H)
106  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: temp !< Pot. temperature (degrees C)
107  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: salt !< Salinity (ppt)
108  type(eos_type), pointer :: eos !< Equation of state
109  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(inout) :: kd_int !< (in) Vertical diffusivity on interfaces (m2/s)
110  !! (out) Modified vertical diffusivity (m2/s)
111 ! Local variables
112  integer :: i, j, k
113  real, dimension( G%ke+1 ) :: n2_1d ! Brunt-Vaisala frequency squared, at interfaces (1/s2)
114  real, dimension( G%ke+1 ) :: kd_1d ! Vertical diffusivity at interfaces (m2/s)
115  real :: gorho, pref, rhok, rhokm1
116 
117  gorho = gv%g_Earth / gv%Rho0
118 
119  n2_1d( 1 ) = 0.
120  n2_1d( g%ke+1 ) = 0.
121  kd_1d( 1 ) = 0.
122  kd_1d( g%ke+1 ) = 0.
123  do j=g%jsc,g%jec ; do i=g%isc,g%iec
124  ! This k-loop calculates external quantities independent of any iterations
125  ! Start at bottom of top level
126  pref = 0. ! Ignore atmospheric pressure
127  do k = 2, g%ke
128  ! Pressure at interface K is incremented by mass of level above
129  pref = pref + gv%g_Earth * gv%Rho0 * h(i,j,k-1) * gv%H_to_m ! Boussinesq approximation!!!! ?????
130  ! Compute Brunt-Vaisala frequency (static stability) on interfaces
131  call calculate_density(temp(i,j,k), salt(i,j,k), pref, rhok, eos)
132  call calculate_density(temp(i,j,k-1), salt(i,j,k-1), pref, rhokm1, eos)
133  n2_1d(k) = gorho * (rhok - rhokm1) / &
134  (0.5*(h(i,j,k-1) + h(i,j,k)) + gv%H_subroundoff) ! Can be negative
135  kd_1d(k) = 0.
136  if (n2_1d(k) < 0.) kd_1d(k) = cs%Kd_convection
137  enddo ! k
138 
139  if (.not. cs%passiveMode) kd_int(i,j,:) = kd_int(i,j,:) + kd_1d(:)
140 
141  if (cs%id_N2 > 0) cs%N2(i,j,:) = n2_1d(:)
142  if (cs%id_Kd_conv > 0) cs%Kd_conv(i,j,:) = kd_1d(:)
143 
144  enddo ; enddo ! j
145 
146  if (cs%id_N2 > 0) call post_data(cs%id_N2, cs%N2, cs%diag)
147  if (cs%id_Kd_conv > 0) call post_data(cs%id_Kd_conv, cs%Kd_conv, cs%diag)
148 
Here is the caller graph for this function:

◆ diffconvection_end()

subroutine, public mom_diffconvection::diffconvection_end ( type(diffconvection_cs), pointer  CS)

Definition at line 153 of file MOM_diffConvection.F90.

Referenced by mom_diabatic_driver::diabatic_driver_end().

153 ! Clear pointers, dealocate memory
154  type(diffconvection_cs), pointer :: cs ! Control structure
155 
156  if (cs%id_N2 > 0) deallocate(cs%N2, cs%diag)
157  if (cs%id_Kd_conv > 0) deallocate(cs%Kd_conv, cs%diag)
158  deallocate(cs)
Here is the caller graph for this function:

◆ diffconvection_init()

logical function, public mom_diffconvection::diffconvection_init ( type(param_file_type), intent(in)  paramFile,
type(ocean_grid_type), intent(in)  G,
type(diag_ctrl), intent(in), target  diag,
type(time_type), intent(in)  Time,
type(diffconvection_cs), pointer  CS 
)
Parameters
[in]paramfileFile parser
[in]gOcean grid
[in]diagDiagnostics
csControl structure

Definition at line 42 of file MOM_diffConvection.F90.

References mom_file_parser::closeparameterblock(), mom_error_handler::mom_error(), mom_file_parser::openparameterblock(), and mom_diag_mediator::register_diag_field().

42 !< Initialize the CVmix KPP module and set up diagnostics
43 !! Returns True if module is to be used, otherwise returns False.
44 
45 ! Arguments
46  type(param_file_type), intent(in) :: paramfile !< File parser
47  type(ocean_grid_type), intent(in) :: g !< Ocean grid
48  type(diag_ctrl), target, intent(in) :: diag !< Diagnostics
49  type(time_type), intent(in) :: time !< Time
50  type(diffconvection_cs), pointer :: cs !< Control structure
51 ! Local variables
52 #include "version_variable.h"
53  character(len=40) :: mdl = 'MOM_diffConvection' ! This module's name.
54 
55  if (associated(cs)) call mom_error(fatal, 'MOM_diffConvection, diffConvection_init: '// &
56  'Control structure has already been initialized')
57  allocate(cs)
58 
59 ! Read parameters
60  call log_version(paramfile, mdl, version, &
61  'This module implements enhanced diffusivity as a\n' // &
62  'function of static stability, N^2.')
63  call get_param(paramfile, mdl, "USE_CONVECTION", diffconvection_init, &
64  "If true, turns on the diffusive convection scheme that\n"// &
65  "increases diapycnal diffusivities at statically unstable\n"// &
66  "interfaces. Relevant parameters are contained in the\n"// &
67  "CONVECTION% parameter block.", &
68  default=.false.)
69 
70  call openparameterblock(paramfile,'CONVECTION')
71  call get_param(paramfile, mdl, 'PASSIVE', cs%passiveMode, &
72  'If True, puts KPP into a passive-diagnostic mode.', &
73  default=.false.)
74  call get_param(paramfile, mdl, 'KD_CONV', cs%Kd_convection, &
75  'DIffusivity used in statically unstable regions of column.', &
76  units='m2/s', default=1.00)
77  call closeparameterblock(paramfile)
78  call get_param(paramfile, mdl, 'DEBUG', cs%debug, default=.false., do_not_log=.true.)
79 
80 ! Forego remainder of initialization if not using this scheme
81  if (.not. diffconvection_init) return
82 
83 ! Register diagnostics
84  cs%diag => diag
85  cs%id_N2 = register_diag_field('ocean_model', 'Conv_N2', diag%axesTi, time, &
86  'Square of Brunt-Vaisala frequency used by diffConvection module', '1/s2')
87  if (cs%id_N2 > 0) allocate( cs%N2( szi_(g), szj_(g), szk_(g)+1 ) )
88  cs%id_Kd_conv = register_diag_field('ocean_model', 'Conv_Kd', diag%axesTi, time, &
89  'Additional diffusivity added by diffConvection module', 'm2/s')
90  if (cs%id_Kd_conv > 0) allocate( cs%Kd_conv( szi_(g), szj_(g), szk_(g)+1 ) )
91 
92  if (cs%id_N2 > 0) cs%N2(:,:,:) = 0.
93  if (cs%id_Kd_conv > 0) cs%Kd_conv(:,:,:) = 0.
94 
Here is the call graph for this function:

Variable Documentation

◆ verbose

logical, parameter mom_diffconvection::verbose = .False.

Definition at line 37 of file MOM_diffConvection.F90.

37 logical, parameter :: verbose = .false.