MOM6
phillips_initialization Module Reference

Detailed Description

Initialization for the "Phillips" channel configuration.

By Robert Hallberg, April 1994 - June 2002

This subroutine initializes the fields for the simulations. The one argument passed to initialize, Time, is set to the current time of the simulation. The fields which are initialized here are: u - Zonal velocity [m s-1]. v - Meridional velocity [m s-1]. h - Layer thickness in m. (Must be positive.) D - Basin depth in m. (Must be positive.) f - The Coriolis parameter [s-1]. g - The reduced gravity at each interface [m s-2] Rlay - Layer potential density (coordinate variable) [kg m-3]. If ENABLE_THERMODYNAMICS is defined: T - Temperature [degC]. S - Salinity [ppt]. If SPONGE is defined: A series of subroutine calls are made to set up the damping rates and reference profiles for all variables that are damped in the sponge. Any user provided tracer code is also first linked through this subroutine.

Forcing-related fields (taux, tauy, buoy, ustar, etc.) are set in MOM_surface_forcing.F90.

These variables are all set in the set of subroutines (in this file) Phillips_initialize_thickness, Phillips_initialize_velocity, Phillips_initialize_topography and Phillips_initialize_sponges that seet up fields that are specific to the Phillips instability test case.

Functions/Subroutines

subroutine, public phillips_initialize_thickness (h, G, GV, US, param_file, just_read_params)
 Initialize the thickness field for the Phillips model test case. More...
 
subroutine, public phillips_initialize_velocity (u, v, G, GV, US, param_file, just_read_params)
 Initialize the velocity fields for the Phillips model test case. More...
 
subroutine, public phillips_initialize_sponges (G, GV, US, tv, param_file, CSp, h)
 Sets up the the inverse restoration time (Idamp), and the values towards which the interface heights and an arbitrary number of tracers should be restored within each sponge for the Phillips model test case. More...
 
real function sech (x)
 sech calculates the hyperbolic secant. More...
 
subroutine, public phillips_initialize_topography (D, G, param_file, max_depth, US)
 Initialize topography. More...
 

Function/Subroutine Documentation

◆ phillips_initialize_sponges()

subroutine, public phillips_initialization::phillips_initialize_sponges ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(thermo_var_ptrs), intent(in)  tv,
type(param_file_type), intent(in)  param_file,
type(sponge_cs), pointer  CSp,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h 
)

Sets up the the inverse restoration time (Idamp), and the values towards which the interface heights and an arbitrary number of tracers should be restored within each sponge for the Phillips model test case.

Parameters
[in]gThe ocean's grid structure.
[in]gvVertical grid structure
[in]usA dimensional unit scaling type
[in]tvA structure containing pointers to any available thermodynamic fields, potential temperature and salinity or mixed layer density. Absent fields have NULL ptrs.
[in]param_fileA structure indicating the open file to parse for model parameter values.
cspA pointer that is set to point to the control structure for the sponge module.
[in]hThickness field [H ~> m or kg m-2].

Definition at line 202 of file Phillips_initialization.F90.

