baroclinic_zone_initialization.F90
1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> Initial conditions for an idealized baroclinic zone
6module baroclinic_zone_initialization
7
8use mom_file_parser, only : get_param, log_version, param_file_type
9use mom_file_parser, only : openparameterblock, closeparameterblock
10use mom_grid, only : ocean_grid_type
11use mom_unit_scaling, only : unit_scale_type
12use mom_verticalgrid, only : verticalgrid_type
13
14implicit none ; private
15
16#include <MOM_memory.h>
17#include "version_variable.h"
18
19! Private (module-wise) parameters
20character(len=40) :: mdl = "baroclinic_zone_initialization" !< This module's name.
21
23
24! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
25! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
26! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
27! vary with the Boussinesq approximation, the Boussinesq variant is given first.
28
29contains
30
31!> Reads the parameters unique to this module
32subroutine bcz_params(G, GV, US, param_file, S_ref, dSdz, delta_S, dSdx, T_ref, dTdz, &
33 delta_T, dTdx, L_zone, just_read)
34 type(ocean_grid_type), intent(in) :: G !< Grid structure
35 type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
36 type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
37 type(param_file_type), intent(in) :: param_file !< Parameter file handle
38 real, intent(out) :: S_ref !< Reference salinity [S ~> ppt]
39 real, intent(out) :: dSdz !< Salinity stratification [S Z-1 ~> ppt m-1]
40 real, intent(out) :: delta_S !< Salinity difference across baroclinic zone [S ~> ppt]
41 real, intent(out) :: dSdx !< Linear salinity gradient, often in [S km-1 ~> ppt km-1]
42 !! or [S degrees_E-1 ~> ppt degrees_E-1], depending on
43 !! the value of G%x_axis_units
44 real, intent(out) :: T_ref !< Reference temperature [C ~> degC]
45 real, intent(out) :: dTdz !< Temperature stratification [C Z-1 ~> degC m-1]
46 real, intent(out) :: delta_T !< Temperature difference across baroclinic zone [C ~> degC]
47 real, intent(out) :: dTdx !< Linear temperature gradient, often in [C km-1 ~> degC km-1]
48 !! or [C degrees_E-1 ~> degC degrees_E-1], depending on
49 !! the value of G%x_axis_units
50 real, intent(out) :: L_zone !< Width of baroclinic zone, often in [km] or [degrees_N],
51 !! depending on the value of G%y_axis_units
52 logical, intent(in) :: just_read !< If true, this call will
53 !! only read parameters without changing h.
54
55 if (.not.just_read) &
56 call log_version(param_file, mdl, version, 'Initialization of an analytic baroclinic zone')
57 call openparameterblock(param_file,'BCZIC')
58 call get_param(param_file, mdl, "S_REF", s_ref, 'Reference salinity', &
59 units='ppt', default=35., scale=us%ppt_to_S, do_not_log=just_read)
60 call get_param(param_file, mdl, "DSDZ", dsdz, 'Salinity stratification', &
61 units='ppt m-1', default=0.0, scale=us%ppt_to_S*us%Z_to_m, do_not_log=just_read)
62 call get_param(param_file, mdl, "DELTA_S",delta_s, 'Salinity difference across baroclinic zone', &
63 units='ppt', default=0.0, scale=us%ppt_to_S, do_not_log=just_read)
64 call get_param(param_file, mdl, "DSDX", dsdx,'Meridional salinity difference', &
65 units='ppt '//trim(g%x_ax_unit_short)//'-1', default=0.0, scale=us%ppt_to_S, do_not_log=just_read)
66 call get_param(param_file, mdl, "T_REF", t_ref, 'Reference temperature', &
67 units='degC', default=10., scale=us%degC_to_C, do_not_log=just_read)
68 call get_param(param_file, mdl, "DTDZ", dtdz, 'Temperature stratification', &
69 units='degC m-1', default=0.0, scale=us%degC_to_C*us%Z_to_m, do_not_log=just_read)
70 call get_param(param_file, mdl, "DELTA_T", delta_t,'Temperature difference across baroclinic zone', &
71 units='degC', default=0.0, scale=us%degC_to_C, do_not_log=just_read)
72 call get_param(param_file, mdl, "DTDX", dtdx,'Meridional temperature difference', &
73 units='degC '//trim(g%x_ax_unit_short)//'-1', default=0.0, scale=us%degC_to_C, do_not_log=just_read)
74 call get_param(param_file, mdl, "L_ZONE", l_zone, 'Width of baroclinic zone', &
75 units=g%y_ax_unit_short, default=0.5*g%len_lat, do_not_log=just_read)
76 call closeparameterblock(param_file)
77
78end subroutine bcz_params
79
80!> Initialization of temperature and salinity with the baroclinic zone initial conditions
81subroutine baroclinic_zone_init_temperature_salinity(T, S, h, depth_tot, G, GV, US, param_file, &
82 just_read)
83 type(ocean_grid_type), intent(in) :: g !< Grid structure
84 type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
85 type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
86 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
87 intent(out) :: t !< Potential temperature [C ~> degC]
88 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
89 intent(out) :: s !< Salinity [S ~> ppt]
90 real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
91 intent(in) :: h !< The model thicknesses [Z ~> m]
92 real, dimension(SZI_(G),SZJ_(G)), &
93 intent(in) :: depth_tot !< The nominal total depth of the ocean [Z ~> m]
94 type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
95 !! to parse for model parameter values.
96 logical, intent(in) :: just_read !< If true, this call will only read
97 !! parameters without changing T & S.
98
99 integer :: i, j, k, is, ie, js, je, nz
100 real :: t_ref, delta_t ! Parameters describing temperature distribution [C ~> degC]
101 real :: dtdz ! Vertical temperature gradients [C Z-1 ~> degC m-1]
102 real :: dtdx ! Zonal temperature gradients [C axis_units-1 ~> degC axis_units-1]
103 real :: s_ref, delta_s ! Parameters describing salinity distribution [S ~> ppt]
104 real :: dsdz ! Vertical salinity gradients [S Z-1 ~> ppt m-1]
105 real :: dsdx ! Zonal salinity gradients [S axis_units-1 ~> ppt axis_units-1]
106 real :: l_zone ! Width of baroclinic zone, often in [km] or [degrees_N], depending
107 ! on the value of G%y_axis_units
108 real :: zc, zi ! Depths in depth units [Z ~> m]
109 real :: x ! X-position relative to the domain center [degrees_E] or [km] or [m]
110 real :: y ! Y-position relative to the domain center [degrees_N] or [km] or [m]
111 real :: fn ! A smooth function based on the position in the baroclinic zone [nondim]
112 real :: xs, xd, yd ! Fractional x- and y-positions relative to the domain extent [nondim]
113 real :: pi ! 3.1415926... calculated as 4*atan(1) [nondim]
114
115 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = gv%ke
116
117 call bcz_params(g, gv, us, param_file, s_ref, dsdz, delta_s, dsdx, t_ref, dtdz, &
118 delta_t, dtdx, l_zone, just_read)
119
120 if (just_read) return ! All run-time parameters have been read, so return.
121
122 t(:,:,:) = 0.
123 s(:,:,:) = 0.
124 pi = 4.*atan(1.)
125
126 do j = g%jsc,g%jec ; do i = g%isc,g%iec
127 zi = -depth_tot(i,j)
128 x = g%geoLonT(i,j) - (g%west_lon + 0.5*g%len_lon) ! Relative to center of domain
129 xd = x / g%len_lon ! -1/2 < xd 1/2
130 y = g%geoLatT(i,j) - (g%south_lat + 0.5*g%len_lat) ! Relative to center of domain
131 yd = y / g%len_lat ! -1/2 < yd 1/2
132 if (l_zone/=0.) then
133 xs = min(1., max(-1., x/l_zone)) ! -1 < ys < 1
134 fn = sin((0.5*pi)*xs)
135 else
136 xs = sign(1., x) ! +/- 1
137 fn = xs
138 endif
139 do k = nz, 1, -1
140 zc = zi + 0.5*h(i,j,k) ! Position of middle of cell
141 zi = zi + h(i,j,k) ! Top interface position
142 t(i,j,k) = t_ref + dtdz * zc & ! Linear temperature stratification
143 + dtdx * x & ! Linear gradient
144 + delta_t * fn ! Smooth fn of width L_zone
145 s(i,j,k) = s_ref + dsdz * zc & ! Linear temperature stratification
146 + dsdx * x & ! Linear gradient
147 + delta_s * fn ! Smooth fn of width L_zone
148 enddo
149 enddo ; enddo
150
151end subroutine baroclinic_zone_init_temperature_salinity
152
153!> \namespace baroclinic_zone_initialization
154!!
155!! \section section_baroclinic_zone Description of the baroclinic zone initial conditions
156!!
157!! yada yada yada
158
159end module baroclinic_zone_initialization