MOM6
user_initialization.F90
Go to the documentation of this file.
2 !***********************************************************************
3 !* GNU General Public License *
4 !* This file is a part of MOM. *
5 !* *
6 !* MOM is free software; you can redistribute it and/or modify it and *
7 !* are expected to follow the terms of the GNU General Public License *
8 !* as published by the Free Software Foundation; either version 2 of *
9 !* the License, or (at your option) any later version. *
10 !* *
11 !* MOM is distributed in the hope that it will be useful, but WITHOUT *
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY *
13 !* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public *
14 !* License for more details. *
15 !* *
16 !* For the full text of the GNU General Public License, *
17 !* write to: Free Software Foundation, Inc., *
18 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
19 !* or see: http://www.gnu.org/licenses/gpl.html *
20 !***********************************************************************
21 
25 use mom_get_input, only : directories
26 use mom_grid, only : ocean_grid_type
27 use mom_io, only : close_file, fieldtype, file_exists
28 use mom_io, only : open_file, read_data, read_axis_data, single_file
29 use mom_io, only : write_field, slasher
38 implicit none ; private
39 
40 #include <MOM_memory.h>
41 
45 
46 logical :: first_call = .true.
47 
48 contains
49 
50 !> Set vertical coordinates.
51 subroutine user_set_coord(Rlay, g_prime, GV, param_file, eqn_of_state)
52  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid
53  !! structure.
54  real, dimension(:), intent(out) :: Rlay !< Layer potential density.
55  real, dimension(:), intent(out) :: g_prime !< The reduced gravity at
56  !! each interface, in m s-2.
57  type(param_file_type), intent(in) :: param_file !< A structure indicating the
58  !! open file to parse for model
59  !! parameter values.
60  type(eos_type), pointer :: eqn_of_state !< Integer that selects the
61  !! equation of state.
62 
63  call mom_error(fatal, &
64  "USER_initialization.F90, USER_set_coord: " // &
65  "Unmodified user routine called - you must edit the routine to use it")
66  rlay(:) = 0.0
67  g_prime(:) = 0.0
68 
69  if (first_call) call write_user_log(param_file)
70 
71 end subroutine user_set_coord
72 
73 !> Initialize topography.
74 subroutine user_initialize_topography(D, G, param_file, max_depth)
75  type(dyn_horgrid_type), intent(in) :: G !< The dynamic horizontal grid type
76  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
77  intent(out) :: D !< Ocean bottom depth in m
78  type(param_file_type), intent(in) :: param_file !< Parameter file structure
79  real, intent(in) :: max_depth !< Maximum depth of model in m
80 
81  call mom_error(fatal, &
82  "USER_initialization.F90, USER_initialize_topography: " // &
83  "Unmodified user routine called - you must edit the routine to use it")
84 
85  d(:,:) = 0.0
86 
87  if (first_call) call write_user_log(param_file)
88 
89 end subroutine user_initialize_topography
90 
91 !> initialize thicknesses.
92 subroutine user_initialize_thickness(h, G, param_file, T, just_read_params)
93  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
94  real, intent(out), dimension(SZI_(G),SZJ_(G),SZK_(G)) :: h !< The thicknesses being
95  !! initialized.
96  type(param_file_type), intent(in) :: param_file !< A structure indicating the
97  !! open file to parse for model
98  !! parameter values.
99  real, intent(in), dimension(SZI_(G),SZJ_(G), SZK_(G)) :: T !< Potential temperature.
100  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
101  !! only read parameters without changing h.
102 
103  logical :: just_read ! If true, just read parameters but set nothing.
104 
105  call mom_error(fatal, &
106  "USER_initialization.F90, USER_initialize_thickness: " // &
107  "Unmodified user routine called - you must edit the routine to use it")
108 
109  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
110 
111  if (just_read) return ! All run-time parameters have been read, so return.
112 
113  h(:,:,1) = 0.0
114 
115  if (first_call) call write_user_log(param_file)
116 
117 end subroutine user_initialize_thickness
118 
119 !> initialize velocities.
120 subroutine user_initialize_velocity(u, v, G, param_file, just_read_params)
121  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
122  real, dimension(SZIB_(G), SZJ_(G), SZK_(G)), intent(out) :: u !< i-component of velocity [m/s]
123  real, dimension(SZI_(G), SZJB_(G), SZK_(G)), intent(out) :: v !< j-component of velocity [m/s]
124  type(param_file_type), intent(in) :: param_file !< A structure indicating the
125  !! open file to parse for model
126  !! parameter values.
127  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
128  !! only read parameters without changing h.
129 
130  logical :: just_read ! If true, just read parameters but set nothing.
131 
132  call mom_error(fatal, &
133  "USER_initialization.F90, USER_initialize_velocity: " // &
134  "Unmodified user routine called - you must edit the routine to use it")
135 
136  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
137 
138  if (just_read) return ! All run-time parameters have been read, so return.
139 
140  u(:,:,1) = 0.0
141  v(:,:,1) = 0.0
142 
143  if (first_call) call write_user_log(param_file)
144 
145 end subroutine user_initialize_velocity
146 
147 !> This function puts the initial layer temperatures and salinities
148 !! into T(:,:,:) and S(:,:,:).
149 subroutine user_init_temperature_salinity(T, S, G, param_file, eqn_of_state, just_read_params)
150  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
151  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: T !< Potential temperature (degC).
152  real, dimension(SZI_(G),SZJ_(G), SZK_(G)), intent(out) :: S !< Salinity (ppt).
153  type(param_file_type), intent(in) :: param_file !< A structure indicating the
154  !! open file to parse for model
155  !! parameter values.
156  type(eos_type), pointer :: eqn_of_state !< Integer that selects the
157  !! equation of state.
158  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
159  !! only read parameters without changing h.
160 
161  logical :: just_read ! If true, just read parameters but set nothing.
162 
163  call mom_error(fatal, &
164  "USER_initialization.F90, USER_init_temperature_salinity: " // &
165  "Unmodified user routine called - you must edit the routine to use it")
166 
167  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
168 
169  if (just_read) return ! All run-time parameters have been read, so return.
170 
171  t(:,:,1) = 0.0
172  s(:,:,1) = 0.0
173 
174  if (first_call) call write_user_log(param_file)
175 
176 end subroutine user_init_temperature_salinity
177 
178 !> Set up the sponges.
179 subroutine user_initialize_sponges(G, use_temperature, tv, param_file, CSp, h)
180  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure.
181  logical, intent(in) :: use_temperature !< Whether to use potential
182  !! temperature.
183  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers
184  !! to any available thermodynamic
185  !! fields, potential temperature and
186  !! salinity or mixed layer density.
187  !! Absent fields have NULL ptrs.
188  type(param_file_type), intent(in) :: param_file !< A structure indicating the
189  !! open file to parse for model
190  !! parameter values.
191  type(sponge_cs), pointer :: CSp !< A pointer to the sponge control
192  !! structure.
193  real, dimension(SZI_(G), SZJ_(G), SZK_(G)), intent(in) :: h !< Layer thicknesses.
194  call mom_error(fatal, &
195  "USER_initialization.F90, USER_initialize_sponges: " // &
196  "Unmodified user routine called - you must edit the routine to use it")
197 
198  if (first_call) call write_user_log(param_file)
199 
200 end subroutine user_initialize_sponges
201 
202 !> This subroutine sets the properties of flow at open boundary conditions.
203 subroutine user_set_obc_data(OBC, tv, G, param_file, tr_Reg)
204  type(ocean_obc_type), pointer :: OBC !< This open boundary condition type specifies
205  !! whether, where, and what open boundary
206  !! conditions are used.
207  type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers to any
208  !! available thermodynamic fields, including potential
209  !! temperature and salinity or mixed layer density. Absent
210  !! fields have NULL ptrs.
211  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
212  type(param_file_type), intent(in) :: param_file !< A structure indicating the
213  !! open file to parse for model
214  !! parameter values.
215  type(tracer_registry_type), pointer :: tr_Reg !< Tracer registry.
216 ! call MOM_error(FATAL, &
217 ! "USER_initialization.F90, USER_set_OBC_data: " // &
218 ! "Unmodified user routine called - you must edit the routine to use it")
219 
220  if (first_call) call write_user_log(param_file)
221 
222 end subroutine user_set_obc_data
223 
224 subroutine user_set_rotation(G, param_file)
225  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
226  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
227  call mom_error(fatal, &
228  "USER_initialization.F90, USER_set_rotation: " // &
229  "Unmodified user routine called - you must edit the routine to use it")
230 
231  if (first_call) call write_user_log(param_file)
232 
233 end subroutine user_set_rotation
234 
235 !> Write output about the parameter values being used.
236 subroutine write_user_log(param_file)
237  type(param_file_type), intent(in) :: param_file !< A structure indicating the
238  !! open file to parse for model
239  !! parameter values.
240 
241 ! This include declares and sets the variable "version".
242 #include "version_variable.h"
243  character(len=40) :: mdl = "user_initialization" ! This module's name.
244 
245  call log_version(param_file, mdl, version)
246  first_call = .false.
247 
248 end subroutine write_user_log
249 
250 !> \namespace user_initialization
251 !!
252 !! By Robert Hallberg, April 1994 - June 2002 *
253 !! *
254 !! This subroutine initializes the fields for the simulations. *
255 !! The one argument passed to initialize, Time, is set to the *
256 !! current time of the simulation. The fields which are initialized *
257 !! here are: *
258 !! u - Zonal velocity in m s-1. *
259 !! v - Meridional velocity in m s-1. *
260 !! h - Layer thickness in m. (Must be positive.) *
261 !! G%bathyT - Basin depth in m. (Must be positive.) *
262 !! G%CoriolisBu - The Coriolis parameter, in s-1. *
263 !! GV%g_prime - The reduced gravity at each interface, in m s-2. *
264 !! GV%Rlay - Layer potential density (coordinate variable), kg m-3. *
265 !! If ENABLE_THERMODYNAMICS is defined: *
266 !! T - Temperature in C. *
267 !! S - Salinity in psu. *
268 !! If BULKMIXEDLAYER is defined: *
269 !! Rml - Mixed layer and buffer layer potential densities in *
270 !! units of kg m-3. *
271 !! If SPONGE is defined: *
272 !! A series of subroutine calls are made to set up the damping *
273 !! rates and reference profiles for all variables that are damped *
274 !! in the sponge. *
275 !! Any user provided tracer code is also first linked through this *
276 !! subroutine. *
277 !! *
278 !! Forcing-related fields (taux, tauy, buoy, ustar, etc.) are set *
279 !! in MOM_surface_forcing.F90. *
280 !! *
281 !! These variables are all set in the set of subroutines (in this *
282 !! file) USER_initialize_bottom_depth, USER_initialize_thickness, *
283 !! USER_initialize_velocity, USER_initialize_temperature_salinity, *
284 !! USER_initialize_mixed_layer_density, USER_initialize_sponges, *
285 !! USER_set_coord, and USER_set_ref_profile. *
286 !! *
287 !! The names of these subroutines should be self-explanatory. They *
288 !! start with "USER_" to indicate that they will likely have to be *
289 !! modified for each simulation to set the initial conditions and *
290 !! boundary conditions. Most of these take two arguments: an integer *
291 !! argument specifying whether the fields are to be calculated *
292 !! internally or read from a NetCDF file; and a string giving the *
293 !! path to that file. If the field is initialized internally, the *
294 !! path is ignored. *
295 !! *
296 !! Macros written all in capital letters are defined in MOM_memory.h. *
297 !! *
298 !! A small fragment of the grid is shown below: *
299 !! *
300 !! j+1 x ^ x ^ x At x: q, CoriolisBu *
301 !! j+1 > o > o > At ^: v, tauy *
302 !! j x ^ x ^ x At >: u, taux *
303 !! j > o > o > At o: h, bathyT, buoy, tr, T, S, Rml, ustar *
304 !! j-1 x ^ x ^ x *
305 !! i-1 i i+1 At x & ^: *
306 !! i i+1 At > & o: *
307 !! *
308 !! The boundaries always run through q grid points (x). *
309 end module user_initialization
subroutine, public read_axis_data(filename, axis_name, var)
Definition: MOM_io.F90:424
subroutine, public user_initialize_sponges(G, use_temperature, tv, param_file, CSp, h)
Set up the sponges.
integer, parameter, public obc_direction_s
Indicates the boundary is an effective southern boundary.
integer, parameter, public obc_direction_w
Indicates the boundary is an effective western boundary.
integer, parameter, public obc_direction_n
Indicates the boundary is an effective northern boundary.
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:19
subroutine, public user_init_temperature_salinity(T, S, G, param_file, eqn_of_state, just_read_params)
This function puts the initial layer temperatures and salinities into T(:,:,:) and S(:...
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
This module contains I/O framework code.
Definition: MOM_io.F90:2
By Robert Hallberg, April 1994 - June 2002 *This subroutine initializes the fields for the simulation...
integer, parameter, public obc_simple
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
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
subroutine, public initialize_sponge(Iresttime, int_height, G, param_file, CS, Iresttime_i_mean, int_height_i_mean)
Definition: MOM_sponge.F90:142
subroutine write_user_log(param_file)
Write output about the parameter values being used.
subroutine, public set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
Definition: MOM_sponge.F90:271
Type to carry basic tracer information.
subroutine, public user_set_obc_data(OBC, tv, G, param_file, tr_Reg)
This subroutine sets the properties of flow at open boundary conditions.
logical function, public is_root_pe()
subroutine, public add_tracer_obc_values(name, Reg, OBC_inflow, OBC_in_u, OBC_in_v)
This subroutine adds open boundary condition concentrations for a tracer that has previously been reg...
subroutine, public mom_mesg(message, verb, all_print)
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
integer, parameter, public obc_direction_e
Indicates the boundary is an effective eastern boundary.
subroutine, public user_initialize_thickness(h, G, param_file, T, just_read_params)
initialize thicknesses.
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
subroutine, public user_set_coord(Rlay, g_prime, GV, param_file, eqn_of_state)
Set vertical coordinates.
integer, parameter, public obc_none
Controls where open boundary conditions are applied.
subroutine, public user_set_rotation(G, param_file)
subroutine, public mom_error(level, message, all_print)
subroutine, public user_initialize_topography(D, G, param_file, max_depth)
Initialize topography.
subroutine, public user_initialize_velocity(u, v, G, param_file, just_read_params)
initialize velocities.
A control structure for the equation of state.
Definition: MOM_EOS.F90:55