FMS  2026.01
Flexible Modeling System
field_manager.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup field_manager_mod field_manager_mod
19 !> @ingroup field_manager
20 !> @brief Reads entries from a field table and stores this
21 !! information along with the type of field it belongs to.
22 !!
23 !> This allows the component models to query the field manager to see if non-default
24 !! methods of operation are desired. In essence the field table is a
25 !! powerful type of namelist. Default values can be provided for all the
26 !! fields through a namelist, individual fields can be modified through
27 !! the field table however.
28 !!
29 !> @author William Cooke
30 !!
31 !! An example of field table entries could be
32 !! <PRE>
33 !! "tracer","atmos_mod","sphum"
34 !!
35 !! "tracer","atmos_mod","sf6"
36 !! "longname","sulf_hex"
37 !! "advection_scheme_horiz","2nd_order"
38 !! "Profile_type","Fixed","surface_value = 0.0E+00"/
39 !!
40 !! "prog_tracers","ocean_mod","age_global"
41 !! horizontal-advection-scheme = mdfl_sweby
42 !! vertical-advection-scheme = mdfl_sweby
43 !! restart_file = ocean_age.res.nc
44 !! </PRE>
45 !!
46 !! The field table consists of entries in the following format.
47 !!
48 !! The first line of an entry should consist of three quoted strings.
49 !!
50 !! The first quoted string will tell the field manager what type of
51 !! field it is.
52 !!
53 !! The second quoted string will tell the field manager which model the
54 !! field is being applied to.
55 !! The supported types at present are
56 !!<PRE>
57 !! "coupler_mod" for the coupler,
58 !! "atmos_mod" for the atmosphere model,
59 !! "ocean_mod" for the ocean model,
60 !! "land_mod" for the land model, and,
61 !! "ice_mod" for the ice model.
62 !!</PRE>
63 !! The third quoted string should be a unique name that can be used as a
64 !! query.
65 !!
66 !! The second and following lines of each entry are called methods in
67 !! this context. Methods can be developed within any module and these
68 !! modules can query the field manager to find any methods that are
69 !! supplied in the field table.
70 !!
71 !! These lines can be coded quite flexibly.
72 !!
73 !! The line can consist of two or three quoted strings or a simple unquoted
74 !! string.
75 !!
76 !! If the line consists two or three quoted strings, then the first string will
77 !! be an identifier that the querying module will ask for.
78 !!
79 !! The second string will be a name that the querying module can use to
80 !! set up values for the module.
81 !!
82 !! The third string, if present, can supply parameters to the calling module that can be
83 !! parsed and used to further modify values.
84 !!
85 !! If the line consists of a simple unquoted string then quotes are not allowed
86 !! in any part of the line.
87 !!
88 !! An entry is ended with a backslash (/) as the final character in a
89 !! row.
90 !!
91 !! Comments can be inserted in the field table by having a # as the
92 !! first character in the line.
93 !!
94 !! In the example above we have three field entries.
95 !!
96 !! The first is a simple declaration of a tracer called "sphum".
97 !!
98 !! The second is for a tracer called "sf6". In this case a field named
99 !! "longname" will be given the value "sulf_hex". A field named
100 !! "advection_scheme_horiz" will be given the value "2nd_order". Finally a field
101 !! name "Profile_type" will be given a child field called "Fixed", and that field
102 !! will be given a field called "surface_value" with a real value of 0.0E+00.
103 !!
104 !! The third entry is an example of a oceanic age tracer. Note that the
105 !! method lines are formatted differently here. This is the flexibility mentioned
106 !! above.
107 !!
108 !! With these formats, a number of restrictions are required.
109 !!
110 !! The following formats are equally valid.
111 !!<PRE>
112 !! "longname","sulf_hex"
113 !! "longname = sulf_hex"
114 !! longname = sulf_hex
115 !!</PRE>
116 !! However the following is not valid.
117 !!<PRE>
118 !! longname = "sulf_hex"
119 !!</PRE>
120 !!
121 !! In the SF6 example above the last line of the entry could be written in the
122 !! following ways.
123 !!<PRE>
124 !! "Profile_type","Fixed","surface_value = 0.0E+00"/
125 !! Profile_type/Fixed/surface_value = 0.0E+00/
126 !!</PRE>
127 !!
128 !! Values supplied with fields are converted to the various types with the
129 !! following assumptions.
130 !!<PRE>
131 !! Real values : These values contain a decimal point or are in exponential format.
132 !! These values only support e or E format for exponentials.
133 !! e.g. 10.0, 1e10 and 1E10 are considered to be real numbers.
134 !!
135 !! Integer values : These values only contain numbers.
136 !! e.g 10 is an integer. 10.0 and 1e10 are not.
137 !!
138 !! Logical values : These values are supplied as one of the following formats.
139 !! T, .T., TRUE, .TRUE.
140 !! t, .t., true, .true.
141 !! F, .F., FALSE, .FALSE.
142 !! f, .f., false, .false.
143 !! These will be converted to T or F in a dump of the field.
144 !!
145 !! Character strings : These values are assumed to be strings if a character
146 !! other than an e (or E) is in the value. Numbers can be suppled in the value.
147 !! If the value does not meet the criteria for a real, integer or logical type,
148 !! it is assumed to be a character type.
149 !!</PRE>
150 !! The entries within the field table can be designed by the individual
151 !! authors of code to allow modification of their routines.
152 !!
153 
154 !> @addtogroup field_manager_mod
155 !> @{
156 module field_manager_mod
157 !TODO this variable can be removed when the legacy table is no longer used
158 #ifndef MAXFIELDS_
159 #define MAXFIELDS_ 250
160 #endif
161 
162 !TODO this variable can be removed when the legacy table is not longer used
163 #ifndef MAXFIELDMETHODS_
164 #define MAXFIELDMETHODS_ 250
165 #endif
166 
167 !
168 ! <CONTACT EMAIL="William.Cooke@noaa.gov"> William Cooke
169 ! </CONTACT>
170 !
171 ! <REVIEWER EMAIL="Richard.Slater@noaa.gov"> Richard D. Slater
172 ! </REVIEWER>
173 !
174 ! <REVIEWER EMAIL="Matthew.Harrison@noaa.gov"> Matthew Harrison
175 ! </REVIEWER>
176 !
177 ! <REVIEWER EMAIL="John.Dunne@noaa.gov"> John P. Dunne
178 ! </REVIEWER>
179 
180 use mpp_mod, only : mpp_error, &
181  fatal, &
182  note, &
183  warning, &
184  mpp_pe, &
185  mpp_root_pe, &
186  stdlog, &
187  stdout, &
188  input_nml_file
189 use fms_mod, only : lowercase, &
192 use fms2_io_mod, only: file_exists, get_instance_filename
193 use platform_mod, only: r4_kind, r8_kind, fms_path_len, fms_file_len
194 #ifdef use_yaml
195 use fm_yaml_mod
196 #endif
197 
198 implicit none
199 private
200 
201 #include<file_version.h>
202 logical :: module_is_initialized = .false.
203 
204 public :: field_manager_init !< (nfields, [table_name]) returns number of fields
205 public :: field_manager_end !< ()
206 public :: find_field_index !< (model, field_name) or (list_path)
207 public :: get_field_info !< (n,fld_type,fld_name,model,num_methods)
208  !! Returns parameters relating to field n.
209 public :: get_field_method !< (n, m, method) Returns the m-th method of field n
210 public :: get_field_methods !< (n, methods) Returns the methods related to field n
211 public :: parse !< (text, label, values) Overloaded function to parse integer,
212  !! real or character. Parse returns the number of values
213  !! decoded (> 1 => an array of values)
214 public :: fm_change_list !< (list) return success
215 public :: fm_change_root !< (list) return success
216 public :: fm_dump_list !< (list [, recursive]) return success
217 public :: fm_exists !< (field) return success
218 public :: fm_get_index !< (field) return index
219 public :: fm_get_current_list !< () return path
220 public :: fm_get_length !< (list) return length
221 public :: fm_get_type !< (field) return string
222 public :: fm_get_value !< (entry, value [, index]) return success !! generic
223 public :: fm_init_loop !< (list, iter)
224 public :: fm_loop_over_list !< (list, name, type, index) return success
225  !! (iter, name, type, index) return success
226 public :: fm_new_list !< (list [, create] [, keep]) return index
227 public :: fm_new_value !< (entry, value [, create] [, index]) return index !! generic
228 public :: fm_reset_loop !< ()
229 public :: fm_return_root !< () return success
230 public :: fm_modify_name !< (oldname, newname) return success
231 public :: fm_query_method !< (name, method_name, method_control) return success and
232  !! name and control strings
233 public :: fm_find_methods !< (list, methods, control) return success and name and
234  !! control strings.
235 public :: fm_copy_list !< (list, suffix, [create]) return index
236 private :: create_field ! (list_p, name) return field pointer
237 private :: dump_list ! (list_p, recursive, depth) return success
238 private :: find_base ! (field, path, base)
239 private :: find_field ! (field, list_p) return field pointer
240 private :: find_head ! (field, head, rest)
241 private :: find_list ! (list, list_p, create) return field pointer
242 private :: get_field ! (field, list_p) return field pointer
243 private :: initialize_module_variables ! ()
244 private :: make_list ! (list_p, name) return field pointer
245 
246 !> The length of a character string representing the field name.
247 integer, parameter, public :: fm_field_name_len = 48
248 !> The length of a character string representing character values for the field.
249 integer, parameter, public :: fm_string_len = 1024
250 !> The length of a character string representing the various types that the values of the field can take.
251 integer, parameter, public :: fm_type_name_len = 8
252 !> Number of models (ATMOS, OCEAN, LAND, ICE, COUPLER).
253 integer, parameter, public :: num_models = 5
254 !> The value returned if a field is not defined.
255 integer, parameter, public :: no_field = -1
256 !> Atmospheric model.
257 integer, parameter, public :: model_atmos = 1
258 !> Ocean model.
259 integer, parameter, public :: model_ocean = 2
260 !> Land model.
261 integer, parameter, public :: model_land = 3
262 !> Ice model.
263 integer, parameter, public :: model_ice = 4
264 !> Ice model.
265 integer, parameter, public :: model_coupler = 5
266 !> Model names, e.g. MODEL_NAMES(MODEL_OCEAN) is 'oceanic'
267 character(len=11), parameter, public, dimension(NUM_MODELS) :: &
268  model_names=(/'atmospheric','oceanic ','land ','ice ','coupler '/)
269 
270 !> @}
271 
272 !> @brief This method_type is a way to allow a component module to alter the parameters it needs
273 !! for various tracers.
274 !!
275 !> In essence this is a way to modify a namelist. A namelist can supply
276 !! default parameters for all tracers. This method will allow the user to modify these
277 !! default parameters for an individual tracer. An example could be that the user wishes to
278 !! use second order advection on a tracer and also use fourth order advection on a second
279 !! tracer within the same model run. The default advection could be second order and the
280 !! field table would then indicate that the second tracer requires fourth order advection.
281 !! This would be parsed by the advection routine.
282 !> @ingroup field_manager_mod
283 type, public :: method_type
284 
285  character(len=fm_string_len) :: method_type !< This string represents a tag that a module
286  !! using this method can key on. Typically this should
287  !! contain some reference to the module that is calling it.
288  character(len=fm_string_len) :: method_name !< This is the name of a method which the module
289  !! can parse and use to assign different default values to
290  !! a field method.
291  character(len=fm_string_len) :: method_control !< This is the string containing parameters that
292  !! the module can use as values for a field method. These should
293  !! override default values within the module.
294 end type
295 
296 !> This method_type is the same as method_type except that the
297 !! method_control string is not present. This is used when you wish to
298 !! change to a scheme within a module but do not need to pass
299 !! parameters. See @ref method_type for member information.
300 !> @ingroup field_manager_mod
301 type, public :: method_type_short
302  character(len=fm_string_len) :: method_type
303  character(len=fm_string_len) :: method_name
304 end type
305 
306 !> This is the same as method_type except that the
307 !! method_control and method_name strings are not present. This is used
308 !! when you wish to change to a scheme within a module but do not need
309 !! to pass parameters.
310 !> @ingroup field_manager_mod
311 type, public :: method_type_very_short
312  character(len=fm_string_len) :: method_type
313 end type
314 
315 !> Iterator over the field manager list
316 !> @ingroup field_manager_mod
317 type, public :: fm_list_iter_type
318  type(field_def), pointer :: ptr => null() !< pointer to the current field
319 end type fm_list_iter_type
320 
321 !> @ingroup field_manager_mod
322 type(method_type), public :: default_method
323 
324 !> @brief Returns an index corresponding to the given field name.
325 !!
326 !> Model number can be given for old method.
327 !! <br>Example usage:
328 !! @code{.F90}
329 !! value=find_field_index( model, field_name )
330 !! value=find_field_index( field_name )
331 !! @endcode
332 !> @ingroup field_manager_mod
334  module procedure find_field_index_old
335  module procedure find_field_index_new
336 end interface
337 
338 !> @brief A function to parse an integer or an array of integers,
339 !! a real or an array of reals, a string or an array of strings.
340 !!
341 !> Parse is an integer function that decodes values from a text string.
342 !! The text string has the form: "label=list" where "label" is an
343 !! arbitrary user defined label describing the values being decoded,
344 !! and "list" is a list of one or more values separated by commas.
345 !! The values may be integer, real, or character.
346 !! Parse returns the number of values decoded.
347 !! <br>Example usage:
348 !! @code{.F90}
349 !! number = parse(text, label, value)
350 !! @endcode
351 !> @ingroup field_manager_mod
352 interface parse
353  module procedure parse_real_r4
354  module procedure parse_real_r8
355  module procedure parse_reals_r4
356  module procedure parse_reals_r8
357  module procedure parse_integer
358  module procedure parse_integers
359  module procedure parse_string
360  module procedure parse_strings
361 end interface
362 
363 !> @brief An overloaded function to assign a value to a field.
364 !!
365 !> Allocate and initialize a new value and return the index.
366 !! If an error condition occurs the parameter NO_FIELD is returned.
367 !!
368 !! If the type of the field is changing (e.g. real values being transformed to
369 !! integers), then any previous values for the field are removed and replaced
370 !! by the value passed in the present call to this function.
371 !!
372 !! If append is present and .true., then index cannot be greater than 0 if
373 !! it is present.
374 !! <br> Example usage:
375 !! @code{.F90}
376 !! field_index= fm_new_value(name, value, [create], [index], [append])
377 !! @endcode
378 !> @ingroup field_manager_mod
379 interface fm_new_value
380  module procedure fm_new_value_integer
381  module procedure fm_new_value_logical
382  module procedure fm_new_value_real_r4
383  module procedure fm_new_value_real_r8
384  module procedure fm_new_value_string
385 end interface
386 
387 !> @brief An overloaded function to find and extract a value for a named field.
388 !!
389 !> Find and extract the value for name. The value may be of type real,
390 !! integer, logical or character. If a single value from an array of values
391 !! is required, an optional index can be supplied.
392 !! Return true for success and false for failure
393 !! <br> Example usage:
394 !! @code{.F90}
395 !! success = fm_get_value(name, value, index)
396 !! @endcode
397 !> @ingroup field_manager_mod
398 interface fm_get_value
399  module procedure fm_get_value_integer
400  module procedure fm_get_value_logical
401  module procedure fm_get_value_real_r4
402  module procedure fm_get_value_real_r8
403  module procedure fm_get_value_string
404 end interface
405 
406 !> @brief A function for looping over a list.
407 !!
408 !> Loop over the list, setting the name, type and index
409 !! of the next field. Return false at the end of the loop.
410 !! <br> Example usage:
411 !! @code{.F90}
412 !! success = fm_loop_over_list(list, name, field_type, index)
413 !! @endcode
414 !> @ingroup field_manager_mod
416  module procedure fm_loop_over_list_new
417  module procedure fm_loop_over_list_old
418 end interface
419 
420 character(len=17), parameter :: module_name = 'field_manager_mod'
421 character(len=33), parameter :: error_header = '==>Error from '//trim(module_name)//': '
422 character(len=35), parameter :: warn_header = '==>Warning from '//trim(module_name)//': '
423 character(len=32), parameter :: note_header = '==>Note from '//trim(module_name)//': '
424 character(len=1), parameter :: comma = ","
425 character(len=1), parameter :: list_sep = '/'
426 !TODO these variable can be removed when the legacy table is no longer used
427 character(len=1), parameter :: comment = '#'
428 character(len=1), parameter :: dquote = '"'
429 character(len=1), parameter :: equal = '='
430 character(len=1), parameter :: squote = "'"
431 !
432 integer, parameter :: null_type = 0
433 integer, parameter :: integer_type = 1
434 integer, parameter :: list_type = 2
435 integer, parameter :: logical_type = 3
436 integer, parameter :: real_type = 4
437 integer, parameter :: string_type = 5
438 integer, parameter :: num_types = 5
439 integer, parameter :: array_increment = 10
440 !TODO these variable can be removed when the legacy table is no longer used
441 integer, parameter :: MAX_FIELDS = maxfields_
442 integer, parameter :: MAX_FIELD_METHODS = maxfieldmethods_
443 !
444 
445 !> @brief Private type for internal use
446 !> @ingroup field_manager_mod
447 type, private :: field_mgr_type
448  character(len=fm_field_name_len) :: field_type
449  character(len=fm_string_len) :: field_name
450  integer :: model, num_methods
451  type(method_type), dimension(:), allocatable :: methods !< methods associated with this field name
452 end type field_mgr_type
453 
454 !TODO These two types: field_names_type and field_names_type_short
455 !! will no longer be needed when the legacy field table is not used
456 !> @brief Private type for internal use
457 !> @ingroup field_manager_mod
458 type, private :: field_names_type
459  character(len=fm_field_name_len) :: fld_type
460  character(len=fm_field_name_len) :: mod_name
461  character(len=fm_string_len) :: fld_name
462 end type field_names_type
463 
464 !> @brief Private type for internal use
465 !> @ingroup field_manager_mod
466 type, private :: field_names_type_short
467  character(len=fm_field_name_len) :: fld_type
468  character(len=fm_field_name_len) :: mod_name
469 end type field_names_type_short
470 
471 !> @brief Private type for internal use
472 !> @ingroup field_manager_mod
473 type, private :: field_def
474  character (len=fm_field_name_len) :: name
475  integer :: index
476  type (field_def), pointer :: parent => null()
477  integer :: field_type
478  integer :: length
479  integer :: array_dim
480  integer :: max_index
481  type (field_def), pointer :: first_field => null()
482  type (field_def), pointer :: last_field => null()
483  integer, allocatable, dimension(:) :: i_value
484  logical, allocatable, dimension(:) :: l_value
485  real(r8_kind), allocatable, dimension(:) :: r_value !< string to real conversion will be done at r8;
486  !! all real values will be stored as r8_kind.
487  character(len=fm_string_len), allocatable, dimension(:) :: s_value
488  type (field_def), pointer :: next => null()
489  type (field_def), pointer :: prev => null()
490 end type field_def
491 
492 !> @addtogroup field_manager_mod
493 !> @{
494 
495 type(field_mgr_type), dimension(:), allocatable, private :: fields !< fields of field_mgr_type
496 
497 character(len=FMS_PATH_LEN) :: loop_list
498 character(len=fm_type_name_len) :: field_type_name(num_types)
499 character(len=fm_field_name_len) :: save_root_name
500 ! The string set is the set of characters.
501 character(len=52) :: set = "ABCDEFGHIJKLMNOPQRSTUVWXYZabcdefghijklmnopqrstuvwxyz"
502 ! If a character in the string being parsed matches a character within
503 ! the string set_nonexp then the string being parsed cannot be a number.
504 character(len=50) :: set_nonexp = "ABCDFGHIJKLMNOPQRSTUVWXYZabcdfghijklmnopqrstuvwxyz"
505 ! If a character in the string being parsed matches a character within
506 ! the string setnum then the string may be a number.
507 character(len=13) :: setnum = "0123456789+-."
508 integer :: num_fields = 0
509 type (field_def), pointer :: loop_list_p => null()
510 type (field_def), pointer :: current_list_p => null()
511 type (field_def), pointer :: root_p => null()
512 type (field_def), pointer :: save_root_parent_p => null()
513 type (field_def), target, save :: root
514 
515 logical :: use_field_table_yaml = .false. !< .True. if using the field_table.yaml,
516  !! .false. if using the legacy field_table
517 
518 namelist /field_manager_nml/ use_field_table_yaml
519 
520 contains
521 
522 !> @brief Routine to initialize the field manager.
523 !!
524 !> This routine reads from a file containing yaml paramaters.
525 !! These yaml parameters contain information on which schemes are
526 !! needed within various modules. The field manager does not
527 !! initialize any of those schemes however. It simply holds the
528 !! information and is queried by the appropriate module.
529 !!
530 !! The routine has two loops. The first loop initializes the my_table object
531 !! and counts the number of fields contained therein. The second loop is the
532 !! main loop that acts on each field in the my_table object, defining a list
533 !! object (in the field_manager definition) from which various fm routines may be
534 !! called, as well as populating the "fields" object and the "methods" objects
535 !! within each field object. The "fields" and "methods" objects are then used
536 !! with the subroutine new_name to append various characteristics to the list
537 !! object. Note that the "fields" and "methods" objects are also used with other
538 !! fm routines in a bit of a parallel system.
539 subroutine field_manager_init(nfields, table_name)
540 integer, intent(out), optional :: nfields !< number of fields
541 character(len=fm_string_len), intent(in), optional :: table_name !< Name of the field table, default
542 
543 if (module_is_initialized) then
544  if(present(nfields)) nfields = num_fields
545  return
546 endif
547 
549 
550 !TODO the use_field_table_yaml namelist can be removed when the legacy table is no longer in used
551 if (use_field_table_yaml) then
552  !Crash if you are not compiling with -Duse_yaml or if the field_table is present
553 #ifndef use_yaml
554  call mpp_error(fatal, "You cannot have use_field_table_yaml=.true. without compiling with -Duse_yaml")
555 #else
556  if (file_exists("field_table")) &
557  call mpp_error(fatal, "You cannot have the legacy field_table if use_field_table_yaml=.true.")
558 
559  call mpp_error(note, "field_manager_init:: You are using the yaml version of the field_table")
560  call mpp_error(note, &
561  &"field_manager_init:: You are using the yaml version of the field_table. &
562 The legacy field_table format will be deprecated in a future release, &
563 please switch to the yaml format.")
564  call read_field_table_yaml(nfields, table_name)
565 #endif
566 else
567  if (file_exists("field_table.yaml")) &
568  call mpp_error(fatal, "You cannot have the yaml field_table if use_field_table_yaml=.false.")
569  call mpp_error(note, "field_manager_init:: You are using the legacy version of the field_table")
570  call mpp_error(note, &
571  &"field_manager_init:: You are using the yaml version of the field_table. &
572 The legacy field_table format will be deprecated in a future release, &
573 please switch to the yaml format.")
574  call read_field_table_legacy(nfields, table_name)
575 endif
576 
577 end subroutine field_manager_init
578 
579 #ifdef use_yaml
580 
581 !> @brief Routine to read and parse the field table yaml
582 subroutine read_field_table_yaml(nfields, table_name)
583 integer, intent(out), optional :: nfields !< number of fields
584 character(len=*), intent(in), optional :: table_name !< Name of the field table file, default is 'field_table.yaml'
585 
586 character(len=FMS_FILE_LEN) :: tbl_name !< field_table yaml file
587 character(len=fm_string_len) :: method_control !< field_table yaml file
588 integer :: h, i, j, k, l, m !< dummy integer buffer
589 type (fmTable_t) :: my_table !< the field table
590 integer :: model !< model assocaited with the current field
591 character(len=FMS_PATH_LEN) :: list_name !< field_manager list name
592 character(len=fm_string_len) :: subparamvalue !< subparam value to be used when defining new name
593 character(len=fm_string_len) :: fm_yaml_null !< useful hack when OG subparam does not contain an equals sign
594 integer :: current_field !< field index within loop
595 integer :: index_list_name !< integer used as check for "no field"
596 integer :: subparamindex !< index to identify whether subparams exist for this field
597 logical :: fm_success !< logical for whether fm_change_list was a success
598 logical :: subparams !< logical whether subparams exist in this iteration
599 
600 character(len=FMS_FILE_LEN) :: filename !< Name of the expected field_table.yaml
601 
602 if (.not.PRESENT(table_name)) then
603  tbl_name = 'field_table.yaml'
604 else
605  tbl_name = trim(table_name)
606 endif
607 
608 call get_instance_filename(tbl_name, filename)
609 if (index(trim(filename), "ens_") .ne. 0) then
610  if (file_exists(filename) .and. file_exists(tbl_name)) &
611  call mpp_error(fatal, "Both "//trim(tbl_name)//" and "//trim(filename)//" exists, pick one!")
612 
613  !< If the end_* file does not exist, revert back to tbl_name
614  !! where every ensemble is using the same yaml
615  if (.not. file_exists(filename)) filename = tbl_name
616 endif
617 
618 if (.not. file_exists(trim(filename))) then
619  if(present(nfields)) nfields = 0
620  return
621 endif
622 
623 ! Construct my_table object
624 call build_fmtable(my_table, trim(filename))
625 
626 do h=1,size(my_table%types)
627  do i=1,size(my_table%types(h)%models)
628  do j=1,size(my_table%types(h)%models(i)%variables)
629  num_fields = num_fields + 1
630  end do
631  end do
632 end do
633 
634 allocate(fields(num_fields))
635 
636 current_field = 0
637 do h=1,size(my_table%types)
638  do i=1,size(my_table%types(h)%models)
639  select case (my_table%types(h)%models(i)%name)
640  case ('coupler_mod')
641  model = model_coupler
642  case ('atmos_mod')
643  model = model_atmos
644  case ('ocean_mod')
645  model = model_ocean
646  case ('land_mod')
647  model = model_land
648  case ('ice_mod')
649  model = model_ice
650  case default
651  call mpp_error(fatal, trim(error_header)//'The model name is unrecognised : &
652  &'//trim(my_table%types(h)%models(i)%name))
653  end select
654  do j=1,size(my_table%types(h)%models(i)%variables)
655  current_field = current_field + 1
656  list_name = list_sep//lowercase(trim(my_table%types(h)%models(i)%name))//list_sep//&
657  lowercase(trim(my_table%types(h)%name))//list_sep//&
658  lowercase(trim(my_table%types(h)%models(i)%variables(j)%name))
659  index_list_name = fm_new_list(list_name, create = .true.)
660  if ( index_list_name == no_field ) &
661  call mpp_error(fatal, trim(error_header)//'Could not set field list for '//trim(list_name))
662  fm_success = fm_change_list(list_name)
663  fields(current_field)%model = model
664  fields(current_field)%field_name = lowercase(trim(my_table%types(h)%models(i)%variables(j)%name))
665  fields(current_field)%field_type = lowercase(trim(my_table%types(h)%name))
666  fields(current_field)%num_methods = size(my_table%types(h)%models(i)%variables(j)%keys)
667  allocate(fields(current_field)%methods(fields(current_field)%num_methods))
668  if(fields(current_field)%num_methods.gt.0) then
669  subparams = (size(my_table%types(h)%models(i)%variables(j)%attributes) .gt. 0)
670  do k=1,size(my_table%types(h)%models(i)%variables(j)%keys)
671  fields(current_field)%methods(k)%method_type = &
672  lowercase(trim(my_table%types(h)%models(i)%variables(j)%keys(k)))
673  fields(current_field)%methods(k)%method_name = &
674  lowercase(trim(my_table%types(h)%models(i)%variables(j)%values(k)))
675  if (.not.subparams) then
676  call new_name(list_name, my_table%types(h)%models(i)%variables(j)%keys(k),&
677  my_table%types(h)%models(i)%variables(j)%values(k) )
678  else
679  subparamindex=-1
680  do l=1,size(my_table%types(h)%models(i)%variables(j)%attributes)
681  if(lowercase(trim(my_table%types(h)%models(i)%variables(j)%attributes(l)%paramname)).eq.&
682  lowercase(trim(fields(current_field)%methods(k)%method_type))) then
683  subparamindex = l
684  exit
685  end if
686  end do
687  if (subparamindex.eq.-1) then
688  call new_name(list_name, my_table%types(h)%models(i)%variables(j)%keys(k),&
689  my_table%types(h)%models(i)%variables(j)%values(k) )
690  else
691  do m=1,size(my_table%types(h)%models(i)%variables(j)%attributes(subparamindex)%keys)
692  method_control = " "
693  subparamvalue = " "
694  if (trim(my_table%types(h)%models(i)%variables(j)%values(k)).eq.'fm_yaml_null') then
695  fm_yaml_null = ''
696  else
697  fm_yaml_null = trim(my_table%types(h)%models(i)%variables(j)%values(k))//'/'
698  end if
699  method_control = trim(my_table%types(h)%models(i)%variables(j)%keys(k))//"/"//&
700  &trim(fm_yaml_null)//&
701  &trim(my_table%types(h)%models(i)%variables(j)%attributes(subparamindex)%keys(m))
702  subparamvalue = trim(my_table%types(h)%models(i)%variables(j)%attributes(subparamindex)%values(m))
703  call new_name(list_name, method_control, subparamvalue)
704  end do
705  end if
706  end if
707  end do
708  end if
709  end do
710  end do
711 end do
712 
713 if (present(nfields)) nfields = num_fields
714 end subroutine read_field_table_yaml
715 
716 !> @brief Subroutine to add new values to list parameters.
717 !!
718 !> This subroutine uses input strings list_name, method_name
719 !! and val_name_in to add new values to the list. Given
720 !! list_name a new list item is created that is named
721 !! method_name and is given the value or values in
722 !! val_name_in. If there is more than 1 value in
723 !! val_name_in, these values should be comma-separated.
724 subroutine new_name_yaml ( list_name, method_name_in , val_name_in)
725 character(len=*), intent(in) :: list_name !< The name of the field that is of interest here.
726 character(len=*), intent(in) :: method_name_in !< The name of the method that values are
727  !! being supplied for.
728 character(len=*), intent(inout) :: val_name_in !< The value or values that will be parsed and
729  !! used as the value when creating a new field or fields.
730 
731 character(len=fm_string_len) :: method_name !< name of method to be attached to new list
732 character(len=fm_string_len) :: val_name !< value name (to be converted to appropriate type)
733 integer, dimension(:), allocatable :: end_val !< end values in comma separated list
734 integer, dimension(:), allocatable :: start_val !< start values in comma separated list
735 integer :: i !< loop index
736 integer :: index_t !< appending index
737 integer :: num_elem !< number of elements in comma list
738 integer :: val_int !< value when converted to integer
739 integer :: val_type !< value type represented as integer for use in select case
740 logical :: append_new !< whether or not to append to existing list structure
741 logical :: val_logic !< value when converted to logical
742 real(r8_kind) :: val_real !< value when converted to real.
743  !! All strings will be converted to r8_kind reals.
744 
745 call strip_front_blanks(val_name_in)
746 method_name = trim(method_name_in)
747 call strip_front_blanks(method_name)
748 
749 index_t = 1
750 num_elem = 1
751 append_new = .false.
752 
753 ! If the array of values being passed in is a comma delimited list then count
754 ! the number of elements.
755 
756 do i = 1, len_trim(val_name_in)
757  if ( val_name_in(i:i) == comma ) then
758  num_elem = num_elem + 1
759  endif
760 enddo
761 
762 allocate(start_val(num_elem))
763 allocate(end_val(num_elem))
764 start_val(1) = 1
765 end_val(:) = len_trim(val_name_in)
766 
767 num_elem = 1
768 do i = 1, len_trim(val_name_in)
769  if ( val_name_in(i:i) == comma ) then
770  end_val(num_elem) = i-1
771  start_val(num_elem+1) = i+1
772  num_elem = num_elem + 1
773  endif
774 enddo
775 
776 do i = 1, num_elem
777 
778  if ( i .gt. 1 .or. index_t .eq. 0 ) then
779  append_new = .true.
780  index_t = 0 ! If append is true then index must be <= 0
781  endif
782  val_type = string_type ! Assume it is a string
783  val_name = val_name_in(start_val(i):end_val(i))
784  call strip_front_blanks(val_name)
785 
786  if ( scan(val_name(1:1), setnum ) > 0 ) then
787  if ( scan(val_name, set_nonexp ) .le. 0 ) then
788  if ( scan(val_name, '.') > 0 .or. scan(val_name, 'e') > 0 .or. scan(val_name, 'E') > 0) then
789  read(val_name, *) val_real
790  val_type = real_type
791  else
792  read(val_name, *) val_int
793  val_type = integer_type
794  endif
795  endif
796  endif
797 
798  if ( len_trim(val_name) == 1 .or. len_trim(val_name) == 3) then
799  if ( val_name == 't' .or. val_name == 'T' .or. val_name == '.t.' .or. val_name == '.T.' ) then
800  val_logic = .true.
801  val_type = logical_type
802  endif
803  if ( val_name == 'f' .or. val_name == 'F' .or. val_name == '.f.' .or. val_name == '.F.' ) then
804  val_logic = .false.
805  val_type = logical_type
806  endif
807  endif
808  if ( trim(lowercase(val_name)) == 'true' .or. trim(lowercase(val_name)) == '.true.' ) then
809  val_logic = .true.
810  val_type = logical_type
811  endif
812  if ( trim(lowercase(val_name)) == 'false' .or. trim(lowercase(val_name)) == '.false.' ) then
813  val_logic = .false.
814  val_type = logical_type
815  endif
816 
817  select case(val_type)
818 
819  case (integer_type)
820  if ( fm_new_value( method_name, val_int, create = .true., index = index_t, append = append_new ) < 0 ) &
821  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
822  ' (I) for '//trim(list_name))
823 
824  case (logical_type)
825  if ( fm_new_value( method_name, val_logic, create = .true., index = index_t, append = append_new) < 0 ) &
826  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
827  ' (L) for '//trim(list_name))
828 
829  case (real_type)
830  if ( fm_new_value( method_name, val_real, create = .true., index = index_t, append = append_new) < 0 ) &
831  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
832  ' (R) for '//trim(list_name))
833 
834  case (string_type)
835  if ( fm_new_value( method_name, val_name, create = .true., index = index_t, append = append_new) < 0 ) &
836  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
837  ' (S) for '//trim(list_name))
838  case default
839  call mpp_error(fatal, trim(error_header)//'Could not find a valid type to set the '//trim(method_name)//&
840  ' for '//trim(list_name))
841 
842  end select
843 
844 enddo
845  deallocate(start_val)
846  deallocate(end_val)
847 
848 end subroutine new_name_yaml
849 #endif
850 
851 !> @brief Routine to read and parse the field table yaml
852 !!
853 !> This routine reads from a file containing formatted strings.
854 !! These formatted strings contain information on which schemes are
855 !! needed within various modules. The field manager does not
856 !! initialize any of those schemes however. It simply holds the
857 !! information and is queried by the appropriate module.
858 subroutine read_field_table_legacy(nfields, table_name)
859 
860 integer, intent(out), optional :: nfields !< number of fields
861 character(len=fm_string_len), intent(in), optional :: table_name !< Name of the field table, default
862  !! is 'field_table'
863 
864 character(len=1024) :: record
865 character(len=fm_string_len) :: control_str
866 character(len=FMS_PATH_LEN) :: list_name
867 character(len=fm_string_len) :: method_name
868 character(len=fm_string_len) :: name_str
869 character(len=fm_string_len) :: type_str
870 character(len=fm_string_len) :: val_name
871 character(len=fm_string_len) :: tbl_name
872 integer :: control_array(MAX_FIELDS,3)
873 integer :: endcont
874 integer :: icount
875 integer :: index_list_name
876 integer :: iunit
877 integer :: l
878 integer :: log_unit
879 integer :: ltrec
880 integer :: m
881 integer :: midcont
882 integer :: model
883 integer :: startcont
884 integer :: io_status
885 logical :: flag_method
886 logical :: fm_success
887 type(field_names_type_short) :: text_names_short
888 type(field_names_type) :: text_names
889 type(method_type_short) :: text_method_short
890 type(method_type) :: text_method
891 type(method_type_very_short) :: text_method_very_short
892 
893 if (.not.PRESENT(table_name)) then
894  tbl_name = 'field_table'
895 else
896  tbl_name = trim(table_name)
897 endif
898 if (.not. file_exists(trim(tbl_name))) then
899  if(present(nfields)) nfields = 0
900  return
901 endif
902 
903 allocate(fields(max_fields))
904 
905 open(newunit=iunit, file=trim(tbl_name), action='READ', iostat=io_status)
906 if(io_status/=0) call mpp_error(fatal, 'field_manager_mod: Error in opening file '//trim(tbl_name))
907 !write_version_number should precede all writes to stdlog from field_manager
908 call write_version_number("FIELD_MANAGER_MOD", version)
909 log_unit = stdlog()
910 do while (.true.)
911  read(iunit,'(a)',end=89,err=99) record
912  write( log_unit,'(a)' )record
913  if (record(1:1) == "#" ) cycle
914  ltrec = len_trim(record)
915  if (ltrec .le. 0 ) cycle ! Blank line
916 
917 
918  icount = 0
919  do l= 1, ltrec
920  if (record(l:l) == '"' ) then
921  icount = icount + 1
922  endif
923  enddo
924  if (icount > 6 ) then
925  call mpp_error(fatal,trim(error_header)//'Too many fields in field table header entry.'//trim(record))
926  endif
927 
928  select case (icount)
929  case (6)
930  read(record,*,end=79,err=79) text_names
931  text_names%fld_type = lowercase(trim(text_names%fld_type))
932  text_names%mod_name = lowercase(trim(text_names%mod_name))
933  text_names%fld_name = lowercase(trim(text_names%fld_name))
934  case(4)
935 ! If there is no control string then the last string can be omitted and there are only 4 '"' in the record.
936  read(record,*,end=79,err=79) text_names_short
937  text_names%fld_type = lowercase(trim(text_names_short%fld_type))
938  text_names%mod_name = lowercase(trim(text_names_short%mod_name))
939  text_names%fld_name = lowercase(trim(text_names_short%mod_name))
940  case(2)
941 ! If there is only the method_type string then the last 2 strings need to be blank and there
942 ! are only 2 '"' in the record.
943  read(record,*,end=79,err=79) text_names_short
944  text_names%fld_type = lowercase(trim(text_names_short%fld_type))
945  text_names%mod_name = lowercase(trim(text_names_short%mod_name))
946  text_names%fld_name = lowercase(trim(text_names_short%mod_name))
947  case default
948 ! There is an unterminated or unquoted string in the field table entry.
949  text_names%fld_type = " "
950  text_names%mod_name = lowercase(trim(record))
951  text_names%fld_name = " "
952  end select
953 
954 ! Create a list with Rick Slaters field manager code
955 
956  list_name = list_sep//trim(text_names%mod_name)//list_sep//trim(text_names%fld_type)//&
957  list_sep//trim(text_names%fld_name)
958 
959  index_list_name = fm_new_list(list_name, create = .true.)
960  if ( index_list_name == no_field ) &
961  call mpp_error(fatal, trim(error_header)//'Could not set field list for '//trim(list_name))
962 
963  fm_success = fm_change_list(list_name)
964  select case (text_names%mod_name)
965  case ('coupler_mod')
966  model = model_coupler
967  case ('atmos_mod')
968  model = model_atmos
969  case ('ocean_mod')
970  model = model_ocean
971  case ('land_mod')
972  model = model_land
973  case ('ice_mod')
974  model = model_ice
975  case default
976  call mpp_error(fatal, trim(error_header)//'The model name is unrecognised : '//trim(text_names%mod_name))
977  end select
978  if (find_field_index(list_name) > 0) then
979  num_fields = num_fields + 1
980 
981  if (num_fields > max_fields) call mpp_error(fatal,trim(error_header)//'max fields exceeded')
982  fields(num_fields)%model = model
983  fields(num_fields)%field_name = lowercase(trim(text_names%fld_name))
984  fields(num_fields)%field_type = lowercase(trim(text_names%fld_type))
985  fields(num_fields)%num_methods = 0
986  allocate(fields(num_fields)%methods(max_field_methods))
987  call check_for_name_duplication
988 
989 ! Check to see that the first line is not the only line
990  if ( record(len_trim(record):len_trim(record)) == list_sep) cycle
991 
992  flag_method = .true.
993  m = 1
994  do while (flag_method)
995  read(iunit,'(a)',end=99,err=99) record
996 ! If the line is blank then fetch the next line.
997  if (len_trim(record) .le. 0) cycle
998 ! If the last character in the line is / then this is the end of the field methods
999  if ( record(len_trim(record):len_trim(record)) == list_sep) then
1000  flag_method = .false.
1001  if (len_trim(record) == 1) cycle
1002  record = record(:len_trim(record)-1) ! Remove the end of field method marker
1003  endif
1004 ! If the line is now blank, after removing the field separator marker, then fetch the next line.
1005  if (len_trim(record) .le. 0) cycle
1006 ! If the first character in the line is # then it is treated as a comment
1007  if (record(1:1) == comment ) cycle
1008 
1009  icount = 0
1010  do l= 1, len_trim(record)
1011  if (record(l:l) == dquote ) then
1012  icount = icount + 1
1013  endif
1014  enddo
1015  if (icount > 6 ) call mpp_error(fatal,trim(error_header)//'Too many fields in field entry.'//trim(record))
1016 
1017  if (.not. fm_change_list( list_name)) &
1018  call mpp_error(fatal, trim(error_header)//'Could not change to '//trim(list_name)//' list')
1019 
1020  select case (icount)
1021  case (6)
1022  read(record,*,end=99,err=99) text_method
1023  fields(num_fields)%methods(m)%method_type = lowercase(trim(text_method%method_type))
1024  fields(num_fields)%methods(m)%method_name = lowercase(trim(text_method%method_name))
1025  fields(num_fields)%methods(m)%method_control = lowercase(trim(text_method%method_control))
1026 
1027  type_str = text_method%method_type
1028  name_str = text_method%method_name
1029  control_str = text_method%method_control
1030 
1031  case(4)
1032 ! If there is no control string then the last string can be omitted and there are only 4 '"' in the record.
1033  read(record,*,end=99,err=99) text_method_short
1034  fields(num_fields)%methods(m)%method_type =&
1035  & lowercase(trim(text_method_short%method_type))
1036  fields(num_fields)%methods(m)%method_name =&
1037  & lowercase(trim(text_method_short%method_name))
1038  fields(num_fields)%methods(m)%method_control = " "
1039 
1040  type_str = text_method_short%method_type
1041  name_str = ""
1042  control_str = text_method_short%method_name
1043 
1044  case(2)
1045 ! If there is only the method_type string then the last 2 strings need to be blank and there
1046 ! are only 2 '"' in the record.
1047  read(record,*,end=99,err=99) text_method_very_short
1048  fields(num_fields)%methods(m)%method_type = lowercase(trim(text_method_very_short%method_type))
1049  fields(num_fields)%methods(m)%method_name = " "
1050  fields(num_fields)%methods(m)%method_control = " "
1051 
1052  type_str = ""
1053  name_str = ""
1054  control_str = text_method_very_short%method_type
1055 
1056  case(0)
1057  read(record,'(A)',end=99,err=99) control_str
1058  type_str = ""
1059  name_str = ""
1060 
1061  case default
1062  call mpp_error(fatal,trim(error_header)//'Unterminated field in field entry.'//trim(record))
1063  end select
1064 
1065 ! This section of code breaks the control string into separate strings.
1066 ! The array control_array contains the following parameters.
1067 ! control_array(:,1) = index within control_str of the first character of the name.
1068 ! control_array(:,2) = index within control_str of the equal sign
1069 ! control_array(:,3) = index within control_str of the last character of the value.
1070 !
1071 ! control_array(:,1) -> control_array(:,2) -1 = name of the parameter.
1072 ! control_array(:,2)+1 -> control_array(:,3) = value of the parameter.
1073 
1074  ltrec= len_trim(control_str)
1075  control_array(:,1) = 1
1076  control_array(:,2:3) = ltrec
1077  icount = 0
1078  do l= 1, ltrec
1079  if (control_str(l:l) == equal ) then
1080  icount = icount + 1
1081  control_array(icount,2) = l ! Middle of string
1082  elseif (control_str(l:l) == comma ) then
1083  if (icount .eq. 0) then
1084  call mpp_error(fatal,trim(error_header) // &
1085  ' Bad format for field entry (comma without equals sign): ''' // &
1086  trim(control_str) // '''')
1087  elseif (icount .gt. max_fields) then
1088  call mpp_error(fatal,trim(error_header) // &
1089  ' Too many fields in field entry: ''' // &
1090  trim(control_str) // '''')
1091  else
1092  control_array(icount,3) = l-1 !End of previous string
1093  control_array(min(max_fields,icount+1),1) = l+1 !Start of next string
1094  endif
1095  endif
1096  enddo
1097 
1098  ! Make sure that we point to the end of the string (minus any trailing comma)
1099  ! for the last set of values. This fixes the case where the last set of values
1100  ! is a comma separated list
1101 
1102  if (control_str(ltrec:ltrec) .ne. comma) then
1103  control_array(max(1,icount),3) = ltrec
1104  endif
1105 
1106  if ( icount == 0 ) then
1107  method_name = type_str
1108  if (len_trim(method_name) > 0 ) then
1109  method_name = trim(method_name)//list_sep// trim(name_str)
1110  else
1111  method_name = trim(name_str)
1112  endif
1113  val_name = control_str
1114 
1115  call new_name(list_name, method_name, val_name )
1116 
1117  else
1118 
1119  do l = 1,icount
1120  startcont = control_array(l,1)
1121  midcont = control_array(l,2)
1122  endcont = control_array(l,3)
1123 
1124  method_name = trim(type_str)
1125  if (len_trim(method_name) > 0 ) then
1126  method_name = trim(method_name)//list_sep// trim(name_str)
1127  else
1128  method_name = trim(name_str)
1129  endif
1130 
1131  if (len_trim(method_name) > 0 ) then
1132  method_name = trim(method_name)//list_sep//&
1133  trim(control_str(startcont:midcont-1))
1134  else
1135  method_name = trim(control_str(startcont:midcont-1))
1136  endif
1137  val_name = trim(control_str(midcont+1:endcont))
1138 
1139  call new_name(list_name, method_name, val_name )
1140  enddo
1141 
1142  endif
1143 
1144  fields(num_fields)%num_methods = fields(num_fields)%num_methods + 1
1145  if (fields(num_fields)%num_methods > max_field_methods) &
1146  call mpp_error(fatal,trim(error_header)//'Maximum number of methods for field exceeded')
1147  m = m + 1
1148  enddo
1149  else
1150  flag_method = .true.
1151  do while (flag_method)
1152  read(iunit,'(A)',end=99,err=99) record
1153  if ( record(len_trim(record):len_trim(record)) == list_sep) then
1154  flag_method = .false.
1155  endif
1156  enddo
1157  endif
1158 79 continue
1159 enddo
1160 
1161 89 continue
1162 close(iunit, iostat=io_status)
1163 if(io_status/=0) call mpp_error(fatal, 'field_manager_mod: Error in closing file '//trim(tbl_name))
1164 
1165 
1166 if(present(nfields)) nfields = num_fields
1167 
1168 default_method%method_type = 'none'
1169 default_method%method_name = 'none'
1170 default_method%method_control = 'none'
1171 return
1172 
1173 99 continue
1174 
1175 call mpp_error(fatal,trim(error_header)//' Error reading field table. Record = '//trim(record))
1176 
1177 end subroutine read_field_table_legacy
1178 
1179 subroutine check_for_name_duplication
1180 integer :: i
1181 
1182 ! Check that name is unique amoung fields of the same field_type and model.
1183 do i=1,num_fields-1
1184  if ( fields(i)%field_type == fields(num_fields)%field_type .and. &
1185  fields(i)%model == fields(num_fields)%model .and. &
1186  fields(i)%field_name == fields(num_fields)%field_name ) then
1187  if (mpp_pe() .eq. mpp_root_pe()) then
1188  call mpp_error(warning,'Error in field_manager_mod. Duplicate field name: Field type='//&
1189  trim(fields(i)%field_type)// &
1190  ', Model='//trim(model_names(fields(i)%model))// &
1191  ', Duplicated name='//trim(fields(i)%field_name))
1192  endif
1193  endif
1194 enddo
1195 
1196 end subroutine check_for_name_duplication
1197 
1198 !> @brief Subroutine to add new values to list parameters.
1199 !!
1200 !> This subroutine uses input strings list_name, method_name
1201 !! and val_name_in to add new values to the list. Given
1202 !! list_name a new list item is created that is named
1203 !! method_name and is given the value or values in
1204 !! val_name_in. If there is more than 1 value in
1205 !! val_name_in, these values should be comma-separated.
1206 subroutine new_name ( list_name, method_name_in , val_name_in)
1207 character(len=*), intent(in) :: list_name !< The name of the field that is of interest here.
1208 character(len=*), intent(in) :: method_name_in !< The name of the method that values are
1209  !! being supplied for.
1210 character(len=*), intent(inout) :: val_name_in !< The value or values that will be parsed and
1211  !! used as the value when creating a new field or fields.
1212 
1213 character(len=fm_string_len) :: method_name
1214 character(len=fm_string_len) :: val_name
1215 integer, dimension(MAX_FIELDS) :: end_val
1216 integer, dimension(MAX_FIELDS) :: start_val
1217 integer :: i
1218 integer :: index_t
1219 integer :: left_br
1220 integer :: num_elem
1221 integer :: out_unit
1222 integer :: right_br
1223 integer :: val_int
1224 integer :: val_type
1225 logical :: append_new
1226 logical :: val_logic
1227 real(r8_kind) :: val_real !< all reals converted from string will be in r8_kind precision
1228 integer :: length
1229 
1230 call strip_front_blanks(val_name_in)
1231 method_name = trim(method_name_in)
1232 call strip_front_blanks(method_name)
1233 
1234 index_t = 1
1235 num_elem = 1
1236 append_new = .false.
1237 start_val(1) = 1
1238 end_val(:) = len_trim(val_name_in)
1239 
1240 ! If the array of values being passed in is a comma delimited list then count
1241 ! the number of elements.
1242 
1243 do i = 1, len_trim(val_name_in)
1244  if ( val_name_in(i:i) == comma ) then
1245  end_val(num_elem) = i-1
1246  start_val(num_elem+1) = i+1
1247  num_elem = num_elem + 1
1248  endif
1249 enddo
1250 
1251 ! Check to see if this is an array element of form array[x] = value
1252 left_br = scan(method_name,'[')
1253 right_br = scan(method_name,']')
1254 if ( num_elem .eq. 1 ) then
1255  if ( left_br > 0 .and. right_br == 0 ) &
1256  call mpp_error(fatal, trim(error_header)//"Left bracket present without right bracket in "//trim(method_name))
1257  if ( left_br== 0 .and. right_br > 0 ) &
1258  call mpp_error(fatal, trim(error_header)//"Right bracket present without left bracket in "//trim(method_name))
1259  if ( left_br > 0 .and. right_br > 0 ) then
1260  if ( scan( method_name(left_br+1:right_br -1), set ) > 0 ) &
1261  call mpp_error(fatal, trim(error_header)//"Using a non-numeric value for index in "//trim(method_name))
1262  read(method_name(left_br+1:right_br -1), *) index_t
1263  method_name = method_name(:left_br -1)
1264  endif
1265 else
1266 ! If there are multiple values then there cannot be a bracket in method_name.
1267  if ( left_br > 0 .or. right_br > 0 ) &
1268  call mpp_error(fatal, &
1269  trim(error_header)//"Using a comma delimited list with an indexed array element in "//trim(method_name))
1270 endif
1271 
1272 do i = 1, num_elem
1273 
1274  if ( i .gt. 1 .or. index_t .eq. 0 ) then
1275  append_new = .true.
1276  index_t = 0 ! If append is true then index must be <= 0
1277  endif
1278  val_type = string_type ! Assume it is a string
1279  val_name = val_name_in(start_val(i):end_val(i))
1280  call strip_front_blanks(val_name)
1281 
1282 !
1283 ! if the string starts and ends with matching single quotes, then this is a string
1284 ! if there are quotes which do not match, then this is an error
1285 !
1286 
1287  length = len_trim(val_name)
1288  if (val_name(1:1) .eq. squote) then
1289 
1290  if (val_name(length:length) .eq. squote) then
1291  val_name = val_name(2:length-1)//repeat(" ",len(val_name)-length+2)
1292  val_type = string_type
1293  elseif (val_name(length:length) .eq. dquote) then
1294  call mpp_error(fatal, trim(error_header) // ' Quotes do not match in ' // trim(val_name) // &
1295  ' for ' // trim(method_name) // ' of ' // trim(list_name))
1296  else
1297  call mpp_error(fatal, trim(error_header) // ' No trailing quote in ' // trim(val_name) // &
1298  ' for ' // trim(method_name) // ' of ' // trim(list_name))
1299  endif
1300 
1301  elseif (val_name(1:1) .eq. dquote .or. val_name(length:length) .eq. dquote) then
1302 
1303  call mpp_error(fatal, trim(error_header) // ' Double quotes not allowed in ' // trim(val_name) // &
1304  ' for ' // trim(method_name) // ' of ' // trim(list_name))
1305 
1306  elseif (val_name(length:length) .eq. squote) then
1307 
1308  call mpp_error(fatal, trim(error_header) // ' No leading quote in ' // trim(val_name) // &
1309  ' for ' // trim(method_name) // ' of ' // trim(list_name))
1310 
1311  else
1312 ! If the string to be parsed is a real then all the characters must be numeric,
1313 ! be a plus/minus, be a decimal point or, for exponentials, be e or E.
1314 
1315 ! If a string is an integer, then all the characters must be numeric.
1316 
1317  if ( scan(val_name(1:1), setnum ) > 0 ) then
1318 
1319 ! If there is a letter in the name it may only be e or E
1320 
1321  if ( scan(val_name, set_nonexp ) .le. 0 ) then
1322 ! It is real if there is a . in the name or the value appears exponential
1323  if ( scan(val_name, '.') > 0 .or. scan(val_name, 'e') > 0 .or. scan(val_name, 'E') > 0) then
1324  read(val_name, *) val_real
1325  val_type = real_type
1326  else
1327  read(val_name, *) val_int
1328  val_type = integer_type
1329  endif
1330  endif
1331 
1332  endif
1333 
1334 ! If val_name is t/T or f/F then this is a logical flag.
1335  if ( len_trim(val_name) == 1 .or. len_trim(val_name) == 3) then
1336  if ( val_name == 't' .or. val_name == 'T' .or. val_name == '.t.' .or. val_name == '.T.' ) then
1337  val_logic = .true.
1338  val_type = logical_type
1339  endif
1340  if ( val_name == 'f' .or. val_name == 'F' .or. val_name == '.f.' .or. val_name == '.F.' ) then
1341  val_logic = .false.
1342  val_type = logical_type
1343  endif
1344  endif
1345  if ( trim(lowercase(val_name)) == 'true' .or. trim(lowercase(val_name)) == '.true.' ) then
1346  val_logic = .true.
1347  val_type = logical_type
1348  endif
1349  if ( trim(lowercase(val_name)) == 'false' .or. trim(lowercase(val_name)) == '.false.' ) then
1350  val_logic = .false.
1351  val_type = logical_type
1352  endif
1353  endif
1354 
1355  select case(val_type)
1356 
1357  case (integer_type)
1358  if ( fm_new_value( method_name, val_int, create = .true., index = index_t, append = append_new ) < 0 ) &
1359  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
1360  ' (I) for '//trim(list_name))
1361 
1362  case (logical_type)
1363  if ( fm_new_value( method_name, val_logic, create = .true., index = index_t, append = append_new) < 0 ) &
1364  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
1365  ' (L) for '//trim(list_name))
1366 
1367  case (real_type)
1368  if ( fm_new_value( method_name, val_real, create = .true., index = index_t, append = append_new) < 0 ) &
1369  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
1370  ' (R) for '//trim(list_name))
1371 
1372  case (string_type)
1373  if ( fm_new_value( method_name, val_name, create = .true., index = index_t, append = append_new) < 0 ) &
1374  call mpp_error(fatal, trim(error_header)//'Could not set "' // trim(val_name) // '" for '//trim(method_name)//&
1375  ' (S) for '//trim(list_name))
1376  case default
1377  call mpp_error(fatal, trim(error_header)//'Could not find a valid type to set the '//trim(method_name)//&
1378  ' for '//trim(list_name))
1379 
1380  end select
1381 
1382 enddo
1383 
1384 end subroutine new_name
1385 
1386 !> @brief Destructor for field manager.
1387 !!
1388 !> This subroutine deallocates allocated variables (if allocated) and
1389 !! changes the initialized flag to false.
1391 integer :: j
1392 
1393 module_is_initialized = .false.
1394 
1395 do j=1,size(fields)
1396  if(allocated(fields(j)%methods)) deallocate(fields(j)%methods)
1397 end do
1398 if(allocated(fields)) deallocate(fields)
1399 
1400 end subroutine field_manager_end
1401 
1402 !> @brief A routine to strip whitespace from the start of character strings.
1403 !!
1404 !> This subroutine removes spaces and tabs from the start of a character string.
1405 subroutine strip_front_blanks(name)
1406 
1407 character(len=*), intent(inout) :: name !< name to remove whitespace from
1408 
1409 name = trim(adjustl(name))
1410 end subroutine strip_front_blanks
1411 
1412 !> @brief Function to return the index of the field
1413 !!
1414 !> This function when passed a model number and a field name will
1415 !! return the index of the field within the field manager. This index
1416 !! can be used to access other information from the field manager.
1417 !! @returns The index of the field corresponding to field_name.
1418 function find_field_index_old(model, field_name)
1419 
1420 integer :: find_field_index_old
1421 integer, intent(in) :: model !< The number indicating which model is used.
1422 character(len=*), intent(in) :: field_name !< The name of the field that an index is being requested for.
1423 
1424 integer :: i
1425 
1427 
1428 do i=1,num_fields
1429  if (fields(i)%model == model .and. fields(i)%field_name == lowercase(field_name)) then
1431  return
1432  endif
1433 enddo
1434 
1435 end function find_field_index_old
1436 
1437 !> @returns index of the field corresponding to field_name
1438 function find_field_index_new(field_name)
1439 
1440 integer :: find_field_index_new
1441 character(len=*), intent(in) :: field_name !< The path to the name of the field that an index is
1442  !! being requested for.
1443 
1445 
1446 find_field_index_new = fm_get_index(field_name)
1447 
1448 end function find_field_index_new
1449 
1450 !> @brief This routine allows access to field information given an index.
1451 !!
1452 !> When passed an index, this routine will return the type of field,
1453 !! the name of the field, the model which the field is associated and
1454 !! the number of methods associated with the field.
1455 !! <br>Example usage:
1456 !! @code{.F90}
1457 !! call get_field_info( n,fld_type,fld_name,model,num_methods )
1458 !! @endcode
1459 subroutine get_field_info(n,fld_type,fld_name,model,num_methods)
1460 integer, intent(in) :: n !< index of field
1461 character (len=*),intent(out) :: fld_type !< field type
1462 character (len=*),intent(out) :: fld_name !< name of the field
1463 integer, intent(out) :: model !< number indicating which model is used
1464 integer, intent(out) :: num_methods !< number of methods
1465 
1466 if (n < 1 .or. n > num_fields) call mpp_error(fatal,trim(error_header)//'Invalid field index')
1467 
1468 fld_type = fields(n)%field_type
1469 fld_name = fields(n)%field_name
1470 model = fields(n)%model
1471 num_methods = fields(n)%num_methods
1472 
1473 end subroutine get_field_info
1474 
1475 !> @brief A routine to get a specified method
1476 !!
1477 !> This routine, when passed a field index and a method index will
1478 !! return the method text associated with the field(n) method(m).
1479 subroutine get_field_method(n,m,method)
1480 
1481 integer, intent(in) :: n !< index of field
1482 integer, intent(in) :: m !< index of method
1483 type(method_type) ,intent(inout) :: method !< the m-th method of field with index n
1484 
1485 if (n < 1 .or. n > num_fields) call mpp_error(fatal,trim(error_header)//'Invalid field index')
1486 if (m < 1 .or. m > fields(n)%num_methods) call mpp_error(fatal,trim(error_header)//'Invalid method index')
1487 
1488  method = fields(n)%methods(m)
1489 
1490 end subroutine get_field_method
1491 
1492 !> @brief A routine to obtain all the methods associated with a field.
1493 !!
1494 !> When passed a field index, this routine will return the text
1495 !! associated with all the methods attached to the field.
1496 subroutine get_field_methods(n,methods)
1497 
1498 integer, intent(in) :: n !< field index
1499 type(method_type),intent(inout) :: methods(:) !< an array of methods for field with index n
1500 
1501  if (n < 1 .or. n > num_fields) &
1502  call mpp_error(fatal,trim(error_header)//'Invalid field index')
1503 
1504  if (size(methods(:)) < fields(n)%num_methods) &
1505  call mpp_error(fatal,trim(error_header)//'Method array too small')
1506 
1507  methods = default_method
1508  methods(1:fields(n)%num_methods) = fields(n)%methods(1:fields(n)%num_methods)
1509 
1510 end subroutine get_field_methods
1511 
1512 !> @returns The number of values that have been decoded. This allows
1513 !! a user to define a large array and fill it partially with
1514 !! values from a list. This should be the size of the value array.
1515 function parse_integers ( text, label, values ) result (parse)
1516 character(len=*), intent(in) :: text !< The text string from which the values will be parsed.
1517 character(len=*), intent(in) :: label !< A label which describes the values being decoded.
1518 integer, intent(out) :: values(:) !< The value or values that have been decoded.
1519 
1520 include 'parse.inc'
1521 end function parse_integers
1522 
1523 function parse_strings ( text, label, values ) result (parse)
1524 character(len=*), intent(in) :: text !< The text string from which the values will be parsed.
1525 character(len=*), intent(in) :: label !< A label which describes the values being decoded.
1526 character(len=*), intent(out) :: values(:) !< The value or values that have been decoded.
1527 
1528 include 'parse.inc'
1529 end function parse_strings
1530 
1531 function parse_integer ( text, label, parse_ival ) result (parse)
1532 character(len=*), intent(in) :: text !< The text string from which the values will be parsed.
1533 character(len=*), intent(in) :: label !< A label which describes the values being decoded.
1534 integer, intent(out) :: parse_ival !< The value or values that have been decoded.
1535 integer :: parse
1536 
1537 integer :: values(1)
1538 
1539  parse = parse_integers( text, label, values )
1540  if (parse > 0) parse_ival = values(1)
1541 end function parse_integer
1542 
1543 function parse_string ( text, label, parse_sval ) result (parse)
1544 character(len=*), intent(in) :: text !< The text string from which the values will be parsed.
1545 character(len=*), intent(in) :: label !< A label which describes the values being decoded.
1546 character(len=*), intent(out) :: parse_sval !< The value or values that have been decoded.
1547 integer :: parse
1548 
1549 character(len=len(parse_sval)) :: values(1)
1550 
1551  parse = parse_strings( text, label, values )
1552  if (parse > 0) parse_sval = values(1)
1553 end function parse_string
1554 
1555 !> @brief A function to create a field as a child of parent_p. This will return
1556 !! a pointer to a field_def type.
1557 !!
1558 !> Allocate and initialize a new field in parent_p list.
1559 !! Return a pointer to the field on success, or a null pointer
1560 !! on failure.
1561 !! <br>Example usage:
1562 !! @code{.F90}
1563 !! list_p => create_field(parent_p, name)
1564 !! @endcode
1565 function create_field(parent_p, name) &
1566  result(list_p)
1567 type (field_def), pointer :: list_p
1568 type (field_def), pointer :: parent_p !< A pointer to the parent of the field that is to be created
1569 character(len=*), intent(in) :: name !< The name of the field that is to be created
1570 
1571 integer :: error, out_unit
1572 ! Check for fatal errors which should never arise
1573 out_unit = stdout()
1574 if (.not. associated(parent_p) .or. name .eq. ' ') then
1575  nullify(list_p)
1576  return
1577 endif
1578 
1579 ! Allocate space for the new list
1580 allocate(list_p, stat = error)
1581 if (error .ne. 0) then
1582  write (out_unit,*) trim(error_header), 'Error ', error, &
1583  ' allocating memory for list ', trim(name)
1584  nullify(list_p)
1585  return
1586 endif
1587 ! Initialize the new field
1588 list_p%name = name
1589 
1590 nullify(list_p%next)
1591 list_p%prev => parent_p%last_field
1592 nullify(list_p%first_field)
1593 nullify(list_p%last_field)
1594 list_p%length = 0
1595 list_p%field_type = null_type
1596 list_p%max_index = 0
1597 list_p%array_dim = 0
1598 if (allocated(list_p%i_value)) deallocate(list_p%i_value)
1599 if (allocated(list_p%l_value)) deallocate(list_p%l_value)
1600 if (allocated(list_p%r_value)) deallocate(list_p%r_value)
1601 if (allocated(list_p%s_value)) deallocate(list_p%s_value)
1602 ! If this is the first field in the parent, then set the pointer
1603 ! to it, otherwise, update the "next" pointer for the last list
1604 if (parent_p%length .le. 0) then
1605  parent_p%first_field => list_p
1606 else
1607  parent_p%last_field%next => list_p
1608 endif
1609 ! Update the pointer for the last list in the parent
1610 parent_p%last_field => list_p
1611 ! Update the length for the parent
1612 parent_p%length = parent_p%length + 1
1613 ! Set the new index as the return value
1614 list_p%index = parent_p%length
1615 ! set the pointer to the parent list
1616 list_p%parent => parent_p
1617 
1618 end function create_field
1619 
1620 !> @brief This is a function that lists the parameters of a field.
1621 !!
1622 !> Given a pointer to a list, this function prints out the fields, and
1623 !! subfields, if recursive is true, associated with the list.
1624 !!
1625 !! This is most likely to be used through fm_dump_list.
1626 !! <br> Example usage:
1627 !! @code{.F90}
1628 !! success = dump_list(list_p, recursive=.true., depth=0)
1629 !! @endcode
1630 logical recursive function dump_list(list_p, recursive, depth, out_unit) result(success)
1631 
1632  type (field_def), pointer :: list_p !< pointer to the field to be printed out
1633  logical, intent(in) :: recursive !< flag to make function recursively print subfields
1634  integer, intent(in) :: depth !< Listing will be padded so that 'depth' spaces appear before
1635  !! the field being printed
1636  integer, intent(in) :: out_unit !< unit number to print to
1637 
1638  integer :: depthp1
1639  integer :: j
1640  character(len=fm_field_name_len) :: num, scratch
1641  type (field_def), pointer :: this_field_p
1642  character(len=depth+fm_field_name_len) :: blank
1643 
1644  blank = ' ' ! initialize blank string
1645 
1646  ! Check for a valid list
1647  success = .false.
1648  if (.not. associated(list_p)) then
1649  return
1650  elseif (list_p%field_type .ne. list_type) then
1651  return
1652  endif
1653 
1654  ! set the default return value
1655  success = .true.
1656 
1657  ! Print the name of this list
1658  write (out_unit,'(a,a,a)') blank(1:depth), trim(list_p%name), list_sep
1659 
1660  ! Increment the indentation depth
1661  ! The following max function is to work around an error in the IBM compiler for len_trim
1662  ! depthp1 = depth + max(len_trim(list_p%name),0) + len_trim(list_sep)
1663  depthp1 = depth + 6
1664 
1665  this_field_p => list_p%first_field
1666 
1667  do while (associated(this_field_p))
1668 
1669  select case(this_field_p%field_type)
1670  case(list_type)
1671  ! If this is a list, then call dump_list
1672  if (recursive) then
1673  ! If recursive is true, then this routine will find and dump sub-fields.
1674  success = dump_list(this_field_p, .true., depthp1, out_unit)
1675  if (.not.success) exit ! quit immediately in case of error
1676  else ! Otherwise it will print out the name of this field.
1677  write (out_unit,'(a,a,a)') blank(1:depthp1), trim(this_field_p%name), list_sep
1678  endif
1679 
1680  case(integer_type)
1681  if (this_field_p%max_index .eq. 0) then
1682  write (out_unit,'(a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = NULL'
1683  elseif (this_field_p%max_index .eq. 1) then
1684  write (scratch,*) this_field_p%i_value(1)
1685  write (out_unit,'(a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = ', &
1686  trim(adjustl(scratch))
1687  else ! Write out the array of values for this field.
1688  do j = 1, this_field_p%max_index
1689  write (scratch,*) this_field_p%i_value(j)
1690  write (num,*) j
1691  write (out_unit,'(a,a,a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), &
1692  '[', trim(adjustl(num)), '] = ', trim(adjustl(scratch))
1693  enddo
1694  endif
1695 
1696  case(logical_type)
1697  if (this_field_p%max_index .eq. 0) then
1698  write (out_unit,'(a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = NULL'
1699  elseif (this_field_p%max_index .eq. 1) then
1700  write (scratch,'(l1)') this_field_p%l_value(1)
1701  write (out_unit,'(a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = ', &
1702  trim(adjustl(scratch))
1703  else ! Write out the array of values for this field.
1704  do j = 1, this_field_p%max_index
1705  write (scratch,'(l1)') this_field_p%l_value(j)
1706  write (num,*) j
1707  write (out_unit,'(a,a,a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), &
1708  '[', trim(adjustl(num)), '] = ', trim(adjustl(scratch))
1709  enddo
1710  endif
1711 
1712  case(real_type)
1713  if (this_field_p%max_index .eq. 0) then
1714  write (out_unit,'(a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = NULL'
1715  elseif (this_field_p%max_index .eq. 1) then
1716  write (scratch,*) this_field_p%r_value(1)
1717  write (out_unit,'(a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = ', &
1718  trim(adjustl(scratch))
1719  else ! Write out the array of values for this field.
1720  do j = 1, this_field_p%max_index
1721  write (scratch,*) this_field_p%r_value(j)
1722  write (num,*) j
1723  write (out_unit,'(a,a,a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), &
1724  '[', trim(adjustl(num)), '] = ', trim(adjustl(scratch))
1725  end do
1726  endif
1727 
1728  case(string_type)
1729  if (this_field_p%max_index .eq. 0) then
1730  write (out_unit,'(a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = NULL'
1731  elseif (this_field_p%max_index .eq. 1) then
1732  write (out_unit,'(a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), ' = ', &
1733  ''''//trim(this_field_p%s_value(1))//''''
1734  else ! Write out the array of values for this field.
1735  do j = 1, this_field_p%max_index
1736  write (num,*) j
1737  write (out_unit,'(a,a,a,a,a,a)') blank(1:depthp1), trim(this_field_p%name), &
1738  '[', trim(adjustl(num)), '] = ', ''''//trim(this_field_p%s_value(j))//''''
1739  enddo
1740  endif
1741 
1742  case default
1743  success = .false.
1744  exit
1745 
1746  end select
1747 
1748  this_field_p => this_field_p%next
1749  enddo
1750 
1751 end function dump_list
1752 
1753 !> @brief A subroutine that splits a listname into a path and a base.
1754 !!
1755 !> Find the base name for a list by splitting the list name into
1756 !! a path and base. The base is the last field within name, while the
1757 !! path is the preceding section of name. The base string can then be
1758 !! used to query for values associated with name.
1759 subroutine find_base(name, path, base)
1760 
1761 character(len=*), intent(in) :: name !< list name for a field
1762 character(len=*), intent(out) :: path !< path of the base field
1763 character(len=*), intent(out) :: base !< A string which can be used to query for values associated with name
1764 
1765 integer :: i
1766 integer :: length
1767 
1768 ! Check for the last occurrence of the list separator in name
1769 ! The following max function is to work around an error in the IBM compiler for len_trim
1770 length = max(len_trim(name),0)
1771 if (length .eq. 0) then
1772 
1773  ! Empty name, so return empty path and base
1774  path = ' '
1775  base = ' '
1776 else
1777  ! Remove trailing list separators
1778  do while (name(length:length) .eq. list_sep)
1779  length = length - 1
1780  if (length .eq. 0) then
1781  exit
1782  endif
1783  enddo
1784  if (length .eq. 0) then
1785 
1786  ! Name only list separators, so return empty path and base
1787  path = ' '
1788  base = ' '
1789  else
1790  ! Check for the last occurrence of the list separator in name
1791  i = index(name(1:length), list_sep, back = .true.)
1792  if (i .eq. 0) then
1793  ! no list separators in the path, so return an empty path
1794  ! and name as the base
1795  path = ' '
1796  base = name(1:length)
1797  else
1798  ! Found a list separator, so return the part up to the last
1799  ! list separator in path, and the remainder in base
1800  path = name(1:i)
1801  base = name(i+1:length)
1802  endif
1803  endif
1804 endif
1805 
1806 end subroutine find_base
1807 
1808 !> @brief Find and return a pointer to the field in the specified
1809 !! list. Return a null pointer on error.
1810 !!
1811 !> Find and return a pointer to the field in the specified
1812 !! list. Return a null pointer on error. Given a pointer to a field,
1813 !! this function searchs for "name" as a sub field.
1814 !> @returns A pointer to the field corresponding to "name" or an unassociated pointer if the field
1815 !! name does not exist.
1816 function find_field(name, this_list_p) &
1817  result(field_p)
1818 type (field_def), pointer :: field_p
1819 character(len=*), intent(in) :: name !< The name of a field that the user wishes to find
1820 type (field_def), pointer :: this_list_p !< A pointer to a list which the user wishes to search
1821  !! for a field "name".
1822 
1823 type (field_def), pointer, save :: temp_p
1824 
1825 
1826 nullify (field_p)
1827 
1828 if (name .eq. '.') then
1829 
1830 ! If the field is '.' then return this list
1831  field_p => this_list_p
1832 elseif (name .eq. '..') then
1833 ! If the field is '..' then return the parent list
1834  field_p => this_list_p%parent
1835 else
1836 ! Loop over each field in this list
1837  temp_p => this_list_p%first_field
1838 
1839  do while (associated(temp_p))
1840 ! If the name matches, then set the return pointer and exit
1841 ! the loop
1842  if (temp_p%name .eq. name) then
1843  field_p => temp_p
1844  exit
1845  endif
1846 
1847  temp_p => temp_p%next
1848 
1849  enddo
1850 endif
1851 
1852 end function find_field
1853 
1854 !> @brief Find the first list for a name by splitting the name into
1855 !! a head and the rest.
1856 !!
1857 !> Find the first list for a name by splitting the name into a head and the
1858 !! rest. The head is the first field within name, while rest is the remaining
1859 !! section of name. The head string can then be used to find other fields that
1860 !! may be associated with name.
1861 subroutine find_head(name, head, rest)
1862 
1863 character(len=*), intent(in) :: name !< The name of a field of interest
1864 character(len=*), intent(out) :: head !< the first field within name
1865 character(len=*), intent(out) :: rest !< the remaining section of name
1866 
1867 integer :: i
1868 ! Check for the first occurrence of the list separator in name
1869 i = index(name, list_sep)
1870 ! Check for additional consecutive list separators and return
1871 ! those also
1872 do while (i .le. len(name))
1873  if (name(i+1:i+1) .eq. list_sep) then
1874  i = i + 1
1875  else
1876  exit
1877  endif
1878 enddo
1879 
1880 if (i .eq. 0) then
1881 ! no list separators in the path, so return an empty head and
1882 ! name as the rest
1883  head = ' '
1884  rest = name
1885 elseif (i .eq. len(name)) then
1886 ! The last character in name is a list separator, so return name
1887 ! as head and an empty rest
1888  head = name
1889  rest = ' '
1890 else
1891 ! Found a list separator, so return the part up to the list
1892 ! separator in head, and the remainder in rest
1893  head = name(1:i)
1894  rest = name(i+1:)
1895 endif
1896 
1897 end subroutine find_head
1898 
1899 !> @brief Find and return a pointer to the specified list, relative to
1900 !! relative_p. Return a null pointer on error.
1901 !!
1902 !> This function, when supplied a pointer to a field and a name of a second
1903 !! field relative to that pointer, will find a list and return the pointer to
1904 !! the second field. If create is .true. and the second field does not exist,
1905 !! it will be created.
1906 !> @returns A pointer to the list to be returned
1907 function find_list(path, relative_p, create) &
1908  result(list_p)
1909 type (field_def), pointer :: list_p
1910 character(len=*), intent(in) :: path !< path to the list of interest
1911 type (field_def), pointer :: relative_p !< pointer to the list to which "path" is relative to
1912 logical, intent(in) :: create !< If the list does not exist, it will be created if set to true
1913 
1914 character(len=FMS_PATH_LEN) :: working_path
1915 character(len=FMS_PATH_LEN) :: rest
1916 character(len=fm_field_name_len) :: this_list
1917 integer :: i, out_unit
1918 type (field_def), pointer, save :: working_path_p
1919 type (field_def), pointer, save :: this_list_p
1920 
1921 out_unit = stdout()
1922 nullify(list_p)
1923 ! If the path is empty, then return the relative list
1924 if (path .eq. ' ') then
1925 
1926  list_p => relative_p
1927 
1928 else
1929 ! If a fully qualified path is given (i.e., starts with the
1930 ! list separator) then do everything relative to root,
1931 ! otherwise, do everything relative to relative list.
1932  if (path(1:1) .eq. list_sep) then
1933  working_path_p => root_p
1934  working_path = path(2:)
1935  else
1936  working_path_p => relative_p
1937  working_path = path
1938  endif
1939 ! Loop over each field in the path
1940  do while (working_path .ne. ' ')
1941 ! Get the first list in the working path
1942  call find_head(working_path, this_list, rest)
1943 ! If the first list is empty, then the 'rest' should hold the
1944 ! final field in the path
1945  if (this_list .eq. ' ') then
1946  this_list = rest
1947  rest = ' '
1948  endif
1949 ! Strip off trailing list separators
1950  i = len_trim(this_list)
1951  do while (i .gt. 0 .and. this_list(i:i) .eq. list_sep)
1952  this_list(i:i) = ' '
1953  i = i - 1
1954  enddo
1955 ! Find a pointer to this field in the working list
1956  this_list_p => find_field(this_list, working_path_p)
1957 
1958  if (.not. associated(this_list_p)) then
1959  if (create) then
1960 ! Create the list if so requested
1961  this_list_p => make_list(working_path_p, this_list)
1962  if (.not. associated(this_list_p)) then
1963  nullify(list_p)
1964  return
1965  endif
1966  else
1967 ! Otherwise, return an error
1968 
1969  nullify(list_p)
1970  return
1971  endif
1972  endif
1973 ! Make sure that the field found is a list, and if so, proceed to
1974 ! the next field in the path, otherwise, return an error
1975  if (this_list_p%field_type .eq. list_type) then
1976  working_path_p => this_list_p
1977  working_path = rest
1978  else
1979  nullify(list_p)
1980  return
1981  endif
1982  enddo
1983  list_p => working_path_p
1984 endif
1985 
1986 end function find_list
1987 
1988 !> @brief Change the current list. Return true on success, false otherwise
1989 !!
1990 !> This function changes the currect list to correspond to the list named name.
1991 !! If the first character of name is the list separator (/) then the list will
1992 !! search for "name" starting from the root of the field tree. Otherwise it
1993 !! will search for name starting from the current list.
1994 !! @return A flag to indicate operation success, true = no errors
1995 function fm_change_list(name) &
1996  result(success)
1997 logical :: success
1998 character(len=*), intent(in) :: name !< name of a list to change to
1999 
2000 type (field_def), pointer, save :: temp_p
2001 ! Initialize the field manager if needed
2002 if (.not. module_is_initialized) then
2004 endif
2005 ! Find the list if path is not empty
2006 temp_p => find_list(name, current_list_p, .false.)
2007 
2008 if (associated(temp_p)) then
2009  current_list_p => temp_p
2010  success = .true.
2011 else
2012  success = .false.
2013 endif
2014 
2015 end function fm_change_list
2016 
2017 !> @brief Change the root list
2018 !!
2019 !> This function changes the root of the field tree to correspond to the
2020 !! field named name. An example of a use of this would be if code is
2021 !! interested in a subset of fields with a common base. This common base
2022 !! could be set using fm_change_root and fields could be referenced using
2023 !! this root.
2024 !!
2025 !! This function should be used in conjunction with fm_return_root.
2026 !! @return A flag to indicate operation success, true = no errors
2027 function fm_change_root(name) &
2028  result(success)
2029 logical :: success
2030 character(len=*), intent(in) :: name !< name of the field which the user wishes to become the root.
2031 
2032 type (field_def), pointer, save :: temp_list_p
2033 integer :: out_unit
2034 ! Initialize the field manager if needed
2035 if (.not. module_is_initialized) then
2037 endif
2038 out_unit = stdout()
2039 ! Must supply a field field name
2040 if (name .eq. ' ') then
2041  success = .false.
2042  return
2043 endif
2044 ! Get a pointer to the list
2045 temp_list_p => find_list(name, current_list_p, .false.)
2046 
2047 if (associated(temp_list_p)) then
2048 ! restore the saved root values if we've already changed root
2049  if (save_root_name .ne. ' ') then
2050  root_p%name = save_root_name
2051  root_p%parent => save_root_parent_p
2052  endif
2053 ! set the pointer for the new root field
2054  root_p => temp_list_p
2055 ! save the new root field's name and parent
2056  save_root_name = root_p%name
2057  save_root_parent_p => root_p%parent
2058 ! set the new root name and parent fields to appropriate values
2059  root_p%name = ' '
2060  nullify(root_p%parent)
2061 ! set the current list to the new root as it likely is not
2062 ! going to be meaningful anymore
2063  current_list_p => root_p
2064  success = .true.
2065 else
2066 ! Couldn't find the list
2067  success = .false.
2068 endif
2069 
2070 end function fm_change_root
2071 
2072 !> @brief A function to list properties associated with a field.
2073 !!
2074 !> This function writes the contents of the field named "name" to stdout.
2075 !! If recursive is present and .true., then this function writes out the
2076 !! contents of any subfields associated with the field named "name".
2077 !! @return A flag to indicate operation success, true = no errors
2078 logical function fm_dump_list(name, recursive, unit) result (success)
2079  character(len=*), intent(in) :: name !< The name of the field for which output is requested.
2080  logical, intent(in), optional :: recursive !< If present and .true., then a recursive listing of
2081  !! fields will be performed.
2082  integer, intent(in), optional :: unit !< file to print to
2083 
2084  logical :: recursive_t
2085  type (field_def), pointer, save :: temp_list_p
2086  integer :: out_unit
2087 
2088  if (present(unit)) then
2089  out_unit = unit
2090  else
2091  out_unit = stdout()
2092  endif
2093 
2094  recursive_t = .false.
2095  if (present(recursive)) recursive_t = recursive
2096  if (.not. module_is_initialized) call initialize_module_variables()
2097 
2098  if (name .eq. ' ') then
2099  ! If list is empty, then dump the current list
2100  temp_list_p => current_list_p
2101  success = .true.
2102  else
2103  ! Get a pointer to the list
2104  temp_list_p => find_list(name, current_list_p, .false.)
2105  if (associated(temp_list_p)) then
2106  success = .true.
2107  else
2108  success = .false.
2109  endif
2110  endif
2111  ! Dump the list
2112  if (success) then
2113  success = dump_list(temp_list_p, recursive_t, 0, out_unit)
2114  endif
2115 end function fm_dump_list
2116 
2117 !> @brief A function to test whether a named field exists.
2118 !!
2119 !> This function determines is a field exists, relative to the current list,
2120 !! and returns true if the list exists, false otherwise.
2121 !! @return A flag to indicate operation success, true = no errors
2122 function fm_exists(name) &
2123  result(success)
2124 logical :: success
2125 character(len=*), intent(in) :: name !< The name of the field that is being queried
2126 
2127 type (field_def), pointer, save :: dummy_p
2128 ! Initialize the field manager if needed
2129 if (.not. module_is_initialized) then
2131 endif
2132 ! Determine whether the field exists
2133 dummy_p => get_field(name, current_list_p)
2134 success = associated(dummy_p)
2135 
2136 end function fm_exists
2137 
2138 !> @brief A function to return the index of a named field.
2139 !!
2140 !> Returns the index for name, returns the parameter NO_FIELD if it does not
2141 !! exist. If the first character of the named field is the list peparator,
2142 !! then the named field will be relative to the root of the field tree.
2143 !! Otherwise the named field will be relative to the current list.
2144 !> @returns index of the named field if it exists, otherwise the parameter NO_FIELD
2145 function fm_get_index(name) &
2146  result(index)
2147 integer :: index
2148 character(len=*), intent(in) :: name !< The name of a field that the user wishes to get an index for
2149 
2150 type (field_def), pointer, save :: temp_field_p
2151 integer :: out_unit
2152 
2153 out_unit = stdout()
2154 ! Initialize the field manager if needed
2155 if (.not. module_is_initialized) then
2157 endif
2158 ! Must supply a field field name
2159 if (name .eq. ' ') then
2160  index = no_field
2161  return
2162 endif
2163 ! Get a pointer to the field
2164 temp_field_p => get_field(name, current_list_p)
2165 if (associated(temp_field_p)) then
2166 ! Set the index
2167  index = temp_field_p%index
2168 else
2169  index = no_field
2170 endif
2171 
2172 end function fm_get_index
2173 
2174 !> @brief A function to return the full path of the current list.
2175 !!
2176 !> This function returns the full path for the current list. A blank
2177 !! path indicates an error condition has occurred.
2178 !> @returns The path corresponding to the current list
2180  result(path)
2181 character(len=FMS_PATH_LEN) :: path
2182 
2183 type (field_def), pointer, save :: temp_list_p
2184 ! Initialize the field manager if needed
2185 if (.not. module_is_initialized) then
2187 endif
2188 ! Set a pointer to the current list and proceed
2189 ! up the tree, filling in the name as we go
2190 temp_list_p => current_list_p
2191 path = ' '
2192 
2193 do while (associated(temp_list_p))
2194 ! Check whether we are at the root field--it is the
2195 ! only field with a blank name
2196  if (temp_list_p%name .eq. ' ') then
2197  exit
2198  endif
2199 ! Append the name to the path
2200  path = list_sep // trim(temp_list_p%name) // path
2201 ! Point to the next field
2202  temp_list_p => temp_list_p%parent
2203 enddo
2204 
2205 if (.not. associated(temp_list_p)) then
2206 ! The pointer is not associated, indicating an error has
2207 ! occurred, so set the path accordingly
2208  path = ' '
2209 elseif (path .eq. ' ') then
2210 ! If path is empty, then the current list must be root,
2211 ! so set path accordingly
2212  path = list_sep
2213 endif
2214 
2215 end function fm_get_current_list
2216 
2217 !> @brief A function to return how many elements are contained within the named
2218 !! list or entry.
2219 !!
2220 !> This function returns the list or entry length for the named list or entry.
2221 !! If the named field or entry does not exist, a value of 0 is returned.
2222 !> @returns The number of elements that the field name has.
2223 function fm_get_length(name) &
2224  result(length)
2225 integer :: length
2226 character(len=*), intent(in) :: name !< The name of a list or entry that the user wishes to get the length of
2227 
2228 type (field_def), pointer, save :: temp_field_p
2229 integer :: out_unit
2230 
2231 out_unit = stdout()
2232 ! Initialize the field manager if needed
2233 if (.not. module_is_initialized) then
2235 endif
2236 ! Must supply a field name
2237 if (name .eq. ' ') then
2238  length = 0
2239  return
2240 endif
2241 ! Get a pointer to the field
2242 temp_field_p => get_field(name, current_list_p)
2243 
2244 if (associated(temp_field_p)) then
2245 ! Set the field length
2246  if (temp_field_p%field_type .eq. list_type) then
2247  length = temp_field_p%length
2248  else
2249  length = temp_field_p%max_index
2250  endif
2251 else
2252  length = 0
2253 endif
2254 
2255 end function fm_get_length
2256 
2257 !> @brief A function to return the type of the named field
2258 !!
2259 !> This function returns the type of the field for name.
2260 !! This indicates whether the named field is a "list" (has children fields),
2261 !! or has values of type "integer", "real", "logical" or "string".
2262 !! If it does not exist it returns a blank string.
2263 !> @returns A string containing the type of the named field
2264 function fm_get_type(name) &
2265  result(name_field_type)
2266 character(len=8) :: name_field_type
2267 character(len=*), intent(in) :: name !< The name of a field that the user wishes to find the type of
2268 
2269 type (field_def), pointer, save :: temp_field_p
2270 integer :: out_unit
2271 
2272 out_unit = stdout()
2273 ! Initialize the field manager if needed
2274 if (.not. module_is_initialized) then
2276 endif
2277 ! Must supply a field name
2278 if (name .eq. ' ') then
2279  name_field_type = ' '
2280  return
2281 endif
2282 ! Get a pointer to the field
2283 temp_field_p => get_field(name, current_list_p)
2284 
2285 if (associated(temp_field_p)) then
2286 ! Set the field type
2287  name_field_type = field_type_name(temp_field_p%field_type)
2288 else
2289  name_field_type = ' '
2290 endif
2291 
2292 end function fm_get_type
2293 
2294 !> @returns A flag to indicate whether the function operated with (false) or without
2295 !! (true) errors.
2296 function fm_get_value_integer(name, get_ival, index) &
2297  result(success)
2298 logical :: success
2299 character(len=*), intent(in) :: name !< The name of a field that the user wishes to get a value for.
2300 integer, intent(out) :: get_ival !< The value associated with the named field.
2301 integer, intent(in), optional :: index !< An optional index to retrieve a single value from an array.
2302 
2303 integer :: index_t
2304 type (field_def), pointer, save :: temp_field_p
2305 integer :: out_unit
2306 
2307 out_unit = stdout()
2308 ! Initialize the field manager if needed
2309 if (.not. module_is_initialized) then
2311 endif
2312 ! Must supply a field field name
2313 if (name .eq. ' ') then
2314  get_ival = 0
2315  success = .false.
2316  return
2317 endif
2318 ! Set index to retrieve
2319 if (present(index)) then
2320  index_t = index
2321 else
2322  index_t = 1
2323 endif
2324 ! Get a pointer to the field
2325 temp_field_p => get_field(name, current_list_p)
2326 
2327 if (associated(temp_field_p)) then
2328 ! check that the field is the correct type
2329  if (temp_field_p%field_type .eq. integer_type) then
2330  if (index_t .lt. 1 .or. index_t .gt. temp_field_p%max_index) then
2331 ! Index is not positive or index too large
2332  get_ival = 0
2333  success = .false.
2334  else
2335 ! extract the value
2336  get_ival = temp_field_p%i_value(index_t)
2337  success = .true.
2338  endif
2339  else
2340 ! Field not corrcet type
2341  get_ival = 0
2342  success = .false.
2343  endif
2344 else
2345  get_ival = 0
2346  success = .false.
2347 endif
2348 
2349 end function fm_get_value_integer
2350 
2351 !> @returns A flag to indicate whether the function operated with (false) or without
2352 !! (true) errors.
2353 function fm_get_value_logical(name, get_lval, index) &
2354  result(success)
2355 logical :: success
2356 character(len=*), intent(in) :: name !< The name of a field that the user wishes to get a value for.
2357 logical, intent(out) :: get_lval !< The value associated with the named field
2358 integer, intent(in), optional :: index !< An optional index to retrieve a single value from an array.
2359 
2360 integer :: index_t
2361 type (field_def), pointer, save :: temp_field_p
2362 integer :: out_unit
2363 
2364 out_unit = stdout()
2365 ! Initialize the field manager if needed
2366 if (.not. module_is_initialized) then
2368 endif
2369 ! Must supply a field field name
2370 if (name .eq. ' ') then
2371  get_lval = .false.
2372  success = .false.
2373  return
2374 endif
2375 ! Set index to retrieve
2376 if (present(index)) then
2377  index_t = index
2378 else
2379  index_t = 1
2380 endif
2381 ! Get a pointer to the field
2382 temp_field_p => get_field(name, current_list_p)
2383 
2384 if (associated(temp_field_p)) then
2385 ! check that the field is the correct type
2386  if (temp_field_p%field_type .eq. logical_type) then
2387 
2388  if (index_t .lt. 1 .or. index_t .gt. temp_field_p%max_index) then
2389 ! Index is not positive or too large
2390  get_lval = .false.
2391  success = .false.
2392  else
2393 ! extract the value
2394  get_lval = temp_field_p%l_value(index_t)
2395  success = .true.
2396  endif
2397  else
2398 ! Field not correct type
2399  get_lval = .false.
2400  success = .false.
2401  endif
2402 else
2403  get_lval = .false.
2404  success = .false.
2405 endif
2406 
2407 end function fm_get_value_logical
2408 
2409 !> @returns A flag to indicate whether the function operated with (false) or without
2410 !! (true) errors.
2411 function fm_get_value_string(name, get_sval, index) &
2412  result(success)
2413 logical :: success
2414 character(len=*), intent(in) :: name !< The name of a field that the user wishes to get a value for.
2415 character(len=*), intent(out) :: get_sval !< The value associated with the named field
2416 integer, intent(in), optional :: index !< An optional index to retrieve a single value from an array.
2417 
2418 integer :: index_t
2419 type (field_def), pointer, save :: temp_field_p
2420 integer :: out_unit
2421 
2422 out_unit = stdout()
2423 ! Initialize the field manager if needed
2424 if (.not. module_is_initialized) then
2426 endif
2427 ! Must supply a field field name
2428 if (name .eq. ' ') then
2429  get_sval = ''
2430  success = .false.
2431  return
2432 endif
2433 ! Set index to retrieve
2434 if (present(index)) then
2435  index_t = index
2436 else
2437  index_t = 1
2438 endif
2439 ! Get a pointer to the field
2440 temp_field_p => get_field(name, current_list_p)
2441 
2442 if (associated(temp_field_p)) then
2443 ! check that the field is the correct type
2444  if (temp_field_p%field_type .eq. string_type) then
2445  if (index_t .lt. 1 .or. index_t .gt. temp_field_p%max_index) then
2446 ! Index is not positive or is too large
2447  get_sval = ''
2448  success = .false.
2449  else
2450 ! extract the value
2451  get_sval = temp_field_p%s_value(index_t)
2452  success = .true.
2453  endif
2454  else
2455 ! Field not correct type
2456  get_sval = ''
2457  success = .false.
2458  endif
2459 else
2460  get_sval = ''
2461  success = .false.
2462 endif
2463 
2464 end function fm_get_value_string
2465 
2466 !> Iterates through the given list
2467 !> @returns A flag to indicate whether the function operated with (FALSE)
2468 !! or without (TRUE) errors
2469 function fm_loop_over_list_old(list, name, field_type, index) &
2470  result(success)
2471 logical :: success
2472 character(len=*), intent(in) :: list !< Name of a list to loop over
2473 character(len=*), intent(out) :: name !< name of a field from list
2474 character(len=fm_type_name_len), intent(out) :: field_type !< type of a list entry
2475 integer, intent(out) :: index !< index of the field within the list
2476 
2477 integer :: out_unit
2478 
2479 out_unit = stdout()
2480 ! Initialize the field manager if needed
2481 if (.not. module_is_initialized) then
2483 endif
2484 
2485 if (list .eq. loop_list .and. associated(loop_list_p)) then
2486 ! We've already started this loop, so continue on
2487  loop_list_p => loop_list_p%next
2488  success = set_list_stuff()
2489 elseif (list .eq. ' ') then
2490 ! If list is empty, then loop over the current list
2491  loop_list = ' '
2492  loop_list_p => current_list_p%first_field
2493  success = set_list_stuff()
2494 else
2495 ! Get a pointer to the list
2496  loop_list = list
2497  loop_list_p => find_list(loop_list, current_list_p, .false.)
2498  if (associated(loop_list_p)) then
2499  loop_list_p => loop_list_p%first_field
2500  success = set_list_stuff()
2501  else
2502  success = .false.
2503  endif
2504 endif
2505 
2506 return
2507 
2508 contains
2509 
2510 !> If the pointer matches to the right list,
2511 !! extract the field information. Used in fm_loop_over_list
2512 !> @returns A flag to indicate whether the function operated with (FALSE)
2513 !! or without (TRUE) errors
2514 function set_list_stuff() &
2515  result(success)
2516 
2517  logical :: success
2518 
2519  if (associated(loop_list_p)) then
2520  name = loop_list_p%name
2521  field_type = field_type_name(loop_list_p%field_type)
2522  index = loop_list_p%index
2523  success = .true.
2524  else
2525  name = ' '
2526  field_type = ' '
2527  index = 0
2528  success = .false.
2529  loop_list = ' '
2530  endif
2531 
2532 end function set_list_stuff
2533 
2534 end function fm_loop_over_list_old
2535 
2536 !> given a name of the list, prepares an iterator over the list content.
2537 !! If the name of the given list is blank, then the current list is used
2538 subroutine fm_init_loop(loop_list, iter)
2539  character(len=*) , intent(in) :: loop_list !< name of the list to iterate over
2540  type(fm_list_iter_type), intent(out) :: iter !< loop iterator
2541 
2542  if (.not.module_is_initialized) call initialize_module_variables
2543 
2544  if (loop_list==' ') then ! looping over current list
2545  iter%ptr => current_list_p%first_field
2546  else
2547  iter%ptr => find_list(loop_list,current_list_p,.false.)
2548  if (associated(iter%ptr)) iter%ptr => iter%ptr%first_field
2549  endif
2550 end subroutine fm_init_loop
2551 
2552 !> given a list iterator, returns information about curren list element
2553 !! and advances the iterator to the next list element. At the end of the
2554 !! list, returns FALSE
2555 function fm_loop_over_list_new(iter, name, field_type, index) &
2556  result(success) ; logical success
2557  type (fm_list_iter_type), intent(inout) :: iter !< list iterator
2558  character(len=*), intent(out) :: name !< name of the current list item
2559  character(len=*), intent(out) :: field_type !< type of the field
2560  integer , intent(out) :: index !< index in the list
2561 
2562  if (.not.module_is_initialized) call initialize_module_variables
2563  if (associated(iter%ptr)) then
2564  name = iter%ptr%name
2565  field_type = field_type_name(iter%ptr%field_type)
2566  index = iter%ptr%index
2567  success = .true.
2568  iter%ptr => iter%ptr%next
2569  else
2570  name = ' '
2571  field_type = ' '
2572  index = 0
2573  success = .false.
2574  endif
2575 end function fm_loop_over_list_new
2576 
2577 !> @brief A function to create a new list
2578 !!
2579 !> Allocate and initialize a new list and return the index of the list. If an
2580 !! error occurs return the parameter NO_FIELD.
2581 !> @return integer index of the newly created list
2582 function fm_new_list(name, create, keep) &
2583  result(index)
2584 integer :: index
2585 character(len=*), intent(in) :: name !< Name of a list that user wishes to create
2586 logical, intent(in), optional :: create !< If present and true, create the list if it does not exist
2587 logical, intent(in), optional :: keep !< If present and true, make this list the current list
2588 
2589 logical :: create_t
2590 logical :: keep_t
2591 character(len=FMS_PATH_LEN) :: path
2592 character(len=fm_field_name_len) :: base
2593 type (field_def), pointer, save :: temp_list_p
2594 integer :: out_unit
2595 
2596 out_unit = stdout()
2597 ! Initialize the field manager if needed
2598 if (.not. module_is_initialized) then
2600 endif
2601 ! Must supply a field list name
2602 if (name .eq. ' ') then
2603  index = no_field
2604  return
2605 endif
2606 ! Check for optional arguments
2607 if (present(create)) then
2608  create_t = create
2609 else
2610  create_t = .false.
2611 endif
2612 
2613 if (present(keep)) then
2614  keep_t = keep
2615 else
2616  keep_t = .false.
2617 endif
2618 ! Get a pointer to the parent list
2619 call find_base(name, path, base)
2620 
2621 temp_list_p => find_list(path, current_list_p, create_t)
2622 
2623 if (associated(temp_list_p)) then
2624 ! Create the list
2625  temp_list_p => make_list(temp_list_p, base)
2626  if (associated(temp_list_p)) then
2627 ! Make this list the current list, if requested
2628  if (keep_t) then
2629  current_list_p => temp_list_p
2630  endif
2631  index = temp_list_p%index
2632  else
2633  index = no_field
2634  endif
2635 else
2636  index = no_field
2637 endif
2638 
2639 end function fm_new_list
2640 
2641 !> @brief Assigns a given value to a given field
2642 !> @returns An index for the named field
2643 function fm_new_value_integer(name, new_ival, create, index, append) &
2644  result(field_index)
2645 integer :: field_index
2646 character(len=*), intent(in) :: name !< The name of a field that the user wishes to create
2647  !! a value for.
2648 integer, intent(in) :: new_ival !< The value that the user wishes to apply to the
2649  !! named field.
2650 logical, intent(in), optional :: create !< If present and .true., then a value for this
2651  !! field will be created.
2652 integer, intent(in), optional :: index !< The index to an array of values that the user
2653  !! wishes to apply a new value.
2654 logical, intent(in), optional :: append !< If present and .true., then append the value to
2655  !! an array of the present values. If present and .true., then index cannot be greater than 0.
2656 
2657 logical :: create_t
2658 integer :: i
2659 integer :: index_t
2660 integer, pointer, dimension(:) :: temp_i_value
2661 character(len=FMS_PATH_LEN) :: path
2662 character(len=fm_field_name_len) :: base
2663 type (field_def), pointer, save :: temp_list_p
2664 type (field_def), pointer, save :: temp_field_p
2665 integer :: out_unit
2666 
2667 out_unit = stdout()
2668 ! Initialize the field manager if needed
2669 if (.not. module_is_initialized) then
2671 endif
2672 ! Must supply a field name
2673 if (name .eq. ' ') then
2674  field_index = no_field
2675  return
2676 endif
2677 ! Check for optional arguments
2678 if (present(create)) then
2679  create_t = create
2680 else
2681  create_t = .false.
2682 endif
2683 ! Check that append is not true and index non-positive
2684 if (present(index) .and. present(append)) then
2685  if (append .and. index .gt. 0) then
2686  field_index = no_field
2687  return
2688  endif
2689 endif
2690 ! Set index to define
2691 if (present(index)) then
2692  index_t = index
2693  if (index_t .lt. 0) then
2694 ! Index is negative
2695  field_index = no_field
2696  return
2697  endif
2698 else
2699  index_t = 1
2700 endif
2701 ! Get a pointer to the parent list
2702 call find_base(name, path, base)
2703 temp_list_p => find_list(path, current_list_p, create_t)
2704 
2705 if (associated(temp_list_p)) then
2706  temp_field_p => find_field(base, temp_list_p)
2707  if (.not. associated(temp_field_p)) then
2708 ! Create the field if it doesn't exist
2709  temp_field_p => create_field(temp_list_p, base)
2710  endif
2711  if (associated(temp_field_p)) then
2712 ! Check if the field_type is the same as previously
2713 ! If not then reset max_index to 0
2714  if (temp_field_p%field_type == real_type ) then
2715  ! promote integer input to real
2716  ! all real field values are stored as r8_kind
2717  field_index = fm_new_value(name, real(new_ival,r8_kind), create, index, append)
2718  return
2719  else if (temp_field_p%field_type /= integer_type ) then
2720  ! slm: why would we reset index? Is it not an error to have a "list" defined
2721  ! with different types in more than one place?
2722  temp_field_p%max_index = 0
2723  endif
2724 ! Assign the type
2725  temp_field_p%field_type = integer_type
2726 ! Set the index if appending
2727  if (present(append)) then
2728  if (append) then
2729  index_t = temp_field_p%max_index + 1
2730  endif
2731  endif
2732 
2733  if (index_t .gt. temp_field_p%max_index + 1) then
2734 ! Index too large
2735  field_index = no_field
2736  return
2737 
2738  elseif (index_t .eq. 0 .and. &
2739  temp_field_p%max_index .gt. 0) then
2740 ! Can't set non-null field to null
2741  field_index = no_field
2742  return
2743 
2744  elseif (.not. allocated(temp_field_p%i_value) .and. &
2745  index_t .gt. 0) then
2746 ! Array undefined, so allocate the array
2747  allocate(temp_field_p%i_value(1))
2748  temp_field_p%max_index = 1
2749  temp_field_p%array_dim = 1
2750  elseif (index_t .gt. temp_field_p%array_dim) then
2751 ! Array is too small, so allocate new array and copy over
2752 ! old values
2753  temp_field_p%array_dim = temp_field_p%array_dim + array_increment
2754  allocate (temp_i_value(temp_field_p%array_dim))
2755  do i = 1, temp_field_p%max_index
2756  temp_i_value(i) = temp_field_p%i_value(i)
2757  enddo
2758  if (allocated(temp_field_p%i_value)) deallocate(temp_field_p%i_value)
2759  temp_field_p%i_value = temp_i_value
2760  temp_field_p%max_index = index_t
2761  endif
2762 ! Assign the value and set the field_index for return
2763 ! for non-null fields (index_t > 0)
2764  if (index_t .gt. 0) then
2765  temp_field_p%i_value(index_t) = new_ival
2766  if (index_t .gt. temp_field_p%max_index) then
2767  temp_field_p%max_index = index_t
2768  endif
2769  endif
2770  field_index = temp_field_p%index
2771 
2772  else
2773  field_index = no_field
2774  endif
2775 else
2776  field_index = no_field
2777 endif
2778 
2779 end function fm_new_value_integer
2780 
2781 !> @brief Assigns a given value to a given field
2782 !> @returns An index for the named field
2783 function fm_new_value_logical(name, new_lval, create, index, append) &
2784  result(field_index)
2785 integer :: field_index
2786 character(len=*), intent(in) :: name !< The name of a field that the user wishes to create
2787  !! a value for.
2788 logical, intent(in) :: new_lval !< The value that the user wishes to apply to the
2789  !! named field.
2790 logical, intent(in), optional :: create !< If present and .true., then a value for this
2791  !! field will be created.
2792 integer, intent(in), optional :: index !< The index to an array of values that the user
2793  !! wishes to apply a new value.
2794 logical, intent(in), optional :: append !< If present and .true., then append the value to
2795  !! an array of the present values. If present and .true., then index cannot be greater than 0.
2796 
2797 character(len=FMS_PATH_LEN) :: path
2798 character(len=fm_field_name_len) :: base
2799 integer :: i
2800 integer :: index_t
2801 logical :: create_t
2802 logical, dimension(:), pointer :: temp_l_value
2803 type (field_def), pointer, save :: temp_list_p
2804 type (field_def), pointer, save :: temp_field_p
2805 integer :: out_unit
2806 
2807 out_unit = stdout()
2808 ! Initialize the field manager if needed
2809 if (.not. module_is_initialized) then
2811 endif
2812 ! Must supply a field name
2813 if (name .eq. ' ') then
2814  field_index = no_field
2815  return
2816 endif
2817 ! Check for optional arguments
2818 if (present(create)) then
2819  create_t = create
2820 else
2821  create_t = .false.
2822 endif
2823 ! Check that append is not true and index greater than 0
2824 if (present(index) .and. present(append)) then
2825  if (append .and. index .gt. 0) then
2826  field_index = no_field
2827  return
2828  endif
2829 endif
2830 ! Set index to define
2831 if (present(index)) then
2832  index_t = index
2833  if (index_t .lt. 0) then
2834 ! Index is negative
2835  field_index = no_field
2836  return
2837  endif
2838 else
2839  index_t = 1
2840 endif
2841 ! Get a pointer to the parent list
2842 call find_base(name, path, base)
2843 temp_list_p => find_list(path, current_list_p, create_t)
2844 
2845 if (associated(temp_list_p)) then
2846  temp_field_p => find_field(base, temp_list_p)
2847  if (.not. associated(temp_field_p)) then
2848 ! Create the field if it doesn't exist
2849  temp_field_p => create_field(temp_list_p, base)
2850  endif
2851  if (associated(temp_field_p)) then
2852 ! Check if the field_type is the same as previously
2853 ! If not then reset max_index to 0
2854  if (temp_field_p%field_type /= logical_type ) then
2855  temp_field_p%max_index = 0
2856  endif
2857 ! Assign the type
2858  temp_field_p%field_type = logical_type
2859 ! Set the index if appending
2860  if (present(append)) then
2861  if (append) then
2862  index_t = temp_field_p%max_index + 1
2863  endif
2864  endif
2865 
2866  if (index_t .gt. temp_field_p%max_index + 1) then
2867 ! Index too large
2868  field_index = no_field
2869  return
2870 
2871  elseif (index_t .eq. 0 .and. &
2872  temp_field_p%max_index .gt. 0) then
2873 ! Can't set non-null field to null
2874  field_index = no_field
2875  return
2876 
2877  elseif (.not. allocated(temp_field_p%l_value) .and. &
2878  index_t .gt. 0) then
2879 ! Array undefined, so allocate the array
2880  allocate(temp_field_p%l_value(1))
2881  temp_field_p%max_index = 1
2882  temp_field_p%array_dim = 1
2883 
2884  elseif (index_t .gt. temp_field_p%array_dim) then
2885 ! Array is too small, so allocate new array and copy over
2886 ! old values
2887  temp_field_p%array_dim = temp_field_p%array_dim + array_increment
2888  allocate (temp_l_value(temp_field_p%array_dim))
2889  do i = 1, temp_field_p%max_index
2890  temp_l_value(i) = temp_field_p%l_value(i)
2891  enddo
2892  if (allocated(temp_field_p%l_value)) deallocate(temp_field_p%l_value)
2893  temp_field_p%l_value = temp_l_value
2894  temp_field_p%max_index = index_t
2895 
2896  endif
2897 ! Assign the value and set the field_index for return
2898 ! for non-null fields (index_t > 0)
2899  if (index_t .gt. 0) then
2900  temp_field_p%l_value(index_t) = new_lval
2901  if (index_t .gt. temp_field_p%max_index) then
2902  temp_field_p%max_index = index_t
2903  endif
2904  endif
2905  field_index = temp_field_p%index
2906  else
2907  field_index = no_field
2908  endif
2909 else
2910  field_index = no_field
2911 endif
2912 
2913 end function fm_new_value_logical
2914 
2915 !> @brief Assigns a given value to a given field
2916 !> @returns An index for the named field
2917 function fm_new_value_string(name, new_sval, create, index, append) &
2918  result(field_index)
2919 integer :: field_index
2920 character(len=*), intent(in) :: name !< The name of a field that the user wishes to create
2921  !! a value for.
2922 character(len=*), intent(in) :: new_sval !< The value that the user wishes to apply to the
2923  !! named field.
2924 logical, intent(in), optional :: create !< If present and .true., then a value for this
2925  !! field will be created.
2926 integer, intent(in), optional :: index !< The index to an array of values that the user
2927  !! wishes to apply a new value.
2928 logical, intent(in), optional :: append !< If present and .true., then append the value to
2929 
2930 character(len=fm_string_len), dimension(:), pointer :: temp_s_value
2931 character(len=FMS_PATH_LEN) :: path
2932 character(len=fm_field_name_len) :: base
2933 integer :: i
2934 integer :: index_t
2935 logical :: create_t
2936 type (field_def), save, pointer :: temp_list_p
2937 type (field_def), save, pointer :: temp_field_p
2938 integer :: out_unit
2939 
2940 out_unit = stdout()
2941 ! Initialize the field manager if needed
2942 if (.not. module_is_initialized) then
2944 endif
2945 ! Must supply a field name
2946 if (name .eq. ' ') then
2947  field_index = no_field
2948  return
2949 endif
2950 ! Check for optional arguments
2951 if (present(create)) then
2952  create_t = create
2953 else
2954  create_t = .false.
2955 endif
2956 ! Check that append is not true and index greater than 0
2957 if (present(index) .and. present(append)) then
2958  if (append .and. index .gt. 0) then
2959  field_index = no_field
2960  return
2961  endif
2962 endif
2963 ! Set index to define
2964 if (present(index)) then
2965  index_t = index
2966  if (index_t .lt. 0) then
2967 ! Index is negative
2968  field_index = no_field
2969  return
2970  endif
2971 else
2972  index_t = 1
2973 endif
2974 ! Get a pointer to the parent list
2975 call find_base(name, path, base)
2976 temp_list_p => find_list(path, current_list_p, create_t)
2977 
2978 if (associated(temp_list_p)) then
2979  temp_field_p => find_field(base, temp_list_p)
2980  if (.not. associated(temp_field_p)) then
2981 ! Create the field if it doesn't exist
2982  temp_field_p => create_field(temp_list_p, base)
2983  endif
2984  if (associated(temp_field_p)) then
2985 ! Check if the field_type is the same as previously
2986 ! If not then reset max_index to 0
2987  if (temp_field_p%field_type /= string_type ) then
2988  temp_field_p%max_index = 0
2989  endif
2990 ! Assign the type
2991  temp_field_p%field_type = string_type
2992 ! Set the index if appending
2993  if (present(append)) then
2994  if (append) then
2995  index_t = temp_field_p%max_index + 1
2996  endif
2997  endif
2998 
2999  if (index_t .gt. temp_field_p%max_index + 1) then
3000 ! Index too large
3001  field_index = no_field
3002  return
3003 
3004  elseif (index_t .eq. 0 .and. &
3005  temp_field_p%max_index .gt. 0) then
3006 ! Can't set non-null field to null
3007  field_index = no_field
3008  return
3009 
3010  elseif (.not.allocated(temp_field_p%s_value) .and. &
3011  index_t .gt. 0) then
3012 ! Array undefined, so allocate the array
3013  allocate(temp_field_p%s_value(1))
3014  temp_field_p%max_index = 1
3015  temp_field_p%array_dim = 1
3016 
3017  elseif (index_t .gt. temp_field_p%array_dim) then
3018 ! Array is too small, so allocate new array and copy over
3019 ! old values
3020  temp_field_p%array_dim = temp_field_p%array_dim + array_increment
3021  allocate (temp_s_value(temp_field_p%array_dim))
3022  do i = 1, temp_field_p%max_index
3023  temp_s_value(i) = temp_field_p%s_value(i)
3024  enddo
3025  if (allocated(temp_field_p%s_value)) deallocate(temp_field_p%s_value)
3026  temp_field_p%s_value = temp_s_value
3027  temp_field_p%max_index = index_t
3028 
3029  endif
3030 ! Assign the value and set the field_index for return
3031 ! for non-null fields (index_t > 0)
3032  if (index_t .gt. 0) then
3033  temp_field_p%s_value(index_t) = new_sval
3034  if (index_t .gt. temp_field_p%max_index) then
3035  temp_field_p%max_index = index_t
3036  endif
3037  endif
3038  field_index = temp_field_p%index
3039  else
3040 ! Error in making the field
3041  field_index = no_field
3042  endif
3043 else
3044 ! Error following the path
3045  field_index = no_field
3046 endif
3047 
3048 end function fm_new_value_string
3049 
3050 
3051 !> Resets the loop variable. For use in conjunction with fm_loop_over_list.
3052 subroutine fm_reset_loop
3053 ! Initialize the field manager if needed
3054 if (.not. module_is_initialized) then
3056 endif
3057 ! Reset the variables
3058 loop_list = ' '
3059 nullify(loop_list_p)
3060 
3061 end subroutine fm_reset_loop
3062 
3063 !> Return the root list to the value at initialization.
3064 !!
3065 !> For use in conjunction with fm_change_root.
3066 !!
3067 !! Users should use this routine before leaving their routine if they
3068 !! previously used fm_change_root.
3069 subroutine fm_return_root
3070 ! Initialize the field manager if needed
3071 if (.not. module_is_initialized) then
3073 endif
3074 ! restore the saved values to the current root
3075 root_p%name = save_root_name
3076 root_p%parent => save_root_parent_p
3077 ! set the pointer to the original root field
3078 root_p => root
3079 ! reset the save root name and parent variables
3080 save_root_name = ' '
3081 nullify(save_root_parent_p)
3082 
3083 end subroutine fm_return_root
3084 
3085 !> Return a pointer to the field if it exists relative to this_list_p,
3086 !! null otherwise
3087 !! @returns A pointer to the field name
3088 function get_field(name, this_list_p) &
3089  result(list_p)
3090 type (field_def), pointer :: list_p
3091 character(len=*), intent(in) :: name !< The name of a list that the user wishes to get information for
3092 type (field_def), pointer :: this_list_p !< A pointer to a list that serves as the base point
3093  !! for searching for name
3094 
3095 character(len=FMS_PATH_LEN) :: path
3096 character(len=fm_field_name_len) :: base
3097 type (field_def), pointer, save :: temp_p
3098 
3099 nullify(list_p)
3100 ! Get the path and base for name
3101 call find_base(name, path, base)
3102 ! Find the list if path is not empty
3103 if (path .ne. ' ') then
3104  temp_p => find_list(path, this_list_p, .false.)
3105  if (associated(temp_p)) then
3106  list_p => find_field(base, temp_p)
3107  else
3108  nullify(list_p)
3109  endif
3110 else
3111  list_p => find_field(base, this_list_p)
3112 endif
3113 
3114 end function get_field
3115 
3116 !> This function allows a user to rename a field without modifying the
3117 !! contents of the field.
3118 !!
3119 !> Function to modify the name of a field.
3120 !! Should be used with caution.
3121 !! @returns A flag to indicate whether the function operated with (FALSE) or
3122 !! without (TRUE) errors.
3123 function fm_modify_name(oldname, newname) &
3124  result(success)
3125 logical :: success
3126 character(len=*), intent(in) :: oldname !< The name of a field that the user wishes to change
3127  !! the name of
3128 character(len=*), intent(in) :: newname !< The name that the user wishes to change the name of
3129  !! the field to.
3130 
3131 character(len=FMS_PATH_LEN) :: path
3132 character(len=fm_field_name_len) :: base
3133 type (field_def), pointer, save :: list_p
3134 type (field_def), pointer, save :: temp_p
3135 ! Get the path and base for name
3136 call find_base(oldname, path, base)
3137 ! Find the list if path is not empty
3138 success = .false.
3139 if (path .ne. ' ') then
3140  temp_p => find_list(path, current_list_p, .false.)
3141  if (associated(temp_p)) then
3142  list_p => find_field(base, temp_p)
3143  if (associated(list_p)) then
3144  list_p%name = newname
3145  success = .true.
3146  endif
3147  else
3148  nullify(list_p)
3149  endif
3150 else
3151  list_p => find_field(base, current_list_p)
3152  if (associated(list_p)) then
3153  list_p%name = newname
3154  success = .true.
3155  endif
3156 endif
3157 
3158 end function fm_modify_name
3159 
3160 !> A function to initialize the values of the pointers. This will remove
3161 !! all fields and reset the field tree to only the root field.
3163  ! Initialize the root field
3164  integer :: io, ierr !< Error codes when reading the namelist
3165  integer :: logunit !< Unit number for the log file
3166 
3167  if (.not. module_is_initialized) then
3168 
3169  read (input_nml_file, nml=field_manager_nml, iostat=io)
3170  ierr = check_nml_error(io,"field_manager_nml")
3171 
3172  logunit = stdlog()
3173  if (mpp_pe() == mpp_root_pe()) write (logunit, nml=field_manager_nml)
3174 
3175  root_p => root
3176 
3177  field_type_name(integer_type) = 'integer'
3178  field_type_name(list_type) = 'list'
3179  field_type_name(logical_type) = 'logical'
3180  field_type_name(real_type) = 'real'
3181  field_type_name(string_type) = 'string'
3182 
3183  root%name = ' '
3184  root%index = 1
3185  root%parent => root_p
3186 
3187  root%field_type = list_type
3188 
3189  root%length = 0
3190  nullify(root%first_field)
3191  nullify(root%last_field)
3192  root%max_index = 0
3193  root%array_dim = 0
3194  if (allocated(root%i_value)) deallocate(root%i_value)
3195  if (allocated(root%l_value)) deallocate(root%l_value)
3196  if (allocated(root%r_value)) deallocate(root%r_value)
3197  if (allocated(root%s_value)) deallocate(root%s_value)
3198 
3199  nullify(root%next)
3200  nullify(root%prev)
3201 
3202  current_list_p => root
3203 
3204  nullify(loop_list_p)
3205  loop_list = ' '
3206 
3207  nullify(save_root_parent_p)
3208  save_root_name = ' '
3209 
3210  module_is_initialized = .true.
3211 
3212 endif
3213 
3214 end subroutine initialize_module_variables
3215 
3216 !> This function creates a new field and returns a pointer to that field.
3217 !!
3218 !> Allocate and initialize a new list in this_list_p list.
3219 !! @return a pointer to the list on success, or a null pointer
3220 !! on failure.
3221 function make_list(this_list_p, name) &
3222  result(list_p)
3223 type (field_def), pointer :: list_p
3224 type (field_def), pointer :: this_list_p !< Base of a list that the user wishes to add a list to
3225 character(len=*), intent(in) :: name !< name of a list that the user wishes to create
3226 
3227 type (field_def), pointer, save :: dummy_p
3228 integer :: out_unit
3229 
3230 out_unit = stdout()
3231 ! Check to see whether there is already a list with
3232 ! this name, and if so, return an error as list names
3233 ! must be unique
3234 dummy_p => find_field(name, this_list_p )
3235 if (associated(dummy_p)) then
3236 ! This list is already specified, return an error
3237  list_p => dummy_p
3238  return
3239 endif
3240 ! Create a field for the new list
3241 nullify(list_p)
3242 list_p => create_field(this_list_p, name)
3243 if (.not. associated(list_p)) then
3244  nullify(list_p)
3245  return
3246 endif
3247 ! Initialize the new list
3248 list_p%length = 0
3249 list_p%field_type = list_type
3250 if (allocated(list_p%i_value)) deallocate(list_p%i_value)
3251 if (allocated(list_p%l_value)) deallocate(list_p%l_value)
3252 if (allocated(list_p%r_value)) deallocate(list_p%r_value)
3253 if (allocated(list_p%s_value)) deallocate(list_p%s_value)
3254 
3255 end function make_list
3256 
3257 !> This is a function that provides the capability to return parameters
3258 !! associated with a field in a pair of strings.
3259 !!
3260 !> Given a name return a list of method names and control strings.
3261 !! This function should return strings similar to those in the field
3262 !! table if a comma delimited format is being used.
3263 !> @return A flag to indicate whether the function operated with (FALSE) or
3264 !! without (TRUE) errors
3265 function fm_query_method(name, method_name, method_control) &
3266  result(success)
3267 logical :: success
3268 character(len=*), intent(in) :: name !< name of a list that the user wishes to change to
3269 character(len=*), intent(out) :: method_name !< name of a parameter associated with the named field
3270 character(len=*), intent(out) :: method_control !< value of parameters associated with the named field
3271 
3272 character(len=FMS_PATH_LEN) :: path
3273 character(len=FMS_PATH_LEN) :: base
3274 character(len=FMS_PATH_LEN) :: name_loc
3275 logical :: recursive_t
3276 type (field_def), pointer, save :: temp_list_p
3277 type (field_def), pointer, save :: temp_value_p
3278 type (field_def), pointer, save :: this_field_p
3279 integer :: out_unit
3280 
3281  out_unit = stdout()
3282  success = .false.
3283  recursive_t = .true.
3284  method_name = " "
3285  method_control = " "
3286 ! Initialize the field manager if needed
3287 if (.not. module_is_initialized) call initialize_module_variables
3288 name_loc = lowercase(name)
3289 call find_base(name_loc, path, base)
3290 
3291  temp_list_p => find_list(name_loc, current_list_p, .false.)
3292 
3293 if (associated(temp_list_p)) then
3294 ! Find the entry values for the list.
3295  success = query_method(temp_list_p, recursive_t, base, method_name, method_control)
3296 else
3297 ! This is not a list but it may be a parameter with a value
3298 ! If so put the parameter value in method_name.
3299 
3300  temp_value_p => find_list(path, current_list_p, .false.)
3301  if (associated(temp_value_p)) then
3302 ! Find the entry values for this item.
3303  this_field_p => temp_value_p%first_field
3304 
3305  do while (associated(this_field_p))
3306  if ( this_field_p%name == base ) then
3307  method_name = this_field_p%s_value(1)
3308  method_control = ""
3309  success = .true.
3310  exit
3311  else
3312  success = .false.
3313  endif
3314  this_field_p => this_field_p%next
3315  enddo
3316 
3317  else
3318 ! Error following the path
3319  success = .false.
3320  endif
3321 endif
3322 
3323 end function fm_query_method
3324 
3325 !> A private function that can recursively recover values for parameters
3326 !! associated with a field.
3327 !> @return A flag to indicate whether the function operated with (FALSE) or
3328 !! without (TRUE) errors
3329 recursive function query_method(list_p, recursive, name, method_name, method_control) &
3330  result(success)
3331 logical :: success
3332 type (field_def), pointer :: list_p !< A pointer to the field that is of interest
3333 logical, intent(in) :: recursive !< A flag to enable recursive searching if true
3334 character(len=*), intent(in) :: name !< name of a list that the user wishes to change to
3335 character(len=*), intent(out) :: method_name !< name of a parameter associated with the named field
3336 character(len=*), intent(out) :: method_control !< value of parameters associated with the named field
3337 
3338 integer :: i
3339 character(len=64) :: scratch
3340 type (field_def), pointer :: this_field_p
3341 integer :: out_unit
3342 
3343 out_unit = stdout()
3344 
3345 ! Check for a valid list
3346 if (.not. associated(list_p) .or. list_p%field_type .ne. list_type) then
3347  success = .false.
3348 else
3349 
3350  ! set the default return value
3351  success = .true.
3352 
3353  this_field_p => list_p%first_field
3354 
3355  do while (associated(this_field_p))
3356  select case(this_field_p%field_type)
3357  case(list_type)
3358  ! If this is a list, then this is the method name
3359  if (recursive) then
3360  if (.not. query_method(this_field_p, .true., this_field_p%name, method_name, method_control)) then
3361  success = .false.
3362  exit
3363  else
3364  method_name = trim(method_name)//trim(this_field_p%name)
3365  endif
3366  endif
3367 
3368  case(integer_type)
3369  write (scratch,*) this_field_p%i_value
3370  call concat_strings(method_control, comma//trim(this_field_p%name)//' = '//trim(adjustl(scratch)))
3371 
3372  case(logical_type)
3373  write (scratch,'(l1)')this_field_p%l_value
3374  call concat_strings(method_control, comma//trim(this_field_p%name)//' = '//trim(adjustl(scratch)))
3375 
3376  case(real_type)
3377  write (scratch,*) this_field_p%r_value
3378  call concat_strings(method_control, comma//trim(this_field_p%name)//' = '//trim(adjustl(scratch)))
3379 
3380  case(string_type)
3381  call concat_strings(method_control, comma//trim(this_field_p%name)//' = '//trim(this_field_p%s_value(1)))
3382  do i = 2, this_field_p%max_index
3383  call concat_strings(method_control, comma//trim(this_field_p%s_value(i)))
3384  enddo
3385 
3386  case default
3387  success = .false.
3388  exit
3389 
3390  end select
3391  this_field_p => this_field_p%next
3392  enddo
3393 endif
3394 
3395 end function query_method
3396 
3397 !> private function: appends str2 to the end of str1, with length check
3398 subroutine concat_strings(str1,str2)
3399  character(*), intent(inout) :: str1
3400  character(*), intent(in) :: str2
3401 
3402  character(64) :: n1,n2 ! for error reporting
3403 
3404  if (len_trim(str1)+len_trim(str2)>len(str1)) then
3405  write(n1,*)len(str1)
3406  write(n2,*)len_trim(str1)+len_trim(str2)
3407  call mpp_error(fatal,'length of output string ('//trim(adjustl(n1))&
3408  //') is not enough for the result of concatenation (len='&
3409  //trim(adjustl(n2))//')')
3410  endif
3411  str1 = trim(str1)//trim(str2)
3412 end subroutine concat_strings
3413 
3414 !> A function that allows the user to copy a field and add a suffix to
3415 !! the name of the new field.
3416 !!
3417 !> Given the name of a pre-existing field and a suffix, this function
3418 !! will create a new field. The name of the new field will be that of
3419 !! the old field with a suffix supplied by the user.
3420 !! @return index of the field that has been created by the copy
3421 function fm_copy_list(list_name, suffix, create ) &
3422  result(index)
3423 integer :: index
3424 character(len=*), intent(in) :: list_name !< name of a field that the user wishes to copy
3425 character(len=*), intent(in) :: suffix !< suffix that will be added to list_name when
3426  !! field is copied
3427 logical, intent(in), optional :: create !< flag to create new list if applicable
3428 
3429 character(len=fm_string_len), dimension(:), allocatable :: control
3430 character(len=fm_string_len), dimension(:), allocatable :: method
3431 character(len=fm_string_len) :: head
3432 character(len=fm_string_len) :: list_name_new
3433 character(len=fm_string_len) :: tail
3434 character(len=fm_string_len) :: val_str
3435 integer :: n
3436 integer :: num_meth
3437 integer :: val_int
3438 logical :: found_methods
3439 logical :: got_value
3440 logical :: recursive_t
3441 logical :: success
3442 logical :: val_logical
3443 real(r8_kind) :: val_real
3444 type (field_def), pointer, save :: temp_field_p
3445 type (field_def), pointer, save :: temp_list_p
3446 integer :: out_unit
3447 
3448 out_unit = stdout()
3449 
3450 
3451 num_meth= 1
3452 list_name_new = trim(list_name)//trim(suffix)
3453  recursive_t = .true.
3454 ! Initialize the field manager if needed
3455 if (.not. module_is_initialized) then
3457 endif
3458 
3459 if (list_name .eq. ' ') then
3460 ! If list is empty, then dump the current list
3461  temp_list_p => current_list_p
3462  success = .true.
3463 else
3464 ! Get a pointer to the list
3465  temp_list_p => find_list(list_name, current_list_p, .false.)
3466  if (associated(temp_list_p)) then
3467  success = .true.
3468  else
3469 ! Error following the path
3470  success = .false.
3471  endif
3472 endif
3473 ! Find the list
3474 if (success) then
3475  found_methods = fm_find_methods(trim(list_name), method, control)
3476  do n = 1, size(method)
3477  if (len_trim(method(n)) > 0 ) then
3478  index = fm_new_list(trim(list_name_new)//list_sep//method(n), create = create)
3479  call find_base(method(n), head, tail)
3480  temp_field_p => find_list(trim(list_name)//list_sep//head,temp_list_p, .false.)
3481  temp_field_p => find_field(tail,temp_field_p)
3482  select case (temp_field_p%field_type)
3483  case (integer_type)
3484  got_value = fm_get_value( trim(list_name)//list_sep//method(n), val_int)
3485  if ( fm_new_value( trim(list_name_new)//list_sep//method(n), val_int, &
3486  create = create, append = .true.) < 0 ) &
3487  call mpp_error(fatal, trim(error_header)//'Could not set the '//trim(method(n))//&
3488  ' for '//trim(list_name)//trim(suffix))
3489 
3490  case (logical_type)
3491  got_value = fm_get_value( trim(list_name)//list_sep//method(n), val_logical)
3492  if ( fm_new_value( trim(list_name_new)//list_sep//method(n), val_logical, &
3493  create = create, append = .true.) < 0 ) &
3494  call mpp_error(fatal, trim(error_header)//'Could not set the '//trim(method(n))//&
3495  ' for '//trim(list_name)//trim(suffix))
3496 
3497  case (real_type)
3498  got_value = fm_get_value( trim(list_name)//list_sep//method(n), val_real)
3499  if ( fm_new_value( trim(list_name_new)//list_sep//method(n), val_real, &
3500  create = create, append = .true.) < 0 ) &
3501  call mpp_error(fatal, trim(error_header)//'Could not set the '//trim(method(n))//&
3502  ' for '//trim(list_name)//trim(suffix))
3503 
3504  case (string_type)
3505  got_value = fm_get_value( trim(list_name)//list_sep//method(n), val_str)
3506  if ( fm_new_value( trim(list_name_new)//list_sep//method(n), val_str, &
3507  create = create, append = .true.) < 0 ) &
3508  call mpp_error(fatal, trim(error_header)//'Could not set the '//trim(method(n))//&
3509  ' for '//trim(list_name)//trim(suffix))
3510  case default
3511  end select
3512 
3513  endif
3514  enddo
3515 endif
3516 
3517 end function fm_copy_list
3518 
3519 !> This function retrieves all the methods associated with a field.
3520 !!
3521 !> This is different from fm_query_method in that this function gets all
3522 !! the methods associated as opposed to 1 method.
3523 !! @return A flag to indicate whether the function operated with (FALSE) or
3524 !! without (TRUE) errors.
3525 function fm_find_methods(list_name, methods, control ) &
3526  result(success)
3527 logical :: success
3528 character(len=*), intent(in) :: list_name !< The name of a list that the user wishes to find methods for
3529 character(len=*), intent(out), dimension(:) :: methods !< An array of the methods associated with list_name
3530 character(len=*), intent(out), dimension(:) :: control !< An array of the parameters associated with methods
3531 
3532 integer :: num_meth
3533 logical :: recursive_t
3534 type (field_def), pointer, save :: temp_list_p
3535 integer :: out_unit
3536 
3537 out_unit = stdout()
3538 num_meth= 1
3539 ! Check whether to do things recursively
3540  recursive_t = .true.
3541 
3542 if (.not. module_is_initialized) then
3544 endif
3545 
3546 if (list_name .eq. ' ') then
3547 ! If list is empty, then dump the current list
3548  temp_list_p => current_list_p
3549  success = .true.
3550 else
3551 ! Get a pointer to the list
3552  temp_list_p => find_list(list_name, current_list_p, .false.)
3553  if (associated(temp_list_p)) then
3554  success = .true.
3555  else
3556 ! Error following the path
3557  success = .false.
3558  endif
3559 endif
3560 ! Find the list
3561 if (success) then
3562  success = find_method(temp_list_p, recursive_t, num_meth, methods, control)
3563 endif
3564 
3565 end function fm_find_methods
3566 
3567 !> Given a field list pointer this function retrieves methods and
3568 !! associated parameters for the field list.
3569 !! @returns A flag to indicate whether the function operated with (FALSE) or
3570 !! without (TRUE) errors.
3571 recursive function find_method(list_p, recursive, num_meth, method, control) &
3572  result(success)
3573 logical :: success
3574 type (field_def), pointer :: list_p !< A pointer to the field of interest
3575 logical, intent(in) :: recursive !< If true, search recursively for fields
3576 integer, intent(inout) :: num_meth !< The number of methods found
3577 character(len=*), intent(out), dimension(:) :: method !< The methods associated with the field pointed to by list_p
3578 character(len=*), intent(out), dimension(:) :: control !< The control parameters for the methods found
3579 
3580 character(len=FMS_PATH_LEN) :: scratch
3581 integer :: i
3582 integer :: n
3583 type (field_def), pointer, save :: this_field_p
3584 integer :: out_unit
3585 
3586 out_unit = stdout()
3587 ! Check for a valid list
3588 if (.not. associated(list_p) .or. list_p%field_type .ne. list_type) then
3589  success = .false.
3590 else
3591 ! set the default return value
3592  success = .true.
3593 
3594  this_field_p => list_p%first_field
3595 
3596  do while (associated(this_field_p))
3597  select case(this_field_p%field_type)
3598  case(list_type)
3599 ! If this is a list, then this is the method name
3600  if ( this_field_p%length > 1) then
3601  do n = num_meth+1, num_meth + this_field_p%length - 1
3602  write (method(n),'(a,a,a,$)') trim(method(num_meth)), &
3603  trim(this_field_p%name), list_sep
3604  enddo
3605  write (method(num_meth),'(a,a,a,$)') trim(method(num_meth)), &
3606  trim(this_field_p%name), list_sep
3607  else
3608  write (method(num_meth),'(a,a,a,$)') trim(method(num_meth)), &
3609  trim(this_field_p%name), list_sep
3610  endif
3611  success = find_method(this_field_p, .true., num_meth, method, control)
3612 
3613  case(integer_type)
3614  write (scratch,*) this_field_p%i_value
3615  call strip_front_blanks(scratch)
3616  write (method(num_meth),'(a,a)') trim(method(num_meth)), &
3617  trim(this_field_p%name)
3618  write (control(num_meth),'(a)') &
3619  trim(scratch)
3620  num_meth = num_meth + 1
3621 
3622 
3623  case(logical_type)
3624 
3625  write (method(num_meth),'(a,a)') trim(method(num_meth)), &
3626  trim(this_field_p%name)
3627  write (control(num_meth),'(l1)') &
3628  this_field_p%l_value
3629  num_meth = num_meth + 1
3630 
3631  case(real_type)
3632 
3633  if(allocated(this_field_p%r_value)) write (scratch,*) this_field_p%r_value
3634  call strip_front_blanks(scratch)
3635  write (method(num_meth),'(a,a)') trim(method(num_meth)), &
3636  trim(this_field_p%name)
3637  write (control(num_meth),'(a)') &
3638  trim(scratch)
3639  num_meth = num_meth + 1
3640 
3641 
3642  case(string_type)
3643  write (method(num_meth),'(a,a)') trim(method(num_meth)), &
3644  trim(this_field_p%name)
3645  write (control(num_meth),'(a)') &
3646  trim(this_field_p%s_value(1))
3647  do i = 2, this_field_p%max_index
3648  write (control(num_meth),'(a,a,$)') comma//trim(this_field_p%s_value(i))
3649  enddo
3650  num_meth = num_meth + 1
3651 
3652 
3653  case default
3654  success = .false.
3655  exit
3656 
3657  end select
3658 
3659  this_field_p => this_field_p%next
3660  enddo
3661 endif
3662 
3663 end function find_method
3664 
3665 #include "field_manager_r4.fh"
3666 #include "field_manager_r8.fh"
3667 
3668 end module field_manager_mod
3669 !> @}
3670 ! close documentation grouping
integer function parse_integer(text, label, parse_ival)
integer function, public fm_copy_list(list_name, suffix, create)
A function that allows the user to copy a field and add a suffix to the name of the new field.
type(field_def) function, pointer, private find_list(path, relative_p, create)
Find and return a pointer to the specified list, relative to relative_p. Return a null pointer on err...
recursive logical function query_method(list_p, recursive, name, method_name, method_control)
A private function that can recursively recover values for parameters associated with a field.
subroutine, private find_base(name, path, base)
A subroutine that splits a listname into a path and a base.
integer function fm_new_value_integer(name, new_ival, create, index, append)
Assigns a given value to a given field.
integer, parameter, public fm_string_len
The length of a character string representing character values for the field.
integer function fm_new_value_logical(name, new_lval, create, index, append)
Assigns a given value to a given field.
logical function set_list_stuff()
If the pointer matches to the right list, extract the field information. Used in fm_loop_over_list.
type(field_def) function, pointer, private create_field(parent_p, name)
A function to create a field as a child of parent_p. This will return a pointer to a field_def type.
subroutine strip_front_blanks(name)
A routine to strip whitespace from the start of character strings.
integer, parameter, public model_land
Land model.
function parse_strings(text, label, values)
subroutine concat_strings(str1, str2)
private function: appends str2 to the end of str1, with length check
integer function find_field_index_new(field_name)
integer function, public fm_get_length(name)
A function to return how many elements are contained within the named list or entry.
logical function fm_get_value_integer(name, get_ival, index)
character(len=8) function, public fm_get_type(name)
A function to return the type of the named field.
subroutine, public fm_reset_loop
Resets the loop variable. For use in conjunction with fm_loop_over_list.
logical function, public fm_exists(name)
A function to test whether a named field exists.
character(len=11), dimension(num_models), parameter, public model_names
Model names, e.g. MODEL_NAMES(MODEL_OCEAN) is 'oceanic'.
logical function, public fm_change_list(name)
Change the current list. Return true on success, false otherwise.
subroutine read_field_table_yaml(nfields, table_name)
Routine to read and parse the field table yaml.
integer function, public fm_new_list(name, create, keep)
A function to create a new list.
type(field_mgr_type), dimension(:), allocatable, private fields
fields of field_mgr_type
subroutine, public get_field_methods(n, methods)
A routine to obtain all the methods associated with a field.
character(len=fms_path_len) function, public fm_get_current_list()
A function to return the full path of the current list.
subroutine, private find_head(name, head, rest)
Find the first list for a name by splitting the name into a head and the rest.
integer function, public fm_get_index(name)
A function to return the index of a named field.
subroutine, public field_manager_end
Destructor for field manager.
integer function find_field_index_old(model, field_name)
Function to return the index of the field.
subroutine new_name_yaml(list_name, method_name_in, val_name_in)
Subroutine to add new values to list parameters.
type(field_def) function, pointer, private get_field(name, this_list_p)
Return a pointer to the field if it exists relative to this_list_p, null otherwise.
subroutine, public field_manager_init(nfields, table_name)
Routine to initialize the field manager.
logical function, public fm_query_method(name, method_name, method_control)
This is a function that provides the capability to return parameters associated with a field in a pai...
integer, parameter, public model_coupler
Ice model.
subroutine, private initialize_module_variables
A function to initialize the values of the pointers. This will remove all fields and reset the field ...
integer, parameter, public fm_field_name_len
The length of a character string representing the field name.
type(field_def) function, pointer, private make_list(this_list_p, name)
This function creates a new field and returns a pointer to that field.
recursive logical function find_method(list_p, recursive, num_meth, method, control)
Given a field list pointer this function retrieves methods and associated parameters for the field li...
subroutine, public get_field_method(n, m, method)
A routine to get a specified method.
subroutine, public get_field_info(n, fld_type, fld_name, model, num_methods)
This routine allows access to field information given an index.
logical function fm_get_value_logical(name, get_lval, index)
integer, parameter, public fm_type_name_len
The length of a character string representing the various types that the values of the field can take...
logical use_field_table_yaml
.True. if using the field_table.yaml, .false. if using the legacy field_table
subroutine read_field_table_legacy(nfields, table_name)
Routine to read and parse the field table yaml.
logical function, public fm_modify_name(oldname, newname)
This function allows a user to rename a field without modifying the contents of the field.
logical function, public fm_dump_list(name, recursive, unit)
A function to list properties associated with a field.
logical recursive function, private dump_list(list_p, recursive, depth, out_unit)
This is a function that lists the parameters of a field.
integer, parameter, public model_ocean
Ocean model.
subroutine, public fm_return_root
Return the root list to the value at initialization.
integer, parameter, public no_field
The value returned if a field is not defined.
subroutine, public fm_init_loop(loop_list, iter)
given a name of the list, prepares an iterator over the list content. If the name of the given list i...
function parse_integers(text, label, values)
logical function fm_loop_over_list_new(iter, name, field_type, index)
given a list iterator, returns information about curren list element and advances the iterator to the...
integer, parameter, public model_ice
Ice model.
integer function fm_new_value_string(name, new_sval, create, index, append)
Assigns a given value to a given field.
logical function, public fm_change_root(name)
Change the root list.
subroutine new_name(list_name, method_name_in, val_name_in)
Subroutine to add new values to list parameters.
integer function parse_string(text, label, parse_sval)
integer, parameter, public model_atmos
Atmospheric model.
integer, parameter, public num_models
Number of models (ATMOS, OCEAN, LAND, ICE, COUPLER).
type(field_def) function, pointer, private find_field(name, this_list_p)
Find and return a pointer to the field in the specified list. Return a null pointer on error.
logical function, public fm_find_methods(list_name, methods, control)
This function retrieves all the methods associated with a field.
logical function fm_loop_over_list_old(list, name, field_type, index)
Iterates through the given list.
logical function fm_get_value_string(name, get_sval, index)
Returns an index corresponding to the given field name.
An overloaded function to find and extract a value for a named field.
A function for looping over a list.
An overloaded function to assign a value to a field.
A function to parse an integer or an array of integers, a real or an array of reals,...
Private type for internal use.
Private type for internal use.
Private type for internal use.
Private type for internal use.
Iterator over the field manager list.
This method_type is a way to allow a component module to alter the parameters it needs for various tr...
This method_type is the same as method_type except that the method_control string is not present....
This is the same as method_type except that the method_control and method_name strings are not presen...
subroutine, public build_fmtable(fmTable, filename)
Subroutine to populate an fmTable by reading a yaml file, given an optional filename.
Definition: fm_yaml.F90:97
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
Definition: fms.F90:523
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:701
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:42
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
Definition: mpp_util.inc:58
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381