MOM6
scm_cvmix_tests Module Reference

Detailed Description

Initial conditions and forcing for the single column model (SCM) CVMix test set.

Data Types

type  scm_cvmix_tests_cs
 Container for surface forcing parameters. More...
 

Functions/Subroutines

subroutine, public scm_cvmix_tests_ts_init (T, S, h, G, GV, US, param_file, just_read_params)
 Initializes temperature and salinity for the SCM CVMix test example. More...
 
subroutine, public scm_cvmix_tests_surface_forcing_init (Time, G, param_file, CS)
 Initializes surface forcing for the CVMix test case suite. More...
 
subroutine, public scm_cvmix_tests_wind_forcing (sfc_state, forces, day, G, US, CS)
 
subroutine, public scm_cvmix_tests_buoyancy_forcing (sfc_state, fluxes, day, G, US, CS)
 

Variables

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

Function/Subroutine Documentation

◆ scm_cvmix_tests_buoyancy_forcing()

subroutine, public scm_cvmix_tests::scm_cvmix_tests_buoyancy_forcing ( type(surface), intent(in)  sfc_state,
type(forcing), intent(inout)  fluxes,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(scm_cvmix_tests_cs), pointer  CS 
)
Parameters
[in]sfc_stateSurface state structure
[in,out]fluxesSurface fluxes structure
[in]dayCurrent model time
[in,out]gGrid structure
[in]usA dimensional unit scaling type
csContainer for SCM parameters

Definition at line 236 of file SCM_CVMix_tests.F90.

237  type(surface), intent(in) :: sfc_state !< Surface state structure
238  type(forcing), intent(inout) :: fluxes !< Surface fluxes structure
239  type(time_type), intent(in) :: day !< Current model time
240  type(ocean_grid_type), intent(inout) :: G !< Grid structure
241  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
242  type(SCM_CVMix_tests_CS), pointer :: CS !< Container for SCM parameters
243 
244  ! Local variables
245  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
246  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
247  real :: PI
248 
249  pi = 4.0*atan(1.0)
250 
251  ! Bounds for loops and memory allocation
252  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
253  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
254  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
255  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
256 
257  if (cs%UseHeatFlux) then
258  ! Note CVMix test inputs give Heat flux in [m K/s]
259  ! therefore must convert to W/m2 by multiplying
260  ! by Rho0*Cp
261  do j=jsq,jeq ; do i=is,ie
262  fluxes%sens(i,j) = cs%surf_HF * cs%Rho0 * fluxes%C_p
263  enddo ; enddo
264  endif
265 
266  if (cs%UseEvaporation) then
267  do j=jsq,jeq ; do i=is,ie
268  ! Note CVMix test inputs give evaporation in [m s-1]
269  ! This therefore must be converted to mass flux in [R Z T-1 ~> kg m-2 s-1]
270  ! by multiplying by density and some unit conversion factors.
271  fluxes%evap(i,j) = cs%surf_evap * cs%Rho0
272  enddo ; enddo
273  endif
274 
275  if (cs%UseDiurnalSW) then
276  do j=jsq,jeq ; do i=is,ie
277  ! Note CVMix test inputs give max sw rad in [m degC/s]
278  ! therefore must convert to W/m2 by multiplying by Rho0*Cp
279  ! Note diurnal cycle peaks at Noon.
280  fluxes%sw(i,j) = cs%Max_sw * max(0.0, cos(2*pi*(time_type_to_real(day)/86400.0 - 0.5))) * cs%RHO0 * fluxes%C_p
281  enddo ; enddo
282  endif
283 

◆ scm_cvmix_tests_surface_forcing_init()

subroutine, public scm_cvmix_tests::scm_cvmix_tests_surface_forcing_init ( type(time_type), intent(in)  Time,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(scm_cvmix_tests_cs), pointer  CS 
)

Initializes surface forcing for the CVMix test case suite.

Parameters
[in]timeModel time
[in]gGrid structure
[in]param_fileInput parameter structure
csParameter container

Definition at line 128 of file SCM_CVMix_tests.F90.

