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