MOM6
MOM_file_parser.F90
1 !> The MOM6 facility to parse input files for runtime parameters
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_coms, only : root_pe, broadcast
7 use mom_error_handler, only : mom_error, fatal, warning, mom_mesg
8 use mom_error_handler, only : is_root_pe, stdlog, stdout
9 use mom_time_manager, only : get_time, time_type, get_ticks_per_second
10 use mom_time_manager, only : set_date, get_date, real_to_time, operator(-), set_time
11 use mom_document, only : doc_param, doc_module, doc_init, doc_end, doc_type
12 use mom_document, only : doc_openblock, doc_closeblock
13 use mom_string_functions, only : left_int, left_ints, slasher
14 use mom_string_functions, only : left_real, left_reals
15 
16 implicit none ; private
17 
18 integer, parameter, public :: max_param_files = 5 !< Maximum number of parameter files.
19 integer, parameter :: input_str_length = 320 !< Maximum line length in parameter file.
20 integer, parameter :: filename_length = 200 !< Maximum number of characters in file names.
21 
22 ! The all_PEs_read option should be eliminated with post-riga shared code.
23 logical :: all_pes_read = .false. !< If true, all PEs read the input files
24  !! TODO: Eliminate this parameter
25 
26 !>@{ Default values for parameters
27 logical, parameter :: report_unused_default = .true.
28 logical, parameter :: unused_params_fatal_default = .false.
29 logical, parameter :: log_to_stdout_default = .false.
30 logical, parameter :: complete_doc_default = .true.
31 logical, parameter :: minimal_doc_default = .true.
32 !>@}
33 
34 !> The valid lines extracted from an input parameter file without comments
35 type, private :: file_data_type ; private
36  integer :: num_lines = 0 !< The number of lines in this type
37  character(len=INPUT_STR_LENGTH), pointer, dimension(:) :: line => null() !< The line content
38  logical, pointer, dimension(:) :: line_used => null() !< If true, the line has been read
39 end type file_data_type
40 
41 !> A link in the list of variables that have already had override warnings issued
42 type :: link_parameter ; private
43  type(link_parameter), pointer :: next => null() !< Facilitates linked list
44  character(len=80) :: name !< Parameter name
45  logical :: hasissuedoverridewarning = .false. !< Has a default value
46 end type link_parameter
47 
48 !> Specify the active parameter block
49 type :: parameter_block ; private
50  character(len=240) :: name = '' !< The active parameter block name
51 end type parameter_block
52 
53 !> A structure that can be parsed to read and document run-time parameters.
54 type, public :: param_file_type ; private
55  integer :: nfiles = 0 !< The number of open files.
56  integer :: iounit(max_param_files) !< The unit numbers of open files.
57  character(len=FILENAME_LENGTH) :: filename(max_param_files) !< The names of the open files.
58  logical :: netcdf_file(max_param_files) !< If true, the input file is in NetCDF.
59  ! This is not yet implemented.
60  type(file_data_type) :: param_data(max_param_files) !< Structures that contain
61  !! the valid data lines from the parameter
62  !! files, enabling all subsequent reads of
63  !! parameter data to occur internally.
64  logical :: report_unused = report_unused_default !< If true, report any
65  !! parameter lines that are not used in the run.
66  logical :: unused_params_fatal = unused_params_fatal_default !< If true, kill
67  !! the run if there are any unused parameters.
68  logical :: log_to_stdout = log_to_stdout_default !< If true, all log
69  !! messages are also sent to stdout.
70  logical :: log_open = .false. !< True if the log file has been opened.
71  integer :: stdout !< The unit number from stdout().
72  integer :: stdlog !< The unit number from stdlog().
73  character(len=240) :: doc_file !< A file where all run-time parameters, their
74  !! settings and defaults are documented.
75  logical :: complete_doc = complete_doc_default !< If true, document all
76  !! run-time parameters.
77  logical :: minimal_doc = minimal_doc_default !< If true, document only those
78  !! run-time parameters that differ from defaults.
79  type(doc_type), pointer :: doc => null() !< A structure that contains information
80  !! related to parameter documentation.
81  type(link_parameter), pointer :: chain => null() !< Facilitates linked list
82  type(parameter_block), pointer :: blockname => null() !< Name of active parameter block
83 end type param_file_type
84 
85 public read_param, open_param_file, close_param_file, log_param, log_version
86 public doc_param, get_param
87 public clearparameterblock, openparameterblock, closeparameterblock
88 
89 !> An overloaded interface to read various types of parameters
90 interface read_param
91  module procedure read_param_int, read_param_real, read_param_logical, &
92  read_param_char, read_param_char_array, read_param_time, &
93  read_param_int_array, read_param_real_array
94 end interface
95 !> An overloaded interface to log the values of various types of parameters
96 interface log_param
97  module procedure log_param_int, log_param_real, log_param_logical, &
98  log_param_char, log_param_time, &
99  log_param_int_array, log_param_real_array
100 end interface
101 !> An overloaded interface to read and log the values of various types of parameters
102 interface get_param
103  module procedure get_param_int, get_param_real, get_param_logical, &
104  get_param_char, get_param_char_array, get_param_time, &
105  get_param_int_array, get_param_real_array
106 end interface
107 
108 !> An overloaded interface to log version information about modules
109 interface log_version
110  module procedure log_version_cs, log_version_plain
111 end interface
112 
113 contains
114 
115 !> Make the contents of a parameter input file availalble in a param_file_type
116 subroutine open_param_file(filename, CS, checkable, component, doc_file_dir)
117  character(len=*), intent(in) :: filename !< An input file name, optionally with the full path
118  type(param_file_type), intent(inout) :: CS !< The control structure for the file_parser module,
119  !! it is also a structure to parse for run-time parameters
120  logical, optional, intent(in) :: checkable !< If this is false, it disables checks of this
121  !! file for unused parameters. The default is True.
122  character(len=*), optional, intent(in) :: component !< If present, this component name is used
123  !! to generate parameter documentation file names; the default is"MOM"
124  character(len=*), optional, intent(in) :: doc_file_dir !< An optional directory in which to write out
125  !! the documentation files. The default is effectively './'.
126 
127  ! Local variables
128  logical :: file_exists, unit_in_use, Netcdf_file, may_check
129  integer :: ios, iounit, strlen, i
130  character(len=240) :: doc_path
131  type(parameter_block), pointer :: block => null()
132 
133  may_check = .true. ; if (present(checkable)) may_check = checkable
134 
135  ! Check for non-blank filename
136  strlen = len_trim(filename)
137  if (strlen == 0) then
138  call mom_error(fatal, "open_param_file: Input file has not been specified.")
139  endif
140 
141  ! Check that this file has not already been opened
142  if (cs%nfiles > 0) then
143  inquire(file=trim(filename), number=iounit)
144  if (iounit /= -1) then
145  do i = 1, cs%nfiles
146  if (cs%iounit(i) == iounit) then
147  if (trim(cs%filename(1)) /= trim(filename)) then
148  call mom_error(fatal, &
149  "open_param_file: internal inconsistency! "//trim(filename)// &
150  " is registered as open but has the wrong unit number!")
151  else
152  call mom_error(warning, &
153  "open_param_file: file "//trim(filename)// &
154  " has already been opened. This should NOT happen!"// &
155  " Did you specify the same file twice in a namelist?")
156  return
157  endif ! filenames
158  endif ! unit numbers
159  enddo ! i
160  endif
161  endif
162 
163  ! Check that the file exists to readstdlog
164  inquire(file=trim(filename), exist=file_exists)
165  if (.not.file_exists) call mom_error(fatal, &
166  "open_param_file: Input file "// trim(filename)//" does not exist.")
167 
168  netcdf_file = .false.
169  if (strlen > 3) then
170  if (filename(strlen-2:strlen) == ".nc") netcdf_file = .true.
171  endif
172 
173  if (netcdf_file) &
174  call mom_error(fatal,"open_param_file: NetCDF files are not yet supported.")
175 
176  if (all_pes_read .or. is_root_pe()) then
177  ! Find an unused unit number.
178  do iounit=10,512
179  INQUIRE(iounit,opened=unit_in_use) ; if (.not.unit_in_use) exit
180  enddo
181  if (iounit >= 512) call mom_error(fatal, &
182  "open_param_file: No unused file unit could be found.")
183 
184  ! Open the parameter file.
185  open(iounit, file=trim(filename), access='SEQUENTIAL', &
186  form='FORMATTED', action='READ', position='REWIND', iostat=ios)
187  if (ios /= 0) call mom_error(fatal, "open_param_file: Error opening "// &
188  trim(filename))
189  else
190  iounit = 1
191  endif
192 
193  ! Store/register the unit and details
194  i = cs%nfiles + 1
195  cs%nfiles = i
196  cs%iounit(i) = iounit
197  cs%filename(i) = filename
198  cs%NetCDF_file(i) = netcdf_file
199  allocate(block) ; block%name = '' ; cs%blockName => block
200 
201  call mom_mesg("open_param_file: "// trim(filename)// &
202  " has been opened successfully.", 5)
203 
204  call populate_param_data(iounit, filename, cs%param_data(i))
205 
206  call read_param(cs,"SEND_LOG_TO_STDOUT",cs%log_to_stdout)
207  call read_param(cs,"REPORT_UNUSED_PARAMS",cs%report_unused)
208  call read_param(cs,"FATAL_UNUSED_PARAMS",cs%unused_params_fatal)
209  cs%doc_file = "MOM_parameter_doc"
210  if (present(component)) cs%doc_file = trim(component)//"_parameter_doc"
211  call read_param(cs,"DOCUMENT_FILE", cs%doc_file)
212  if (.not.may_check) then
213  cs%report_unused = .false.
214  cs%unused_params_fatal = .false.
215  endif
216 
217  ! Open the log file.
218  cs%stdlog = stdlog() ; cs%stdout = stdout()
219  cs%log_open = (stdlog() > 0)
220 
221  doc_path = cs%doc_file
222  if (len_trim(cs%doc_file) > 0) then
223  cs%complete_doc = complete_doc_default
224  call read_param(cs, "COMPLETE_DOCUMENTATION", cs%complete_doc)
225  cs%minimal_doc = minimal_doc_default
226  call read_param(cs, "MINIMAL_DOCUMENTATION", cs%minimal_doc)
227  if (present(doc_file_dir)) then ; if (len_trim(doc_file_dir) > 0) then
228  doc_path = trim(slasher(doc_file_dir))//trim(cs%doc_file)
229  endif ; endif
230  else
231  cs%complete_doc = .false.
232  cs%minimal_doc = .false.
233  endif
234  call doc_init(doc_path, cs%doc, minimal=cs%minimal_doc, complete=cs%complete_doc, &
235  layout=cs%complete_doc, debugging=cs%complete_doc)
236 
237 end subroutine open_param_file
238 
239 !> Close any open input files and deallocate memory associated with this param_file_type.
240 !! To use this type again, open_param_file would have to be called again.
241 subroutine close_param_file(CS, quiet_close, component)
242  type(param_file_type), intent(inout) :: CS !< The control structure for the file_parser module,
243  !! it is also a structure to parse for run-time parameters
244  logical, optional, intent(in) :: quiet_close !< if present and true, do not do any
245  !! logging with this call.
246  character(len=*), optional, intent(in) :: component !< If present, this component name is used
247  !! to generate parameter documentation file names
248  ! Local variables
249  logical :: all_default
250  character(len=128) :: docfile_default
251  character(len=40) :: mdl ! This module's name.
252  ! This include declares and sets the variable "version".
253 # include "version_variable.h"
254  integer :: i, n, num_unused
255 
256  if (present(quiet_close)) then ; if (quiet_close) then
257  do i = 1, cs%nfiles
258  if (all_pes_read .or. is_root_pe()) close(cs%iounit(i))
259  call mom_mesg("close_param_file: "// trim(cs%filename(i))// &
260  " has been closed successfully.", 5)
261  cs%iounit(i) = -1
262  cs%filename(i) = ''
263  cs%NetCDF_file(i) = .false.
264  deallocate (cs%param_data(i)%line)
265  deallocate (cs%param_data(i)%line_used)
266  enddo
267  cs%log_open = .false.
268  call doc_end(cs%doc)
269  return
270  endif ; endif
271 
272  ! Log the parameters for the parser.
273  docfile_default = "MOM_parameter_doc"
274  if (present(component)) docfile_default = trim(component)//"_parameter_doc"
275 
276  all_default = (cs%log_to_stdout .eqv. log_to_stdout_default)
277  all_default = all_default .and. (trim(cs%doc_file) == trim(docfile_default))
278  if (len_trim(cs%doc_file) > 0) then
279  all_default = all_default .and. (cs%complete_doc .eqv. complete_doc_default)
280  all_default = all_default .and. (cs%minimal_doc .eqv. minimal_doc_default)
281  endif
282 
283  mdl = "MOM_file_parser"
284  call log_version(cs, mdl, version, "", debugging=.true., log_to_all=.true., all_default=all_default)
285  call log_param(cs, mdl, "SEND_LOG_TO_STDOUT", cs%log_to_stdout, &
286  "If true, all log messages are also sent to stdout.", &
287  default=log_to_stdout_default)
288  call log_param(cs, mdl, "REPORT_UNUSED_PARAMS", cs%report_unused, &
289  "If true, report any parameter lines that are not used "//&
290  "in the run.", default=report_unused_default, &
291  debuggingparam=.true.)
292  call log_param(cs, mdl, "FATAL_UNUSED_PARAMS", cs%unused_params_fatal, &
293  "If true, kill the run if there are any unused "//&
294  "parameters.", default=unused_params_fatal_default, &
295  debuggingparam=.true.)
296  call log_param(cs, mdl, "DOCUMENT_FILE", cs%doc_file, &
297  "The basename for files where run-time parameters, their "//&
298  "settings, units and defaults are documented. Blank will "//&
299  "disable all parameter documentation.", default=docfile_default)
300  if (len_trim(cs%doc_file) > 0) then
301  call log_param(cs, mdl, "COMPLETE_DOCUMENTATION", cs%complete_doc, &
302  "If true, all run-time parameters are "//&
303  "documented in "//trim(cs%doc_file)//&
304  ".all .", default=complete_doc_default)
305  call log_param(cs, mdl, "MINIMAL_DOCUMENTATION", cs%minimal_doc, &
306  "If true, non-default run-time parameters are "//&
307  "documented in "//trim(cs%doc_file)//&
308  ".short .", default=minimal_doc_default)
309  endif
310 
311  num_unused = 0
312  do i = 1, cs%nfiles
313  if (is_root_pe() .and. (cs%report_unused .or. &
314  cs%unused_params_fatal)) then
315  ! Check for unused lines.
316  do n=1,cs%param_data(i)%num_lines
317  if (.not.cs%param_data(i)%line_used(n)) then
318  num_unused = num_unused + 1
319  if (cs%report_unused) &
320  call mom_error(warning, "Unused line in "//trim(cs%filename(i))// &
321  " : "//trim(cs%param_data(i)%line(n)))
322  endif
323  enddo
324  endif
325 
326  if (all_pes_read .or. is_root_pe()) close(cs%iounit(i))
327  call mom_mesg("close_param_file: "// trim(cs%filename(i))// &
328  " has been closed successfully.", 5)
329  cs%iounit(i) = -1
330  cs%filename(i) = ''
331  cs%NetCDF_file(i) = .false.
332  deallocate (cs%param_data(i)%line)
333  deallocate (cs%param_data(i)%line_used)
334  enddo
335 
336  if (is_root_pe() .and. (num_unused>0) .and. cs%unused_params_fatal) &
337  call mom_error(fatal, "Run stopped because of unused parameter lines.")
338 
339  cs%log_open = .false.
340  call doc_end(cs%doc)
341 
342 end subroutine close_param_file
343 
344 !> Read the contents of a parameter input file, and store the contents in a
345 !! file_data_type after removing comments and simplifying white space
346 subroutine populate_param_data(iounit, filename, param_data)
347  integer, intent(in) :: iounit !< The IO unit number that is open for filename
348  character(len=*), intent(in) :: filename !< An input file name, optionally with the full path
349  type(file_data_type), intent(inout) :: param_data !< A list of the input lines that set parameters
350  !! after comments have been stripped out.
351 
352  ! Local variables
353  character(len=INPUT_STR_LENGTH) :: line
354  integer :: num_lines
355  logical :: inMultiLineComment
356 
357  ! Find the number of keyword lines in a parameter file
358  ! Allocate the space to hold the lines in param_data%line
359  ! Populate param_data%line with the keyword lines from parameter file
360 
361  if (iounit <= 0) return
362 
363  if (all_pes_read .or. is_root_pe()) then
364  ! rewind the parameter file
365  rewind(iounit)
366 
367  ! count the number of valid entries in the parameter file
368  num_lines = 0
369  inmultilinecomment = .false.
370  do while(.true.)
371  read(iounit, '(a)', end=8, err=9) line
372  line = replacetabs(line)
373  if (inmultilinecomment) then
374  if (closemultilinecomment(line)) inmultilinecomment=.false.
375  else
376  if (lastnoncommentnonblank(line)>0) num_lines = num_lines + 1
377  if (openmultilinecomment(line)) inmultilinecomment=.true.
378  endif
379  enddo ! while (.true.)
380  8 continue ! get here when read() reaches EOF
381 
382  if (inmultilinecomment .and. is_root_pe()) &
383  call mom_error(fatal, 'MOM_file_parser : A C-style multi-line comment '// &
384  '(/* ... */) was not closed before the end of '//trim(filename))
385 
386  ! allocate space to hold contents of the parameter file
387  param_data%num_lines = num_lines
388  endif ! (is_root_pe())
389 
390  ! Broadcast the number of valid entries in parameter file
391  if (.not. all_pes_read) then
392  call broadcast(param_data%num_lines, root_pe())
393  endif
394 
395  ! Set up the space for storing the actual lines.
396  num_lines = param_data%num_lines
397  allocate (param_data%line(num_lines))
398  allocate (param_data%line_used(num_lines))
399  param_data%line(:) = ' '
400  param_data%line_used(:) = .false.
401 
402  ! Read the actual lines.
403  if (all_pes_read .or. is_root_pe()) then
404  ! rewind the parameter file
405  rewind(iounit)
406 
407  ! Populate param_data%line
408  num_lines = 0
409  do while(.true.)
410  read(iounit, '(a)', end=18, err=9) line
411  line = replacetabs(line)
412  if (inmultilinecomment) then
413  if (closemultilinecomment(line)) inmultilinecomment=.false.
414  else
415  if (lastnoncommentnonblank(line)>0) then
416  line = removecomments(line)
417  line = simplifywhitespace(line(:len_trim(line)))
418  num_lines = num_lines + 1
419  param_data%line(num_lines) = line
420  endif
421  if (openmultilinecomment(line)) inmultilinecomment=.true.
422  endif
423  enddo ! while (.true.)
424 18 continue ! get here when read() reaches EOF
425 
426  if (num_lines /= param_data%num_lines) &
427  call mom_error(fatal, 'MOM_file_parser : Found different number of '// &
428  'valid lines on second reading of '//trim(filename))
429  endif ! (is_root_pe())
430 
431  ! Broadcast the populated array param_data%line
432  if (.not. all_pes_read) then
433  call broadcast(param_data%line, input_str_length, root_pe())
434  endif
435 
436  return
437 
438 9 call mom_error(fatal, "MOM_file_parser : "//&
439  "Error while reading file "//trim(filename))
440 
441 end subroutine populate_param_data
442 
443 
444 !> Return True if a /* appears on this line without a closing */
445 function openmultilinecomment(string)
446  character(len=*), intent(in) :: string !< The input string to process
447  logical :: openMultiLineComment
448 
449  ! Local variables
450  integer :: icom, last
451 
452  openmultilinecomment = .false.
453  last = lastnoncommentindex(string)+1
454  icom = index(string(last:), "/*")
455  if (icom > 0) then
456  openmultilinecomment=.true.
457  last = last+icom+1
458  endif
459  icom = index(string(last:), "*/") ; if (icom > 0) openmultilinecomment=.false.
460 end function openmultilinecomment
461 
462 !> Return True if a */ appears on this line
463 function closemultilinecomment(string)
464  character(len=*), intent(in) :: string !< The input string to process
465  logical :: closeMultiLineComment
466 ! True if a */ appears on this line
467  closemultilinecomment = .false.
468  if (index(string, "*/")>0) closemultilinecomment=.true.
469 end function closemultilinecomment
470 
471 !> Find position of last character before any comments, As marked by "!", "//", or "/*"
472 !! following F90, C++, or C syntax
473 function lastnoncommentindex(string)
474  character(len=*), intent(in) :: string !< The input string to process
475  integer :: lastNonCommentIndex
476 
477  ! Local variables
478  integer :: icom, last
479 
480  ! This subroutine is the only place where a comment needs to be defined
481  last = len_trim(string)
482  icom = index(string(:last), "!") ; if (icom > 0) last = icom-1 ! F90 style
483  icom = index(string(:last), "//") ; if (icom > 0) last = icom-1 ! C++ style
484  icom = index(string(:last), "/*") ; if (icom > 0) last = icom-1 ! C style
485  lastnoncommentindex = last
486 end function lastnoncommentindex
487 
488 !> Find position of last non-blank character before any comments
489 function lastnoncommentnonblank(string)
490  character(len=*), intent(in) :: string !< The input string to process
491  integer :: lastNonCommentNonBlank
492 
493  lastnoncommentnonblank = len_trim(string(:lastnoncommentindex(string))) ! Ignore remaining trailing blanks
494 end function lastnoncommentnonblank
495 
496 !> Returns a string with tabs replaced by a blank
497 function replacetabs(string)
498  character(len=*), intent(in) :: string !< The input string to process
499  character(len=len(string)) :: replaceTabs
500 
501  integer :: i
502 
503  do i=1, len(string)
504  if (string(i:i)==achar(9)) then
505  replacetabs(i:i)=" "
506  else
507  replacetabs(i:i)=string(i:i)
508  endif
509  enddo
510 end function replacetabs
511 
512 !> Trims comments and leading blanks from string
513 function removecomments(string)
514  character(len=*), intent(in) :: string !< The input string to process
515  character(len=len(string)) :: removeComments
516 
517  integer :: last
518 
519  removecomments=repeat(" ",len(string))
520  last = lastnoncommentnonblank(string)
521  removecomments(:last)=adjustl(string(:last)) ! Copy only the non-comment part of string
522 end function removecomments
523 
524 !> Constructs a string with all repeated whitespace replaced with single blanks
525 !! and insert white space where it helps delineate tokens (e.g. around =)
526 function simplifywhitespace(string)
527  character(len=*), intent(in) :: string !< A string to modify to simpify white space
528  character(len=len(string)+16) :: simplifyWhiteSpace
529 
530  ! Local variables
531  integer :: i,j
532  logical :: nonBlank = .false., insidestring = .false.
533  character(len=1) :: quoteChar=" "
534 
535  nonblank = .false.; insidestring = .false. ! NOTE: For some reason this line is needed??
536  i=0
537  simplifywhitespace=repeat(" ",len(string)+16)
538  do j=1,len_trim(string)
539  if (insidestring) then ! Do not change formatting inside strings
540  i=i+1
541  simplifywhitespace(i:i)=string(j:j)
542  if (string(j:j)==quotechar) insidestring=.false. ! End of string
543  else ! The following is outside of string delimiters
544  if (string(j:j)==" " .or. string(j:j)==achar(9)) then ! Space or tab
545  if (nonblank) then ! Only copy a blank if the preceeding character was non-blank
546  i=i+1
547  simplifywhitespace(i:i)=" " ! Not string(j:j) so that tabs are replace by blanks
548  nonblank=.false.
549  endif
550  elseif (string(j:j)=='"' .or. string(j:j)=="'") then ! Start a sting
551  i=i+1
552  simplifywhitespace(i:i)=string(j:j)
553  insidestring=.true.
554  quotechar=string(j:j) ! Keep copy of starting quote
555  nonblank=.true. ! For exit from string
556  elseif (string(j:j)=='=') then
557  ! Insert spaces if this character is "=" so that line contains " = "
558  if (nonblank) then
559  i=i+1
560  simplifywhitespace(i:i)=" "
561  endif
562  i=i+2
563  simplifywhitespace(i-1:i)=string(j:j)//" "
564  nonblank=.false.
565  else ! All other characters
566  i=i+1
567  simplifywhitespace(i:i)=string(j:j)
568  nonblank=.true.
569  endif
570  endif ! if (insideString)
571  enddo ! j
572  if (insidestring) then ! A missing close quote should be flagged
573  if (is_root_pe()) call mom_error(fatal, &
574  "There is a mismatched quote in the parameter file line: "// &
575  trim(string))
576  endif
577 end function simplifywhitespace
578 
579 !> This subroutine reads the value of an integer model parameter from a parameter file.
580 subroutine read_param_int(CS, varname, value, fail_if_missing)
581  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
582  !! it is also a structure to parse for run-time parameters
583  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
584  integer, intent(inout) :: value !< The value of the parameter that may be
585  !! read from the parameter file
586  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
587  !! if this variable is not found in the parameter file
588  ! Local variables
589  character(len=INPUT_STR_LENGTH) :: value_string(1)
590  logical :: found, defined
591 
592  call get_variable_line(cs, varname, found, defined, value_string)
593  if (found .and. defined .and. (len_trim(value_string(1)) > 0)) then
594  read(value_string(1),*,err = 1001) value
595  else
596  if (present(fail_if_missing)) then ; if (fail_if_missing) then
597  if (.not.found) then
598  call mom_error(fatal,'read_param_int: Unable to find variable '//trim(varname)// &
599  ' in any input files.')
600  else
601  call mom_error(fatal,'read_param_int: Variable '//trim(varname)// &
602  ' found but not set in input files.')
603  endif
604  endif ; endif
605  endif
606  return
607  1001 call mom_error(fatal,'read_param_int: read error for integer variable '//trim(varname)// &
608  ' parsing "'//trim(value_string(1))//'"')
609 end subroutine read_param_int
610 
611 !> This subroutine reads the values of an array of integer model parameters from a parameter file.
612 subroutine read_param_int_array(CS, varname, value, fail_if_missing)
613  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
614  !! it is also a structure to parse for run-time parameters
615  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
616  integer, dimension(:), intent(inout) :: value !< The value of the parameter that may be
617  !! read from the parameter file
618  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
619  !! if this variable is not found in the parameter file
620  ! Local variables
621  character(len=INPUT_STR_LENGTH) :: value_string(1)
622  logical :: found, defined
623 
624  call get_variable_line(cs, varname, found, defined, value_string)
625  if (found .and. defined .and. (len_trim(value_string(1)) > 0)) then
626  read(value_string(1),*,end=991,err=1002) value
627  991 return
628  else
629  if (present(fail_if_missing)) then ; if (fail_if_missing) then
630  if (.not.found) then
631  call mom_error(fatal,'read_param_int_array: Unable to find variable '//trim(varname)// &
632  ' in any input files.')
633  else
634  call mom_error(fatal,'read_param_int_array: Variable '//trim(varname)// &
635  ' found but not set in input files.')
636  endif
637  endif ; endif
638  endif
639  return
640  1002 call mom_error(fatal,'read_param_int_array: read error for integer array '//trim(varname)// &
641  ' parsing "'//trim(value_string(1))//'"')
642 end subroutine read_param_int_array
643 
644 !> This subroutine reads the value of a real model parameter from a parameter file.
645 subroutine read_param_real(CS, varname, value, fail_if_missing, scale)
646  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
647  !! it is also a structure to parse for run-time parameters
648  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
649  real, intent(inout) :: value !< The value of the parameter that may be
650  !! read from the parameter file
651  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
652  !! if this variable is not found in the parameter file
653  real, optional, intent(in) :: scale !< A scaling factor that the parameter is multiplied
654  !! by before it is returned.
655 
656  ! Local variables
657  character(len=INPUT_STR_LENGTH) :: value_string(1)
658  logical :: found, defined
659 
660  call get_variable_line(cs, varname, found, defined, value_string)
661  if (found .and. defined .and. (len_trim(value_string(1)) > 0)) then
662  read(value_string(1),*,err=1003) value
663  if (present(scale)) value = scale*value
664  else
665  if (present(fail_if_missing)) then ; if (fail_if_missing) then
666  if (.not.found) then
667  call mom_error(fatal,'read_param_real: Unable to find variable '//trim(varname)// &
668  ' in any input files.')
669  else
670  call mom_error(fatal,'read_param_real: Variable '//trim(varname)// &
671  ' found but not set in input files.')
672  endif
673  endif ; endif
674  endif
675  return
676  1003 call mom_error(fatal,'read_param_real: read error for real variable '//trim(varname)// &
677  ' parsing "'//trim(value_string(1))//'"')
678 end subroutine read_param_real
679 
680 !> This subroutine reads the values of an array of real model parameters from a parameter file.
681 subroutine read_param_real_array(CS, varname, value, fail_if_missing, scale)
682  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
683  !! it is also a structure to parse for run-time parameters
684  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
685  real, dimension(:), intent(inout) :: value !< The value of the parameter that may be
686  !! read from the parameter file
687  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
688  !! if this variable is not found in the parameter file
689  real, optional, intent(in) :: scale !< A scaling factor that the parameter is multiplied
690  !! by before it is returned.
691 
692  ! Local variables
693  character(len=INPUT_STR_LENGTH) :: value_string(1)
694  logical :: found, defined
695 
696  call get_variable_line(cs, varname, found, defined, value_string)
697  if (found .and. defined .and. (len_trim(value_string(1)) > 0)) then
698  read(value_string(1),*,end=991,err=1004) value
699 991 continue
700  if (present(scale)) value(:) = scale*value(:)
701  return
702  else
703  if (present(fail_if_missing)) then ; if (fail_if_missing) then
704  if (.not.found) then
705  call mom_error(fatal,'read_param_real_array: Unable to find variable '//trim(varname)// &
706  ' in any input files.')
707  else
708  call mom_error(fatal,'read_param_real_array: Variable '//trim(varname)// &
709  ' found but not set in input files.')
710  endif
711  endif ; endif
712  endif
713  return
714  1004 call mom_error(fatal,'read_param_real_array: read error for real array '//trim(varname)// &
715  ' parsing "'//trim(value_string(1))//'"')
716 end subroutine read_param_real_array
717 
718 !> This subroutine reads the value of a character string model parameter from a parameter file.
719 subroutine read_param_char(CS, varname, value, fail_if_missing)
720  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
721  !! it is also a structure to parse for run-time parameters
722  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
723  character(len=*), intent(inout) :: value !< The value of the parameter that may be
724  !! read from the parameter file
725  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
726  !! if this variable is not found in the parameter file
727  ! Local variables
728  character(len=INPUT_STR_LENGTH) :: value_string(1)
729  logical :: found, defined
730 
731  call get_variable_line(cs, varname, found, defined, value_string)
732  if (found) then
733  value = trim(strip_quotes(value_string(1)))
734  elseif (present(fail_if_missing)) then ; if (fail_if_missing) then
735  call mom_error(fatal,'Unable to find variable '//trim(varname)// &
736  ' in any input files.')
737  endif ; endif
738 
739 end subroutine read_param_char
740 
741 !> This subroutine reads the values of an array of character string model parameters from a parameter file.
742 subroutine read_param_char_array(CS, varname, value, fail_if_missing)
743  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
744  !! it is also a structure to parse for run-time parameters
745  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
746  character(len=*), dimension(:), intent(inout) :: value !< The value of the parameter that may be
747  !! read from the parameter file
748  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
749  !! if this variable is not found in the parameter file
750 
751  ! Local variables
752  character(len=INPUT_STR_LENGTH) :: value_string(1), loc_string
753  logical :: found, defined
754  integer :: i, i_out
755 
756  call get_variable_line(cs, varname, found, defined, value_string)
757  if (found) then
758  loc_string = trim(value_string(1))
759  i = index(loc_string,",")
760  i_out = 1
761  do while(i>0)
762  value(i_out) = trim(strip_quotes(loc_string(:i-1)))
763  i_out = i_out+1
764  loc_string = trim(adjustl(loc_string(i+1:)))
765  i = index(loc_string,",")
766  enddo
767  if (len_trim(loc_string)>0) then
768  value(i_out) = trim(strip_quotes(adjustl(loc_string)))
769  i_out = i_out+1
770  endif
771  do i=i_out,SIZE(value) ; value(i) = " " ; enddo
772  elseif (present(fail_if_missing)) then ; if (fail_if_missing) then
773  call mom_error(fatal,'Unable to find variable '//trim(varname)// &
774  ' in any input files.')
775  endif ; endif
776 
777 end subroutine read_param_char_array
778 
779 !> This subroutine reads the value of a logical model parameter from a parameter file.
780 subroutine read_param_logical(CS, varname, value, fail_if_missing)
781  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
782  !! it is also a structure to parse for run-time parameters
783  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
784  logical, intent(inout) :: value !< The value of the parameter that may be
785  !! read from the parameter file
786  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
787  !! if this variable is not found in the parameter file
788 
789  ! Local variables
790  character(len=INPUT_STR_LENGTH) :: value_string(1)
791  logical :: found, defined
792 
793  call get_variable_line(cs, varname, found, defined, value_string, paramislogical=.true.)
794  if (found) then
795  value = defined
796  elseif (present(fail_if_missing)) then ; if (fail_if_missing) then
797  call mom_error(fatal,'Unable to find variable '//trim(varname)// &
798  ' in any input files.')
799  endif ; endif
800 end subroutine read_param_logical
801 
802 !> This subroutine reads the value of a time_type model parameter from a parameter file.
803 subroutine read_param_time(CS, varname, value, timeunit, fail_if_missing, date_format)
804  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
805  !! it is also a structure to parse for run-time parameters
806  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
807  type(time_type), intent(inout) :: value !< The value of the parameter that may be
808  !! read from the parameter file
809  real, optional, intent(in) :: timeunit !< The number of seconds in a time unit for real-number input.
810  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
811  !! if this variable is not found in the parameter file
812  logical, optional, intent(out) :: date_format !< If present, this indicates whether this
813  !! parameter was read in a date format, so that it can
814  !! later be logged in the same format.
815 
816  ! Local variables
817  character(len=INPUT_STR_LENGTH) :: value_string(1)
818  character(len=240) :: err_msg
819  logical :: found, defined
820  real :: real_time, time_unit
821  integer :: vals(7)
822 
823  if (present(date_format)) date_format = .false.
824 
825  call get_variable_line(cs, varname, found, defined, value_string)
826  if (found .and. defined .and. (len_trim(value_string(1)) > 0)) then
827  ! Determine whether value string should be parsed for a real number
828  ! or a date, in either a string format or a comma-delimited list of values.
829  if ((index(value_string(1),'-') > 0) .and. &
830  (index(value_string(1),'-',back=.true.) > index(value_string(1),'-'))) then
831  ! There are two dashes, so this must be a date format.
832  value = set_date(value_string(1), err_msg=err_msg)
833  if (len_trim(err_msg) > 0) call mom_error(fatal,'read_param_time: '//&
834  trim(err_msg)//' in integer list read error for time-type variable '//&
835  trim(varname)// ' parsing "'//trim(value_string(1))//'"')
836  if (present(date_format)) date_format = .true.
837  elseif (index(value_string(1),',') > 0) then
838  ! Initialize vals with an invalid date.
839  vals(:) = (/ -999, -999, -999, 0, 0, 0, 0 /)
840  read(value_string(1),*,end=995,err=1005) vals
841  995 continue
842  if ((vals(1) < 0) .or. (vals(2) < 0) .or. (vals(3) < 0)) &
843  call mom_error(fatal,'read_param_time: integer list read error for time-type variable '//&
844  trim(varname)// ' parsing "'//trim(value_string(1))//'"')
845  value = set_date(vals(1), vals(2), vals(3), vals(4), vals(5), vals(6), &
846  vals(7), err_msg=err_msg)
847  if (len_trim(err_msg) > 0) call mom_error(fatal,'read_param_time: '//&
848  trim(err_msg)//' in integer list read error for time-type variable '//&
849  trim(varname)// ' parsing "'//trim(value_string(1))//'"')
850  if (present(date_format)) date_format = .true.
851  else
852  time_unit = 1.0 ; if (present(timeunit)) time_unit = timeunit
853  read( value_string(1), *) real_time
854  value = real_to_time(real_time*time_unit)
855  endif
856  else
857  if (present(fail_if_missing)) then ; if (fail_if_missing) then
858  if (.not.found) then
859  call mom_error(fatal,'Unable to find variable '//trim(varname)// &
860  ' in any input files.')
861  else
862  call mom_error(fatal,'Variable '//trim(varname)// &
863  ' found but not set in input files.')
864  endif
865  endif ; endif
866  endif
867  return
868  1005 call mom_error(fatal,'read_param_time: read error for time-type variable '//&
869  trim(varname)// ' parsing "'//trim(value_string(1))//'"')
870 end subroutine read_param_time
871 
872 !> This function removes single and double quotes from a character string
873 function strip_quotes(val_str)
874  character(len=*) :: val_str !< The character string to work on
875  character(len=INPUT_STR_LENGTH) :: strip_quotes
876  ! Local variables
877  integer :: i
878  strip_quotes = val_str
879  i = index(strip_quotes,achar(34)) ! Double quote
880  do while (i>0)
881  if (i > 1) then ; strip_quotes = strip_quotes(:i-1)//strip_quotes(i+1:)
882  else ; strip_quotes = strip_quotes(2:) ; endif
883  i = index(strip_quotes,achar(34)) ! Double quote
884  enddo
885  i = index(strip_quotes,achar(39)) ! Single quote
886  do while (i>0)
887  if (i > 1) then ; strip_quotes = strip_quotes(:i-1)//strip_quotes(i+1:)
888  else ; strip_quotes = strip_quotes(2:) ; endif
889  i = index(strip_quotes,achar(39)) ! Single quote
890  enddo
891 end function strip_quotes
892 
893 !> This subtoutine extracts the contents of lines in the param_file_type that refer to
894 !! a named parameter. The value_string that is returned must be interepreted in a way
895 !! that depends on the type of this variable.
896 subroutine get_variable_line(CS, varname, found, defined, value_string, paramIsLogical)
897  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
898  !! it is also a structure to parse for run-time parameters
899  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
900  logical, intent(out) :: found !< If true, this parameter has been found in CS
901  logical, intent(out) :: defined !< If true, this parameter is set (or true) in the CS
902  character(len=*), intent(out) :: value_string(:) !< A string that encodes the new value
903  logical, optional, intent(in) :: paramIsLogical !< If true, this is a logical parameter
904  !! that can be simply defined without parsing a value_string.
905 
906  ! Local variables
907  character(len=INPUT_STR_LENGTH) :: val_str, lname, origLine
908  character(len=INPUT_STR_LENGTH) :: line, continuationBuffer, blockName
909  character(len=FILENAME_LENGTH) :: filename
910  integer :: is, id, isd, isu, ise, iso, verbose, ipf
911  integer :: last, last1, ival, oval, max_vals, count, contBufSize
912  character(len=52) :: set
913  logical :: found_override, found_equals
914  logical :: found_define, found_undef
915  logical :: force_cycle, defined_in_line, continuedLine
916  logical :: variableKindIsLogical, valueIsSame
917  logical :: inWrongBlock, fullPathParameter
918  logical, parameter :: requireNamedClose = .false.
919  set = "ABCDEFGHIJKLMNOPQRSTUVWXYZabcdefghijklmnopqrstuvwxyz"
920  continuationbuffer = repeat(" ",input_str_length)
921  contbufsize = 0
922  verbose = 1
923 
924  variablekindislogical=.false.
925  if (present(paramislogical)) variablekindislogical = paramislogical
926 
927  ! Find the first instance (if any) where the named variable is found, and
928  ! return variables indicating whether this variable is defined and the string
929  ! that contains the value of this variable.
930  found = .false.
931  oval = 0; ival = 0
932  max_vals = SIZE(value_string)
933  do is=1,max_vals ; value_string(is) = " " ; enddo
934 
935  paramfile_loop: do ipf = 1, cs%nfiles
936  filename = cs%filename(ipf)
937  continuedline = .false.
938  blockname = ''
939 
940  ! Scan through each line of the file
941  do count = 1, cs%param_data(ipf)%num_lines
942  line = cs%param_data(ipf)%line(count)
943  last = len_trim(line)
944 
945  last1 = max(1,last)
946  ! Check if line ends in continuation character (either & or \)
947  ! Note achar(92) is a backslash
948  if (line(last1:last1) == achar(92).or.line(last1:last1) == "&") then
949  continuationbuffer(contbufsize+1:contbufsize+len_trim(line))=line(:last-1)
950  contbufsize=contbufsize + len_trim(line)-1
951  continuedline = .true.
952  if (count==cs%param_data(ipf)%num_lines .and. is_root_pe()) &
953  call mom_error(fatal, "MOM_file_parser : the last line"// &
954  " of the file ends in a continuation character but"// &
955  " there are no more lines to read. "// &
956  " Line: '"//trim(line(:last))//"'"//&
957  " in file "//trim(filename)//".")
958  cycle ! cycle inorder to append the next line of the file
959  elseif (continuedline) then
960  ! If we reached this point then this is the end of line continuation
961  continuationbuffer(contbufsize+1:contbufsize+len_trim(line))=line(:last)
962  line = continuationbuffer
963  continuationbuffer=repeat(" ",input_str_length) ! Clear for next use
964  contbufsize = 0
965  continuedline = .false.
966  last = len_trim(line)
967  endif
968 
969  origline = trim(line) ! Keep original for error messages
970 
971  ! Check for '#override' at start of line
972  found_override = .false.; found_define = .false.; found_undef = .false.
973  iso = index(line(:last), "#override " )!; if (is > 0) found_override = .true.
974  if (iso>1) call mom_error(fatal, "MOM_file_parser : #override was found "// &
975  " but was not the first keyword."// &
976  " Line: '"//trim(line(:last))//"'"//&
977  " in file "//trim(filename)//".")
978  if (iso==1) then
979  found_override = .true.
980  if (index(line(:last), "#override define ")==1) found_define = .true.
981  if (index(line(:last), "#override undef ")==1) found_undef = .true.
982  line = trim(adjustl(line(iso+10:last))); last = len_trim(line)
983  endif
984 
985  ! Check for start of fortran namelist, ie. '&namelist'
986  if (index(line(:last),'&')==1) then
987  iso=index(line(:last),' ')
988  if (iso>0) then ! possibly simething else on this line
989  blockname = pushblocklevel(blockname,line(2:iso-1))
990  line=trim(adjustl(line(iso:last)))
991  last=len_trim(line)
992  if (last==0) cycle ! nothing else on this line
993  else ! just the namelist on this line
994  if (len_trim(blockname)>0) then
995  blockname = trim(blockname) // '%' //trim(line(2:last))
996  else
997  blockname = trim(line(2:last))
998  endif
999  call flag_line_as_read(cs%param_data(ipf)%line_used,count)
1000  cycle
1001  endif
1002  endif
1003 
1004  ! Newer form of parameter block, block%, %block or block%param or
1005  iso=index(line(:last),'%')
1006  fullpathparameter = .false.
1007  if (iso==1) then ! % is first character means this is a close
1008  if (len_trim(blockname)==0 .and. is_root_pe()) call mom_error(fatal, &
1009  'get_variable_line: An extra close block was encountered. Line="'// &
1010  trim(line(:last))//'"' )
1011  if (last>1 .and. trim(blockname)/=trim(line(2:last)) .and. is_root_pe()) &
1012  call mom_error(fatal, 'get_variable_line: A named close for a parameter'// &
1013  ' block did not match the open block. Line="'//trim(line(:last))//'"' )
1014  if (last==1 .and. requirenamedclose) & ! line = '%' is a generic (unnamed) close
1015  call mom_error(fatal, 'get_variable_line: A named close for a parameter'// &
1016  ' block is required but found "%". Block="'//trim(blockname)//'"' )
1017  blockname = popblocklevel(blockname)
1018  call flag_line_as_read(cs%param_data(ipf)%line_used,count)
1019  elseif (iso==last) then ! This is a new block if % is last character
1020  blockname = pushblocklevel(blockname, line(:iso-1))
1021  call flag_line_as_read(cs%param_data(ipf)%line_used,count)
1022  else ! This is of the form block%parameter = ... (full path parameter)
1023  iso=index(line(:last),'%',.true.)
1024  ! Check that the parameter block names on the line matches the state set by the caller
1025  if (iso>0 .and. trim(cs%blockName%name)==trim(line(:iso-1))) then
1026  fullpathparameter = .true.
1027  line = trim(line(iso+1:last)) ! Strip away the block name for subsequent processing
1028  last = len_trim(line)
1029  endif
1030  endif
1031 
1032  ! We should only interpret this line if this block is the active block
1033  inwrongblock = .false.
1034  if (len_trim(blockname)>0) then ! In a namelist block in file
1035  if (trim(cs%blockName%name)/=trim(blockname)) inwrongblock = .true. ! Not in the required block
1036  endif
1037  if (len_trim(cs%blockName%name)>0) then ! In a namelist block in the model
1038  if (trim(cs%blockName%name)/=trim(blockname)) inwrongblock = .true. ! Not in the required block
1039  endif
1040 
1041  ! Check for termination of a fortran namelist (with a '/')
1042  if (line(last:last)=='/') then
1043  if (len_trim(blockname)==0 .and. is_root_pe()) call mom_error(fatal, &
1044  'get_variable_line: An extra namelist/block end was encountered. Line="'// &
1045  trim(line(:last))//'"' )
1046  blockname = popblocklevel(blockname)
1047  last = last - 1 ! Ignore the termination character from here on
1048  endif
1049  if (inwrongblock .and. .not. fullpathparameter) then
1050  if (index(" "//line(:last+1), " "//trim(varname)//" ")>0) &
1051  call mom_error(warning,"MOM_file_parser : "//trim(varname)// &
1052  ' found outside of block '//trim(cs%blockName%name)//'%. Ignoring.')
1053  cycle
1054  endif
1055 
1056  ! Determine whether this line mentions the named parameter or not
1057  if (index(" "//line(:last)//" ", " "//trim(varname)//" ") == 0) cycle
1058 
1059  ! Detect keywords
1060  found_equals = .false.
1061  isd = index(line(:last), "define" )!; if (isd > 0) found_define = .true.
1062  isu = index(line(:last), "undef" )!; if (isu > 0) found_undef = .true.
1063  ise = index(line(:last), " = " ); if (ise > 1) found_equals = .true.
1064  if (index(line(:last), "#define ")==1) found_define = .true.
1065  if (index(line(:last), "#undef ")==1) found_undef = .true.
1066 
1067  ! Check for missing, mutually exclusive or incomplete keywords
1068  if (is_root_pe()) then
1069  if (.not. (found_define .or. found_undef .or. found_equals)) &
1070  call mom_error(fatal, "MOM_file_parser : the parameter name '"// &
1071  trim(varname)//"' was found without define or undef."// &
1072  " Line: '"//trim(line(:last))//"'"//&
1073  " in file "//trim(filename)//".")
1074  if (found_define .and. found_undef) call mom_error(fatal, &
1075  "MOM_file_parser : Both 'undef' and 'define' occur."// &
1076  " Line: '"//trim(line(:last))//"'"//&
1077  " in file "//trim(filename)//".")
1078  if (found_equals .and. (found_define .or. found_undef)) &
1079  call mom_error(fatal, &
1080  "MOM_file_parser : Both 'a=b' and 'undef/define' syntax occur."// &
1081  " Line: '"//trim(line(:last))//"'"//&
1082  " in file "//trim(filename)//".")
1083  if (found_override .and. .not. (found_define .or. found_undef .or. found_equals)) &
1084  call mom_error(fatal, "MOM_file_parser : override was found "// &
1085  " without a define or undef."// &
1086  " Line: '"//trim(line(:last))//"'"//&
1087  " in file "//trim(filename)//".")
1088  endif
1089 
1090  ! Interpret the line and collect values, if any
1091  if (found_define) then
1092  ! Move starting pointer to first letter of defined name.
1093  is = isd + 5 + scan(line(isd+6:last), set)
1094 
1095  id = scan(line(is:last), ' ') ! Find space between name and value
1096  if ( id == 0 ) then
1097  ! There is no space so the name is simply being defined.
1098  lname = trim(line(is:last))
1099  if (trim(lname) /= trim(varname)) cycle
1100  val_str = " "
1101  else
1102  ! There is a string or number after the name.
1103  lname = trim(line(is:is+id-1))
1104  if (trim(lname) /= trim(varname)) cycle
1105  val_str = trim(adjustl(line(is+id:last)))
1106  endif
1107  found = .true. ; defined_in_line = .true.
1108  elseif (found_undef) then
1109  ! Move starting pointer to first letter of undefined name.
1110  is = isu + 4 + scan(line(isu+5:last), set)
1111 
1112  id = scan(line(is:last), ' ') ! Find the first space after the name.
1113  if (id > 0) last = is + id - 1
1114  lname = trim(line(is:last))
1115  if (trim(lname) /= trim(varname)) cycle
1116  val_str = " "
1117  found = .true. ; defined_in_line = .false.
1118  elseif (found_equals) then
1119  ! Move starting pointer to first letter of defined name.
1120  is = scan(line(1:ise), set)
1121  lname = trim(line(is:ise-1))
1122  if (trim(lname) /= trim(varname)) cycle
1123  val_str = trim(adjustl(line(ise+3:last)))
1124  if (variablekindislogical) then ! Special handling for logicals
1125  read(val_str(:len_trim(val_str)),*) defined_in_line
1126  else
1127  defined_in_line = .true.
1128  endif
1129  found = .true.
1130  else
1131  call mom_error(fatal, "MOM_file_parser (non-root PE?): the parameter name '"// &
1132  trim(varname)//"' was found without an assignment, define or undef."// &
1133  " Line: '"//trim(line(:last))//"'"//" in file "//trim(filename)//".")
1134  endif
1135 
1136  ! This line has now been used.
1137  call flag_line_as_read(cs%param_data(ipf)%line_used,count)
1138 
1139  ! Detect inconsistencies
1140  force_cycle = .false.
1141  valueissame = (trim(val_str) == trim(value_string(max_vals)))
1142  if (found_override .and. (oval >= max_vals)) then
1143  if (is_root_pe()) then
1144  if ((defined_in_line .neqv. defined) .or. .not. valueissame) then
1145  call mom_error(fatal,"MOM_file_parser : "//trim(varname)// &
1146  " found with multiple inconsistent overrides."// &
1147  " Line A: '"//trim(value_string(max_vals))//"'"//&
1148  " Line B: '"//trim(line(:last))//"'"//&
1149  " in file "//trim(filename)//" caused the model failure.")
1150  else
1151  call mom_error(warning,"MOM_file_parser : "//trim(varname)// &
1152  " over-ridden more times than is permitted."// &
1153  " Line: '"//trim(line(:last))//"'"//&
1154  " in file "//trim(filename)//" is being ignored.")
1155  endif
1156  endif
1157  force_cycle = .true.
1158  endif
1159  if (.not.found_override .and. (oval > 0)) then
1160  if (is_root_pe()) &
1161  call mom_error(warning,"MOM_file_parser : "//trim(varname)// &
1162  " has already been over-ridden."// &
1163  " Line: '"//trim(line(:last))//"'"//&
1164  " in file "//trim(filename)//" is being ignored.")
1165  force_cycle = .true.
1166  endif
1167  if (.not.found_override .and. (ival >= max_vals)) then
1168  if (is_root_pe()) then
1169  if ((defined_in_line .neqv. defined) .or. .not. valueissame) then
1170  call mom_error(fatal,"MOM_file_parser : "//trim(varname)// &
1171  " found with multiple inconsistent definitions."// &
1172  " Line A: '"//trim(value_string(max_vals))//"'"//&
1173  " Line B: '"//trim(line(:last))//"'"//&
1174  " in file "//trim(filename)//" caused the model failure.")
1175  else
1176  call mom_error(warning,"MOM_file_parser : "//trim(varname)// &
1177  " occurs more times than is permitted."// &
1178  " Line: '"//trim(line(:last))//"'"//&
1179  " in file "//trim(filename)//" is being ignored.")
1180  endif
1181  endif
1182  force_cycle = .true.
1183  endif
1184  if (force_cycle) cycle
1185 
1186  ! Store new values
1187  if (found_override) then
1188  oval = oval + 1
1189  value_string(oval) = trim(val_str)
1190  defined = defined_in_line
1191  if (verbose > 0 .and. ival > 0 .and. is_root_pe() .and. &
1192  .not. overridewarninghasbeenissued(cs%chain, trim(varname)) ) &
1193  call mom_error(warning,"MOM_file_parser : "//trim(varname)// &
1194  " over-ridden. Line: '"//trim(line(:last))//"'"//&
1195  " in file "//trim(filename)//".")
1196  else ! (.not. found_overide)
1197  ival = ival + 1
1198  value_string(ival) = trim(val_str)
1199  defined = defined_in_line
1200  if (verbose > 1 .and. is_root_pe()) &
1201  call mom_error(warning,"MOM_file_parser : "//trim(varname)// &
1202  " set. Line: '"//trim(line(:last))//"'"//&
1203  " in file "//trim(filename)//".")
1204  endif
1205 
1206  enddo ! CS%param_data(ipf)%num_lines
1207 
1208  if (len_trim(blockname)>0 .and. is_root_pe()) call mom_error(fatal, &
1209  'A namelist/parameter block was not closed. Last open block appears '// &
1210  'to be "'//trim(blockname)//'".')
1211 
1212  enddo paramfile_loop
1213 
1214 end subroutine get_variable_line
1215 
1216 !> Record that a line has been used to set a parameter
1217 subroutine flag_line_as_read(line_used, count)
1218  logical, dimension(:), pointer :: line_used !< A structure indicating which lines have been read
1219  integer, intent(in) :: count !< The parameter on this line number has been read
1220  line_used(count) = .true.
1221 end subroutine flag_line_as_read
1222 
1223 !> Returns true if an override warning has been issued for the variable varName
1224 function overridewarninghasbeenissued(chain, varName)
1225  type(link_parameter), pointer :: chain !< The linked list of variables that have already had
1226  !! override warnings issued
1227  character(len=*), intent(in) :: varName !< The name of the variable being queried for warnings
1228  logical :: overrideWarningHasBeenIssued
1229  ! Local variables
1230  type(link_parameter), pointer :: newLink => null(), this => null()
1231 
1232  overridewarninghasbeenissued = .false.
1233  this => chain
1234  do while( associated(this) )
1235  if (trim(varname) == trim(this%name)) then
1236  overridewarninghasbeenissued = .true.
1237  return
1238  endif
1239  this => this%next
1240  enddo
1241  allocate(newlink)
1242  newlink%name = trim(varname)
1243  newlink%hasIssuedOverrideWarning = .true.
1244  newlink%next => chain
1245  chain => newlink
1246 end function overridewarninghasbeenissued
1247 
1248 ! The following subroutines write out to a log file.
1249 
1250 !> Log the version of a module to a log file and/or stdout, and/or to the
1251 !! parameter documentation file.
1252 subroutine log_version_cs(CS, modulename, version, desc, log_to_all, all_default, layout, debugging)
1253  type(param_file_type), intent(in) :: CS !< File parser type
1254  character(len=*), intent(in) :: modulename !< Name of calling module
1255  character(len=*), intent(in) :: version !< Version string of module
1256  character(len=*), optional, intent(in) :: desc !< Module description
1257  logical, optional, intent(in) :: log_to_all !< If present and true, log this parameter to the
1258  !! ..._doc.all files, even if this module also has layout
1259  !! or debugging parameters.
1260  logical, optional, intent(in) :: all_default !< If true, all parameters take their default values.
1261  logical, optional, intent(in) :: layout !< If present and true, this module has layout parameters.
1262  logical, optional, intent(in) :: debugging !< If present and true, this module has debugging parameters.
1263  ! Local variables
1264  character(len=240) :: mesg
1265 
1266  mesg = trim(modulename)//": "//trim(version)
1267  if (is_root_pe()) then
1268  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1269  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1270  endif
1271 
1272  if (present(desc)) call doc_module(cs%doc, modulename, desc, log_to_all, all_default, layout, debugging)
1273 
1274 end subroutine log_version_cs
1275 
1276 !> Log the version of a module to a log file and/or stdout.
1277 subroutine log_version_plain(modulename, version)
1278  character(len=*), intent(in) :: modulename !< Name of calling module
1279  character(len=*), intent(in) :: version !< Version string of module
1280  ! Local variables
1281  character(len=240) :: mesg
1282 
1283  mesg = trim(modulename)//": "//trim(version)
1284  if (is_root_pe()) then
1285  write(stdlog(),'(a)') trim(mesg)
1286  endif
1287 
1288 end subroutine log_version_plain
1289 
1290 !> Log the name and value of an integer model parameter in documentation files.
1291 subroutine log_param_int(CS, modulename, varname, value, desc, units, &
1292  default, layoutParam, debuggingParam, like_default)
1293  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1294  !! it is also a structure to parse for run-time parameters
1295  character(len=*), intent(in) :: modulename !< The name of the module using this parameter
1296  character(len=*), intent(in) :: varname !< The name of the parameter to log
1297  integer, intent(in) :: value !< The value of the parameter to log
1298  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1299  !! present, this parameter is not written to a doc file
1300  character(len=*), optional, intent(in) :: units !< The units of this parameter
1301  integer, optional, intent(in) :: default !< The default value of the parameter
1302  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1303  !! logged in the layout parameter file
1304  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1305  !! logged in the debugging parameter file
1306  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1307  !! though it has the default value, even if there is no default.
1308 
1309  character(len=240) :: mesg, myunits
1310 
1311  write(mesg, '(" ",a," ",a,": ",a)') trim(modulename), trim(varname), trim(left_int(value))
1312  if (is_root_pe()) then
1313  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1314  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1315  endif
1316 
1317  myunits=" "; if (present(units)) write(myunits(1:240),'(A)') trim(units)
1318  if (present(desc)) &
1319  call doc_param(cs%doc, varname, desc, myunits, value, default, &
1320  layoutparam=layoutparam, debuggingparam=debuggingparam, like_default=like_default)
1321 
1322 end subroutine log_param_int
1323 
1324 !> Log the name and values of an array of integer model parameter in documentation files.
1325 subroutine log_param_int_array(CS, modulename, varname, value, desc, &
1326  units, default, layoutParam, debuggingParam, like_default)
1327  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1328  !! it is also a structure to parse for run-time parameters
1329  character(len=*), intent(in) :: modulename !< The name of the module using this parameter
1330  character(len=*), intent(in) :: varname !< The name of the parameter to log
1331  integer, dimension(:), intent(in) :: value !< The value of the parameter to log
1332  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1333  !! present, this parameter is not written to a doc file
1334  character(len=*), optional, intent(in) :: units !< The units of this parameter
1335  integer, optional, intent(in) :: default !< The default value of the parameter
1336  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1337  !! logged in the layout parameter file
1338  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1339  !! logged in the debugging parameter file
1340  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1341  !! though it has the default value, even if there is no default.
1342 
1343  character(len=1320) :: mesg
1344  character(len=240) :: myunits
1345 
1346  write(mesg, '(" ",a," ",a,": ",A)') trim(modulename), trim(varname), trim(left_ints(value))
1347  if (is_root_pe()) then
1348  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1349  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1350  endif
1351 
1352  myunits=" "; if (present(units)) write(myunits(1:240),'(A)') trim(units)
1353  if (present(desc)) &
1354  call doc_param(cs%doc, varname, desc, myunits, value, default, &
1355  layoutparam=layoutparam, debuggingparam=debuggingparam, like_default=like_default)
1356 
1357 end subroutine log_param_int_array
1358 
1359 !> Log the name and value of a real model parameter in documentation files.
1360 subroutine log_param_real(CS, modulename, varname, value, desc, units, &
1361  default, debuggingParam, like_default)
1362  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1363  !! it is also a structure to parse for run-time parameters
1364  character(len=*), intent(in) :: modulename !< The name of the calling module
1365  character(len=*), intent(in) :: varname !< The name of the parameter to log
1366  real, intent(in) :: value !< The value of the parameter to log
1367  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1368  !! present, this parameter is not written to a doc file
1369  character(len=*), optional, intent(in) :: units !< The units of this parameter
1370  real, optional, intent(in) :: default !< The default value of the parameter
1371  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1372  !! logged in the debugging parameter file
1373  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1374  !! though it has the default value, even if there is no default.
1375 
1376  character(len=240) :: mesg, myunits
1377 
1378  write(mesg, '(" ",a," ",a,": ",a)') &
1379  trim(modulename), trim(varname), trim(left_real(value))
1380  if (is_root_pe()) then
1381  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1382  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1383  endif
1384 
1385  myunits="not defined"; if (present(units)) write(myunits(1:240),'(A)') trim(units)
1386  if (present(desc)) &
1387  call doc_param(cs%doc, varname, desc, myunits, value, default, &
1388  debuggingparam=debuggingparam, like_default=like_default)
1389 
1390 end subroutine log_param_real
1391 
1392 !> Log the name and values of an array of real model parameter in documentation files.
1393 subroutine log_param_real_array(CS, modulename, varname, value, desc, &
1394  units, default, debuggingParam, like_default)
1395  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1396  !! it is also a structure to parse for run-time parameters
1397  character(len=*), intent(in) :: modulename !< The name of the calling module
1398  character(len=*), intent(in) :: varname !< The name of the parameter to log
1399  real, dimension(:), intent(in) :: value !< The value of the parameter to log
1400  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1401  !! present, this parameter is not written to a doc file
1402  character(len=*), optional, intent(in) :: units !< The units of this parameter
1403  real, optional, intent(in) :: default !< The default value of the parameter
1404  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1405  !! logged in the debugging parameter file
1406  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1407  !! though it has the default value, even if there is no default.
1408 
1409  character(len=1320) :: mesg
1410  character(len=240) :: myunits
1411 
1412  !write(mesg, '(" ",a," ",a,": ",ES19.12,99(",",ES19.12))') &
1413  !write(mesg, '(" ",a," ",a,": ",G,99(",",G))') &
1414  ! trim(modulename), trim(varname), value
1415  write(mesg, '(" ",a," ",a,": ",a)') &
1416  trim(modulename), trim(varname), trim(left_reals(value))
1417  if (is_root_pe()) then
1418  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1419  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1420  endif
1421 
1422  myunits="not defined"; if (present(units)) write(myunits(1:240),'(A)') trim(units)
1423  if (present(desc)) &
1424  call doc_param(cs%doc, varname, desc, myunits, value, default, &
1425  debuggingparam=debuggingparam, like_default=like_default)
1426 
1427 end subroutine log_param_real_array
1428 
1429 !> Log the name and value of a logical model parameter in documentation files.
1430 subroutine log_param_logical(CS, modulename, varname, value, desc, &
1431  units, default, layoutParam, debuggingParam, like_default)
1432  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1433  !! it is also a structure to parse for run-time parameters
1434  character(len=*), intent(in) :: modulename !< The name of the calling module
1435  character(len=*), intent(in) :: varname !< The name of the parameter to log
1436  logical, intent(in) :: value !< The value of the parameter to log
1437  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1438  !! present, this parameter is not written to a doc file
1439  character(len=*), optional, intent(in) :: units !< The units of this parameter
1440  logical, optional, intent(in) :: default !< The default value of the parameter
1441  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1442  !! logged in the layout parameter file
1443  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1444  !! logged in the debugging parameter file
1445  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1446  !! though it has the default value, even if there is no default.
1447 
1448  character(len=240) :: mesg, myunits
1449 
1450  if (value) then
1451  write(mesg, '(" ",a," ",a,": True")') trim(modulename), trim(varname)
1452  else
1453  write(mesg, '(" ",a," ",a,": False")') trim(modulename), trim(varname)
1454  endif
1455  if (is_root_pe()) then
1456  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1457  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1458  endif
1459 
1460  myunits="Boolean"; if (present(units)) write(myunits(1:240),'(A)') trim(units)
1461  if (present(desc)) &
1462  call doc_param(cs%doc, varname, desc, myunits, value, default, &
1463  layoutparam=layoutparam, debuggingparam=debuggingparam, like_default=like_default)
1464 
1465 end subroutine log_param_logical
1466 
1467 !> Log the name and value of a character string model parameter in documentation files.
1468 subroutine log_param_char(CS, modulename, varname, value, desc, units, &
1469  default, layoutParam, debuggingParam, like_default)
1470  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1471  !! it is also a structure to parse for run-time parameters
1472  character(len=*), intent(in) :: modulename !< The name of the calling module
1473  character(len=*), intent(in) :: varname !< The name of the parameter to log
1474  character(len=*), intent(in) :: value !< The value of the parameter to log
1475  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1476  !! present, this parameter is not written to a doc file
1477  character(len=*), optional, intent(in) :: units !< The units of this parameter
1478  character(len=*), optional, intent(in) :: default !< The default value of the parameter
1479  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1480  !! logged in the layout parameter file
1481  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1482  !! logged in the debugging parameter file
1483  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1484  !! though it has the default value, even if there is no default.
1485 
1486  character(len=1024) :: mesg, myunits
1487 
1488  write(mesg, '(" ",a," ",a,": ",a)') &
1489  trim(modulename), trim(varname), trim(value)
1490  if (is_root_pe()) then
1491  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1492  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1493  endif
1494 
1495  myunits=" "; if (present(units)) write(myunits(1:1024),'(A)') trim(units)
1496  if (present(desc)) &
1497  call doc_param(cs%doc, varname, desc, myunits, value, default, &
1498  layoutparam=layoutparam, debuggingparam=debuggingparam, like_default=like_default)
1499 
1500 end subroutine log_param_char
1501 
1502 !> This subroutine writes the value of a time-type parameter to a log file,
1503 !! along with its name and the module it came from.
1504 subroutine log_param_time(CS, modulename, varname, value, desc, units, &
1505  default, timeunit, layoutParam, debuggingParam, log_date, like_default)
1506  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1507  !! it is also a structure to parse for run-time parameters
1508  character(len=*), intent(in) :: modulename !< The name of the calling module
1509  character(len=*), intent(in) :: varname !< The name of the parameter to log
1510  type(time_type), intent(in) :: value !< The value of the parameter to log
1511  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1512  !! present, this parameter is not written to a doc file
1513  character(len=*), optional, intent(in) :: units !< The units of this parameter
1514  type(time_type), optional, intent(in) :: default !< The default value of the parameter
1515  real, optional, intent(in) :: timeunit !< The number of seconds in a time unit for
1516  !! real-number output.
1517  logical, optional, intent(in) :: log_date !< If true, log the time_type in date format.
1518  !! If missing the default is false.
1519  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1520  !! logged in the layout parameter file
1521  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1522  !! logged in the debugging parameter file
1523  logical, optional, intent(in) :: like_default !< If present and true, log this parameter as
1524  !! though it has the default value, even if there is no default.
1525 
1526  ! Local variables
1527  real :: real_time, real_default
1528  logical :: use_timeunit, date_format
1529  character(len=240) :: mesg, myunits
1530  character(len=80) :: date_string, default_string
1531  integer :: days, secs, ticks, ticks_per_sec
1532 
1533  use_timeunit = .false.
1534  date_format = .false. ; if (present(log_date)) date_format = log_date
1535 
1536  call get_time(value, secs, days, ticks)
1537 
1538  if (ticks == 0) then
1539  write(mesg, '(" ",a," ",a," (Time): ",i0,":",i0)') trim(modulename), &
1540  trim(varname), days, secs
1541  else
1542  write(mesg, '(" ",a," ",a," (Time): ",i0,":",i0,":",i0)') trim(modulename), &
1543  trim(varname), days, secs, ticks
1544  endif
1545  if (is_root_pe()) then
1546  if (cs%log_open) write(cs%stdlog,'(a)') trim(mesg)
1547  if (cs%log_to_stdout) write(cs%stdout,'(a)') trim(mesg)
1548  endif
1549 
1550  if (present(desc)) then
1551  if (present(timeunit)) use_timeunit = (timeunit > 0.0)
1552  if (date_format) then
1553  myunits='[date]'
1554 
1555  date_string = convert_date_to_string(value)
1556  if (present(default)) then
1557  default_string = convert_date_to_string(default)
1558  call doc_param(cs%doc, varname, desc, myunits, date_string, &
1559  default=default_string, layoutparam=layoutparam, &
1560  debuggingparam=debuggingparam, like_default=like_default)
1561  else
1562  call doc_param(cs%doc, varname, desc, myunits, date_string, &
1563  layoutparam=layoutparam, debuggingparam=debuggingparam, like_default=like_default)
1564  endif
1565  elseif (use_timeunit) then
1566  if (present(units)) then
1567  write(myunits(1:240),'(A)') trim(units)
1568  else
1569  if (abs(timeunit-1.0) < 0.01) then ; myunits = "seconds"
1570  elseif (abs(timeunit-3600.0) < 1.0) then ; myunits = "hours"
1571  elseif (abs(timeunit-86400.0) < 1.0) then ; myunits = "days"
1572  elseif (abs(timeunit-3.1e7) < 1.0e6) then ; myunits = "years"
1573  else ; write(myunits,'(es8.2," sec")') timeunit ; endif
1574  endif
1575  real_time = (86400.0/timeunit)*days + secs/timeunit
1576  if (ticks > 0) real_time = real_time + &
1577  real(ticks) / (timeunit*get_ticks_per_second())
1578  if (present(default)) then
1579  call get_time(default, secs, days, ticks)
1580  real_default = (86400.0/timeunit)*days + secs/timeunit
1581  if (ticks > 0) real_default = real_default + &
1582  real(ticks) / (timeunit*get_ticks_per_second())
1583  call doc_param(cs%doc, varname, desc, myunits, real_time, real_default, like_default=like_default)
1584  else
1585  call doc_param(cs%doc, varname, desc, myunits, real_time, like_default=like_default)
1586  endif
1587  else
1588  call doc_param(cs%doc, varname, desc, value, default, units=units, like_default=like_default)
1589  endif
1590  endif
1591 
1592 end subroutine log_param_time
1593 
1594 !> This function converts a date into a string, valid with ticks and for dates up to year 99,999,999
1595 function convert_date_to_string(date) result(date_string)
1596  type(time_type), intent(in) :: date !< The date to be translated into a string.
1597  character(len=40) :: date_string !< A date string in a format like YYYY-MM-DD HH:MM:SS.sss
1598 
1599  ! Local variables
1600  character(len=40) :: sub_string
1601  real :: real_secs
1602  integer :: yrs, mons, days, hours, mins, secs, ticks, ticks_per_sec
1603 
1604  call get_date(date, yrs, mons, days, hours, mins, secs, ticks)
1605  write (date_string, '(i8.4)') yrs
1606  write (sub_string, '("-", i2.2, "-", I2.2, " ", i2.2, ":", i2.2, ":")') &
1607  mons, days, hours, mins
1608  date_string = trim(adjustl(date_string)) // trim(sub_string)
1609  if (ticks > 0) then
1610  ticks_per_sec = get_ticks_per_second()
1611  real_secs = secs + ticks/ticks_per_sec
1612  if (ticks_per_sec <= 100) then
1613  write (sub_string, '(F7.3)') real_secs
1614  else
1615  write (sub_string, '(F10.6)') real_secs
1616  endif
1617  else
1618  write (sub_string, '(i2.2)') secs
1619  endif
1620  date_string = trim(date_string) // trim(adjustl(sub_string))
1621 
1622 end function convert_date_to_string
1623 
1624 !> This subroutine reads the value of an integer model parameter from a parameter file
1625 !! and logs it in documentation files.
1626 subroutine get_param_int(CS, modulename, varname, value, desc, units, &
1627  default, fail_if_missing, do_not_read, do_not_log, &
1628  static_value, layoutParam, debuggingParam)
1629  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1630  !! it is also a structure to parse for run-time parameters
1631  character(len=*), intent(in) :: modulename !< The name of the calling module
1632  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1633  integer, intent(inout) :: value !< The value of the parameter that may be
1634  !! read from the parameter file and logged
1635  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1636  !! present, this parameter is not written to a doc file
1637  character(len=*), optional, intent(in) :: units !< The units of this parameter
1638  integer, optional, intent(in) :: default !< The default value of the parameter
1639  integer, optional, intent(in) :: static_value !< If this parameter is static, it takes
1640  !! this value, which can be compared for consistency with
1641  !! what is in the parameter file.
1642  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1643  !! if this variable is not found in the parameter file
1644  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1645  !! value for this parameter, although it might be logged.
1646  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1647  !! parameter to the documentation files
1648  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1649  !! logged in the layout parameter file
1650  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1651  !! logged in the debugging parameter file
1652 
1653  logical :: do_read, do_log
1654 
1655  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1656  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1657 
1658  if (do_read) then
1659  if (present(default)) value = default
1660  if (present(static_value)) value = static_value
1661  call read_param_int(cs, varname, value, fail_if_missing)
1662  endif
1663 
1664  if (do_log) then
1665  call log_param_int(cs, modulename, varname, value, desc, units, &
1666  default, layoutparam, debuggingparam)
1667  endif
1668 
1669 end subroutine get_param_int
1670 
1671 !> This subroutine reads the values of an array of integer model parameters from a parameter file
1672 !! and logs them in documentation files.
1673 subroutine get_param_int_array(CS, modulename, varname, value, desc, units, &
1674  default, fail_if_missing, do_not_read, do_not_log, &
1675  static_value, layoutParam, debuggingParam)
1676  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1677  !! it is also a structure to parse for run-time parameters
1678  character(len=*), intent(in) :: modulename !< The name of the calling module
1679  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1680  integer, dimension(:), intent(inout) :: value !< The value of the parameter that may be reset
1681  !! from the parameter file
1682  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1683  !! present, this parameter is not written to a doc file
1684  character(len=*), optional, intent(in) :: units !< The units of this parameter
1685  integer, optional, intent(in) :: default !< The default value of the parameter
1686  integer, optional, intent(in) :: static_value !< If this parameter is static, it takes
1687  !! this value, which can be compared for consistency with
1688  !! what is in the parameter file.
1689  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1690  !! if this variable is not found in the parameter file
1691  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1692  !! value for this parameter, although it might be logged.
1693  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1694  !! parameter to the documentation files
1695  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1696  !! logged in the layout parameter file
1697  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1698  !! logged in the debugging parameter file
1699 
1700  logical :: do_read, do_log
1701 
1702  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1703  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1704 
1705  if (do_read) then
1706  if (present(default)) then ; value(:) = default ; endif
1707  if (present(static_value)) then ; value(:) = static_value ; endif
1708  call read_param_int_array(cs, varname, value, fail_if_missing)
1709  endif
1710 
1711  if (do_log) then
1712  call log_param_int_array(cs, modulename, varname, value, desc, &
1713  units, default, layoutparam, debuggingparam)
1714  endif
1715 
1716 end subroutine get_param_int_array
1717 
1718 !> This subroutine reads the value of a real model parameter from a parameter file
1719 !! and logs it in documentation files.
1720 subroutine get_param_real(CS, modulename, varname, value, desc, units, &
1721  default, fail_if_missing, do_not_read, do_not_log, &
1722  static_value, debuggingParam, scale, unscaled)
1723  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1724  !! it is also a structure to parse for run-time parameters
1725  character(len=*), intent(in) :: modulename !< The name of the calling module
1726  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1727  real, intent(inout) :: value !< The value of the parameter that may be
1728  !! read from the parameter file and logged
1729  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1730  !! present, this parameter is not written to a doc file
1731  character(len=*), optional, intent(in) :: units !< The units of this parameter
1732  real, optional, intent(in) :: default !< The default value of the parameter
1733  real, optional, intent(in) :: static_value !< If this parameter is static, it takes
1734  !! this value, which can be compared for consistency with
1735  !! what is in the parameter file.
1736  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1737  !! if this variable is not found in the parameter file
1738  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1739  !! value for this parameter, although it might be logged.
1740  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1741  !! parameter to the documentation files
1742  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1743  !! logged in the debugging parameter file
1744  real, optional, intent(in) :: scale !< A scaling factor that the parameter is
1745  !! multiplied by before it is returned.
1746  real, optional, intent(out) :: unscaled !< The value of the parameter that would be
1747  !! returned without any multiplication by a scaling factor.
1748 
1749  logical :: do_read, do_log
1750 
1751  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1752  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1753 
1754  if (do_read) then
1755  if (present(default)) value = default
1756  if (present(static_value)) value = static_value
1757  call read_param_real(cs, varname, value, fail_if_missing)
1758  endif
1759 
1760  if (do_log) then
1761  call log_param_real(cs, modulename, varname, value, desc, units, &
1762  default, debuggingparam)
1763  endif
1764 
1765  if (present(unscaled)) unscaled = value
1766  if (present(scale)) value = scale*value
1767 
1768 end subroutine get_param_real
1769 
1770 !> This subroutine reads the values of an array of real model parameters from a parameter file
1771 !! and logs them in documentation files.
1772 subroutine get_param_real_array(CS, modulename, varname, value, desc, units, &
1773  default, fail_if_missing, do_not_read, do_not_log, debuggingParam, &
1774  static_value, scale, unscaled)
1775  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1776  !! it is also a structure to parse for run-time parameters
1777  character(len=*), intent(in) :: modulename !< The name of the calling module
1778  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1779  real, dimension(:), intent(inout) :: value !< The value of the parameter that may be
1780  !! read from the parameter file and logged
1781  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1782  !! present, this parameter is not written to a doc file
1783  character(len=*), optional, intent(in) :: units !< The units of this parameter
1784  real, optional, intent(in) :: default !< The default value of the parameter
1785  real, optional, intent(in) :: static_value !< If this parameter is static, it takes
1786  !! this value, which can be compared for consistency with
1787  !! what is in the parameter file.
1788  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1789  !! if this variable is not found in the parameter file
1790  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1791  !! value for this parameter, although it might be logged.
1792  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1793  !! parameter to the documentation files
1794  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1795  !! logged in the debugging parameter file
1796  real, optional, intent(in) :: scale !< A scaling factor that the parameter is
1797  !! multiplied by before it is returned.
1798  real, dimension(:), optional, intent(out) :: unscaled !< The value of the parameter that would be
1799  !! returned without any multiplication by a scaling factor.
1800 
1801  logical :: do_read, do_log
1802 
1803  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1804  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1805 
1806  if (do_read) then
1807  if (present(default)) then ; value(:) = default ; endif
1808  if (present(static_value)) then ; value(:) = static_value ; endif
1809  call read_param_real_array(cs, varname, value, fail_if_missing)
1810  endif
1811 
1812  if (do_log) then
1813  call log_param_real_array(cs, modulename, varname, value, desc, &
1814  units, default, debuggingparam)
1815  endif
1816 
1817  if (present(unscaled)) unscaled(:) = value(:)
1818  if (present(scale)) value(:) = scale*value(:)
1819 
1820 end subroutine get_param_real_array
1821 
1822 !> This subroutine reads the value of a character string model parameter from a parameter file
1823 !! and logs it in documentation files.
1824 subroutine get_param_char(CS, modulename, varname, value, desc, units, &
1825  default, fail_if_missing, do_not_read, do_not_log, &
1826  static_value, layoutParam, debuggingParam)
1827  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1828  !! it is also a structure to parse for run-time parameters
1829  character(len=*), intent(in) :: modulename !< The name of the calling module
1830  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1831  character(len=*), intent(inout) :: value !< The value of the parameter that may be
1832  !! read from the parameter file and logged
1833  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1834  !! present, this parameter is not written to a doc file
1835  character(len=*), optional, intent(in) :: units !< The units of this parameter
1836  character(len=*), optional, intent(in) :: default !< The default value of the parameter
1837  character(len=*), optional, intent(in) :: static_value !< If this parameter is static, it takes
1838  !! this value, which can be compared for consistency with
1839  !! what is in the parameter file.
1840  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1841  !! if this variable is not found in the parameter file
1842  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1843  !! value for this parameter, although it might be logged.
1844  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1845  !! parameter to the documentation files
1846  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1847  !! logged in the layout parameter file
1848  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1849  !! logged in the debugging parameter file
1850 
1851  logical :: do_read, do_log
1852 
1853  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1854  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1855 
1856  if (do_read) then
1857  if (present(default)) value = default
1858  if (present(static_value)) value = static_value
1859  call read_param_char(cs, varname, value, fail_if_missing)
1860  endif
1861 
1862  if (do_log) then
1863  call log_param_char(cs, modulename, varname, value, desc, units, &
1864  default, layoutparam, debuggingparam)
1865  endif
1866 
1867 end subroutine get_param_char
1868 
1869 !> This subroutine reads the values of an array of character string model parameters
1870 !! from a parameter file and logs them in documentation files.
1871 subroutine get_param_char_array(CS, modulename, varname, value, desc, units, &
1872  default, fail_if_missing, do_not_read, do_not_log, static_value)
1873  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1874  !! it is also a structure to parse for run-time parameters
1875  character(len=*), intent(in) :: modulename !< The name of the calling module
1876  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1877  character(len=*), dimension(:), intent(inout) :: value !< The value of the parameter that may be
1878  !! read from the parameter file and logged
1879  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1880  !! present, this parameter is not written to a doc file
1881  character(len=*), optional, intent(in) :: units !< The units of this parameter
1882  character(len=*), optional, intent(in) :: default !< The default value of the parameter
1883  character(len=*), optional, intent(in) :: static_value !< If this parameter is static, it takes
1884  !! this value, which can be compared for consistency with
1885  !! what is in the parameter file.
1886  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1887  !! if this variable is not found in the parameter file
1888  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1889  !! value for this parameter, although it might be logged.
1890  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1891  !! parameter to the documentation files
1892 
1893  ! Local variables
1894  logical :: do_read, do_log
1895  integer :: i, len_tot, len_val
1896  character(len=1024) :: cat_val
1897 
1898  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1899  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1900 
1901  if (do_read) then
1902  if (present(default)) then ; value(:) = default ; endif
1903  if (present(static_value)) then ; value(:) = static_value ; endif
1904  call read_param_char_array(cs, varname, value, fail_if_missing)
1905  endif
1906 
1907  if (do_log) then
1908  cat_val = trim(value(1)); len_tot = len_trim(value(1))
1909  do i=2,size(value)
1910  len_val = len_trim(value(i))
1911  if ((len_val > 0) .and. (len_tot + len_val + 2 < 240)) then
1912  cat_val = trim(cat_val)//achar(34)// ", "//achar(34)//trim(value(i))
1913  len_tot = len_tot + len_val
1914  endif
1915  enddo
1916  call log_param_char(cs, modulename, varname, cat_val, desc, &
1917  units, default)
1918  endif
1919 
1920 end subroutine get_param_char_array
1921 
1922 !> This subroutine reads the value of a logical model parameter from a parameter file
1923 !! and logs it in documentation files.
1924 subroutine get_param_logical(CS, modulename, varname, value, desc, units, &
1925  default, fail_if_missing, do_not_read, do_not_log, &
1926  static_value, layoutParam, debuggingParam)
1927  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1928  !! it is also a structure to parse for run-time parameters
1929  character(len=*), intent(in) :: modulename !< The name of the calling module
1930  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1931  logical, intent(inout) :: value !< The value of the parameter that may be
1932  !! read from the parameter file and logged
1933  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1934  !! present, this parameter is not written to a doc file
1935  character(len=*), optional, intent(in) :: units !< The units of this parameter
1936  logical, optional, intent(in) :: default !< The default value of the parameter
1937  logical, optional, intent(in) :: static_value !< If this parameter is static, it takes
1938  !! this value, which can be compared for consistency with
1939  !! what is in the parameter file.
1940  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1941  !! if this variable is not found in the parameter file
1942  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1943  !! value for this parameter, although it might be logged.
1944  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1945  !! parameter to the documentation files
1946  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1947  !! logged in the layout parameter file
1948  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1949  !! logged in the debugging parameter file
1950 
1951  logical :: do_read, do_log
1952 
1953  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
1954  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
1955 
1956  if (do_read) then
1957  if (present(default)) value = default
1958  if (present(static_value)) value = static_value
1959  call read_param_logical(cs, varname, value, fail_if_missing)
1960  endif
1961 
1962  if (do_log) then
1963  call log_param_logical(cs, modulename, varname, value, desc, &
1964  units, default, layoutparam, debuggingparam)
1965  endif
1966 
1967 end subroutine get_param_logical
1968 
1969 !> This subroutine reads the value of a time-type model parameter from a parameter file
1970 !! and logs it in documentation files.
1971 subroutine get_param_time(CS, modulename, varname, value, desc, units, &
1972  default, fail_if_missing, do_not_read, do_not_log, &
1973  timeunit, static_value, layoutParam, debuggingParam, &
1974  log_as_date)
1975  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
1976  !! it is also a structure to parse for run-time parameters
1977  character(len=*), intent(in) :: modulename !< The name of the calling module
1978  character(len=*), intent(in) :: varname !< The case-sensitive name of the parameter to read
1979  type(time_type), intent(inout) :: value !< The value of the parameter that may be
1980  !! read from the parameter file and logged
1981  character(len=*), optional, intent(in) :: desc !< A description of this variable; if not
1982  !! present, this parameter is not written to a doc file
1983  character(len=*), optional, intent(in) :: units !< The units of this parameter
1984  type(time_type), optional, intent(in) :: default !< The default value of the parameter
1985  type(time_type), optional, intent(in) :: static_value !< If this parameter is static, it takes
1986  !! this value, which can be compared for consistency with
1987  !! what is in the parameter file.
1988  logical, optional, intent(in) :: fail_if_missing !< If present and true, a fatal error occurs
1989  !! if this variable is not found in the parameter file
1990  logical, optional, intent(in) :: do_not_read !< If present and true, do not read a
1991  !! value for this parameter, although it might be logged.
1992  logical, optional, intent(in) :: do_not_log !< If present and true, do not log this
1993  !! parameter to the documentation files
1994  real, optional, intent(in) :: timeunit !< The number of seconds in a time unit for
1995  !! real-number input to be translated to a time.
1996  logical, optional, intent(in) :: layoutParam !< If present and true, this parameter is
1997  !! logged in the layout parameter file
1998  logical, optional, intent(in) :: debuggingParam !< If present and true, this parameter is
1999  !! logged in the debugging parameter file
2000  logical, optional, intent(in) :: log_as_date !< If true, log the time_type in date
2001  !! format. The default is false.
2002 
2003  logical :: do_read, do_log, date_format, log_date
2004 
2005  do_read = .true. ; if (present(do_not_read)) do_read = .not.do_not_read
2006  do_log = .true. ; if (present(do_not_log)) do_log = .not.do_not_log
2007  log_date = .false.
2008 
2009  if (do_read) then
2010  if (present(default)) value = default
2011  if (present(static_value)) value = static_value
2012  call read_param_time(cs, varname, value, timeunit, fail_if_missing, date_format=log_date)
2013  endif
2014 
2015  if (do_log) then
2016  if (present(log_as_date)) log_date = log_as_date
2017  call log_param_time(cs, modulename, varname, value, desc, units, default, &
2018  timeunit, layoutparam=layoutparam, &
2019  debuggingparam=debuggingparam, log_date=log_date)
2020  endif
2021 
2022 end subroutine get_param_time
2023 
2024 ! -----------------------------------------------------------------------------
2025 
2026 !> Resets the parameter block name to blank
2027 subroutine clearparameterblock(CS)
2028  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
2029  !! it is also a structure to parse for run-time parameters
2030 
2031  type(parameter_block), pointer :: block => null()
2032  if (associated(cs%blockName)) then
2033  block => cs%blockName
2034  block%name = ''
2035  else
2036  if (is_root_pe()) call mom_error(fatal, &
2037  'clearParameterBlock: A clear was attempted before allocation.')
2038  endif
2039 end subroutine clearparameterblock
2040 
2041 !> Tags blockName onto the end of the active parameter block name
2042 subroutine openparameterblock(CS,blockName,desc)
2043  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
2044  !! it is also a structure to parse for run-time parameters
2045  character(len=*), intent(in) :: blockName !< The name of a parameter block being added
2046  character(len=*), optional, intent(in) :: desc !< A description of the parameter block being added
2047 
2048  type(parameter_block), pointer :: block => null()
2049  if (associated(cs%blockName)) then
2050  block => cs%blockName
2051  block%name = pushblocklevel(block%name,blockname)
2052  call doc_openblock(cs%doc,block%name,desc)
2053  else
2054  if (is_root_pe()) call mom_error(fatal, &
2055  'openParameterBlock: A push was attempted before allocation.')
2056  endif
2057 end subroutine openparameterblock
2058 
2059 !> Remove the lowest level of recursion from the active block name
2060 subroutine closeparameterblock(CS)
2061  type(param_file_type), intent(in) :: CS !< The control structure for the file_parser module,
2062  !! it is also a structure to parse for run-time parameters
2063 
2064  type(parameter_block), pointer :: block => null()
2065 
2066  if (associated(cs%blockName)) then
2067  block => cs%blockName
2068  if (is_root_pe().and.len_trim(block%name)==0) call mom_error(fatal, &
2069  'closeParameterBlock: A pop was attempted on an empty stack. ("'//&
2070  trim(block%name)//'")')
2071  call doc_closeblock(cs%doc,block%name)
2072  else
2073  if (is_root_pe()) call mom_error(fatal, &
2074  'closeParameterBlock: A pop was attempted before allocation.')
2075  endif
2076  block%name = popblocklevel(block%name)
2077 end subroutine closeparameterblock
2078 
2079 !> Extends block name (deeper level of parameter block)
2080 function pushblocklevel(oldblockName,newBlockName)
2081  character(len=*), intent(in) :: oldBlockName !< A sequence of hierarchical parameter block names
2082  character(len=*), intent(in) :: newBlockName !< A new block name to add to the end of the sequence
2083  character(len=len(oldBlockName)+40) :: pushBlockLevel
2084 
2085  if (len_trim(oldblockname)>0) then
2086  pushblocklevel=trim(oldblockname)//'%'//trim(newblockname)
2087  else
2088  pushblocklevel=trim(newblockname)
2089  endif
2090 end function pushblocklevel
2091 
2092 !> Truncates block name (shallower level of parameter block)
2093 function popblocklevel(oldblockName)
2094  character(len=*), intent(in) :: oldBlockName !< A sequence of hierarchical parameter block names
2095  character(len=len(oldBlockName)+40) :: popBlockLevel
2096 
2097  integer :: i
2098  i = index(trim(oldblockname), '%', .true.)
2099  if (i>1) then
2100  popblocklevel = trim(oldblockname(1:i-1))
2101  elseif (i==0) then
2102  popblocklevel = ''
2103  else ! i==1
2104  if (is_root_pe()) call mom_error(fatal, &
2105  'popBlockLevel: A pop was attempted leaving an empty block name.')
2106  endif
2107 end function popblocklevel
2108 
2109 !> \namespace mom_file_parser
2110 !!
2111 !! By Robert Hallberg and Alistair Adcroft, updated 9/2013.
2112 !!
2113 !! The subroutines here parse a set of input files for the value
2114 !! a named parameter and sets that parameter at run time. Currently
2115 !! these files use use one of several formats:
2116 !! \#define VAR ! To set the logical VAR to true.
2117 !! VAR = True ! To set the logical VAR to true.
2118 !! \#undef VAR ! To set the logical VAR to false.
2119 !! VAR = False ! To set the logical VAR to false.
2120 !! \#define VAR 999 ! To set the real or integer VAR to 999.
2121 !! VAR = 999 ! To set the real or integer VAR to 999.
2122 !! \#override VAR = 888 ! To override a previously set value.
2123 !! VAR = 1.1, 2.2, 3.3 ! To set an array of real values.
2124  ! Note that in the comments above, dOxygen translates \# to # .
2125 !!
2126 !! In addition, when set by the get_param interface, the values of
2127 !! parameters are automatically logged, along with defaults, units,
2128 !! and a description. It is an error for a variable to be overridden
2129 !! more than once, and MOM6 has a facility to check for unused lines
2130 !! to set variables, which may indicate miss-spelled or archaic
2131 !! parameters. Parameter names are case-specific, and lines may use
2132 !! a F90 or C++ style comment, starting with ! or //.
2133 
2134 end module mom_file_parser
Wraps the FMS time manager functions.
The subroutines here provide hooks for document generation functions at various levels of granularity...
Definition: MOM_document.F90:3
A structure that can be parsed to read and document run-time parameters.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
Document parameter values.
Routines for error handling and I/O management.
Specify the active parameter block.
An overloaded interface to log version information about modules.
An overloaded interface to read various types of parameters.
A structure that controls where the documentation occurs, its veborsity and formatting.
Handy functions for manipulating strings.
The valid lines extracted from an input parameter file without comments.
An overloaded interface to read and log the values of various types of parameters.