MOM6
MOM_shared_initialization.F90
1 !> Code that initializes fixed aspects of the model grid, such as horizontal
2 !! grid metrics, topography and Coriolis, and can be shared between components.
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
7 use mom_coms, only : max_across_pes, reproducing_sum
8 use mom_domains, only : pass_var, pass_vector, sum_across_pes, broadcast
9 use mom_domains, only : root_pe, to_all, scalar_pair, cgrid_ne, agrid
11 use mom_error_handler, only : mom_mesg, mom_error, fatal, warning, is_root_pe
12 use mom_error_handler, only : calltree_enter, calltree_leave, calltree_waypoint
14 use mom_io, only : close_file, create_file, fieldtype, file_exists
15 use mom_io, only : mom_read_data, mom_read_vector, single_file, multiple
16 use mom_io, only : slasher, vardesc, write_field, var_desc
17 use mom_string_functions, only : uppercase
19 
20 use netcdf
21 
22 implicit none ; private
23 
24 public mom_shared_init_init
25 public mom_initialize_rotation, mom_calculate_grad_coriolis
26 public initialize_topography_from_file, apply_topography_edits_from_file
27 public initialize_topography_named, limit_topography, diagnosemaximumdepth
28 public set_rotation_planetary, set_rotation_beta_plane, initialize_grid_rotation_angle
29 public reset_face_lengths_named, reset_face_lengths_file, reset_face_lengths_list
30 public read_face_length_list, set_velocity_depth_max, set_velocity_depth_min
31 public compute_global_grid_integrals, write_ocean_geometry_file
32 
33 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
34 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
35 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
36 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
37 
38 contains
39 
40 ! -----------------------------------------------------------------------------
41 !> MOM_shared_init_init just writes the code version.
42 subroutine mom_shared_init_init(PF)
43  type(param_file_type), intent(in) :: pf !< A structure indicating the open file
44  !! to parse for model parameter values.
45 
46  character(len=40) :: mdl = "MOM_shared_initialization" ! This module's name.
47 
48 ! This include declares and sets the variable "version".
49 #include "version_variable.h"
50  call log_version(pf, mdl, version, &
51  "Sharable code to initialize time-invariant fields, like bathymetry and Coriolis parameters.")
52 
53 end subroutine mom_shared_init_init
54 ! -----------------------------------------------------------------------------
55 
56 !> MOM_initialize_rotation makes the appropriate call to set up the Coriolis parameter.
57 subroutine mom_initialize_rotation(f, G, PF, US)
58  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
59  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), intent(out) :: f !< The Coriolis parameter [T-1 ~> s-1]
60  type(param_file_type), intent(in) :: pf !< Parameter file structure
61  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
62 
63 ! This subroutine makes the appropriate call to set up the Coriolis parameter.
64 ! This is a separate subroutine so that it can be made public and shared with
65 ! the ice-sheet code or other components.
66 ! Set up the Coriolis parameter, f, either analytically or from file.
67  character(len=40) :: mdl = "MOM_initialize_rotation" ! This subroutine's name.
68  character(len=200) :: config
69 
70  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
71  call get_param(pf, mdl, "ROTATION", config, &
72  "This specifies how the Coriolis parameter is specified: \n"//&
73  " \t 2omegasinlat - Use twice the planetary rotation rate \n"//&
74  " \t\t times the sine of latitude.\n"//&
75  " \t betaplane - Use a beta-plane or f-plane.\n"//&
76  " \t USER - call a user modified routine.", &
77  default="2omegasinlat")
78  select case (trim(config))
79  case ("2omegasinlat"); call set_rotation_planetary(f, g, pf, us)
80  case ("beta"); call set_rotation_beta_plane(f, g, pf, us)
81  case ("betaplane"); call set_rotation_beta_plane(f, g, pf, us)
82  !case ("nonrotating") ! Note from AJA: Missing case?
83  case default ; call mom_error(fatal,"MOM_initialize: "// &
84  "Unrecognized rotation setup "//trim(config))
85  end select
86  call calltree_leave(trim(mdl)//'()')
87 end subroutine mom_initialize_rotation
88 
89 !> Calculates the components of grad f (Coriolis parameter)
90 subroutine mom_calculate_grad_coriolis(dF_dx, dF_dy, G, US)
91  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid type
92  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
93  intent(out) :: df_dx !< x-component of grad f [T-1 m-1 ~> s-1 m-1]
94  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
95  intent(out) :: df_dy !< y-component of grad f [T-1 m-1 ~> s-1 m-1]
96  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
97  ! Local variables
98  integer :: i,j
99  real :: f1, f2
100 
101  if ((lbound(g%CoriolisBu,1) > g%isc-1) .or. &
102  (lbound(g%CoriolisBu,2) > g%isc-1)) then
103  ! The gradient of the Coriolis parameter can not be calculated with this grid.
104  df_dx(:,:) = 0.0 ; df_dy(:,:) = 0.0
105  return
106  endif
107 
108  do j=g%jsc, g%jec ; do i=g%isc, g%iec
109  f1 = 0.5*( g%CoriolisBu(i,j) + g%CoriolisBu(i,j-1) )
110  f2 = 0.5*( g%CoriolisBu(i-1,j) + g%CoriolisBu(i-1,j-1) )
111  df_dx(i,j) = g%IdxT(i,j) * ( f1 - f2 )
112  f1 = 0.5*( g%CoriolisBu(i,j) + g%CoriolisBu(i-1,j) )
113  f2 = 0.5*( g%CoriolisBu(i,j-1) + g%CoriolisBu(i-1,j-1) )
114  df_dy(i,j) = g%IdyT(i,j) * ( f1 - f2 )
115  enddo ; enddo
116  call pass_vector(df_dx, df_dy, g%Domain, stagger=agrid)
117 end subroutine mom_calculate_grad_coriolis
118 
119 !> Return the global maximum ocean bottom depth in the same units as the input depth.
120 function diagnosemaximumdepth(D, G)
121  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
122  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
123  intent(in) :: d !< Ocean bottom depth in m or Z
124  real :: diagnosemaximumdepth !< The global maximum ocean bottom depth in m or Z
125  ! Local variables
126  integer :: i,j
127  diagnosemaximumdepth = d(g%isc,g%jsc)
128  do j=g%jsc, g%jec ; do i=g%isc, g%iec
129  diagnosemaximumdepth = max(diagnosemaximumdepth,d(i,j))
130  enddo ; enddo
131  call max_across_pes(diagnosemaximumdepth)
132 end function diagnosemaximumdepth
133 
134 
135 !> Read gridded depths from file
136 subroutine initialize_topography_from_file(D, G, param_file, US)
137  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
138  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
139  intent(out) :: d !< Ocean bottom depth in m or Z if US is present
140  type(param_file_type), intent(in) :: param_file !< Parameter file structure
141  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
142  ! Local variables
143  real :: m_to_z ! A dimensional rescaling factor.
144  character(len=200) :: filename, topo_file, inputdir ! Strings for file/path
145  character(len=200) :: topo_varname ! Variable name in file
146  character(len=40) :: mdl = "initialize_topography_from_file" ! This subroutine's name.
147 
148  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
149 
150  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
151 
152  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
153  inputdir = slasher(inputdir)
154  call get_param(param_file, mdl, "TOPO_FILE", topo_file, &
155  "The file from which the bathymetry is read.", &
156  default="topog.nc")
157  call get_param(param_file, mdl, "TOPO_VARNAME", topo_varname, &
158  "The name of the bathymetry variable in TOPO_FILE.", &
159  default="depth")
160 
161  filename = trim(inputdir)//trim(topo_file)
162  call log_param(param_file, mdl, "INPUTDIR/TOPO_FILE", filename)
163 
164  if (.not.file_exists(filename, g%Domain)) call mom_error(fatal, &
165  " initialize_topography_from_file: Unable to open "//trim(filename))
166 
167  d(:,:) = -9.e30*m_to_z ! Initializing to a very large negative depth (tall mountains) everywhere
168  ! before reading from a file should do nothing. However, in the instance of
169  ! masked-out PEs, halo regions are not updated when a processor does not
170  ! exist. We need to ensure the depth in masked-out PEs appears to be that
171  ! of land so this line does that in the halo regions. For non-masked PEs
172  ! the halo region is filled properly with a later pass_var().
173  call mom_read_data(filename, trim(topo_varname), d, g%Domain, scale=m_to_z)
174 
175  call apply_topography_edits_from_file(d, g, param_file, us)
176 
177  call calltree_leave(trim(mdl)//'()')
178 end subroutine initialize_topography_from_file
179 
180 !> Applies a list of topography overrides read from a netcdf file
181 subroutine apply_topography_edits_from_file(D, G, param_file, US)
182  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
183  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
184  intent(inout) :: d !< Ocean bottom depth in m or Z if US is present
185  type(param_file_type), intent(in) :: param_file !< Parameter file structure
186  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
187 
188  ! Local variables
189  real :: m_to_z ! A dimensional rescaling factor.
190  character(len=200) :: topo_edits_file, inputdir ! Strings for file/path
191  character(len=40) :: mdl = "apply_topography_edits_from_file" ! This subroutine's name.
192  integer :: n_edits, n, ashape(5), i, j, ncid, id, ncstatus, iid, jid, zid
193  integer, dimension(:), allocatable :: ig, jg
194  real, dimension(:), allocatable :: new_depth
195 
196  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
197 
198  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
199 
200  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
201  inputdir = slasher(inputdir)
202  call get_param(param_file, mdl, "TOPO_EDITS_FILE", topo_edits_file, &
203  "The file from which to read a list of i,j,z topography overrides.", &
204  default="")
205 
206  if (len_trim(topo_edits_file)==0) return
207 
208  topo_edits_file = trim(inputdir)//trim(topo_edits_file)
209  if (.not.file_exists(topo_edits_file, g%Domain)) call mom_error(fatal, &
210  'initialize_topography_from_file: Unable to open '//trim(topo_edits_file))
211 
212  ncstatus = nf90_open(trim(topo_edits_file), nf90_nowrite, ncid)
213  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
214  'Failed to open '//trim(topo_edits_file))
215 
216  ! Get nEdits
217  ncstatus = nf90_inq_dimid(ncid, 'nEdits', id)
218  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
219  'Failed to inq_dimid nEdits for '//trim(topo_edits_file))
220  ncstatus = nf90_inquire_dimension(ncid, id, len=n_edits)
221  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
222  'Failed to inquire_dimension nEdits for '//trim(topo_edits_file))
223 
224  ! Read ni
225  ncstatus = nf90_inq_varid(ncid, 'ni', id)
226  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
227  'Failed to inq_varid ni for '//trim(topo_edits_file))
228  ncstatus = nf90_get_var(ncid, id, i)
229  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
230  'Failed to get_var ni for '//trim(topo_edits_file))
231  if (i /= g%ieg) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
232  'Incompatible i-dimension of grid in '//trim(topo_edits_file))
233 
234  ! Read nj
235  ncstatus = nf90_inq_varid(ncid, 'nj', id)
236  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
237  'Failed to inq_varid nj for '//trim(topo_edits_file))
238  ncstatus = nf90_get_var(ncid, id, j)
239  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
240  'Failed to get_var nj for '//trim(topo_edits_file))
241  if (j /= g%jeg) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
242  'Incompatible j-dimension of grid in '//trim(topo_edits_file))
243 
244  ! Read iEdit
245  ncstatus = nf90_inq_varid(ncid, 'iEdit', id)
246  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
247  'Failed to inq_varid iEdit for '//trim(topo_edits_file))
248  allocate(ig(n_edits))
249  ncstatus = nf90_get_var(ncid, id, ig)
250  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
251  'Failed to get_var iEdit for '//trim(topo_edits_file))
252 
253  ! Read jEdit
254  ncstatus = nf90_inq_varid(ncid, 'jEdit', id)
255  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
256  'Failed to inq_varid jEdit for '//trim(topo_edits_file))
257  allocate(jg(n_edits))
258  ncstatus = nf90_get_var(ncid, id, jg)
259  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
260  'Failed to get_var jEdit for '//trim(topo_edits_file))
261 
262  ! Read zEdit
263  ncstatus = nf90_inq_varid(ncid, 'zEdit', id)
264  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
265  'Failed to inq_varid zEdit for '//trim(topo_edits_file))
266  allocate(new_depth(n_edits))
267  ncstatus = nf90_get_var(ncid, id, new_depth)
268  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
269  'Failed to get_var zEdit for '//trim(topo_edits_file))
270 
271  ! Close file
272  ncstatus = nf90_close(ncid)
273  if (ncstatus /= nf90_noerr) call mom_error(fatal, 'apply_topography_edits_from_file: '//&
274  'Failed to close '//trim(topo_edits_file))
275 
276  do n = 1, n_edits
277  i = ig(n) - g%isd_global + 2 ! +1 for python indexing and +1 for ig-isd_global+1
278  j = jg(n) - g%jsd_global + 2
279  if (i>=g%isc .and. i<=g%iec .and. j>=g%jsc .and. j<=g%jec) then
280  if (new_depth(n)/=0.) then
281  write(*,'(a,3i5,f8.2,a,f8.2,2i4)') &
282  'Ocean topography edit: ',n,ig(n),jg(n),d(i,j)/m_to_z,'->',abs(new_depth(n)),i,j
283  d(i,j) = abs(m_to_z*new_depth(n)) ! Allows for height-file edits (i.e. converts negatives)
284  else
285  call mom_error(fatal, ' apply_topography_edits_from_file: '//&
286  "A zero depth edit would change the land mask and is not allowed in"//trim(topo_edits_file))
287  endif
288  endif
289  enddo
290 
291  deallocate( ig, jg, new_depth )
292 
293  call calltree_leave(trim(mdl)//'()')
294 end subroutine apply_topography_edits_from_file
295 
296 !> initialize the bathymetry based on one of several named idealized configurations
297 subroutine initialize_topography_named(D, G, param_file, topog_config, max_depth, US)
298  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
299  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
300  intent(out) :: d !< Ocean bottom depth in m or Z if US is present
301  type(param_file_type), intent(in) :: param_file !< Parameter file structure
302  character(len=*), intent(in) :: topog_config !< The name of an idealized
303  !! topographic configuration
304  real, intent(in) :: max_depth !< Maximum depth of model in the units of D
305  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
306 
307  ! This subroutine places the bottom depth in m into D(:,:), shaped according to the named config.
308 
309  ! Local variables
310  real :: m_to_z ! A dimensional rescaling factor.
311  real :: min_depth ! The minimum depth [Z ~> m].
312  real :: pi ! 3.1415926... calculated as 4*atan(1)
313  real :: d0 ! A constant to make the maximum basin depth MAXIMUM_DEPTH.
314  real :: expdecay ! A decay scale of associated with the sloping boundaries [m].
315  real :: dedge ! The depth [Z ~> m], at the basin edge
316 ! real :: south_lat, west_lon, len_lon, len_lat, Rad_earth
317  integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
318  character(len=40) :: mdl = "initialize_topography_named" ! This subroutine's name.
319  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
320  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
321 
322  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
323  call mom_mesg(" MOM_shared_initialization.F90, initialize_topography_named: "//&
324  "TOPO_CONFIG = "//trim(topog_config), 5)
325 
326  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
327 
328  call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, &
329  "The minimum depth of the ocean.", units="m", default=0.0, scale=m_to_z)
330  if (max_depth<=0.) call mom_error(fatal,"initialize_topography_named: "// &
331  "MAXIMUM_DEPTH has a non-sensical value! Was it set?")
332 
333  if (trim(topog_config) /= "flat") then
334  call get_param(param_file, mdl, "EDGE_DEPTH", dedge, &
335  "The depth at the edge of one of the named topographies.", &
336  units="m", default=100.0, scale=m_to_z)
337 ! call get_param(param_file, mdl, "SOUTHLAT", south_lat, &
338 ! "The southern latitude of the domain.", units="degrees", &
339 ! fail_if_missing=.true.)
340 ! call get_param(param_file, mdl, "LENLAT", len_lat, &
341 ! "The latitudinal length of the domain.", units="degrees", &
342 ! fail_if_missing=.true.)
343 ! call get_param(param_file, mdl, "WESTLON", west_lon, &
344 ! "The western longitude of the domain.", units="degrees", &
345 ! default=0.0)
346 ! call get_param(param_file, mdl, "LENLON", len_lon, &
347 ! "The longitudinal length of the domain.", units="degrees", &
348 ! fail_if_missing=.true.)
349 ! call get_param(param_file, mdl, "RAD_EARTH", Rad_Earth, &
350 ! "The radius of the Earth.", units="m", default=6.378e6)
351  call get_param(param_file, mdl, "TOPOG_SLOPE_SCALE", expdecay, &
352  "The exponential decay scale used in defining some of "//&
353  "the named topographies.", units="m", default=400000.0)
354  endif
355 
356 
357  pi = 4.0*atan(1.0)
358 
359  if (trim(topog_config) == "flat") then
360  do i=is,ie ; do j=js,je ; d(i,j) = max_depth ; enddo ; enddo
361  elseif (trim(topog_config) == "spoon") then
362  d0 = (max_depth - dedge) / &
363  ((1.0 - exp(-0.5*g%len_lat*g%Rad_earth*pi/(180.0 *expdecay))) * &
364  (1.0 - exp(-0.5*g%len_lat*g%Rad_earth*pi/(180.0 *expdecay))))
365  do i=is,ie ; do j=js,je
366  ! This sets a bowl shaped (sort of) bottom topography, with a !
367  ! maximum depth of max_depth. !
368  d(i,j) = dedge + d0 * &
369  (sin(pi * (g%geoLonT(i,j) - (g%west_lon)) / g%len_lon) * &
370  (1.0 - exp((g%geoLatT(i,j) - (g%south_lat+g%len_lat))*g%Rad_earth*pi / &
371  (180.0*expdecay)) ))
372  enddo ; enddo
373  elseif (trim(topog_config) == "bowl") then
374  d0 = (max_depth - dedge) / &
375  ((1.0 - exp(-0.5*g%len_lat*g%Rad_earth*pi/(180.0 *expdecay))) * &
376  (1.0 - exp(-0.5*g%len_lat*g%Rad_earth*pi/(180.0 *expdecay))))
377 
378  ! This sets a bowl shaped (sort of) bottom topography, with a
379  ! maximum depth of max_depth.
380  do i=is,ie ; do j=js,je
381  d(i,j) = dedge + d0 * &
382  (sin(pi * (g%geoLonT(i,j) - g%west_lon) / g%len_lon) * &
383  ((1.0 - exp(-(g%geoLatT(i,j) - g%south_lat)*g%Rad_Earth*pi/ &
384  (180.0*expdecay))) * &
385  (1.0 - exp((g%geoLatT(i,j) - (g%south_lat+g%len_lat))* &
386  g%Rad_Earth*pi/(180.0*expdecay)))))
387  enddo ; enddo
388  elseif (trim(topog_config) == "halfpipe") then
389  d0 = max_depth - dedge
390  do i=is,ie ; do j=js,je
391  d(i,j) = dedge + d0 * abs(sin(pi*(g%geoLatT(i,j) - g%south_lat)/g%len_lat))
392  enddo ; enddo
393  else
394  call mom_error(fatal,"initialize_topography_named: "// &
395  "Unrecognized topography name "//trim(topog_config))
396  endif
397 
398  ! This is here just for safety. Hopefully it doesn't do anything.
399  do i=is,ie ; do j=js,je
400  if (d(i,j) > max_depth) d(i,j) = max_depth
401  if (d(i,j) < min_depth) d(i,j) = 0.5*min_depth
402  enddo ; enddo
403 
404  call calltree_leave(trim(mdl)//'()')
405 end subroutine initialize_topography_named
406 ! -----------------------------------------------------------------------------
407 
408 ! -----------------------------------------------------------------------------
409 !> limit_topography ensures that min_depth < D(x,y) < max_depth
410 subroutine limit_topography(D, G, param_file, max_depth, US)
411  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
412  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
413  intent(inout) :: d !< Ocean bottom depth in m or Z if US is present
414  type(param_file_type), intent(in) :: param_file !< Parameter file structure
415  real, intent(in) :: max_depth !< Maximum depth of model in the units of D
416  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
417 
418  ! Local variables
419  real :: m_to_z ! A dimensional rescaling factor.
420  integer :: i, j
421  character(len=40) :: mdl = "limit_topography" ! This subroutine's name.
422  real :: min_depth, mask_depth
423 
424  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
425 
426  m_to_z = 1.0 ; if (present(us)) m_to_z = us%m_to_Z
427 
428  call get_param(param_file, mdl, "MINIMUM_DEPTH", min_depth, &
429  "If MASKING_DEPTH is unspecified, then anything shallower than "//&
430  "MINIMUM_DEPTH is assumed to be land and all fluxes are masked out. "//&
431  "If MASKING_DEPTH is specified, then all depths shallower than "//&
432  "MINIMUM_DEPTH but deeper than MASKING_DEPTH are rounded to MINIMUM_DEPTH.", &
433  units="m", default=0.0, scale=m_to_z)
434  call get_param(param_file, mdl, "MASKING_DEPTH", mask_depth, &
435  "The depth below which to mask the ocean as land.", &
436  units="m", default=-9999.0, scale=m_to_z, do_not_log=.true.)
437 
438 ! Make sure that min_depth < D(x,y) < max_depth
439  if (mask_depth < -9990.*m_to_z) then
440  do j=g%jsd,g%jed ; do i=g%isd,g%ied
441  d(i,j) = min( max( d(i,j), 0.5*min_depth ), max_depth )
442  enddo ; enddo
443  else
444  do j=g%jsd,g%jed ; do i=g%isd,g%ied
445  if (d(i,j)>0.) then
446  d(i,j) = min( max( d(i,j), min_depth ), max_depth )
447  else
448  d(i,j) = 0.
449  endif
450  enddo ; enddo
451  endif
452 
453  call calltree_leave(trim(mdl)//'()')
454 end subroutine limit_topography
455 ! -----------------------------------------------------------------------------
456 
457 ! -----------------------------------------------------------------------------
458 !> This subroutine sets up the Coriolis parameter for a sphere
459 subroutine set_rotation_planetary(f, G, param_file, US)
460  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid
461  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
462  intent(out) :: f !< Coriolis parameter (vertical component) [T-1 ~> s-1]
463  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
464  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
465 
466 ! This subroutine sets up the Coriolis parameter for a sphere
467  character(len=30) :: mdl = "set_rotation_planetary" ! This subroutine's name.
468  integer :: i, j
469  real :: pi
470  real :: omega ! The planetary rotation rate [T-1 ~> s-1]
471  real :: t_to_s ! A time unit conversion factor
472 
473  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
474 
475  t_to_s = 1.0 ; if (present(us)) t_to_s = us%T_to_s
476 
477  call get_param(param_file, "set_rotation_planetary", "OMEGA", omega, &
478  "The rotation rate of the earth.", units="s-1", &
479  default=7.2921e-5, scale=t_to_s)
480  pi = 4.0*atan(1.0)
481 
482  do i=g%IsdB,g%IedB ; do j=g%JsdB,g%JedB
483  f(i,j) = ( 2.0 * omega ) * sin( ( pi * g%geoLatBu(i,j) ) / 180.)
484  enddo ; enddo
485 
486  call calltree_leave(trim(mdl)//'()')
487 end subroutine set_rotation_planetary
488 ! -----------------------------------------------------------------------------
489 
490 ! -----------------------------------------------------------------------------
491 !> This subroutine sets up the Coriolis parameter for a beta-plane or f-plane
492 subroutine set_rotation_beta_plane(f, G, param_file, US)
493  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid
494  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
495  intent(out) :: f !< Coriolis parameter (vertical component) [T-1 ~> s-1]
496  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
497  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
498 
499 ! This subroutine sets up the Coriolis parameter for a beta-plane
500  integer :: i, j
501  real :: f_0 ! The reference value of the Coriolis parameter [T-1 ~> s-1]
502  real :: beta ! The meridional gradient of the Coriolis parameter [T-1 m-1 ~> s-1 m-1]
503  real :: y_scl, rad_earth
504  real :: t_to_s ! A time unit conversion factor
505  real :: pi
506  character(len=40) :: mdl = "set_rotation_beta_plane" ! This subroutine's name.
507  character(len=200) :: axis_units
508 
509  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
510 
511  t_to_s = 1.0 ; if (present(us)) t_to_s = us%T_to_s
512 
513  call get_param(param_file, mdl, "F_0", f_0, &
514  "The reference value of the Coriolis parameter with the "//&
515  "betaplane option.", units="s-1", default=0.0, scale=t_to_s)
516  call get_param(param_file, mdl, "BETA", beta, &
517  "The northward gradient of the Coriolis parameter with "//&
518  "the betaplane option.", units="m-1 s-1", default=0.0, scale=t_to_s)
519  call get_param(param_file, mdl, "AXIS_UNITS", axis_units, default="degrees")
520 
521  pi = 4.0*atan(1.0)
522  select case (axis_units(1:1))
523  case ("d")
524  call get_param(param_file, mdl, "RAD_EARTH", rad_earth, &
525  "The radius of the Earth.", units="m", default=6.378e6)
526  y_scl = rad_earth/pi
527  case ("k"); y_scl = 1.e3
528  case ("m"); y_scl = 1.
529  case ("c"); y_scl = 1.e-2
530  case default ; call mom_error(fatal, &
531  " set_rotation_beta_plane: unknown AXIS_UNITS = "//trim(axis_units))
532  end select
533 
534  do i=g%IsdB,g%IedB ; do j=g%JsdB,g%JedB
535  f(i,j) = f_0 + beta * ( g%geoLatBu(i,j) * y_scl )
536  enddo ; enddo
537 
538  call calltree_leave(trim(mdl)//'()')
539 end subroutine set_rotation_beta_plane
540 
541 !> initialize_grid_rotation_angle initializes the arrays with the sine and
542 !! cosine of the angle between logical north on the grid and true north.
543 subroutine initialize_grid_rotation_angle(G, PF)
544  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
545  type(param_file_type), intent(in) :: pf !< A structure indicating the open file
546  !! to parse for model parameter values.
547 
548  real :: angle, lon_scale
549  real :: len_lon ! The periodic range of longitudes, usually 360 degrees.
550  real :: pi_720deg ! One quarter the conversion factor from degrees to radians.
551  real :: lonb(2,2) ! The longitude of a point, shifted to have about the same value.
552  character(len=40) :: mdl = "initialize_grid_rotation_angle" ! This subroutine's name.
553  logical :: use_bugs
554  integer :: i, j, m, n
555 
556  call get_param(pf, mdl, "GRID_ROTATION_ANGLE_BUGS", use_bugs, &
557  "If true, use an older algorithm to calculate the sine and "//&
558  "cosines needed rotate between grid-oriented directions and "//&
559  "true north and east. Differences arise at the tripolar fold.", &
560  default=.true.)
561 
562  if (use_bugs) then
563  do j=g%jsc,g%jec ; do i=g%isc,g%iec
564  lon_scale = cos((g%geoLatBu(i-1,j-1) + g%geoLatBu(i,j-1 ) + &
565  g%geoLatBu(i-1,j) + g%geoLatBu(i,j)) * atan(1.0)/180)
566  angle = atan2((g%geoLonBu(i-1,j) + g%geoLonBu(i,j) - &
567  g%geoLonBu(i-1,j-1) - g%geoLonBu(i,j-1))*lon_scale, &
568  g%geoLatBu(i-1,j) + g%geoLatBu(i,j) - &
569  g%geoLatBu(i-1,j-1) - g%geoLatBu(i,j-1) )
570  g%sin_rot(i,j) = sin(angle) ! angle is the clockwise angle from lat/lon to ocean
571  g%cos_rot(i,j) = cos(angle) ! grid (e.g. angle of ocean "north" from true north)
572  enddo ; enddo
573 
574  ! This is not right at a tripolar or cubed-sphere fold.
575  call pass_var(g%cos_rot, g%Domain)
576  call pass_var(g%sin_rot, g%Domain)
577  else
578  pi_720deg = atan(1.0) / 180.0
579  len_lon = 360.0 ; if (g%len_lon > 0.0) len_lon = g%len_lon
580  do j=g%jsc,g%jec ; do i=g%isc,g%iec
581  do n=1,2 ; do m=1,2
582  lonb(m,n) = modulo_around_point(g%geoLonBu(i+m-2,j+n-2), g%geoLonT(i,j), len_lon)
583  enddo ; enddo
584  lon_scale = cos(pi_720deg*((g%geoLatBu(i-1,j-1) + g%geoLatBu(i,j)) + &
585  (g%geoLatBu(i,j-1) + g%geoLatBu(i-1,j)) ) )
586  angle = atan2(lon_scale*((lonb(1,2) - lonb(2,1)) + (lonb(2,2) - lonb(1,1))), &
587  (g%geoLatBu(i-1,j) - g%geoLatBu(i,j-1)) + &
588  (g%geoLatBu(i,j) - g%geoLatBu(i-1,j-1)) )
589  g%sin_rot(i,j) = sin(angle) ! angle is the clockwise angle from lat/lon to ocean
590  g%cos_rot(i,j) = cos(angle) ! grid (e.g. angle of ocean "north" from true north)
591  enddo ; enddo
592 
593  call pass_vector(g%cos_rot, g%sin_rot, g%Domain, stagger=agrid)
594  endif
595 
596 end subroutine initialize_grid_rotation_angle
597 
598 ! -----------------------------------------------------------------------------
599 !> Return the modulo value of x in an interval [xc-(Lx/2) xc+(Lx/2)]
600 !! If Lx<=0, then it returns x without applying modulo arithmetic.
601 function modulo_around_point(x, xc, Lx) result(x_mod)
602  real, intent(in) :: x !< Value to which to apply modulo arithmetic
603  real, intent(in) :: xc !< Center of modulo range
604  real, intent(in) :: lx !< Modulo range width
605  real :: x_mod !< x shifted by an integer multiple of Lx to be close to xc.
606 
607  if (lx > 0.0) then
608  x_mod = modulo(x - (xc - 0.5*lx), lx) + (xc - 0.5*lx)
609  else
610  x_mod = x
611  endif
612 end function modulo_around_point
613 
614 ! -----------------------------------------------------------------------------
615 !> This subroutine sets the open face lengths at selected points to restrict
616 !! passages to their observed widths based on a named set of sizes.
617 subroutine reset_face_lengths_named(G, param_file, name, US)
618  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
619  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
620  character(len=*), intent(in) :: name !< The name for the set of face lengths. Only "global_1deg"
621  !! is currently implemented.
622  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
623 
624  ! Local variables
625  character(len=256) :: mesg ! Message for error messages.
626  real :: dx_2 = -1.0, dy_2 = -1.0
627  real :: pi_180
628  integer :: option = -1
629  integer :: i, j, isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
630  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
631  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
632  pi_180 = (4.0*atan(1.0))/180.0
633 
634  select case ( trim(name) )
635  case ("global_1deg") ; option = 1 ; dx_2 = 0.5*1.0
636  case default ; call mom_error(fatal, "reset_face_lengths_named: "//&
637  "Unrecognized channel configuration name "//trim(name))
638  end select
639 
640  if (option==1) then ! 1-degree settings.
641  do j=jsd,jed ; do i=isdb,iedb ! Change any u-face lengths within this loop.
642  dy_2 = dx_2 * g%dyCu(i,j)*g%IdxCu(i,j) * cos(pi_180 * g%geoLatCu(i,j))
643 
644  if ((abs(g%geoLatCu(i,j)-35.5) < dy_2) .and. (g%geoLonCu(i,j) < -4.5) .and. &
645  (g%geoLonCu(i,j) > -6.5)) &
646  g%dy_Cu(i,j) = g%mask2dCu(i,j)*12000.0 ! Gibraltar
647 
648  if ((abs(g%geoLatCu(i,j)-12.5) < dy_2) .and. (abs(g%geoLonCu(i,j)-43.0) < dx_2)) &
649  g%dy_Cu(i,j) = g%mask2dCu(i,j)*10000.0 ! Red Sea
650 
651  if ((abs(g%geoLatCu(i,j)-40.5) < dy_2) .and. (abs(g%geoLonCu(i,j)-26.0) < dx_2)) &
652  g%dy_Cu(i,j) = g%mask2dCu(i,j)*5000.0 ! Dardanelles
653 
654  if ((abs(g%geoLatCu(i,j)-41.5) < dy_2) .and. (abs(g%geoLonCu(i,j)+220.0) < dx_2)) &
655  g%dy_Cu(i,j) = g%mask2dCu(i,j)*35000.0 ! Tsugaru strait at 140.0e
656 
657  if ((abs(g%geoLatCu(i,j)-45.5) < dy_2) .and. (abs(g%geoLonCu(i,j)+217.5) < 0.9)) &
658  g%dy_Cu(i,j) = g%mask2dCu(i,j)*15000.0 ! Betw Hokkaido and Sakhalin at 217&218 = 142e
659 
660 
661  ! Greater care needs to be taken in the tripolar region.
662  if ((abs(g%geoLatCu(i,j)-80.84) < 0.2) .and. (abs(g%geoLonCu(i,j)+64.9) < 0.8)) &
663  g%dy_Cu(i,j) = g%mask2dCu(i,j)*38000.0 ! Smith Sound in Canadian Arch - tripolar region
664 
665  enddo ; enddo
666 
667  do j=jsdb,jedb ; do i=isd,ied ! Change any v-face lengths within this loop.
668  dy_2 = dx_2 * g%dyCv(i,j)*g%IdxCv(i,j) * cos(pi_180 * g%geoLatCv(i,j))
669  if ((abs(g%geoLatCv(i,j)-41.0) < dy_2) .and. (abs(g%geoLonCv(i,j)-28.5) < dx_2)) &
670  g%dx_Cv(i,j) = g%mask2dCv(i,j)*2500.0 ! Bosporus - should be 1000.0 m wide.
671 
672  if ((abs(g%geoLatCv(i,j)-13.0) < dy_2) .and. (abs(g%geoLonCv(i,j)-42.5) < dx_2)) &
673  g%dx_Cv(i,j) = g%mask2dCv(i,j)*10000.0 ! Red Sea
674 
675  if ((abs(g%geoLatCv(i,j)+2.8) < 0.8) .and. (abs(g%geoLonCv(i,j)+241.5) < dx_2)) &
676  g%dx_Cv(i,j) = g%mask2dCv(i,j)*40000.0 ! Makassar Straits at 241.5 W = 118.5 E
677 
678  if ((abs(g%geoLatCv(i,j)-0.56) < 0.5) .and. (abs(g%geoLonCv(i,j)+240.5) < dx_2)) &
679  g%dx_Cv(i,j) = g%mask2dCv(i,j)*80000.0 ! entry to Makassar Straits at 240.5 W = 119.5 E
680 
681  if ((abs(g%geoLatCv(i,j)-0.19) < 0.5) .and. (abs(g%geoLonCv(i,j)+230.5) < dx_2)) &
682  g%dx_Cv(i,j) = g%mask2dCv(i,j)*25000.0 ! Channel betw N Guinea and Halmahara 230.5 W = 129.5 E
683 
684  if ((abs(g%geoLatCv(i,j)-0.19) < 0.5) .and. (abs(g%geoLonCv(i,j)+229.5) < dx_2)) &
685  g%dx_Cv(i,j) = g%mask2dCv(i,j)*25000.0 ! Channel betw N Guinea and Halmahara 229.5 W = 130.5 E
686 
687  if ((abs(g%geoLatCv(i,j)-0.0) < 0.25) .and. (abs(g%geoLonCv(i,j)+228.5) < dx_2)) &
688  g%dx_Cv(i,j) = g%mask2dCv(i,j)*25000.0 ! Channel betw N Guinea and Halmahara 228.5 W = 131.5 E
689 
690  if ((abs(g%geoLatCv(i,j)+8.5) < 0.5) .and. (abs(g%geoLonCv(i,j)+244.5) < dx_2)) &
691  g%dx_Cv(i,j) = g%mask2dCv(i,j)*20000.0 ! Lombok Straits at 244.5 W = 115.5 E
692 
693  if ((abs(g%geoLatCv(i,j)+8.5) < 0.5) .and. (abs(g%geoLonCv(i,j)+235.5) < dx_2)) &
694  g%dx_Cv(i,j) = g%mask2dCv(i,j)*20000.0 ! Timor Straits at 235.5 W = 124.5 E
695 
696  if ((abs(g%geoLatCv(i,j)-52.5) < dy_2) .and. (abs(g%geoLonCv(i,j)+218.5) < dx_2)) &
697  g%dx_Cv(i,j) = g%mask2dCv(i,j)*2500.0 ! Russia and Sakhalin Straits at 218.5 W = 141.5 E
698 
699  ! Greater care needs to be taken in the tripolar region.
700  if ((abs(g%geoLatCv(i,j)-76.8) < 0.06) .and. (abs(g%geoLonCv(i,j)+88.7) < dx_2)) &
701  g%dx_Cv(i,j) = g%mask2dCv(i,j)*8400.0 ! Jones Sound in Canadian Arch - tripolar region
702 
703  enddo ; enddo
704  endif
705 
706  ! These checks apply regardless of the chosen option.
707 
708  do j=jsd,jed ; do i=isdb,iedb
709  if (g%dy_Cu(i,j) > g%dyCu(i,j)) then
710  write(mesg,'("dy_Cu of ",ES11.4," exceeds unrestricted width of ",ES11.4,&
711  &" by ",ES11.4," at lon/lat of ", ES11.4, ES11.4)') &
712  g%dy_Cu(i,j), g%dyCu(i,j), g%dy_Cu(i,j)-g%dyCu(i,j), &
713  g%geoLonCu(i,j), g%geoLatCu(i,j)
714  call mom_error(fatal,"reset_face_lengths_named "//mesg)
715  endif
716  g%areaCu(i,j) = g%dxCu(i,j)*g%dy_Cu(i,j)
717  g%IareaCu(i,j) = 0.0
718  if (g%areaCu(i,j) > 0.0) g%IareaCu(i,j) = g%mask2dCu(i,j) / g%areaCu(i,j)
719  enddo ; enddo
720 
721  do j=jsdb,jedb ; do i=isd,ied
722  if (g%dx_Cv(i,j) > g%dxCv(i,j)) then
723  write(mesg,'("dx_Cv of ",ES11.4," exceeds unrestricted width of ",ES11.4,&
724  &" by ",ES11.4, " at lon/lat of ", ES11.4, ES11.4)') &
725  g%dx_Cv(i,j), g%dxCv(i,j), g%dx_Cv(i,j)-g%dxCv(i,j), &
726  g%geoLonCv(i,j), g%geoLatCv(i,j)
727 
728  call mom_error(fatal,"reset_face_lengths_named "//mesg)
729  endif
730  g%areaCv(i,j) = g%dyCv(i,j)*g%dx_Cv(i,j)
731  g%IareaCv(i,j) = 0.0
732  if (g%areaCv(i,j) > 0.0) g%IareaCv(i,j) = g%mask2dCv(i,j) / g%areaCv(i,j)
733  enddo ; enddo
734 
735 end subroutine reset_face_lengths_named
736 ! -----------------------------------------------------------------------------
737 
738 ! -----------------------------------------------------------------------------
739 !> This subroutine sets the open face lengths at selected points to restrict
740 !! passages to their observed widths from a arrays read from a file.
741 subroutine reset_face_lengths_file(G, param_file, US)
742  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
743  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
744  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
745 
746  ! Local variables
747  character(len=40) :: mdl = "reset_face_lengths_file" ! This subroutine's name.
748  character(len=256) :: mesg ! Message for error messages.
749  character(len=200) :: filename, chan_file, inputdir ! Strings for file/path
750  integer :: i, j, isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
751  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
752  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
753  ! These checks apply regardless of the chosen option.
754 
755  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
756 
757  call get_param(param_file, mdl, "CHANNEL_WIDTH_FILE", chan_file, &
758  "The file from which the list of narrowed channels is read.", &
759  default="ocean_geometry.nc")
760  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
761  inputdir = slasher(inputdir)
762  filename = trim(inputdir)//trim(chan_file)
763  call log_param(param_file, mdl, "INPUTDIR/CHANNEL_WIDTH_FILE", filename)
764 
765  if (is_root_pe()) then ; if (.not.file_exists(filename)) &
766  call mom_error(fatal," reset_face_lengths_file: Unable to open "//&
767  trim(filename))
768  endif
769 
770  call mom_read_vector(filename, "dyCuo", "dxCvo", g%dy_Cu, g%dx_Cv, g%Domain)
771  call pass_vector(g%dy_Cu, g%dx_Cv, g%Domain, to_all+scalar_pair, cgrid_ne)
772 
773  do j=jsd,jed ; do i=isdb,iedb
774  if (g%dy_Cu(i,j) > g%dyCu(i,j)) then
775  write(mesg,'("dy_Cu of ",ES11.4," exceeds unrestricted width of ",ES11.4,&
776  &" by ",ES11.4," at lon/lat of ", ES11.4, ES11.4)') &
777  g%dy_Cu(i,j), g%dyCu(i,j), g%dy_Cu(i,j)-g%dyCu(i,j), &
778  g%geoLonCu(i,j), g%geoLatCu(i,j)
779  call mom_error(fatal,"reset_face_lengths_file "//mesg)
780  endif
781  g%areaCu(i,j) = g%dxCu(i,j)*g%dy_Cu(i,j)
782  g%IareaCu(i,j) = 0.0
783  if (g%areaCu(i,j) > 0.0) g%IareaCu(i,j) = g%mask2dCu(i,j) / g%areaCu(i,j)
784  enddo ; enddo
785 
786  do j=jsdb,jedb ; do i=isd,ied
787  if (g%dx_Cv(i,j) > g%dxCv(i,j)) then
788  write(mesg,'("dx_Cv of ",ES11.4," exceeds unrestricted width of ",ES11.4,&
789  &" by ",ES11.4, " at lon/lat of ", ES11.4, ES11.4)') &
790  g%dx_Cv(i,j), g%dxCv(i,j), g%dx_Cv(i,j)-g%dxCv(i,j), &
791  g%geoLonCv(i,j), g%geoLatCv(i,j)
792 
793  call mom_error(fatal,"reset_face_lengths_file "//mesg)
794  endif
795  g%areaCv(i,j) = g%dyCv(i,j)*g%dx_Cv(i,j)
796  g%IareaCv(i,j) = 0.0
797  if (g%areaCv(i,j) > 0.0) g%IareaCv(i,j) = g%mask2dCv(i,j) / g%areaCv(i,j)
798  enddo ; enddo
799 
800  call calltree_leave(trim(mdl)//'()')
801 end subroutine reset_face_lengths_file
802 ! -----------------------------------------------------------------------------
803 
804 ! -----------------------------------------------------------------------------
805 !> This subroutine sets the open face lengths at selected points to restrict
806 !! passages to their observed widths from a list read from a file.
807 subroutine reset_face_lengths_list(G, param_file, US)
808  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
809  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
810  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
811 
812  ! Local variables
813  character(len=120), pointer, dimension(:) :: lines => null()
814  character(len=120) :: line
815  character(len=200) :: filename, chan_file, inputdir, mesg ! Strings for file/path
816  character(len=40) :: mdl = "reset_face_lengths_list" ! This subroutine's name.
817  real, pointer, dimension(:,:) :: &
818  u_lat => null(), u_lon => null(), v_lat => null(), v_lon => null()
819  real, pointer, dimension(:) :: &
820  u_width => null(), v_width => null()
821  real :: lat, lon ! The latitude and longitude of a point.
822  real :: len_lon ! The periodic range of longitudes, usually 360 degrees.
823  real :: len_lat ! The range of latitudes, usually 180 degrees.
824  real :: lon_p, lon_m ! The longitude of a point shifted by 360 degrees.
825  logical :: check_360 ! If true, check for longitudes that are shifted by
826  ! +/- 360 degrees from the specified range of values.
827  logical :: found_u, found_v
828  logical :: unit_in_use
829  integer :: ios, iounit, isu, isv
830  integer :: last, num_lines, nl_read, ln, npt, u_pt, v_pt
831  integer :: i, j, isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
832  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
833  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
834 
835  call calltree_enter(trim(mdl)//"(), MOM_shared_initialization.F90")
836 
837  call get_param(param_file, mdl, "CHANNEL_LIST_FILE", chan_file, &
838  "The file from which the list of narrowed channels is read.", &
839  default="MOM_channel_list")
840  call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
841  inputdir = slasher(inputdir)
842  filename = trim(inputdir)//trim(chan_file)
843  call log_param(param_file, mdl, "INPUTDIR/CHANNEL_LIST_FILE", filename)
844  call get_param(param_file, mdl, "CHANNEL_LIST_360_LON_CHECK", check_360, &
845  "If true, the channel configuration list works for any "//&
846  "longitudes in the range of -360 to 360.", default=.true.)
847 
848  if (is_root_pe()) then
849  ! Open the input file.
850  if (.not.file_exists(filename)) call mom_error(fatal, &
851  " reset_face_lengths_list: Unable to open "//trim(filename))
852 
853  ! Find an unused unit number.
854  do iounit=10,512
855  INQUIRE(iounit,opened=unit_in_use) ; if (.not.unit_in_use) exit
856  enddo
857  if (iounit >= 512) call mom_error(fatal, &
858  "reset_face_lengths_list: No unused file unit could be found.")
859 
860  ! Open the parameter file.
861  open(iounit, file=trim(filename), access='SEQUENTIAL', &
862  form='FORMATTED', action='READ', position='REWIND', iostat=ios)
863  if (ios /= 0) call mom_error(fatal, &
864  "reset_face_lengths_list: Error opening "//trim(filename))
865 
866  ! Count the number of u_width and v_width entries.
867  call read_face_length_list(iounit, filename, num_lines, lines)
868  endif
869 
870  len_lon = 360.0 ; if (g%len_lon > 0.0) len_lon = g%len_lon
871  len_lat = 180.0 ; if (g%len_lat > 0.0) len_lat = g%len_lat
872  ! Broadcast the number of lines and allocate the required space.
873  call broadcast(num_lines, root_pe())
874  u_pt = 0 ; v_pt = 0
875  if (num_lines > 0) then
876  allocate (lines(num_lines))
877  if (num_lines > 0) then
878  allocate(u_lat(2,num_lines)) ; u_lat(:,:) = -1e34
879  allocate(u_lon(2,num_lines)) ; u_lon(:,:) = -1e34
880  allocate(u_width(num_lines)) ; u_width(:) = -1e34
881 
882  allocate(v_lat(2,num_lines)) ; v_lat(:,:) = -1e34
883  allocate(v_lon(2,num_lines)) ; v_lon(:,:) = -1e34
884  allocate(v_width(num_lines)) ; v_width(:) = -1e34
885  endif
886 
887  ! Actually read the lines.
888  if (is_root_pe()) then
889  call read_face_length_list(iounit, filename, nl_read, lines)
890  if (nl_read /= num_lines) &
891  call mom_error(fatal, 'reset_face_lengths_list : Found different '// &
892  'number of valid lines on second reading of '//trim(filename))
893  close(iounit) ; iounit = -1
894  endif
895 
896  ! Broadcast the lines.
897  call broadcast(lines, 120, root_pe())
898 
899  ! Populate the u_width, etc., data.
900  do ln=1,num_lines
901  line = lines(ln)
902  ! Detect keywords
903  found_u = .false.; found_v = .false.
904  isu = index(uppercase(line), "U_WIDTH" ); if (isu > 0) found_u = .true.
905  isv = index(uppercase(line), "V_WIDTH" ); if (isv > 0) found_v = .true.
906 
907  ! Store and check the relevant values.
908  if (found_u) then
909  u_pt = u_pt + 1
910  read(line(isu+8:),*) u_lon(1:2,u_pt), u_lat(1:2,u_pt), u_width(u_pt)
911  if (is_root_pe()) then
912  if (check_360) then
913  if ((abs(u_lon(1,u_pt)) > len_lon) .or. (abs(u_lon(2,u_pt)) > len_lon)) &
914  call mom_error(warning, "reset_face_lengths_list : Out-of-bounds "//&
915  "u-longitude found when reading line "//trim(line)//" from file "//&
916  trim(filename))
917  if ((abs(u_lat(1,u_pt)) > len_lat) .or. (abs(u_lat(2,u_pt)) > len_lat)) &
918  call mom_error(warning, "reset_face_lengths_list : Out-of-bounds "//&
919  "u-latitude found when reading line "//trim(line)//" from file "//&
920  trim(filename))
921  endif
922  if (u_lat(1,u_pt) > u_lat(2,u_pt)) &
923  call mom_error(warning, "reset_face_lengths_list : Out-of-order "//&
924  "u-face latitudes found when reading line "//trim(line)//" from file "//&
925  trim(filename))
926  if (u_lon(1,u_pt) > u_lon(2,u_pt)) &
927  call mom_error(warning, "reset_face_lengths_list : Out-of-order "//&
928  "u-face longitudes found when reading line "//trim(line)//" from file "//&
929  trim(filename))
930  if (u_width(u_pt) < 0.0) &
931  call mom_error(warning, "reset_face_lengths_list : Negative "//&
932  "u-width found when reading line "//trim(line)//" from file "//&
933  trim(filename))
934  endif
935  elseif (found_v) then
936  v_pt = v_pt + 1
937  read(line(isv+8:),*) v_lon(1:2,v_pt), v_lat(1:2,v_pt), v_width(v_pt)
938  if (is_root_pe()) then
939  if (check_360) then
940  if ((abs(v_lon(1,v_pt)) > len_lon) .or. (abs(v_lon(2,v_pt)) > len_lon)) &
941  call mom_error(warning, "reset_face_lengths_list : Out-of-bounds "//&
942  "v-longitude found when reading line "//trim(line)//" from file "//&
943  trim(filename))
944  if ((abs(v_lat(1,v_pt)) > len_lat) .or. (abs(v_lat(2,v_pt)) > len_lat)) &
945  call mom_error(warning, "reset_face_lengths_list : Out-of-bounds "//&
946  "v-latitude found when reading line "//trim(line)//" from file "//&
947  trim(filename))
948  endif
949  if (v_lat(1,v_pt) > v_lat(2,v_pt)) &
950  call mom_error(warning, "reset_face_lengths_list : Out-of-order "//&
951  "v-face latitudes found when reading line "//trim(line)//" from file "//&
952  trim(filename))
953  if (v_lon(1,v_pt) > v_lon(2,v_pt)) &
954  call mom_error(warning, "reset_face_lengths_list : Out-of-order "//&
955  "v-face longitudes found when reading line "//trim(line)//" from file "//&
956  trim(filename))
957  if (v_width(v_pt) < 0.0) &
958  call mom_error(warning, "reset_face_lengths_list : Negative "//&
959  "v-width found when reading line "//trim(line)//" from file "//&
960  trim(filename))
961  endif
962  endif
963  enddo
964 
965  deallocate(lines)
966  endif
967 
968  do j=jsd,jed ; do i=isdb,iedb
969  lat = g%geoLatCu(i,j) ; lon = g%geoLonCu(i,j)
970  if (check_360) then ; lon_p = lon+len_lon ; lon_m = lon-len_lon
971  else ; lon_p = lon ; lon_m = lon ; endif
972 
973  do npt=1,u_pt
974  if (((lat >= u_lat(1,npt)) .and. (lat <= u_lat(2,npt))) .and. &
975  (((lon >= u_lon(1,npt)) .and. (lon <= u_lon(2,npt))) .or. &
976  ((lon_p >= u_lon(1,npt)) .and. (lon_p <= u_lon(2,npt))) .or. &
977  ((lon_m >= u_lon(1,npt)) .and. (lon_m <= u_lon(2,npt)))) ) then
978 
979  g%dy_Cu(i,j) = g%mask2dCu(i,j) * min(g%dyCu(i,j), max(u_width(npt), 0.0))
980  if (j>=g%jsc .and. j<=g%jec .and. i>=g%isc .and. i<=g%iec) then ! Limit messages/checking to compute domain
981  if ( g%mask2dCu(i,j) == 0.0 ) then
982  write(*,'(A,2F8.2,A,4F8.2,A)') "read_face_lengths_list : G%mask2dCu=0 at ",lat,lon," (",&
983  u_lat(1,npt), u_lat(2,npt), u_lon(1,npt), u_lon(2,npt),") so grid metric is unmodified."
984  else
985  write(*,'(A,2F8.2,A,4F8.2,A5,F9.2,A1)') &
986  "read_face_lengths_list : Modifying dy_Cu gridpoint at ",lat,lon," (",&
987  u_lat(1,npt), u_lat(2,npt), u_lon(1,npt), u_lon(2,npt),") to ",g%dy_Cu(i,j),"m"
988  endif
989  endif
990  endif
991  enddo
992 
993  g%areaCu(i,j) = g%dxCu(i,j)*g%dy_Cu(i,j)
994  g%IareaCu(i,j) = 0.0
995  if (g%areaCu(i,j) > 0.0) g%IareaCu(i,j) = g%mask2dCu(i,j) / g%areaCu(i,j)
996  enddo ; enddo
997 
998  do j=jsdb,jedb ; do i=isd,ied
999  lat = g%geoLatCv(i,j) ; lon = g%geoLonCv(i,j)
1000  if (check_360) then ; lon_p = lon+len_lon ; lon_m = lon-len_lon
1001  else ; lon_p = lon ; lon_m = lon ; endif
1002 
1003  do npt=1,v_pt
1004  if (((lat >= v_lat(1,npt)) .and. (lat <= v_lat(2,npt))) .and. &
1005  (((lon >= v_lon(1,npt)) .and. (lon <= v_lon(2,npt))) .or. &
1006  ((lon_p >= v_lon(1,npt)) .and. (lon_p <= v_lon(2,npt))) .or. &
1007  ((lon_m >= v_lon(1,npt)) .and. (lon_m <= v_lon(2,npt)))) ) then
1008  g%dx_Cv(i,j) = g%mask2dCv(i,j) * min(g%dxCv(i,j), max(v_width(npt), 0.0))
1009  if (i>=g%isc .and. i<=g%iec .and. j>=g%jsc .and. j<=g%jec) then ! Limit messages/checking to compute domain
1010  if ( g%mask2dCv(i,j) == 0.0 ) then
1011  write(*,'(A,2F8.2,A,4F8.2,A)') "read_face_lengths_list : G%mask2dCv=0 at ",lat,lon," (",&
1012  v_lat(1,npt), v_lat(2,npt), v_lon(1,npt), v_lon(2,npt),") so grid metric is unmodified."
1013  else
1014  write(*,'(A,2F8.2,A,4F8.2,A5,F9.2,A1)') &
1015  "read_face_lengths_list : Modifying dx_Cv gridpoint at ",lat,lon," (",&
1016  v_lat(1,npt), v_lat(2,npt), v_lon(1,npt), v_lon(2,npt),") to ",g%dx_Cv(i,j),"m"
1017  endif
1018  endif
1019  endif
1020  enddo
1021 
1022  g%areaCv(i,j) = g%dyCv(i,j)*g%dx_Cv(i,j)
1023  g%IareaCv(i,j) = 0.0
1024  if (g%areaCv(i,j) > 0.0) g%IareaCv(i,j) = g%mask2dCv(i,j) / g%areaCv(i,j)
1025  enddo ; enddo
1026 
1027  if (num_lines > 0) then
1028  deallocate(u_lat) ; deallocate(u_lon) ; deallocate(u_width)
1029  deallocate(v_lat) ; deallocate(v_lon) ; deallocate(v_width)
1030  endif
1031 
1032  call calltree_leave(trim(mdl)//'()')
1033 end subroutine reset_face_lengths_list
1034 ! -----------------------------------------------------------------------------
1035 
1036 ! -----------------------------------------------------------------------------
1037 !> This subroutine reads and counts the non-blank lines in the face length list file, after removing comments.
1038 subroutine read_face_length_list(iounit, filename, num_lines, lines)
1039  integer, intent(in) :: iounit !< An open I/O unit number for the file
1040  character(len=*), intent(in) :: filename !< The name of the face-length file to read
1041  integer, intent(out) :: num_lines !< The number of non-blank lines in the file
1042  character(len=120), dimension(:), pointer :: lines !< The non-blank lines, after removing comments
1043 
1044  ! This subroutine reads and counts the non-blank lines in the face length
1045  ! list file, after removing comments.
1046  character(len=120) :: line, line_up
1047  logical :: found_u, found_v
1048  integer :: isu, isv, icom, verbose
1049  integer :: last
1050 
1051  num_lines = 0
1052 
1053  if (iounit <= 0) return
1054  rewind(iounit)
1055  do while(.true.)
1056  read(iounit, '(a)', end=8, err=9) line
1057  last = len_trim(line)
1058  ! Eliminate either F90 or C comments from the line.
1059  icom = index(line(:last), "!") ; if (icom > 0) last = icom-1
1060  icom = index(line(:last), "/*") ; if (icom > 0) last = icom-1
1061  if (last < 1) cycle
1062 
1063  ! Detect keywords
1064  line_up = uppercase(line)
1065  found_u = .false.; found_v = .false.
1066  isu = index(line_up(:last), "U_WIDTH" ); if (isu > 0) found_u = .true.
1067  isv = index(line_up(:last), "V_WIDTH" ); if (isv > 0) found_v = .true.
1068 
1069  if (found_u .and. found_v) call mom_error(fatal, &
1070  "read_face_length_list : both U_WIDTH and V_WIDTH found when "//&
1071  "reading the line "//trim(line(:last))//" in file "//trim(filename))
1072  if (found_u .or. found_v) then
1073  num_lines = num_lines + 1
1074  if (associated(lines)) then
1075  lines(num_lines) = line(1:last)
1076  endif
1077  endif
1078  enddo ! while (.true.)
1079 
1080 8 continue
1081  return
1082 
1083 9 call mom_error(fatal, "read_face_length_list : "//&
1084  "Error while reading file "//trim(filename))
1085 
1086 end subroutine read_face_length_list
1087 ! -----------------------------------------------------------------------------
1088 
1089 ! -----------------------------------------------------------------------------
1090 !> Set the bathymetry at velocity points to be the maximum of the depths at the
1091 !! neighoring tracer points.
1092 subroutine set_velocity_depth_max(G)
1093  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
1094  ! This subroutine sets the 4 bottom depths at velocity points to be the
1095  ! maximum of the adjacent depths.
1096  integer :: i, j
1097 
1098  do i=g%isd,g%ied-1 ; do j=g%jsd,g%jed
1099  g%Dblock_u(i,j) = g%mask2dCu(i,j) * max(g%bathyT(i,j), g%bathyT(i+1,j))
1100  g%Dopen_u(i,j) = g%Dblock_u(i,j)
1101  enddo ; enddo
1102  do i=g%isd,g%ied ; do j=g%jsd,g%jed-1
1103  g%Dblock_v(i,j) = g%mask2dCv(i,j) * max(g%bathyT(i,j), g%bathyT(i,j+1))
1104  g%Dopen_v(i,j) = g%Dblock_v(i,j)
1105  enddo ; enddo
1106 end subroutine set_velocity_depth_max
1107 ! -----------------------------------------------------------------------------
1108 
1109 ! -----------------------------------------------------------------------------
1110 !> Set the bathymetry at velocity points to be the minimum of the depths at the
1111 !! neighoring tracer points.
1112 subroutine set_velocity_depth_min(G)
1113  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
1114  ! This subroutine sets the 4 bottom depths at velocity points to be the
1115  ! minimum of the adjacent depths.
1116  integer :: i, j
1117 
1118  do i=g%isd,g%ied-1 ; do j=g%jsd,g%jed
1119  g%Dblock_u(i,j) = g%mask2dCu(i,j) * min(g%bathyT(i,j), g%bathyT(i+1,j))
1120  g%Dopen_u(i,j) = g%Dblock_u(i,j)
1121  enddo ; enddo
1122  do i=g%isd,g%ied ; do j=g%jsd,g%jed-1
1123  g%Dblock_v(i,j) = g%mask2dCv(i,j) * min(g%bathyT(i,j), g%bathyT(i,j+1))
1124  g%Dopen_v(i,j) = g%Dblock_v(i,j)
1125  enddo ; enddo
1126 end subroutine set_velocity_depth_min
1127 ! -----------------------------------------------------------------------------
1128 
1129 ! -----------------------------------------------------------------------------
1130 !> Pre-compute global integrals of grid quantities (like masked ocean area) for
1131 !! later use in reporting diagnostics
1132 subroutine compute_global_grid_integrals(G)
1133  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
1134 
1135  ! Local variables
1136  real, dimension(G%isc:G%iec, G%jsc:G%jec) :: tmpforsumming
1137  integer :: i,j
1138 
1139  tmpforsumming(:,:) = 0.
1140  g%areaT_global = 0.0 ; g%IareaT_global = 0.0
1141  do j=g%jsc,g%jec ; do i=g%isc,g%iec
1142  tmpforsumming(i,j) = g%areaT(i,j) * g%mask2dT(i,j)
1143  enddo ; enddo
1144  g%areaT_global = reproducing_sum(tmpforsumming)
1145 
1146  if (g%areaT_global == 0.0) &
1147  call mom_error(fatal, "compute_global_grid_integrals: "//&
1148  "zero ocean area (check topography?)")
1149 
1150  g%IareaT_global = 1. / g%areaT_global
1151 end subroutine compute_global_grid_integrals
1152 ! -----------------------------------------------------------------------------
1153 
1154 ! -----------------------------------------------------------------------------
1155 !> Write out a file describing the topography, Coriolis parameter, grid locations
1156 !! and various other fixed fields from the grid.
1157 subroutine write_ocean_geometry_file(G, param_file, directory, geom_file, US)
1158  type(dyn_horgrid_type), intent(inout) :: g !< The dynamic horizontal grid
1159  type(param_file_type), intent(in) :: param_file !< Parameter file structure
1160  character(len=*), intent(in) :: directory !< The directory into which to place the geometry file.
1161  character(len=*), optional, intent(in) :: geom_file !< If present, the name of the geometry file
1162  !! (otherwise the file is "ocean_geometry")
1163  type(unit_scale_type), optional, intent(in) :: us !< A dimensional unit scaling type
1164 
1165  ! Local variables.
1166  character(len=240) :: filepath
1167  character(len=40) :: mdl = "write_ocean_geometry_file"
1168  integer, parameter :: nflds=23
1169  type(vardesc) :: vars(nflds)
1170  type(fieldtype) :: fields(nflds)
1171  real :: z_to_m_scale ! A unit conversion factor from Z to m.
1172  real :: s_to_t_scale ! A unit conversion factor from T-1 to s-1.
1173  integer :: unit
1174  integer :: file_threading
1175  integer :: nflds_used
1176  integer :: i, j, is, ie, js, je, isq, ieq, jsq, jeq
1177  integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
1178  logical :: multiple_files
1179  real, dimension(G%isd :G%ied ,G%jsd :G%jed ) :: out_h
1180  real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB) :: out_q
1181  real, dimension(G%IsdB:G%IedB,G%jsd :G%jed ) :: out_u
1182  real, dimension(G%isd :G%ied ,G%JsdB:G%JedB) :: out_v
1183 
1184  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1185  isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
1186  isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1187  isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
1188 
1189  z_to_m_scale = 1.0 ; if (present(us)) z_to_m_scale = us%Z_to_m
1190  s_to_t_scale = 1.0 ; if (present(us)) s_to_t_scale = us%s_to_T
1191 
1192 ! vardesc is a structure defined in MOM_io.F90. The elements of
1193 ! this structure, in order, are:
1194 ! (1) the variable name for the NetCDF file
1195 ! (2) the variable's long name
1196 ! (3) a character indicating the horizontal grid, which may be '1' (column),
1197 ! 'h', 'q', 'u', or 'v', for the corresponding C-grid variable
1198 ! (4) a character indicating the vertical grid, which may be 'L' (layer),
1199 ! 'i' (interface), or '1' (no vertical location)
1200 ! (5) a character indicating the time levels of the field, which may be
1201 ! 's' (snap-shot), 'p' (periodic), or '1' (no time variation)
1202 ! (6) the variable's units
1203  vars(1) = var_desc("geolatb","degree","latitude at corner (Bu) points",'q','1','1')
1204  vars(2) = var_desc("geolonb","degree","longitude at corner (Bu) points",'q','1','1')
1205  vars(3) = var_desc("geolat","degree", "latitude at tracer (T) points", 'h','1','1')
1206  vars(4) = var_desc("geolon","degree","longitude at tracer (T) points",'h','1','1')
1207  vars(5) = var_desc("D","meter","Basin Depth",'h','1','1')
1208  vars(6) = var_desc("f","s-1","Coriolis Parameter",'q','1','1')
1209  vars(7) = var_desc("dxCv","m","Zonal grid spacing at v points",'v','1','1')
1210  vars(8) = var_desc("dyCu","m","Meridional grid spacing at u points",'u','1','1')
1211  vars(9) = var_desc("dxCu","m","Zonal grid spacing at u points",'u','1','1')
1212  vars(10)= var_desc("dyCv","m","Meridional grid spacing at v points",'v','1','1')
1213  vars(11)= var_desc("dxT","m","Zonal grid spacing at h points",'h','1','1')
1214  vars(12)= var_desc("dyT","m","Meridional grid spacing at h points",'h','1','1')
1215  vars(13)= var_desc("dxBu","m","Zonal grid spacing at q points",'q','1','1')
1216  vars(14)= var_desc("dyBu","m","Meridional grid spacing at q points",'q','1','1')
1217  vars(15)= var_desc("Ah","m2","Area of h cells",'h','1','1')
1218  vars(16)= var_desc("Aq","m2","Area of q cells",'q','1','1')
1219 
1220  vars(17)= var_desc("dxCvo","m","Open zonal grid spacing at v points",'v','1','1')
1221  vars(18)= var_desc("dyCuo","m","Open meridional grid spacing at u points",'u','1','1')
1222  vars(19)= var_desc("wet", "nondim", "land or ocean?", 'h','1','1')
1223 
1224  vars(20) = var_desc("Dblock_u","m","Blocked depth at u points",'u','1','1')
1225  vars(21) = var_desc("Dopen_u","m","Open depth at u points",'u','1','1')
1226  vars(22) = var_desc("Dblock_v","m","Blocked depth at v points",'v','1','1')
1227  vars(23) = var_desc("Dopen_v","m","Open depth at v points",'v','1','1')
1228 
1229 
1230  nflds_used = 19 ; if (g%bathymetry_at_vel) nflds_used = 23
1231 
1232  if (present(geom_file)) then
1233  filepath = trim(directory) // trim(geom_file)
1234  else
1235  filepath = trim(directory) // "ocean_geometry"
1236  endif
1237 
1238  out_h(:,:) = 0.0
1239  out_u(:,:) = 0.0
1240  out_v(:,:) = 0.0
1241  out_q(:,:) = 0.0
1242 
1243  call get_param(param_file, mdl, "PARALLEL_RESTARTFILES", multiple_files, &
1244  "If true, each processor writes its own restart file, "//&
1245  "otherwise a single restart file is generated", &
1246  default=.false.)
1247  file_threading = single_file
1248  if (multiple_files) file_threading = multiple
1249 
1250  call create_file(unit, trim(filepath), vars, nflds_used, fields, &
1251  file_threading, dg=g)
1252 
1253  do j=jsq,jeq; do i=isq,ieq; out_q(i,j) = g%geoLatBu(i,j); enddo ; enddo
1254  call write_field(unit, fields(1), g%Domain%mpp_domain, out_q)
1255  do j=jsq,jeq; do i=isq,ieq; out_q(i,j) = g%geoLonBu(i,j); enddo ; enddo
1256  call write_field(unit, fields(2), g%Domain%mpp_domain, out_q)
1257  call write_field(unit, fields(3), g%Domain%mpp_domain, g%geoLatT)
1258  call write_field(unit, fields(4), g%Domain%mpp_domain, g%geoLonT)
1259 
1260  do j=js,je ; do i=is,ie ; out_h(i,j) = z_to_m_scale*g%bathyT(i,j) ; enddo ; enddo
1261  call write_field(unit, fields(5), g%Domain%mpp_domain, out_h)
1262  do j=jsq,jeq ; do i=isq,ieq ; out_q(i,j) = s_to_t_scale*g%CoriolisBu(i,j) ; enddo ; enddo
1263  call write_field(unit, fields(6), g%Domain%mpp_domain, out_q)
1264 
1265  ! I think that all of these copies are holdovers from a much earlier
1266  ! ancestor code in which many of the metrics were macros that could have
1267  ! had reduced dimensions, and that they are no longer needed in MOM6. -RWH
1268  do j=jsq,jeq ; do i=is,ie ; out_v(i,j) = g%dxCv(i,j) ; enddo ; enddo
1269  call write_field(unit, fields(7), g%Domain%mpp_domain, out_v)
1270  do j=js,je ; do i=isq,ieq ; out_u(i,j) = g%dyCu(i,j) ; enddo ; enddo
1271  call write_field(unit, fields(8), g%Domain%mpp_domain, out_u)
1272 
1273  do j=js,je ; do i=isq,ieq ; out_u(i,j) = g%dxCu(i,j) ; enddo ; enddo
1274  call write_field(unit, fields(9), g%Domain%mpp_domain, out_u)
1275  do j=jsq,jeq ; do i=is,ie ; out_v(i,j) = g%dyCv(i,j) ; enddo ; enddo
1276  call write_field(unit, fields(10), g%Domain%mpp_domain, out_v)
1277 
1278  do j=js,je ; do i=is,ie ; out_h(i,j) = g%dxT(i,j); enddo ; enddo
1279  call write_field(unit, fields(11), g%Domain%mpp_domain, out_h)
1280  do j=js,je ; do i=is,ie ; out_h(i,j) = g%dyT(i,j) ; enddo ; enddo
1281  call write_field(unit, fields(12), g%Domain%mpp_domain, out_h)
1282 
1283  do j=jsq,jeq ; do i=isq,ieq ; out_q(i,j) = g%dxBu(i,j) ; enddo ; enddo
1284  call write_field(unit, fields(13), g%Domain%mpp_domain, out_q)
1285  do j=jsq,jeq ; do i=isq,ieq ; out_q(i,j) = g%dyBu(i,j) ; enddo ; enddo
1286  call write_field(unit, fields(14), g%Domain%mpp_domain, out_q)
1287 
1288  do j=js,je ; do i=is,ie ; out_h(i,j) = g%areaT(i,j) ; enddo ; enddo
1289  call write_field(unit, fields(15), g%Domain%mpp_domain, out_h)
1290  do j=jsq,jeq ; do i=isq,ieq ; out_q(i,j) = g%areaBu(i,j) ; enddo ; enddo
1291  call write_field(unit, fields(16), g%Domain%mpp_domain, out_q)
1292 
1293  call write_field(unit, fields(17), g%Domain%mpp_domain, g%dx_Cv)
1294  call write_field(unit, fields(18), g%Domain%mpp_domain, g%dy_Cu)
1295  call write_field(unit, fields(19), g%Domain%mpp_domain, g%mask2dT)
1296 
1297  if (g%bathymetry_at_vel) then
1298  do j=js,je ; do i=isq,ieq ; out_u(i,j) = z_to_m_scale*g%Dblock_u(i,j) ; enddo ; enddo
1299  call write_field(unit, fields(20), g%Domain%mpp_domain, out_u)
1300  do j=js,je ; do i=isq,ieq ; out_u(i,j) = z_to_m_scale*g%Dopen_u(i,j) ; enddo ; enddo
1301  call write_field(unit, fields(21), g%Domain%mpp_domain, out_u)
1302  do j=jsq,jeq ; do i=is,ie ; out_v(i,j) = z_to_m_scale*g%Dblock_v(i,j) ; enddo ; enddo
1303  call write_field(unit, fields(22), g%Domain%mpp_domain, out_v)
1304  do j=jsq,jeq ; do i=is,ie ; out_v(i,j) = z_to_m_scale*g%Dopen_v(i,j) ; enddo ; enddo
1305  call write_field(unit, fields(23), g%Domain%mpp_domain, out_v)
1306  endif
1307 
1308  call close_file(unit)
1309 
1310 end subroutine write_ocean_geometry_file
1311 
1312 end module mom_shared_initialization
mom_shared_initialization
Code that initializes fixed aspects of the model grid, such as horizontal grid metrics,...
Definition: MOM_shared_initialization.F90:3
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_dyn_horgrid
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Definition: MOM_dyn_horgrid.F90:3
mom_string_functions
Handy functions for manipulating strings.
Definition: MOM_string_functions.F90:2
mom_domains::pass_var
Do a halo update on an array.
Definition: MOM_domains.F90:49
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_io
This module contains I/O framework code.
Definition: MOM_io.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
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
mom_domains
Describes the decomposed MOM domain and has routines for communications across PEs.
Definition: MOM_domains.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_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_io::vardesc
Type for describing a variable, typically a tracer.
Definition: MOM_io.F90:53
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_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_dyn_horgrid::dyn_horgrid_type
Describes the horizontal ocean grid with only dynamic memory arrays.
Definition: MOM_dyn_horgrid.F90:22