129  type(time_type), intent(in) :: Time !< Model time
130  type(ocean_grid_type), intent(in) :: G !< Grid structure
131  type(param_file_type), intent(in) :: param_file !< Input parameter structure
132  type(SCM_CVMix_tests_CS), pointer :: CS !< Parameter container
133 
134 
135  ! This include declares and sets the variable "version".
136 # include "version_variable.h"
137  type(unit_scale_type), pointer :: US => null() !< A dimensional unit scaling type
138 
139  us => g%US
140 
141  if (associated(cs)) then
142  call mom_error(fatal, "SCM_CVMix_tests_surface_forcing_init called with an associated "// &
143  "control structure.")
144  return
145  endif
146  allocate(cs)
147 
148  ! Read all relevant parameters and write them to the model log.
149  call log_version(param_file, mdl, version, "")
150  call get_param(param_file, mdl, "SCM_USE_WIND_STRESS", &
151  cs%UseWindStress, "Wind Stress switch "// &
152  "used in the SCM CVMix surface forcing.", &
153  units='', default=.false.)
154  call get_param(param_file, mdl, "SCM_USE_HEAT_FLUX", &
155  cs%UseHeatFlux, "Heat flux switch "// &
156  "used in the SCM CVMix test surface forcing.", &
157  units='', default=.false.)
158  call get_param(param_file, mdl, "SCM_USE_EVAPORATION", &
159  cs%UseEvaporation, "Evaporation switch "// &
160  "used in the SCM CVMix test surface forcing.", &
161  units='', default=.false.)
162  call get_param(param_file, mdl, "SCM_USE_DIURNAL_SW", &
163  cs%UseDiurnalSW, "Diurnal sw radation switch "// &
164  "used in the SCM CVMix test surface forcing.", &
165  units='', default=.false.)
166  if (cs%UseWindStress) then
167  call get_param(param_file, mdl, "SCM_TAU_X", &
168  cs%tau_x, "Constant X-dir wind stress "// &
169  "used in the SCM CVMix test surface forcing.", &
170  units='N/m2', scale=us%kg_m2s_to_RZ_T*us%m_s_to_L_T, fail_if_missing=.true.)
171  call get_param(param_file, mdl, "SCM_TAU_Y", &
172  cs%tau_y, "Constant y-dir wind stress "// &
173  "used in the SCM CVMix test surface forcing.", &
174  units='N/m2', scale=us%kg_m2s_to_RZ_T*us%m_s_to_L_T, fail_if_missing=.true.)
175  endif
176  if (cs%UseHeatFlux) then
177  call get_param(param_file, mdl, "SCM_HEAT_FLUX", &
178  cs%surf_HF, "Constant surface heat flux "// &
179  "used in the SCM CVMix test surface forcing.", &
180  units='m K/s', scale=us%m_to_Z*us%T_to_s, fail_if_missing=.true.)
181  endif
182  if (cs%UseEvaporation) then
183  call get_param(param_file, mdl, "SCM_EVAPORATION", &
184  cs%surf_evap, "Constant surface evaporation "// &
185  "used in the SCM CVMix test surface forcing.", &
186  units='m/s', scale=us%m_to_Z*us%T_to_s, fail_if_missing=.true.)
187  endif
188  if (cs%UseDiurnalSW) then
189  call get_param(param_file, mdl, "SCM_DIURNAL_SW_MAX", &
190  cs%Max_sw, "Maximum diurnal sw radiation "// &
191  "used in the SCM CVMix test surface forcing.", &
192  units='m K/s', scale=us%m_to_Z*us%T_to_s, fail_if_missing=.true.)
193  endif
194  call get_param(param_file, mdl, "RHO_0", cs%Rho0, &
195  "The mean ocean density used with BOUSSINESQ true to "//&
196  "calculate accelerations and the mass for conservation "//&
197  "properties, or with BOUSSINSEQ false to convert some "//&
198  "parameters from vertical units of m to kg m-2.", &
199  units="kg m-3", default=1035.0, scale=us%kg_m3_to_R)
200 

◆ scm_cvmix_tests_ts_init()

subroutine, public scm_cvmix_tests::scm_cvmix_tests_ts_init ( real, dimension(szi_(g),szj_(g),szk_(gv)), intent(out)  T,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(out)  S,
real, dimension(szi_(g),szj_(g),szk_(gv)), intent(in)  h,
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,
logical, intent(in), optional  just_read_params 
)

Initializes temperature and salinity for the SCM CVMix test example.

Parameters
[in]gGrid structure
[in]gvVertical grid structure
[out]tPotential temperature [degC]
[out]sSalinity [psu]
[in]hLayer thickness [H ~> m or kg m-2]
[in]usA dimensional unit scaling type
[in]param_fileInput parameter structure
[in]just_read_paramsIf present and true, this call will only read parameters without changing h.

Definition at line 55 of file SCM_CVMix_tests.F90.

