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