MOM6
SCM_CVMix_tests.F90
1 !> Initial conditions and forcing for the single column model (SCM) CVMix test set.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_domains, only : pass_var, pass_vector, to_all
7 use mom_error_handler, only : mom_error, fatal
10 use mom_grid, only : ocean_grid_type
14 use mom_time_manager, only : time_type, operator(+), operator(/), time_type_to_real
17 
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
22 public scm_cvmix_tests_ts_init
23 public scm_cvmix_tests_surface_forcing_init
24 public scm_cvmix_tests_wind_forcing
25 public scm_cvmix_tests_buoyancy_forcing
26 public scm_cvmix_tests_cs
27 
28 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
29 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
30 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
31 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
32 
33 !> Container for surface forcing parameters
34 type scm_cvmix_tests_cs ; private
35  logical :: usewindstress !< True to use wind stress
36  logical :: useheatflux !< True to use heat flux
37  logical :: useevaporation !< True to use evaporation
38  logical :: usediurnalsw !< True to use diurnal sw radiation
39  real :: tau_x !< (Constant) Wind stress, X [Pa]
40  real :: tau_y !< (Constant) Wind stress, Y [Pa]
41  real :: surf_hf !< (Constant) Heat flux [m degC s-1]
42  real :: surf_evap !< (Constant) Evaporation rate [m s-1]
43  real :: max_sw !< maximum of diurnal sw radiation [m degC s-1]
44  real,public :: rho0 !< reference density copied for easy passing [kg m-3]
45 end type
46 
47 ! This include declares and sets the variable "version".
48 #include "version_variable.h"
49 
50 character(len=40) :: mdl = "SCM_CVMix_tests" !< This module's name.
51 
52 contains
53 
54 !> Initializes temperature and salinity for the SCM CVMix test example
55 subroutine scm_cvmix_tests_ts_init(T, S, h, G, GV, US, param_file, just_read_params)
56  real, dimension(NIMEM_,NJMEM_, NKMEM_), intent(out) :: t !< Potential temperature [degC]
57  real, dimension(NIMEM_,NJMEM_, NKMEM_), intent(out) :: s !< Salinity [psu]
58  real, dimension(NIMEM_,NJMEM_, NKMEM_), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
59  type(ocean_grid_type), intent(in) :: g !< Grid structure
60  type(verticalgrid_type), intent(in) :: gv !< Vertical grid structure
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 will
64  !! 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 
125 end subroutine scm_cvmix_tests_ts_init
126 
127 !> Initializes surface forcing for the CVMix test case suite
128 subroutine scm_cvmix_tests_surface_forcing_init(Time, G, param_file, CS)
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 ! This include declares and sets the variable "version".
135 #include "version_variable.h"
136 
137  if (associated(cs)) then
138  call mom_error(fatal, "SCM_CVMix_tests_surface_forcing_init called with an associated "// &
139  "control structure.")
140  return
141  endif
142  allocate(cs)
143 
144  ! Read all relevant parameters and write them to the model log.
145  call log_version(param_file, mdl, version, "")
146  call get_param(param_file, mdl, "SCM_USE_WIND_STRESS", &
147  cs%UseWindStress, "Wind Stress switch "// &
148  "used in the SCM CVMix surface forcing.", &
149  units='', default=.false.)
150  call get_param(param_file, mdl, "SCM_USE_HEAT_FLUX", &
151  cs%UseHeatFlux, "Heat flux switch "// &
152  "used in the SCM CVMix test surface forcing.", &
153  units='', default=.false.)
154  call get_param(param_file, mdl, "SCM_USE_EVAPORATION", &
155  cs%UseEvaporation, "Evaporation switch "// &
156  "used in the SCM CVMix test surface forcing.", &
157  units='', default=.false.)
158  call get_param(param_file, mdl, "SCM_USE_DIURNAL_SW", &
159  cs%UseDiurnalSW, "Diurnal sw radation switch "// &
160  "used in the SCM CVMix test surface forcing.", &
161  units='', default=.false.)
162  if (cs%UseWindStress) then
163  call get_param(param_file, mdl, "SCM_TAU_X", &
164  cs%tau_x, "Constant X-dir wind stress "// &
165  "used in the SCM CVMix test surface forcing.", &
166  units='N/m2', fail_if_missing=.true.)
167  call get_param(param_file, mdl, "SCM_TAU_Y", &
168  cs%tau_y, "Constant y-dir wind stress "// &
169  "used in the SCM CVMix test surface forcing.", &
170  units='N/m2', fail_if_missing=.true.)
171  endif
172  if (cs%UseHeatFlux) then
173  call get_param(param_file, mdl, "SCM_HEAT_FLUX", &
174  cs%surf_HF, "Constant surface heat flux "// &
175  "used in the SCM CVMix test surface forcing.", &
176  units='m K/s', fail_if_missing=.true.)
177  endif
178  if (cs%UseEvaporation) then
179  call get_param(param_file, mdl, "SCM_EVAPORATION", &
180  cs%surf_evap, "Constant surface evaporation "// &
181  "used in the SCM CVMix test surface forcing.", &
182  units='m/s', fail_if_missing=.true.)
183  endif
184  if (cs%UseDiurnalSW) then
185  call get_param(param_file, mdl, "SCM_DIURNAL_SW_MAX", &
186  cs%Max_sw, "Maximum diurnal sw radiation "// &
187  "used in the SCM CVMix test surface forcing.", &
188  units='m K/s', fail_if_missing=.true.)
189  endif
190 
191 end subroutine scm_cvmix_tests_surface_forcing_init
192 
193 subroutine scm_cvmix_tests_wind_forcing(state, forces, day, G, US, CS)
194  type(surface), intent(in) :: state !< Surface state structure
195  type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
196  type(time_type), intent(in) :: day !< Time in days
197  type(ocean_grid_type), intent(inout) :: g !< Grid structure
198  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
199  type(scm_cvmix_tests_cs), pointer :: cs !< Container for SCM parameters
200  ! Local variables
201  integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq
202  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
203  real :: mag_tau
204  ! Bounds for loops and memory allocation
205  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
206  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
207  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
208  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
209 
210  do j=js,je ; do i=isq,ieq
211  forces%taux(i,j) = cs%tau_x
212  enddo ; enddo
213  do j=jsq,jeq ; do i=is,ie
214  forces%tauy(i,j) = cs%tau_y
215  enddo ; enddo
216  call pass_vector(forces%taux, forces%tauy, g%Domain, to_all)
217 
218 
219  mag_tau = sqrt(cs%tau_x*cs%tau_x + cs%tau_y*cs%tau_y)
220  if (associated(forces%ustar)) then ; do j=js,je ; do i=is,ie
221  forces%ustar(i,j) = us%m_to_Z*us%T_to_s * sqrt( mag_tau / cs%Rho0 )
222  enddo ; enddo ; endif
223 
224 end subroutine scm_cvmix_tests_wind_forcing
225 
226 
227 subroutine scm_cvmix_tests_buoyancy_forcing(state, fluxes, day, G, CS)
228  type(surface), intent(in) :: state !< Surface state structure
229  type(forcing), intent(inout) :: fluxes !< Surface fluxes structure
230  type(time_type), intent(in) :: day !< Current model time
231  type(ocean_grid_type), intent(inout) :: g !< Grid structure
232  type(scm_cvmix_tests_cs), pointer :: cs !< Container for SCM parameters
233 
234  ! Local variables
235  integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq
236  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
237  real :: pi
238 
239  pi = 4.0*atan(1.0)
240 
241  ! Bounds for loops and memory allocation
242  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
243  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
244  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
245  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
246 
247  if (cs%UseHeatFlux) then
248  ! Note CVMix test inputs give Heat flux in [m K/s]
249  ! therefore must convert to W/m2 by multiplying
250  ! by Rho0*Cp
251  do j=jsq,jeq ; do i=is,ie
252  fluxes%sens(i,j) = cs%surf_HF * cs%Rho0 * fluxes%C_p
253  enddo ; enddo
254  endif
255 
256  if (cs%UseEvaporation) then
257  do j=jsq,jeq ; do i=is,ie
258  ! Note CVMix test inputs give evaporation in [m s-1]
259  ! This therefore must be converted to mass flux
260  ! by multiplying by density
261  fluxes%evap(i,j) = cs%surf_evap * cs%Rho0
262  enddo ; enddo
263  endif
264 
265  if (cs%UseDiurnalSW) then
266  do j=jsq,jeq ; do i=is,ie
267  ! Note CVMix test inputs give max sw rad in [m K/s]
268  ! therefore must convert to W/m2 by multiplying
269  ! by Rho0*Cp
270  ! Note diurnal cycle peaks at Noon.
271  fluxes%sw(i,j) = cs%Max_sw * max(0.0,cos(2*pi* &
272  (time_type_to_real(day)/86400.-0.5))) * cs%RHO0 * fluxes%C_p
273  enddo ; enddo
274  endif
275 
276 end subroutine scm_cvmix_tests_buoyancy_forcing
277 
278 end module scm_cvmix_tests
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_forcing_type::mech_forcing
Structure that contains pointers to the mechanical forcing at the surface used to drive the liquid oc...
Definition: MOM_forcing_type.F90:185
mom_variables::surface
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Definition: MOM_variables.F90:38
mom_safe_alloc
Convenience functions for safely allocating memory without accidentally reallocating pointer and caus...
Definition: MOM_safe_alloc.F90:3
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:82
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
scm_cvmix_tests
Initial conditions and forcing for the single column model (SCM) CVMix test set.
Definition: SCM_CVMix_tests.F90:2
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
mom_forcing_type
This module implements boundary forcing for MOM6.
Definition: MOM_forcing_type.F90:2
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_forcing_type::forcing
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Definition: MOM_forcing_type.F90:49
scm_cvmix_tests::scm_cvmix_tests_cs
Container for surface forcing parameters.
Definition: SCM_CVMix_tests.F90:34
mom_safe_alloc::safe_alloc_ptr
Allocate a pointer to a 1-d, 2-d or 3-d array.
Definition: MOM_safe_alloc.F90:12
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25