202  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
203  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure
204  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
205  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers
206  !! to any available thermodynamic
207  !! fields, potential temperature and
208  !! salinity or mixed layer density.
209  !! Absent fields have NULL ptrs.
210  type(param_file_type), intent(in) :: param_file !< A structure indicating the
211  !! open file to parse for model
212  !! parameter values.
213  type(sponge_CS), pointer :: CSp !< A pointer that is set to point to
214  !! the control structure for the
215  !! sponge module.
216  real, intent(in), dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h !< Thickness field [H ~> m or kg m-2].
217 
218  ! Local variables
219  real :: eta0(SZK_(G)+1) ! The 1-d nominal positions of the interfaces.
220  real :: eta(SZI_(G),SZJ_(G),SZK_(G)+1) ! A temporary array for eta [Z ~> m].
221  real :: temp(SZI_(G),SZJ_(G),SZK_(G)) ! A temporary array for other variables.
222  real :: Idamp(SZI_(G),SZJ_(G)) ! The inverse damping rate [s-1].
223  real :: eta_im(SZJ_(G),SZK_(G)+1) ! A temporary array for zonal-mean eta [Z ~> m].
224  real :: Idamp_im(SZJ_(G)) ! The inverse zonal-mean damping rate [s-1].
225  real :: damp_rate ! The inverse zonal-mean damping rate [s-1].
226  real :: jet_width ! The width of the zonal mean jet, in km.
227  real :: jet_height ! The interface height scale associated with the zonal-mean jet [Z ~> m].
228  real :: y_2 ! The y-position relative to the channel center, in km.
229  real :: half_strat ! The fractional depth where the straficiation is centered [nondim].
230  real :: half_depth ! The depth where the stratification is centered [Z ~> m].
231  character(len=40) :: mdl = "Phillips_initialize_sponges" ! This subroutine's name.
232 
233  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
234  logical, save :: first_call = .true.
235 
236  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
237  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
238 
239  eta(:,:,:) = 0.0 ; temp(:,:,:) = 0.0 ; idamp(:,:) = 0.0
240  eta_im(:,:) = 0.0 ; idamp_im(:) = 0.0
241 
242  if (first_call) call log_version(param_file, mdl, version)
243  first_call = .false.
244  call get_param(param_file, mdl, "HALF_STRAT_DEPTH", half_strat, &
245  "The fractional depth where the stratificaiton is centered.", &
246  units="nondim", default = 0.5)
247  call get_param(param_file, mdl, "SPONGE_RATE", damp_rate, &
248  "The rate at which the zonal-mean sponges damp.", units="s-1", &
249  default = 1.0/(10.0*86400.0))
250 
251  call get_param(param_file, mdl, "JET_WIDTH", jet_width, &
252  "The width of the zonal-mean jet.", units="km", &
253  fail_if_missing=.true.)
254  call get_param(param_file, mdl, "JET_HEIGHT", jet_height, &
255  "The interface height scale associated with the "//&
256  "zonal-mean jet.", units="m", scale=us%m_to_Z, &
257  fail_if_missing=.true.)
258 
259  half_depth = g%max_depth*half_strat
260  eta0(1) = 0.0 ; eta0(nz+1) = -g%max_depth
261  do k=2,1+nz/2 ; eta0(k) = -half_depth*(2.0*(k-1)/real(nz)) ; enddo
262  do k=2+nz/2,nz+1
263  eta0(k) = -g%max_depth - 2.0*(g%max_depth-half_depth) * ((k-(nz+1))/real(nz))
264  enddo
265 
266  do j=js,je
267  idamp_im(j) = damp_rate
268  eta_im(j,1) = 0.0 ; eta_im(j,nz+1) = -g%max_depth
269  enddo
270  do k=2,nz ; do j=js,je
271  y_2 = g%geoLatT(is,j) - g%south_lat - 0.5*g%len_lat
272  eta_im(j,k) = eta0(k) + jet_height * tanh(y_2 / jet_width)
273 ! jet_height * atan(y_2 / jet_width)
274  if (eta_im(j,k) > 0.0) eta_im(j,k) = 0.0
275  if (eta_im(j,k) < -g%max_depth) eta_im(j,k) = -g%max_depth
276  enddo ; enddo
277 
278  call initialize_sponge(idamp, eta, g, param_file, csp, gv, idamp_im, eta_im)
279 

References mom_sponge::initialize_sponge().

Referenced by mom_state_initialization::mom_initialize_state().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ phillips_initialize_thickness()

subroutine, public phillips_initialization::phillips_initialize_thickness ( real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(out)  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 
)

Initialize the thickness field for the Phillips model test case.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[out]hThe thickness that is being initialized [H ~> m or kg m-2]
[in]param_fileA structure indicating the open file to parse for model parameter values.
[in]just_read_paramsIf present and true, this call will only read parameters without changing h.

Definition at line 39 of file Phillips_initialization.F90.

