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