56  type(ocean_grid_type), intent(in) :: G !< Grid structure
57  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure
58  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: T !< Potential temperature [degC]
59  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: S !< Salinity [psu]
60  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
61  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
62  type(param_file_type), intent(in) :: param_file !< Input parameter structure
63  logical, optional, intent(in) :: just_read_params !< If present and true, this call
64  !! will only read parameters without changing h.
65  ! Local variables
66  real :: UpperLayerTempMLD !< Upper layer Temp MLD thickness [Z ~> m].
67  real :: UpperLayerSaltMLD !< Upper layer Salt MLD thickness [Z ~> m].
68  real :: UpperLayerTemp !< Upper layer temperature (SST if thickness 0) [degC]
69  real :: UpperLayerSalt !< Upper layer salinity (SSS if thickness 0) [ppt]
70  real :: LowerLayerTemp !< Temp at top of lower layer [degC]
71  real :: LowerLayerSalt !< Salt at top of lower layer [ppt]
72  real :: LowerLayerdTdz !< Temp gradient in lower layer [degC / Z ~> degC m-1].
73  real :: LowerLayerdSdz !< Salt gradient in lower layer [ppt / Z ~> ppt m-1].
74  real :: LowerLayerMinTemp !< Minimum temperature in lower layer [degC]
75  real :: zC, DZ, top, bottom ! Depths and thicknesses [Z ~> m].
76  logical :: just_read ! If true, just read parameters but set nothing.
77  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
78 
79  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
80  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
81 
82 
83  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
84 
85  if (.not.just_read) call log_version(param_file, mdl, version)
86  call get_param(param_file, mdl, "SCM_TEMP_MLD", upperlayertempmld, &
87  'Initial temp mixed layer depth', &
88  units='m', default=0.0, scale=us%m_to_Z, do_not_log=just_read)
89  call get_param(param_file, mdl, "SCM_SALT_MLD", upperlayersaltmld, &
90  'Initial salt mixed layer depth', &
91  units='m', default=0.0, scale=us%m_to_Z, do_not_log=just_read)
92  call get_param(param_file, mdl, "SCM_L1_SALT", upperlayersalt, &
93  'Layer 2 surface salinity', units='1e-3', default=35.0, do_not_log=just_read)
94  call get_param(param_file, mdl, "SCM_L1_TEMP", upperlayertemp, &
95  'Layer 1 surface temperature', units='C', default=20.0, do_not_log=just_read)
96  call get_param(param_file, mdl, "SCM_L2_SALT", lowerlayersalt, &
97  'Layer 2 surface salinity', units='1e-3', default=35.0, do_not_log=just_read)
98  call get_param(param_file, mdl, "SCM_L2_TEMP", lowerlayertemp, &
99  'Layer 2 surface temperature', units='C', default=20.0, do_not_log=just_read)
100  call get_param(param_file, mdl, "SCM_L2_DTDZ", lowerlayerdtdz, &
101  'Initial temperature stratification in layer 2', &
102  units='C/m', default=0.0, scale=us%Z_to_m, do_not_log=just_read)
103  call get_param(param_file, mdl, "SCM_L2_DSDZ", lowerlayerdsdz, &
104  'Initial salinity stratification in layer 2', &
105  units='PPT/m', default=0.0, scale=us%Z_to_m, do_not_log=just_read)
106  call get_param(param_file, mdl, "SCM_L2_MINTEMP",lowerlayermintemp, &
107  'Layer 2 minimum temperature', units='C', default=4.0, do_not_log=just_read)
108 
109  if (just_read) return ! All run-time parameters have been read, so return.
110 
111  do j=js,je ; do i=is,ie
112  top = 0. ! Reference to surface
113  bottom = 0.
114  do k=1,nz
115  bottom = bottom - h(i,j,k)*gv%H_to_Z ! Interface below layer [Z ~> m]
116  zc = 0.5*( top + bottom ) ! Z of middle of layer [Z ~> m]
117  dz = min(0., zc + upperlayertempmld)
118  t(i,j,k) = max(lowerlayermintemp,lowerlayertemp + lowerlayerdtdz * dz)
119  dz = min(0., zc + upperlayersaltmld)
120  s(i,j,k) = lowerlayersalt + lowerlayerdsdz * dz
121  top = bottom
122  enddo ! k
123  enddo ; enddo
124 

◆ scm_cvmix_tests_wind_forcing()

subroutine, public scm_cvmix_tests::scm_cvmix_tests_wind_forcing ( type(surface), intent(in)  sfc_state,
type(mech_forcing), intent(inout)  forces,
type(time_type), intent(in)  day,
type(ocean_grid_type), intent(inout)  G,
type(unit_scale_type), intent(in)  US,
type(scm_cvmix_tests_cs), pointer  CS 
)
Parameters
[in]sfc_stateSurface state structure
[in,out]forcesA structure with the driving mechanical forces
[in]dayTime in days
[in,out]gGrid structure
[in]usA dimensional unit scaling type
csContainer for SCM parameters

Definition at line 203 of file SCM_CVMix_tests.F90.

204  type(surface), intent(in) :: sfc_state !< Surface state structure
205  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
206  type(time_type), intent(in) :: day !< Time in days
207  type(ocean_grid_type), intent(inout) :: G !< Grid structure
208  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
209  type(SCM_CVMix_tests_CS), pointer :: CS !< Container for SCM parameters
210  ! Local variables
211  integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq
212  integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
213  real :: mag_tau
214  ! Bounds for loops and memory allocation
215  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
216  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
217  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
218  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
219 
220  do j=js,je ; do i=isq,ieq
221  forces%taux(i,j) = cs%tau_x
222  enddo ; enddo
223  do j=jsq,jeq ; do i=is,ie
224  forces%tauy(i,j) = cs%tau_y
225  enddo ; enddo
226  call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
227 
228  mag_tau = sqrt(cs%tau_x*cs%tau_x + cs%tau_y*cs%tau_y)
229  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
230  forces%ustar(i,j) = sqrt( us%L_to_Z * mag_tau / (cs%Rho0) )
231  enddo ; enddo ; endif
232