MOM6
neverland_initialization Module Reference

Detailed Description

Initialization for the "Neverland" configuration.

Functions/Subroutines

subroutine, public neverland_initialize_topography (D, G, param_file, max_depth)
 This subroutine sets up the Neverland test case topography. More...
 
real function cosbell (x, L)
 Returns the value of a cosine-bell function evaluated at x/L. More...
 
real function spike (x, L)
 Returns the value of a sin-spike function evaluated at x/L. More...
 
subroutine, public neverland_initialize_thickness (h, G, GV, US, param_file, eqn_of_state, P_ref)
 This subroutine initializes layer thicknesses for the Neverland test case, by finding the depths of interfaces in a specified latitude-dependent temperature profile with an exponentially decaying thermocline on top of a linear stratification. More...
 

Function/Subroutine Documentation

◆ cosbell()

real function neverland_initialization::cosbell ( real, intent(in)  x,
real, intent(in)  L 
)
private

Returns the value of a cosine-bell function evaluated at x/L.

Parameters
[in]xnon-dimensional position
[in]lnon-dimensional width

Definition at line 90 of file Neverland_initialization.F90.

90  real , intent(in) :: x !< non-dimensional position
91  real , intent(in) :: L !< non-dimensional width
92  real :: PI !< 3.1415926... calculated as 4*atan(1)
93 
94  pi = 4.0*atan(1.0)
95  cosbell = 0.5 * (1 + cos(pi*min(abs(x/l),1.0)))

Referenced by neverland_initialize_topography().

Here is the caller graph for this function:

◆ neverland_initialize_thickness()

subroutine, public neverland_initialization::neverland_initialize_thickness ( real, dimension(szi_(g),szj_(g),szk_(gv)), intent(out)  h,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state,
real, intent(in)  P_ref 
)

This subroutine initializes layer thicknesses for the Neverland test case, by finding the depths of interfaces in a specified latitude-dependent temperature profile with an exponentially decaying thermocline on top of a linear stratification.

Parameters
[in]gThe ocean's grid structure.
[in]gvThe ocean's vertical grid structure.
[in]usA dimensional unit scaling type
[out]hThe thickness that is being initialized [H ~> m or kg m-2].
[in]param_fileA structure indicating the open file to parse for model parameter values.
eqn_of_stateinteger that selects the equation of state.
[in]p_refThe coordinate-density reference pressure [Pa].

Definition at line 114 of file Neverland_initialization.F90.

114  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
115  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
116  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
117  real, intent(out), dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h !< The thickness that is being
118  !! initialized [H ~> m or kg m-2].
119  type(param_file_type), intent(in) :: param_file !< A structure indicating the open
120  !! file to parse for model
121  !! parameter values.
122  type(EOS_type), pointer :: eqn_of_state !< integer that selects the
123  !! equation of state.
124  real, intent(in) :: P_Ref !< The coordinate-density
125  !! reference pressure [Pa].
126  ! Local variables
127  real :: e0(SZK_(G)+1) ! The resting interface heights, in depth units [Z ~> m],
128  ! usually negative because it is positive upward.
129  real, dimension(SZK_(G)) :: h_profile ! Vector of initial thickness profile [Z ~> m].
130  real :: e_interface ! Current interface position [Z ~> m].
131  real :: x,y,r1,r2 ! x,y and radial coordinates for computation of initial pert.
132  real :: pert_amp ! Amplitude of perturbations measured in Angstrom_H
133  real :: h_noise ! Amplitude of noise to scale h by
134  real :: noise ! Noise
135  type(randomNumberStream) :: rns ! Random numbers for stochastic tidal parameterization
136  character(len=40) :: mdl = "Neverland_initialize_thickness" ! This subroutine's name.
137  integer :: i, j, k, k1, is, ie, js, je, nz, itt
138 
139  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
140 
141  call mom_mesg(" Neverland_initialization.F90, Neverland_initialize_thickness: setting thickness", 5)
142  call get_param(param_file, mdl, "INIT_THICKNESS_PROFILE", h_profile, &
143  "Profile of initial layer thicknesses.", units="m", scale=us%m_to_Z, &
144  fail_if_missing=.true.)
145  call get_param(param_file, mdl, "NL_THICKNESS_PERT_AMP", pert_amp, &
146  "Amplitude of finite scale perturbations as fraction of depth.", &
147  units="nondim", default=0.)
148  call get_param(param_file, mdl, "NL_THICKNESS_NOISE_AMP", h_noise, &
149  "Amplitude of noise to scale layer by.", units="nondim", default=0.)
150 
151  ! e0 is the notional position of interfaces
152  e0(1) = 0. ! The surface
153  do k=1,nz
154  e0(k+1) = e0(k) - h_profile(k)
155  enddo
156 
157  do j=js,je ; do i=is,ie
158  e_interface = -g%bathyT(i,j)
159  do k=nz,2,-1
160  h(i,j,k) = gv%Z_to_H * (e0(k) - e_interface) ! Nominal thickness
161  x=(g%geoLonT(i,j)-g%west_lon)/g%len_lon
162  y=(g%geoLatT(i,j)-g%south_lat)/g%len_lat
163  r1=sqrt((x-0.7)**2+(y-0.2)**2)
164  r2=sqrt((x-0.3)**2+(y-0.25)**2)
165  h(i,j,k) = h(i,j,k) + pert_amp * (e0(k) - e0(nz+1)) * gv%Z_to_H * &
166  (spike(r1,0.15)-spike(r2,0.15)) ! Prescribed perturbation
167  if (h_noise /= 0.) then
168  rns = initializerandomnumberstream( int( 4096*(x + (y+1.)) ) )
169  call getrandomnumbers(rns, noise) ! x will be in (0,1)
170  noise = h_noise * 2. * ( noise - 0.5 ) ! range -h_noise to h_noise
171  h(i,j,k) = ( 1. + noise ) * h(i,j,k)
172  endif
173  h(i,j,k) = max( gv%Angstrom_H, h(i,j,k) ) ! Limit to non-negative
174  e_interface = e_interface + gv%H_to_Z * h(i,j,k) ! Actual position of upper interface
175  enddo
176  h(i,j,1) = gv%Z_to_H * (e0(1) - e_interface) ! Nominal thickness
177  h(i,j,1) = max( gv%Angstrom_H, h(i,j,1) ) ! Limit to non-negative
178  enddo ; enddo
179 

