14 use regrid_consts,
only : coordinatemode, default_coordinate_mode
18 implicit none ;
private 20 #include <MOM_memory.h> 23 character(len=40) :: mdl =
"soliton_initialization" 25 public soliton_initialize_thickness
26 public soliton_initialize_velocity
31 subroutine soliton_initialize_thickness(h, G, GV, US)
35 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
38 integer :: i, j, k, is, ie, js, je, nz
40 real :: val1, val2, val3, val4
41 character(len=40) :: verticalcoordinate
43 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
45 call mom_mesg(
"soliton_initialization.F90, soliton_initialize_thickness: setting thickness")
47 x0 = 2.0*g%len_lon/3.0
50 val2 = us%m_to_Z * 0.771*(val1*val1)
52 do j = g%jsc,g%jec ;
do i = g%isc,g%iec
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))
62 end subroutine soliton_initialize_thickness
66 subroutine soliton_initialize_velocity(u, v, h, G, US)
68 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(out) :: u
69 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(out) :: v
70 real,
dimension(SZI_(G),SZJ_(G), SZK_(G)),
intent(in) :: h
80 integer :: i, j, k, is, ie, js, je, nz
82 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
84 x0 = 2.0*g%len_lon/3.0
87 val2 = us%m_s_to_L_T * 0.771*(val1*val1)
92 do j = g%jsc,g%jec ;
do i = g%isc-1,g%iec+1
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
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)
101 do j = g%jsc-1,g%jec+1 ;
do i = g%isc,g%iec
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
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)
111 end subroutine soliton_initialize_velocity
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
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...
Provides subroutines for quantities specific to the equation of state.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
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.