39  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
40  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
41  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
42  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
43  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2]
44  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
45  !! to parse for model parameter values.
46  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
47  !! only read parameters without changing h.
48 
49  real :: eta0(SZK_(G)+1) ! The 1-d nominal positions of the interfaces [Z ~> m]
50  real :: eta_im(SZJ_(G),SZK_(G)+1) ! A temporary array for zonal-mean eta [Z ~> m]
51  real :: eta1D(SZK_(G)+1) ! Interface height relative to the sea surface, positive upward [Z ~> m]
52  real :: jet_width ! The width of the zonal-mean jet [km]
53  real :: jet_height ! The interface height scale associated with the zonal-mean jet [Z ~> m]
54  real :: y_2 ! The y-position relative to the center of the domain [km]
55  real :: half_strat ! The fractional depth where the stratification is centered [nondim]
56  real :: half_depth ! The depth where the stratification is centered [Z ~> m]
57  logical :: just_read ! If true, just read parameters but set nothing.
58  character(len=40) :: mdl = "Phillips_initialize_thickness" ! This subroutine's name.
59  integer :: i, j, k, is, ie, js, je, isd, ied, jsd, jed, nz
60 
61  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
62  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
63 
64  eta_im(:,:) = 0.0
65 
66  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
67 
68  if (.not.just_read) call log_version(param_file, mdl, version)
69  call get_param(param_file, mdl, "HALF_STRAT_DEPTH", half_strat, &
70  "The fractional depth where the stratification is centered.", &
71  units="nondim", default = 0.5, do_not_log=just_read)
72  call get_param(param_file, mdl, "JET_WIDTH", jet_width, &
73  "The width of the zonal-mean jet.", units="km", &
74  fail_if_missing=.not.just_read, do_not_log=just_read)
75  call get_param(param_file, mdl, "JET_HEIGHT", jet_height, &
76  "The interface height scale associated with the "//&
77  "zonal-mean jet.", units="m", scale=us%m_to_Z, &
78  fail_if_missing=.not.just_read, do_not_log=just_read)
79 
80  if (just_read) return ! All run-time parameters have been read, so return.
81 
82  half_depth = g%max_depth*half_strat
83  eta0(1) = 0.0 ; eta0(nz+1) = -g%max_depth
84  do k=2,1+nz/2 ; eta0(k) = -half_depth*(2.0*(k-1)/real(nz)) ; enddo
85  do k=2+nz/2,nz+1
86  eta0(k) = -g%max_depth - 2.0*(g%max_depth-half_depth) * ((k-(nz+1))/real(nz))
87  enddo
88 
89  do j=js,je
90  eta_im(j,1) = 0.0 ; eta_im(j,nz+1) = -g%max_depth
91  enddo
92  do k=2,nz ; do j=js,je
93  y_2 = g%geoLatT(is,j) - g%south_lat - 0.5*g%len_lat
94  eta_im(j,k) = eta0(k) + jet_height * tanh(y_2 / jet_width)
95  ! or ... + jet_height * atan(y_2 / jet_width)
96  if (eta_im(j,k) > 0.0) eta_im(j,k) = 0.0
97  if (eta_im(j,k) < -g%max_depth) eta_im(j,k) = -g%max_depth
98  enddo ; enddo
99 
100  do j=js,je ; do i=is,ie
101  ! This sets the initial thickness in [H ~> m or kg m-2] of the layers. The
102  ! thicknesses are set to insure that: 1. each layer is at least an Angstrom thick, and
103  ! 2. the interfaces are where they should be based on the resting depths and interface
104  ! height perturbations, as long at this doesn't interfere with 1.
105  eta1d(nz+1) = -g%bathyT(i,j)
106  do k=nz,1,-1
107  eta1d(k) = eta_im(j,k)
108  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z)) then
109  eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
110  h(i,j,k) = gv%Angstrom_H
111  else
112  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
113  endif
114  enddo
115  enddo ; enddo
116 

◆ phillips_initialize_topography()