References mom_error_handler::mom_mesg(), and spike().

Referenced by mom_state_initialization::mom_initialize_state().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ neverland_initialize_topography()

subroutine, public neverland_initialization::neverland_initialize_topography ( real, dimension(g%isd:g%ied,g%jsd:g%jed), intent(out)  D,
type(dyn_horgrid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
real, intent(in)  max_depth 
)

This subroutine sets up the Neverland test case topography.

Parameters
[in]gThe dynamic horizontal grid type
[out]dOcean bottom depth in the units of depth_max
[in]param_fileParameter file structure
[in]max_depthMaximum ocean depth in arbitrary units

Definition at line 36 of file Neverland_initialization.F90.

36  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
37  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
38  intent(out) :: D !< Ocean bottom depth in the units of depth_max
39  type(param_file_type), intent(in) :: param_file !< Parameter file structure
40  real, intent(in) :: max_depth !< Maximum ocean depth in arbitrary units
41 
42  ! Local variables
43  real :: PI ! 3.1415926... calculated as 4*atan(1)
44  real :: D0 ! A constant to make the maximum !
45  ! basin depth MAXIMUM_DEPTH. !
46  real :: x, y
47  ! This include declares and sets the variable "version".
48 # include "version_variable.h"
49  character(len=40) :: mdl = "Neverland_initialize_topography" ! This subroutine's name.
50  integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
51  real :: nl_roughness_amp, nl_top_amp
52  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
53  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
54 
55  call mom_mesg(" Neverland_initialization.F90, Neverland_initialize_topography: setting topography", 5)
56 
57  call log_version(param_file, mdl, version, "")
58  call get_param(param_file, mdl, "NL_ROUGHNESS_AMP", nl_roughness_amp, &
59  "Amplitude of wavy signal in bathymetry.", default=0.05)
60  call get_param(param_file, mdl, "NL_CONTINENT_AMP", nl_top_amp, &
61  "Scale factor for topography - 0.0 for no continents.", default=1.0)
62 
63  pi = 4.0*atan(1.0)
64 
65 ! Calculate the depth of the bottom.
66  do j=js,je ; do i=is,ie
67  x = (g%geoLonT(i,j)-g%west_lon) / g%len_lon
68  y =( g%geoLatT(i,j)-g%south_lat) / g%len_lat
69 ! This sets topography that has a reentrant channel to the south.
70  d(i,j) = 1.0 - 1.1 * spike(y-1,0.12) - 1.1 * spike(y,0.12) - & !< The great northern wall and Antarctica
71  nl_top_amp*( &
72  (1.2 * spike(x,0.2) + 1.2 * spike(x-1.0,0.2)) * spike(min(0.0,y-0.3),0.2) & !< South America
73  + 1.2 * spike(x-0.5,0.2) * spike(min(0.0,y-0.55),0.2) & !< Africa
74  + 1.2 * (spike(x,0.12) + spike(x-1,0.12)) * spike(max(0.0,y-0.06),0.12) & !< Antarctic Peninsula
75  + 0.1 * (cosbell(x,0.1) + cosbell(x-1,0.1)) & !< Drake Passage ridge
76  + 0.5 * cosbell(x-0.16,0.05) * (cosbell(y-0.18,0.13)**0.4) & !< Scotia Arc East
77  + 0.4 * (cosbell(x-0.09,0.08)**0.4) * cosbell(y-0.26,0.05) & !< Scotia Arc North
78  + 0.4 * (cosbell(x-0.08,0.08)**0.4) * cosbell(y-0.1,0.05)) & !< Scotia Arc South
79  - nl_roughness_amp * cos(14*pi*x) * sin(14*pi*y) & !< roughness
80  - nl_roughness_amp * cos(20*pi*x) * cos(20*pi*y) !< roughness
81  if (d(i,j) < 0.0) d(i,j) = 0.0
82  d(i,j) = d(i,j) * max_depth
83  enddo ; enddo
84 

References cosbell(), mom_error_handler::mom_mesg(), and spike().

Referenced by mom_fixed_initialization::mom_initialize_topography().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ spike()

real function neverland_initialization::spike ( real, intent(in)  x,
real, intent(in)  L 
)
private

Returns the value of a sin-spike function evaluated at x/L.

Parameters
[in]xnon-dimensional position
[in]lnon-dimensional width

Definition at line 100 of file Neverland_initialization.F90.

100 
101  real , intent(in) :: x !< non-dimensional position
102  real , intent(in) :: L !< non-dimensional width
103  real :: PI !< 3.1415926... calculated as 4*atan(1)
104 
105  pi = 4.0*atan(1.0)
106  spike = (1 - sin(pi*min(abs(x/l),0.5)))

Referenced by neverland_initialize_thickness(), and neverland_initialize_topography().

Here is the caller graph for this function: