MOM6
soliton_initialization.F90
Go to the documentation of this file.
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 
8 use mom_get_input, only : directories
9 use mom_grid, only : ocean_grid_type
10 use mom_io, only : close_file, fieldtype, file_exists
11 use mom_io, only : open_file, read_data, read_axis_data, single_file
12 use mom_io, only : write_field, slasher, vardesc
19 
20 implicit none ; private
21 
22 #include <MOM_memory.h>
23 
24 ! Private (module-wise) parameters
25 character(len=40) :: mdl = "soliton_initialization" !< This module's name.
26 
29 
30 contains
31 
32 !> Initialization of thicknesses in Equatorial Rossby soliton test
33 subroutine soliton_initialize_thickness(h, G)
34  type(ocean_grid_type), intent(in) :: G !< Grid structure
35  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: h !< Thickness
36 
37  integer :: i, j, k, is, ie, js, je, nz
38  real :: x, y, x0, y0
39  real :: val1, val2, val3, val4
40  character(len=40) :: verticalCoordinate
41 
42  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
43 
44  call mom_mesg("soliton_initialization.F90, soliton_initialize_thickness: setting thickness")
45 
46  x0 = 2.0*g%len_lon/3.0
47  y0 = 0.0
48  val1 = 0.395
49  val2 = 0.771*(val1*val1)
50 
51  do j = g%jsc,g%jec ; do i = g%isc,g%iec
52  do k = 1, nz
53  x = g%geoLonT(i,j)-x0
54  y = g%geoLatT(i,j)-y0
55  val3 = exp(-val1*x)
56  val4 = val2*((2.0*val3/(1.0+(val3*val3)))**2)
57  h(i,j,k) = 0.25*val4*(6.0*y*y+3.0)* &
58  exp(-0.5*y*y)
59  enddo
60  end do ; end do
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]
69  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), intent(out) :: v !< j-component of velocity [m/s]
70  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(in) :: h !< Thickness [H]
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)* &
93  exp(-0.5*y*y)
94  enddo
95  enddo ; enddo
96  do j = g%jsc-1,g%jec+1 ; do i = g%isc,g%iec
97  do k = 1, nz
98  x = 0.5*(g%geoLonT(i,j+1)+g%geoLonT(i,j))-x0
99  y = 0.5*(g%geoLatT(i,j+1)+g%geoLatT(i,j))-y0
100  val3 = exp(-val1*x)
101  val4 = val2*((2.0*val3/(1.0+(val3*val3)))**2)
102  v(i,j,k) = 2.0*val4*y*(-2.0*val1*tanh(val1*x))* &
103  exp(-0.5*y*y)
104  enddo
105  enddo ; enddo
106 
107 end subroutine soliton_initialize_velocity
108 
109 
110 !> \namespace soliton_initialization
111 !!
112 !! \section section_soliton Description of the equatorial Rossby soliton initial
113 !! conditions
114 !!
115 
116 end module soliton_initialization
integer, parameter regridding_layer
Layer mode.
subroutine, public read_axis_data(filename, axis_name, var)
Definition: MOM_io.F90:424
integer function coordinatemode(string)
Parse a string parameter specifying the coordinate mode and return the appropriate enumerated integer...
integer, parameter regridding_sigma
Sigma coordinates.
subroutine, public soliton_initialize_thickness(h, G)
Initialization of thicknesses in Equatorial Rossby soliton test.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:45
Provides the ocean grid type.
Definition: MOM_grid.F90:2
character(len= *), parameter default_coordinate_mode
Default coordinate mode.
This module contains I/O framework code.
Definition: MOM_io.F90:2
Initial conditions for the Equatorial Rossby soliton test (Boyd).
subroutine, public calculate_density_derivs(T, S, pressure, drho_dT, drho_dS, start, npts, EOS)
Calls the appropriate subroutine to calculate density derivatives for 1-D array inputs.
Definition: MOM_EOS.F90:214
logical function, public is_root_pe()
subroutine, public soliton_initialize_velocity(u, v, h, G)
Initialization of u and v in the equatorial Rossby soliton test.
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:51
integer, parameter regridding_zstar
z* coordinates
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
Contains constants for interpreting input parameters that control regridding.
subroutine, public mom_error(level, message, all_print)
integer, parameter regridding_rho
Target interface densities.
character(len=40) mdl
This module&#39;s name.
A control structure for the equation of state.
Definition: MOM_EOS.F90:55