MOM6
MOM_state_initialization.F90 File Reference
#include <MOM_memory.h>
Include dependency graph for MOM_state_initialization.F90:

Go to the source code of this file.

Modules

module  mom_state_initialization
 Initialize state variables, u, v, h, T and S.
 

Functions/Subroutines

subroutine, public mom_state_initialization::mom_initialize_state (u, v, h, tv, Time, G, GV, PF, dirs, restart_CS, ALE_CSp, tracer_Reg, sponge_CSp, ALE_sponge_CSp, OBC, Time_in)
 
subroutine mom_state_initialization::initialize_thickness_from_file (h, G, GV, param_file, file_has_thickness, just_read_params)
 This subroutine reads the layer thicknesses or interface heights from a file. More...
 
subroutine mom_state_initialization::adjustetatofitbathymetry (G, GV, eta, h)
 Adjust interface heights to fit the bathymetry and diagnose layer thickness. If the bottom most interface is below the topography then the bottom-most layers are contracted to GVAngstrom_z. If the bottom most interface is above the topography then the entire column is dilated (expanded) to fill the void. More...
 
subroutine mom_state_initialization::initialize_thickness_uniform (h, G, GV, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_thickness_search
 
subroutine mom_state_initialization::convert_thickness (h, G, GV, tv)
 
subroutine mom_state_initialization::depress_surface (h, G, GV, param_file, tv, just_read_params)
 
subroutine mom_state_initialization::trim_for_ice (PF, G, GV, ALE_CSp, tv, h, just_read_params)
 Adjust the layer thicknesses by cutting away the top of each model column at the depth where the hydrostatic pressure matches an imposed surface pressure read from file. More...
 
subroutine mom_state_initialization::cut_off_column_top (nk, tv, Rho0, G_earth, depth, min_thickness, T, T_t, T_b, S, S_t, S_b, p_surf, h, remap_CS)
 Adjust the layer thicknesses by cutting away the top at the depth where the hydrostatic pressure matches p_surf. More...
 
subroutine mom_state_initialization::initialize_velocity_from_file (u, v, G, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_velocity_zero (u, v, G, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_velocity_uniform (u, v, G, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_velocity_circular (u, v, G, param_file, just_read_params)
 
real function my_psi (ig, jg)
 
subroutine mom_state_initialization::initialize_temp_salt_from_file (T, S, G, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_temp_salt_from_profile (T, S, G, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_temp_salt_fit (T, S, G, GV, param_file, eqn_of_state, P_Ref, just_read_params)
 
subroutine mom_state_initialization::initialize_temp_salt_linear (T, S, G, param_file, just_read_params)
 
subroutine mom_state_initialization::initialize_sponges_file (G, GV, use_temperature, tv, param_file, CSp)
 This subroutine sets the inverse restoration time (Idamp), and the values towards which the interface heights and an arbitrary number of tracers should be restored within each sponge. The interface height is always subject to damping, and must always be the first registered field. More...
 
subroutine mom_state_initialization::set_velocity_depth_max (G)
 
subroutine mom_state_initialization::compute_global_grid_integrals (G)
 
subroutine mom_state_initialization::set_velocity_depth_min (G)
 
subroutine mom_state_initialization::mom_temp_salt_initialize_from_z (h, tv, G, GV, PF, just_read_params)
 This subroutine determines the isopycnal or other coordinate interfaces and layer potential temperatures and salinities directly from a z-space file on a latitude-longitude grid. More...
 
subroutine mom_state_initialization::mom_state_init_tests (G, GV, tv)
 Run simple unit tests. More...
 

Variables

character(len=40) mom_state_initialization::mdl = "MOM_state_initialization"
 

Function/Subroutine Documentation

◆ my_psi()

real function initialize_velocity_circular::my_psi ( integer  ig,
integer  jg 
)
private

Definition at line 1317 of file MOM_state_initialization.F90.

Referenced by mom_state_initialization::initialize_velocity_circular().

1317  integer :: ig, jg
1318  real :: x, y, r
1319  x = 2.0*(g%geoLonBu(ig,jg)-g%west_lon)/g%len_lon-1.0 ! -1<x<1
1320  y = 2.0*(g%geoLatBu(ig,jg)-g%south_lat)/g%len_lat-1.0 ! -1<y<1
1321  r = sqrt( x**2 + y**2 ) ! Circulat stream fn nis fn of radius only
1322  r = min(1.0,r) ! Flatten stream function in corners of box
1323  my_psi = 0.5*(1.0 - cos(dpi*r))
1324  my_psi = my_psi * (circular_max_u*g%len_lon*1e3/dpi) ! len_lon is in km
real function my_psi(ig, jg)
Here is the caller graph for this function: