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)
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 [m s-1]
69  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: v !< j-component of velocity [m s-1]
70  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness [H ~> m or kg m-2]
71 
72  real :: x, y, x0, y0
73  real :: val1, val2, val3, val4
74  integer :: i, j, k, is, ie, js, je, nz
75 
76  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
77 
78  x0 = 2.0*g%len_lon/3.0
79  y0 = 0.0
80  val1 = 0.395
81  val2 = 0.771*(val1*val1)
82 
83  v(:,:,:) = 0.0
84  u(:,:,:) = 0.0
85 
86  do j = g%jsc,g%jec ; do i = g%isc-1,g%iec+1
87  do k = 1, nz
88  x = 0.5*(g%geoLonT(i+1,j)+g%geoLonT(i,j))-x0
89  y = 0.5*(g%geoLatT(i+1,j)+g%geoLatT(i,j))-y0
90  val3 = exp(-val1*x)
91  val4 = val2*((2.0*val3/(1.0+(val3*val3)))**2)
92  u(i,j,k) = 0.25*val4*(6.0*y*y-9.0) * exp(-0.5*y*y)
93  enddo
94  enddo ; enddo
95  do j = g%jsc-1,g%jec+1 ; do i = g%isc,g%iec
96  do k = 1, nz
97  x = 0.5*(g%geoLonT(i,j+1)+g%geoLonT(i,j))-x0
98  y = 0.5*(g%geoLatT(i,j+1)+g%geoLatT(i,j))-y0
99  val3 = exp(-val1*x)
100  val4 = val2*((2.0*val3/(1.0+(val3*val3)))**2)
101  v(i,j,k) = 2.0*val4*y*(-2.0*val1*tanh(val1*x)) * exp(-0.5*y*y)
102  enddo
103  enddo ; enddo
104 
105 end subroutine soliton_initialize_velocity
106 
107 
108 !> \namespace soliton_initialization
109 !!
110 !! \section section_soliton Description of the equatorial Rossby soliton initial
111 !! conditions
112 !!
113 
114 end module soliton_initialization
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_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
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
soliton_initialization
Initial conditions for the Equatorial Rossby soliton test (Boyd).
Definition: soliton_initialization.F90:2
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
regrid_consts
Contains constants for interpreting input parameters that control regridding.
Definition: regrid_consts.F90:2
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:86
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_eos::calculate_density_derivs
Calculate the derivatives of density with temperature and salinity from T, S, and P.
Definition: MOM_EOS.F90:70
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
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:60