MOM6
lock_exchange_initialization.F90
1 !> Initialization of the "lock exchange" experiment.
2 !! lock_exchange = A 2-d density driven hydraulic exchange flow.
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
7 use mom_error_handler, only : mom_mesg, mom_error, fatal, is_root_pe
9 use mom_get_input, only : directories
10 use mom_grid, only : ocean_grid_type
15 
16 implicit none ; private
17 
18 #include <MOM_memory.h>
19 
20 public lock_exchange_initialize_thickness
21 
22 contains
23 
24 !> This subroutine initializes layer thicknesses for the lock_exchange experiment.
25 ! -----------------------------------------------------------------------------
26 subroutine lock_exchange_initialize_thickness(h, G, GV, US, param_file, just_read_params)
27  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure.
28  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
29  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
30  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
31  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
32  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
33  !! to parse for model parameter values.
34  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
35  !! only read parameters without changing h.
36 
37  real :: e0(szk_(gv)) ! The resting interface heights [Z ~> m], usually
38  ! negative because it is positive upward.
39  real :: e_pert(szk_(gv)) ! Interface height perturbations, positive upward [Z ~> m].
40  real :: eta1d(szk_(gv)+1)! Interface height relative to the sea surface
41  ! positive upward [Z ~> m].
42  real :: front_displacement ! Vertical displacement acrodd front
43  real :: thermocline_thickness ! Thickness of stratified region
44  logical :: just_read ! If true, just read parameters but set nothing.
45 ! This include declares and sets the variable "version".
46 #include "version_variable.h"
47  character(len=40) :: mdl = "lock_exchange_initialize_thickness" ! This subroutine's name.
48  integer :: i, j, k, is, ie, js, je, nz
49 
50  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
51 
52  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
53 
54  if (.not.just_read) &
55  call mom_mesg(" lock_exchange_initialization.F90, lock_exchange_initialize_thickness: setting thickness", 5)
56 
57  if (.not.just_read) call log_version(param_file, mdl, version, "")
58  call get_param(param_file, mdl, "FRONT_DISPLACEMENT", front_displacement, &
59  "The vertical displacement of interfaces across the front. "//&
60  "A value larger in magnitude that MAX_DEPTH is truncated,", &
61  units="m", fail_if_missing=.not.just_read, do_not_log=just_read, scale=us%m_to_Z)
62  call get_param(param_file, mdl, "THERMOCLINE_THICKNESS", thermocline_thickness, &
63  "The thickness of the thermocline in the lock exchange "//&
64  "experiment. A value of zero creates a two layer system "//&
65  "with vanished layers in between the two inflated layers.", &
66  default=0., units="m", do_not_log=just_read, scale=us%m_to_Z)
67 
68  if (just_read) return ! All run-time parameters have been read, so return.
69 
70  do j=g%jsc,g%jec ; do i=g%isc,g%iec
71  do k=2,nz
72  eta1d(k) = -0.5 * g%max_depth & ! Middle of column
73  - thermocline_thickness * ( (real(k-1))/real(nz) -0.5 ) ! Stratification
74  if (g%geoLonT(i,j)-g%west_lon < 0.5 * g%len_lon) then
75  eta1d(k) = eta1d(k) + 0.5 * front_displacement
76  elseif (g%geoLonT(i,j)-g%west_lon > 0.5 * g%len_lon) then
77  eta1d(k) = eta1d(k) - 0.5 * front_displacement
78  endif
79  enddo
80  eta1d(nz+1) = -g%max_depth ! Force bottom interface to bottom
81  do k=nz,2,-1 ! Make sure interfaces increase upwards
82  eta1d(k) = max( eta1d(k), eta1d(k+1) + gv%Angstrom_Z )
83  enddo
84  eta1d(1) = 0. ! Force bottom interface to bottom
85  do k=2,nz ! Make sure interfaces decrease downwards
86  eta1d(k) = min( eta1d(k), eta1d(k-1) - gv%Angstrom_Z )
87  enddo
88  do k=nz,1,-1
89  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
90  enddo
91  enddo ; enddo
92 
93 end subroutine lock_exchange_initialize_thickness
94 ! -----------------------------------------------------------------------------
95 
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_get_input::directories
Container for paths and parameter file names.
Definition: MOM_get_input.F90:20
mom_tracer_registry
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Definition: MOM_tracer_registry.F90:5
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_get_input
Reads the only Fortran name list needed to boot-strap the model.
Definition: MOM_get_input.F90:6
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
lock_exchange_initialization
Initialization of the "lock exchange" experiment. lock_exchange = A 2-d density driven hydraulic exch...
Definition: lock_exchange_initialization.F90:3
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
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_tracer_registry::tracer_registry_type
Type to carry basic tracer information.
Definition: MOM_tracer_registry.F90:122
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_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