MOM6
baroclinic_zone_initialization.F90
Go to the documentation of this file.
1 !> Initial conditions for an idealized baroclinic zone
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
8 use mom_grid, only : ocean_grid_type
9 
10 implicit none ; private
11 
12 #include <MOM_memory.h>
13 #include "version_variable.h"
14 
15 ! Private (module-wise) parameters
16 character(len=40) :: mdl = "baroclinic_zone_initialization" !< This module's name.
17 
19 
20 contains
21 
22 !> Reads the parameters unique to this module
23 subroutine bcz_params(G, param_file, S_ref, dSdz, delta_S, dSdx, T_ref, dTdz, &
24  delta_T, dTdx, L_zone, just_read_params)
25  type(ocean_grid_type), intent(in) :: G !< Grid structure
26  type(param_file_type), intent(in) :: param_file !< Parameter file handle
27  real, intent(out) :: S_ref !< Reference salinity (ppt)
28  real, intent(out) :: dSdz !< Salinity stratification (ppt/m)
29  real, intent(out) :: delta_S !< Salinity difference across baroclinic zone (ppt)
30  real, intent(out) :: dSdx !< Linear salinity gradient (ppt/m)
31  real, intent(out) :: T_ref !< Reference temperature (ppt)
32  real, intent(out) :: dTdz !< Temperature stratification (ppt/m)
33  real, intent(out) :: delta_T !< Temperature difference across baroclinic zone (ppt)
34  real, intent(out) :: dTdx !< Linear temperature gradient (ppt/m)
35  real, intent(out) :: L_zone !< Width of baroclinic zone (m)
36  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
37  !! only read parameters without changing h.
38  logical :: just_read ! If true, just read parameters but set nothing.
39 
40  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
41 
42  if (.not.just_read) &
43  call log_version(param_file, mdl, version, 'Initialization of an analytic baroclinic zone')
44  call openparameterblock(param_file,'BCZIC')
45  call get_param(param_file, mdl, "S_REF", s_ref, 'Reference salinity', units='ppt', &
46  default=35., do_not_log=just_read)
47  call get_param(param_file, mdl,"DSDZ",dsdz,'Salinity stratification',units='ppt/m', &
48  default=0.0, do_not_log=just_read)
49  call get_param(param_file, mdl,"DELTA_S",delta_s,'Salinity difference across baroclinic zone', &
50  units='ppt', default=0.0, do_not_log=just_read)
51  call get_param(param_file, mdl,"DSDX",dsdx,'Meridional salinity difference', &
52  units='ppt/'//trim(g%x_axis_units), default=0.0, do_not_log=just_read)
53  call get_param(param_file, mdl,"T_REF",t_ref,'Reference temperature',units='C', &
54  default=10., do_not_log=just_read)
55  call get_param(param_file, mdl,"DTDZ",dtdz,'Temperature stratification',units='C/m', &
56  default=0.0, do_not_log=just_read)
57  call get_param(param_file, mdl,"DELTA_T",delta_t,'Temperature difference across baroclinic zone', &
58  units='C', default=0.0, do_not_log=just_read)
59  call get_param(param_file, mdl,"DTDX",dtdx,'Meridional temperature difference', &
60  units='C/'//trim(g%x_axis_units), default=0.0, do_not_log=just_read)
61  call get_param(param_file, mdl,"L_ZONE",l_zone,'Width of baroclinic zone', &
62  units=g%x_axis_units, default=0.5*g%len_lat, do_not_log=just_read)
63  call closeparameterblock(param_file)
64 
65 end subroutine bcz_params
66 
67 !> Initialization of temperature and salinity with the baroclinic zone initial conditions
68 subroutine baroclinic_zone_init_temperature_salinity(T, S, h, G, param_file, &
69  just_read_params)
70  type(ocean_grid_type), intent(in) :: G !< Grid structure
71  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: T !< Potential temperature [deg C]
72  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: S !< Salinity [ppt]
73  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness
74  type(param_file_type), intent(in) :: param_file !< Parameter file handle
75  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
76  !! only read parameters without changing h.
77 
78  integer :: i, j, k, is, ie, js, je, nz
79  real :: T_ref, dTdz, dTdx, delta_T ! Parameters describing temperature distribution
80  real :: S_ref, dSdz, dSdx, delta_S ! Parameters describing salinity distribution
81  real :: L_zone ! Width of baroclinic zone
82  real :: zc, zi, x, xd, xs, y, yd, fn
83  real :: PI ! 3.1415926... calculated as 4*atan(1)
84  logical :: just_read ! If true, just read parameters but set nothing.
85 
86  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
87  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
88 
89  call bcz_params(g, param_file, s_ref, dsdz, delta_s, dsdx, t_ref, dtdz, delta_t, dtdx, l_zone, just_read_params)
90 
91  if (just_read) return ! All run-time parameters have been read, so return.
92 
93  t(:,:,:) = 0.
94  s(:,:,:) = 0.
95  pi = 4.*atan(1.)
96 
97  do j = g%jsc,g%jec ; do i = g%isc,g%iec
98  zi = -g%bathyT(i,j)
99  x = g%geoLonT(i,j) - (g%west_lon + 0.5*g%len_lon) ! Relative to center of domain
100  xd = x / g%len_lon ! -1/2 < xd 1/2
101  y = g%geoLatT(i,j) - (g%south_lat + 0.5*g%len_lat) ! Relative to center of domain
102  yd = y / g%len_lat ! -1/2 < yd 1/2
103  if (l_zone/=0.) then
104  xs = min(1., max(-1., x/l_zone)) ! -1 < ys < 1
105  fn = sin((0.5*pi)*xs)
106  else
107  xs = sign(1., x) ! +/- 1
108  fn = xs
109  endif
110  do k = nz, 1, -1
111  zc = zi + 0.5*h(i,j,k) ! Position of middle of cell
112  zi = zi + h(i,j,k) ! Top interface position
113  t(i,j,k) = t_ref + dtdz * zc & ! Linear temperature stratification
114  + dtdx * x & ! Linear gradient
115  + delta_t * fn ! Smooth fn of width L_zone
116  s(i,j,k) = s_ref + dsdz * zc & ! Linear temperature stratification
117  + dsdx * x & ! Linear gradient
118  + delta_s * fn ! Smooth fn of width L_zone
119  enddo
120  end do ; end do
121 
123 
124 !> \namespace baroclinic_zone_initialization
125 !!
126 !! \section section_baroclinic_zone Description of the baroclinic zone initial conditions
127 !!
128 !! yada yada yada
129 
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Provides the ocean grid type.
Definition: MOM_grid.F90:2
subroutine bcz_params(G, param_file, S_ref, dSdz, delta_S, dSdx, T_ref, dTdz, delta_T, dTdx, L_zone, just_read_params)
Reads the parameters unique to this module.
Initial conditions for an idealized baroclinic zone.
subroutine, public closeparameterblock(CS)
subroutine, public openparameterblock(CS, blockName, desc)
subroutine, public baroclinic_zone_init_temperature_salinity(T, S, h, G, param_file, just_read_params)
Initialization of temperature and salinity with the baroclinic zone initial conditions.
character(len=40) mdl
This module&#39;s name.