MOM_restart.F90

1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> The MOM6 facility for reading and writing restart files, and querying what has been read.
6module mom_restart
7
8use, intrinsic :: iso_fortran_env, only : int64
9use mom_array_transform, only : rotate_array, rotate_vector, rotate_array_pair
10use mom_checksums, only : chksum => field_checksum
11use mom_domains, only : pe_here, num_pes, agrid, bgrid_ne, cgrid_ne
12use mom_error_handler, only : mom_error, mom_mesg, fatal, warning, note, is_root_pe, mom_get_verbosity
13use mom_file_parser, only : get_param, log_param, log_version, param_file_type
14use mom_grid, only : ocean_grid_type
15use mom_io, only : create_mom_file, file_exists
16use mom_io, only : mom_infra_file, mom_field
17use mom_io, only : mom_read_data, read_data, mom_write_field, field_exists
19use mom_io, only : multiple, readonly_file, single_file
20use mom_io, only : center, corner, north_face, east_face
23use mom_time_manager, only : time_type, time_type_to_real, real_to_time
24use mom_time_manager, only : days_in_month, get_date, set_date
26
27implicit none ; private
28
29public restart_init, restart_end, restore_state, register_restart_field
30public copy_restart_var, copy_restart_vector
31public save_restart, query_initialized, set_initialized, only_read_from_restarts
34public register_restart_field_as_obsolete, register_restart_pair
35public lock_check
36
37! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
38! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
39! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
40! vary with the Boussinesq approximation, the Boussinesq variant is given first.
41! The functions in this module work with variables with arbitrary units, in which case the
42! arbitrary rescaled units are indicated with [A ~> a], while the unscaled units are just [a].
43
44!> A type for making arrays of pointers to 4-d arrays
45type p4d
46 real, dimension(:,:,:,:), pointer :: p => null() !< A pointer to a 4d array in arbitrary rescaled units [A ~> a]
47end type p4d
48
49!> A type for making arrays of pointers to 3-d arrays
50type p3d
51 real, dimension(:,:,:), pointer :: p => null() !< A pointer to a 3d array in arbitrary rescaled units [A ~> a]
52end type p3d
53
54!> A type for making arrays of pointers to 2-d arrays
55type p2d
56 real, dimension(:,:), pointer :: p => null() !< A pointer to a 2d array in arbitrary rescaled units [A ~> a]
57end type p2d
58
59!> A type for making arrays of pointers to 1-d arrays
60type p1d
61 real, dimension(:), pointer :: p => null() !< A pointer to a 1d array in arbitrary rescaled units [A ~> a]
62end type p1d
63
64!> A type for making arrays of pointers to scalars
65type p0d
66 real, pointer :: p => null() !< A pointer to a scalar in arbitrary rescaled units [A ~> a]
67end type p0d
68
69!> A structure with information about a single restart field
71 type(vardesc) :: vars !< Description of a field that is to be read from or written
72 !! to the restart file.
73 logical :: mand_var !< If .true. the run will abort if this field is not successfully
74 !! read from the restart file.
75 logical :: initialized !< .true. if this field has been read from the restart file.
76 character(len=32) :: var_name !< A name by which a variable may be queried.
77 real :: conv = 1.0 !< A factor by which a restart field should be multiplied before it
78 !! is written to a restart file, usually to convert it to MKS or
79 !! other standard units [a A-1 ~> 1]. When read, the restart field
80 !! is multiplied by the reciprocal of this factor.
81end type field_restart
82
83!> A structure to store information about restart fields that are no longer used
85 character(len=32) :: field_name !< Name of restart field that is no longer in use
86 character(len=32) :: replacement_name !< Name of replacement restart field, if applicable
87end type obsolete_restart
88
89!> A restart registry and the control structure for restarts
90type, public :: mom_restart_cs ; private
91 logical :: initialized = .false. !< True if this control structure has been initialized.
92 logical :: restart !< restart is set to .true. if the run has been started from a full restart
93 !! file. Otherwise some fields must be initialized approximately.
94 integer :: novars = 0 !< The number of restart fields that have been registered.
95 integer :: num_obsolete_vars = 0 !< The number of obsolete restart fields that have been registered.
96 logical :: parallel_restartfiles !< If true, the IO layout is used to group processors that write
97 !! to the same restart file or each processor writes its own
98 !! (numbered) restart file. If false, a single restart file is
99 !! generated after internally combining output from all PEs.
100 logical :: new_run !< If true, the input filenames and restart file existence will
101 !! result in a new run that is not initialized from restart files.
102 logical :: new_run_set = .false. !< If true, new_run has been determined for this restart_CS.
103 logical :: checksum_required !< If true, require the restart checksums to match and error out otherwise.
104 !! Users may want to avoid this comparison if for example the restarts are
105 !! made from a run with a different mask_table than the current run,
106 !! in which case the checksums will not match and cause crash.
107 logical :: symmetric_checksums !< If true, do the restart checksums on all the edge points for
108 !! a non-reentrant grid. Setting this to true requires that
109 !! SYMMETRIC_MEMORY_ is defined at compile time.
110 logical :: unsigned_zeros !< If true, convert any negative zeros that would be written to
111 !! the restart file into ordinary unsigned zeros. This does not
112 !! change answers, but it can be helpful in comparing restart
113 !! files after grid rotation, for example.
114 logical :: reentrant_x !< If true, the domain is reentrant in the x-direction. This is only
115 !! used here to determine the extent of the restart checksums.
116 logical :: reentrant_y !< If true, the domain is reentrant in the y-direction. This is only
117 !! used here to determine the extent of the restart checksums.
118 character(len=240) :: restartfile !< The name or name root for MOM restart files.
119 integer :: turns !< Number of quarter turns from input to model domain
120 logical :: locked = .false. !< If true this registry has been locked and no further restart
121 !! fields can be added without explicitly unlocking the registry.
122
123 !> An array of descriptions of the registered fields
124 type(field_restart), pointer :: restart_field(:) => null()
125
126 !> An array of obsolete restart fields
127 type(obsolete_restart), pointer :: restart_obsolete(:) => null()
128
129 !>@{ Pointers to the fields that have been registered for restarts
130 type(p0d), pointer :: var_ptr0d(:) => null()
131 type(p1d), pointer :: var_ptr1d(:) => null()
132 type(p2d), pointer :: var_ptr2d(:) => null()
133 type(p3d), pointer :: var_ptr3d(:) => null()
134 type(p4d), pointer :: var_ptr4d(:) => null()
135 !>@}
136 integer :: max_fields !< The maximum number of restart fields
137end type mom_restart_cs
138
139!> Register fields for restarts
140interface register_restart_field
146end interface
147
148!> Register a pair of restart fields whose rotations map onto each other
149interface register_restart_pair
150 module procedure register_restart_pair_ptr2d
151 module procedure register_restart_pair_ptr3d
152 module procedure register_restart_pair_ptr4d
153end interface register_restart_pair
154
155!> Indicate whether a field has been read from a restart file
156interface query_initialized
157 module procedure query_initialized_name
163end interface
164
165!> Specify that a field has been initialized, even if it was not read from a restart file
166interface set_initialized
170end interface
171
172!> Copy the restart variable with the specified name into an array, perhaps after rotation
173interface copy_restart_var
174 module procedure copy_restart_var_3d
175end interface copy_restart_var
176
177!> Copy the restart vector component variables with the specified names into a pair of arrays,
178!! perhaps after rotation
179interface copy_restart_vector
180 module procedure copy_restart_vector_3d
181end interface copy_restart_vector
182
183!> Read optional variables from restart files.
184interface only_read_from_restarts
185 module procedure only_read_restart_field_4d
186 module procedure only_read_restart_field_3d
187 module procedure only_read_restart_field_2d
188! module procedure only_read_restart_field_1d
189! module procedure only_read_restart_field_0d
190 module procedure only_read_restart_pair_3d
191end interface
192
193contains
194
195!> Register a restart field as obsolete
196subroutine register_restart_field_as_obsolete(field_name, replacement_name, CS)
197 character(*), intent(in) :: field_name !< Name of restart field that is no longer in use
198 character(*), intent(in) :: replacement_name !< Name of replacement restart field, if applicable
199 type(mom_restart_cs), intent(inout) :: cs !< MOM restart control struct
200
201 cs%num_obsolete_vars = cs%num_obsolete_vars+1
202 cs%restart_obsolete(cs%num_obsolete_vars)%field_name = field_name
203 cs%restart_obsolete(cs%num_obsolete_vars)%replacement_name = replacement_name
205
206!> Register a 3-d field for restarts, providing the metadata in a structure
207subroutine register_restart_field_ptr3d(f_ptr, var_desc, mandatory, CS, conversion)
208 real, dimension(:,:,:), &
209 target, intent(in) :: f_ptr !< A pointer to the field to be read or written
210 !! in arbitrary rescaled units [A ~> a]
211 type(vardesc), intent(in) :: var_desc !< A structure with metadata about this variable
212 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
213 !! successfully read from the restart file.
214 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
215 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
216 !! before it is written [a A-1 ~> 1], 1 by default.
217
218 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
219 "register_restart_field: Module must be initialized before it is used.")
220
221 call lock_check(cs, var_desc)
222
223 cs%novars = cs%novars+1
224 if (cs%novars > cs%max_fields) return ! This is an error that will be reported
225 ! once the total number of fields is known.
226
227 cs%restart_field(cs%novars)%vars = var_desc
228 cs%restart_field(cs%novars)%mand_var = mandatory
229 cs%restart_field(cs%novars)%initialized = .false.
230 cs%restart_field(cs%novars)%conv = 1.0
231 if (present(conversion)) cs%restart_field(cs%novars)%conv = conversion
232 call query_vardesc(cs%restart_field(cs%novars)%vars, &
233 name=cs%restart_field(cs%novars)%var_name, &
234 caller="register_restart_field_ptr3d")
235
236 cs%var_ptr3d(cs%novars)%p => f_ptr
237 cs%var_ptr4d(cs%novars)%p => null()
238 cs%var_ptr2d(cs%novars)%p => null()
239 cs%var_ptr1d(cs%novars)%p => null()
240 cs%var_ptr0d(cs%novars)%p => null()
241
242end subroutine register_restart_field_ptr3d
243
244!> Register a 4-d field for restarts, providing the metadata in a structure
245subroutine register_restart_field_ptr4d(f_ptr, var_desc, mandatory, CS, conversion)
246 real, dimension(:,:,:,:), &
247 target, intent(in) :: f_ptr !< A pointer to the field to be read or written
248 !! in arbitrary rescaled units [A ~> a]
249 type(vardesc), intent(in) :: var_desc !< A structure with metadata about this variable
250 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
251 !! successfully read from the restart file.
252 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
253 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
254 !! before it is written [a A-1 ~> 1], 1 by default.
255
256 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
257 "register_restart_field: Module must be initialized before it is used.")
258
259 call lock_check(cs, var_desc)
260
261 cs%novars = cs%novars+1
262 if (cs%novars > cs%max_fields) return ! This is an error that will be reported
263 ! once the total number of fields is known.
264
265 cs%restart_field(cs%novars)%vars = var_desc
266 cs%restart_field(cs%novars)%mand_var = mandatory
267 cs%restart_field(cs%novars)%initialized = .false.
268 cs%restart_field(cs%novars)%conv = 1.0
269 if (present(conversion)) cs%restart_field(cs%novars)%conv = conversion
270 call query_vardesc(cs%restart_field(cs%novars)%vars, &
271 name=cs%restart_field(cs%novars)%var_name, &
272 caller="register_restart_field_ptr4d")
273
274 cs%var_ptr4d(cs%novars)%p => f_ptr
275 cs%var_ptr3d(cs%novars)%p => null()
276 cs%var_ptr2d(cs%novars)%p => null()
277 cs%var_ptr1d(cs%novars)%p => null()
278 cs%var_ptr0d(cs%novars)%p => null()
279
280end subroutine register_restart_field_ptr4d
281
282!> Register a 2-d field for restarts, providing the metadata in a structure
283subroutine register_restart_field_ptr2d(f_ptr, var_desc, mandatory, CS, conversion)
284 real, dimension(:,:), &
285 target, intent(in) :: f_ptr !< A pointer to the field to be read or written
286 !! in arbitrary rescaled units [A ~> a]
287 type(vardesc), intent(in) :: var_desc !< A structure with metadata about this variable
288 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
289 !! successfully read from the restart file.
290 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
291 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
292 !! before it is written [a A-1 ~> 1], 1 by default.
293
294 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
295 "register_restart_field: Module must be initialized before it is used.")
296
297 call lock_check(cs, var_desc)
298
299 cs%novars = cs%novars+1
300 if (cs%novars > cs%max_fields) return ! This is an error that will be reported
301 ! once the total number of fields is known.
302
303 cs%restart_field(cs%novars)%vars = var_desc
304 cs%restart_field(cs%novars)%mand_var = mandatory
305 cs%restart_field(cs%novars)%initialized = .false.
306 cs%restart_field(cs%novars)%conv = 1.0
307 if (present(conversion)) cs%restart_field(cs%novars)%conv = conversion
308 call query_vardesc(cs%restart_field(cs%novars)%vars, &
309 name=cs%restart_field(cs%novars)%var_name, &
310 caller="register_restart_field_ptr2d")
311
312 cs%var_ptr2d(cs%novars)%p => f_ptr
313 cs%var_ptr4d(cs%novars)%p => null()
314 cs%var_ptr3d(cs%novars)%p => null()
315 cs%var_ptr1d(cs%novars)%p => null()
316 cs%var_ptr0d(cs%novars)%p => null()
317
318end subroutine register_restart_field_ptr2d
319
320!> Register a 1-d field for restarts, providing the metadata in a structure
321subroutine register_restart_field_ptr1d(f_ptr, var_desc, mandatory, CS, conversion)
322 real, dimension(:), target, intent(in) :: f_ptr !< A pointer to the field to be read or written
323 !! in arbitrary rescaled units [A ~> a]
324 type(vardesc), intent(in) :: var_desc !< A structure with metadata about this variable
325 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
326 !! successfully read from the restart file.
327 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
328 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
329 !! before it is written [a A-1 ~> 1], 1 by default.
330
331 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
332 "register_restart_field: Module must be initialized before it is used.")
333
334 call lock_check(cs, var_desc)
335
336 cs%novars = cs%novars+1
337 if (cs%novars > cs%max_fields) return ! This is an error that will be reported
338 ! once the total number of fields is known.
339
340 cs%restart_field(cs%novars)%vars = var_desc
341 cs%restart_field(cs%novars)%mand_var = mandatory
342 cs%restart_field(cs%novars)%initialized = .false.
343 cs%restart_field(cs%novars)%conv = 1.0
344 if (present(conversion)) cs%restart_field(cs%novars)%conv = conversion
345 call query_vardesc(cs%restart_field(cs%novars)%vars, &
346 name=cs%restart_field(cs%novars)%var_name, &
347 caller="register_restart_field_ptr1d")
348
349 cs%var_ptr1d(cs%novars)%p => f_ptr
350 cs%var_ptr4d(cs%novars)%p => null()
351 cs%var_ptr3d(cs%novars)%p => null()
352 cs%var_ptr2d(cs%novars)%p => null()
353 cs%var_ptr0d(cs%novars)%p => null()
354
355end subroutine register_restart_field_ptr1d
356
357!> Register a 0-d field for restarts, providing the metadata in a structure
358subroutine register_restart_field_ptr0d(f_ptr, var_desc, mandatory, CS, conversion)
359 real, target, intent(in) :: f_ptr !< A pointer to the field to be read or written
360 !! in arbitrary rescaled units [A ~> a]
361 type(vardesc), intent(in) :: var_desc !< A structure with metadata about this variable
362 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
363 !! successfully read from the restart file.
364 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
365 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
366 !! before it is written [a A-1 ~> 1], 1 by default.
367
368 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
369 "register_restart_field: Module must be initialized before it is used.")
370
371 call lock_check(cs, var_desc)
372
373 cs%novars = cs%novars+1
374 if (cs%novars > cs%max_fields) return ! This is an error that will be reported
375 ! once the total number of fields is known.
376
377 cs%restart_field(cs%novars)%vars = var_desc
378 cs%restart_field(cs%novars)%mand_var = mandatory
379 cs%restart_field(cs%novars)%initialized = .false.
380 cs%restart_field(cs%novars)%conv = 1.0
381 if (present(conversion)) cs%restart_field(cs%novars)%conv = conversion
382 call query_vardesc(cs%restart_field(cs%novars)%vars, &
383 name=cs%restart_field(cs%novars)%var_name, &
384 caller="register_restart_field_ptr0d")
385
386 cs%var_ptr0d(cs%novars)%p => f_ptr
387 cs%var_ptr4d(cs%novars)%p => null()
388 cs%var_ptr3d(cs%novars)%p => null()
389 cs%var_ptr2d(cs%novars)%p => null()
390 cs%var_ptr1d(cs%novars)%p => null()
391
392end subroutine register_restart_field_ptr0d
393
394
395!> Register a pair of rotationally equivalent 2d restart fields
396subroutine register_restart_pair_ptr2d(a_ptr, b_ptr, a_desc, b_desc, &
397 mandatory, CS, conversion, scalar_pair)
398 real, dimension(:,:), target, intent(in) :: a_ptr !< First field pointer
399 !! in arbitrary rescaled units [A ~> a]
400 real, dimension(:,:), target, intent(in) :: b_ptr !< Second field pointer
401 !! in arbitrary rescaled units [A ~> a]
402 type(vardesc), intent(in) :: a_desc !< First field descriptor
403 type(vardesc), intent(in) :: b_desc !< Second field descriptor
404 logical, intent(in) :: mandatory !< If true, abort if field is missing
405 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control structure
406 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
407 !! before it is written [a A-1 ~> 1], 1 by default.
408 logical, optional, intent(in) :: scalar_pair !< If true, the arrays describe a pair of
409 !! scalars, instead of vector components
410 !! whose signs change when rotated
411
412 ! Local variables
413 real :: a_conv, b_conv ! Factors to multipy the a- and b-components by before they are written,
414 ! including sign changes to account for grid rotation [a A-1 ~> 1]
415
416 call lock_check(cs, a_desc)
417 call set_conversion_pair(a_conv, b_conv, cs%turns, conversion, scalar_pair)
418
419 if (modulo(cs%turns, 2) == 0) then ! This is the usual case.
420 call register_restart_field(a_ptr, a_desc, mandatory, cs, conversion=a_conv)
421 call register_restart_field(b_ptr, b_desc, mandatory, cs, conversion=b_conv)
422 else
423 call register_restart_field(b_ptr, a_desc, mandatory, cs, conversion=a_conv)
424 call register_restart_field(a_ptr, b_desc, mandatory, cs, conversion=b_conv)
425 endif
426end subroutine register_restart_pair_ptr2d
427
428
429!> Register a pair of rotationally equivalent 3d restart fields
430subroutine register_restart_pair_ptr3d(a_ptr, b_ptr, a_desc, b_desc, &
431 mandatory, CS, conversion, scalar_pair)
432 real, dimension(:,:,:), target, intent(in) :: a_ptr !< First field pointer
433 !! in arbitrary rescaled units [A ~> a]
434 real, dimension(:,:,:), target, intent(in) :: b_ptr !< Second field pointer
435 !! in arbitrary rescaled units [A ~> a]
436 type(vardesc), intent(in) :: a_desc !< First field descriptor
437 type(vardesc), intent(in) :: b_desc !< Second field descriptor
438 logical, intent(in) :: mandatory !< If true, abort if field is missing
439 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control structure
440 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
441 !! before it is written [a A-1 ~> 1], 1 by default.
442 logical, optional, intent(in) :: scalar_pair !< If true, the arrays describe a pair of
443 !! scalars, instead of vector components
444 !! whose signs change when rotated
445
446 ! Local variables
447 real :: a_conv, b_conv ! Factors to multipy the a- and b-components by before they are written,
448 ! including sign changes to account for grid rotation [a A-1 ~> 1]
449
450 call lock_check(cs, a_desc)
451 call set_conversion_pair(a_conv, b_conv, cs%turns, conversion, scalar_pair)
452
453 if (modulo(cs%turns, 2) == 0) then ! This is the usual case.
454 call register_restart_field(a_ptr, a_desc, mandatory, cs, conversion=a_conv)
455 call register_restart_field(b_ptr, b_desc, mandatory, cs, conversion=b_conv)
456 else
457 call register_restart_field(b_ptr, a_desc, mandatory, cs, conversion=a_conv)
458 call register_restart_field(a_ptr, b_desc, mandatory, cs, conversion=b_conv)
459 endif
460end subroutine register_restart_pair_ptr3d
461
462
463!> Register a pair of rotationally equivalent 2d restart fields
464subroutine register_restart_pair_ptr4d(a_ptr, b_ptr, a_desc, b_desc, &
465 mandatory, CS, conversion, scalar_pair)
466 real, dimension(:,:,:,:), target, intent(in) :: a_ptr !< First field pointer
467 !! in arbitrary rescaled units [A ~> a]
468 real, dimension(:,:,:,:), target, intent(in) :: b_ptr !< Second field pointer
469 !! in arbitrary rescaled units [A ~> a]
470 type(vardesc), intent(in) :: a_desc !< First field descriptor
471 type(vardesc), intent(in) :: b_desc !< Second field descriptor
472 logical, intent(in) :: mandatory !< If true, abort if field is missing
473 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control structure
474 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
475 !! before it is written [a A-1 ~> 1], 1 by default.
476 logical, optional, intent(in) :: scalar_pair !< If true, the arrays describe a pair of
477 !! scalars, instead of vector components
478 !! whose signs change when rotated
479
480 ! Local variables
481 real :: a_conv, b_conv ! Factors to multipy the a- and b-components by before they are written,
482 ! including sign changes to account for grid rotation [a A-1 ~> 1]
483
484 call lock_check(cs, a_desc)
485 call set_conversion_pair(a_conv, b_conv, cs%turns, conversion, scalar_pair)
486
487 if (modulo(cs%turns, 2) == 0) then ! This is the usual case.
488 call register_restart_field(a_ptr, a_desc, mandatory, cs, conversion=a_conv)
489 call register_restart_field(b_ptr, b_desc, mandatory, cs, conversion=b_conv)
490 else
491 call register_restart_field(b_ptr, a_desc, mandatory, cs, conversion=a_conv)
492 call register_restart_field(a_ptr, b_desc, mandatory, cs, conversion=b_conv)
493 endif
494end subroutine register_restart_pair_ptr4d
495
496!> Set a pair of factors to multiply by the components of a vector when writing
497!! that include any sign changes needed to account for grid rotation.
498subroutine set_conversion_pair(u_conv, v_conv, turns, conversion, scalar_pair)
499 real, intent(out) :: u_conv !< A factor to multiply the u-component of a vector by before it is
500 !! written, including sign changes due to grid rotation [a A-1 ~> 1]
501 real, intent(out) :: v_conv !< A factor to multiply the u-component of a vector by before it is
502 !! written, including sign changes due to grid rotation [a A-1 ~> 1]
503 integer, intent(in) :: turns !< Number of quarter turns from input to model domain
504 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
505 !! before it is written [a A-1 ~> 1], 1 by default.
506 logical, optional, intent(in) :: scalar_pair !< If true, the arrays describe a pair of scalars,
507 !! instead of vector components whose signs change when rotated
508
509 ! Local variables
510 integer :: q_turns
511 logical :: scalars
512
513 u_conv = 1.0 ; v_conv = 1.0
514 if (present(conversion)) then
515 u_conv = conversion ; v_conv = conversion
516 endif
517
518 scalars = .false. ; if (present(scalar_pair)) scalars = scalar_pair
519 if (scalars) return
520
521 q_turns = modulo(turns, 4)
522 if (q_turns == 1) then
523 v_conv = -1.0*v_conv
524 elseif (q_turns == 2) then
525 u_conv = -1.0*u_conv ; v_conv = -1.0*v_conv
526 elseif (q_turns == 3) then
527 u_conv = -1.0*u_conv
528 endif
529
530end subroutine set_conversion_pair
531
532
533! The following provide alternate interfaces to register restarts.
534
535!> Register a 4-d field for restarts, providing the metadata as individual arguments
536subroutine register_restart_field_4d(f_ptr, name, mandatory, CS, longname, units, conversion, &
537 hor_grid, z_grid, t_grid, extra_axes)
538 real, dimension(:,:,:,:), &
539 target, intent(in) :: f_ptr !< A pointer to the field to be read or written
540 !! in arbitrary rescaled units [A ~> a]
541 character(len=*), intent(in) :: name !< variable name to be used in the restart file
542 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
543 !! successfully read from the restart file.
544 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
545 character(len=*), optional, intent(in) :: longname !< variable long name
546 character(len=*), optional, intent(in) :: units !< variable units
547 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
548 !! before it is written [a A-1 ~> 1], 1 by default.
549 character(len=*), optional, intent(in) :: hor_grid !< variable horizontal staggering, 'h' if absent
550 character(len=*), optional, intent(in) :: z_grid !< variable vertical staggering, 'L' if absent
551 character(len=*), optional, intent(in) :: t_grid !< time description: s, p, or 1, 's' if absent
552 type(axis_info), dimension(:), &
553 optional, intent(in) :: extra_axes !< dimensions other than space-time
554
555 type(vardesc) :: vd
556 character(len=32), dimension(:), allocatable :: dim_names
557 integer :: n, n_extradims
558
559 ! first 2 dimensions in dim_names are reserved for i,j
560 ! so extra_dimensions are shifted to index 3.
561 ! this is designed not to break the behavior in SIS2
562 ! (see register_restart_field_4d in SIS_restart.F90)
563 if (present(extra_axes)) then
564 n_extradims = size(extra_axes)
565 allocate(dim_names(n_extradims+2))
566 dim_names(1) = ""
567 dim_names(2) = ""
568 do n=3,n_extradims+2
569 dim_names(n) = extra_axes(n-2)%name
570 enddo
571 endif
572
573 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart: " // &
574 "register_restart_field_4d: Module must be initialized before "//&
575 "it is used to register "//trim(name))
576
577 call lock_check(cs, name=name)
578
579 if (present(extra_axes)) then
580 vd = var_desc(name, units=units, longname=longname, hor_grid=hor_grid, &
581 z_grid=z_grid, t_grid=t_grid, dim_names=dim_names, extra_axes=extra_axes)
582 else
583 vd = var_desc(name, units=units, longname=longname, hor_grid=hor_grid, &
584 z_grid=z_grid, t_grid=t_grid)
585 endif
586
587 call register_restart_field_ptr4d(f_ptr, vd, mandatory, cs, conversion)
588
589end subroutine register_restart_field_4d
590
591!> Register a 3-d field for restarts, providing the metadata as individual arguments
592subroutine register_restart_field_3d(f_ptr, name, mandatory, CS, longname, units, conversion, &
593 hor_grid, z_grid, t_grid, extra_axes)
594 real, dimension(:,:,:), &
595 target, intent(in) :: f_ptr !< A pointer to the field to be read or written
596 !! in arbitrary rescaled units [A ~> a]
597 character(len=*), intent(in) :: name !< variable name to be used in the restart file
598 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
599 !! successfully read from the restart file.
600 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
601 character(len=*), optional, intent(in) :: longname !< variable long name
602 character(len=*), optional, intent(in) :: units !< variable units
603 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
604 !! before it is written [a A-1 ~> 1], 1 by default.
605 character(len=*), optional, intent(in) :: hor_grid !< variable horizontal staggering, 'h' if absent
606 character(len=*), optional, intent(in) :: z_grid !< variable vertical staggering, 'L' if absent
607 character(len=*), optional, intent(in) :: t_grid !< time description: s, p, or 1, 's' if absent
608 type(axis_info), dimension(:), &
609 optional, intent(in) :: extra_axes !< dimensions other than space-time
610
611 type(vardesc) :: vd
612 character(len=32), dimension(:), allocatable :: dim_names
613 integer :: n, n_extradims
614
615 ! first 2 dimensions in dim_names are reserved for i,j
616 ! so extra_dimensions are shifted to index 3.
617 ! this is designed not to break the behavior in SIS2
618 ! (see register_restart_field_4d in SIS_restart.F90)
619 if (present(extra_axes)) then
620 n_extradims = size(extra_axes)
621 allocate(dim_names(n_extradims+2))
622 dim_names(1) = ""
623 dim_names(2) = ""
624 do n=3,n_extradims+2
625 dim_names(n) = extra_axes(n-2)%name
626 enddo
627 endif
628
629 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart: " // &
630 "register_restart_field_3d: Module must be initialized before "//&
631 "it is used to register "//trim(name))
632
633 call lock_check(cs, name=name)
634
635 if (present(extra_axes)) then
636 vd = var_desc(name, units=units, longname=longname, hor_grid=hor_grid, &
637 z_grid=z_grid, t_grid=t_grid, dim_names=dim_names, extra_axes=extra_axes)
638 else
639 vd = var_desc(name, units=units, longname=longname, hor_grid=hor_grid, &
640 z_grid=z_grid, t_grid=t_grid)
641 endif
642
643 call register_restart_field_ptr3d(f_ptr, vd, mandatory, cs, conversion)
644
645end subroutine register_restart_field_3d
646
647!> Register a 2-d field for restarts, providing the metadata as individual arguments
648subroutine register_restart_field_2d(f_ptr, name, mandatory, CS, longname, units, conversion, &
649 hor_grid, z_grid, t_grid)
650 real, dimension(:,:), &
651 target, intent(in) :: f_ptr !< A pointer to the field to be read or written
652 !! in arbitrary rescaled units [A ~> a]
653 character(len=*), intent(in) :: name !< variable name to be used in the restart file
654 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
655 !! successfully read from the restart file.
656 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
657 character(len=*), optional, intent(in) :: longname !< variable long name
658 character(len=*), optional, intent(in) :: units !< variable units
659 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
660 !! before it is written [a A-1 ~> 1], 1 by default.
661 character(len=*), optional, intent(in) :: hor_grid !< variable horizontal staggering, 'h' if absent
662 character(len=*), optional, intent(in) :: z_grid !< variable vertical staggering, '1' if absent
663 character(len=*), optional, intent(in) :: t_grid !< time description: s, p, or 1, 's' if absent
664
665 type(vardesc) :: vd
666 character(len=8) :: Zgrid
667
668 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart: " // &
669 "register_restart_field_2d: Module must be initialized before "//&
670 "it is used to register "//trim(name))
671
672 zgrid = '1' ; if (present(z_grid)) zgrid = z_grid
673
674 call lock_check(cs, name=name)
675
676 vd = var_desc(name, units=units, longname=longname, hor_grid=hor_grid, &
677 z_grid=zgrid, t_grid=t_grid)
678
679 call register_restart_field_ptr2d(f_ptr, vd, mandatory, cs, conversion)
680
681end subroutine register_restart_field_2d
682
683!> Register a 1-d field for restarts, providing the metadata as individual arguments
684subroutine register_restart_field_1d(f_ptr, name, mandatory, CS, longname, units, conversion, &
685 hor_grid, z_grid, t_grid)
686 real, dimension(:), target, intent(in) :: f_ptr !< A pointer to the field to be read or written
687 !! in arbitrary rescaled units [A ~> a]
688 character(len=*), intent(in) :: name !< variable name to be used in the restart file
689 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
690 !! successfully read from the restart file.
691 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
692 character(len=*), optional, intent(in) :: longname !< variable long name
693 character(len=*), optional, intent(in) :: units !< variable units
694 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
695 !! before it is written [a A-1 ~> 1], 1 by default.
696 character(len=*), optional, intent(in) :: hor_grid !< variable horizontal staggering, '1' if absent
697 character(len=*), optional, intent(in) :: z_grid !< variable vertical staggering, 'L' if absent
698 character(len=*), optional, intent(in) :: t_grid !< time description: s, p, or 1, 's' if absent
699
700 type(vardesc) :: vd
701 character(len=8) :: hgrid
702
703 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart: " // &
704 "register_restart_field_3d: Module must be initialized before "//&
705 "it is used to register "//trim(name))
706
707 hgrid = '1' ; if (present(hor_grid)) hgrid = hor_grid
708
709 call lock_check(cs, name=name)
710
711 vd = var_desc(name, units=units, longname=longname, hor_grid=hgrid, &
712 z_grid=z_grid, t_grid=t_grid)
713
714 call register_restart_field_ptr1d(f_ptr, vd, mandatory, cs, conversion)
715
716end subroutine register_restart_field_1d
717
718!> Register a 0-d field for restarts, providing the metadata as individual arguments
719subroutine register_restart_field_0d(f_ptr, name, mandatory, CS, longname, units, conversion, &
720 t_grid)
721 real, target, intent(in) :: f_ptr !< A pointer to the field to be read or written
722 !! in arbitrary rescaled units [A ~> a]
723 character(len=*), intent(in) :: name !< variable name to be used in the restart file
724 logical, intent(in) :: mandatory !< If true, the run will abort if this field is not
725 !! successfully read from the restart file.
726 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
727 character(len=*), optional, intent(in) :: longname !< variable long name
728 character(len=*), optional, intent(in) :: units !< variable units
729 real, optional, intent(in) :: conversion !< A factor to multiply a restart field by
730 !! before it is written [a A-1 ~> 1], 1 by default.
731 character(len=*), optional, intent(in) :: t_grid !< time description: s, p, or 1, 's' if absent
732
733 type(vardesc) :: vd
734
735 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart: " // &
736 "register_restart_field_0d: Module must be initialized before "//&
737 "it is used to register "//trim(name))
738
739 call lock_check(cs, name=name)
740
741 vd = var_desc(name, units=units, longname=longname, hor_grid='1', &
742 z_grid='1', t_grid=t_grid)
743
744 call register_restart_field_ptr0d(f_ptr, vd, mandatory, cs, conversion)
745
746end subroutine register_restart_field_0d
747
748
749!> query_initialized_name determines whether a named field has been successfully
750!! read from a restart file or has otherwise been recorded as being initialized.
751function query_initialized_name(name, CS) result(query_initialized)
752 character(len=*), intent(in) :: name !< The name of the field that is being queried
753 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
754 logical :: query_initialized
755
756 integer :: m, n
757
758 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
759 "query_initialized: Module must be initialized before it is used.")
760
761 if (cs%novars > cs%max_fields) call restart_error(cs)
762
763 query_initialized = .false.
764 n = cs%novars+1
765 do m=1,cs%novars
766 if (trim(name) == cs%restart_field(m)%var_name) then
767 if (cs%restart_field(m)%initialized) query_initialized = .true.
768 n = m ; exit
769 endif
770 enddo
771 if ((n==cs%novars+1) .and. (is_root_pe())) &
772 call mom_error(note,"MOM_restart: Unknown restart variable "//name// &
773 " queried for initialization.")
774
775 if ((is_root_pe()) .and. query_initialized) &
776 call mom_error(note,"MOM_restart: "//name// &
777 " initialization confirmed by name.")
778
779end function query_initialized_name
780
781!> Indicate whether the field pointed to by f_ptr has been initialized from a restart file.
782function query_initialized_0d(f_ptr, CS) result(query_initialized)
783 real, target, intent(in) :: f_ptr !< A pointer to the field that is being queried [arbitrary]
784 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
785 logical :: query_initialized
786
787 integer :: m, n
788
789 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
790 "query_initialized: Module must be initialized before it is used.")
791
792 if (cs%novars > cs%max_fields) call restart_error(cs)
793
794 query_initialized = .false.
795 n = cs%novars+1
796 do m=1,cs%novars
797 if (associated(cs%var_ptr0d(m)%p,f_ptr)) then
798 if (cs%restart_field(m)%initialized) query_initialized = .true.
799 n = m ; exit
800 endif
801 enddo
802
803end function query_initialized_0d
804
805!> Indicate whether the field pointed to by f_ptr has been initialized from a restart file.
806function query_initialized_1d(f_ptr, CS) result(query_initialized)
807 real, dimension(:), target, intent(in) :: f_ptr !< A pointer to the field that is being queried [arbitrary]
808 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
809 logical :: query_initialized
810
811 integer :: m, n
812
813 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
814 "query_initialized: Module must be initialized before it is used.")
815
816 if (cs%novars > cs%max_fields) call restart_error(cs)
817
818 query_initialized = .false.
819 n = cs%novars+1
820 do m=1,cs%novars
821 if (associated(cs%var_ptr1d(m)%p,f_ptr)) then
822 if (cs%restart_field(m)%initialized) query_initialized = .true.
823 n = m ; exit
824 endif
825 enddo
826
827end function query_initialized_1d
828
829!> Indicate whether the field pointed to by f_ptr has been initialized from a restart file.
830function query_initialized_2d(f_ptr, CS) result(query_initialized)
831 real, dimension(:,:), &
832 target, intent(in) :: f_ptr !< A pointer to the field that is being queried [arbitrary]
833 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
834 logical :: query_initialized
835
836 integer :: m, n
837
838 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
839 "query_initialized: Module must be initialized before it is used.")
840
841 if (cs%novars > cs%max_fields) call restart_error(cs)
842
843 query_initialized = .false.
844 n = cs%novars+1
845 do m=1,cs%novars
846 if (associated(cs%var_ptr2d(m)%p,f_ptr)) then
847 if (cs%restart_field(m)%initialized) query_initialized = .true.
848 n = m ; exit
849 endif
850 enddo
851
852end function query_initialized_2d
853
854!> Indicate whether the field pointed to by f_ptr has been initialized from a restart file.
855function query_initialized_3d(f_ptr, CS) result(query_initialized)
856 real, dimension(:,:,:), &
857 target, intent(in) :: f_ptr !< A pointer to the field that is being queried [arbitrary]
858 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
859 logical :: query_initialized
860
861 integer :: m, n
862
863 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
864 "query_initialized: Module must be initialized before it is used.")
865
866 if (cs%novars > cs%max_fields) call restart_error(cs)
867
868 query_initialized = .false.
869 n = cs%novars+1
870 do m=1,cs%novars
871 if (associated(cs%var_ptr3d(m)%p,f_ptr)) then
872 if (cs%restart_field(m)%initialized) query_initialized = .true.
873 n = m ; exit
874 endif
875 enddo
876
877end function query_initialized_3d
878
879!> Indicate whether the field pointed to by f_ptr has been initialized from a restart file.
880function query_initialized_4d(f_ptr, CS) result(query_initialized)
881 real, dimension(:,:,:,:), &
882 target, intent(in) :: f_ptr !< A pointer to the field that is being queried [arbitrary]
883 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
884 logical :: query_initialized
885
886 integer :: m, n
887
888 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
889 "query_initialized: Module must be initialized before it is used.")
890
891 if (cs%novars > cs%max_fields) call restart_error(cs)
892
893 query_initialized = .false.
894 n = cs%novars+1
895 do m=1,cs%novars
896 if (associated(cs%var_ptr4d(m)%p,f_ptr)) then
897 if (cs%restart_field(m)%initialized) query_initialized = .true.
898 n = m ; exit
899 endif
900 enddo
901
902end function query_initialized_4d
903
904!> Indicate whether the field stored in f_ptr or with the specified variable
905!! name has been initialized from a restart file.
906function query_initialized_0d_name(f_ptr, name, CS) result(query_initialized)
907 real, target, intent(in) :: f_ptr !< The field that is being queried [arbitrary]
908 character(len=*), intent(in) :: name !< The name of the field that is being queried
909 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
910 logical :: query_initialized
911
912 integer :: m, n
913
914 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
915 "query_initialized: Module must be initialized before it is used.")
916
917 if (cs%novars > cs%max_fields) call restart_error(cs)
918
919 query_initialized = .false.
920 n = cs%novars+1
921 do m=1,cs%novars
922 if (associated(cs%var_ptr0d(m)%p,f_ptr)) then
923 if (cs%restart_field(m)%initialized) query_initialized = .true.
924 n = m ; exit
925 endif
926 enddo
927 if (n==cs%novars+1) then
928 if (is_root_pe()) &
929 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
930 "probably because of the suspect comparison of pointers by ASSOCIATED.")
931 query_initialized = query_initialized_name(name, cs)
932 endif
933
934end function query_initialized_0d_name
935
936!> Indicate whether the field stored in f_ptr or with the specified variable
937!! name has been initialized from a restart file.
938function query_initialized_1d_name(f_ptr, name, CS) result(query_initialized)
939 real, dimension(:), &
940 target, intent(in) :: f_ptr !< The field that is being queried [arbitrary]
941 character(len=*), intent(in) :: name !< The name of the field that is being queried
942 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
943 logical :: query_initialized
944
945 integer :: m, n
946
947 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
948 "query_initialized: Module must be initialized before it is used.")
949
950 if (cs%novars > cs%max_fields) call restart_error(cs)
951
952 query_initialized = .false.
953 n = cs%novars+1
954 do m=1,cs%novars
955 if (associated(cs%var_ptr1d(m)%p,f_ptr)) then
956 if (cs%restart_field(m)%initialized) query_initialized = .true.
957 n = m ; exit
958 endif
959 enddo
960 if (n==cs%novars+1) then
961 if (is_root_pe()) &
962 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
963 "probably because of the suspect comparison of pointers by ASSOCIATED.")
964 query_initialized = query_initialized_name(name, cs)
965 endif
966
967end function query_initialized_1d_name
968
969!> Indicate whether the field stored in f_ptr or with the specified variable
970!! name has been initialized from a restart file.
971function query_initialized_2d_name(f_ptr, name, CS) result(query_initialized)
972 real, dimension(:,:), &
973 target, intent(in) :: f_ptr !< The field that is being queried [arbitrary]
974 character(len=*), intent(in) :: name !< The name of the field that is being queried
975 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
976 logical :: query_initialized
977
978 integer :: m, n
979
980 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
981 "query_initialized: Module must be initialized before it is used.")
982
983 if (cs%novars > cs%max_fields) call restart_error(cs)
984
985 query_initialized = .false.
986 n = cs%novars+1
987 do m=1,cs%novars
988 if (associated(cs%var_ptr2d(m)%p,f_ptr)) then
989 if (cs%restart_field(m)%initialized) query_initialized = .true.
990 n = m ; exit
991 endif
992 enddo
993 if (n==cs%novars+1) then
994 if (is_root_pe()) &
995 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
996 "probably because of the suspect comparison of pointers by ASSOCIATED.")
997 query_initialized = query_initialized_name(name, cs)
998 endif
999
1000end function query_initialized_2d_name
1001
1002!> Indicate whether the field stored in f_ptr or with the specified variable
1003!! name has been initialized from a restart file.
1004function query_initialized_3d_name(f_ptr, name, CS) result(query_initialized)
1005 real, dimension(:,:,:), &
1006 target, intent(in) :: f_ptr !< The field that is being queried [arbitrary]
1007 character(len=*), intent(in) :: name !< The name of the field that is being queried
1008 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
1009 logical :: query_initialized
1010
1011 integer :: m, n
1012
1013 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1014 "query_initialized: Module must be initialized before it is used.")
1015
1016 if (cs%novars > cs%max_fields) call restart_error(cs)
1017
1018 query_initialized = .false.
1019 n = cs%novars+1
1020 do m=1,cs%novars
1021 if (associated(cs%var_ptr3d(m)%p,f_ptr)) then
1022 if (cs%restart_field(m)%initialized) query_initialized = .true.
1023 n = m ; exit
1024 endif
1025 enddo
1026 if (n==cs%novars+1) then
1027 if (is_root_pe()) &
1028 call mom_error(note, "MOM_restart: Unable to find "//name//" queried by pointer, "//&
1029 "possibly because of the suspect comparison of pointers by ASSOCIATED.")
1030 query_initialized = query_initialized_name(name, cs)
1031 endif
1032
1033end function query_initialized_3d_name
1034
1035!> Indicate whether the field stored in f_ptr or with the specified variable
1036!! name has been initialized from a restart file.
1037function query_initialized_4d_name(f_ptr, name, CS) result(query_initialized)
1038 real, dimension(:,:,:,:), &
1039 target, intent(in) :: f_ptr !< The field that is being queried [arbitrary]
1040 character(len=*), intent(in) :: name !< The name of the field that is being queried
1041 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
1042 logical :: query_initialized
1043
1044 integer :: m, n
1045
1046 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1047 "query_initialized: Module must be initialized before it is used.")
1048
1049 if (cs%novars > cs%max_fields) call restart_error(cs)
1050
1051 query_initialized = .false.
1052 n = cs%novars+1
1053 do m=1,cs%novars
1054 if (associated(cs%var_ptr4d(m)%p,f_ptr)) then
1055 if (cs%restart_field(m)%initialized) query_initialized = .true.
1056 n = m ; exit
1057 endif
1058 enddo
1059 if (n==cs%novars+1) then
1060 if (is_root_pe()) &
1061 call mom_error(note, "MOM_restart: Unable to find "//name//" queried by pointer, "//&
1062 "possibly because of the suspect comparison of pointers by ASSOCIATED.")
1063 query_initialized = query_initialized_name(name, cs)
1064 endif
1065
1066end function query_initialized_4d_name
1067
1068!> set_initialized_name records that a named field has been initialized.
1069subroutine set_initialized_name(name, CS)
1070 character(len=*), intent(in) :: name !< The name of the field that is being set
1071 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
1072
1073 integer :: m
1074
1075 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1076 "set_initialized: Module must be initialized before it is used.")
1077
1078 do m=1,cs%novars ; if (trim(name) == trim(cs%restart_field(m)%var_name)) then
1079 cs%restart_field(m)%initialized = .true. ; exit
1080 endif ; enddo
1081
1082 if ((m==cs%novars+1) .and. (is_root_pe())) &
1083 call mom_error(note,"MOM_restart: Unknown restart variable "//name// &
1084 " used in set_initialized call.")
1085
1086end subroutine set_initialized_name
1087
1088!> Record that the array in f_ptr with the given name has been initialized.
1089subroutine set_initialized_0d_name(f_ptr, name, CS)
1090 real, target, intent(in) :: f_ptr !< The variable that has been initialized [arbitrary]
1091 character(len=*), intent(in) :: name !< The name of the field that has been initialized
1092 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
1093
1094 integer :: m
1095
1096 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1097 "set_initialized: Module must be initialized before it is used.")
1098
1099 do m=1,cs%novars ; if (associated(cs%var_ptr0d(m)%p,f_ptr)) then
1100 cs%restart_field(m)%initialized = .true. ; exit
1101 endif ; enddo
1102
1103 if (m==cs%novars+1) then
1104 if (is_root_pe()) &
1105 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
1106 "probably because of the suspect comparison of pointers by ASSOCIATED.")
1107 call set_initialized_name(name, cs)
1108 endif
1109
1110end subroutine set_initialized_0d_name
1111
1112!> Record that the array in f_ptr with the given name has been initialized.
1113subroutine set_initialized_1d_name(f_ptr, name, CS)
1114 real, dimension(:), &
1115 target, intent(in) :: f_ptr !< The array that has been initialized [arbitrary]
1116 character(len=*), intent(in) :: name !< The name of the field that has been initialized
1117 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
1118
1119 integer :: m
1120
1121 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1122 "set_initialized: Module must be initialized before it is used.")
1123
1124 do m=1,cs%novars ; if (associated(cs%var_ptr1d(m)%p,f_ptr)) then
1125 cs%restart_field(m)%initialized = .true. ; exit
1126 endif ; enddo
1127
1128 if (m==cs%novars+1) then
1129 if (is_root_pe()) &
1130 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
1131 "probably because of the suspect comparison of pointers by ASSOCIATED.")
1132 call set_initialized_name(name, cs)
1133 endif
1134
1135end subroutine set_initialized_1d_name
1136
1137!> Record that the array in f_ptr with the given name has been initialized.
1138subroutine set_initialized_2d_name(f_ptr, name, CS)
1139 real, dimension(:,:), &
1140 target, intent(in) :: f_ptr !< The array that has been initialized [arbitrary]
1141 character(len=*), intent(in) :: name !< The name of the field that has been initialized
1142 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
1143
1144 integer :: m
1145
1146 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1147 "set_initialized: Module must be initialized before it is used.")
1148
1149 do m=1,cs%novars ; if (associated(cs%var_ptr2d(m)%p,f_ptr)) then
1150 cs%restart_field(m)%initialized = .true. ; exit
1151 endif ; enddo
1152
1153 if (m==cs%novars+1) then
1154 if (is_root_pe()) &
1155 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
1156 "probably because of the suspect comparison of pointers by ASSOCIATED.")
1157 call set_initialized_name(name, cs)
1158 endif
1159
1160end subroutine set_initialized_2d_name
1161
1162!> Record that the array in f_ptr with the given name has been initialized.
1163subroutine set_initialized_3d_name(f_ptr, name, CS)
1164 real, dimension(:,:,:), &
1165 target, intent(in) :: f_ptr !< The array that has been initialized [arbitrary]
1166 character(len=*), intent(in) :: name !< The name of the field that has been initialized
1167 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
1168
1169 integer :: m
1170
1171 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1172 "set_initialized: Module must be initialized before it is used.")
1173
1174 do m=1,cs%novars ; if (associated(cs%var_ptr3d(m)%p,f_ptr)) then
1175 cs%restart_field(m)%initialized = .true. ; exit
1176 endif ; enddo
1177
1178 if (m==cs%novars+1) then
1179 if (is_root_pe()) &
1180 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
1181 "probably because of the suspect comparison of pointers by ASSOCIATED.")
1182 call set_initialized_name(name, cs)
1183 endif
1184
1185end subroutine set_initialized_3d_name
1186
1187!> Record that the array in f_ptr with the given name has been initialized.
1188subroutine set_initialized_4d_name(f_ptr, name, CS)
1189 real, dimension(:,:,:,:), &
1190 target, intent(in) :: f_ptr !< The array that has been initialized [arbitrary]
1191 character(len=*), intent(in) :: name !< The name of the field that has been initialized
1192 type(mom_restart_cs), intent(inout) :: CS !< MOM restart control struct
1193
1194 integer :: m
1195
1196 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1197 "set_initialized: Module must be initialized before it is used.")
1198
1199 do m=1,cs%novars ; if (associated(cs%var_ptr4d(m)%p,f_ptr)) then
1200 cs%restart_field(m)%initialized = .true. ; exit
1201 endif ; enddo
1202
1203 if (m==cs%novars+1) then
1204 if (is_root_pe()) &
1205 call mom_error(note,"MOM_restart: Unable to find "//name//" queried by pointer, "//&
1206 "probably because of the suspect comparison of pointers by ASSOCIATED.")
1207 call set_initialized_name(name, cs)
1208 endif
1209
1210end subroutine set_initialized_4d_name
1211
1212
1213!====================== only_read_from_restarts variants =======================
1214
1215!> Try to read a named 4-d field from the restart files
1216subroutine only_read_restart_field_4d(varname, f_ptr, G, CS, position, filename, directory, success, scale)
1217 character(len=*), intent(in) :: varname !< The variable name to be used in the restart file
1218 real, dimension(:,:,:,:), intent(inout) :: f_ptr !< The array for the field to be read
1219 !! in arbitrary rescaled units [A ~> a]
1220 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1221 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
1222 integer, optional, intent(in) :: position !< A coded integer indicating the horizontal
1223 !! position of this variable
1224 character(len=*), optional, intent(in) :: filename !< The list of restart file names or a single
1225 !! character 'r' to read automatically named files
1226 character(len=*), optional, intent(in) :: directory !< The directory in which to seek restart files.
1227 logical, optional, intent(out) :: success !< True if the field was read successfully
1228 real, optional, intent(in) :: scale !< A factor by which the field will be scaled
1229 !! [A a-1 ~> 1] to convert from the units in
1230 !! the file to the internal units of this field
1231
1232 ! Local variables
1233 character(len=:), allocatable :: file_path ! The full path to the file with the variable
1234 logical :: found ! True if the variable was found.
1235 logical :: is_global ! True if the variable is in a global file.
1236
1237 found = find_var_in_restart_files(varname, g, cs, file_path, filename, directory, is_global)
1238
1239 if (found) then
1240 call mom_read_data(file_path, varname, f_ptr, g%domain, timelevel=1, position=position, &
1241 scale=scale, global_file=is_global)
1242 endif
1243 if (present(success)) success = found
1244
1245end subroutine only_read_restart_field_4d
1246
1247!> Try to read a named 3-d field from the restart files
1248subroutine only_read_restart_field_3d(varname, f_ptr, G, CS, position, filename, directory, success, scale)
1249 character(len=*), intent(in) :: varname !< The variable name to be used in the restart file
1250 real, dimension(:,:,:), intent(inout) :: f_ptr !< The array for the field to be read
1251 !! in arbitrary rescaled units [A ~> a]
1252 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1253 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
1254 integer, optional, intent(in) :: position !< A coded integer indicating the horizontal
1255 !! position of this variable
1256 character(len=*), optional, intent(in) :: filename !< The list of restart file names or a single
1257 !! character 'r' to read automatically named files
1258 character(len=*), optional, intent(in) :: directory !< The directory in which to seek restart files.
1259 logical, optional, intent(out) :: success !< True if the field was read successfully
1260 real, optional, intent(in) :: scale !< A factor by which the field will be scaled
1261 !! [A a-1 ~> 1] to convert from the units in
1262 !! the file to the internal units of this field
1263
1264 ! Local variables
1265 character(len=:), allocatable :: file_path ! The full path to the file with the variable
1266 logical :: found ! True if the variable was found.
1267 logical :: is_global ! True if the variable is in a global file.
1268
1269 found = find_var_in_restart_files(varname, g, cs, file_path, filename, directory, is_global)
1270
1271 if (found) then
1272 call mom_read_data(file_path, varname, f_ptr, g%domain, timelevel=1, position=position, &
1273 scale=scale, global_file=is_global)
1274 endif
1275 if (present(success)) success = found
1276
1277end subroutine only_read_restart_field_3d
1278
1279!> Try to read a named 2-d field from the restart files
1280subroutine only_read_restart_field_2d(varname, f_ptr, G, CS, position, filename, directory, success, scale)
1281 character(len=*), intent(in) :: varname !< The variable name to be used in the restart file
1282 real, dimension(:,:), intent(inout) :: f_ptr !< The array for the field to be read
1283 !! in arbitrary rescaled units [A ~> a]
1284 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1285 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
1286 integer, optional, intent(in) :: position !< A coded integer indicating the horizontal
1287 !! position of this variable
1288 character(len=*), optional, intent(in) :: filename !< The list of restart file names or a single
1289 !! character 'r' to read automatically named files
1290 character(len=*), optional, intent(in) :: directory !< The directory in which to seek restart files.
1291 logical, optional, intent(out) :: success !< True if the field was read successfully
1292 real, optional, intent(in) :: scale !< A factor by which the field will be scaled
1293 !! [A a-1 ~> 1] to convert from the units in
1294 !! the file to the internal units of this field
1295
1296 ! Local variables
1297 character(len=:), allocatable :: file_path ! The full path to the file with the variable
1298 logical :: found ! True if the variable was found.
1299 logical :: is_global ! True if the variable is in a global file.
1300
1301 found = find_var_in_restart_files(varname, g, cs, file_path, filename, directory, is_global)
1302
1303 if (found) then
1304 call mom_read_data(file_path, varname, f_ptr, g%domain, timelevel=1, position=position, &
1305 scale=scale, global_file=is_global)
1306 endif
1307 if (present(success)) success = found
1308
1309end subroutine only_read_restart_field_2d
1310
1311
1312!> Try to read a named 3-d field from the restart files
1313subroutine only_read_restart_pair_3d(a_ptr, b_ptr, a_name, b_name, G, CS, &
1314 stagger, filename, directory, success, scale)
1315 real, dimension(:,:,:), intent(inout) :: a_ptr !< The array for the first field to be read
1316 !! in arbitrary rescaled units [A ~> a]
1317 real, dimension(:,:,:), intent(inout) :: b_ptr !< The array for the second field to be read
1318 !! in arbitrary rescaled units [A ~> a]
1319 character(len=*), intent(in) :: a_name !< The first variable name to be used in the restart file
1320 character(len=*), intent(in) :: b_name !< The second variable name to be used in the restart file
1321 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
1322 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
1323 integer, optional, intent(in) :: stagger !< A coded integer indicating the horizontal
1324 !! position of this pair of variables
1325 character(len=*), optional, intent(in) :: filename !< The list of restart file names or a single
1326 !! character 'r' to read automatically named files
1327 character(len=*), optional, intent(in) :: directory !< The directory in which to seek restart files.
1328 logical, optional, intent(out) :: success !< True if the field was read successfully
1329 real, optional, intent(in) :: scale !< A factor by which the fields will be scaled
1330 !! [A a-1 ~> 1] to convert from the units in
1331 !! the file to the internal units of this field
1332
1333 ! Local variables
1334 character(len=:), allocatable :: file_path_a ! The full path to the file with the first variable
1335 character(len=:), allocatable :: file_path_b ! The full path to the file with the second variable
1336 integer :: a_pos, b_pos ! A coded position for the two variables.
1337 logical :: a_found, b_found ! True if the variables were found.
1338 logical :: global_a, global_b ! True if the variables are in global files.
1339
1340 a_found = find_var_in_restart_files(a_name, g, cs, file_path_a, filename, directory, global_a)
1341 b_found = find_var_in_restart_files(b_name, g, cs, file_path_b, filename, directory, global_b)
1342
1343 a_pos = east_face ; b_pos = north_face
1344 if (present(stagger)) then ; select case (stagger)
1345 case (agrid) ; a_pos = center ; b_pos = center
1346 case (bgrid_ne) ; a_pos = corner ; b_pos = corner
1347 case (cgrid_ne) ; a_pos = east_face ; b_pos = north_face
1348 case default ; a_pos = east_face ; b_pos = north_face
1349 end select ; endif
1350
1351 if (a_found .and. b_found) then
1352 call mom_read_data(file_path_a, a_name, a_ptr, g%domain, timelevel=1, position=a_pos, &
1353 scale=scale, global_file=global_b, file_may_be_4d=.true.)
1354 call mom_read_data(file_path_b, b_name, b_ptr, g%domain, timelevel=1, position=b_pos, &
1355 scale=scale, global_file=global_b, file_may_be_4d=.true.)
1356 endif
1357 if (present(success)) success = (a_found .and. b_found)
1358
1359end subroutine only_read_restart_pair_3d
1360
1361!> Return an indication of whether the named variable is in the restart files, and provide the full path
1362!! to the restart file in which a variable is found.
1363function find_var_in_restart_files(varname, G, CS, file_path, filename, directory, is_global) result (found)
1364 character(len=*), intent(in) :: varname !< The variable name to be used in the restart file
1365 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
1366 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
1367 character(len=:), allocatable, intent(out) :: file_path !< The full path to the file in which the
1368 !! variable is found
1369 character(len=*), optional, intent(in) :: filename !< The list of restart file names or a single
1370 !! character 'r' to read automatically named files
1371 character(len=*), optional, intent(in) :: directory !< The directory in which to seek restart files.
1372 logical, optional, intent(out) :: is_global !< True if the file is global.
1373 logical :: found !< True if the named variable was found in the restart files.
1374
1375 ! Local variables
1376 character(len=240), allocatable, dimension(:) :: file_paths ! The possible file names.
1377 character(len=:), allocatable :: dir ! The directory to read from.
1378 character(len=:), allocatable :: fname ! The list of file names.
1379 logical, allocatable, dimension(:) :: global_file ! True if the file is global
1380 integer :: n, num_files
1381
1382 dir = "./INPUT/" ; if (present(directory)) dir = trim(directory)
1383
1384 ! Set the default return values.
1385 found = .false.
1386 file_path = ""
1387 if (present(is_global)) is_global = .false.
1388
1389 fname = 'r'
1390 if (present(filename)) then
1391 if (.not.((len_trim(filename) == 1) .and. (filename(1:1) == 'F'))) fname = filename
1392 endif
1393
1394 num_files = get_num_restart_files(fname, dir, g, cs)
1395 if (num_files == 0) return
1396 allocate(file_paths(num_files), global_file(num_files))
1397 num_files = open_restart_units(fname, dir, g, cs, file_paths=file_paths, global_files=global_file)
1398
1399 do n=1,num_files ; if (field_exists(file_paths(n), varname, mom_domain=g%domain)) then
1400 found = .true.
1401 file_path = file_paths(n)
1402 if (present(is_global)) is_global = global_file(n)
1403 exit
1404 endif ; enddo
1405
1406 deallocate(file_paths, global_file)
1407
1408end function find_var_in_restart_files
1409
1410!====================== end of the only_read_from_restarts variants =======================
1411
1412
1413!> Copy the restart variable with the specified name into a 3-d array, perhaps after rotation
1414subroutine copy_restart_var_3d(var, name, CS, unrotate)
1415 real, dimension(:,:,:), intent(inout) :: var !< The field that is being copied [arbitrary]
1416 character(len=*), intent(in) :: name !< The name of the field that is being copied
1417 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
1418 logical, optional, intent(in) :: unrotate !< If present and true, the output is on an unrotated grid.
1419
1420 logical :: keep_rotation
1421 character(len=256) :: size_msg !< The array sizes
1422 integer :: m, n
1423
1424 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1425 "query_initialized: Module must be initialized before it is used.")
1426
1427 if (cs%novars > cs%max_fields) call restart_error(cs)
1428
1429 keep_rotation = .true. ; if (present(unrotate)) keep_rotation = .not.unrotate
1430
1431 n = cs%novars+1
1432 do m=1,cs%novars
1433 if (trim(name) == cs%restart_field(m)%var_name) then
1434 if (.not.associated(cs%var_ptr3d(m)%p)) then
1435 call mom_error(fatal, "MOM_restart: copy_restart_var(_3d) "//&
1436 "attempted to copy restart variable "//name//" with the wrong rank.")
1437 elseif (cs%restart_field(m)%initialized) then
1438 if (cs%turns == 0 .or. keep_rotation) then
1439 if ( size_mismatch_3d(var, cs%var_ptr3d(m)%p, cs%turns, size_msg) ) &
1440 call mom_error(fatal, "MOM_restart: copy_restart_var(_3d) "//&
1441 "attempted to copy restart variable "//name//" with the wrong sizes, "//trim(size_msg))
1442
1443 var(:,:,:) = cs%var_ptr3d(m)%p(:,:,:)
1444 else
1445 call rotate_array(cs%var_ptr3d(m)%p, -cs%turns, var)
1446 endif
1447 else
1448 call mom_error(note, "MOM_restart: copy_restart_var(_3d) "//&
1449 "attempted to copy uninitialized restart variable "//name//".")
1450 endif
1451 n = m ; exit
1452 endif
1453 enddo
1454 if ((n==cs%novars+1) .and. (is_root_pe())) &
1455 call mom_error(note, "MOM_restart: copy_restart_var(_3d) "//&
1456 "attempted to copy unknown restart variable "//name//".")
1457
1458end subroutine copy_restart_var_3d
1459
1460
1461!> Copy the restart vector component variables with the specified names into a pair
1462!! of 3-d arrays, perhaps after rotation
1463subroutine copy_restart_vector_3d(u_var, v_var, u_name, v_name, CS, unrotate, scalar_pair)
1464 real, dimension(:,:,:), intent(inout) :: u_var !< The u-component of the field that is being copied [arbitrary]
1465 real, dimension(:,:,:), intent(inout) :: v_var !< The u-component of the field that is being copied [arbitrary]
1466 character(len=*), intent(in) :: u_name !< The name of the u-component of the field that is being copied
1467 character(len=*), intent(in) :: v_name !< The name of the v-component of the field that is being copied
1468 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
1469 logical, optional, intent(in) :: unrotate !< If present and true, the output is on an unrotated grid.
1470 logical, optional, intent(in) :: scalar_pair !< If true, the arrays describe a pair of
1471 !! scalars, instead of vector components
1472 !! whose signs change when rotated
1473
1474 logical :: keep_rotation, scalars
1475 character(len=256) :: size_msg !< The array sizes
1476 integer :: m, n_u, n_v
1477
1478 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1479 "query_initialized: Module must be initialized before it is used.")
1480
1481 if (cs%novars > cs%max_fields) call restart_error(cs)
1482
1483 keep_rotation = .true. ; if (present(unrotate)) keep_rotation = .not.unrotate
1484
1485 n_u = cs%novars+1 ; n_v = cs%novars+1
1486 do m=1,cs%novars
1487 if (trim(u_name) == cs%restart_field(m)%var_name) then
1488 if (.not.associated(cs%var_ptr3d(m)%p)) then
1489 call mom_error(fatal, "MOM_restart: copy_restart_vector(_3d) "//&
1490 "attempted to copy restart variable "//trim(u_name)//" with the wrong rank.")
1491 elseif (cs%restart_field(m)%initialized) then
1492 n_u = m
1493 else
1494 call mom_error(note, "MOM_restart: copy_restart_vector(_3d) "//&
1495 "attempted to copy uninitialized restart variable "//trim(u_name)//".")
1496 n_u = -1
1497 endif
1498 endif
1499 if (trim(v_name) == cs%restart_field(m)%var_name) then
1500 if (.not.associated(cs%var_ptr3d(m)%p)) then
1501 call mom_error(fatal, "MOM_restart: copy_restart_vector(_3d) "//&
1502 "attempted to copy restart variable "//trim(v_name)//" with the wrong rank.")
1503 elseif (cs%restart_field(m)%initialized) then
1504 n_v = m
1505 else
1506 call mom_error(note, "MOM_restart: copy_restart_vector(_3d) "//&
1507 "attempted to copy uninitialized restart variable "//trim(v_name)//".")
1508 n_v = -1
1509 endif
1510 endif
1511 enddo
1512 if ((n_u==cs%novars+1) .and. (is_root_pe())) &
1513 call mom_error(note, "MOM_restart: copy_restart_vector(_3d) "//&
1514 "attempted to copy unknown restart variable "//trim(u_name)//".")
1515 if ((n_v==cs%novars+1) .and. (is_root_pe())) &
1516 call mom_error(note, "MOM_restart: copy_restart_vector(_3d) "//&
1517 "attempted to copy unknown restart variable "//trim(v_name)//".")
1518
1519 if ((n_u>0) .and. (n_u<=cs%novars) .and. (n_v>0) .and. (n_v<=cs%novars)) then
1520 ! Now actually update the vector.
1521 if ( size_mismatch_3d(u_var, cs%var_ptr3d(n_u)%p, cs%turns, size_msg) ) &
1522 call mom_error(fatal, "MOM_restart: copy_restart_vector(_3d) "//&
1523 "attempted to copy restart variable "//trim(u_name)//" with the wrong sizes, "//trim(size_msg))
1524 if ( size_mismatch_3d(v_var, cs%var_ptr3d(n_v)%p, cs%turns, size_msg) ) &
1525 call mom_error(fatal, "MOM_restart: copy_restart_vector(_3d) "//&
1526 "attempted to copy restart variable "//trim(v_name)//" with the wrong sizes, "//trim(size_msg))
1527
1528 if (cs%turns == 0 .or. keep_rotation) then
1529 u_var(:,:,:) = cs%var_ptr3d(n_u)%p(:,:,:)
1530 v_var(:,:,:) = cs%var_ptr3d(n_v)%p(:,:,:)
1531 else
1532 scalars = .false. ; if (present(scalar_pair)) scalars = scalar_pair
1533 if ((modulo(cs%turns, 2) == 0) .and. scalars) then
1534 call rotate_array_pair(cs%var_ptr3d(n_u)%p, cs%var_ptr3d(n_v)%p, -cs%turns, u_var, v_var)
1535 elseif (modulo(cs%turns, 2) == 0) then
1536 call rotate_vector(cs%var_ptr3d(n_u)%p, cs%var_ptr3d(n_v)%p, -cs%turns, u_var, v_var)
1537 elseif (scalars) then ! This is less common
1538 call rotate_array_pair(cs%var_ptr3d(n_v)%p, cs%var_ptr3d(n_u)%p, -cs%turns, u_var, v_var)
1539 else
1540 call rotate_vector(cs%var_ptr3d(n_v)%p, cs%var_ptr3d(n_u)%p, -cs%turns, u_var, v_var)
1541 endif
1542 endif
1543 endif
1544
1545end subroutine copy_restart_vector_3d
1546
1547!> Indicate if two 3-d arrays are not of the same size after rotation is considered.
1548logical function size_mismatch_3d(var_a, var_b, turns, size_msg)
1549 real, intent(in) :: var_a(:,:,:) !< The first field being compared
1550 real, intent(in) :: var_b(:,:,:) !< The second field being compared
1551 integer, intent(in) :: turns !< Number of quarter turns from input to model domain
1552 character(len=256), intent(out) :: size_msg !< The array sizes
1553
1554 if (modulo(turns, 2) == 0) then
1555 size_mismatch_3d = ( (size(var_a,1) /= size(var_b,1)) .or. &
1556 (size(var_a,2) /= size(var_b,2)) .or. &
1557 (size(var_a,3) /= size(var_b,3)) )
1558 else
1559 size_mismatch_3d = ( (size(var_a,1) /= size(var_b,2)) .or. &
1560 (size(var_a,2) /= size(var_b,1)) .or. &
1561 (size(var_a,3) /= size(var_b,3)) )
1562 endif
1563 write(size_msg, '(3(1x,I0), " vs ", 3(1x,I0))') size(var_a,1), size(var_a,2), size(var_a,3), &
1564 size(var_b,1), size(var_b,2), size(var_b,3)
1565end function size_mismatch_3d
1566
1567
1568!> save_restart saves all registered variables to restart files.
1569subroutine save_restart(directory, time, G, CS, time_stamped, filename, GV, num_rest_files, write_IC)
1570 character(len=*), intent(in) :: directory !< The directory where the restart files
1571 !! are to be written
1572 type(time_type), intent(in) :: time !< The current model time
1573 type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure as seen from the driver.
1574 type(mom_restart_cs), intent(inout) :: cs !< MOM restart control struct
1575 logical, optional, intent(in) :: time_stamped !< If present and true, add time-stamp
1576 !! to the restart file names
1577 character(len=*), optional, intent(in) :: filename !< A filename that overrides the name in CS%restartfile
1578 type(verticalgrid_type), &
1579 optional, intent(in) :: gv !< The ocean's vertical grid structure
1580 integer, optional, intent(out) :: num_rest_files !< number of restart files written
1581 logical, optional, intent(in) :: write_ic !< If present and true, initial conditions
1582 !! are being written
1583
1584 ! Local variables
1585 type(vardesc) :: vars(cs%max_fields) ! Descriptions of the fields that
1586 ! are to be read from the restart file.
1587 type(mom_field) :: fields(cs%max_fields) ! Opaque types containing metadata describing
1588 ! each variable that will be written.
1589 character(len=512) :: restartpath ! The restart file path (dir/file).
1590 character(len=256) :: restartname ! The restart file name (no dir).
1591 character(len=8) :: suffix ! A suffix (like _2) that is appended
1592 ! to the name of files after the first.
1593 integer(kind=int64) :: var_sz, size_in_file ! The size in bytes of each variable
1594 ! and the variables already in a file.
1595 integer(kind=int64), parameter :: max_file_size = 4294967292_int64 ! The maximum size in bytes for the
1596 ! starting position of each variable in a file's record,
1597 ! based on the use of NetCDF 3.6 or later. For earlier
1598 ! versions of NetCDF, the value was 2147483647_int64.
1599 integer :: start_var, next_var ! The starting variables of the
1600 ! current and next files.
1601 type(mom_infra_file) :: io_handle ! The I/O handle of the open fileset
1602 integer :: m, nz, na
1603 integer :: num_files ! The number of restart files that will be used.
1604 integer :: seconds, days, year, month, hour, minute
1605 character(len=8) :: z_grid, t_grid ! Variable grid info.
1606 integer :: pos ! A coded integer indicating the horizontal staggering of a variable
1607 real :: conv ! Shorthand for the conversion factor [a A-1 ~> 1]
1608 real :: restart_time ! The model time at whic the restart file is being written [days]
1609 character(len=32) :: filename_appendix = '' ! Appendix to filename for ensemble runs
1610 integer :: length ! The length of a text string.
1611 character(len=256) :: mesg, var_name
1612 integer(kind=int64) :: check_val(cs%max_fields,1)
1613 logical :: verbose
1614 integer :: isl, iel, jsl, jel
1615 integer :: turns ! Number of quarter turns from input to model domain
1616 integer, parameter :: nmax_extradims = 5
1617 type(axis_info), dimension(:), allocatable :: extra_axes
1618
1619 turns = cs%turns
1620
1621 allocate (extra_axes(nmax_extradims))
1622
1623 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1624 "save_restart: Module must be initialized before it is used.")
1625
1626 if (cs%novars > cs%max_fields) call restart_error(cs)
1627 verbose = (is_root_pe() .and. (mom_get_verbosity() >= 7))
1628
1629 ! With parallel read & write, it is possible to disable the following...
1630 num_files = 0
1631 next_var = 0
1632 nz = 1 ; if (present(gv)) nz = gv%ke
1633
1634 restart_time = time_type_to_real(time) / 86400.0
1635
1636 restartname = trim(cs%restartfile)
1637 if (present(filename)) restartname = trim(filename)
1638 if (PRESENT(time_stamped)) then ; if (time_stamped) then
1639 call get_date(time, year, month, days, hour, minute, seconds)
1640 ! Compute the year-day, because I don't like months. - RWH
1641 do m=1,month-1
1642 days = days + days_in_month(set_date(year, m, 2, 0, 0, 0))
1643 enddo
1644 seconds = seconds + 60*minute + 3600*hour
1645 if (year <= 9999) then
1646 write(restartname,'("_Y",I4.4,"_D",I3.3,"_S",I5.5)') year, days, seconds
1647 elseif (year <= 99999) then
1648 write(restartname,'("_Y",I5.5,"_D",I3.3,"_S",I5.5)') year, days, seconds
1649 else
1650 write(restartname,'("_Y",I10.10,"_D",I3.3,"_S",I5.5)') year, days, seconds
1651 endif
1652 restartname = trim(cs%restartfile)//trim(restartname)
1653 endif ; endif
1654
1655 ! Determine if there is a filename_appendix (used for ensemble runs).
1656 call get_filename_appendix(filename_appendix)
1657 if (len_trim(filename_appendix) > 0) then
1658 length = len_trim(restartname)
1659 if (restartname(length-2:length) == '.nc') then
1660 restartname = restartname(1:length-3)//'.'//trim(filename_appendix)//'.nc'
1661 else
1662 restartname = restartname(1:length) //'.'//trim(filename_appendix)
1663 endif
1664 endif
1665
1666 next_var = 1
1667 do while (next_var <= cs%novars )
1668 start_var = next_var
1669 size_in_file = 8*(2*g%Domain%niglobal+2*g%Domain%njglobal+2*nz+1000)
1670
1671 do m=start_var,cs%novars
1672 call query_vardesc(cs%restart_field(m)%vars, position=pos, &
1673 z_grid=z_grid, t_grid=t_grid, caller="save_restart", &
1674 extra_axes=extra_axes)
1675
1676 var_sz = get_variable_byte_size(pos, z_grid, t_grid, g, nz)
1677 ! factor in size of extra axes, or multiply by 1
1678 do na=1,nmax_extradims
1679 var_sz = var_sz*extra_axes(na)%ax_size
1680 enddo
1681
1682 if ((m==start_var) .OR. (size_in_file < max_file_size-var_sz)) then
1683 size_in_file = size_in_file + var_sz
1684 else ; exit
1685 endif
1686
1687 enddo
1688 next_var = m
1689
1690 restartpath = trim(directory) // trim(restartname)
1691
1692 write(suffix,'("_",I0)') num_files
1693
1694 length = len_trim(restartpath)
1695 if (length < 3) then ! This case is very uncommon but this test avoids segmentation-faults.
1696 if (num_files > 0) restartpath = trim(restartpath) // suffix
1697 restartpath = trim(restartpath)//".nc"
1698 elseif (restartpath(length-2:length) == ".nc") then
1699 if (num_files > 0) restartpath = restartpath(1:length-3)//trim(suffix)//".nc"
1700 else
1701 if (num_files > 0) restartpath = trim(restartpath) // suffix
1702 restartpath = trim(restartpath)//".nc"
1703 endif
1704
1705 do m=start_var,next_var-1
1706 vars(m-start_var+1) = cs%restart_field(m)%vars
1707 enddo
1708 call query_vardesc(vars(1), t_grid=t_grid, position=pos, caller="save_restart")
1709 t_grid = adjustl(t_grid)
1710 if (t_grid(1:1) /= 'p') &
1711 call modify_vardesc(vars(1), t_grid='s', caller="save_restart")
1712
1713 !Prepare the checksum of the restart fields to be written to restart files
1714 do m=start_var,next_var-1
1715
1716 call query_vardesc(vars(m), position=pos, name=var_name, caller="save_restart")
1717 if (modulo(turns, 2) == 0) then
1718 call get_checksum_loop_ranges(g, cs, pos, isl, iel, jsl, jel)
1719 else ! Note that G is always the unrotated grid as it is seen by the driver level.
1720 call get_checksum_loop_ranges(g, cs, pos, jsl, jel, isl, iel)
1721 endif
1722 if (verbose) then
1723 if (pos == center) then
1724 write(mesg, '(" is in CENTER position, checksum range",4(1x,I0))') isl, iel, jsl, jel
1725 elseif (pos == corner) then
1726 write(mesg, '(" is in CORNER position, checksum range",4(1x,I0))') isl, iel, jsl, jel
1727 elseif (pos == north_face) then
1728 write(mesg, '(" is in NORTH_FACE position, checksum range",4(1x,I0))') isl, iel, jsl, jel
1729 elseif (pos == east_face) then
1730 write(mesg, '(" is in EAST_FACE position, checksum range",4(1x,I0))') isl, iel, jsl, jel
1731 else
1732 write(mesg, '(" is in another position, ",I0,", checksum range",4(1x,I0))') pos, isl, iel, jsl, jel
1733 endif
1734 call mom_mesg(trim(var_name)//mesg)
1735 endif
1736
1737 conv = cs%restart_field(m)%conv
1738 if (associated(cs%var_ptr3d(m)%p)) then
1739 check_val(m-start_var+1,1) = chksum(cs%var_ptr3d(m)%p(isl:iel,jsl:jel,:), turns=-turns, unscale=conv)
1740 elseif (associated(cs%var_ptr2d(m)%p)) then
1741 check_val(m-start_var+1,1) = chksum(cs%var_ptr2d(m)%p(isl:iel,jsl:jel), turns=-turns, unscale=conv)
1742 elseif (associated(cs%var_ptr4d(m)%p)) then
1743 check_val(m-start_var+1,1) = chksum(cs%var_ptr4d(m)%p(isl:iel,jsl:jel,:,:), turns=-turns, unscale=conv)
1744 elseif (associated(cs%var_ptr1d(m)%p)) then
1745 check_val(m-start_var+1,1) = chksum(cs%var_ptr1d(m)%p(:), unscale=conv)
1746 elseif (associated(cs%var_ptr0d(m)%p)) then
1747 check_val(m-start_var+1,1) = chksum(cs%var_ptr0d(m)%p, pelist=(/pe_here()/), unscale=conv)
1748 endif
1749 enddo
1750
1751 if (cs%parallel_restartfiles) then
1752 call create_mom_file(io_handle, trim(restartpath), vars, next_var-start_var, &
1753 fields, multiple, g=g, gv=gv, checksums=check_val, extra_axes=extra_axes)
1754 else
1755 call create_mom_file(io_handle, trim(restartpath), vars, next_var-start_var, &
1756 fields, single_file, g=g, gv=gv, checksums=check_val, extra_axes=extra_axes)
1757 endif
1758
1759 do m=start_var,next_var-1
1760 if (associated(cs%var_ptr3d(m)%p)) then
1761 call mom_write_field(io_handle, fields(m-start_var+1), g%Domain, cs%var_ptr3d(m)%p, &
1762 restart_time, unscale=cs%restart_field(m)%conv, turns=-turns, &
1763 zero_zeros=cs%unsigned_zeros)
1764 elseif (associated(cs%var_ptr2d(m)%p)) then
1765 call mom_write_field(io_handle, fields(m-start_var+1), g%Domain, cs%var_ptr2d(m)%p, &
1766 restart_time, unscale=cs%restart_field(m)%conv, turns=-turns, &
1767 zero_zeros=cs%unsigned_zeros)
1768 elseif (associated(cs%var_ptr4d(m)%p)) then
1769 call mom_write_field(io_handle, fields(m-start_var+1), g%Domain, cs%var_ptr4d(m)%p, &
1770 restart_time, unscale=cs%restart_field(m)%conv, turns=-turns, &
1771 zero_zeros=cs%unsigned_zeros)
1772 elseif (associated(cs%var_ptr1d(m)%p)) then
1773 call mom_write_field(io_handle, fields(m-start_var+1), cs%var_ptr1d(m)%p, &
1774 restart_time, unscale=cs%restart_field(m)%conv, &
1775 zero_zeros=cs%unsigned_zeros)
1776 elseif (associated(cs%var_ptr0d(m)%p)) then
1777 call mom_write_field(io_handle, fields(m-start_var+1), cs%var_ptr0d(m)%p, &
1778 restart_time, unscale=cs%restart_field(m)%conv, &
1779 zero_zeros=cs%unsigned_zeros)
1780 endif
1781 enddo
1782
1783 call io_handle%close()
1784
1785 num_files = num_files+1
1786
1787 enddo
1788
1789 if (present(num_rest_files)) num_rest_files = num_files
1790
1791end subroutine save_restart
1792
1793!> restore_state reads the model state from previously generated files. All
1794!! restart variables are read from the first file in the input filename list
1795!! in which they are found.
1796subroutine restore_state(filename, directory, day, G, CS)
1797 character(len=*), intent(in) :: filename !< The list of restart file names or a single
1798 !! character 'r' to read automatically named files
1799 character(len=*), intent(in) :: directory !< The directory in which to find restart files
1800 type(time_type), intent(out) :: day !< The time of the restarted run
1801 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
1802 type(mom_restart_cs), intent(inout) :: cs !< MOM restart control struct
1803
1804 ! Local variables
1805 real :: scale ! A scaling factor for reading a field [A a-1 ~> 1] to convert
1806 ! from the units in the file to the internal units of this field
1807 real :: conv ! The output conversion factor for writing a field [a A-1 ~> 1]
1808 character(len=512) :: mesg ! A message for warnings.
1809 character(len=80) :: varname ! A variable's name.
1810 integer :: num_file ! The number of files (restart files and others
1811 ! explicitly in filename) that are open.
1812 integer :: i, n, m, missing_fields
1813 integer :: isl, iel, jsl, jel
1814 integer :: nvar, ntime, pos
1815
1816 type(mom_infra_file) :: io_handles(cs%max_fields) ! The I/O units of all open files.
1817 character(len=200) :: unit_path(cs%max_fields) ! The file names.
1818 logical :: unit_is_global(cs%max_fields) ! True if the file is global.
1819
1820 real :: t1, t2 ! Two times from the start of different files [days].
1821 real, allocatable :: time_vals(:) ! Times from a file extracted with getl_file_times [days]
1822 type(mom_field), allocatable :: fields(:)
1823 logical :: is_there_a_checksum ! Is there a valid checksum that should be checked.
1824 integer(kind=int64) :: checksum_file ! The checksum value recorded in the input file.
1825 integer(kind=int64) :: checksum_data ! The checksum value for the data that was read in.
1826
1827 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
1828 "restore_state: Module must be initialized before it is used.")
1829
1830 if (cs%novars > cs%max_fields) call restart_error(cs)
1831
1832 ! Get NetCDF ids for all of the restart files.
1833 if ((len_trim(filename) == 1) .and. (filename(1:1) == 'F')) then
1834 num_file = open_restart_units('r', directory, g, cs, io_handles=io_handles, &
1835 file_paths=unit_path, global_files=unit_is_global)
1836 else
1837 num_file = open_restart_units(filename, directory, g, cs, io_handles=io_handles, &
1838 file_paths=unit_path, global_files=unit_is_global)
1839 endif
1840
1841 if (num_file == 0) then
1842 write(mesg,'("Unable to find any restart files specified by ",A," in directory ",A,".")') &
1843 trim(filename), trim(directory)
1844 call mom_error(fatal,"MOM_restart: "//mesg)
1845 endif
1846
1847 ! Get the time from the first file in the list that has one.
1848 do n=1,num_file
1849 call io_handles(n)%get_file_times(time_vals, ntime)
1850 if (ntime < 1) cycle
1851
1852 t1 = time_vals(1)
1853 deallocate(time_vals)
1854
1855 day = real_to_time(t1*86400.0)
1856 exit
1857 enddo
1858
1859 if (n>num_file) call mom_error(warning, "MOM_restart: No times found in restart files.")
1860
1861 ! Check the remaining files for different times and issue a warning
1862 ! if they differ from the first time.
1863 do m = n+1,num_file
1864 call io_handles(n)%get_file_times(time_vals, ntime)
1865 if (ntime < 1) cycle
1866
1867 t2 = time_vals(1)
1868 deallocate(time_vals)
1869
1870 if (t1 /= t2 .and. is_root_pe()) then
1871 write(mesg,'("WARNING: Restart file ",I0," has time ",F10.4,"whereas &
1872 &simulation is restarted at ",F10.4," (differing by ",F10.4,").")') &
1873 m, t1, t2, t1-t2
1874 call mom_error(warning, "MOM_restart: "//mesg)
1875 endif
1876 enddo
1877
1878 ! Read each variable from the first file in which it is found.
1879 do n=1,num_file
1880 call io_handles(n)%get_file_info(nvar=nvar)
1881
1882 allocate(fields(nvar))
1883 call io_handles(n)%get_file_fields(fields(1:nvar))
1884
1885 do m=1, nvar
1886 call io_handles(n)%get_field_atts(fields(m), name=varname)
1887 do i=1,cs%num_obsolete_vars
1888 if (adjustl(lowercase(trim(varname))) == adjustl(lowercase(trim(cs%restart_obsolete(i)%field_name)))) then
1889 call mom_error(fatal, "MOM_restart restore_state: Attempting to use obsolete restart field "//&
1890 trim(varname)//" - the new corresponding restart field is "//&
1891 trim(cs%restart_obsolete(i)%replacement_name))
1892 endif
1893 enddo
1894 enddo
1895
1896 missing_fields = 0
1897
1898 do m=1,cs%novars
1899 if (cs%restart_field(m)%initialized) cycle
1900 call query_vardesc(cs%restart_field(m)%vars, position=pos, caller="restore_state")
1901 conv = cs%restart_field(m)%conv
1902 if (conv == 0.0) then ; scale = 1.0 ; else ; scale = 1.0 / conv ; endif
1903
1904 if (modulo(cs%turns, 2) == 0) then
1905 call get_checksum_loop_ranges(g, cs, pos, isl, iel, jsl, jel)
1906 else ! Note that G is always the unrotated grid as it is used during initialization.
1907 call get_checksum_loop_ranges(g, cs, pos, jsl, jel, isl, iel)
1908 endif
1909 do i=1, nvar
1910 call io_handles(n)%get_field_atts(fields(i), name=varname)
1911 if (lowercase(trim(varname)) == lowercase(trim(cs%restart_field(m)%var_name))) then
1912 checksum_data = -1
1913 if (cs%checksum_required) then
1914 call io_handles(n)%read_field_chksum(fields(i), checksum_file, is_there_a_checksum)
1915 else
1916 checksum_file = -1
1917 is_there_a_checksum = .false. ! Do not need to do data checksumming.
1918 endif
1919
1920 if (associated(cs%var_ptr1d(m)%p)) then
1921 ! Read a 1d array, which should be invariant to domain decomposition.
1922 call mom_read_data(unit_path(n), varname, cs%var_ptr1d(m)%p, &
1923 timelevel=1, scale=scale, mom_domain=g%Domain)
1924 if (is_there_a_checksum) checksum_data = chksum(cs%var_ptr1d(m)%p(:), unscale=conv)
1925 elseif (associated(cs%var_ptr0d(m)%p)) then ! Read a scalar...
1926 call mom_read_data(unit_path(n), varname, cs%var_ptr0d(m)%p, &
1927 timelevel=1, scale=scale, mom_domain=g%Domain)
1928 if (is_there_a_checksum) checksum_data = chksum(cs%var_ptr0d(m)%p, pelist=(/pe_here()/), unscale=conv)
1929 elseif (associated(cs%var_ptr2d(m)%p)) then ! Read a 2d array.
1930 if (pos /= 0) then
1931 call mom_read_data(unit_path(n), varname, cs%var_ptr2d(m)%p, &
1932 g%Domain, timelevel=1, position=pos, scale=scale, turns=cs%turns)
1933 else ! This array is not domain-decomposed. This variant may be under-tested.
1934 call mom_error(fatal, &
1935 "MOM_restart does not support 2-d arrays without domain decomposition.")
1936 ! call read_data(unit_path(n), varname, CS%var_ptr2d(m)%p,no_domain=.true., timelevel=1)
1937 endif
1938 if (is_there_a_checksum) checksum_data = chksum(cs%var_ptr2d(m)%p(isl:iel,jsl:jel), unscale=conv)
1939 elseif (associated(cs%var_ptr3d(m)%p)) then ! Read a 3d array.
1940 if (pos /= 0) then
1941 call mom_read_data(unit_path(n), varname, cs%var_ptr3d(m)%p, &
1942 g%Domain, timelevel=1, position=pos, scale=scale, turns=cs%turns)
1943 else ! This array is not domain-decomposed. This variant may be under-tested.
1944 call mom_error(fatal, &
1945 "MOM_restart does not support 3-d arrays without domain decomposition.")
1946 ! call read_data(unit_path(n), varname, CS%var_ptr3d(m)%p, no_domain=.true., timelevel=1)
1947 endif
1948 if (is_there_a_checksum) checksum_data = chksum(cs%var_ptr3d(m)%p(isl:iel,jsl:jel,:), unscale=conv)
1949 elseif (associated(cs%var_ptr4d(m)%p)) then ! Read a 4d array.
1950 if (pos /= 0) then
1951 call mom_read_data(unit_path(n), varname, cs%var_ptr4d(m)%p, &
1952 g%Domain, timelevel=1, position=pos, scale=scale, &
1953 global_file=unit_is_global(n), turns=cs%turns)
1954 else ! This array is not domain-decomposed. This variant may be under-tested.
1955 call mom_error(fatal, &
1956 "MOM_restart does not support 4-d arrays without domain decomposition.")
1957 ! call read_data(unit_path(n), varname, CS%var_ptr4d(m)%p, no_domain=.true., timelevel=1)
1958 endif
1959 if (is_there_a_checksum) checksum_data = chksum(cs%var_ptr4d(m)%p(isl:iel,jsl:jel,:,:), unscale=conv)
1960 else
1961 call mom_error(fatal, "MOM_restart restore_state: No pointers set for "//trim(varname))
1962 endif
1963
1964 if (is_root_pe() .and. is_there_a_checksum .and. (checksum_file /= checksum_data)) then
1965 write (mesg,'(a,Z16,a,Z16,a)') "Checksum of input field "// trim(varname)//" ",checksum_data,&
1966 " does not match value ", checksum_file, &
1967 " stored in "//trim(unit_path(n)//"." )
1968 call mom_error(fatal, "MOM_restart(restore_state): "//trim(mesg) )
1969 endif
1970
1971 cs%restart_field(m)%initialized = .true.
1972 exit ! Start search for next restart variable.
1973 endif
1974 enddo
1975 if (i>nvar) missing_fields = missing_fields+1
1976 enddo
1977
1978 deallocate(fields)
1979 if (missing_fields == 0) exit
1980 enddo
1981
1982 do n=1,num_file
1983 call io_handles(n)%close()
1984 enddo
1985
1986 ! Check whether any mandatory fields have not been found.
1987 cs%restart = .true.
1988 do m=1,cs%novars
1989 if (.not.(cs%restart_field(m)%initialized)) then
1990 cs%restart = .false.
1991 if (cs%restart_field(m)%mand_var) then
1992 call mom_error(fatal,"MOM_restart: Unable to find mandatory variable " &
1993 //trim(cs%restart_field(m)%var_name)//" in restart files.")
1994 endif
1995 endif
1996 enddo
1997
1998 ! Lock the restart registry so that no further variables can be registered.
1999 cs%locked = .true.
2000
2001end subroutine restore_state
2002
2003
2004
2005!> restart_files_exist determines whether any restart files exist.
2006function restart_files_exist(filename, directory, G, CS)
2007 character(len=*), intent(in) :: filename !< The list of restart file names or a single
2008 !! character 'r' to read automatically named files
2009 character(len=*), intent(in) :: directory !< The directory in which to find restart files
2010 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
2011 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
2012 logical :: restart_files_exist !< The function result, which indicates whether
2013 !! any of the explicitly or automatically named
2014 !! restart files exist in directory
2015 integer :: num_files
2016
2017 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
2018 "restart_files_exist: Module must be initialized before it is used.")
2019
2020 if ((len_trim(filename) == 1) .and. (filename(1:1) == 'F')) then
2021 num_files = get_num_restart_files('r', directory, g, cs)
2022 else
2023 num_files = get_num_restart_files(filename, directory, g, cs)
2024 endif
2025 restart_files_exist = (num_files > 0)
2026end function restart_files_exist
2027
2028!> determine_is_new_run determines from the value of filename and the existence
2029!! automatically named restart files in directory whether this would be a new,
2030!! and as a side effect stores this information in CS.
2031function determine_is_new_run(filename, directory, G, CS) result(is_new_run)
2032 character(len=*), intent(in) :: filename !< The list of restart file names or a single
2033 !! character 'r' to read automatically named files
2034 character(len=*), intent(in) :: directory !< The directory in which to find restart files
2035 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
2036 type(mom_restart_cs), intent(inout) :: cs !< MOM restart control struct
2037 logical :: is_new_run !< The function result, which indicates whether
2038 !! this is a new run, based on the value of
2039 !! filename and whether restart files exist
2040
2041 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
2042 "determine_is_new_run: Module must be initialized before it is used.")
2043
2044 if (len_trim(filename) > 1) then
2045 cs%new_run = .false.
2046 elseif (len_trim(filename) == 0) then
2047 cs%new_run = .true.
2048 elseif (filename(1:1) == 'n') then
2049 cs%new_run = .true.
2050 elseif (filename(1:1) == 'F') then
2051 cs%new_run = (get_num_restart_files('r', directory, g, cs) == 0)
2052 else
2053 cs%new_run = .false.
2054 endif
2055
2056 cs%new_run_set = .true.
2057 is_new_run = cs%new_run
2058end function determine_is_new_run
2059
2060!> is_new_run returns whether this is going to be a new run based on the
2061!! information stored in CS by a previous call to determine_is_new_run.
2062function is_new_run(CS)
2063 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
2064
2065 logical :: is_new_run !< The function result, which had been stored in CS during
2066 !! a previous call to determine_is_new_run
2067
2068 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
2069 "is_new_run: Module must be initialized before it is used.")
2070
2071 if (.not.cs%new_run_set) call mom_error(fatal, "MOM_restart " // &
2072 "determine_is_new_run must be called for a restart file before is_new_run.")
2073
2074 is_new_run = cs%new_run
2075end function is_new_run
2076
2077!> open_restart_units determines the number of existing restart files and optionally opens
2078!! them and returns unit ids, paths and whether the files are global or spatially decomposed.
2079function open_restart_units(filename, directory, G, CS, IO_handles, file_paths, &
2080 global_files) result(num_files)
2081 character(len=*), intent(in) :: filename !< The list of restart file names or a single
2082 !! character 'r' to read automatically named files
2083 character(len=*), intent(in) :: directory !< The directory in which to find restart files
2084 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
2085 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
2086
2087 type(mom_infra_file), dimension(:), &
2088 optional, intent(out) :: io_handles !< The I/O handles of all opened files
2089 character(len=*), dimension(:), &
2090 optional, intent(out) :: file_paths !< The full paths to open files
2091 logical, dimension(:), &
2092 optional, intent(out) :: global_files !< True if a file is global
2093
2094 integer :: num_files !< The number of files (both automatically named restart
2095 !! files and others explicitly in filename) that have been opened.
2096
2097 ! Local variables
2098 character(len=256) :: filepath ! The path (dir/file) to the file being opened.
2099 character(len=256) :: fname ! The name of the current file.
2100 character(len=8) :: suffix ! A suffix (like "_2") that is added to any
2101 ! additional restart files.
2102 integer :: num_restart ! The number of restart files that have already
2103 ! been opened using their numbered suffix.
2104 integer :: start_char ! The location of the starting character in the
2105 ! current file name.
2106 integer :: nf ! The number of files that have been found so far
2107 integer :: m, length
2108 logical :: still_looking ! If true, the code is still looking for automatically named files
2109 logical :: fexists ! True if a file has been found
2110 character(len=32) :: filename_appendix = '' ! Filename appendix for ensemble runs
2111 character(len=80) :: restartname
2112
2113 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
2114 "open_restart_units: Module must be initialized before it is used.")
2115
2116 ! Get NetCDF ids for all of the restart files.
2117 num_restart = 0 ; nf = 0 ; start_char = 1
2118 do while (start_char <= len_trim(filename) )
2119 do m=start_char,len_trim(filename)
2120 if (filename(m:m) == ' ') exit
2121 enddo
2122 fname = filename(start_char:m-1)
2123 start_char = m
2124 do while (start_char <= len_trim(filename))
2125 if (filename(start_char:start_char) == ' ') then
2126 start_char = start_char + 1
2127 else
2128 exit
2129 endif
2130 enddo
2131
2132 if ((fname(1:1)=='r') .and. ( len_trim(fname) == 1)) then
2133 still_looking = (num_restart <= 0) ! Avoid going through the file list twice.
2134 do while (still_looking)
2135 restartname = trim(cs%restartfile)
2136
2137 ! Determine if there is a filename_appendix (used for ensemble runs).
2138 call get_filename_appendix(filename_appendix)
2139 if (len_trim(filename_appendix) > 0) then
2140 length = len_trim(restartname)
2141 if (restartname(length-2:length) == '.nc') then
2142 restartname = restartname(1:length-3)//'.'//trim(filename_appendix)//'.nc'
2143 else
2144 restartname = restartname(1:length) //'.'//trim(filename_appendix)
2145 endif
2146 endif
2147 filepath = trim(directory) // trim(restartname)
2148
2149 write(suffix,'("_",I0)') num_restart
2150 if (num_restart > 0) filepath = trim(filepath) // suffix
2151
2152 filepath = trim(filepath)//".nc"
2153
2154 num_restart = num_restart + 1
2155 ! Look for a global netCDF file.
2156 inquire(file=filepath, exist=fexists)
2157 if (fexists) then
2158 nf = nf + 1
2159 if (present(io_handles)) &
2160 call io_handles(nf)%open(trim(filepath), readonly_file, &
2161 mom_domain=g%Domain, threading=multiple, fileset=single_file)
2162 if (present(global_files)) global_files(nf) = .true.
2163 if (present(file_paths)) file_paths(nf) = filepath
2164 elseif (cs%parallel_restartfiles) then
2165 ! Look for decomposed files using the I/O Layout.
2166 fexists = file_exists(filepath, g%Domain)
2167 if (fexists) then
2168 nf = nf + 1
2169 if (present(io_handles)) &
2170 call io_handles(nf)%open(trim(filepath), readonly_file, mom_domain=g%Domain)
2171 if (present(global_files)) global_files(nf) = .false.
2172 if (present(file_paths)) file_paths(nf) = filepath
2173 endif
2174 endif
2175
2176 if (fexists) then
2177 if (is_root_pe() .and. (present(io_handles))) &
2178 call mom_error(note, "MOM_restart: MOM run restarted using : "//trim(filepath))
2179 else
2180 still_looking = .false. ; exit
2181 endif
2182 enddo ! while (still_looking) loop
2183 else
2184 filepath = trim(directory)//trim(fname)
2185 inquire(file=filepath, exist=fexists)
2186 if (.not. fexists) filepath = trim(filepath)//".nc"
2187
2188 inquire(file=filepath, exist=fexists)
2189 if (fexists) then
2190 nf = nf + 1
2191 if (present(io_handles)) &
2192 call io_handles(nf)%open(trim(filepath), readonly_file, &
2193 mom_domain=g%Domain, threading=multiple, fileset=single_file)
2194 if (present(global_files)) global_files(nf) = .true.
2195 if (present(file_paths)) file_paths(nf) = filepath
2196 if (is_root_pe() .and. (present(io_handles))) &
2197 call mom_error(note, "MOM_restart: MOM run restarted using : "//trim(filepath))
2198 else
2199 if (present(io_handles)) &
2200 call mom_error(warning, "MOM_restart: Unable to find restart file : "//trim(filepath))
2201 endif
2202
2203 endif
2204 enddo ! while (start_char < len_trim(filename)) loop
2205 num_files = nf
2206
2207end function open_restart_units
2208
2209!> get_num_restart_files returns the number of existing restart files that match the provided
2210!! directory structure and other information stored in the control structure and optionally
2211!! also provides the full paths to these files.
2212function get_num_restart_files(filenames, directory, G, CS, file_paths) result(num_files)
2213 character(len=*), intent(in) :: filenames !< The list of restart file names or a single
2214 !! character 'r' to read automatically named files
2215 character(len=*), intent(in) :: directory !< The directory in which to find restart files
2216 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
2217 type(mom_restart_cs), intent(in) :: cs !< MOM restart control struct
2218 character(len=*), dimension(:), &
2219 optional, intent(out) :: file_paths !< The full paths to the restart files.
2220
2221 integer :: num_files !< The function result, the number of files (both automatically named
2222 !! restart files and others explicitly in filename) that have been opened
2223
2224 if (.not.cs%initialized) call mom_error(fatal, "MOM_restart " // &
2225 "get_num_restart_files: Module must be initialized before it is used.")
2226
2227 ! This call uses open_restart_units without the optional arguments needed to actually
2228 ! open the files to determine the number of restart files.
2229 num_files = open_restart_units(filenames, directory, g, cs, file_paths=file_paths)
2230
2231end function get_num_restart_files
2232
2233
2234!> Initialize this module and set up a restart control structure.
2235subroutine restart_init(param_file, CS, restart_root)
2236 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
2237 type(mom_restart_cs), pointer :: cs !< A pointer to a MOM_restart_CS object that is allocated here
2238 character(len=*), optional, &
2239 intent(in) :: restart_root !< A filename root that overrides the value
2240 !! set by RESTARTFILE to enable the use of this module by
2241 !! other components than MOM.
2242
2243 logical :: rotate_index
2244
2245 ! This include declares and sets the variable "version".
2246# include "version_variable.h"
2247 character(len=40) :: mdl = "MOM_restart" ! This module's name.
2248 logical :: all_default ! If true, all parameters are using their default values.
2249
2250 if (associated(cs)) then
2251 call mom_error(warning, "restart_init called with an associated control structure.")
2252 return
2253 endif
2254 allocate(cs)
2255
2256 cs%initialized = .true.
2257
2258 ! Determine whether all paramters are set to their default values.
2259 call get_param(param_file, mdl, "PARALLEL_RESTARTFILES", cs%parallel_restartfiles, &
2260 default=.false., do_not_log=.true.)
2261 call get_param(param_file, mdl, "MAX_FIELDS", cs%max_fields, default=100, do_not_log=.true.)
2262 call get_param(param_file, mdl, "RESTART_CHECKSUMS_REQUIRED", cs%checksum_required, &
2263 default=.true., do_not_log=.true.)
2264 call get_param(param_file, mdl, "RESTART_SYMMETRIC_CHECKSUMS", cs%symmetric_checksums, &
2265 default=.false., do_not_log=.true.)
2266 call get_param(param_file, mdl, "RESTART_UNSIGNED_ZEROS", cs%unsigned_zeros, &
2267 default=.false., do_not_log=.true.)
2268 all_default = ((.not.cs%parallel_restartfiles) .and. (cs%max_fields == 100) .and. &
2269 (cs%checksum_required) .and. (.not.cs%symmetric_checksums) .and. (.not.cs%unsigned_zeros))
2270 if (.not.present(restart_root)) then
2271 call get_param(param_file, mdl, "RESTARTFILE", cs%restartfile, &
2272 default="MOM.res", do_not_log=.true.)
2273 all_default = (all_default .and. (trim(cs%restartfile) == trim("MOM.res")))
2274 endif
2275
2276 ! Read all relevant parameters and write them to the model log.
2277 call log_version(param_file, mdl, version, "", all_default=all_default)
2278 call get_param(param_file, mdl, "PARALLEL_RESTARTFILES", cs%parallel_restartfiles, &
2279 "If true, the IO layout is used to group processors that write to the same "//&
2280 "restart file or each processor writes its own (numbered) restart file. "//&
2281 "If false, a single restart file is generated combining output from all PEs.", &
2282 default=.false.)
2283
2284 if (present(restart_root)) then
2285 cs%restartfile = restart_root
2286 call log_param(param_file, mdl, "RESTARTFILE from argument", cs%restartfile)
2287 else
2288 call get_param(param_file, mdl, "RESTARTFILE", cs%restartfile, &
2289 "The name-root of the restart file.", default="MOM.res")
2290 endif
2291 call get_param(param_file, mdl, "MAX_FIELDS", cs%max_fields, &
2292 "The maximum number of restart fields that can be used.", &
2293 default=100)
2294 call get_param(param_file, mdl, "RESTART_CHECKSUMS_REQUIRED", cs%checksum_required, &
2295 "If true, require the restart checksums to match and error out otherwise. "//&
2296 "Users may want to avoid this comparison if for example the restarts are "//&
2297 "made from a run with a different mask_table than the current run, "//&
2298 "in which case the checksums will not match and cause crash.",&
2299 default=.true.)
2300 call get_param(param_file, mdl, "RESTART_SYMMETRIC_CHECKSUMS", cs%symmetric_checksums, &
2301 "If true, do the restart checksums on all the edge points for a non-reentrant "//&
2302 "grid. This requires that SYMMETRIC_MEMORY_ is defined at compile time.", &
2303 default=.false.)
2304 call get_param(param_file, mdl, "RESTART_UNSIGNED_ZEROS", cs%unsigned_zeros, &
2305 "If true, convert any negative zeros that would be written to the restart file "//&
2306 "into ordinary unsigned zeros. This does not change answers, but it can be "//&
2307 "helpful in comparing restart files after grid rotation, for example.", &
2308 default=.false.)
2309 call get_param(param_file, mdl, "REENTRANT_X", cs%reentrant_x, &
2310 "If true, the domain is zonally reentrant.", default=.true., do_not_log=.true.)
2311 call get_param(param_file, mdl, "REENTRANT_Y", cs%reentrant_y, &
2312 "If true, the domain is meridionally reentrant.", default=.false., do_not_log=.true.)
2313
2314 ! Maybe not the best place to do this?
2315 call get_param(param_file, mdl, "ROTATE_INDEX", rotate_index, &
2316 default=.false., do_not_log=.true.)
2317
2318 cs%turns = 0
2319 if (rotate_index) then
2320 call get_param(param_file, mdl, "INDEX_TURNS", cs%turns, &
2321 default=1, do_not_log=.true.)
2322 endif
2323
2324 allocate(cs%restart_field(cs%max_fields))
2325 allocate(cs%restart_obsolete(cs%max_fields))
2326 allocate(cs%var_ptr0d(cs%max_fields))
2327 allocate(cs%var_ptr1d(cs%max_fields))
2328 allocate(cs%var_ptr2d(cs%max_fields))
2329 allocate(cs%var_ptr3d(cs%max_fields))
2330 allocate(cs%var_ptr4d(cs%max_fields))
2331
2332 cs%locked = .false.
2333
2334end subroutine restart_init
2335
2336!> Issue an error message if the restart_registry is locked.
2337subroutine lock_check(CS, var_desc, name)
2338 type(mom_restart_cs), intent(in) :: cs !< A MOM_restart_CS object (intent in)
2339 type(vardesc), optional, intent(in) :: var_desc !< A structure with metadata about this variable
2340 character(len=*), optional, intent(in) :: name !< variable name to be used in the restart file
2341
2342 character(len=256) :: var_name ! A variable name.
2343
2344 if (cs%locked) then
2345 if (present(var_desc)) then
2346 call query_vardesc(var_desc, name=var_name)
2347 call mom_error(fatal, "Attempted to register "//trim(var_name)//" but the restart registry is locked.")
2348 elseif (present(name)) then
2349 call mom_error(fatal, "Attempted to register "//trim(name)//" but the restart registry is locked.")
2350 else
2351 call mom_error(fatal, "Attempted to register a variable but the restart registry is locked.")
2352 endif
2353 endif
2354
2355end subroutine lock_check
2356
2357!> Lock the restart registry so that an error is issued if any further restart variables are registered.
2358subroutine restart_registry_lock(CS, unlocked)
2359 type(mom_restart_cs), intent(inout) :: cs !< A MOM_restart_CS object (intent inout)
2360 logical, optional, intent(in) :: unlocked !< If present and true, unlock the registry
2361
2362 cs%locked = .true.
2363 if (present(unlocked)) cs%locked = .not.unlocked
2364end subroutine restart_registry_lock
2365
2366!> Indicate that all variables have now been registered and lock the registry.
2367subroutine restart_init_end(CS)
2368 type(mom_restart_cs), pointer :: cs !< A pointer to a MOM_restart_CS object
2369
2370 if (associated(cs)) then
2371 cs%locked = .true.
2372
2373 if (cs%novars == 0) call restart_end(cs)
2374 endif
2375
2376end subroutine restart_init_end
2377
2378!> Deallocate memory associated with a MOM_restart_CS variable.
2379subroutine restart_end(CS)
2380 type(mom_restart_cs), pointer :: cs !< A pointer to a MOM_restart_CS object
2381
2382 if (associated(cs%restart_field)) deallocate(cs%restart_field)
2383 if (associated(cs%restart_obsolete)) deallocate(cs%restart_obsolete)
2384 if (associated(cs%var_ptr0d)) deallocate(cs%var_ptr0d)
2385 if (associated(cs%var_ptr1d)) deallocate(cs%var_ptr1d)
2386 if (associated(cs%var_ptr2d)) deallocate(cs%var_ptr2d)
2387 if (associated(cs%var_ptr3d)) deallocate(cs%var_ptr3d)
2388 if (associated(cs%var_ptr4d)) deallocate(cs%var_ptr4d)
2389 deallocate(cs)
2390
2391end subroutine restart_end
2392
2393subroutine restart_error(CS)
2394 type(mom_restart_cs), intent(in) :: CS !< MOM restart control struct
2395
2396 character(len=16) :: num ! String for error messages
2397
2398 if (cs%novars > cs%max_fields) then
2399 write(num,'(I0)') cs%novars
2400 call mom_error(fatal,"MOM_restart: Too many fields registered for " // &
2401 "restart. Set MAX_FIELDS to be at least "//trim(num)//" in the MOM input file.")
2402 else
2403 call mom_error(fatal,"MOM_restart: Unspecified fatal error.")
2404 endif
2405end subroutine restart_error
2406
2407!> Return bounds for computing checksums to store in restart files
2408subroutine get_checksum_loop_ranges(G, CS, pos, isL, ieL, jsL, jeL)
2409 type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
2410 type(mom_restart_cs), intent(in) :: CS !< MOM restart control structure
2411 integer, intent(in) :: pos !< A coded integer indicating the horizontal staggering
2412 !! of a variable
2413 integer, intent(out) :: isL !< i-start for checksum
2414 integer, intent(out) :: ieL !< i-end for checksum
2415 integer, intent(out) :: jsL !< j-start for checksum
2416 integer, intent(out) :: jeL !< j-end for checksum
2417
2418 ! Regular non-symmetric compute domain
2419 isl = g%isc-g%isd+1
2420 iel = g%iec-g%isd+1
2421 jsl = g%jsc-g%jsd+1
2422 jel = g%jec-g%jsd+1
2423
2424 ! Expand range east or south for symmetric arrays
2425 if (cs%symmetric_checksums) then
2426 if (.not.g%symmetric) call mom_error(fatal, &
2427 "Setting SYMMETRIC_RESTART_CHECKSUMS to true only works with symmetric memory allocation, "//&
2428 "which is specified at compile time by defining the cpp macro SYMMETRIC_MEMORY_.")
2429
2430 if (((pos == east_face) .or. (pos == corner)) .and. (.not.cs%reentrant_x)) then ! For u-, q-points only
2431 if (g%isc+g%idg_offset == 1) isl = isl - 1 ! Include western edge in checksums only for western PEs
2432 endif
2433 if (((pos == north_face) .or. (pos == corner)) .and. (.not.cs%reentrant_y)) then ! For v-, q-points only
2434 if (g%jsc+g%jdg_offset == 1) jsl = jsl - 1 ! Include southern edge in checksums only for southern PEs
2435 endif
2436 endif
2437
2438end subroutine get_checksum_loop_ranges
2439
2440!> get the size of a variable in bytes
2441function get_variable_byte_size(pos, z_grid, t_grid, G, num_z) result(var_sz)
2442 integer, intent(in) :: pos !< An integer indicating the horizontal staggering position
2443 character(len=8), intent(in) :: z_grid !< The vertical grid string to interpret
2444 character(len=8), intent(in) :: t_grid !< A time string to interpret
2445 type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
2446 integer, intent(in) :: num_z !< The number of vertical layers in the grid
2447 integer(kind=int64) :: var_sz !< The function result, the size in bytes of a variable
2448
2449 ! Local variables
2450 integer :: var_periods ! The number of entries in a time-periodic axis
2451 character(len=8) :: t_grid_read, t_grid_tmp ! Modified versions of t_grid
2452
2453 if (pos == 0) then
2454 var_sz = 8
2455 else ! This may be an overestimate, as it is based on symmetric-memory corner points.
2456 var_sz = 8*(g%Domain%niglobal+1)*(g%Domain%njglobal+1)
2457 endif
2458
2459 select case (trim(z_grid))
2460 case ('L') ; var_sz = var_sz * num_z
2461 case ('i') ; var_sz = var_sz * (num_z+1)
2462 end select
2463
2464 t_grid_tmp = adjustl(t_grid)
2465 if (t_grid_tmp(1:1) == 'p') then
2466 if (len_trim(t_grid_tmp(2:8)) > 0) then
2467 var_periods = -1
2468 t_grid_read = adjustl(t_grid_tmp(2:8))
2469 read(t_grid_read,*) var_periods
2470 if (var_periods > 1) var_sz = var_sz * var_periods
2471 endif
2472 endif
2473
2474end function get_variable_byte_size
2475
2476end module mom_restart