subroutine, public phillips_initialization::phillips_initialize_topography ( real, dimension(g%isd:g%ied,g%jsd:g%jed), intent(out)  D,
type(dyn_horgrid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
real, intent(in)  max_depth,
type(unit_scale_type), intent(in), optional  US 
)

Initialize topography.

Parameters
[in]gThe dynamic horizontal grid type
[out]dOcean bottom depth in m or Z if US is present
[in]param_fileParameter file structure
[in]max_depthMaximum model depth in the units of D
[in]usA dimensional unit scaling type

Definition at line 297 of file Phillips_initialization.F90.

297  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
298  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
299  intent(out) :: D !< Ocean bottom depth in m or Z if US is present
300  type(param_file_type), intent(in) :: param_file !< Parameter file structure
301  real, intent(in) :: max_depth !< Maximum model depth in the units of D
302  type(unit_scale_type), optional, intent(in) :: US !< A dimensional unit scaling type
303 
304  ! Local variables
305  real :: m_to_Z ! A dimensional rescaling factor.
306  real :: PI, Htop, Wtop, Ltop, offset, dist
307  real :: x1, x2, x3, x4, y1, y2
308  integer :: i,j,is,ie,js,je
309  character(len=40) :: mdl = "Phillips_initialize_topography" ! This subroutine's name.
310 
311  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
312 
313  pi = 4.0*atan(1.0)
314  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
315 
316  call get_param(param_file, mdl, "PHILLIPS_HTOP", htop, &
317  "The maximum height of the topography.", units="m", scale=m_to_z, &
318  fail_if_missing=.true.)
319 ! Htop=0.375*max_depth ! max height of topog. above max_depth
320  wtop = 0.5*g%len_lat ! meridional width of drake and mount
321  ltop = 0.25*g%len_lon ! zonal width of topographic features
322  offset = 0.1*g%len_lat ! meridional offset from center
323  dist = 0.333*g%len_lon ! distance between drake and mount
324  ! should be longer than Ltop/2
325 
326  y1=g%south_lat+0.5*g%len_lat+offset-0.5*wtop; y2=y1+wtop
327  x1=g%west_lon+0.1*g%len_lon; x2=x1+ltop; x3=x1+dist; x4=x3+3.0/2.0*ltop
328 
329  do i=is,ie ; do j=js,je
330  d(i,j)=0.0
331  if (g%geoLonT(i,j)>x1 .and. g%geoLonT(i,j)<x2) then
332  d(i,j) = htop*sin(pi*(g%geoLonT(i,j)-x1)/(x2-x1))**2
333  if (g%geoLatT(i,j)>y1 .and. g%geoLatT(i,j)<y2) then
334  d(i,j) = d(i,j)*(1-sin(pi*(g%geoLatT(i,j)-y1)/(y2-y1))**2)
335  endif
336  elseif (g%geoLonT(i,j)>x3 .and. g%geoLonT(i,j)<x4 .and. &
337  g%geoLatT(i,j)>y1 .and. g%geoLatT(i,j)<y2) then
338  d(i,j) = 2.0/3.0*htop*sin(pi*(g%geoLonT(i,j)-x3)/(x4-x3))**2 &
339  *sin(pi*(g%geoLatT(i,j)-y1)/(y2-y1))**2
340  endif
341  d(i,j) = max_depth - d(i,j)
342  enddo ; enddo
343 

Referenced by mom_fixed_initialization::mom_initialize_topography().

Here is the caller graph for this function:

◆ phillips_initialize_velocity()

subroutine, public phillips_initialization::phillips_initialize_velocity ( real, dimension( g %isdb: g %iedb, g %jsd: g %jed, gv %ke), intent(out)  u,
real, dimension( g %isd: g %ied, g %jsdb: g %jedb, gv %ke), intent(out)  v,
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 
)

Initialize the velocity fields for the Phillips model test case.

Parameters
[in]gGrid structure
[in]gvVertical grid structure
[out]ui-component of velocity [L T-1 ~> m s-1]
[out]vj-component of velocity [L T-1 ~> m s-1]
[in]usA dimensional unit scaling type
[in]param_fileA structure indicating the open file to parse for modelparameter values.
[in]just_read_paramsIf present and true, this call will only read parameters without changing h.

Definition at line 121 of file Phillips_initialization.F90.

121  type(ocean_grid_type), intent(in) :: G !< Grid structure
122  type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure
123  real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
124  intent(out) :: u !< i-component of velocity [L T-1 ~> m s-1]
125  real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
126  intent(out) :: v !< j-component of velocity [L T-1 ~> m s-1]
127  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
128  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to
129  !! parse for modelparameter values.
130  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
131  !! only read parameters without changing h.
132 
133  real :: jet_width ! The width of the zonal-mean jet [km]
134  real :: jet_height ! The interface height scale associated with the zonal-mean jet [Z ~> m]
135  real :: x_2 ! The x-position relative to the center of the domain [nondim]
136  real :: y_2 ! The y-position relative to the center of the domain [km] or [nondim]
137  real :: velocity_amplitude ! The amplitude of velocity perturbations [L T-1 ~> m s-1]
138  real :: pi ! The ratio of the circumference of a circle to its diameter [nondim]
139  integer :: i, j, k, is, ie, js, je, nz, m
140  logical :: just_read ! If true, just read parameters but set nothing.
141  character(len=40) :: mdl = "Phillips_initialize_velocity" ! This subroutine's name.
142  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
143 
144  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
145 
146  if (.not.just_read) call log_version(param_file, mdl, version)
147  call get_param(param_file, mdl, "VELOCITY_IC_PERTURB_AMP", velocity_amplitude, &
148  "The magnitude of the initial velocity perturbation.", &
149  units="m s-1", default=0.001, scale=us%m_s_to_L_T, do_not_log=just_read)
150  call get_param(param_file, mdl, "JET_WIDTH", jet_width, &
151  "The width of the zonal-mean jet.", units="km", &
152  fail_if_missing=.not.just_read, do_not_log=just_read)
153  call get_param(param_file, mdl, "JET_HEIGHT", jet_height, &
154  "The interface height scale associated with the "//&
155  "zonal-mean jet.", units="m", scale=us%m_to_Z, &
156  fail_if_missing=.not.just_read, do_not_log=just_read)
157 
158  if (just_read) return ! All run-time parameters have been read, so return.
159 
160  u(:,:,:) = 0.0
161  v(:,:,:) = 0.0
162 
163  pi = 4.0*atan(1.0)
164 
165  ! Use thermal wind shear to give a geostrophically balanced flow.
166  do k=nz-1,1 ; do j=js,je ; do i=is-1,ie
167  y_2 = g%geoLatCu(i,j) - g%south_lat - 0.5*g%len_lat
168 ! This uses d/d y_2 atan(y_2 / jet_width)
169 ! u(I,j,k) = u(I,j,k+1) + (1e-3 * jet_height / &
170 ! (US%m_to_L*jet_width * (1.0 + (y_2 / jet_width)**2))) * &
171 ! (2.0 * GV%g_prime(K+1) / (G%CoriolisBu(I,J) + G%CoriolisBu(I,J-1)))
172 ! This uses d/d y_2 tanh(y_2 / jet_width)
173  u(i,j,k) = u(i,j,k+1) + (1e-3 * (jet_height / (us%m_to_L*jet_width)) * &
174  (sech(y_2 / jet_width))**2 ) * &
175  (2.0 * gv%g_prime(k+1) / (g%CoriolisBu(i,j) + g%CoriolisBu(i,j-1)))
176  enddo ; enddo ; enddo
177 
178  do k=1,nz ; do j=js,je ; do i=is-1,ie
179  y_2 = (g%geoLatCu(i,j) - g%south_lat - 0.5*g%len_lat) / g%len_lat
180  x_2 = (g%geoLonCu(i,j) - g%west_lon - 0.5*g%len_lon) / g%len_lon
181  if (g%geoLonCu(i,j) == g%west_lon) then
182  ! This modification is required so that the perturbations are identical for
183  ! symmetric and non-symmetric memory. It is exactly equivalent to
184  ! taking the longitude at the eastern edge of the domain, so that x_2 ~= 0.5.
185  x_2 = ((g%west_lon + g%len_lon*real(g%ieg-(g%isg-1))/real(g%Domain%niglobal)) - &
186  g%west_lon - 0.5*g%len_lon) / g%len_lon
187  endif
188  u(i,j,k) = u(i,j,k) + velocity_amplitude * ((real(k)-0.5)/real(nz)) * &
189  (0.5 - abs(2.0*x_2) + 0.1*abs(cos(10.0*pi*x_2)) - abs(sin(5.0*pi*y_2)))
190  do m=1,10
191  u(i,j,k) = u(i,j,k) + 0.2*velocity_amplitude * ((real(k)-0.5)/real(nz)) * &
192  cos(2.0*m*pi*x_2 + 2*m) * cos(6.0*pi*y_2)
193  enddo
194  enddo ; enddo ; enddo
195 

References sech().

Here is the call graph for this function:

◆ sech()

real function phillips_initialization::sech ( real, intent(in)  x)
private

sech calculates the hyperbolic secant.

Parameters
[in]xInput value.
Returns
Result.

Definition at line 284 of file Phillips_initialization.F90.

284  real, intent(in) :: x !< Input value.
285  real :: sech !< Result.
286 
287  ! This is here to prevent overflows or underflows.
288  if (abs(x) > 228.) then
289  sech = 0.0
290  else
291  sech = 2.0 / (exp(x) + exp(-x))
292  endif

Referenced by phillips_initialize_velocity().

Here is the caller graph for this function: