MOM6
MOM_state_initialization.F90
1 !> Initialization functions for state variables, u, v, h, T and S.
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_debugging, only : hchksum, qchksum, uvchksum
7 use mom_coms, only : max_across_pes, min_across_pes, reproducing_sum
8 use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
9 use mom_cpu_clock, only : clock_routine, clock_loop
10 use mom_domains, only : pass_var, pass_vector, sum_across_pes, broadcast
11 use mom_domains, only : root_pe, to_all, scalar_pair, cgrid_ne, agrid
12 use mom_eos, only : find_depth_of_pressure_in_cell
13 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
14 use mom_error_handler, only : calltree_enter, calltree_leave, calltree_waypoint
16 use mom_file_parser, only : log_version
17 use mom_get_input, only : directories
18 use mom_grid, only : ocean_grid_type, ispointincell
20 use mom_io, only : file_exists
22 use mom_io, only : slasher
23 use mom_open_boundary, only : ocean_obc_type, open_boundary_init
24 use mom_open_boundary, only : obc_none, obc_simple
25 use mom_open_boundary, only : open_boundary_query
26 use mom_open_boundary, only : set_tracer_data
27 use mom_open_boundary, only : open_boundary_test_extern_h
28 use mom_open_boundary, only : fill_temp_salt_segments
29 !use MOM_open_boundary, only : set_3D_OBC_data
30 use mom_grid_initialize, only : initialize_masks, set_grid_metrics
31 use mom_restart, only : restore_state, determine_is_new_run, mom_restart_cs
32 use mom_sponge, only : set_up_sponge_field, set_up_sponge_ml_density
33 use mom_sponge, only : initialize_sponge, sponge_cs
35 use mom_ale_sponge, only : ale_sponge_cs
36 use mom_string_functions, only : uppercase, lowercase
37 use mom_time_manager, only : time_type
41 use mom_verticalgrid, only : setverticalgridaxes, verticalgrid_type
42 use mom_ale, only : pressure_gradient_plm
44 use mom_eos, only : int_specific_vol_dp, convert_temp_salt_for_teos10
45 use user_initialization, only : user_initialize_thickness, user_initialize_velocity
46 use user_initialization, only : user_init_temperature_salinity
47 use user_initialization, only : user_set_obc_data
48 use user_initialization, only : user_initialize_sponges
49 use dome_initialization, only : dome_initialize_thickness
50 use dome_initialization, only : dome_set_obc_data
51 use dome_initialization, only : dome_initialize_sponges
52 use isomip_initialization, only : isomip_initialize_thickness
53 use isomip_initialization, only : isomip_initialize_sponges
54 use isomip_initialization, only : isomip_initialize_temperature_salinity
55 use rgc_initialization, only : rgc_initialize_sponges
56 use baroclinic_zone_initialization, only : baroclinic_zone_init_temperature_salinity
57 use benchmark_initialization, only : benchmark_initialize_thickness
58 use benchmark_initialization, only : benchmark_init_temperature_salinity
59 use neverland_initialization, only : neverland_initialize_thickness
60 use circle_obcs_initialization, only : circle_obcs_initialize_thickness
61 use lock_exchange_initialization, only : lock_exchange_initialize_thickness
62 use external_gwave_initialization, only : external_gwave_initialize_thickness
63 use dome2d_initialization, only : dome2d_initialize_thickness
64 use dome2d_initialization, only : dome2d_initialize_temperature_salinity
65 use dome2d_initialization, only : dome2d_initialize_sponges
66 use adjustment_initialization, only : adjustment_initialize_thickness
67 use adjustment_initialization, only : adjustment_initialize_temperature_salinity
68 use sloshing_initialization, only : sloshing_initialize_thickness
69 use sloshing_initialization, only : sloshing_initialize_temperature_salinity
70 use seamount_initialization, only : seamount_initialize_thickness
71 use seamount_initialization, only : seamount_initialize_temperature_salinity
72 use dumbbell_initialization, only : dumbbell_initialize_thickness
73 use dumbbell_initialization, only : dumbbell_initialize_temperature_salinity
74 use phillips_initialization, only : phillips_initialize_thickness
75 use phillips_initialization, only : phillips_initialize_velocity
76 use phillips_initialization, only : phillips_initialize_sponges
77 use rossby_front_2d_initialization, only : rossby_front_initialize_thickness
78 use rossby_front_2d_initialization, only : rossby_front_initialize_temperature_salinity
79 use rossby_front_2d_initialization, only : rossby_front_initialize_velocity
80 use scm_cvmix_tests, only: scm_cvmix_tests_ts_init
81 use dyed_channel_initialization, only : dyed_channel_set_obc_tracer_data
82 use dyed_obcs_initialization, only : dyed_obcs_set_obc_data
83 use supercritical_initialization, only : supercritical_set_obc_data
84 use soliton_initialization, only : soliton_initialize_velocity
85 use soliton_initialization, only : soliton_initialize_thickness
86 use bfb_initialization, only : bfb_initialize_sponges_southonly
87 use dense_water_initialization, only : dense_water_initialize_ts
88 use dense_water_initialization, only : dense_water_initialize_sponges
89 use dumbbell_initialization, only : dumbbell_initialize_sponges
90 
91 use midas_vertmap, only : find_interfaces, tracer_z_init
92 use midas_vertmap, only : determine_temperature
93 
94 use mom_ale, only : ale_initregridding, ale_cs, ale_initthicknesstocoord
95 use mom_ale, only : ale_remap_scalar, ale_build_grid, ale_regrid_accelerated
96 use mom_regridding, only : regridding_cs, set_regrid_params, getcoordinateresolution
97 use mom_regridding, only : regridding_main
98 use mom_remapping, only : remapping_cs, initialize_remapping
99 use mom_remapping, only : remapping_core_h
101 
102 use fms_io_mod, only : field_size
103 
104 implicit none ; private
105 
106 #include <MOM_memory.h>
107 
108 public mom_initialize_state
109 
110 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
111 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
112 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
113 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
114 
115 character(len=40) :: mdl = "MOM_state_initialization" !< This module's name.
116 
117 contains
118 
119 !> Initialize temporally evolving fields, either as initial
120 !! conditions or by reading them from a restart (or saves) file.
121 subroutine mom_initialize_state(u, v, h, tv, Time, G, GV, US, PF, dirs, &
122  restart_CS, ALE_CSp, tracer_Reg, sponge_CSp, &
123  ALE_sponge_CSp, OBC, Time_in)
124  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure.
125  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure.
126  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
127  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
128  intent(out) :: u !< The zonal velocity that is being
129  !! initialized [m s-1]
130  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
131  intent(out) :: v !< The meridional velocity that is being
132  !! initialized [m s-1]
133  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
134  intent(out) :: h !< Layer thicknesses [H ~> m or kg m-2]
135  type(thermo_var_ptrs), intent(inout) :: tv !< A structure pointing to various thermodynamic
136  !! variables
137  type(time_type), intent(inout) :: time !< Time at the start of the run segment.
138  type(param_file_type), intent(in) :: pf !< A structure indicating the open file to parse
139  !! for model parameter values.
140  type(directories), intent(in) :: dirs !< A structure containing several relevant
141  !! directory paths.
142  type(mom_restart_cs), pointer :: restart_cs !< A pointer to the restart control
143  !! structure.
144  type(ale_cs), pointer :: ale_csp !< The ALE control structure for remapping
145  type(tracer_registry_type), pointer :: tracer_reg !< A pointer to the tracer registry
146  type(sponge_cs), pointer :: sponge_csp !< The layerwise sponge control structure.
147  type(ale_sponge_cs), pointer :: ale_sponge_csp !< The ALE sponge control structure.
148  type(ocean_obc_type), pointer :: obc !< The open boundary condition control structure.
149  type(time_type), optional, intent(in) :: time_in !< Time at the start of the run segment.
150  !! Time_in overrides any value set for Time.
151  ! Local variables
152  character(len=200) :: filename ! The name of an input file.
153  character(len=200) :: filename2 ! The name of an input files.
154  character(len=200) :: inputdir ! The directory where NetCDF input files are.
155  character(len=200) :: config
156  real :: h_rescale ! A rescaling factor for thicknesses from the representation in
157  ! a restart file to the internal representation in this run.
158  real :: dt ! The baroclinic dynamics timestep for this run [s].
159  logical :: from_z_file, useale
160  logical :: new_sim
161  integer :: write_geom
162  logical :: use_temperature, use_sponge, use_obc
163  logical :: use_eos ! If true, density is calculated from T & S using an equation of state.
164  logical :: depress_sfc ! If true, remove the mass that would be displaced
165  ! by a large surface pressure by squeezing the column.
166  logical :: trim_ic_for_p_surf ! If true, remove the mass that would be displaced
167  ! by a large surface pressure, such as with an ice sheet.
168  logical :: regrid_accelerate
169  integer :: regrid_iterations
170 ! logical :: Analytic_FV_PGF, obsol_test
171  logical :: convert
172  logical :: just_read ! If true, only read the parameters because this
173  ! is a run from a restart file; this option
174  ! allows the use of Fatal unused parameters.
175  type(eos_type), pointer :: eos => null()
176  logical :: debug ! If true, write debugging output.
177  logical :: debug_obc ! If true, do debugging calls related to OBCs.
178  logical :: debug_layers = .false.
179  character(len=80) :: mesg
180 ! This include declares and sets the variable "version".
181 #include "version_variable.h"
182  integer :: i, j, k, is, ie, js, je, isq, ieq, jsq, jeq, nz
183  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
184 
185  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
186  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
187  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
188  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
189 
190  call calltree_enter("MOM_initialize_state(), MOM_state_initialization.F90")
191  call log_version(pf, mdl, version, "")
192  call get_param(pf, mdl, "DEBUG", debug, default=.false.)
193  call get_param(pf, mdl, "DEBUG_OBC", debug_obc, default=.false.)
194 
195  new_sim = determine_is_new_run(dirs%input_filename, dirs%restart_input_dir, &
196  g, restart_cs)
197  just_read = .not.new_sim
198 
199  call get_param(pf, mdl, "INPUTDIR", inputdir, &
200  "The directory in which input files are found.", default=".")
201  inputdir = slasher(inputdir)
202 
203  use_temperature = associated(tv%T)
204  useale = associated(ale_csp)
205  use_eos = associated(tv%eqn_of_state)
206  use_obc = associated(obc)
207  if (use_eos) eos => tv%eqn_of_state
208 
209  !====================================================================
210  ! Initialize temporally evolving fields, either as initial
211  ! conditions or by reading them from a restart (or saves) file.
212  !====================================================================
213 
214  if (new_sim) then
215  call mom_mesg("Run initialized internally.", 3)
216 
217  if (present(time_in)) time = time_in
218  ! Otherwise leave Time at its input value.
219 
220  ! This initialization should not be needed. Certainly restricting it
221  ! to the computational domain helps detect possible uninitialized
222  ! data in halos which should be covered by the pass_var(h) later.
223  !do k = 1, nz; do j = js, je; do i = is, ie
224  ! h(i,j,k) = 0.
225  !enddo
226  endif
227 
228  ! The remaining initialization calls are done, regardless of whether the
229  ! fields are actually initialized here (if just_read=.false.) or whether it
230  ! is just to make sure that all valid parameters are read to enable the
231  ! detection of unused parameters.
232  call get_param(pf, mdl, "INIT_LAYERS_FROM_Z_FILE", from_z_file, &
233  "If true, initialize the layer thicknesses, temperatures, "//&
234  "and salinities from a Z-space file on a latitude-longitude "//&
235  "grid.", default=.false., do_not_log=just_read)
236 
237  if (from_z_file) then
238  ! Initialize thickness and T/S from z-coordinate data in a file.
239  if (.NOT.use_temperature) call mom_error(fatal,"MOM_initialize_state : "//&
240  "use_temperature must be true if INIT_LAYERS_FROM_Z_FILE is true")
241 
242  call mom_temp_salt_initialize_from_z(h, tv, g, gv, us, pf, just_read_params=just_read)
243 
244  else
245  ! Initialize thickness, h.
246  call get_param(pf, mdl, "THICKNESS_CONFIG", config, &
247  "A string that determines how the initial layer "//&
248  "thicknesses are specified for a new run: \n"//&
249  " \t file - read interface heights from the file specified \n"//&
250  " \t thickness_file - read thicknesses from the file specified \n"//&
251  " \t\t by (THICKNESS_FILE).\n"//&
252  " \t coord - determined by ALE coordinate.\n"//&
253  " \t uniform - uniform thickness layers evenly distributed \n"//&
254  " \t\t between the surface and MAXIMUM_DEPTH. \n"//&
255  " \t list - read a list of positive interface depths. \n"//&
256  " \t DOME - use a slope and channel configuration for the \n"//&
257  " \t\t DOME sill-overflow test case. \n"//&
258  " \t ISOMIP - use a configuration for the \n"//&
259  " \t\t ISOMIP test case. \n"//&
260  " \t benchmark - use the benchmark test case thicknesses. \n"//&
261  " \t Neverland - use the Neverland test case thicknesses. \n"//&
262  " \t search - search a density profile for the interface \n"//&
263  " \t\t densities. This is not yet implemented. \n"//&
264  " \t circle_obcs - the circle_obcs test case is used. \n"//&
265  " \t DOME2D - 2D version of DOME initialization. \n"//&
266  " \t adjustment2d - 2D lock exchange thickness ICs. \n"//&
267  " \t sloshing - sloshing gravity thickness ICs. \n"//&
268  " \t seamount - no motion test with seamount ICs. \n"//&
269  " \t dumbbell - sloshing channel ICs. \n"//&
270  " \t soliton - Equatorial Rossby soliton. \n"//&
271  " \t rossby_front - a mixed layer front in thermal wind balance.\n"//&
272  " \t USER - call a user modified routine.", &
273  fail_if_missing=new_sim, do_not_log=just_read)
274  select case (trim(config))
275  case ("file")
276  call initialize_thickness_from_file(h, g, gv, us, pf, .false., just_read_params=just_read)
277  case ("thickness_file")
278  call initialize_thickness_from_file(h, g, gv, us, pf, .true., just_read_params=just_read)
279  case ("coord")
280  if (new_sim .and. useale) then
281  call ale_initthicknesstocoord( ale_csp, g, gv, h )
282  elseif (new_sim) then
283  call mom_error(fatal, "MOM_initialize_state: USE_REGRIDDING must be True "//&
284  "for THICKNESS_CONFIG of 'coord'")
285  endif
286  case ("uniform"); call initialize_thickness_uniform(h, g, gv, pf, &
287  just_read_params=just_read)
288  case ("list"); call initialize_thickness_list(h, g, gv, us, pf, &
289  just_read_params=just_read)
290  case ("DOME"); call dome_initialize_thickness(h, g, gv, pf, &
291  just_read_params=just_read)
292  case ("ISOMIP"); call isomip_initialize_thickness(h, g, gv, us, pf, tv, &
293  just_read_params=just_read)
294  case ("benchmark"); call benchmark_initialize_thickness(h, g, gv, us, pf, &
295  tv%eqn_of_state, tv%P_Ref, just_read_params=just_read)
296  case ("Neverland"); call neverland_initialize_thickness(h, g, gv, us, pf, &
297  tv%eqn_of_state, tv%P_Ref)
298  case ("search"); call initialize_thickness_search
299  case ("circle_obcs"); call circle_obcs_initialize_thickness(h, g, gv, pf, &
300  just_read_params=just_read)
301  case ("lock_exchange"); call lock_exchange_initialize_thickness(h, g, gv, us, &
302  pf, just_read_params=just_read)
303  case ("external_gwave"); call external_gwave_initialize_thickness(h, g, gv, us, &
304  pf, just_read_params=just_read)
305  case ("DOME2D"); call dome2d_initialize_thickness(h, g, gv, us, pf, &
306  just_read_params=just_read)
307  case ("adjustment2d"); call adjustment_initialize_thickness(h, g, gv, us, &
308  pf, just_read_params=just_read)
309  case ("sloshing"); call sloshing_initialize_thickness(h, g, gv, us, pf, &
310  just_read_params=just_read)
311  case ("seamount"); call seamount_initialize_thickness(h, g, gv, us, pf, &
312  just_read_params=just_read)
313  case ("dumbbell"); call dumbbell_initialize_thickness(h, g, gv, us, pf, &
314  just_read_params=just_read)
315  case ("soliton"); call soliton_initialize_thickness(h, g, gv, us)
316  case ("phillips"); call phillips_initialize_thickness(h, g, gv, us, pf, &
317  just_read_params=just_read)
318  case ("rossby_front"); call rossby_front_initialize_thickness(h, g, gv, &
319  pf, just_read_params=just_read)
320  case ("USER"); call user_initialize_thickness(h, g, gv, pf, &
321  just_read_params=just_read)
322  case default ; call mom_error(fatal, "MOM_initialize_state: "//&
323  "Unrecognized layer thickness configuration "//trim(config))
324  end select
325 
326  ! Initialize temperature and salinity (T and S).
327  if ( use_temperature ) then
328  call get_param(pf, mdl, "TS_CONFIG", config, &
329  "A string that determines how the initial tempertures "//&
330  "and salinities are specified for a new run: \n"//&
331  " \t file - read velocities from the file specified \n"//&
332  " \t\t by (TS_FILE). \n"//&
333  " \t fit - find the temperatures that are consistent with \n"//&
334  " \t\t the layer densities and salinity S_REF. \n"//&
335  " \t TS_profile - use temperature and salinity profiles \n"//&
336  " \t\t (read from TS_FILE) to set layer densities. \n"//&
337  " \t benchmark - use the benchmark test case T & S. \n"//&
338  " \t linear - linear in logical layer space. \n"//&
339  " \t DOME2D - 2D DOME initialization. \n"//&
340  " \t ISOMIP - ISOMIP initialization. \n"//&
341  " \t adjustment2d - 2d lock exchange T/S ICs. \n"//&
342  " \t sloshing - sloshing mode T/S ICs. \n"//&
343  " \t seamount - no motion test with seamount ICs. \n"//&
344  " \t dumbbell - sloshing channel ICs. \n"//&
345  " \t rossby_front - a mixed layer front in thermal wind balance.\n"//&
346  " \t SCM_CVMix_tests - used in the SCM CVMix tests.\n"//&
347  " \t USER - call a user modified routine.", &
348  fail_if_missing=new_sim, do_not_log=just_read)
349 ! " \t baroclinic_zone - an analytic baroclinic zone. \n"//&
350  select case (trim(config))
351  case ("fit"); call initialize_temp_salt_fit(tv%T, tv%S, g, gv, pf, &
352  eos, tv%P_Ref, just_read_params=just_read)
353  case ("file"); call initialize_temp_salt_from_file(tv%T, tv%S, g, &
354  pf, just_read_params=just_read)
355  case ("benchmark"); call benchmark_init_temperature_salinity(tv%T, tv%S, &
356  g, gv, pf, eos, tv%P_Ref, just_read_params=just_read)
357  case ("TS_profile") ; call initialize_temp_salt_from_profile(tv%T, tv%S, &
358  g, pf, just_read_params=just_read)
359  case ("linear"); call initialize_temp_salt_linear(tv%T, tv%S, g, pf, &
360  just_read_params=just_read)
361  case ("DOME2D"); call dome2d_initialize_temperature_salinity ( tv%T, &
362  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
363  case ("ISOMIP"); call isomip_initialize_temperature_salinity ( tv%T, &
364  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
365  case ("adjustment2d"); call adjustment_initialize_temperature_salinity ( tv%T, &
366  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
367  case ("baroclinic_zone"); call baroclinic_zone_init_temperature_salinity( tv%T, &
368  tv%S, h, g, gv, us, pf, just_read_params=just_read)
369  case ("sloshing"); call sloshing_initialize_temperature_salinity(tv%T, &
370  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
371  case ("seamount"); call seamount_initialize_temperature_salinity(tv%T, &
372  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
373  case ("dumbbell"); call dumbbell_initialize_temperature_salinity(tv%T, &
374  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
375  case ("rossby_front"); call rossby_front_initialize_temperature_salinity ( tv%T, &
376  tv%S, h, g, gv, pf, eos, just_read_params=just_read)
377  case ("SCM_CVMix_tests"); call scm_cvmix_tests_ts_init(tv%T, tv%S, h, &
378  g, gv, us, pf, just_read_params=just_read)
379  case ("dense"); call dense_water_initialize_ts(g, gv, pf, eos, tv%T, tv%S, &
380  h, just_read_params=just_read)
381  case ("USER"); call user_init_temperature_salinity(tv%T, tv%S, g, pf, eos, &
382  just_read_params=just_read)
383  case default ; call mom_error(fatal, "MOM_initialize_state: "//&
384  "Unrecognized Temp & salt configuration "//trim(config))
385  end select
386  endif
387  endif ! not from_Z_file.
388  if (use_temperature .and. use_obc) &
389  call fill_temp_salt_segments(g, obc, tv)
390 
391  ! The thicknesses in halo points might be needed to initialize the velocities.
392  if (new_sim) call pass_var(h, g%Domain)
393 
394  ! Initialize velocity components, u and v
395  call get_param(pf, mdl, "VELOCITY_CONFIG", config, &
396  "A string that determines how the initial velocities "//&
397  "are specified for a new run: \n"//&
398  " \t file - read velocities from the file specified \n"//&
399  " \t\t by (VELOCITY_FILE). \n"//&
400  " \t zero - the fluid is initially at rest. \n"//&
401  " \t uniform - the flow is uniform (determined by\n"//&
402  " \t\t parameters INITIAL_U_CONST and INITIAL_V_CONST).\n"//&
403  " \t rossby_front - a mixed layer front in thermal wind balance.\n"//&
404  " \t soliton - Equatorial Rossby soliton.\n"//&
405  " \t USER - call a user modified routine.", default="zero", &
406  do_not_log=just_read)
407  select case (trim(config))
408  case ("file"); call initialize_velocity_from_file(u, v, g, pf, &
409  just_read_params=just_read)
410  case ("zero"); call initialize_velocity_zero(u, v, g, pf, &
411  just_read_params=just_read)
412  case ("uniform"); call initialize_velocity_uniform(u, v, g, pf, &
413  just_read_params=just_read)
414  case ("circular"); call initialize_velocity_circular(u, v, g, pf, &
415  just_read_params=just_read)
416  case ("phillips"); call phillips_initialize_velocity(u, v, g, gv, us, pf, &
417  just_read_params=just_read)
418  case ("rossby_front"); call rossby_front_initialize_velocity(u, v, h, &
419  g, gv, us, pf, just_read_params=just_read)
420  case ("soliton"); call soliton_initialize_velocity(u, v, h, g)
421  case ("USER"); call user_initialize_velocity(u, v, g, pf, &
422  just_read_params=just_read)
423  case default ; call mom_error(fatal, "MOM_initialize_state: "//&
424  "Unrecognized velocity configuration "//trim(config))
425  end select
426 
427  if (new_sim) call pass_vector(u, v, g%Domain)
428  if (debug .and. new_sim) then
429  call uvchksum("MOM_initialize_state [uv]", u, v, g%HI, haloshift=1)
430  endif
431 
432  ! Optionally convert the thicknesses from m to kg m-2. This is particularly
433  ! useful in a non-Boussinesq model.
434  call get_param(pf, mdl, "CONVERT_THICKNESS_UNITS", convert, &
435  "If true, convert the thickness initial conditions from "//&
436  "units of m to kg m-2 or vice versa, depending on whether "//&
437  "BOUSSINESQ is defined. This does not apply if a restart "//&
438  "file is read.", default=.not.gv%Boussinesq, do_not_log=just_read)
439 
440  if (new_sim .and. convert .and. .not.gv%Boussinesq) &
441  ! Convert thicknesses from geomtric distances to mass-per-unit-area.
442  call convert_thickness(h, g, gv, us, tv)
443 
444  ! Remove the mass that would be displaced by an ice shelf or inverse barometer.
445  call get_param(pf, mdl, "DEPRESS_INITIAL_SURFACE", depress_sfc, &
446  "If true, depress the initial surface to avoid huge "//&
447  "tsunamis when a large surface pressure is applied.", &
448  default=.false., do_not_log=just_read)
449  call get_param(pf, mdl, "TRIM_IC_FOR_P_SURF", trim_ic_for_p_surf, &
450  "If true, cuts way the top of the column for initial conditions "//&
451  "at the depth where the hydrostatic pressure matches the imposed "//&
452  "surface pressure which is read from file.", default=.false., &
453  do_not_log=just_read)
454  if (depress_sfc .and. trim_ic_for_p_surf) call mom_error(fatal, "MOM_initialize_state: "//&
455  "DEPRESS_INITIAL_SURFACE and TRIM_IC_FOR_P_SURF are exclusive and cannot both be True")
456  if (new_sim .and. debug .and. (depress_sfc .or. trim_ic_for_p_surf)) &
457  call hchksum(h, "Pre-depress: h ", g%HI, haloshift=1, scale=gv%H_to_m)
458  if (depress_sfc) call depress_surface(h, g, gv, us, pf, tv, just_read_params=just_read)
459  if (trim_ic_for_p_surf) call trim_for_ice(pf, g, gv, us, ale_csp, tv, h, just_read_params=just_read)
460 
461  ! Perhaps we want to run the regridding coordinate generator for multiple
462  ! iterations here so the initial grid is consistent with the coordinate
463  if (useale) then
464  call get_param(pf, mdl, "REGRID_ACCELERATE_INIT", regrid_accelerate, &
465  "If true, runs REGRID_ACCELERATE_ITERATIONS iterations of the regridding "//&
466  "algorithm to push the initial grid to be consistent with the initial "//&
467  "condition. Useful only for state-based and iterative coordinates.", &
468  default=.false., do_not_log=just_read)
469  if (regrid_accelerate) then
470  call get_param(pf, mdl, "REGRID_ACCELERATE_ITERATIONS", regrid_iterations, &
471  "The number of regridding iterations to perform to generate "//&
472  "an initial grid that is consistent with the initial conditions.", &
473  default=1, do_not_log=just_read)
474 
475  call get_param(pf, mdl, "DT", dt, "Timestep", fail_if_missing=.true.)
476 
477  if (new_sim .and. debug) &
478  call hchksum(h, "Pre-ALE_regrid: h ", g%HI, haloshift=1, scale=gv%H_to_m)
479  call ale_regrid_accelerated(ale_csp, g, gv, h, tv, regrid_iterations, u, v, tracer_reg, &
480  dt=dt, initial=.true.)
481  endif
482  endif
483  ! This is the end of the block of code that might have initialized fields
484  ! internally at the start of a new run.
485 
486  if (.not.new_sim) then ! This block restores the state from a restart file.
487  ! This line calls a subroutine that reads the initial conditions
488  ! from a previously generated file.
489  call restore_state(dirs%input_filename, dirs%restart_input_dir, time, &
490  g, restart_cs)
491  if (present(time_in)) time = time_in
492  if ((gv%m_to_H_restart /= 0.0) .and. (gv%m_to_H_restart /= gv%m_to_H)) then
493  h_rescale = gv%m_to_H / gv%m_to_H_restart
494  do k=1,nz ; do j=js,je ; do i=is,ie ; h(i,j,k) = h_rescale * h(i,j,k) ; enddo ; enddo ; enddo
495  endif
496  endif
497 
498  if ( use_temperature ) then
499  call pass_var(tv%T, g%Domain, complete=.false.)
500  call pass_var(tv%S, g%Domain, complete=.false.)
501  endif
502  call pass_var(h, g%Domain)
503 
504  if (debug) then
505  call hchksum(h, "MOM_initialize_state: h ", g%HI, haloshift=1, scale=gv%H_to_m)
506  if ( use_temperature ) call hchksum(tv%T, "MOM_initialize_state: T ", g%HI, haloshift=1)
507  if ( use_temperature ) call hchksum(tv%S, "MOM_initialize_state: S ", g%HI, haloshift=1)
508  if ( use_temperature .and. debug_layers) then ; do k=1,nz
509  write(mesg,'("MOM_IS: T[",I2,"]")') k
510  call hchksum(tv%T(:,:,k), mesg, g%HI, haloshift=1)
511  write(mesg,'("MOM_IS: S[",I2,"]")') k
512  call hchksum(tv%S(:,:,k), mesg, g%HI, haloshift=1)
513  enddo ; endif
514  endif
515 
516  call get_param(pf, mdl, "SPONGE", use_sponge, &
517  "If true, sponges may be applied anywhere in the domain. "//&
518  "The exact location and properties of those sponges are "//&
519  "specified via SPONGE_CONFIG.", default=.false.)
520  if ( use_sponge ) then
521  call get_param(pf, mdl, "SPONGE_CONFIG", config, &
522  "A string that sets how the sponges are configured: \n"//&
523  " \t file - read sponge properties from the file \n"//&
524  " \t\t specified by (SPONGE_FILE).\n"//&
525  " \t ISOMIP - apply ale sponge in the ISOMIP case \n"//&
526  " \t RGC - apply sponge in the rotating_gravity_current case \n"//&
527  " \t DOME - use a slope and channel configuration for the \n"//&
528  " \t\t DOME sill-overflow test case. \n"//&
529  " \t BFB - Sponge at the southern boundary of the domain\n"//&
530  " \t\t for buoyancy-forced basin case.\n"//&
531  " \t USER - call a user modified routine.", default="file")
532  select case (trim(config))
533  case ("DOME"); call dome_initialize_sponges(g, gv, us, tv, pf, sponge_csp)
534  case ("DOME2D"); call dome2d_initialize_sponges(g, gv, tv, pf, useale, &
535  sponge_csp, ale_sponge_csp)
536  case ("ISOMIP"); call isomip_initialize_sponges(g, gv, us, tv, pf, useale, &
537  sponge_csp, ale_sponge_csp)
538  case("RGC"); call rgc_initialize_sponges(g, gv, tv, u, v, pf, useale, &
539  sponge_csp, ale_sponge_csp)
540  case ("USER"); call user_initialize_sponges(g, gv, use_temperature, tv, pf, sponge_csp, h)
541  case ("BFB"); call bfb_initialize_sponges_southonly(g, gv, us, use_temperature, tv, pf, &
542  sponge_csp, h)
543  case ("DUMBBELL"); call dumbbell_initialize_sponges(g, gv, us, tv, pf, useale, &
544  sponge_csp, ale_sponge_csp)
545  case ("phillips"); call phillips_initialize_sponges(g, gv, us, tv, pf, sponge_csp, h)
546  case ("dense"); call dense_water_initialize_sponges(g, gv, tv, pf, useale, &
547  sponge_csp, ale_sponge_csp)
548  case ("file"); call initialize_sponges_file(g, gv, us, use_temperature, tv, pf, &
549  sponge_csp, ale_sponge_csp, time)
550  case default ; call mom_error(fatal, "MOM_initialize_state: "//&
551  "Unrecognized sponge configuration "//trim(config))
552  end select
553  endif
554 
555  ! Reads OBC parameters not pertaining to the location of the boundaries
556  call open_boundary_init(g, pf, obc)
557 
558  ! This controls user code for setting open boundary data
559  if (associated(obc)) then
560  call get_param(pf, mdl, "OBC_USER_CONFIG", config, &
561  "A string that sets how the user code is invoked to set open boundary data: \n"//&
562  " DOME - specified inflow on northern boundary\n"//&
563  " dyed_channel - supercritical with dye on the inflow boundary\n"//&
564  " dyed_obcs - circle_obcs with dyes on the open boundaries\n"//&
565  " Kelvin - barotropic Kelvin wave forcing on the western boundary\n"//&
566  " shelfwave - Flather with shelf wave forcing on western boundary\n"//&
567  " supercritical - now only needed here for the allocations\n"//&
568  " tidal_bay - Flather with tidal forcing on eastern boundary\n"//&
569  " USER - user specified", default="none")
570  if (trim(config) == "DOME") then
571  call dome_set_obc_data(obc, tv, g, gv, us, pf, tracer_reg)
572  elseif (trim(config) == "dyed_channel") then
573  call dyed_channel_set_obc_tracer_data(obc, g, gv, pf, tracer_reg)
574  obc%update_OBC = .true.
575  elseif (trim(config) == "dyed_obcs") then
576  call dyed_obcs_set_obc_data(obc, g, gv, pf, tracer_reg)
577  elseif (trim(config) == "Kelvin") then
578  obc%update_OBC = .true.
579  elseif (trim(config) == "shelfwave") then
580  obc%update_OBC = .true.
581  elseif (lowercase(trim(config)) == "supercritical") then
582  call supercritical_set_obc_data(obc, g, pf)
583  elseif (trim(config) == "tidal_bay") then
584  obc%update_OBC = .true.
585  elseif (trim(config) == "USER") then
586  call user_set_obc_data(obc, tv, g, pf, tracer_reg)
587  elseif (.not. trim(config) == "none") then
588  call mom_error(fatal, "The open boundary conditions specified by "//&
589  "OBC_USER_CONFIG = "//trim(config)//" have not been fully implemented.")
590  endif
591  if (open_boundary_query(obc, apply_open_obc=.true.)) then
592  call set_tracer_data(obc, tv, h, g, pf, tracer_reg)
593  endif
594  endif
595 ! if (open_boundary_query(OBC, apply_nudged_OBC=.true.)) then
596 ! call set_3D_OBC_data(OBC, tv, h, G, PF, tracer_Reg)
597 ! endif
598  ! Still need a way to specify the boundary values
599  if (debug.and.associated(obc)) then
600  call hchksum(g%mask2dT, 'MOM_initialize_state: mask2dT ', g%HI)
601  call uvchksum('MOM_initialize_state: mask2dC[uv]', g%mask2dCu, &
602  g%mask2dCv, g%HI)
603  call qchksum(g%mask2dBu, 'MOM_initialize_state: mask2dBu ', g%HI)
604  endif
605 
606  if (debug_obc) call open_boundary_test_extern_h(g, obc, h)
607  call calltree_leave('MOM_initialize_state()')
608 
609 end subroutine mom_initialize_state
610 
611 !> Reads the layer thicknesses or interface heights from a file.
612 subroutine initialize_thickness_from_file(h, G, GV, US, param_file, file_has_thickness, &
613  just_read_params)
614  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
615  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
616  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
617  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
618  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
619  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
620  !! to parse for model parameter values.
621  logical, intent(in) :: file_has_thickness !< If true, this file contains layer
622  !! thicknesses; otherwise it contains
623  !! interface heights.
624  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
625  !! only read parameters without changing h.
626 
627  ! Local variables
628  real :: eta(SZI_(G),SZJ_(G),SZK_(G)+1) ! Interface heights, in depth units.
629  integer :: inconsistent = 0
630  logical :: correct_thickness
631  logical :: just_read ! If true, just read parameters but set nothing.
632  character(len=40) :: mdl = "initialize_thickness_from_file" ! This subroutine's name.
633  character(len=200) :: filename, thickness_file, inputdir, mesg ! Strings for file/path
634  integer :: i, j, k, is, ie, js, je, nz
635 
636  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
637 
638  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
639 
640  if (.not.just_read) &
641  call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
642 
643  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".", do_not_log=just_read)
644  inputdir = slasher(inputdir)
645  call get_param(param_file, mdl, "THICKNESS_FILE", thickness_file, &
646  "The name of the thickness file.", &
647  fail_if_missing=.not.just_read, do_not_log=just_read)
648 
649  filename = trim(inputdir)//trim(thickness_file)
650  if (.not.just_read) call log_param(param_file, mdl, "INPUTDIR/THICKNESS_FILE", filename)
651 
652  if ((.not.just_read) .and. (.not.file_exists(filename, g%Domain))) call mom_error(fatal, &
653  " initialize_thickness_from_file: Unable to open "//trim(filename))
654 
655  if (file_has_thickness) then
656  !### Consider adding a parameter to use to rescale h.
657  if (just_read) return ! All run-time parameters have been read, so return.
658  call mom_read_data(filename, "h", h(:,:,:), g%Domain, scale=gv%m_to_H)
659  else
660  call get_param(param_file, mdl, "ADJUST_THICKNESS", correct_thickness, &
661  "If true, all mass below the bottom removed if the "//&
662  "topography is shallower than the thickness input file "//&
663  "would indicate.", default=.false., do_not_log=just_read)
664  if (just_read) return ! All run-time parameters have been read, so return.
665 
666  call mom_read_data(filename, "eta", eta(:,:,:), g%Domain, scale=us%m_to_Z)
667 
668  if (correct_thickness) then
669  call adjustetatofitbathymetry(g, gv, us, eta, h)
670  else
671  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
672  if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z)) then
673  eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
674  h(i,j,k) = gv%Angstrom_H
675  else
676  h(i,j,k) = gv%Z_to_H * (eta(i,j,k) - eta(i,j,k+1))
677  endif
678  enddo ; enddo ; enddo
679 
680  do j=js,je ; do i=is,ie
681  if (abs(eta(i,j,nz+1) + g%bathyT(i,j)) > 1.0*us%m_to_Z) &
682  inconsistent = inconsistent + 1
683  enddo ; enddo
684  call sum_across_pes(inconsistent)
685 
686  if ((inconsistent > 0) .and. (is_root_pe())) then
687  write(mesg,'("Thickness initial conditions are inconsistent ",'// &
688  '"with topography in ",I8," places.")') inconsistent
689  call mom_error(warning, mesg)
690  endif
691  endif
692 
693  endif
694  call calltree_leave(trim(mdl)//'()')
695 end subroutine initialize_thickness_from_file
696 
697 !> Adjust interface heights to fit the bathymetry and diagnose layer thickness.
698 !!
699 !! If the bottom most interface is below the topography then the bottom-most
700 !! layers are contracted to GV%Angstrom_m.
701 !! If the bottom most interface is above the topography then the entire column
702 !! is dilated (expanded) to fill the void.
703 !! @remark{There is a (hard-wired) "tolerance" parameter such that the
704 !! criteria for adjustment must equal or exceed 10cm.}
705 subroutine adjustetatofitbathymetry(G, GV, US, eta, h)
706  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
707  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
708  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
709  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), intent(inout) :: eta !< Interface heights [Z ~> m].
710  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
711  ! Local variables
712  integer :: i, j, k, is, ie, js, je, nz, contractions, dilations
713  real :: hTolerance = 0.1 !< Tolerance to exceed adjustment criteria [Z ~> m]
714  real :: hTmp, eTmp, dilate
715  character(len=100) :: mesg
716 
717  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
718  htolerance = 0.1*us%m_to_Z
719 
720  contractions = 0
721  do j=js,je ; do i=is,ie
722  if (-eta(i,j,nz+1) > g%bathyT(i,j) + htolerance) then
723  eta(i,j,nz+1) = -g%bathyT(i,j)
724  contractions = contractions + 1
725  endif
726  enddo ; enddo
727  call sum_across_pes(contractions)
728  if ((contractions > 0) .and. (is_root_pe())) then
729  write(mesg,'("Thickness initial conditions were contracted ",'// &
730  '"to fit topography in ",I8," places.")') contractions
731  call mom_error(warning, 'adjustEtaToFitBathymetry: '//mesg)
732  endif
733 
734  ! To preserve previous answers in non-Boussinesq cases, delay converting
735  ! thicknesses to units of H until the end of this routine.
736  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
737  ! Collapse layers to thinnest possible if the thickness less than
738  ! the thinnest possible (or negative).
739  if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z)) then
740  eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
741  h(i,j,k) = gv%Angstrom_Z
742  else
743  h(i,j,k) = (eta(i,j,k) - eta(i,j,k+1))
744  endif
745  enddo ; enddo ; enddo
746 
747  dilations = 0
748  do j=js,je ; do i=is,ie
749  ! The whole column is dilated to accommodate deeper topography than
750  ! the bathymetry would indicate.
751  ! This should be... if ((G%mask2dt(i,j)*(eta(i,j,1)-eta(i,j,nz+1)) > 0.0) .and. &
752  if (-eta(i,j,nz+1) < g%bathyT(i,j) - htolerance) then
753  dilations = dilations + 1
754  if (eta(i,j,1) <= eta(i,j,nz+1)) then
755  do k=1,nz ; h(i,j,k) = (eta(i,j,1) + g%bathyT(i,j)) / real(nz) ; enddo
756  else
757  dilate = (eta(i,j,1) + g%bathyT(i,j)) / (eta(i,j,1) - eta(i,j,nz+1))
758  do k=1,nz ; h(i,j,k) = h(i,j,k) * dilate ; enddo
759  endif
760  do k=nz,2,-1 ; eta(i,j,k) = eta(i,j,k+1) + h(i,j,k) ; enddo
761  endif
762  enddo ; enddo
763 
764  ! Now convert thicknesses to units of H.
765  do k=1,nz ; do j=js,je ; do i=is,ie
766  h(i,j,k) = h(i,j,k)*gv%Z_to_H
767  enddo ; enddo ; enddo
768 
769  call sum_across_pes(dilations)
770  if ((dilations > 0) .and. (is_root_pe())) then
771  write(mesg,'("Thickness initial conditions were dilated ",'// &
772  '"to fit topography in ",I8," places.")') dilations
773  call mom_error(warning, 'adjustEtaToFitBathymetry: '//mesg)
774  endif
775 
776 end subroutine adjustetatofitbathymetry
777 
778 !> Initializes thickness to be uniform
779 subroutine initialize_thickness_uniform(h, G, GV, param_file, just_read_params)
780  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
781  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
782  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
783  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
784  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
785  !! to parse for model parameter values.
786  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
787  !! only read parameters without changing h.
788  ! Local variables
789  character(len=40) :: mdl = "initialize_thickness_uniform" ! This subroutine's name.
790  real :: e0(SZK_(G)+1) ! The resting interface heights, in depth units, usually
791  ! negative because it is positive upward.
792  real :: eta1D(SZK_(G)+1)! Interface height relative to the sea surface
793  ! positive upward, in depth units.
794  logical :: just_read ! If true, just read parameters but set nothing.
795  integer :: i, j, k, is, ie, js, je, nz
796 
797  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
798 
799  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
800 
801  if (just_read) return ! This subroutine has no run-time parameters.
802 
803  call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
804 
805  if (g%max_depth<=0.) call mom_error(fatal,"initialize_thickness_uniform: "// &
806  "MAXIMUM_DEPTH has a non-sensical value! Was it set?")
807 
808  do k=1,nz
809  e0(k) = -g%max_depth * real(k-1) / real(nz)
810  enddo
811 
812  do j=js,je ; do i=is,ie
813  ! This sets the initial thickness (in m) of the layers. The
814  ! thicknesses are set to insure that: 1. each layer is at least an
815  ! Angstrom thick, and 2. the interfaces are where they should be
816  ! based on the resting depths and interface height perturbations,
817  ! as long at this doesn't interfere with 1.
818  eta1d(nz+1) = -g%bathyT(i,j)
819  do k=nz,1,-1
820  eta1d(k) = e0(k)
821  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z)) then
822  eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
823  h(i,j,k) = gv%Angstrom_H
824  else
825  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
826  endif
827  enddo
828  enddo ; enddo
829 
830  call calltree_leave(trim(mdl)//'()')
831 end subroutine initialize_thickness_uniform
832 
833 !> Initialize thickness from a 1D list
834 subroutine initialize_thickness_list(h, G, GV, US, param_file, just_read_params)
835  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
836  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
837  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
838  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
839  intent(out) :: h !< The thickness that is being initialized [H ~> m or kg m-2].
840  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file
841  !! to parse for model parameter values.
842  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
843  !! only read parameters without changing h.
844  ! Local variables
845  character(len=40) :: mdl = "initialize_thickness_list" ! This subroutine's name.
846  real :: e0(SZK_(G)+1) ! The resting interface heights, in depth units [Z ~> m],
847  ! usually negative because it is positive upward.
848  real :: eta1D(SZK_(G)+1)! Interface height relative to the sea surface
849  ! positive upward, in depth units [Z ~> m].
850  logical :: just_read ! If true, just read parameters but set nothing.
851  character(len=200) :: filename, eta_file, inputdir ! Strings for file/path
852  character(len=72) :: eta_var
853  integer :: i, j, k, is, ie, js, je, nz
854 
855  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
856 
857  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
858 
859  call get_param(param_file, mdl, "INTERFACE_IC_FILE", eta_file, &
860  "The file from which horizontal mean initial conditions "//&
861  "for interface depths can be read.", fail_if_missing=.true.)
862  call get_param(param_file, mdl, "INTERFACE_IC_VAR", eta_var, &
863  "The variable name for horizontal mean initial conditions "//&
864  "for interface depths relative to mean sea level.", &
865  default="eta")
866 
867  if (just_read) return
868 
869  call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
870 
871  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
872  filename = trim(slasher(inputdir))//trim(eta_file)
873  call log_param(param_file, mdl, "INPUTDIR/INTERFACE_IC_FILE", filename)
874 
875  e0(:) = 0.0
876  call mom_read_data(filename, eta_var, e0(:), scale=us%m_to_Z)
877 
878  if ((abs(e0(1)) - 0.0) > 0.001) then
879  ! This list probably starts with the interior interface, so shift it up.
880  do k=nz+1,2,-1 ; e0(k) = e0(k-1) ; enddo
881  e0(1) = 0.0
882  endif
883 
884  if (e0(2) > e0(1)) then ! Switch to the convention for interface heights increasing upward.
885  do k=1,nz ; e0(k) = -e0(k) ; enddo
886  endif
887 
888  do j=js,je ; do i=is,ie
889  ! This sets the initial thickness (in m) of the layers. The
890  ! thicknesses are set to insure that: 1. each layer is at least an
891  ! Angstrom thick, and 2. the interfaces are where they should be
892  ! based on the resting depths and interface height perturbations,
893  ! as long at this doesn't interfere with 1.
894  eta1d(nz+1) = -g%bathyT(i,j)
895  do k=nz,1,-1
896  eta1d(k) = e0(k)
897  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z)) then
898  eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
899  h(i,j,k) = gv%Angstrom_H
900  else
901  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
902  endif
903  enddo
904  enddo ; enddo
905 
906  call calltree_leave(trim(mdl)//'()')
907 end subroutine initialize_thickness_list
908 
909 !> Search density space for location of layers (not implemented!)
910 subroutine initialize_thickness_search
911  call mom_error(fatal," MOM_state_initialization.F90, initialize_thickness_search: NOT IMPLEMENTED")
912 end subroutine initialize_thickness_search
913 
914 !> Converts thickness from geometric to pressure units
915 subroutine convert_thickness(h, G, GV, US, tv)
916  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
917  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
918  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
919  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
920  intent(inout) :: h !< Input geometric layer thicknesses being converted
921  !! to layer pressure [H ~> m or kg m-2].
922  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various
923  !! thermodynamic variables
924  ! Local variables
925  real, dimension(SZI_(G),SZJ_(G)) :: &
926  p_top, p_bot
927  real :: dz_geo(SZI_(G),SZJ_(G)) ! The change in geopotential height
928  ! across a layer [m2 s-2].
929  real :: rho(SZI_(G))
930  real :: I_gEarth
931  real :: Hm_rho_to_Pa ! A conversion factor from the input geometric thicknesses times the
932  ! layer densities into Pa [Pa m3 H-1 kg-1 ~> s-2 m2 or s-2 m5 kg-1].
933  logical :: Boussinesq
934  integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
935  integer :: itt, max_itt
936 
937  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
938  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
939  max_itt = 10
940  boussinesq = gv%Boussinesq
941  i_gearth = 1.0 / (gv%mks_g_Earth)
942  hm_rho_to_pa = gv%mks_g_Earth * gv%H_to_m ! = GV%H_to_Pa / GV%Rho0
943 
944  if (boussinesq) then
945  call mom_error(fatal,"Not yet converting thickness with Boussinesq approx.")
946  else
947  if (associated(tv%eqn_of_state)) then
948  do j=jsq,jeq+1 ; do i=isq,ieq+1
949  p_bot(i,j) = 0.0 ; p_top(i,j) = 0.0
950  enddo ; enddo
951  do k=1,nz
952  do j=js,je
953  do i=is,ie ; p_top(i,j) = p_bot(i,j) ; enddo
954  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_top(:,j), rho, &
955  is, ie-is+1, tv%eqn_of_state)
956  do i=is,ie
957  p_bot(i,j) = p_top(i,j) + hm_rho_to_pa * (h(i,j,k) * rho(i))
958  enddo
959  enddo
960 
961  do itt=1,max_itt
962  call int_specific_vol_dp(tv%T(:,:,k), tv%S(:,:,k), p_top, p_bot, &
963  0.0, g%HI, tv%eqn_of_state, dz_geo)
964  if (itt < max_itt) then ; do j=js,je
965  call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_bot(:,j), rho, &
966  is, ie-is+1, tv%eqn_of_state)
967  ! Use Newton's method to correct the bottom value.
968  ! The hydrostatic equation is linear to such a
969  ! high degree that no bounds-checking is needed.
970  do i=is,ie
971  p_bot(i,j) = p_bot(i,j) + rho(i) * &
972  (hm_rho_to_pa*h(i,j,k) - dz_geo(i,j))
973  enddo
974  enddo ; endif
975  enddo
976 
977  do j=js,je ; do i=is,ie
978  h(i,j,k) = (p_bot(i,j) - p_top(i,j)) * gv%kg_m2_to_H * i_gearth
979  enddo ; enddo
980  enddo
981  else
982  do k=1,nz ; do j=js,je ; do i=is,ie
983  h(i,j,k) = (h(i,j,k) * gv%Rlay(k)) * hm_rho_to_pa * gv%kg_m2_to_H**2
984  ! This is mathematically equivalent to
985  ! h(i,j,k) = h(i,j,k) * (GV%Rlay(k) / GV%Rho0)
986  enddo ; enddo ; enddo
987  endif
988  endif
989 
990 end subroutine convert_thickness
991 
992 !> Depress the sea-surface based on an initial condition file
993 subroutine depress_surface(h, G, GV, US, param_file, tv, just_read_params)
994  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
995  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
996  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
997  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
998  intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
999  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1000  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic variables
1001  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1002  !! only read parameters without changing h.
1003  ! Local variables
1004  real, dimension(SZI_(G),SZJ_(G)) :: &
1005  eta_sfc ! The free surface height that the model should use [m].
1006  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
1007  eta ! The free surface height that the model should use [m].
1008  real :: dilate ! A ratio by which layers are dilated [nondim].
1009  real :: scale_factor ! A scaling factor for the eta_sfc values that are read
1010  ! in, which can be used to change units, for example.
1011  character(len=40) :: mdl = "depress_surface" ! This subroutine's name.
1012  character(len=200) :: inputdir, eta_srf_file ! Strings for file/path
1013  character(len=200) :: filename, eta_srf_var ! Strings for file/path
1014  logical :: just_read ! If true, just read parameters but set nothing.
1015  integer :: i, j, k, is, ie, js, je, nz
1016  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1017 
1018  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1019 
1020  ! Read the surface height (or pressure) from a file.
1021 
1022  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1023  inputdir = slasher(inputdir)
1024  call get_param(param_file, mdl, "SURFACE_HEIGHT_IC_FILE", eta_srf_file,&
1025  "The initial condition file for the surface height.", &
1026  fail_if_missing=.not.just_read, do_not_log=just_read)
1027  call get_param(param_file, mdl, "SURFACE_HEIGHT_IC_VAR", eta_srf_var, &
1028  "The initial condition variable for the surface height.",&
1029  default="SSH", do_not_log=just_read)
1030  filename = trim(inputdir)//trim(eta_srf_file)
1031  if (.not.just_read) &
1032  call log_param(param_file, mdl, "INPUTDIR/SURFACE_HEIGHT_IC_FILE", filename)
1033 
1034  call get_param(param_file, mdl, "SURFACE_HEIGHT_IC_SCALE", scale_factor, &
1035  "A scaling factor to convert SURFACE_HEIGHT_IC_VAR into "//&
1036  "units of m", units="variable", default=1.0, do_not_log=just_read)
1037 
1038  if (just_read) return ! All run-time parameters have been read, so return.
1039 
1040  call mom_read_data(filename, eta_srf_var, eta_sfc, g%Domain, scale=scale_factor)
1041 
1042  ! Convert thicknesses to interface heights.
1043  call find_eta(h, tv, g, gv, us, eta, eta_to_m=1.0)
1044 
1045  do j=js,je ; do i=is,ie ; if (g%mask2dT(i,j) > 0.0) then
1046 ! if (eta_sfc(i,j) < eta(i,j,nz+1)) then
1047  ! Issue a warning?
1048 ! endif
1049  if (eta_sfc(i,j) > eta(i,j,1)) then
1050  ! Dilate the water column to agree, but only up to 10-fold.
1051  if (eta_sfc(i,j) - eta(i,j,nz+1) > 10.0*(eta(i,j,1) - eta(i,j,nz+1))) then
1052  dilate = 10.0
1053  call mom_error(warning, "Free surface height dilation attempted "//&
1054  "to exceed 10-fold.", all_print=.true.)
1055  else
1056  dilate = (eta_sfc(i,j) - eta(i,j,nz+1)) / (eta(i,j,1) - eta(i,j,nz+1))
1057  endif
1058  do k=1,nz ; h(i,j,k) = h(i,j,k) * dilate ; enddo
1059  elseif (eta(i,j,1) > eta_sfc(i,j)) then
1060  ! Remove any mass that is above the target free surface.
1061  do k=1,nz
1062  if (eta(i,j,k) <= eta_sfc(i,j)) exit
1063  if (eta(i,j,k+1) >= eta_sfc(i,j)) then
1064  h(i,j,k) = gv%Angstrom_H
1065  else
1066  h(i,j,k) = max(gv%Angstrom_H, h(i,j,k) * &
1067  (eta_sfc(i,j) - eta(i,j,k+1)) / (eta(i,j,k) - eta(i,j,k+1)) )
1068  endif
1069  enddo
1070  endif
1071  endif ; enddo ; enddo
1072 
1073 end subroutine depress_surface
1074 
1075 !> Adjust the layer thicknesses by cutting away the top of each model column at the depth
1076 !! where the hydrostatic pressure matches an imposed surface pressure read from file.
1077 subroutine trim_for_ice(PF, G, GV, US, ALE_CSp, tv, h, just_read_params)
1078  type(param_file_type), intent(in) :: PF !< Parameter file structure
1079  type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
1080  type(verticalgrid_type), intent(in) :: GV !< Vertical grid structure
1081  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1082  type(ale_cs), pointer :: ALE_CSp !< ALE control structure
1083  type(thermo_var_ptrs), intent(inout) :: tv !< Thermodynamics structure
1084  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1085  intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
1086  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1087  !! only read parameters without changing h.
1088  ! Local variables
1089  character(len=200) :: mdl = "trim_for_ice"
1090  real, dimension(SZI_(G),SZJ_(G)) :: p_surf ! Imposed pressure on ocean at surface [Pa]
1091  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: S_t, S_b ! Top and bottom edge values for reconstructions
1092  real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: T_t, T_b ! of salinity and temperature within each layer.
1093  character(len=200) :: inputdir, filename, p_surf_file, p_surf_var ! Strings for file/path
1094  real :: scale_factor ! A file-dependent scaling vactor for the input pressurs.
1095  real :: min_thickness ! The minimum layer thickness, recast into Z units.
1096  integer :: i, j, k
1097  logical :: just_read ! If true, just read parameters but set nothing.
1098  logical :: use_remapping ! If true, remap the initial conditions.
1099  type(remapping_cs), pointer :: remap_CS => null()
1100 
1101  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1102 
1103  call get_param(pf, mdl, "SURFACE_PRESSURE_FILE", p_surf_file, &
1104  "The initial condition file for the surface height.", &
1105  fail_if_missing=.not.just_read, do_not_log=just_read)
1106  call get_param(pf, mdl, "SURFACE_PRESSURE_VAR", p_surf_var, &
1107  "The initial condition variable for the surface height.", &
1108  units="kg m-2", default="", do_not_log=just_read)
1109  call get_param(pf, mdl, "INPUTDIR", inputdir, default=".", do_not_log=.true.)
1110  filename = trim(slasher(inputdir))//trim(p_surf_file)
1111  if (.not.just_read) call log_param(pf, mdl, "!INPUTDIR/SURFACE_HEIGHT_IC_FILE", filename)
1112 
1113  call get_param(pf, mdl, "SURFACE_PRESSURE_SCALE", scale_factor, &
1114  "A scaling factor to convert SURFACE_PRESSURE_VAR from "//&
1115  "file SURFACE_PRESSURE_FILE into a surface pressure.", &
1116  units="file dependent", default=1., do_not_log=just_read)
1117  call get_param(pf, mdl, "MIN_THICKNESS", min_thickness, 'Minimum layer thickness', &
1118  units='m', default=1.e-3, do_not_log=just_read, scale=us%m_to_Z)
1119  call get_param(pf, mdl, "TRIMMING_USES_REMAPPING", use_remapping, &
1120  'When trimming the column, also remap T and S.', &
1121  default=.false., do_not_log=just_read)
1122 
1123  if (just_read) return ! All run-time parameters have been read, so return.
1124 
1125  call mom_read_data(filename, p_surf_var, p_surf, g%Domain, scale=scale_factor)
1126 
1127  if (use_remapping) then
1128  allocate(remap_cs)
1129  call initialize_remapping(remap_cs, 'PLM', boundary_extrapolation=.true.)
1130  endif
1131 
1132  ! Find edge values of T and S used in reconstructions
1133  if ( associated(ale_csp) ) then ! This should only be associated if we are in ALE mode
1134  call pressure_gradient_plm(ale_csp, s_t, s_b, t_t, t_b, g, gv, tv, h, .true.)
1135  else
1136 ! call MOM_error(FATAL, "trim_for_ice: Does not work without ALE mode")
1137  do k=1,g%ke ; do j=g%jsc,g%jec ; do i=g%isc,g%iec
1138  t_t(i,j,k) = tv%T(i,j,k) ; t_b(i,j,k) = tv%T(i,j,k)
1139  s_t(i,j,k) = tv%S(i,j,k) ; s_b(i,j,k) = tv%S(i,j,k)
1140  enddo ; enddo ; enddo
1141  endif
1142 
1143  do j=g%jsc,g%jec ; do i=g%isc,g%iec
1144  call cut_off_column_top(gv%ke, tv, gv, gv%mks_g_Earth*us%Z_to_m, g%bathyT(i,j), &
1145  min_thickness, tv%T(i,j,:), t_t(i,j,:), t_b(i,j,:), &
1146  tv%S(i,j,:), s_t(i,j,:), s_b(i,j,:), p_surf(i,j), h(i,j,:), remap_cs, &
1147  z_tol=1.0e-5*us%m_to_Z)
1148  enddo ; enddo
1149 
1150 end subroutine trim_for_ice
1151 
1152 
1153 !> Adjust the layer thicknesses by removing the top of the water column above the
1154 !! depth where the hydrostatic pressure matches p_surf
1155 subroutine cut_off_column_top(nk, tv, GV, G_earth, depth, min_thickness, &
1156  T, T_t, T_b, S, S_t, S_b, p_surf, h, remap_CS, z_tol)
1157  integer, intent(in) :: nk !< Number of layers
1158  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
1159  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1160  real, intent(in) :: G_earth !< Gravitational acceleration [m2 Z-1 s-2 ~> m s-2]
1161  real, intent(in) :: depth !< Depth of ocean column [Z ~> m].
1162  real, intent(in) :: min_thickness !< Smallest thickness allowed [Z ~> m].
1163  real, dimension(nk), intent(inout) :: T !< Layer mean temperature [degC]
1164  real, dimension(nk), intent(in) :: T_t !< Temperature at top of layer [degC]
1165  real, dimension(nk), intent(in) :: T_b !< Temperature at bottom of layer [degC]
1166  real, dimension(nk), intent(inout) :: S !< Layer mean salinity [ppt]
1167  real, dimension(nk), intent(in) :: S_t !< Salinity at top of layer [ppt]
1168  real, dimension(nk), intent(in) :: S_b !< Salinity at bottom of layer [ppt]
1169  real, intent(in) :: p_surf !< Imposed pressure on ocean at surface [Pa]
1170  real, dimension(nk), intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
1171  type(remapping_cs), pointer :: remap_CS !< Remapping structure for remapping T and S,
1172  !! if associated
1173  real, optional, intent(in) :: z_tol !< The tolerance with which to find the depth
1174  !! matching the specified pressure [Z ~> m].
1175 
1176  ! Local variables
1177  real, dimension(nk+1) :: e ! Top and bottom edge values for reconstructions
1178  real, dimension(nk) :: h0, S0, T0, h1, S1, T1
1179  real :: P_t, P_b, z_out, e_top
1180  integer :: k
1181 
1182  ! Calculate original interface positions
1183  e(nk+1) = -depth
1184  do k=nk,1,-1
1185  e(k) = e(k+1) + gv%H_to_Z*h(k)
1186  h0(k) = h(nk+1-k) ! Keep a copy to use in remapping
1187  enddo
1188 
1189  p_t = 0.
1190  e_top = e(1)
1191  do k=1,nk
1192  call find_depth_of_pressure_in_cell(t_t(k), t_b(k), s_t(k), s_b(k), e(k), e(k+1), &
1193  p_t, p_surf, gv%Rho0, g_earth, tv%eqn_of_state, &
1194  p_b, z_out, z_tol=z_tol)
1195  if (z_out>=e(k)) then
1196  ! Imposed pressure was less that pressure at top of cell
1197  exit
1198  elseif (z_out<=e(k+1)) then
1199  ! Imposed pressure was greater than pressure at bottom of cell
1200  e_top = e(k+1)
1201  else
1202  ! Imposed pressure was fell between pressures at top and bottom of cell
1203  e_top = z_out
1204  exit
1205  endif
1206  p_t = p_b
1207  enddo
1208  if (e_top<e(1)) then
1209  ! Clip layers from the top down, if at all
1210  do k=1,nk
1211  if (e(k) > e_top) then
1212  ! Original e(K) is too high
1213  e(k) = e_top
1214  e_top = e_top - min_thickness ! Next interface must be at least this deep
1215  endif
1216  ! This layer needs trimming
1217  h(k) = gv%Z_to_H * max( min_thickness, e(k) - e(k+1) )
1218  if (e(k) < e_top) exit ! No need to go further
1219  enddo
1220  endif
1221 
1222  ! Now we need to remap but remapping assumes the surface is at the
1223  ! same place in the two columns so we turn the column upside down.
1224  if (associated(remap_cs)) then
1225  do k=1,nk
1226  s0(k) = s(nk+1-k)
1227  t0(k) = t(nk+1-k)
1228  h1(k) = h(nk+1-k)
1229  enddo
1230  call remapping_core_h(remap_cs, nk, h0, t0, nk, h1, t1, 1.0e-30*gv%m_to_H, 1.0e-10*gv%m_to_H)
1231  call remapping_core_h(remap_cs, nk, h0, s0, nk, h1, s1, 1.0e-30*gv%m_to_H, 1.0e-10*gv%m_to_H)
1232  do k=1,nk
1233  s(k) = s1(nk+1-k)
1234  t(k) = t1(nk+1-k)
1235  enddo
1236  endif
1237 
1238 end subroutine cut_off_column_top
1239 
1240 !> Initialize horizontal velocity components from file
1241 subroutine initialize_velocity_from_file(u, v, G, param_file, just_read_params)
1242  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1243  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1244  intent(out) :: u !< The zonal velocity that is being initialized [m s-1]
1245  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1246  intent(out) :: v !< The meridional velocity that is being initialized [m s-1]
1247  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to
1248  !! parse for modelparameter values.
1249  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1250  !! only read parameters without changing h.
1251  ! Local variables
1252  character(len=40) :: mdl = "initialize_velocity_from_file" ! This subroutine's name.
1253  character(len=200) :: filename,velocity_file,inputdir ! Strings for file/path
1254  logical :: just_read ! If true, just read parameters but set nothing.
1255 
1256  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1257 
1258  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
1259 
1260  call get_param(param_file, mdl, "VELOCITY_FILE", velocity_file, &
1261  "The name of the velocity initial condition file.", &
1262  fail_if_missing=.not.just_read, do_not_log=just_read)
1263  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1264  inputdir = slasher(inputdir)
1265 
1266  if (just_read) return ! All run-time parameters have been read, so return.
1267 
1268  filename = trim(inputdir)//trim(velocity_file)
1269  call log_param(param_file, mdl, "INPUTDIR/VELOCITY_FILE", filename)
1270 
1271  if (.not.file_exists(filename, g%Domain)) call mom_error(fatal, &
1272  " initialize_velocity_from_file: Unable to open "//trim(filename))
1273 
1274  ! Read the velocities from a netcdf file.
1275  call mom_read_vector(filename, "u", "v", u(:,:,:), v(:,:,:),g%Domain)
1276 
1277  call calltree_leave(trim(mdl)//'()')
1278 end subroutine initialize_velocity_from_file
1279 
1280 !> Initialize horizontal velocity components to zero.
1281 subroutine initialize_velocity_zero(u, v, G, param_file, just_read_params)
1282  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1283  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1284  intent(out) :: u !< The zonal velocity that is being initialized [m s-1]
1285  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1286  intent(out) :: v !< The meridional velocity that is being initialized [m s-1]
1287  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to
1288  !! parse for modelparameter values.
1289  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1290  !! only read parameters without changing h.
1291  ! Local variables
1292  character(len=200) :: mdl = "initialize_velocity_zero" ! This subroutine's name.
1293  logical :: just_read ! If true, just read parameters but set nothing.
1294  integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1295  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1296  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1297 
1298  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1299 
1300  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
1301 
1302  if (just_read) return ! All run-time parameters have been read, so return.
1303 
1304  do k=1,nz ; do j=js,je ; do i=isq,ieq
1305  u(i,j,k) = 0.0
1306  enddo ; enddo ; enddo
1307  do k=1,nz ; do j=jsq,jeq ; do i=is,ie
1308  v(i,j,k) = 0.0
1309  enddo ; enddo ; enddo
1310 
1311  call calltree_leave(trim(mdl)//'()')
1312 end subroutine initialize_velocity_zero
1313 
1314 !> Sets the initial velocity components to uniform
1315 subroutine initialize_velocity_uniform(u, v, G, param_file, just_read_params)
1316  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1317  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1318  intent(out) :: u !< The zonal velocity that is being initialized [m s-1]
1319  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1320  intent(out) :: v !< The meridional velocity that is being initialized [m s-1]
1321  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to
1322  !! parse for modelparameter values.
1323  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1324  !! only read parameters without changing h.
1325  ! Local variables
1326  integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1327  real :: initial_u_const, initial_v_const
1328  logical :: just_read ! If true, just read parameters but set nothing.
1329  character(len=200) :: mdl = "initialize_velocity_uniform" ! This subroutine's name.
1330  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1331  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1332 
1333  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1334 
1335  call get_param(param_file, mdl, "INITIAL_U_CONST", initial_u_const, &
1336  "A initial uniform value for the zonal flow.", &
1337  units="m s-1", fail_if_missing=.not.just_read, do_not_log=just_read)
1338  call get_param(param_file, mdl, "INITIAL_V_CONST", initial_v_const, &
1339  "A initial uniform value for the meridional flow.", &
1340  units="m s-1", fail_if_missing=.not.just_read, do_not_log=just_read)
1341 
1342  if (just_read) return ! All run-time parameters have been read, so return.
1343 
1344  do k=1,nz ; do j=js,je ; do i=isq,ieq
1345  u(i,j,k) = initial_u_const
1346  enddo ; enddo ; enddo
1347  do k=1,nz ; do j=jsq,jeq ; do i=is,ie
1348  v(i,j,k) = initial_v_const
1349  enddo ; enddo ; enddo
1350 
1351 end subroutine initialize_velocity_uniform
1352 
1353 !> Sets the initial velocity components to be circular with
1354 !! no flow at edges of domain and center.
1355 subroutine initialize_velocity_circular(u, v, G, param_file, just_read_params)
1356  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1357  real, dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
1358  intent(out) :: u !< The zonal velocity that is being initialized [m s-1]
1359  real, dimension(SZI_(G),SZJB_(G),SZK_(G)), &
1360  intent(out) :: v !< The meridional velocity that is being initialized [m s-1]
1361  type(param_file_type), intent(in) :: param_file !< A structure indicating the open file to
1362  !! parse for model parameter values.
1363  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1364  !! only read parameters without changing h.
1365  ! Local variables
1366  character(len=200) :: mdl = "initialize_velocity_circular"
1367  real :: circular_max_u
1368  real :: dpi, psi1, psi2
1369  logical :: just_read ! If true, just read parameters but set nothing.
1370  integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
1371  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1372  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1373 
1374  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1375 
1376  call get_param(param_file, mdl, "CIRCULAR_MAX_U", circular_max_u, &
1377  "The amplitude of zonal flow from which to scale the "// &
1378  "circular stream function [m s-1].", &
1379  units="m s-1", default=0., do_not_log=just_read)
1380 
1381  if (just_read) return ! All run-time parameters have been read, so return.
1382 
1383  dpi=acos(0.0)*2.0 ! pi
1384 
1385  do k=1,nz ; do j=js,je ; do i=isq,ieq
1386  psi1 = my_psi(i,j)
1387  psi2 = my_psi(i,j-1)
1388  u(i,j,k) = (psi1-psi2)/g%dy_Cu(i,j)! *(circular_max_u*G%len_lon/(2.0*dpi))
1389  enddo ; enddo ; enddo
1390  do k=1,nz ; do j=jsq,jeq ; do i=is,ie
1391  psi1 = my_psi(i,j)
1392  psi2 = my_psi(i-1,j)
1393  v(i,j,k) = (psi2-psi1)/g%dx_Cv(i,j)! *(circular_max_u*G%len_lon/(2.0*dpi))
1394  enddo ; enddo ; enddo
1395 
1396  contains
1397 
1398  !> Returns the value of a circular stream function at (ig,jg)
1399  real function my_psi(ig,jg)
1400  integer :: ig !< Global i-index
1401  integer :: jg !< Global j-index
1402  ! Local variables
1403  real :: x, y, r
1404 
1405  x = 2.0*(g%geoLonBu(ig,jg)-g%west_lon)/g%len_lon-1.0 ! -1<x<1
1406  y = 2.0*(g%geoLatBu(ig,jg)-g%south_lat)/g%len_lat-1.0 ! -1<y<1
1407  r = sqrt( x**2 + y**2 ) ! Circulat stream fn nis fn of radius only
1408  r = min(1.0,r) ! Flatten stream function in corners of box
1409  my_psi = 0.5*(1.0 - cos(dpi*r))
1410  my_psi = my_psi * (circular_max_u*g%len_lon*1e3/dpi) ! len_lon is in km
1411  end function my_psi
1412 
1413 end subroutine initialize_velocity_circular
1414 
1415 !> Initializes temperature and salinity from file
1416 subroutine initialize_temp_salt_from_file(T, S, G, param_file, just_read_params)
1417  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1418  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: T !< The potential temperature that is
1419  !! being initialized [degC]
1420  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: S !< The salinity that is
1421  !! being initialized [ppt]
1422  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1423  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1424  !! only read parameters without changing h.
1425  ! Local variables
1426  logical :: just_read ! If true, just read parameters but set nothing.
1427  character(len=200) :: filename, salt_filename ! Full paths to input files
1428  character(len=200) :: ts_file, salt_file, inputdir ! Strings for file/path
1429  character(len=40) :: mdl = "initialize_temp_salt_from_file"
1430  character(len=64) :: temp_var, salt_var ! Temperature and salinity names in files
1431 
1432  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1433 
1434  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
1435 
1436  call get_param(param_file, mdl, "TS_FILE", ts_file, &
1437  "The initial condition file for temperature.", &
1438  fail_if_missing=.not.just_read, do_not_log=just_read)
1439  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1440  inputdir = slasher(inputdir)
1441 
1442  filename = trim(inputdir)//trim(ts_file)
1443  if (.not.just_read) call log_param(param_file, mdl, "INPUTDIR/TS_FILE", filename)
1444  call get_param(param_file, mdl, "TEMP_IC_VAR", temp_var, &
1445  "The initial condition variable for potential temperature.", &
1446  default="PTEMP", do_not_log=just_read)
1447  call get_param(param_file, mdl, "SALT_IC_VAR", salt_var, &
1448  "The initial condition variable for salinity.", &
1449  default="SALT", do_not_log=just_read)
1450  call get_param(param_file, mdl, "SALT_FILE", salt_file, &
1451  "The initial condition file for salinity.", &
1452  default=trim(ts_file), do_not_log=just_read)
1453 
1454  if (just_read) return ! All run-time parameters have been read, so return.
1455 
1456  if (.not.file_exists(filename, g%Domain)) call mom_error(fatal, &
1457  " initialize_temp_salt_from_file: Unable to open "//trim(filename))
1458 
1459  ! Read the temperatures and salinities from netcdf files.
1460  call mom_read_data(filename, temp_var, t(:,:,:), g%Domain)
1461 
1462  salt_filename = trim(inputdir)//trim(salt_file)
1463  if (.not.file_exists(salt_filename, g%Domain)) call mom_error(fatal, &
1464  " initialize_temp_salt_from_file: Unable to open "//trim(salt_filename))
1465 
1466  call mom_read_data(salt_filename, salt_var, s(:,:,:), g%Domain)
1467 
1468  call calltree_leave(trim(mdl)//'()')
1469 end subroutine initialize_temp_salt_from_file
1470 
1471 !> Initializes temperature and salinity from a 1D profile
1472 subroutine initialize_temp_salt_from_profile(T, S, G, param_file, just_read_params)
1473  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1474  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: T !< The potential temperature that is
1475  !! being initialized [degC]
1476  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: S !< The salinity that is
1477  !! being initialized [ppt]
1478  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
1479  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1480  !! only read parameters without changing h.
1481  ! Local variables
1482  real, dimension(SZK_(G)) :: T0, S0
1483  integer :: i, j, k
1484  logical :: just_read ! If true, just read parameters but set nothing.
1485  character(len=200) :: filename, ts_file, inputdir ! Strings for file/path
1486  character(len=40) :: mdl = "initialize_temp_salt_from_profile"
1487 
1488  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1489 
1490  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
1491 
1492  call get_param(param_file, mdl, "TS_FILE", ts_file, &
1493  "The file with the reference profiles for temperature "//&
1494  "and salinity.", fail_if_missing=.not.just_read, do_not_log=just_read)
1495 
1496  if (just_read) return ! All run-time parameters have been read, so return.
1497 
1498  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1499  inputdir = slasher(inputdir)
1500  filename = trim(inputdir)//trim(ts_file)
1501  call log_param(param_file, mdl, "INPUTDIR/TS_FILE", filename)
1502  if (.not.file_exists(filename)) call mom_error(fatal, &
1503  " initialize_temp_salt_from_profile: Unable to open "//trim(filename))
1504 
1505  ! Read the temperatures and salinities from a netcdf file.
1506  call mom_read_data(filename, "PTEMP", t0(:))
1507  call mom_read_data(filename, "SALT", s0(:))
1508 
1509  do k=1,g%ke ; do j=g%jsc,g%jec ; do i=g%isc,g%iec
1510  t(i,j,k) = t0(k) ; s(i,j,k) = s0(k)
1511  enddo ; enddo ; enddo
1512 
1513  call calltree_leave(trim(mdl)//'()')
1514 end subroutine initialize_temp_salt_from_profile
1515 
1516 !> Initializes temperature and salinity by fitting to density
1517 subroutine initialize_temp_salt_fit(T, S, G, GV, param_file, eqn_of_state, P_Ref, just_read_params)
1518  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1519  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1520  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: T !< The potential temperature that is
1521  !! being initialized [degC].
1522  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: S !< The salinity that is being
1523  !! initialized [ppt].
1524  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time
1525  !! parameters.
1526  type(eos_type), pointer :: eqn_of_state !< Integer that selects the equatio of state.
1527  real, intent(in) :: P_Ref !< The coordinate-density reference pressure
1528  !! [Pa].
1529  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1530  !! only read parameters without changing h.
1531  ! Local variables
1532  real :: T0(SZK_(G)) ! Layer potential temperatures [degC]
1533  real :: S0(SZK_(G)) ! Layer salinities [degC]
1534  real :: T_Ref ! Reference Temperature [degC]
1535  real :: S_Ref ! Reference Salinity [ppt]
1536  real :: pres(SZK_(G)) ! An array of the reference pressure [Pa].
1537  real :: drho_dT(SZK_(G)) ! Derivative of density with temperature [kg m-3 degC-1].
1538  real :: drho_dS(SZK_(G)) ! Derivative of density with salinity [kg m-3 ppt-1].
1539  real :: rho_guess(SZK_(G)) ! Potential density at T0 & S0 [kg m-3].
1540  logical :: fit_salin ! If true, accept the prescribed temperature and fit the salinity.
1541  logical :: just_read ! If true, just read parameters but set nothing.
1542  character(len=40) :: mdl = "initialize_temp_salt_fit" ! This subroutine's name.
1543  integer :: i, j, k, itt, nz
1544  nz = g%ke
1545 
1546  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1547 
1548  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
1549 
1550  call get_param(param_file, mdl, "T_REF", t_ref, &
1551  "A reference temperature used in initialization.", &
1552  units="degC", fail_if_missing=.not.just_read, do_not_log=just_read)
1553  call get_param(param_file, mdl, "S_REF", s_ref, &
1554  "A reference salinity used in initialization.", units="PSU", &
1555  default=35.0, do_not_log=just_read)
1556  call get_param(param_file, mdl, "FIT_SALINITY", fit_salin, &
1557  "If true, accept the prescribed temperature and fit the "//&
1558  "salinity; otherwise take salinity and fit temperature.", &
1559  default=.false., do_not_log=just_read)
1560 
1561  if (just_read) return ! All run-time parameters have been read, so return.
1562 
1563  do k=1,nz
1564  pres(k) = p_ref ; s0(k) = s_ref
1565  t0(k) = t_ref
1566  enddo
1567 
1568  call calculate_density(t0(1),s0(1),pres(1),rho_guess(1),eqn_of_state)
1569  call calculate_density_derivs(t0,s0,pres,drho_dt,drho_ds,1,1,eqn_of_state)
1570 
1571  if (fit_salin) then
1572  ! A first guess of the layers' temperatures.
1573  do k=nz,1,-1
1574  s0(k) = max(0.0, s0(1) + (gv%Rlay(k) - rho_guess(1)) / drho_ds(1))
1575  enddo
1576  ! Refine the guesses for each layer.
1577  do itt=1,6
1578  call calculate_density(t0,s0,pres,rho_guess,1,nz,eqn_of_state)
1579  call calculate_density_derivs(t0,s0,pres,drho_dt,drho_ds,1,nz,eqn_of_state)
1580  do k=1,nz
1581  s0(k) = max(0.0, s0(k) + (gv%Rlay(k) - rho_guess(k)) / drho_ds(k))
1582  enddo
1583  enddo
1584  else
1585  ! A first guess of the layers' temperatures.
1586  do k=nz,1,-1
1587  t0(k) = t0(1) + (gv%Rlay(k) - rho_guess(1)) / drho_dt(1)
1588  enddo
1589  do itt=1,6
1590  call calculate_density(t0,s0,pres,rho_guess,1,nz,eqn_of_state)
1591  call calculate_density_derivs(t0,s0,pres,drho_dt,drho_ds,1,nz,eqn_of_state)
1592  do k=1,nz
1593  t0(k) = t0(k) + (gv%Rlay(k) - rho_guess(k)) / drho_dt(k)
1594  enddo
1595  enddo
1596  endif
1597 
1598  do k=1,nz ; do j=g%jsd,g%jed ; do i=g%isd,g%ied
1599  t(i,j,k) = t0(k) ; s(i,j,k) = s0(k)
1600  enddo ; enddo ; enddo
1601 
1602  call calltree_leave(trim(mdl)//'()')
1603 end subroutine initialize_temp_salt_fit
1604 
1605 !> Initializes T and S with linear profiles according to reference surface
1606 !! layer salinity and temperature and a specified range.
1607 !!
1608 !! \remark Note that the linear distribution is set up with respect to the layer
1609 !! number, not the physical position).
1610 subroutine initialize_temp_salt_linear(T, S, G, param_file, just_read_params)
1611  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1612  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: T !< The potential temperature that is
1613  !! being initialized [degC]
1614  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(out) :: S !< The salinity that is
1615  !! being initialized [ppt]
1616  type(param_file_type), intent(in) :: param_file !< A structure to parse for
1617  !! run-time parameters
1618  logical, optional, intent(in) :: just_read_params !< If present and true,
1619  !! this call will only read
1620  !! parameters without
1621  !! changing h.
1622 
1623  integer :: k
1624  real :: delta_S, delta_T
1625  real :: S_top, T_top ! Reference salinity and temerature within surface layer
1626  real :: S_range, T_range ! Range of salinities and temperatures over the vertical
1627  real :: delta
1628  logical :: just_read ! If true, just read parameters but set nothing.
1629  character(len=40) :: mdl = "initialize_temp_salt_linear" ! This subroutine's name.
1630 
1631  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
1632 
1633  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
1634  call get_param(param_file, mdl, "T_TOP", t_top, &
1635  "Initial temperature of the top surface.", &
1636  units="degC", fail_if_missing=.not.just_read, do_not_log=just_read)
1637  call get_param(param_file, mdl, "T_RANGE", t_range, &
1638  "Initial temperature difference (top-bottom).", &
1639  units="degC", fail_if_missing=.not.just_read, do_not_log=just_read)
1640  call get_param(param_file, mdl, "S_TOP", s_top, &
1641  "Initial salinity of the top surface.", &
1642  units="PSU", fail_if_missing=.not.just_read, do_not_log=just_read)
1643  call get_param(param_file, mdl, "S_RANGE", s_range, &
1644  "Initial salinity difference (top-bottom).", &
1645  units="PSU", fail_if_missing=.not.just_read, do_not_log=just_read)
1646 
1647  if (just_read) return ! All run-time parameters have been read, so return.
1648 
1649  ! Prescribe salinity
1650 ! delta_S = S_range / ( G%ke - 1.0 )
1651 ! S(:,:,1) = S_top
1652 ! do k = 2,G%ke
1653 ! S(:,:,k) = S(:,:,k-1) + delta_S
1654 ! enddo
1655  do k = 1,g%ke
1656  s(:,:,k) = s_top - s_range*((real(k)-0.5)/real(g%ke))
1657  t(:,:,k) = t_top - t_range*((real(k)-0.5)/real(g%ke))
1658  enddo
1659 
1660  ! Prescribe temperature
1661 ! delta_T = T_range / ( G%ke - 1.0 )
1662 ! T(:,:,1) = T_top
1663 ! do k = 2,G%ke
1664 ! T(:,:,k) = T(:,:,k-1) + delta_T
1665 ! enddo
1666 ! delta = 1
1667 ! T(:,:,G%ke/2 - (delta-1):G%ke/2 + delta) = 1.0
1668 
1669  call calltree_leave(trim(mdl)//'()')
1670 end subroutine initialize_temp_salt_linear
1671 
1672 !> This subroutine sets the inverse restoration time (Idamp), and
1673 !! the values towards which the interface heights and an arbitrary
1674 !! number of tracers should be restored within each sponge. The
1675 !! interface height is always subject to damping, and must always be
1676 !! the first registered field.
1677 subroutine initialize_sponges_file(G, GV, US, use_temperature, tv, param_file, CSp, ALE_CSp, Time)
1678  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
1679  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
1680  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1681  logical, intent(in) :: use_temperature !< If true, T & S are state variables.
1682  type(thermo_var_ptrs), intent(in) :: tv !< A structure pointing to various thermodynamic
1683  !! variables.
1684  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters.
1685  type(sponge_cs), pointer :: CSp !< A pointer that is set to point to the control
1686  !! structure for this module (in layered mode).
1687  type(ale_sponge_cs), pointer :: ALE_CSp !< A pointer that is set to point to the control
1688  !! structure for this module (in ALE mode).
1689  type(time_type), intent(in) :: Time !< Time at the start of the run segment. Time_in
1690  !! overrides any value set for Time.
1691  ! Local variables
1692  real, allocatable, dimension(:,:,:) :: eta ! The target interface heights [Z ~> m].
1693  real, allocatable, dimension(:,:,:) :: h ! The target interface thicknesses [H ~> m or kg m-2].
1694 
1695  real, dimension (SZI_(G),SZJ_(G),SZK_(G)) :: &
1696  tmp, tmp2 ! A temporary array for tracers.
1697  real, dimension (SZI_(G),SZJ_(G)) :: &
1698  tmp_2d ! A temporary array for tracers.
1699 
1700  real :: Idamp(SZI_(G),SZJ_(G)) ! The inverse damping rate [s-1].
1701  real :: pres(SZI_(G)) ! An array of the reference pressure [Pa].
1702 
1703  integer :: i, j, k, is, ie, js, je, nz
1704  integer :: isd, ied, jsd, jed
1705  integer, dimension(4) :: siz
1706  integer :: nz_data ! The size of the sponge source grid
1707  character(len=40) :: potemp_var, salin_var, Idamp_var, eta_var
1708  character(len=40) :: mdl = "initialize_sponges_file"
1709  character(len=200) :: damping_file, state_file ! Strings for filenames
1710  character(len=200) :: filename, inputdir ! Strings for file/path and path.
1711 
1712  logical :: use_ALE ! True if ALE is being used, False if in layered mode
1713  logical :: new_sponges ! True if using the newer sponges which do not
1714  ! need to reside on the model horizontal grid.
1715 
1716  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
1717  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1718 
1719  pres(:) = 0.0 ; eta(:,:,:) = 0.0 ; tmp(:,:,:) = 0.0 ; idamp(:,:) = 0.0
1720 
1721  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
1722  inputdir = slasher(inputdir)
1723  call get_param(param_file, mdl, "SPONGE_DAMPING_FILE", damping_file, &
1724  "The name of the file with the sponge damping rates.", &
1725  fail_if_missing=.true.)
1726  call get_param(param_file, mdl, "SPONGE_STATE_FILE", state_file, &
1727  "The name of the file with the state to damp toward.", &
1728  default=damping_file)
1729  call get_param(param_file, mdl, "SPONGE_PTEMP_VAR", potemp_var, &
1730  "The name of the potential temperature variable in "//&
1731  "SPONGE_STATE_FILE.", default="PTEMP")
1732  call get_param(param_file, mdl, "SPONGE_SALT_VAR", salin_var, &
1733  "The name of the salinity variable in "//&
1734  "SPONGE_STATE_FILE.", default="SALT")
1735  call get_param(param_file, mdl, "SPONGE_ETA_VAR", eta_var, &
1736  "The name of the interface height variable in "//&
1737  "SPONGE_STATE_FILE.", default="ETA")
1738  call get_param(param_file, mdl, "SPONGE_IDAMP_VAR", idamp_var, &
1739  "The name of the inverse damping rate variable in "//&
1740  "SPONGE_DAMPING_FILE.", default="IDAMP")
1741  call get_param(param_file, mdl, "USE_REGRIDDING", use_ale, do_not_log = .true.)
1742 
1743  call get_param(param_file, mdl, "NEW_SPONGES", new_sponges, &
1744  "Set True if using the newer sponging code which "//&
1745  "performs on-the-fly regridding in lat-lon-time.",&
1746  "of sponge restoring data.", default=.false.)
1747 
1748 ! if (use_ALE) then
1749 ! call get_param(param_file, mdl, "SPONGE_RESTORE_ETA", restore_eta, &
1750 ! "If true, then restore the interface positions towards "//&
1751 ! "target values (in ALE mode)", default = .false.)
1752 ! endif
1753 
1754  filename = trim(inputdir)//trim(damping_file)
1755  call log_param(param_file, mdl, "INPUTDIR/SPONGE_DAMPING_FILE", filename)
1756  if (.not.file_exists(filename, g%Domain)) &
1757  call mom_error(fatal, " initialize_sponges: Unable to open "//trim(filename))
1758 
1759  if (new_sponges .and. .not. use_ale) &
1760  call mom_error(fatal, " initialize_sponges: Newer sponges are currently unavailable in layered mode ")
1761 
1762  call mom_read_data(filename, "Idamp", idamp(:,:), g%Domain)
1763 
1764  ! Now register all of the fields which are damped in the sponge.
1765  ! By default, momentum is advected vertically within the sponge, but
1766  ! momentum is typically not damped within the sponge.
1767 
1768  filename = trim(inputdir)//trim(state_file)
1769  call log_param(param_file, mdl, "INPUTDIR/SPONGE_STATE_FILE", filename)
1770  if (.not.file_exists(filename, g%Domain)) &
1771  call mom_error(fatal, " initialize_sponges: Unable to open "//trim(filename))
1772 
1773  ! The first call to set_up_sponge_field is for the interface heights if in layered mode.!
1774 
1775  if (.not. use_ale) then
1776  allocate(eta(isd:ied,jsd:jed,nz+1))
1777  call mom_read_data(filename, eta_var, eta(:,:,:), g%Domain, scale=us%m_to_Z)
1778 
1779  do j=js,je ; do i=is,ie
1780  eta(i,j,nz+1) = -g%bathyT(i,j)
1781  enddo ; enddo
1782  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
1783  if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z)) &
1784  eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
1785  enddo ; enddo ; enddo
1786  ! Set the inverse damping rates so that the model will know where to
1787  ! apply the sponges, along with the interface heights.
1788  call initialize_sponge(idamp, eta, g, param_file, csp, gv)
1789  deallocate(eta)
1790  elseif (.not. new_sponges) then ! ALE mode
1791 
1792  call field_size(filename,eta_var,siz,no_domain=.true.)
1793  if (siz(1) /= g%ieg-g%isg+1 .or. siz(2) /= g%jeg-g%jsg+1) &
1794  call mom_error(fatal,"initialize_sponge_file: Array size mismatch for sponge data.")
1795 
1796 ! ALE_CSp%time_dependent_target = .false.
1797 ! if (siz(4) > 1) ALE_CSp%time_dependent_target = .true.
1798  nz_data = siz(3)-1
1799  allocate(eta(isd:ied,jsd:jed,nz_data+1))
1800  allocate(h(isd:ied,jsd:jed,nz_data))
1801 
1802  call mom_read_data(filename, eta_var, eta(:,:,:), g%Domain, scale=us%m_to_Z)
1803 
1804  do j=js,je ; do i=is,ie
1805  eta(i,j,nz+1) = -g%bathyT(i,j)
1806  enddo ; enddo
1807 
1808  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
1809  if (eta(i,j,k) < (eta(i,j,k+1) + gv%Angstrom_Z)) &
1810  eta(i,j,k) = eta(i,j,k+1) + gv%Angstrom_Z
1811  enddo ; enddo ; enddo
1812  do k=1,nz; do j=js,je ; do i=is,ie
1813  h(i,j,k) = gv%Z_to_H*(eta(i,j,k)-eta(i,j,k+1))
1814  enddo ; enddo ; enddo
1815  call initialize_ale_sponge(idamp, g, param_file, ale_csp, h, nz_data)
1816  deallocate(eta)
1817  deallocate(h)
1818  else
1819  ! Initialize sponges without supplying sponge grid
1820  call initialize_ale_sponge(idamp, g, param_file, ale_csp)
1821  endif
1822 
1823  ! Now register all of the tracer fields which are damped in the
1824  ! sponge. By default, momentum is advected vertically within the
1825  ! sponge, but momentum is typically not damped within the sponge.
1826 
1827  if ( gv%nkml>0 .and. .not. new_sponges) then
1828  ! This call to set_up_sponge_ML_density registers the target values of the
1829  ! mixed layer density, which is used in determining which layers can be
1830  ! inflated without causing static instabilities.
1831  do i=is-1,ie ; pres(i) = tv%P_Ref ; enddo
1832 
1833  call mom_read_data(filename, potemp_var, tmp(:,:,:), g%Domain)
1834  call mom_read_data(filename, salin_var, tmp2(:,:,:), g%Domain)
1835 
1836  do j=js,je
1837  call calculate_density(tmp(:,j,1), tmp2(:,j,1), pres, tmp_2d(:,j), &
1838  is, ie-is+1, tv%eqn_of_state)
1839  enddo
1840 
1841  call set_up_sponge_ml_density(tmp_2d, g, csp)
1842  endif
1843 
1844  ! The remaining calls to set_up_sponge_field can be in any order.
1845  if ( use_temperature .and. .not. new_sponges) then
1846  call mom_read_data(filename, potemp_var, tmp(:,:,:), g%Domain)
1847  call set_up_sponge_field(tmp, tv%T, g, nz, csp)
1848  call mom_read_data(filename, salin_var, tmp(:,:,:), g%Domain)
1849  call set_up_sponge_field(tmp, tv%S, g, nz, csp)
1850  elseif (use_temperature) then
1851  call set_up_ale_sponge_field(filename, potemp_var, time, g, gv, tv%T, ale_csp)
1852  call set_up_ale_sponge_field(filename, salin_var, time, g, gv, tv%S, ale_csp)
1853  endif
1854 
1855 end subroutine initialize_sponges_file
1856 
1857 !> This subroutine sets the 4 bottom depths at velocity points to be the
1858 !! maximum of the adjacent depths.
1859 subroutine set_velocity_depth_max(G)
1860  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1861  ! Local variables
1862  integer :: i, j
1863 
1864  do i=g%isd,g%ied-1 ; do j=g%jsd,g%jed
1865  g%Dblock_u(i,j) = g%mask2dCu(i,j) * max(g%bathyT(i,j), g%bathyT(i+1,j))
1866  g%Dopen_u(i,j) = g%Dblock_u(i,j)
1867  enddo ; enddo
1868  do i=g%isd,g%ied ; do j=g%jsd,g%jed-1
1869  g%Dblock_v(i,j) = g%mask2dCv(i,j) * max(g%bathyT(i,j), g%bathyT(i,j+1))
1870  g%Dopen_v(i,j) = g%Dblock_v(i,j)
1871  enddo ; enddo
1872 end subroutine set_velocity_depth_max
1873 
1874 !> Subroutine to pre-compute global integrals of grid quantities for
1875 !! later use in reporting diagnostics
1876 subroutine compute_global_grid_integrals(G)
1877  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1878  ! Local variables
1879  real, dimension(G%isc:G%iec, G%jsc:G%jec) :: tmpForSumming
1880  integer :: i,j
1881 
1882  tmpforsumming(:,:) = 0.
1883  g%areaT_global = 0.0 ; g%IareaT_global = 0.0
1884  do j=g%jsc,g%jec ; do i=g%isc,g%iec
1885  tmpforsumming(i,j) = g%areaT(i,j) * g%mask2dT(i,j)
1886  enddo ; enddo
1887  g%areaT_global = reproducing_sum(tmpforsumming)
1888  g%IareaT_global = 1. / g%areaT_global
1889 end subroutine compute_global_grid_integrals
1890 
1891 !> This subroutine sets the 4 bottom depths at velocity points to be the
1892 !! minimum of the adjacent depths.
1893 subroutine set_velocity_depth_min(G)
1894  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1895  ! Local variables
1896  integer :: i, j
1897 
1898  do i=g%isd,g%ied-1 ; do j=g%jsd,g%jed
1899  g%Dblock_u(i,j) = g%mask2dCu(i,j) * min(g%bathyT(i,j), g%bathyT(i+1,j))
1900  g%Dopen_u(i,j) = g%Dblock_u(i,j)
1901  enddo ; enddo
1902  do i=g%isd,g%ied ; do j=g%jsd,g%jed-1
1903  g%Dblock_v(i,j) = g%mask2dCv(i,j) * min(g%bathyT(i,j), g%bathyT(i,j+1))
1904  g%Dopen_v(i,j) = g%Dblock_v(i,j)
1905  enddo ; enddo
1906 end subroutine set_velocity_depth_min
1907 
1908 !> This subroutine determines the isopycnal or other coordinate interfaces and
1909 !! layer potential temperatures and salinities directly from a z-space file on
1910 !! a latitude-longitude grid.
1911 subroutine mom_temp_salt_initialize_from_z(h, tv, G, GV, US, PF, just_read_params)
1912  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
1913  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
1914  intent(out) :: h !< Layer thicknesses being initialized [H ~> m or kg m-2]
1915  type(thermo_var_ptrs), intent(inout) :: tv !< A structure pointing to various thermodynamic
1916  !! variables including temperature and salinity
1917  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure
1918  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
1919  type(param_file_type), intent(in) :: PF !< A structure indicating the open file
1920  !! to parse for model parameter values.
1921  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
1922  !! only read parameters without changing h.
1923 
1924  ! Local variables
1925  character(len=200) :: filename !< The name of an input file containing temperature
1926  !! and salinity in z-space; also used for ice shelf area.
1927  character(len=200) :: tfilename !< The name of an input file containing only temperature
1928  !! in z-space.
1929  character(len=200) :: sfilename !< The name of an input file containing only salinity
1930  !! in z-space.
1931  character(len=200) :: shelf_file !< The name of an input file used for ice shelf area.
1932  character(len=200) :: inputdir !! The directory where NetCDF input filesare.
1933  character(len=200) :: mesg, area_varname, ice_shelf_file
1934 
1935  type(eos_type), pointer :: eos => null()
1936  type(thermo_var_ptrs) :: tv_loc ! A temporary thermo_var container
1937  type(verticalgrid_type) :: GV_loc ! A temporary vertical grid structure
1938  ! This include declares and sets the variable "version".
1939 # include "version_variable.h"
1940  character(len=40) :: mdl = "MOM_initialize_layers_from_Z" ! This module's name.
1941 
1942  integer :: is, ie, js, je, nz ! compute domain indices
1943  integer :: isc,iec,jsc,jec ! global compute domain indices
1944  integer :: isg, ieg, jsg, jeg ! global extent
1945  integer :: isd, ied, jsd, jed ! data domain indices
1946 
1947  integer :: i, j, k, ks, np, ni, nj
1948  integer :: idbg, jdbg
1949  integer :: nkml, nkbl ! number of mixed and buffer layers
1950 
1951  integer :: kd, inconsistent
1952  integer :: nkd ! number of levels to use for regridding input arrays
1953  real :: eps_Z ! A negligibly thin layer thickness [Z ~> m].
1954  real :: PI_180 ! for conversion from degrees to radians
1955 
1956  real, dimension(:,:), pointer :: shelf_area => null()
1957  real :: min_depth ! The minimum depth [Z ~> m].
1958  real :: dilate
1959  real :: missing_value_temp, missing_value_salt
1960  logical :: correct_thickness
1961  character(len=40) :: potemp_var, salin_var
1962  character(len=8) :: laynum
1963 
1964  integer, parameter :: niter=10 ! number of iterations for t/s adjustment to layer density
1965  logical :: just_read ! If true, just read parameters but set nothing.
1966  logical :: adjust_temperature = .true. ! fit t/s to target densities
1967  real, parameter :: missing_value = -1.e20
1968  real, parameter :: temp_land_fill = 0.0, salt_land_fill = 35.0
1969  logical :: reentrant_x, tripolar_n,dbg
1970  logical :: debug = .false. ! manually set this to true for verbose output
1971 
1972  ! data arrays
1973  real, dimension(:), allocatable :: z_edges_in, z_in, Rb
1974  real, dimension(:,:,:), allocatable, target :: temp_z, salt_z, mask_z
1975  real, dimension(:,:,:), allocatable :: rho_z
1976  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: zi ! Interface heights [Z ~> m].
1977  real, dimension(SZI_(G),SZJ_(G)) :: nlevs
1978  real, dimension(SZI_(G)) :: press ! Pressures [Pa].
1979 
1980  ! Local variables for ALE remapping
1981  real, dimension(:), allocatable :: hTarget ! Target thicknesses [Z ~> m].
1982  real, dimension(:,:), allocatable :: area_shelf_h
1983  real, dimension(:,:), allocatable, target :: frac_shelf_h
1984  real, dimension(:,:,:), allocatable, target :: tmpT1dIn, tmpS1dIn
1985  real, dimension(:,:,:), allocatable :: tmp_mask_in
1986  real, dimension(:,:,:), allocatable :: h1 ! Thicknesses [H ~> m or kg m-2].
1987  real, dimension(:,:,:), allocatable :: dz_interface ! Change in position of interface due to regridding
1988  real :: zTopOfCell, zBottomOfCell ! Heights in Z units [Z ~> m].
1989  type(regridding_cs) :: regridCS ! Regridding parameters and work arrays
1990  type(remapping_cs) :: remapCS ! Remapping parameters and work arrays
1991 
1992  logical :: homogenize, useALEremapping, remap_full_column, remap_general, remap_old_alg
1993  logical :: use_ice_shelf
1994  character(len=10) :: remappingScheme
1995  real :: tempAvg, saltAvg
1996  integer :: nPoints, ans
1997  integer :: id_clock_routine, id_clock_read, id_clock_interp, id_clock_fill, id_clock_ALE
1998 
1999  id_clock_routine = cpu_clock_id('(Initialize from Z)', grain=clock_routine)
2000  id_clock_ale = cpu_clock_id('(Initialize from Z) ALE', grain=clock_loop)
2001 
2002  call cpu_clock_begin(id_clock_routine)
2003 
2004  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
2005  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
2006  isg = g%isg ; ieg = g%ieg ; jsg = g%jsg ; jeg = g%jeg
2007 
2008  pi_180=atan(1.0)/45.
2009 
2010  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
2011 
2012  if (.not.just_read) call calltree_enter(trim(mdl)//"(), MOM_state_initialization.F90")
2013  if (.not.just_read) call log_version(pf, mdl, version, "")
2014 
2015  inputdir = "." ; call get_param(pf, mdl, "INPUTDIR", inputdir)
2016  inputdir = slasher(inputdir)
2017 
2018  eos => tv%eqn_of_state
2019 
2020 ! call mpp_get_compute_domain(G%domain%mpp_domain,isc,iec,jsc,jec)
2021 
2022  reentrant_x = .false. ; call get_param(pf, mdl, "REENTRANT_X", reentrant_x,default=.true.)
2023  tripolar_n = .false. ; call get_param(pf, mdl, "TRIPOLAR_N", tripolar_n, default=.false.)
2024  call get_param(pf, mdl, "MINIMUM_DEPTH", min_depth, default=0.0, scale=us%m_to_Z)
2025 
2026  call get_param(pf, mdl, "NKML",nkml,default=0)
2027  call get_param(pf, mdl, "NKBL",nkbl,default=0)
2028 
2029  call get_param(pf, mdl, "TEMP_SALT_Z_INIT_FILE",filename, &
2030  "The name of the z-space input file used to initialize "//&
2031  "temperatures (T) and salinities (S). If T and S are not "//&
2032  "in the same file, TEMP_Z_INIT_FILE and SALT_Z_INIT_FILE "//&
2033  "must be set.",default="temp_salt_z.nc",do_not_log=just_read)
2034  call get_param(pf, mdl, "TEMP_Z_INIT_FILE",tfilename, &
2035  "The name of the z-space input file used to initialize "//&
2036  "temperatures, only.", default=trim(filename),do_not_log=just_read)
2037  call get_param(pf, mdl, "SALT_Z_INIT_FILE",sfilename, &
2038  "The name of the z-space input file used to initialize "//&
2039  "temperatures, only.", default=trim(filename),do_not_log=just_read)
2040  filename = trim(inputdir)//trim(filename)
2041  tfilename = trim(inputdir)//trim(tfilename)
2042  sfilename = trim(inputdir)//trim(sfilename)
2043  call get_param(pf, mdl, "Z_INIT_FILE_PTEMP_VAR", potemp_var, &
2044  "The name of the potential temperature variable in "//&
2045  "TEMP_Z_INIT_FILE.", default="ptemp",do_not_log=just_read)
2046  call get_param(pf, mdl, "Z_INIT_FILE_SALT_VAR", salin_var, &
2047  "The name of the salinity variable in "//&
2048  "SALT_Z_INIT_FILE.", default="salt",do_not_log=just_read)
2049  call get_param(pf, mdl, "Z_INIT_HOMOGENIZE", homogenize, &
2050  "If True, then horizontally homogenize the interpolated "//&
2051  "initial conditions.", default=.false., do_not_log=just_read)
2052  call get_param(pf, mdl, "Z_INIT_ALE_REMAPPING", usealeremapping, &
2053  "If True, then remap straight to model coordinate from file.", &
2054  default=.false., do_not_log=just_read)
2055  call get_param(pf, mdl, "Z_INIT_REMAPPING_SCHEME", remappingscheme, &
2056  "The remapping scheme to use if using Z_INIT_ALE_REMAPPING "//&
2057  "is True.", default="PPM_IH4", do_not_log=just_read)
2058  call get_param(pf, mdl, "Z_INIT_REMAP_GENERAL", remap_general, &
2059  "If false, only initializes to z* coordinates. "//&
2060  "If true, allows initialization directly to general coordinates.",&
2061  default=.false., do_not_log=just_read)
2062  call get_param(pf, mdl, "Z_INIT_REMAP_FULL_COLUMN", remap_full_column, &
2063  "If false, only reconstructs profiles for valid data points. "//&
2064  "If true, inserts vanished layers below the valid data.", &
2065  default=remap_general, do_not_log=just_read)
2066  call get_param(pf, mdl, "Z_INIT_REMAP_OLD_ALG", remap_old_alg, &
2067  "If false, uses the preferred remapping algorithm for initialization. "//&
2068  "If true, use an older, less robust algorithm for remapping.", &
2069  default=.true., do_not_log=just_read)
2070  call get_param(pf, mdl, "ICE_SHELF", use_ice_shelf, default=.false.)
2071  if (use_ice_shelf) then
2072  call get_param(pf, mdl, "ICE_THICKNESS_FILE", ice_shelf_file, &
2073  "The file from which the ice bathymetry and area are read.", &
2074  fail_if_missing=.not.just_read, do_not_log=just_read)
2075  shelf_file = trim(inputdir)//trim(ice_shelf_file)
2076  if (.not.just_read) call log_param(pf, mdl, "INPUTDIR/THICKNESS_FILE", shelf_file)
2077  call get_param(pf, mdl, "ICE_AREA_VARNAME", area_varname, &
2078  "The name of the area variable in ICE_THICKNESS_FILE.", &
2079  fail_if_missing=.not.just_read, do_not_log=just_read)
2080  endif
2081  if (.not.usealeremapping) then
2082  call get_param(pf, mdl, "ADJUST_THICKNESS", correct_thickness, &
2083  "If true, all mass below the bottom removed if the "//&
2084  "topography is shallower than the thickness input file "//&
2085  "would indicate.", default=.false., do_not_log=just_read)
2086 
2087  call get_param(pf, mdl, "FIT_TO_TARGET_DENSITY_IC", adjust_temperature, &
2088  "If true, all the interior layers are adjusted to "//&
2089  "their target densities using mostly temperature "//&
2090  "This approach can be problematic, particularly in the "//&
2091  "high latitudes.", default=.true., do_not_log=just_read)
2092  endif
2093  if (just_read) then
2094  call cpu_clock_end(id_clock_routine)
2095  return ! All run-time parameters have been read, so return.
2096  endif
2097 
2098  !### Change this to GV%Angstrom_Z
2099  eps_z = 1.0e-10*us%m_to_Z
2100 
2101  ! Read input grid coordinates for temperature and salinity field
2102  ! in z-coordinate dataset. The file is REQUIRED to contain the
2103  ! following:
2104  !
2105  ! dimension variables:
2106  ! lon (degrees_E), lat (degrees_N), depth(meters)
2107  ! variables:
2108  ! ptemp(lon,lat,depth) : degC, potential temperature
2109  ! salt (lon,lat,depth) : ppt, salinity
2110  !
2111  ! The first record will be read if there are multiple time levels.
2112  ! The observation grid MUST tile the model grid. If the model grid extends
2113  ! to the North/South Pole past the limits of the input data, they are extrapolated using the average
2114  ! value at the northernmost/southernmost latitude.
2115 
2116  call horiz_interp_and_extrap_tracer(tfilename, potemp_var, 1.0, 1, &
2117  g, temp_z, mask_z, z_in, z_edges_in, missing_value_temp, reentrant_x, &
2118  tripolar_n, homogenize, m_to_z=us%m_to_Z)
2119 
2120  call horiz_interp_and_extrap_tracer(sfilename, salin_var, 1.0, 1, &
2121  g, salt_z, mask_z, z_in, z_edges_in, missing_value_salt, reentrant_x, &
2122  tripolar_n, homogenize, m_to_z=us%m_to_Z)
2123 
2124  kd = size(z_in,1)
2125 
2126  ! Convert the sign convention of Z_edges_in.
2127  do k=1,size(z_edges_in,1) ; z_edges_in(k) = -z_edges_in(k) ; enddo
2128 
2129  allocate(rho_z(isd:ied,jsd:jed,kd))
2130  allocate(area_shelf_h(isd:ied,jsd:jed))
2131  allocate(frac_shelf_h(isd:ied,jsd:jed))
2132 
2133  press(:) = tv%p_ref
2134 
2135  ! Convert T&S to Absolute Salinity and Conservative Temperature if using TEOS10 or NEMO
2136  call convert_temp_salt_for_teos10(temp_z, salt_z, press, g, kd, mask_z, eos)
2137 
2138  do k=1,kd ; do j=js,je
2139  call calculate_density(temp_z(:,j,k), salt_z(:,j,k), press, rho_z(:,j,k), is, ie, eos)
2140  enddo ; enddo
2141 
2142  call pass_var(temp_z,g%Domain)
2143  call pass_var(salt_z,g%Domain)
2144  call pass_var(mask_z,g%Domain)
2145  call pass_var(rho_z,g%Domain)
2146 
2147  ! This is needed for building an ALE grid under ice shelves
2148  if (use_ice_shelf) then
2149  if (.not.file_exists(shelf_file, g%Domain)) call mom_error(fatal, &
2150  "MOM_temp_salt_initialize_from_Z: Unable to open shelf file "//trim(shelf_file))
2151 
2152  call mom_read_data(shelf_file, trim(area_varname), area_shelf_h, g%Domain)
2153 
2154  ! Initialize frac_shelf_h with zeros (open water everywhere)
2155  frac_shelf_h(:,:) = 0.0
2156  ! Compute fractional ice shelf coverage of h
2157  do j=jsd,jed ; do i=isd,ied
2158  if (g%areaT(i,j) > 0.0) &
2159  frac_shelf_h(i,j) = area_shelf_h(i,j) / g%areaT(i,j)
2160  enddo ; enddo
2161  ! Pass to the pointer for use as an argument to regridding_main
2162  shelf_area => frac_shelf_h
2163 
2164  endif
2165 
2166  ! Done with horizontal interpolation.
2167  ! Now remap to model coordinates
2168  if (usealeremapping) then
2169  call cpu_clock_begin(id_clock_ale)
2170  nkd = max(gv%ke, kd)
2171  ! The regridding tools (grid generation) are coded to work on model arrays of the same
2172  ! vertical shape. We need to re-write the regridding if the model has fewer layers
2173  ! than the data. -AJA
2174 ! if (kd>nz) call MOM_error(FATAL,"MOM_initialize_state, MOM_temp_salt_initialize_from_Z(): "//&
2175 ! "Data has more levels than the model - this has not been coded yet!")
2176  ! Build the source grid and copy data onto model-shaped arrays with vanished layers
2177  allocate( tmp_mask_in(isd:ied,jsd:jed,nkd) ) ; tmp_mask_in(:,:,:) = 0.
2178  allocate( h1(isd:ied,jsd:jed,nkd) ) ; h1(:,:,:) = 0.
2179  allocate( tmpt1din(isd:ied,jsd:jed,nkd) ) ; tmpt1din(:,:,:) = 0.
2180  allocate( tmps1din(isd:ied,jsd:jed,nkd) ) ; tmps1din(:,:,:) = 0.
2181  do j = js, je ; do i = is, ie
2182  if (g%mask2dT(i,j)>0.) then
2183  ztopofcell = 0. ; zbottomofcell = 0.
2184  tmp_mask_in(i,j,1:kd) = mask_z(i,j,:)
2185  do k = 1, nkd
2186  if (tmp_mask_in(i,j,k)>0. .and. k<=kd) then
2187  zbottomofcell = max( z_edges_in(k+1), -g%bathyT(i,j) )
2188  tmpt1din(i,j,k) = temp_z(i,j,k)
2189  tmps1din(i,j,k) = salt_z(i,j,k)
2190  elseif (k>1) then
2191  zbottomofcell = -g%bathyT(i,j)
2192  tmpt1din(i,j,k) = tmpt1din(i,j,k-1)
2193  tmps1din(i,j,k) = tmps1din(i,j,k-1)
2194  else ! This next block should only ever be reached over land
2195  tmpt1din(i,j,k) = -99.9
2196  tmps1din(i,j,k) = -99.9
2197  endif
2198  h1(i,j,k) = gv%Z_to_H * (ztopofcell - zbottomofcell)
2199  ztopofcell = zbottomofcell ! Bottom becomes top for next value of k
2200  enddo
2201  h1(i,j,kd) = h1(i,j,kd) + gv%Z_to_H * max(0., ztopofcell + g%bathyT(i,j) )
2202  ! The max here is in case the data data is shallower than model
2203  endif ! mask2dT
2204  enddo ; enddo
2205  deallocate( tmp_mask_in )
2206  call pass_var(h1, g%Domain)
2207  call pass_var(tmpt1din, g%Domain)
2208  call pass_var(tmps1din, g%Domain)
2209 
2210  ! Build the target grid (and set the model thickness to it)
2211  ! This call can be more general but is hard-coded for z* coordinates... ????
2212  call ale_initregridding( gv, us, g%max_depth, pf, mdl, regridcs ) ! sets regridCS
2213 
2214  if (.not. remap_general) then
2215  ! This is the old way of initializing to z* coordinates only
2216  allocate( htarget(nz) )
2217  htarget = getcoordinateresolution( regridcs )
2218  do j = js, je ; do i = is, ie
2219  h(i,j,:) = 0.
2220  if (g%mask2dT(i,j)>0.) then
2221  ! Build the target grid combining hTarget and topography
2222  ztopofcell = 0. ; zbottomofcell = 0.
2223  do k = 1, nz
2224  zbottomofcell = max( ztopofcell - htarget(k), -g%bathyT(i,j) )
2225  h(i,j,k) = gv%Z_to_H * (ztopofcell - zbottomofcell)
2226  ztopofcell = zbottomofcell ! Bottom becomes top for next value of k
2227  enddo
2228  else
2229  h(i,j,:) = 0.
2230  endif ! mask2dT
2231  enddo ; enddo
2232  call pass_var(h, g%Domain)
2233  deallocate( htarget )
2234  endif
2235 
2236  ! Now remap from source grid to target grid
2237  call initialize_remapping( remapcs, remappingscheme, boundary_extrapolation=.false. ) ! Reconstruction parameters
2238  if (remap_general) then
2239  call set_regrid_params( regridcs, min_thickness=0. )
2240  tv_loc = tv
2241  tv_loc%T => tmpt1din
2242  tv_loc%S => tmps1din
2243  gv_loc = gv
2244  gv_loc%ke = nkd
2245  allocate( dz_interface(isd:ied,jsd:jed,nkd+1) ) ! Need for argument to regridding_main() but is not used
2246  if (use_ice_shelf) then
2247  call regridding_main( remapcs, regridcs, g, gv_loc, h1, tv_loc, h, dz_interface, shelf_area )
2248  else
2249  call regridding_main( remapcs, regridcs, g, gv_loc, h1, tv_loc, h, dz_interface )
2250  endif
2251  deallocate( dz_interface )
2252  endif
2253  call ale_remap_scalar(remapcs, g, gv, nkd, h1, tmpt1din, h, tv%T, all_cells=remap_full_column, &
2254  old_remap=remap_old_alg )
2255  call ale_remap_scalar(remapcs, g, gv, nkd, h1, tmps1din, h, tv%S, all_cells=remap_full_column, &
2256  old_remap=remap_old_alg )
2257  deallocate( h1 )
2258  deallocate( tmpt1din )
2259  deallocate( tmps1din )
2260 
2261  call cpu_clock_end(id_clock_ale)
2262 
2263  else ! remap to isopycnal layer space
2264 
2265  ! Next find interface positions using local arrays
2266  ! nlevs contains the number of valid data points in each column
2267  nlevs = sum(mask_z,dim=3)
2268 
2269  ! Rb contains the layer interface densities
2270  allocate(rb(nz+1))
2271  do k=2,nz ; rb(k)=0.5*(gv%Rlay(k-1)+gv%Rlay(k)) ; enddo
2272  rb(1) = 0.0 ; rb(nz+1) = 2.0*gv%Rlay(nz) - gv%Rlay(nz-1)
2273 
2274  zi(is:ie,js:je,:) = find_interfaces(rho_z(is:ie,js:je,:), z_in, rb, g%bathyT(is:ie,js:je), &
2275  nlevs(is:ie,js:je), nkml, nkbl, min_depth, eps_z=eps_z)
2276 
2277  if (correct_thickness) then
2278  call adjustetatofitbathymetry(g, gv, us, zi, h)
2279  else
2280  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
2281  if (zi(i,j,k) < (zi(i,j,k+1) + gv%Angstrom_Z)) then
2282  zi(i,j,k) = zi(i,j,k+1) + gv%Angstrom_Z
2283  h(i,j,k) = gv%Angstrom_H
2284  else
2285  h(i,j,k) = gv%Z_to_H * (zi(i,j,k) - zi(i,j,k+1))
2286  endif
2287  enddo ; enddo ; enddo
2288  inconsistent=0
2289  do j=js,je ; do i=is,ie
2290  if (abs(zi(i,j,nz+1) + g%bathyT(i,j)) > 1.0*us%m_to_Z) &
2291  inconsistent = inconsistent + 1
2292  enddo ; enddo
2293  call sum_across_pes(inconsistent)
2294 
2295  if ((inconsistent > 0) .and. (is_root_pe())) then
2296  write(mesg, '("Thickness initial conditions are inconsistent ",'// &
2297  '"with topography in ",I5," places.")') inconsistent
2298  call mom_error(warning, mesg)
2299  endif
2300  endif
2301 
2302  tv%T(is:ie,js:je,:) = tracer_z_init(temp_z(is:ie,js:je,:), z_edges_in, zi(is:ie,js:je,:), &
2303  nkml, nkbl, missing_value, g%mask2dT(is:ie,js:je), nz, &
2304  nlevs(is:ie,js:je),dbg,idbg,jdbg, eps_z=eps_z)
2305  tv%S(is:ie,js:je,:) = tracer_z_init(salt_z(is:ie,js:je,:), z_edges_in, zi(is:ie,js:je,:), &
2306  nkml, nkbl, missing_value, g%mask2dT(is:ie,js:je), nz, &
2307  nlevs(is:ie,js:je), eps_z=eps_z)
2308 
2309  do k=1,nz
2310  npoints = 0 ; tempavg = 0. ; saltavg = 0.
2311  do j=js,je ; do i=is,ie ; if (g%mask2dT(i,j) >= 1.0) then
2312  npoints = npoints + 1
2313  tempavg = tempavg + tv%T(i,j,k)
2314  saltavg = saltavg + tv%S(i,j,k)
2315  endif ; enddo ; enddo
2316 
2317  ! Horizontally homogenize data to produce perfectly "flat" initial conditions
2318  if (homogenize) then
2319  call sum_across_pes(npoints)
2320  call sum_across_pes(tempavg)
2321  call sum_across_pes(saltavg)
2322  if (npoints>0) then
2323  tempavg = tempavg / real(npoints)
2324  saltavg = saltavg / real(npoints)
2325  endif
2326  tv%T(:,:,k) = tempavg
2327  tv%S(:,:,k) = saltavg
2328  endif
2329  enddo
2330 
2331  endif ! useALEremapping
2332 
2333  ! Fill land values
2334  do k=1,nz ; do j=js,je ; do i=is,ie
2335  if (tv%T(i,j,k) == missing_value) then
2336  tv%T(i,j,k) = temp_land_fill
2337  tv%S(i,j,k) = salt_land_fill
2338  endif
2339  enddo ; enddo ; enddo
2340 
2341  ! Finally adjust to target density
2342  ks = max(0,nkml)+max(0,nkbl)+1
2343 
2344  if (adjust_temperature .and. .not. usealeremapping) then
2345  call determine_temperature(tv%T(is:ie,js:je,:), tv%S(is:ie,js:je,:), &
2346  gv%Rlay(1:nz), tv%p_ref, niter, missing_value, h(is:ie,js:je,:), ks, eos)
2347 
2348  endif
2349 
2350  deallocate(z_in, z_edges_in, temp_z, salt_z, mask_z)
2351  deallocate(rho_z, area_shelf_h, frac_shelf_h)
2352 
2353  call pass_var(h, g%Domain)
2354  call pass_var(tv%T, g%Domain)
2355  call pass_var(tv%S, g%Domain)
2356 
2357  call calltree_leave(trim(mdl)//'()')
2358  call cpu_clock_end(id_clock_routine)
2359 
2360 end subroutine mom_temp_salt_initialize_from_z
2361 
2362 !> Run simple unit tests
2363 subroutine mom_state_init_tests(G, GV, US, tv)
2364  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
2365  type(verticalgrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
2366  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
2367  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure.
2368  ! Local variables
2369  integer, parameter :: nk=5
2370  real, dimension(nk) :: T, T_t, T_b, S, S_t, S_b, rho, h, z
2371  real, dimension(nk+1) :: e
2372  integer :: k
2373  real :: P_tot, P_t, P_b, z_out
2374  type(remapping_cs), pointer :: remap_CS => null()
2375 
2376  do k = 1, nk
2377  h(k) = 100.
2378  enddo
2379  e(1) = 0.
2380  do k = 1, nk
2381  e(k+1) = e(k) - h(k)
2382  enddo
2383  p_tot = 0.
2384  do k = 1, nk
2385  z(k) = 0.5 * ( e(k) + e(k+1) )
2386  t_t(k) = 20.+(0./500.)*e(k)
2387  t(k) = 20.+(0./500.)*z(k)
2388  t_b(k) = 20.+(0./500.)*e(k+1)
2389  s_t(k) = 35.-(0./500.)*e(k)
2390  s(k) = 35.+(0./500.)*z(k)
2391  s_b(k) = 35.-(0./500.)*e(k+1)
2392  call calculate_density(0.5*(t_t(k)+t_b(k)), 0.5*(s_t(k)+s_b(k)), -gv%Rho0*gv%mks_g_Earth*z(k), &
2393  rho(k), tv%eqn_of_state)
2394  p_tot = p_tot + gv%mks_g_Earth * rho(k) * h(k)
2395  enddo
2396 
2397  p_t = 0.
2398  do k = 1, nk
2399  call find_depth_of_pressure_in_cell(t_t(k), t_b(k), s_t(k), s_b(k), e(k), e(k+1), &
2400  p_t, 0.5*p_tot, gv%Rho0, gv%mks_g_Earth, tv%eqn_of_state, p_b, z_out)
2401  write(0,*) k,p_t,p_b,0.5*p_tot,e(k),e(k+1),z_out
2402  p_t = p_b
2403  enddo
2404  write(0,*) p_b,p_tot
2405 
2406  write(0,*) ''
2407  write(0,*) ' ==================================================================== '
2408  write(0,*) ''
2409  write(0,*) h
2410  call cut_off_column_top(nk, tv, gv, gv%mks_g_Earth, -e(nk+1), gv%Angstrom_H, &
2411  t, t_t, t_b, s, s_t, s_b, 0.5*p_tot, h, remap_cs)
2412  write(0,*) h
2413 
2414 end subroutine mom_state_init_tests
2415 
2416 end module mom_state_initialization
dome2d_initialization
Initialization of the 2D DOME experiment with density water initialized on a coastal shelf.
Definition: DOME2d_initialization.F90:2
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
rossby_front_2d_initialization
Initial conditions for the 2D Rossby front test.
Definition: Rossby_front_2d_initialization.F90:2
mom_state_initialization
Initialization functions for state variables, u, v, h, T and S.
Definition: MOM_state_initialization.F90:2
mom_regridding::regridding_cs
Regridding control structure.
Definition: MOM_regridding.F90:45
circle_obcs_initialization
Configures the model for the "circle_obcs" experiment which tests Open Boundary Conditions radiating ...
Definition: circle_obcs_initialization.F90:3
mom_horizontal_regridding::horiz_interp_and_extrap_tracer
Extrapolate and interpolate data.
Definition: MOM_horizontal_regridding.F90:52
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
external_gwave_initialization
Initialization for the "external gravity wave wave" configuration.
Definition: external_gwave_initialization.F90:2
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_ale_sponge::initialize_ale_sponge
Ddetermine the number of points which are within sponges in this computational domain.
Definition: MOM_ALE_sponge.F90:49
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_ale_sponge::ale_sponge_cs
ALE sponge control structure.
Definition: MOM_ALE_sponge.F90:85
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
mom_remapping::remapping_cs
Container for remapping parameters.
Definition: MOM_remapping.F90:22
mom_tracer_registry
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Definition: MOM_tracer_registry.F90:5
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_ale_sponge
This module contains the routines used to apply sponge layers when using the ALE mode.
Definition: MOM_ALE_sponge.F90:11
mom_horizontal_regridding
Horizontal interpolation.
Definition: MOM_horizontal_regridding.F90:2
mom_io
This module contains I/O framework code.
Definition: MOM_io.F90:2
mom_restart::mom_restart_cs
A restart registry and the control structure for restarts.
Definition: MOM_restart.F90:72
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
mom_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
neverland_initialization
Initialization for the "Neverland" configuration.
Definition: Neverland_initialization.F90:2
mom_remapping
Provides column-wise vertical remapping functions.
Definition: MOM_remapping.F90:2
scm_cvmix_tests
Initial conditions and forcing for the single column model (SCM) CVMix test set.
Definition: SCM_CVMix_tests.F90:2
dumbbell_initialization
Configures the model for the idealized dumbbell test case.
Definition: dumbbell_initialization.F90:2
dyed_obcs_initialization
Dyed open boundary conditions.
Definition: dyed_obcs_initialization.F90:2
mom_domains::pass_vector
Do a halo update on a pair of arrays representing the two components of a vector.
Definition: MOM_domains.F90:54
lock_exchange_initialization
Initialization of the "lock exchange" experiment. lock_exchange = A 2-d density driven hydraulic exch...
Definition: lock_exchange_initialization.F90:3
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:86
dense_water_initialization
Initialization routines for the dense water formation and overflow experiment.
Definition: dense_water_initialization.F90:3
mom_ale_sponge::set_up_ale_sponge_field
Store the reference profile at h points for a variable.
Definition: MOM_ALE_sponge.F90:34
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_restart
The MOM6 facility for reading and writing restart files, and querying what has been read.
Definition: MOM_restart.F90:2
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.F90:2
mom_ale
This module contains the main regridding routines.
Definition: MOM_ALE.F90:9
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
benchmark_initialization
Initialization for the "bench mark" configuration.
Definition: benchmark_initialization.F90:2
mom_io::mom_read_vector
Read a pair of data fields representing the two components of a vector from a file.
Definition: MOM_io.F90:82
mom_io::mom_read_data
Read a data field from a file.
Definition: MOM_io.F90:74
mom_open_boundary
Controls where open boundary conditions are applied.
Definition: MOM_open_boundary.F90:2
dome_initialization
Configures the model for the "DOME" experiment. DOME = Dynamics of Overflows and Mixing Experiment.
Definition: DOME_initialization.F90:3
mom_cpu_clock
Wraps the MPP cpu clock functions.
Definition: MOM_cpu_clock.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_sponge
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
sloshing_initialization
Initialization for the "sloshing" internal waves configuration.
Definition: sloshing_initialization.F90:2
bfb_initialization
Initialization of the boundary-forced-basing configuration.
Definition: BFB_initialization.F90:2
mom_ale::ale_cs
ALE control structure.
Definition: MOM_ALE.F90:63
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_tracer_registry::tracer_registry_type
Type to carry basic tracer information.
Definition: MOM_tracer_registry.F90:122
mom_interface_heights::find_eta
Calculates the heights of sruface or all interfaces from layer thicknesses.
Definition: MOM_interface_heights.F90:21
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_open_boundary::ocean_obc_type
Open-boundary data.
Definition: MOM_open_boundary.F90:186
mom_sponge::sponge_cs
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:40
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_grid_initialize
Initializes horizontal grid.
Definition: MOM_grid_initialize.F90:2
mom_debugging
Provides checksumming functions for debugging.
Definition: MOM_debugging.F90:7
mom_io::file_exists
Indicate whether a file exists, perhaps with domain decomposition.
Definition: MOM_io.F90:68
mom_coms::reproducing_sum
Find an accurate and order-invariant sum of distributed 2d or 3d fields.
Definition: MOM_coms.F90:54
mom_regridding
Generates vertical grids as part of the ALE algorithm.
Definition: MOM_regridding.F90:2
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
dyed_channel_initialization
Initialization for the dyed_channel configuration.
Definition: dyed_channel_initialization.F90:2
seamount_initialization
Configures the model for the idealized seamount test case.
Definition: seamount_initialization.F90:2
supercritical_initialization
The "super critical" configuration.
Definition: supercritical_initialization.F90:2
phillips_initialization
Initialization for the "Phillips" channel configuration.
Definition: Phillips_initialization.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_interface_heights
Functions for calculating interface heights, including free surface height.
Definition: MOM_interface_heights.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
isomip_initialization
Configures the ISOMIP test case.
Definition: ISOMIP_initialization.F90:2
midas_vertmap
Routines for initialization callable from MOM6 or Python (MIDAS)
Definition: midas_vertmap.F90:2
baroclinic_zone_initialization
Initial conditions for an idealized baroclinic zone.
Definition: baroclinic_zone_initialization.F90:2
adjustment_initialization
Configures the model for the geostrophic adjustment test case.
Definition: adjustment_initialization.F90:2
user_initialization
A template of a user to code up customized initial conditions.
Definition: user_initialization.F90:2
mom_eos::calculate_density
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:60
mom_file_parser::read_param
An overloaded interface to read various types of parameters.
Definition: MOM_file_parser.F90:90