MOM6
soliton_initialization.F90
1 !> Initial conditions for the Equatorial Rossby soliton test (Boyd).
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_mesg, mom_error, fatal, is_root_pe
8 use mom_get_input, only : directories
9 use mom_grid, only : ocean_grid_type
14 use regrid_consts, only : coordinatemode, default_coordinate_mode
15 use regrid_consts, only : regridding_layer, regridding_zstar
16 use regrid_consts, only : regridding_rho, regridding_sigma
17 
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
22 ! Private (module-wise) parameters
23 character(len=40) :: mdl = "soliton_initialization" !< This module's name.
24 
25 public soliton_initialize_thickness
26 public soliton_initialize_velocity
27 
28 contains
29 
30 !> Initialization of thicknesses in Equatorial Rossby soliton test
31 subroutine soliton_initialize_thickness(h, G, GV, US)
32  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
33  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
34  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
35  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
36  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
37 
38  integer :: i, j, k, is, ie, js, je, nz
39  real :: x, y, x0, y0
40  real :: val1, val2, val3, val4
41  character(len=40) :: verticalCoordinate
42 
43  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
44 
45  call mom_mesg("soliton_initialization.F90, soliton_initialize_thickness: setting thickness")
46 
47  x0 = 2.0*g%len_lon/3.0
48  y0 = 0.0
49  val1 = 0.395
50  val2 = us%m_to_Z * 0.771*(val1*val1)
51 
52  do j = g%jsc,g%jec ; do i = g%isc,g%iec
53  do k = 1, nz
54  x = g%geoLonT(i,j)-x0
55  y = g%geoLatT(i,j)-y0
56  val3 = exp(-val1*x)
57  val4 = val2 * ( 2.0*val3 / (1.0 + (val3*val3)) )**2
58  h(i,j,k) = gv%Z_to_H * (0.25*val4*(6.0*y*y + 3.0) * exp(-0.5*y*y) + g%bathyT(i,j))
59  enddo
60  enddo ; enddo
61 
62 end subroutine soliton_initialize_thickness
63 
64 
65 !> Initialization of u and v in the equatorial Rossby soliton test
66 subroutine soliton_initialize_velocity(u, v, h, G, US)
67  type(ocean_grid_type), intent(in) :: G !< Grid structure
68  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), intent(out) :: u !< i-component of velocity [L T-1 ~> m s-1]
69  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: v !< j-component of velocity [L T-1 ~> m s-1]
70  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness [H ~> m or kg m-2]
71  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
72 
73  ! Local variables
74  real :: x, x0 ! Positions in the same units as geoLonT.
75  real :: y, y0 ! Positions in the same units as geoLatT.
76  real :: val1 ! A zonal decay scale in the inverse of the units of geoLonT.
77  real :: val2 ! An overall velocity amplitude [L T-1 ~> m s-1]
78  real :: val3 ! A decay factor [nondim]
79  real :: val4 ! The local velocity amplitude [L T-1 ~> m s-1]
80  integer :: i, j, k, is, ie, js, je, nz
81 
82  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
83 
84  x0 = 2.0*g%len_lon/3.0
85  y0 = 0.0
86  val1 = 0.395
87  val2 = us%m_s_to_L_T * 0.771*(val1*val1)
88 
89  v(:,:,:) = 0.0
90  u(:,:,:) = 0.0
91 
92  do j = g%jsc,g%jec ; do i = g%isc-1,g%iec+1
93  do k = 1, nz
94  x = 0.5*(g%geoLonT(i+1,j)+g%geoLonT(i,j))-x0
95  y = 0.5*(g%geoLatT(i+1,j)+g%geoLatT(i,j))-y0
96  val3 = exp(-val1*x)
97  val4 = val2*((2.0*val3/(1.0+(val3*val3)))**2)
98  u(i,j,k) = 0.25*val4*(6.0*y*y-9.0) * exp(-0.5*y*y)
99  enddo
100  enddo ; enddo
101  do j = g%jsc-1,g%jec+1 ; do i = g%isc,g%iec
102  do k = 1, nz
103  x = 0.5*(g%geoLonT(i,j+1)+g%geoLonT(i,j))-x0
104  y = 0.5*(g%geoLatT(i,j+1)+g%geoLatT(i,j))-y0
105  val3 = exp(-val1*x)
106  val4 = val2*((2.0*val3/(1.0+(val3*val3)))**2)
107  v(i,j,k) = 2.0*val4*y*(-2.0*val1*tanh(val1*x)) * exp(-0.5*y*y)
108  enddo
109  enddo ; enddo
110 
111 end subroutine soliton_initialize_velocity
112 
113 
114 !> \namespace soliton_initialization
115 !!
116 !! \section section_soliton Description of the equatorial Rossby soliton initial
117 !! conditions
118 !!
119 
120 end module soliton_initialization
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Reads the only Fortran name list needed to boot-strap the model.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Definition: MOM_grid.F90:2
The MOM6 facility to parse input files for runtime parameters.
Initial conditions for the Equatorial Rossby soliton test (Boyd).
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Definition: MOM_EOS.F90:81
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Container for paths and parameter file names.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Contains constants for interpreting input parameters that control regridding.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108