FV3 Bundle
fm_util.F90
Go to the documentation of this file.
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 
20 module fm_util_mod !{
21 !
22 !<CONTACT EMAIL="Richard.Slater@noaa.gov"> Richard D. Slater
23 !</CONTACT>
24 !
25 !<REVIEWER EMAIL="John.Dunne@noaa.gov"> John P. Dunne
26 !</REVIEWER>
27 !
28 !<OVERVIEW>
29 ! Utility routines for the field manager
30 !</OVERVIEW>
31 !
32 !<DESCRIPTION>
33 ! This module provides utility routines for the field manager.
34 ! Basically, it provides for error catching, reporting and
35 ! termination while interfacing with the field manager.
36 !</DESCRIPTION>
37 !
38 ! <INFO>
39 ! </INFO>
40 !
41 
47 use fms_mod, only: fatal, stdout
48 use mpp_mod, only: mpp_error
49 
50 implicit none
51 
52 private
53 
57 public fm_util_set_caller
63 public fm_util_get_length
66 public fm_util_get_real
67 public fm_util_get_string
72 public fm_util_set_value
81 !public fm_util_get_index
84 
85 !
86 ! Public variables
87 !
88 
89 character(len=128), public :: fm_util_default_caller = ' '
90 
91 !
92 ! private parameters
93 !
94 
95 character(len=48), parameter :: mod_name = 'fm_util_mod'
96 
97 !
98 ! Private variables
99 !
100 
101 character(len=128) :: save_default_caller = ' '
102 character(len=128) :: default_good_name_list = ' '
103 character(len=128) :: save_default_good_name_list = ' '
104 logical :: default_no_overwrite = .false.
105 logical :: save_default_no_overwrite = .false.
106 character(len=fm_path_name_len) :: save_current_list
107 character(len=fm_path_name_len) :: save_path
108 character(len=fm_path_name_len) :: save_name
109 ! Include variable "version" to be written to log file.
110 #include<file_version.h>
111 
112 !
113 ! Interface definitions for overloaded routines
114 !
115 
116 !interface fm_util_get_value !{
117  !module procedure fm_util_get_value_integer
118  !module procedure fm_util_get_value_logical
119  !module procedure fm_util_get_value_real
120  !module procedure fm_util_get_value_string
121  !module procedure fm_util_get_value_integer_array
122  !module procedure fm_util_get_value_logical_array
123  !module procedure fm_util_get_value_real_array
124  !module procedure fm_util_get_value_string_array
125 !end interface !}
126 
127 interface fm_util_set_value !{
128  module procedure fm_util_set_value_integer_array
129  module procedure fm_util_set_value_logical_array
130  module procedure fm_util_set_value_real_array
131  module procedure fm_util_set_value_string_array
132  module procedure fm_util_set_value_integer
133  module procedure fm_util_set_value_logical
134  module procedure fm_util_set_value_real
135  module procedure fm_util_set_value_string
136 end interface !}
137 
138 !interface fm_util_get_index !{
139  !module procedure fm_util_get_index_list
140  !module procedure fm_util_get_index_string
141 !end interface !}
142 
143 
144 contains
145 
146 
147 !#######################################################################
148 ! <SUBROUTINE NAME="fm_util_set_caller">
149 !
150 ! <DESCRIPTION>
151 ! Set the default value for the optional "caller" variable used in many of these
152 ! subroutines. If the argument is blank, then set the default to blank, otherwise
153 ! the deault will have brackets placed around the argument.
154 !
155 ! </DESCRIPTION>
156 !
157 
158 subroutine fm_util_set_caller(caller) !{
160 implicit none
161 
162 !
163 ! arguments
164 !
165 
166 character(len=*), intent(in) :: caller
167 
168 !
169 ! Local variables
170 !
171 
172 !
173 ! save the default caller string
174 !
175 
177 
178 !
179 ! set the default caller string
180 !
181 
182 if (caller .eq. ' ') then !{
184 else !}{
185  fm_util_default_caller = '[' // trim(caller) // ']'
186 endif !}
187 
188 return
189 
190 end subroutine fm_util_set_caller !}
191 ! </SUBROUTINE> NAME="fm_util_set_caller"
192 
193 
194 !#######################################################################
195 ! <SUBROUTINE NAME="fm_util_reset_caller">
196 !
197 ! <DESCRIPTION>
198 ! Reset the default value for the optional "caller" variable used in many of these
199 ! subroutines to blank.
200 !
201 ! </DESCRIPTION>
202 !
203 
204 subroutine fm_util_reset_caller !{
206 implicit none
207 
208 !
209 ! arguments
210 !
211 
212 !
213 ! Local variables
214 !
215 
216 !
217 ! reset the default caller string
218 !
219 
222 
223 return
224 
225 end subroutine fm_util_reset_caller !}
226 ! </SUBROUTINE> NAME="fm_util_reset_caller"
227 
228 
229 !#######################################################################
230 ! <SUBROUTINE NAME="fm_util_set_good_name_list">
231 !
232 ! <DESCRIPTION>
233 ! Set the default value for the optional "good_name_list" variable used in many of these
234 ! subroutines.
235 !
236 ! </DESCRIPTION>
237 !
238 
239 subroutine fm_util_set_good_name_list(good_name_list) !{
241 implicit none
242 
243 !
244 ! arguments
245 !
246 
247 character(len=*), intent(in) :: good_name_list
248 
249 !
250 ! Local variables
251 !
252 
253 !
254 ! save the default good_name_list string
255 !
256 
258 
259 !
260 ! set the default good_name_list string
261 !
262 
263 default_good_name_list = good_name_list
264 
265 return
266 
267 end subroutine fm_util_set_good_name_list !}
268 ! </SUBROUTINE> NAME="fm_util_set_good_name_list"
269 
270 
271 !#######################################################################
272 ! <SUBROUTINE NAME="fm_util_reset_good_name_list">
273 !
274 ! <DESCRIPTION>
275 ! Reset the default value for the optional "good_name_list" variable used in many of these
276 ! subroutines to the saved value.
277 !
278 ! </DESCRIPTION>
279 !
280 
281 subroutine fm_util_reset_good_name_list !{
283 implicit none
284 
285 !
286 ! arguments
287 !
288 
289 !
290 ! Local variables
291 !
292 
293 !
294 ! reset the default good_name_list string
295 !
296 
299 
300 return
301 
302 end subroutine fm_util_reset_good_name_list !}
303 ! </SUBROUTINE> NAME="fm_util_reset_good_name_list"
304 
305 
306 !#######################################################################
307 ! <SUBROUTINE NAME="fm_util_set_no_overwrite">
308 !
309 ! <DESCRIPTION>
310 ! Set the default value for the optional "no_overwrite" variable used in some of these
311 ! subroutines.
312 !
313 ! </DESCRIPTION>
314 !
315 
316 subroutine fm_util_set_no_overwrite(no_overwrite) !{
318 implicit none
319 
320 !
321 ! arguments
322 !
323 
324 logical, intent(in) :: no_overwrite
325 
326 !
327 ! Local variables
328 !
329 
330 !
331 ! save the default no_overwrite string
332 !
333 
335 
336 !
337 ! set the default no_overwrite value
338 !
339 
340 default_no_overwrite = no_overwrite
341 
342 return
343 
344 end subroutine fm_util_set_no_overwrite !}
345 ! </SUBROUTINE> NAME="fm_util_set_no_overwrite"
346 
347 
348 !#######################################################################
349 ! <SUBROUTINE NAME="fm_util_reset_no_overwrite">
350 !
351 ! <DESCRIPTION>
352 ! Reset the default value for the optional "no_overwrite" variable used in some of these
353 ! subroutines to false.
354 !
355 ! </DESCRIPTION>
356 !
357 
358 subroutine fm_util_reset_no_overwrite !{
360 implicit none
361 
362 !
363 ! arguments
364 !
365 
366 !
367 ! Local variables
368 !
369 
370 !
371 ! reset the default no_overwrite value
372 !
373 
376 
377 return
378 
379 end subroutine fm_util_reset_no_overwrite !}
380 ! </SUBROUTINE> NAME="fm_util_reset_no_overwrite"
381 
382 
383 !#######################################################################
384 ! <SUBROUTINE NAME="fm_util_check_for_bad_fields">
385 !
386 ! <DESCRIPTION>
387 ! Check for unrecognized fields in a list
388 !
389 ! </DESCRIPTION>
390 !
391 
392 subroutine fm_util_check_for_bad_fields(list, good_fields, caller) !{
394 implicit none
395 
396 !
397 ! arguments
398 !
399 
400 character(len=*), intent(in) :: list
401 character(len=*), intent(in), dimension(:) :: good_fields
402 character(len=*), intent(in), optional :: caller
403 
404 !
405 ! Local parameters
406 !
407 
408 character(len=48), parameter :: sub_name = 'fm_util_check_for_bad_fields'
409 
410 !
411 ! Local variables
412 !
413 
414 logical :: fm_success
415 integer :: i
416 integer :: ind
417 integer :: list_length
418 integer :: good_length
419 character(len=fm_type_name_len) :: typ
420 character(len=fm_field_name_len) :: name
421 logical :: found
422 character(len=256) :: error_header
423 character(len=256) :: warn_header
424 character(len=256) :: note_header
425 character(len=128) :: caller_str
426 integer :: out_unit
427 
428 out_unit = stdout()
429 
430 !
431 ! set the caller string and headers
432 !
433 
434 if (present(caller)) then !{
435  caller_str = '[' // trim(caller) // ']'
436 else !}{
437  caller_str = fm_util_default_caller
438 endif !}
439 
440 error_header = '==>Error from ' // trim(mod_name) // &
441  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
442 warn_header = '==>Warning from ' // trim(mod_name) // &
443  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
444 note_header = '==>Note from ' // trim(mod_name) // &
445  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
446 
447 !
448 ! check that a list is given (fatal if not)
449 !
450 
451 if (list .eq. ' ') then !{
452  write (out_unit,*) trim(error_header) // ' Empty list given'
453  call mpp_error(fatal, trim(error_header) // ' Empty list given')
454 endif !}
455 
456 !
457 ! Check that we have been given a list
458 !
459 
460 if (fm_get_type(list) .ne. 'list') then !{
461  write (out_unit,*) trim(error_header) // ' Not given a list: ' // trim(list)
462  call mpp_error(fatal, trim(error_header) // ' Not given a list: ' // trim(list))
463 endif !}
464 
465 !
466 ! Get the list length
467 !
468 
469 list_length = fm_get_length(list)
470 if (list_length .lt. 0) then !{
471  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(list))
472 endif !}
473 
474 !
475 ! Get the number of good fields
476 !
477 
478 good_length = size(good_fields)
479 
480 if (list_length .lt. good_length) then !{
481 
482 !
483 ! If the list length is less than the number of good fields this is an error
484 ! as the list should be fully populated and we'll check which extra fields
485 ! are given in good_fields
486 !
487 
488  write (out_unit,*) trim(error_header), ' List length < number of good fields (', &
489  list_length, ' < ', good_length, ') in list ', trim(list)
490 
491  write (out_unit,*)
492  write (out_unit,*) 'The list contains the following fields:'
493  fm_success= fm_dump_list(list, .false.)
494  write (out_unit,*)
495  write (out_unit,*) 'The supposed list of good fields is:'
496  do i = 1, good_length !{
497  if (fm_exists(trim(list) // '/' // good_fields(i))) then !{
498  write (out_unit,*) 'List field: "', trim(good_fields(i)), '"'
499  else !}{
500  write (out_unit,*) 'EXTRA good field: "', trim(good_fields(i)), '"'
501  endif !}
502  enddo !} i
503  write (out_unit,*)
504 
505  call mpp_error(fatal, trim(error_header) // &
506  ' List length < number of good fields for list: ' // trim(list))
507 
508 elseif (list_length .gt. good_length) then !}{
509 
510 !
511 ! If the list length is greater than the number of good fields this is an error
512 ! as the there should not be any more fields than those given in the good fields list
513 ! and we'll check which extra fields are given in the list
514 !
515 
516  write (out_unit,*) trim(warn_header), 'List length > number of good fields (', &
517  list_length, ' > ', good_length, ') in list ', trim(list)
518 
519  write (out_unit,*) trim(error_header), ' Start of list of fields'
520  do while (fm_loop_over_list(list, name, typ, ind)) !{
521  found = .false.
522  do i = 1, good_length !{
523  found = found .or. (name .eq. good_fields(i))
524  enddo !} i
525  if (found) then !{
526  write (out_unit,*) 'Good list field: "', trim(name), '"'
527  else !}{
528  write (out_unit,*) 'EXTRA list field: "', trim(name), '"'
529  endif !}
530  enddo !}
531  write (out_unit,*) trim(error_header), ' End of list of fields'
532 
533  call mpp_error(fatal, trim(error_header) // &
534  ' List length > number of good fields for list: ' // trim(list))
535 
536 endif !}
537 
538 !
539 ! If the list length equals the number of good fields then all is good
540 !
541 
542 return
543 
544 end subroutine fm_util_check_for_bad_fields !}
545 ! </SUBROUTINE> NAME="fm_util_check_for_bad_fields"
546 
547 
548 !#######################################################################
549 ! <FUNCTION NAME="fm_util_get_length">
550 !
551 ! <DESCRIPTION>
552 ! Get the length of an element of the Field Manager tree
553 ! </DESCRIPTION>
554 !
555 function fm_util_get_length(name, caller) &
556  result(field_length) !{
558 implicit none
559 
560 !
561 ! Return type
562 !
563 
564 integer :: field_length
565 
566 !
567 ! arguments
568 !
569 
570 character(len=*), intent(in) :: name
571 character(len=*), intent(in), optional :: caller
572 
573 !
574 ! Local parameters
575 !
576 
577 character(len=48), parameter :: sub_name = 'fm_util_get_length'
578 
579 !
580 ! Local variables
581 !
582 
583 character(len=256) :: error_header
584 character(len=256) :: warn_header
585 character(len=256) :: note_header
586 character(len=128) :: caller_str
587 
588 !
589 ! set the caller string and headers
590 !
591 
592 if (present(caller)) then !{
593  caller_str = '[' // trim(caller) // ']'
594 else !}{
595  caller_str = fm_util_default_caller
596 endif !}
597 
598 error_header = '==>Error from ' // trim(mod_name) // &
599  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
600 warn_header = '==>Warning from ' // trim(mod_name) // &
601  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
602 note_header = '==>Note from ' // trim(mod_name) // &
603  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
604 
605 !
606 ! check that a name is given (fatal if not)
607 !
608 
609 if (name .eq. ' ') then !{
610  call mpp_error(fatal, trim(error_header) // ' Empty name given')
611 endif !}
612 
613 !
614 ! Get the field's length
615 !
616 
617 field_length = fm_get_length(name)
618 if (field_length .lt. 0) then !{
619  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
620 endif !}
621 
622 return
623 
624 end function fm_util_get_length !}
625 ! </FUNCTION> NAME="fm_util_get_length"
626 
627 
628 !#######################################################################
629 ! <FUNCTION NAME="fm_util_get_index_string">
630 !
631 ! <DESCRIPTION>
632 ! Get the index of an element of a string in the Field Manager tree
633 ! </DESCRIPTION>
634 !
635 function fm_util_get_index_string(name, string, caller) &
636  result(fm_index) !{
638 implicit none
639 
640 !
641 ! Return type
642 !
643 
644 integer :: fm_index
645 
646 !
647 ! arguments
648 !
649 
650 character(len=*), intent(in) :: name
651 character(len=*), intent(in) :: string
652 character(len=*), intent(in), optional :: caller
653 
654 !
655 ! Local parameters
656 !
657 
658 character(len=48), parameter :: sub_name = 'fm_util_get_index_string'
659 
660 !
661 ! Local variables
662 !
663 
664 character(len=256) :: error_header
665 character(len=256) :: warn_header
666 character(len=256) :: note_header
667 character(len=128) :: caller_str
668 character(len=32) :: index_str
669 character(len=fm_type_name_len) :: fm_type
670 character(len=fm_string_len) :: fm_string
671 integer :: i
672 integer :: length
673 
674 !
675 ! set the caller string and headers
676 !
677 
678 if (present(caller)) then !{
679  caller_str = '[' // trim(caller) // ']'
680 else !}{
681  caller_str = fm_util_default_caller
682 endif !}
683 
684 error_header = '==>Error from ' // trim(mod_name) // &
685  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
686 warn_header = '==>Warning from ' // trim(mod_name) // &
687  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
688 note_header = '==>Note from ' // trim(mod_name) // &
689  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
690 
691 !
692 ! check that a name is given (fatal if not)
693 !
694 
695 if (name .eq. ' ') then !{
696  call mpp_error(fatal, trim(error_header) // ' Empty name given')
697 endif !}
698 
699 !
700 ! Check the field's type and get the index
701 !
702 
703 fm_index = 0
704 fm_type = fm_get_type(name)
705 if (fm_type .eq. 'string') then !{
706  length = fm_get_length(name)
707  if (length .lt. 0) then !{
708  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
709  endif !}
710  if (length .gt. 0) then !{
711  do i = 1, length !{
712  if (.not. fm_get_value(name, fm_string, index = i)) then !{
713  write (index_str,*) '(', i, ')'
714  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name) // trim(index_str))
715  endif !}
716  if (fm_string .eq. string) then !{
717  fm_index = i
718  exit
719  endif !}
720  enddo !} i
721  endif !}
722 elseif (fm_type .eq. ' ') then !}{
723  call mpp_error(fatal, trim(error_header) // ' Array does not exist: ' // trim(name))
724 else !}{
725  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
726 endif !}
727 
728 !if (fm_index .eq. 0) then !{
729  !call mpp_error(FATAL, trim(error_header) // ' "' // trim(string) // '" does not exist in ' // trim(name))
730 !endif !}
731 
732 return
733 
734 end function fm_util_get_index_string !}
735 ! </FUNCTION> NAME="fm_util_get_index_string"
736 
737 
738 !#######################################################################
739 ! <FUNCTION NAME="fm_util_get_index_list">
740 !
741 ! <DESCRIPTION>
742 ! Get the length of an element of the Field Manager tree
743 ! </DESCRIPTION>
744 !
745 function fm_util_get_index_list(name, caller) &
746  result(fm_index) !{
748 implicit none
749 
750 !
751 ! Return type
752 !
753 
754 integer :: fm_index
755 
756 !
757 ! arguments
758 !
759 
760 character(len=*), intent(in) :: name
761 character(len=*), intent(in), optional :: caller
762 
763 !
764 ! Local parameters
765 !
766 
767 character(len=48), parameter :: sub_name = 'fm_util_get_index_list'
768 
769 !
770 ! Local variables
771 !
772 
773 character(len=256) :: error_header
774 character(len=256) :: warn_header
775 character(len=256) :: note_header
776 character(len=128) :: caller_str
777 character(len=fm_type_name_len) :: fm_type
778 
779 !
780 ! set the caller string and headers
781 !
782 
783 if (present(caller)) then !{
784  caller_str = '[' // trim(caller) // ']'
785 else !}{
786  caller_str = fm_util_default_caller
787 endif !}
788 
789 error_header = '==>Error from ' // trim(mod_name) // &
790  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
791 warn_header = '==>Warning from ' // trim(mod_name) // &
792  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
793 note_header = '==>Note from ' // trim(mod_name) // &
794  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
795 
796 !
797 ! check that a name is given (fatal if not)
798 !
799 
800 if (name .eq. ' ') then !{
801  call mpp_error(fatal, trim(error_header) // ' Empty name given')
802 endif !}
803 
804 !
805 ! Check the field's type and get the index
806 !
807 
808 fm_index = 0
809 fm_type = fm_get_type(name)
810 if (fm_type .eq. 'list') then !{
811  fm_index = fm_get_index(name)
812  if (fm_index .le. 0) then !{
813  call mpp_error(fatal, trim(error_header) // ' List does not exist: ' // trim(name))
814  endif !}
815 elseif (fm_type .eq. ' ') then !}{
816  call mpp_error(fatal, trim(error_header) // ' List does not exist: ' // trim(name))
817 else !}{
818  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
819 endif !}
820 
821 
822 return
823 
824 end function fm_util_get_index_list !}
825 ! </FUNCTION> NAME="fm_util_get_index_list"
826 
827 
828 !#######################################################################
829 ! <FUNCTION NAME="fm_util_get_integer_array">
830 !
831 ! <DESCRIPTION>
832 ! Get an integer value from the Field Manager tree.
833 ! </DESCRIPTION>
834 !
835 function fm_util_get_integer_array(name, caller) &
836  result(array) !{
838 implicit none
839 
840 !
841 ! Return type
842 !
843 
844 integer, pointer, dimension(:) :: array
845 
846 !
847 ! arguments
848 !
849 
850 character(len=*), intent(in) :: name
851 character(len=*), intent(in), optional :: caller
852 
853 !
854 ! Local parameters
855 !
856 
857 character(len=48), parameter :: sub_name = 'fm_util_get_integer_array'
858 
859 !
860 ! Local variables
861 !
862 
863 character(len=256) :: error_header
864 character(len=256) :: warn_header
865 character(len=256) :: note_header
866 character(len=128) :: caller_str
867 character(len=32) :: index_str
868 character(len=fm_type_name_len) :: fm_type
869 integer :: i
870 integer :: length
871 
872 nullify(array)
873 
874 !
875 ! set the caller string and headers
876 !
877 
878 if (present(caller)) then !{
879  caller_str = '[' // trim(caller) // ']'
880 else !}{
881  caller_str = fm_util_default_caller
882 endif !}
883 
884 error_header = '==>Error from ' // trim(mod_name) // &
885  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
886 warn_header = '==>Warning from ' // trim(mod_name) // &
887  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
888 note_header = '==>Note from ' // trim(mod_name) // &
889  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
890 
891 !
892 ! check that a name is given (fatal if not)
893 !
894 
895 if (name .eq. ' ') then !{
896  call mpp_error(fatal, trim(error_header) // ' Empty name given')
897 endif !}
898 
899 fm_type = fm_get_type(name)
900 if (fm_type .eq. 'integer') then !{
901  length = fm_get_length(name)
902  if (length .lt. 0) then !{
903  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
904  endif !}
905  if (length .gt. 0) then !{
906  allocate(array(length))
907  do i = 1, length !{
908  if (.not. fm_get_value(name, array(i), index = i)) then !{
909  write (index_str,*) '(', i, ')'
910  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name) // trim(index_str))
911  endif !}
912  enddo !} i
913  endif !}
914 elseif (fm_type .eq. ' ') then !}{
915  call mpp_error(fatal, trim(error_header) // ' Array does not exist: ' // trim(name))
916 else !}{
917  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
918 endif !}
919 
920 return
921 
922 end function fm_util_get_integer_array !}
923 ! </FUNCTION> NAME="fm_util_get_integer_array"
924 
925 
926 !#######################################################################
927 ! <FUNCTION NAME="fm_util_get_logical_array">
928 !
929 ! <DESCRIPTION>
930 ! Get a logical value from the Field Manager tree.
931 ! </DESCRIPTION>
932 !
933 function fm_util_get_logical_array(name, caller) &
934  result(array) !{
936 implicit none
937 
938 !
939 ! Return type
940 !
941 
942 logical, pointer, dimension(:) :: array
943 
944 !
945 ! arguments
946 !
947 
948 character(len=*), intent(in) :: name
949 character(len=*), intent(in), optional :: caller
950 
951 !
952 ! Local parameters
953 !
954 
955 character(len=48), parameter :: sub_name = 'fm_util_get_logical_array'
956 
957 !
958 ! Local variables
959 !
960 
961 character(len=256) :: error_header
962 character(len=256) :: warn_header
963 character(len=256) :: note_header
964 character(len=128) :: caller_str
965 character(len=32) :: index_str
966 character(len=fm_type_name_len) :: fm_type
967 integer :: i
968 integer :: length
969 
970 nullify(array)
971 
972 !
973 ! set the caller string and headers
974 !
975 
976 if (present(caller)) then !{
977  caller_str = '[' // trim(caller) // ']'
978 else !}{
979  caller_str = fm_util_default_caller
980 endif !}
981 
982 error_header = '==>Error from ' // trim(mod_name) // &
983  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
984 warn_header = '==>Warning from ' // trim(mod_name) // &
985  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
986 note_header = '==>Note from ' // trim(mod_name) // &
987  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
988 
989 !
990 ! check that a name is given (fatal if not)
991 !
992 
993 if (name .eq. ' ') then !{
994  call mpp_error(fatal, trim(error_header) // ' Empty name given')
995 endif !}
996 
997 fm_type = fm_get_type(name)
998 if (fm_type .eq. 'logical') then !{
999  length = fm_get_length(name)
1000  if (length .lt. 0) then !{
1001  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1002  endif !}
1003  if (length .gt. 0) then !{
1004  allocate(array(length))
1005  do i = 1, length !{
1006  if (.not. fm_get_value(name, array(i), index = i)) then !{
1007  write (index_str,*) '(', i, ')'
1008  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name) // trim(index_str))
1009  endif !}
1010  enddo !} i
1011  endif !}
1012 elseif (fm_type .eq. ' ') then !}{
1013  call mpp_error(fatal, trim(error_header) // ' Array does not exist: ' // trim(name))
1014 else !}{
1015  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1016 endif !}
1017 
1018 return
1019 
1020 end function fm_util_get_logical_array !}
1021 ! </FUNCTION> NAME="fm_util_get_logical_array"
1022 
1023 
1024 !#######################################################################
1025 ! <FUNCTION NAME="fm_util_get_real_array">
1026 !
1027 ! <DESCRIPTION>
1028 ! Get a real value from the Field Manager tree.
1029 ! </DESCRIPTION>
1030 !
1031 function fm_util_get_real_array(name, caller) &
1032  result(array) !{
1034 implicit none
1035 
1036 !
1037 ! Return type
1038 !
1039 
1040 real, pointer, dimension(:) :: array
1041 
1042 !
1043 ! arguments
1044 !
1045 
1046 character(len=*), intent(in) :: name
1047 character(len=*), intent(in), optional :: caller
1048 
1049 !
1050 ! Local parameters
1051 !
1052 
1053 character(len=48), parameter :: sub_name = 'fm_util_get_real_array'
1054 
1055 !
1056 ! Local variables
1057 !
1058 
1059 character(len=256) :: error_header
1060 character(len=256) :: warn_header
1061 character(len=256) :: note_header
1062 character(len=128) :: caller_str
1063 character(len=32) :: index_str
1064 character(len=fm_type_name_len) :: fm_type
1065 integer :: i
1066 integer :: length
1067 
1068 nullify(array)
1069 
1070 !
1071 ! set the caller string and headers
1072 !
1073 
1074 if (present(caller)) then !{
1075  caller_str = '[' // trim(caller) // ']'
1076 else !}{
1077  caller_str = fm_util_default_caller
1078 endif !}
1079 
1080 error_header = '==>Error from ' // trim(mod_name) // &
1081  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1082 warn_header = '==>Warning from ' // trim(mod_name) // &
1083  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1084 note_header = '==>Note from ' // trim(mod_name) // &
1085  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1086 
1087 !
1088 ! check that a name is given (fatal if not)
1089 !
1090 
1091 if (name .eq. ' ') then !{
1092  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1093 endif !}
1094 
1095 fm_type = fm_get_type(name)
1096 if (fm_type .eq. 'real') then !{
1097  length = fm_get_length(name)
1098  if (length .lt. 0) then !{
1099  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1100  endif !}
1101  if (length .gt. 0) then !{
1102  allocate(array(length))
1103  do i = 1, length !{
1104  if (.not. fm_get_value(name, array(i), index = i)) then !{
1105  write (index_str,*) '(', i, ')'
1106  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name) // trim(index_str))
1107  endif !}
1108  enddo !} i
1109  endif !}
1110 elseif (fm_type .eq. ' ') then !}{
1111  call mpp_error(fatal, trim(error_header) // ' Array does not exist: ' // trim(name))
1112 else !}{
1113  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1114 endif !}
1115 
1116 return
1117 
1118 end function fm_util_get_real_array !}
1119 ! </FUNCTION> NAME="fm_util_get_real_array"
1120 
1121 
1122 !#######################################################################
1123 ! <FUNCTION NAME="fm_util_get_string_array">
1124 !
1125 ! <DESCRIPTION>
1126 ! Get a string value from the Field Manager tree.
1127 ! </DESCRIPTION>
1128 !
1129 function fm_util_get_string_array(name, caller) &
1130  result(array) !{
1132 implicit none
1133 
1134 !
1135 ! Return type
1136 !
1137 
1138 character(len=fm_string_len), pointer, dimension(:) :: array
1139 
1140 !
1141 ! arguments
1142 !
1143 
1144 character(len=*), intent(in) :: name
1145 character(len=*), intent(in), optional :: caller
1146 
1147 !
1148 ! Local parameters
1149 !
1150 
1151 character(len=48), parameter :: sub_name = 'fm_util_get_string_array'
1152 
1153 !
1154 ! Local variables
1155 !
1156 
1157 character(len=256) :: error_header
1158 character(len=256) :: warn_header
1159 character(len=256) :: note_header
1160 character(len=128) :: caller_str
1161 character(len=32) :: index_str
1162 character(len=fm_type_name_len) :: fm_type
1163 integer :: i
1164 integer :: length
1165 
1166 nullify(array)
1167 
1168 !
1169 ! set the caller string and headers
1170 !
1171 
1172 if (present(caller)) then !{
1173  caller_str = '[' // trim(caller) // ']'
1174 else !}{
1175  caller_str = fm_util_default_caller
1176 endif !}
1177 
1178 error_header = '==>Error from ' // trim(mod_name) // &
1179  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1180 warn_header = '==>Warning from ' // trim(mod_name) // &
1181  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1182 note_header = '==>Note from ' // trim(mod_name) // &
1183  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1184 
1185 !
1186 ! check that a name is given (fatal if not)
1187 !
1188 
1189 if (name .eq. ' ') then !{
1190  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1191 endif !}
1192 
1193 fm_type = fm_get_type(name)
1194 if (fm_type .eq. 'string') then !{
1195  length = fm_get_length(name)
1196  if (length .lt. 0) then !{
1197  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1198  endif !}
1199  if (length .gt. 0) then !{
1200  allocate(array(length))
1201  do i = 1, length !{
1202  if (.not. fm_get_value(name, array(i), index = i)) then !{
1203  write (index_str,*) '(', i, ')'
1204  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name) // trim(index_str))
1205  endif !}
1206  enddo !} i
1207  endif !}
1208 elseif (fm_type .eq. ' ') then !}{
1209  call mpp_error(fatal, trim(error_header) // ' Array does not exist: ' // trim(name))
1210 else !}{
1211  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1212 endif !}
1213 
1214 return
1215 
1216 end function fm_util_get_string_array !}
1217 ! </FUNCTION> NAME="fm_util_get_string_array"
1218 
1219 
1220 !#######################################################################
1221 ! <FUNCTION NAME="fm_util_get_integer">
1222 !
1223 ! <DESCRIPTION>
1224 ! Get an integer value from the Field Manager tree.
1225 ! </DESCRIPTION>
1226 !
1227 function fm_util_get_integer(name, caller, index, default_value, scalar) &
1228  result(value) !{
1230 implicit none
1231 
1232 !
1233 ! Return type
1234 !
1235 
1236 integer :: value
1237 
1238 !
1239 ! arguments
1240 !
1241 
1242 character(len=*), intent(in) :: name
1243 character(len=*), intent(in), optional :: caller
1244 integer, intent(in), optional :: index
1245 integer, intent(in), optional :: default_value
1246 logical, intent(in), optional :: scalar
1247 
1248 !
1249 ! Local parameters
1250 !
1251 
1252 character(len=48), parameter :: sub_name = 'fm_util_get_integer'
1253 
1254 !
1255 ! Local variables
1256 !
1257 
1258 character(len=256) :: error_header
1259 character(len=256) :: warn_header
1260 character(len=256) :: note_header
1261 character(len=128) :: caller_str
1262 integer :: index_t
1263 character(len=fm_type_name_len) :: fm_type
1264 integer :: field_length
1265 
1266 !
1267 ! set the caller string and headers
1268 !
1269 
1270 if (present(caller)) then !{
1271  caller_str = '[' // trim(caller) // ']'
1272 else !}{
1273  caller_str = fm_util_default_caller
1274 endif !}
1275 
1276 error_header = '==>Error from ' // trim(mod_name) // &
1277  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1278 warn_header = '==>Warning from ' // trim(mod_name) // &
1279  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1280 note_header = '==>Note from ' // trim(mod_name) // &
1281  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1282 
1283 !
1284 ! check that a name is given (fatal if not)
1285 !
1286 
1287 if (name .eq. ' ') then !{
1288  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1289 endif !}
1290 
1291 !
1292 ! Check whether we require a scalar (length=1) and return
1293 ! an error if we do, and it isn't
1294 !
1295 
1296 if (present(scalar)) then !{
1297  if (scalar) then !{
1298  field_length = fm_get_length(name)
1299  if (field_length .lt. 0) then !{
1300  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1301  elseif (field_length .gt. 1) then !}{
1302  call mpp_error(fatal, trim(error_header) // trim(name) // ' not scalar')
1303  endif !}
1304  endif !}
1305 endif !}
1306 
1307 !
1308 ! set the index
1309 !
1310 
1311 if (present(index)) then !{
1312  index_t = index
1313  if (index .le. 0) then !{
1314  call mpp_error(fatal, trim(error_header) // ' Index not positive')
1315  endif !}
1316 else !}{
1317  index_t = 1
1318 endif !}
1319 
1320 fm_type = fm_get_type(name)
1321 if (fm_type .eq. 'integer') then !{
1322  if (.not. fm_get_value(name, value, index = index_t)) then !{
1323  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name))
1324  endif !}
1325 elseif (fm_type .eq. ' ' .and. present(default_value)) then !}{
1326  value = default_value
1327 elseif (fm_type .eq. ' ') then !}{
1328  call mpp_error(fatal, trim(error_header) // ' Field does not exist: ' // trim(name))
1329 else !}{
1330  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1331 endif !}
1332 
1333 return
1334 
1335 end function fm_util_get_integer !}
1336 ! </FUNCTION> NAME="fm_util_get_integer"
1337 
1338 
1339 !#######################################################################
1340 ! <FUNCTION NAME="fm_util_get_logical">
1341 !
1342 ! <DESCRIPTION>
1343 ! Get a logical value from the Field Manager tree.
1344 ! </DESCRIPTION>
1345 !
1346 function fm_util_get_logical(name, caller, index, default_value, scalar) &
1347  result(value) !{
1349 implicit none
1350 
1351 !
1352 ! Return type
1353 !
1354 
1355 logical :: value
1356 
1357 !
1358 ! arguments
1359 !
1360 
1361 character(len=*), intent(in) :: name
1362 character(len=*), intent(in), optional :: caller
1363 integer, intent(in), optional :: index
1364 logical, intent(in), optional :: default_value
1365 logical, intent(in), optional :: scalar
1366 
1367 !
1368 ! Local parameters
1369 !
1370 
1371 character(len=48), parameter :: sub_name = 'fm_util_get_logical'
1372 
1373 !
1374 ! Local variables
1375 !
1376 
1377 character(len=256) :: error_header
1378 character(len=256) :: warn_header
1379 character(len=256) :: note_header
1380 character(len=128) :: caller_str
1381 integer :: index_t
1382 character(len=fm_type_name_len) :: fm_type
1383 integer :: field_length
1384 
1385 !
1386 ! set the caller string and headers
1387 !
1388 
1389 if (present(caller)) then !{
1390  caller_str = '[' // trim(caller) // ']'
1391 else !}{
1392  caller_str = fm_util_default_caller
1393 endif !}
1394 
1395 error_header = '==>Error from ' // trim(mod_name) // &
1396  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1397 warn_header = '==>Warning from ' // trim(mod_name) // &
1398  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1399 note_header = '==>Note from ' // trim(mod_name) // &
1400  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1401 
1402 !
1403 ! check that a name is given (fatal if not)
1404 !
1405 
1406 if (name .eq. ' ') then !{
1407  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1408 endif !}
1409 
1410 !
1411 ! Check whether we require a scalar (length=1) and return
1412 ! an error if we do, and it isn't
1413 !
1414 
1415 if (present(scalar)) then !{
1416  if (scalar) then !{
1417  field_length = fm_get_length(name)
1418  if (field_length .lt. 0) then !{
1419  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1420  elseif (field_length .gt. 1) then !}{
1421  call mpp_error(fatal, trim(error_header) // trim(name) // ' not scalar')
1422  endif !}
1423  endif !}
1424 endif !}
1425 
1426 !
1427 ! set the index
1428 !
1429 
1430 if (present(index)) then !{
1431  index_t = index
1432  if (index .le. 0) then !{
1433  call mpp_error(fatal, trim(error_header) // ' Index not positive')
1434  endif !}
1435 else !}{
1436  index_t = 1
1437 endif !}
1438 
1439 fm_type = fm_get_type(name)
1440 if (fm_type .eq. 'logical') then !{
1441  if (.not. fm_get_value(name, value, index = index_t)) then !{
1442  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name))
1443  endif !}
1444 elseif (fm_type .eq. ' ' .and. present(default_value)) then !}{
1445  value = default_value
1446 elseif (fm_type .eq. ' ') then !}{
1447  call mpp_error(fatal, trim(error_header) // ' Field does not exist: ' // trim(name))
1448 else !}{
1449  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1450 endif !}
1451 
1452 return
1453 
1454 end function fm_util_get_logical !}
1455 ! </FUNCTION> NAME="fm_util_get_logical"
1456 
1457 
1458 !#######################################################################
1459 ! <FUNCTION NAME="fm_util_get_real">
1460 !
1461 ! <DESCRIPTION>
1462 ! Get a real value from the Field Manager tree.
1463 ! </DESCRIPTION>
1464 !
1465 function fm_util_get_real(name, caller, index, default_value, scalar) &
1466  result(value) !{
1468 implicit none
1469 
1470 !
1471 ! Return type
1472 !
1473 
1474 real :: value
1475 
1476 !
1477 ! arguments
1478 !
1479 
1480 character(len=*), intent(in) :: name
1481 character(len=*), intent(in), optional :: caller
1482 integer, intent(in), optional :: index
1483 real, intent(in), optional :: default_value
1484 logical, intent(in), optional :: scalar
1485 
1486 !
1487 ! Local parameters
1488 !
1489 
1490 character(len=48), parameter :: sub_name = 'fm_util_get_real'
1491 
1492 !
1493 ! Local variables
1494 !
1495 
1496 character(len=256) :: error_header
1497 character(len=256) :: warn_header
1498 character(len=256) :: note_header
1499 character(len=128) :: caller_str
1500 integer :: index_t
1501 character(len=fm_type_name_len) :: fm_type
1502 integer :: field_length
1503 integer :: ivalue
1504 
1505 !
1506 ! set the caller string and headers
1507 !
1508 
1509 if (present(caller)) then !{
1510  caller_str = '[' // trim(caller) // ']'
1511 else !}{
1512  caller_str = fm_util_default_caller
1513 endif !}
1514 
1515 error_header = '==>Error from ' // trim(mod_name) // &
1516  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1517 warn_header = '==>Warning from ' // trim(mod_name) // &
1518  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1519 note_header = '==>Note from ' // trim(mod_name) // &
1520  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1521 
1522 !
1523 ! check that a name is given (fatal if not)
1524 !
1525 
1526 if (name .eq. ' ') then !{
1527  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1528 endif !}
1529 
1530 !
1531 ! Check whether we require a scalar (length=1) and return
1532 ! an error if we do, and it isn't
1533 !
1534 
1535 if (present(scalar)) then !{
1536  if (scalar) then !{
1537  field_length = fm_get_length(name)
1538  if (field_length .lt. 0) then !{
1539  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1540  elseif (field_length .gt. 1) then !}{
1541  call mpp_error(fatal, trim(error_header) // trim(name) // ' not scalar')
1542  endif !}
1543  endif !}
1544 endif !}
1545 
1546 !
1547 ! set the index
1548 !
1549 
1550 if (present(index)) then !{
1551  index_t = index
1552  if (index .le. 0) then !{
1553  call mpp_error(fatal, trim(error_header) // ' Index not positive')
1554  endif !}
1555 else !}{
1556  index_t = 1
1557 endif !}
1558 
1559 fm_type = fm_get_type(name)
1560 if (fm_type .eq. 'real') then !{
1561  if (.not. fm_get_value(name, value, index = index_t)) then !{
1562  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name))
1563  endif !}
1564 else if (fm_type .eq. 'integer') then
1565  if (.not. fm_get_value(name, ivalue, index = index_t)) then
1566  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name))
1567  endif
1568  value = ivalue
1569 elseif (fm_type .eq. ' ' .and. present(default_value)) then !}{
1570  value = default_value
1571 elseif (fm_type .eq. ' ') then !}{
1572  call mpp_error(fatal, trim(error_header) // ' Field does not exist: ' // trim(name))
1573 else !}{
1574  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1575 endif !}
1576 
1577 return
1578 
1579 end function fm_util_get_real !}
1580 ! </FUNCTION> NAME="fm_util_get_real"
1581 
1582 
1583 !#######################################################################
1584 ! <FUNCTION NAME="fm_util_get_string">
1585 !
1586 ! <DESCRIPTION>
1587 ! Get a string value from the Field Manager tree.
1588 ! </DESCRIPTION>
1589 !
1590 function fm_util_get_string(name, caller, index, default_value, scalar) &
1591  result(value) !{
1593 implicit none
1594 
1595 !
1596 ! Return type
1597 !
1598 
1599 character(len=fm_string_len) :: value
1600 
1601 !
1602 ! arguments
1603 !
1604 
1605 character(len=*), intent(in) :: name
1606 character(len=*), intent(in), optional :: caller
1607 integer, intent(in), optional :: index
1608 character(len=*), intent(in), optional :: default_value
1609 logical, intent(in), optional :: scalar
1610 
1611 !
1612 ! Local parameters
1613 !
1614 
1615 character(len=48), parameter :: sub_name = 'fm_util_get_string'
1616 
1617 !
1618 ! Local variables
1619 !
1620 
1621 character(len=256) :: error_header
1622 character(len=256) :: warn_header
1623 character(len=256) :: note_header
1624 character(len=128) :: caller_str
1625 integer :: index_t
1626 character(len=fm_type_name_len) :: fm_type
1627 integer :: field_length
1628 
1629 !
1630 ! set the caller string and headers
1631 !
1632 
1633 if (present(caller)) then !{
1634  caller_str = '[' // trim(caller) // ']'
1635 else !}{
1636  caller_str = fm_util_default_caller
1637 endif !}
1638 
1639 error_header = '==>Error from ' // trim(mod_name) // &
1640  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1641 warn_header = '==>Warning from ' // trim(mod_name) // &
1642  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1643 note_header = '==>Note from ' // trim(mod_name) // &
1644  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1645 
1646 !
1647 ! check that a name is given (fatal if not)
1648 !
1649 
1650 if (name .eq. ' ') then !{
1651  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1652 endif !}
1653 
1654 !
1655 ! Check whether we require a scalar (length=1) and return
1656 ! an error if we do, and it isn't
1657 !
1658 
1659 if (present(scalar)) then !{
1660  if (scalar) then !{
1661  field_length = fm_get_length(name)
1662  if (field_length .lt. 0) then !{
1663  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1664  elseif (field_length .gt. 1) then !}{
1665  call mpp_error(fatal, trim(error_header) // trim(name) // ' not scalar')
1666  endif !}
1667  endif !}
1668 endif !}
1669 
1670 !
1671 ! set the index
1672 !
1673 
1674 if (present(index)) then !{
1675  index_t = index
1676  if (index .le. 0) then !{
1677  call mpp_error(fatal, trim(error_header) // ' Index not positive')
1678  endif !}
1679 else !}{
1680  index_t = 1
1681 endif !}
1682 
1683 fm_type = fm_get_type(name)
1684 if (fm_type .eq. 'string') then !{
1685  if (.not. fm_get_value(name, value, index = index_t)) then !{
1686  call mpp_error(fatal, trim(error_header) // ' Problem getting ' // trim(name))
1687  endif !}
1688 elseif (fm_type .eq. ' ' .and. present(default_value)) then !}{
1689  value = default_value
1690 elseif (fm_type .eq. ' ') then !}{
1691  call mpp_error(fatal, trim(error_header) // ' Field does not exist: ' // trim(name))
1692 else !}{
1693  call mpp_error(fatal, trim(error_header) // ' Wrong type for ' // trim(name) // ', found (' // trim(fm_type) // ')')
1694 endif !}
1695 
1696 return
1697 
1698 end function fm_util_get_string !}
1699 ! </FUNCTION> NAME="fm_util_get_string"
1700 
1701 
1702 !#######################################################################
1703 ! <SUBROUTINE NAME="fm_util_set_value_integer_array">
1704 !
1705 ! <DESCRIPTION>
1706 ! Set an integer array in the Field Manager tree.
1707 ! </DESCRIPTION>
1708 !
1709 
1710 subroutine fm_util_set_value_integer_array(name, value, length, caller, no_overwrite, good_name_list) !{
1712 implicit none
1713 
1714 !
1715 ! arguments
1716 !
1717 
1718 character(len=*), intent(in) :: name
1719 integer, intent(in) :: length
1720 integer, intent(in) :: value(length)
1721 character(len=*), intent(in), optional :: caller
1722 logical, intent(in), optional :: no_overwrite
1723 character(len=fm_path_name_len), intent(in), optional :: good_name_list
1724 
1725 !
1726 ! Local parameters
1727 !
1728 
1729 character(len=48), parameter :: sub_name = 'fm_util_set_value_integer_array'
1730 
1731 !
1732 ! Local variables
1733 !
1734 
1735 character(len=256) :: error_header
1736 character(len=256) :: warn_header
1737 character(len=256) :: note_header
1738 character(len=128) :: caller_str
1739 character(len=32) :: str_error
1740 integer :: field_index
1741 integer :: field_length
1742 integer :: n
1743 logical :: no_overwrite_use
1744 character(len=fm_path_name_len) :: good_name_list_use
1745 logical :: add_name
1746 
1747 !
1748 ! set the caller string and headers
1749 !
1750 
1751 if (present(caller)) then !{
1752  caller_str = '[' // trim(caller) // ']'
1753 else !}{
1754  caller_str = fm_util_default_caller
1755 endif !}
1756 
1757 error_header = '==>Error from ' // trim(mod_name) // &
1758  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1759 warn_header = '==>Warning from ' // trim(mod_name) // &
1760  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1761 note_header = '==>Note from ' // trim(mod_name) // &
1762  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1763 
1764 !
1765 ! check that a name is given (fatal if not)
1766 !
1767 
1768 if (name .eq. ' ') then !{
1769  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1770 endif !}
1771 
1772 !
1773 ! check that the length is non-negative
1774 !
1775 
1776 if (length .lt. 0) then !{
1777  call mpp_error(fatal, trim(error_header) // ' Negative array length')
1778 endif !}
1779 
1780 !
1781 ! check for whether to overwrite existing values
1782 !
1783 
1784 if (present(no_overwrite)) then !{
1785  no_overwrite_use = no_overwrite
1786 else !}{
1787  no_overwrite_use = default_no_overwrite
1788 endif !}
1789 
1790 !
1791 ! check for whether to save the name in a list
1792 !
1793 
1794 if (present(good_name_list)) then !{
1795  good_name_list_use = good_name_list
1796 else !}{
1797  good_name_list_use = default_good_name_list
1798 endif !}
1799 
1800 !
1801 ! write the data array
1802 !
1803 
1804 if (length .eq. 0) then !{
1805  if (.not. (no_overwrite_use .and. fm_exists(name))) then !{
1806  field_index = fm_new_value(name, 0, index = 0)
1807  if (field_index .le. 0) then !{
1808  write (str_error,*) ' with length = ', length
1809  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
1810  endif !}
1811  endif !}
1812 else !}{
1813  if (no_overwrite_use .and. fm_exists(name)) then !{
1814  field_length = fm_get_length(name)
1815  if (field_length .lt. 0) then !{
1816  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1817  endif !}
1818  do n = field_length + 1, length !{
1819  field_index = fm_new_value(name, value(n), index = n)
1820  if (field_index .le. 0) then !{
1821  write (str_error,*) ' with index = ', n
1822  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
1823  endif !}
1824  enddo !} n
1825  else !}{
1826  field_index = fm_new_value(name, value(1))
1827  if (field_index .le. 0) then !{
1828  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name))
1829  endif !}
1830  do n = 2, length !{
1831  field_index = fm_new_value(name, value(n), index = n)
1832  if (field_index .le. 0) then !{
1833  write (str_error,*) ' with index = ', n
1834  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
1835  endif !}
1836  enddo !} n
1837  endif !}
1838 endif !}
1839 
1840 !
1841 ! Add the variable name to the list of good names, to be used
1842 ! later for a consistency check
1843 !
1844 
1845 if (good_name_list_use .ne. ' ') then !{
1846  if (fm_exists(good_name_list_use)) then !{
1847  add_name = fm_util_get_index_string(good_name_list_use, name, &
1848  caller = caller_str) .le. 0 ! true if name does not exist in string array
1849  else !}{
1850  add_name = .true. ! always add to new list
1851  endif !}
1852  if (add_name .and. fm_exists(name)) then !{
1853  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
1854  call mpp_error(fatal, trim(error_header) // &
1855  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
1856  endif !}
1857  endif !}
1858 endif !}
1859 
1860 return
1861 
1862 end subroutine fm_util_set_value_integer_array !}
1863 ! </SUBROUTINE> NAME="fm_util_set_value_integer_array"
1864 
1865 
1866 !#######################################################################
1867 ! <SUBROUTINE NAME="fm_util_set_value_logical_array">
1868 !
1869 ! <DESCRIPTION>
1870 ! Set a logical array in the Field Manager tree.
1871 ! </DESCRIPTION>
1872 !
1873 
1874 subroutine fm_util_set_value_logical_array(name, value, length, caller, no_overwrite, good_name_list) !{
1876 implicit none
1877 
1878 !
1879 ! arguments
1880 !
1881 
1882 character(len=*), intent(in) :: name
1883 integer, intent(in) :: length
1884 logical, intent(in) :: value(length)
1885 character(len=*), intent(in), optional :: caller
1886 logical, intent(in), optional :: no_overwrite
1887 character(len=fm_path_name_len), intent(in), optional :: good_name_list
1888 
1889 !
1890 ! Local parameters
1891 !
1892 
1893 character(len=48), parameter :: sub_name = 'fm_util_set_value_logical_array'
1894 
1895 !
1896 ! Local variables
1897 !
1898 
1899 character(len=256) :: error_header
1900 character(len=256) :: warn_header
1901 character(len=256) :: note_header
1902 character(len=128) :: caller_str
1903 character(len=32) :: str_error
1904 integer :: field_index
1905 integer :: field_length
1906 integer :: n
1907 logical :: no_overwrite_use
1908 character(len=fm_path_name_len) :: good_name_list_use
1909 logical :: add_name
1910 
1911 !
1912 ! set the caller string and headers
1913 !
1914 
1915 if (present(caller)) then !{
1916  caller_str = '[' // trim(caller) // ']'
1917 else !}{
1918  caller_str = fm_util_default_caller
1919 endif !}
1920 
1921 error_header = '==>Error from ' // trim(mod_name) // &
1922  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1923 warn_header = '==>Warning from ' // trim(mod_name) // &
1924  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1925 note_header = '==>Note from ' // trim(mod_name) // &
1926  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
1927 
1928 !
1929 ! check that a name is given (fatal if not)
1930 !
1931 
1932 if (name .eq. ' ') then !{
1933  call mpp_error(fatal, trim(error_header) // ' Empty name given')
1934 endif !}
1935 
1936 !
1937 ! check that the length is non-negative
1938 !
1939 
1940 if (length .lt. 0) then !{
1941  call mpp_error(fatal, trim(error_header) // ' Negative array length')
1942 endif !}
1943 
1944 !
1945 ! check for whether to overwrite existing values
1946 !
1947 
1948 if (present(no_overwrite)) then !{
1949  no_overwrite_use = no_overwrite
1950 else !}{
1951  no_overwrite_use = default_no_overwrite
1952 endif !}
1953 
1954 !
1955 ! check for whether to save the name in a list
1956 !
1957 
1958 if (present(good_name_list)) then !{
1959  good_name_list_use = good_name_list
1960 else !}{
1961  good_name_list_use = default_good_name_list
1962 endif !}
1963 
1964 !
1965 ! write the data array
1966 !
1967 
1968 if (length .eq. 0) then !{
1969  if (.not. (no_overwrite_use .and. fm_exists(name))) then !{
1970  field_index = fm_new_value(name, .false., index = 0)
1971  if (field_index .le. 0) then !{
1972  write (str_error,*) ' with length = ', length
1973  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
1974  endif !}
1975  endif !}
1976 else !}{
1977  if (no_overwrite_use .and. fm_exists(name)) then !{
1978  field_length = fm_get_length(name)
1979  if (field_length .lt. 0) then !{
1980  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
1981  endif !}
1982  do n = field_length + 1, length !{
1983  field_index = fm_new_value(name, value(n), index = n)
1984  if (field_index .le. 0) then !{
1985  write (str_error,*) ' with index = ', n
1986  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
1987  endif !}
1988  enddo !} n
1989  else !}{
1990  field_index = fm_new_value(name, value(1))
1991  if (field_index .le. 0) then !{
1992  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name))
1993  endif !}
1994  do n = 2, length !{
1995  field_index = fm_new_value(name, value(n), index = n)
1996  if (field_index .le. 0) then !{
1997  write (str_error,*) ' with index = ', n
1998  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
1999  endif !}
2000  enddo !} n
2001  endif !}
2002 endif !}
2003 
2004 !
2005 ! Add the variable name to the list of good names, to be used
2006 ! later for a consistency check
2007 !
2008 
2009 if (good_name_list_use .ne. ' ') then !{
2010  if (fm_exists(good_name_list_use)) then !{
2011  add_name = fm_util_get_index_string(good_name_list_use, name, &
2012  caller = caller_str) .le. 0 ! true if name does not exist in string array
2013  else !}{
2014  add_name = .true. ! always add to new list
2015  endif !}
2016  if (add_name .and. fm_exists(name)) then !{
2017  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
2018  call mpp_error(fatal, trim(error_header) // &
2019  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
2020  endif !}
2021  endif !}
2022 endif !}
2023 
2024 return
2025 
2026 end subroutine fm_util_set_value_logical_array !}
2027 ! </SUBROUTINE> NAME="fm_util_set_value_logical_array"
2028 
2029 
2030 !#######################################################################
2031 ! <SUBROUTINE NAME="fm_util_set_value_real_array">
2032 !
2033 ! <DESCRIPTION>
2034 ! Set a real array in the Field Manager tree.
2035 ! </DESCRIPTION>
2036 !
2037 
2038 subroutine fm_util_set_value_real_array(name, value, length, caller, no_overwrite, good_name_list) !{
2040 implicit none
2041 
2042 !
2043 ! arguments
2044 !
2045 
2046 character(len=*), intent(in) :: name
2047 integer, intent(in) :: length
2048 real, intent(in) :: value(length)
2049 character(len=*), intent(in), optional :: caller
2050 logical, intent(in), optional :: no_overwrite
2051 character(len=fm_path_name_len), intent(in), optional :: good_name_list
2052 
2053 !
2054 ! Local parameters
2055 !
2056 
2057 character(len=48), parameter :: sub_name = 'fm_util_set_value_real_array'
2058 
2059 !
2060 ! Local variables
2061 !
2062 
2063 character(len=256) :: error_header
2064 character(len=256) :: warn_header
2065 character(len=256) :: note_header
2066 character(len=128) :: caller_str
2067 character(len=32) :: str_error
2068 integer :: field_index
2069 integer :: field_length
2070 integer :: n
2071 logical :: no_overwrite_use
2072 character(len=fm_path_name_len) :: good_name_list_use
2073 logical :: add_name
2074 
2075 !
2076 ! set the caller string and headers
2077 !
2078 
2079 if (present(caller)) then !{
2080  caller_str = '[' // trim(caller) // ']'
2081 else !}{
2082  caller_str = fm_util_default_caller
2083 endif !}
2084 
2085 error_header = '==>Error from ' // trim(mod_name) // &
2086  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2087 warn_header = '==>Warning from ' // trim(mod_name) // &
2088  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2089 note_header = '==>Note from ' // trim(mod_name) // &
2090  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2091 
2092 !
2093 ! check that a name is given (fatal if not)
2094 !
2095 
2096 if (name .eq. ' ') then !{
2097  call mpp_error(fatal, trim(error_header) // ' Empty name given')
2098 endif !}
2099 
2100 !
2101 ! check that the length is non-negative
2102 !
2103 
2104 if (length .lt. 0) then !{
2105  call mpp_error(fatal, trim(error_header) // ' Negative array length')
2106 endif !}
2107 
2108 !
2109 ! check for whether to overwrite existing values
2110 !
2111 
2112 if (present(no_overwrite)) then !{
2113  no_overwrite_use = no_overwrite
2114 else !}{
2115  no_overwrite_use = default_no_overwrite
2116 endif !}
2117 
2118 !
2119 ! check for whether to save the name in a list
2120 !
2121 
2122 if (present(good_name_list)) then !{
2123  good_name_list_use = good_name_list
2124 else !}{
2125  good_name_list_use = default_good_name_list
2126 endif !}
2127 
2128 !
2129 ! write the data array
2130 !
2131 
2132 if (length .eq. 0) then !{
2133  if (.not. (no_overwrite_use .and. fm_exists(name))) then !{
2134  field_index = fm_new_value(name, 0.0, index = 0)
2135  if (field_index .le. 0) then !{
2136  write (str_error,*) ' with length = ', length
2137  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2138  endif !}
2139  endif !}
2140 else !}{
2141  if (no_overwrite_use .and. fm_exists(name)) then !{
2142  field_length = fm_get_length(name)
2143  if (field_length .lt. 0) then !{
2144  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
2145  endif !}
2146  do n = field_length + 1, length !{
2147  field_index = fm_new_value(name, value(n), index = n)
2148  if (field_index .le. 0) then !{
2149  write (str_error,*) ' with index = ', n
2150  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2151  endif !}
2152  enddo !} n
2153  else !}{
2154  field_index = fm_new_value(name, value(1))
2155  if (field_index .le. 0) then !{
2156  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name))
2157  endif !}
2158  do n = 2, length !{
2159  field_index = fm_new_value(name, value(n), index = n)
2160  if (field_index .le. 0) then !{
2161  write (str_error,*) ' with index = ', n
2162  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2163  endif !}
2164  enddo !} n
2165  endif !}
2166 endif !}
2167 
2168 !
2169 ! Add the variable name to the list of good names, to be used
2170 ! later for a consistency check
2171 !
2172 
2173 if (good_name_list_use .ne. ' ') then !{
2174  if (fm_exists(good_name_list_use)) then !{
2175  add_name = fm_util_get_index_string(good_name_list_use, name, &
2176  caller = caller_str) .le. 0 ! true if name does not exist in string array
2177  else !}{
2178  add_name = .true. ! always add to new list
2179  endif !}
2180  if (add_name .and. fm_exists(name)) then !{
2181  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
2182  call mpp_error(fatal, trim(error_header) // &
2183  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
2184  endif !}
2185  endif !}
2186 endif !}
2187 
2188 return
2189 
2190 end subroutine fm_util_set_value_real_array !}
2191 ! </SUBROUTINE> NAME="fm_util_set_value_real_array"
2192 
2193 
2194 !#######################################################################
2195 ! <SUBROUTINE NAME="fm_util_set_value_string_array">
2196 !
2197 ! <DESCRIPTION>
2198 ! Set a string array in the Field Manager tree.
2199 ! </DESCRIPTION>
2200 !
2201 
2202 subroutine fm_util_set_value_string_array(name, value, length, caller, no_overwrite, good_name_list) !{
2204 implicit none
2205 
2206 !
2207 ! arguments
2208 !
2209 
2210 character(len=*), intent(in) :: name
2211 integer, intent(in) :: length
2212 character(len=*), intent(in) :: value(length)
2213 character(len=*), intent(in), optional :: caller
2214 logical, intent(in), optional :: no_overwrite
2215 character(len=fm_path_name_len), intent(in), optional :: good_name_list
2216 
2217 !
2218 ! Local parameters
2219 !
2220 
2221 character(len=48), parameter :: sub_name = 'fm_util_set_value_string_array'
2222 
2223 !
2224 ! Local variables
2225 !
2226 
2227 character(len=256) :: error_header
2228 character(len=256) :: warn_header
2229 character(len=256) :: note_header
2230 character(len=128) :: caller_str
2231 character(len=32) :: str_error
2232 integer :: field_index
2233 integer :: field_length
2234 integer :: n
2235 logical :: no_overwrite_use
2236 character(len=fm_path_name_len) :: good_name_list_use
2237 logical :: add_name
2238 
2239 !
2240 ! set the caller string and headers
2241 !
2242 
2243 if (present(caller)) then !{
2244  caller_str = '[' // trim(caller) // ']'
2245 else !}{
2246  caller_str = fm_util_default_caller
2247 endif !}
2248 
2249 error_header = '==>Error from ' // trim(mod_name) // &
2250  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2251 warn_header = '==>Warning from ' // trim(mod_name) // &
2252  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2253 note_header = '==>Note from ' // trim(mod_name) // &
2254  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2255 
2256 !
2257 ! check that a name is given (fatal if not)
2258 !
2259 
2260 if (name .eq. ' ') then !{
2261  call mpp_error(fatal, trim(error_header) // ' Empty name given')
2262 endif !}
2263 
2264 !
2265 ! check that the length is non-negative
2266 !
2267 
2268 if (length .lt. 0) then !{
2269  call mpp_error(fatal, trim(error_header) // ' Negative array length')
2270 endif !}
2271 
2272 !
2273 ! check for whether to overwrite existing values
2274 !
2275 
2276 if (present(no_overwrite)) then !{
2277  no_overwrite_use = no_overwrite
2278 else !}{
2279  no_overwrite_use = default_no_overwrite
2280 endif !}
2281 
2282 !
2283 ! check for whether to save the name in a list
2284 !
2285 
2286 if (present(good_name_list)) then !{
2287  good_name_list_use = good_name_list
2288 else !}{
2289  good_name_list_use = default_good_name_list
2290 endif !}
2291 
2292 !
2293 ! write the data array
2294 !
2295 
2296 if (length .eq. 0) then !{
2297  if (.not. (no_overwrite_use .and. fm_exists(name))) then !{
2298  field_index = fm_new_value(name, ' ', index = 0)
2299  if (field_index .le. 0) then !{
2300  write (str_error,*) ' with length = ', length
2301  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2302  endif !}
2303  endif !}
2304 else !}{
2305  if (no_overwrite_use .and. fm_exists(name)) then !{
2306  field_length = fm_get_length(name)
2307  if (field_length .lt. 0) then !{
2308  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
2309  endif !}
2310  do n = field_length + 1, length !{
2311  field_index = fm_new_value(name, value(n), index = n)
2312  if (field_index .le. 0) then !{
2313  write (str_error,*) ' with index = ', n
2314  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2315  endif !}
2316  enddo !} n
2317  else !}{
2318  field_index = fm_new_value(name, value(1))
2319  if (field_index .le. 0) then !{
2320  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name))
2321  endif !}
2322  do n = 2, length !{
2323  field_index = fm_new_value(name, value(n), index = n)
2324  if (field_index .le. 0) then !{
2325  write (str_error,*) ' with index = ', n
2326  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2327  endif !}
2328  enddo !} n
2329  endif !}
2330 endif !}
2331 
2332 !
2333 ! Add the variable name to the list of good names, to be used
2334 ! later for a consistency check
2335 !
2336 
2337 if (good_name_list_use .ne. ' ') then !{
2338  if (fm_exists(good_name_list_use)) then !{
2339  add_name = fm_util_get_index_string(good_name_list_use, name, &
2340  caller = caller_str) .le. 0 ! true if name does not exist in string array
2341  else !}{
2342  add_name = .true. ! always add to new list
2343  endif !}
2344  if (add_name .and. fm_exists(name)) then !{
2345  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
2346  call mpp_error(fatal, trim(error_header) // &
2347  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
2348  endif !}
2349  endif !}
2350 endif !}
2351 
2352 return
2353 
2354 end subroutine fm_util_set_value_string_array !}
2355 ! </SUBROUTINE> NAME="fm_util_set_value_string_array"
2356 
2357 
2358 !#######################################################################
2359 ! <SUBROUTINE NAME="fm_util_set_value_integer">
2360 !
2361 ! <DESCRIPTION>
2362 ! Set an integer value in the Field Manager tree.
2363 ! </DESCRIPTION>
2364 !
2365 
2366 subroutine fm_util_set_value_integer(name, value, caller, index, append, no_create, &
2367  no_overwrite, good_name_list) !{
2369 implicit none
2370 
2371 !
2372 ! arguments
2373 !
2374 
2375 character(len=*), intent(in) :: name
2376 integer, intent(in) :: value
2377 character(len=*), intent(in), optional :: caller
2378 integer, intent(in), optional :: index
2379 logical, intent(in), optional :: append
2380 logical, intent(in), optional :: no_create
2381 logical, intent(in), optional :: no_overwrite
2382 character(len=*), intent(in), optional :: good_name_list
2383 
2384 !
2385 ! Local parameters
2386 !
2387 
2388 character(len=48), parameter :: sub_name = 'fm_util_set_value_integer'
2389 
2390 !
2391 ! Local variables
2392 !
2393 
2394 character(len=256) :: error_header
2395 character(len=256) :: warn_header
2396 character(len=256) :: note_header
2397 character(len=128) :: caller_str
2398 character(len=32) :: str_error
2399 integer :: field_index
2400 logical :: no_overwrite_use
2401 integer :: field_length
2402 character(len=fm_path_name_len) :: good_name_list_use
2403 logical :: create
2404 logical :: add_name
2405 
2406 !
2407 ! set the caller string and headers
2408 !
2409 
2410 if (present(caller)) then !{
2411  caller_str = '[' // trim(caller) // ']'
2412 else !}{
2413  caller_str = fm_util_default_caller
2414 endif !}
2415 
2416 error_header = '==>Error from ' // trim(mod_name) // &
2417  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2418 warn_header = '==>Warning from ' // trim(mod_name) // &
2419  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2420 note_header = '==>Note from ' // trim(mod_name) // &
2421  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2422 
2423 !
2424 ! check that a name is given (fatal if not)
2425 !
2426 
2427 if (name .eq. ' ') then !{
2428  call mpp_error(fatal, trim(error_header) // ' Empty name given')
2429 endif !}
2430 
2431 !
2432 ! check that append and index are not both given
2433 !
2434 
2435 if (present(index) .and. present(append)) then !{
2436  call mpp_error(fatal, trim(error_header) // ' Append and index both given as arguments')
2437 endif !}
2438 
2439 !
2440 ! check for whether to overwrite existing values
2441 !
2442 
2443 if (present(no_overwrite)) then !{
2444  no_overwrite_use = no_overwrite
2445 else !}{
2446  no_overwrite_use = default_no_overwrite
2447 endif !}
2448 
2449 !
2450 ! check for whether to save the name in a list
2451 !
2452 
2453 if (present(good_name_list)) then !{
2454  good_name_list_use = good_name_list
2455 else !}{
2456  good_name_list_use = default_good_name_list
2457 endif !}
2458 
2459 if (present(no_create)) then !{
2460  create = .not. no_create
2461  if (no_create .and. (present(append) .or. present(index))) then !{
2462  call mpp_error(fatal, trim(error_header) // ' append or index are present when no_create is true for ' // trim(name))
2463  endif !}
2464 else !}{
2465  create = .true.
2466 endif !}
2467 
2468 if (present(index)) then !{
2469  if (fm_exists(name)) then !{
2470  field_length = fm_get_length(name)
2471  if (field_length .lt. 0) then !{
2472  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
2473  endif !}
2474  if (.not. (no_overwrite_use .and. field_length .ge. index)) then !{
2475  field_index = fm_new_value(name, value, index = index)
2476  if (field_index .le. 0) then !{
2477  write (str_error,*) ' with index = ', index
2478  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name) // trim(str_error))
2479  endif !}
2480  endif !}
2481  else !}{
2482  field_index = fm_new_value(name, value, index = index)
2483  if (field_index .le. 0) then !{
2484  write (str_error,*) ' with index = ', index
2485  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2486  endif !}
2487  endif !}
2488 elseif (present(append)) then !}{
2489  field_index = fm_new_value(name, value, append = append)
2490  if (field_index .le. 0) then !{
2491  write (str_error,*) ' with append = ', append
2492  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2493  endif !}
2494 else !}{
2495  if (fm_exists(name)) then !{
2496  if (.not. no_overwrite_use) then !{
2497  field_index = fm_new_value(name, value)
2498  if (field_index .le. 0) then !{
2499  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name))
2500  endif !}
2501  endif !}
2502  elseif (create) then !}{
2503  field_index = fm_new_value(name, value)
2504  if (field_index .le. 0) then !{
2505  call mpp_error(fatal, trim(error_header) // ' Problem creating ' // trim(name))
2506  endif !}
2507  endif !}
2508 endif !}
2509 
2510 !
2511 ! Add the variable name to the list of good names, to be used
2512 ! later for a consistency check, unless the field did not exist and we did not create it
2513 !
2514 
2515 if (good_name_list_use .ne. ' ') then !{
2516  if (fm_exists(good_name_list_use)) then !{
2517  add_name = fm_util_get_index_string(good_name_list_use, name, &
2518  caller = caller_str) .le. 0 ! true if name does not exist in string array
2519  else !}{
2520  add_name = .true. ! always add to new list
2521  endif !}
2522  if (add_name .and. fm_exists(name)) then !{
2523  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
2524  call mpp_error(fatal, trim(error_header) // &
2525  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
2526  endif !}
2527  endif !}
2528 endif !}
2529 
2530 return
2531 
2532 end subroutine fm_util_set_value_integer !}
2533 ! </SUBROUTINE> NAME="fm_util_set_value_integer"
2534 
2535 
2536 !#######################################################################
2537 ! <SUBROUTINE NAME="fm_util_set_value_logical">
2538 !
2539 ! <DESCRIPTION>
2540 ! Set a logical value in the Field Manager tree.
2541 ! </DESCRIPTION>
2542 !
2543 
2544 subroutine fm_util_set_value_logical(name, value, caller, index, append, no_create, &
2545  no_overwrite, good_name_list) !{
2547 implicit none
2548 
2549 !
2550 ! arguments
2551 !
2552 
2553 character(len=*), intent(in) :: name
2554 logical, intent(in) :: value
2555 character(len=*), intent(in), optional :: caller
2556 integer, intent(in), optional :: index
2557 logical, intent(in), optional :: append
2558 logical, intent(in), optional :: no_create
2559 logical, intent(in), optional :: no_overwrite
2560 character(len=*), intent(in), optional :: good_name_list
2561 
2562 !
2563 ! Local parameters
2564 !
2565 
2566 character(len=48), parameter :: sub_name = 'fm_util_set_value_logical'
2567 
2568 !
2569 ! Local variables
2570 !
2571 
2572 character(len=256) :: error_header
2573 character(len=256) :: warn_header
2574 character(len=256) :: note_header
2575 character(len=128) :: caller_str
2576 character(len=32) :: str_error
2577 integer :: field_index
2578 logical :: no_overwrite_use
2579 integer :: field_length
2580 character(len=fm_path_name_len) :: good_name_list_use
2581 logical :: create
2582 logical :: add_name
2583 
2584 !
2585 ! set the caller string and headers
2586 !
2587 
2588 if (present(caller)) then !{
2589  caller_str = '[' // trim(caller) // ']'
2590 else !}{
2591  caller_str = fm_util_default_caller
2592 endif !}
2593 
2594 error_header = '==>Error from ' // trim(mod_name) // &
2595  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2596 warn_header = '==>Warning from ' // trim(mod_name) // &
2597  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2598 note_header = '==>Note from ' // trim(mod_name) // &
2599  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2600 
2601 !
2602 ! check that a name is given (fatal if not)
2603 !
2604 
2605 if (name .eq. ' ') then !{
2606  call mpp_error(fatal, trim(error_header) // ' Empty name given')
2607 endif !}
2608 
2609 !
2610 ! check that append and index are not both given
2611 !
2612 
2613 if (present(index) .and. present(append)) then !{
2614  call mpp_error(fatal, trim(error_header) // ' Append and index both given as arguments')
2615 endif !}
2616 
2617 !
2618 ! check for whether to overwrite existing values
2619 !
2620 
2621 if (present(no_overwrite)) then !{
2622  no_overwrite_use = no_overwrite
2623 else !}{
2624  no_overwrite_use = default_no_overwrite
2625 endif !}
2626 
2627 !
2628 ! check for whether to save the name in a list
2629 !
2630 
2631 if (present(good_name_list)) then !{
2632  good_name_list_use = good_name_list
2633 else !}{
2634  good_name_list_use = default_good_name_list
2635 endif !}
2636 
2637 if (present(no_create)) then !{
2638  create = .not. no_create
2639  if (no_create .and. (present(append) .or. present(index))) then !{
2640  call mpp_error(fatal, trim(error_header) // ' append or index are present when no_create is true for ' // trim(name))
2641  endif !}
2642 else !}{
2643  create = .true.
2644 endif !}
2645 
2646 if (present(index)) then !{
2647  if (fm_exists(name)) then !{
2648  field_length = fm_get_length(name)
2649  if (field_length .lt. 0) then !{
2650  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
2651  endif !}
2652  if (.not. (no_overwrite_use .and. field_length .ge. index)) then !{
2653  field_index = fm_new_value(name, value, index = index)
2654  if (field_index .le. 0) then !{
2655  write (str_error,*) ' with index = ', index
2656  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name) // trim(str_error))
2657  endif !}
2658  endif !}
2659  else !}{
2660  field_index = fm_new_value(name, value, index = index)
2661  if (field_index .le. 0) then !{
2662  write (str_error,*) ' with index = ', index
2663  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2664  endif !}
2665  endif !}
2666 elseif (present(append)) then !}{
2667  field_index = fm_new_value(name, value, append = append)
2668  if (field_index .le. 0) then !{
2669  write (str_error,*) ' with append = ', append
2670  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2671  endif !}
2672 else !}{
2673  if (fm_exists(name)) then !{
2674  if (.not. no_overwrite_use) then !{
2675  field_index = fm_new_value(name, value)
2676  if (field_index .le. 0) then !{
2677  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name))
2678  endif !}
2679  endif !}
2680  elseif (create) then !}{
2681  field_index = fm_new_value(name, value)
2682  if (field_index .le. 0) then !{
2683  call mpp_error(fatal, trim(error_header) // ' Problem creating ' // trim(name))
2684  endif !}
2685  endif !}
2686 endif !}
2687 
2688 !
2689 ! Add the variable name to the list of good names, to be used
2690 ! later for a consistency check, unless the field did not exist and we did not create it
2691 !
2692 
2693 if (good_name_list_use .ne. ' ') then !{
2694  if (fm_exists(good_name_list_use)) then !{
2695  add_name = fm_util_get_index_string(good_name_list_use, name, &
2696  caller = caller_str) .le. 0 ! true if name does not exist in string array
2697  else !}{
2698  add_name = .true. ! always add to new list
2699  endif !}
2700  if (add_name .and. fm_exists(name)) then !{
2701  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
2702  call mpp_error(fatal, trim(error_header) // &
2703  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
2704  endif !}
2705  endif !}
2706 endif !}
2707 
2708 return
2709 
2710 end subroutine fm_util_set_value_logical !}
2711 ! </SUBROUTINE> NAME="fm_util_set_value_logical"
2712 
2713 
2714 !#######################################################################
2715 ! <SUBROUTINE NAME="fm_util_set_value_real">
2716 !
2717 ! <DESCRIPTION>
2718 ! Set a real value in the Field Manager tree.
2719 ! </DESCRIPTION>
2720 !
2721 
2722 subroutine fm_util_set_value_real(name, value, caller, index, append, no_create, &
2723  no_overwrite, good_name_list) !{
2725 implicit none
2726 
2727 !
2728 ! arguments
2729 !
2730 
2731 character(len=*), intent(in) :: name
2732 real, intent(in) :: value
2733 character(len=*), intent(in), optional :: caller
2734 integer, intent(in), optional :: index
2735 logical, intent(in), optional :: append
2736 logical, intent(in), optional :: no_create
2737 logical, intent(in), optional :: no_overwrite
2738 character(len=*), intent(in), optional :: good_name_list
2739 
2740 !
2741 ! Local parameters
2742 !
2743 
2744 character(len=48), parameter :: sub_name = 'fm_util_set_value_real'
2745 
2746 !
2747 ! Local variables
2748 !
2749 
2750 character(len=256) :: error_header
2751 character(len=256) :: warn_header
2752 character(len=256) :: note_header
2753 character(len=128) :: caller_str
2754 character(len=32) :: str_error
2755 integer :: field_index
2756 logical :: no_overwrite_use
2757 integer :: field_length
2758 character(len=fm_path_name_len) :: good_name_list_use
2759 logical :: create
2760 logical :: add_name
2761 
2762 !
2763 ! set the caller string and headers
2764 !
2765 
2766 if (present(caller)) then !{
2767  caller_str = '[' // trim(caller) // ']'
2768 else !}{
2769  caller_str = fm_util_default_caller
2770 endif !}
2771 
2772 error_header = '==>Error from ' // trim(mod_name) // &
2773  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2774 warn_header = '==>Warning from ' // trim(mod_name) // &
2775  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2776 note_header = '==>Note from ' // trim(mod_name) // &
2777  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2778 
2779 !
2780 ! check that a name is given (fatal if not)
2781 !
2782 
2783 if (name .eq. ' ') then !{
2784  call mpp_error(fatal, trim(error_header) // ' Empty name given')
2785 endif !}
2786 
2787 !
2788 ! check that append and index are not both given
2789 !
2790 
2791 if (present(index) .and. present(append)) then !{
2792  call mpp_error(fatal, trim(error_header) // ' Append and index both given as arguments')
2793 endif !}
2794 
2795 !
2796 ! check for whether to overwrite existing values
2797 !
2798 
2799 if (present(no_overwrite)) then !{
2800  no_overwrite_use = no_overwrite
2801 else !}{
2802  no_overwrite_use = default_no_overwrite
2803 endif !}
2804 
2805 !
2806 ! check for whether to save the name in a list
2807 !
2808 
2809 if (present(good_name_list)) then !{
2810  good_name_list_use = good_name_list
2811 else !}{
2812  good_name_list_use = default_good_name_list
2813 endif !}
2814 
2815 if (present(no_create)) then !{
2816  create = .not. no_create
2817  if (no_create .and. (present(append) .or. present(index))) then !{
2818  call mpp_error(fatal, trim(error_header) // ' append or index are present when no_create is true for ' // trim(name))
2819  endif !}
2820 else !}{
2821  create = .true.
2822 endif !}
2823 
2824 if (present(index)) then !{
2825  if (fm_exists(name)) then !{
2826  field_length = fm_get_length(name)
2827  if (field_length .lt. 0) then !{
2828  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
2829  endif !}
2830  if (.not. (no_overwrite_use .and. field_length .ge. index)) then !{
2831  field_index = fm_new_value(name, value, index = index)
2832  if (field_index .le. 0) then !{
2833  write (str_error,*) ' with index = ', index
2834  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name) // trim(str_error))
2835  endif !}
2836  endif !}
2837  else !}{
2838  field_index = fm_new_value(name, value, index = index)
2839  if (field_index .le. 0) then !{
2840  write (str_error,*) ' with index = ', index
2841  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2842  endif !}
2843  endif !}
2844 elseif (present(append)) then !}{
2845  field_index = fm_new_value(name, value, append = append)
2846  if (field_index .le. 0) then !{
2847  write (str_error,*) ' with append = ', append
2848  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
2849  endif !}
2850 else !}{
2851  if (fm_exists(name)) then !{
2852  if (.not. no_overwrite_use) then !{
2853  field_index = fm_new_value(name, value)
2854  if (field_index .le. 0) then !{
2855  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name))
2856  endif !}
2857  endif !}
2858  elseif (create) then !}{
2859  field_index = fm_new_value(name, value)
2860  if (field_index .le. 0) then !{
2861  call mpp_error(fatal, trim(error_header) // ' Problem creating ' // trim(name))
2862  endif !}
2863  endif !}
2864 endif !}
2865 
2866 !
2867 ! Add the variable name to the list of good names, to be used
2868 ! later for a consistency check, unless the field did not exist and we did not create it
2869 !
2870 
2871 if (good_name_list_use .ne. ' ') then !{
2872  if (fm_exists(good_name_list_use)) then !{
2873  add_name = fm_util_get_index_string(good_name_list_use, name, &
2874  caller = caller_str) .le. 0 ! true if name does not exist in string array
2875  else !}{
2876  add_name = .true. ! always add to new list
2877  endif !}
2878  if (add_name .and. fm_exists(name)) then !{
2879  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
2880  call mpp_error(fatal, trim(error_header) // &
2881  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
2882  endif !}
2883  endif !}
2884 endif !}
2885 
2886 return
2887 
2888 end subroutine fm_util_set_value_real !}
2889 ! </SUBROUTINE> NAME="fm_util_set_value_real"
2890 
2891 
2892 !#######################################################################
2893 ! <SUBROUTINE NAME="fm_util_set_value_string">
2894 !
2895 ! <DESCRIPTION>
2896 ! Set a string value in the Field Manager tree.
2897 ! </DESCRIPTION>
2898 !
2899 
2900 subroutine fm_util_set_value_string(name, value, caller, index, append, no_create, &
2901  no_overwrite, good_name_list) !{
2903 implicit none
2904 
2905 !
2906 ! arguments
2907 !
2908 
2909 character(len=*), intent(in) :: name
2910 character(len=*), intent(in) :: value
2911 character(len=*), intent(in), optional :: caller
2912 integer, intent(in), optional :: index
2913 logical, intent(in), optional :: append
2914 logical, intent(in), optional :: no_create
2915 logical, intent(in), optional :: no_overwrite
2916 character(len=*), intent(in), optional :: good_name_list
2917 
2918 !
2919 ! Local parameters
2920 !
2921 
2922 character(len=48), parameter :: sub_name = 'fm_util_set_value_string'
2923 
2924 !
2925 ! Local variables
2926 !
2927 
2928 character(len=256) :: error_header
2929 character(len=256) :: warn_header
2930 character(len=256) :: note_header
2931 character(len=128) :: caller_str
2932 character(len=32) :: str_error
2933 integer :: field_index
2934 logical :: no_overwrite_use
2935 integer :: field_length
2936 character(len=fm_path_name_len) :: good_name_list_use
2937 logical :: create
2938 logical :: add_name
2939 
2940 !
2941 ! set the caller string and headers
2942 !
2943 
2944 if (present(caller)) then !{
2945  caller_str = '[' // trim(caller) // ']'
2946 else !}{
2947  caller_str = fm_util_default_caller
2948 endif !}
2949 
2950 error_header = '==>Error from ' // trim(mod_name) // &
2951  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2952 warn_header = '==>Warning from ' // trim(mod_name) // &
2953  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2954 note_header = '==>Note from ' // trim(mod_name) // &
2955  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
2956 
2957 !
2958 ! check that a name is given (fatal if not)
2959 !
2960 
2961 if (name .eq. ' ') then !{
2962  call mpp_error(fatal, trim(error_header) // ' Empty name given')
2963 endif !}
2964 
2965 !
2966 ! check that append and index are not both given
2967 !
2968 
2969 if (present(index) .and. present(append)) then !{
2970  call mpp_error(fatal, trim(error_header) // ' Append and index both given as arguments')
2971 endif !}
2972 
2973 !
2974 ! check for whether to overwrite existing values
2975 !
2976 
2977 if (present(no_overwrite)) then !{
2978  no_overwrite_use = no_overwrite
2979 else !}{
2980  no_overwrite_use = default_no_overwrite
2981 endif !}
2982 
2983 !
2984 ! check for whether to save the name in a list
2985 !
2986 
2987 if (present(good_name_list)) then !{
2988  good_name_list_use = good_name_list
2989 else !}{
2990  good_name_list_use = default_good_name_list
2991 endif !}
2992 
2993 if (present(no_create)) then !{
2994  create = .not. no_create
2995  if (no_create .and. (present(append) .or. present(index))) then !{
2996  call mpp_error(fatal, trim(error_header) // ' append or index are present when no_create is true for ' // trim(name))
2997  endif !}
2998 else !}{
2999  create = .true.
3000 endif !}
3001 
3002 if (present(index)) then !{
3003  if (fm_exists(name)) then !{
3004  field_length = fm_get_length(name)
3005  if (field_length .lt. 0) then !{
3006  call mpp_error(fatal, trim(error_header) // ' Problem getting length of ' // trim(name))
3007  endif !}
3008  if (.not. (no_overwrite_use .and. field_length .ge. index)) then !{
3009  field_index = fm_new_value(name, value, index = index)
3010  if (field_index .le. 0) then !{
3011  write (str_error,*) ' with index = ', index
3012  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name) // trim(str_error))
3013  endif !}
3014  endif !}
3015  else !}{
3016  field_index = fm_new_value(name, value, index = index)
3017  if (field_index .le. 0) then !{
3018  write (str_error,*) ' with index = ', index
3019  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
3020  endif !}
3021  endif !}
3022 elseif (present(append)) then !}{
3023  field_index = fm_new_value(name, value, append = append)
3024  if (field_index .le. 0) then !{
3025  write (str_error,*) ' with append = ', append
3026  call mpp_error(fatal, trim(error_header) // ' Problem setting ' // trim(name) // trim(str_error))
3027  endif !}
3028 else !}{
3029  if (fm_exists(name)) then !{
3030  if (.not. no_overwrite_use) then !{
3031  field_index = fm_new_value(name, value)
3032  if (field_index .le. 0) then !{
3033  call mpp_error(fatal, trim(error_header) // ' Problem overwriting ' // trim(name))
3034  endif !}
3035  endif !}
3036  elseif (create) then !}{
3037  field_index = fm_new_value(name, value)
3038  if (field_index .le. 0) then !{
3039  call mpp_error(fatal, trim(error_header) // ' Problem creating ' // trim(name))
3040  endif !}
3041  endif !}
3042 endif !}
3043 
3044 !
3045 ! Add the variable name to the list of good names, to be used
3046 ! later for a consistency check, unless the field did not exist and we did not create it
3047 !
3048 
3049 if (good_name_list_use .ne. ' ') then !{
3050  if (fm_exists(good_name_list_use)) then !{
3051  add_name = fm_util_get_index_string(good_name_list_use, name, &
3052  caller = caller_str) .le. 0 ! true if name does not exist in string array
3053  else !}{
3054  add_name = .true. ! always add to new list
3055  endif !}
3056  if (add_name .and. fm_exists(name)) then !{
3057  if (fm_new_value(good_name_list_use, name, append = .true., create = .true.) .le. 0) then !{
3058  call mpp_error(fatal, trim(error_header) // &
3059  ' Could not add ' // trim(name) // ' to "' // trim(good_name_list_use) // '" list')
3060  endif !}
3061  endif !}
3062 endif !}
3063 
3064 return
3065 
3066 end subroutine fm_util_set_value_string !}
3067 ! </SUBROUTINE> NAME="fm_util_set_value_string"
3068 
3069 
3070 !#######################################################################
3071 ! <SUBROUTINE NAME="fm_util_start_namelist">
3072 !
3073 ! <DESCRIPTION>
3074 ! Start processing a namelist
3075 ! </DESCRIPTION>
3076 !
3077 subroutine fm_util_start_namelist(path, name, caller, no_overwrite, check) !{
3079 implicit none
3080 
3081 !
3082 ! arguments
3083 !
3084 
3085 character(len=*), intent(in) :: path
3086 character(len=*), intent(in) :: name
3087 character(len=*), intent(in), optional :: caller
3088 logical, intent(in), optional :: no_overwrite
3089 logical, intent(in), optional :: check
3090 
3091 !
3092 ! Local parameters
3093 !
3094 
3095 character(len=48), parameter :: sub_name = 'fm_util_start_namelist'
3096 
3097 !
3098 ! Local variables
3099 !
3100 
3101 integer :: namelist_index
3102 character(len=fm_path_name_len) :: path_name
3103 character(len=256) :: error_header
3104 character(len=256) :: warn_header
3105 character(len=256) :: note_header
3106 character(len=128) :: caller_str
3107 integer :: out_unit
3108 
3109 out_unit = stdout()
3110 
3111 !
3112 ! set the caller string and headers
3113 !
3114 
3115 if (present(caller)) then !{
3116  caller_str = '[' // trim(caller) // ']'
3117 else !}{
3118  caller_str = fm_util_default_caller
3119 endif !}
3120 
3121 error_header = '==>Error from ' // trim(mod_name) // &
3122  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
3123 warn_header = '==>Warning from ' // trim(mod_name) // &
3124  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
3125 note_header = '==>Note from ' // trim(mod_name) // &
3126  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
3127 
3128 !
3129 ! check that a name is given (fatal if not)
3130 !
3131 
3132 if (name .eq. ' ') then !{
3133  call mpp_error(fatal, trim(error_header) // ' Empty name given')
3134 endif !}
3135 
3136 !
3137 ! Concatenate the path and name
3138 !
3139 
3140 if (path .eq. ' ') then !{
3141  path_name = name
3142 else !}{
3143  path_name = trim(path) // '/' // name
3144 endif !}
3145 save_path = path
3146 save_name = name
3147 
3148 !
3149 ! set the default caller string, if desired
3150 !
3151 
3152 if (present(caller)) then !{
3153  call fm_util_set_caller(caller)
3154 else !}{
3156 endif !}
3157 
3158 !
3159 ! set the default no_overwrite flag, if desired
3160 !
3161 
3162 if (present(no_overwrite)) then !{
3163  call fm_util_set_no_overwrite(no_overwrite)
3164 else !}{
3166 endif !}
3167 
3168 !
3169 ! set the default good_name_list string, if desired
3170 !
3171 
3172 if (present(check)) then !{
3173  if (check) then !{
3174  call fm_util_set_good_name_list('/ocean_mod/GOOD/namelists/' // trim(path_name) // '/good_list')
3175  else !}{
3177  endif !}
3178 else !}{
3180 endif !}
3181 
3182 !
3183 ! Process the namelist
3184 !
3185 
3186 write (out_unit,*)
3187 write (out_unit,*) trim(note_header), ' Processing namelist ', trim(path_name)
3188 
3189 !
3190 ! Check whether the namelist already exists. If so, then use that one
3191 !
3192 
3193 namelist_index = fm_get_index('/ocean_mod/namelists/' // trim(path_name))
3194 if (namelist_index .gt. 0) then !{
3195 
3196  !write (out_unit,*) trim(note_header), ' Namelist already set with index ', namelist_index
3197 
3198 else !}{
3199 
3200 !
3201 ! Set a new namelist and get its index
3202 !
3203 
3204  namelist_index = fm_new_list('/ocean_mod/namelists/' // trim(path_name), create = .true.)
3205  if (namelist_index .le. 0) then !{
3206  call mpp_error(fatal, trim(error_header) // ' Could not set namelist ' // trim(path_name))
3207  endif !}
3208 
3209 endif !}
3210 
3211 !
3212 ! Add the namelist name to the list of good namelists, to be used
3213 ! later for a consistency check
3214 !
3215 
3216 if (fm_new_value('/ocean_mod/GOOD/namelists/' // trim(path) // '/good_values', &
3217  name, append = .true., create = .true.) .le. 0) then !{
3218  call mpp_error(fatal, trim(error_header) // &
3219  ' Could not add ' // trim(name) // ' to "' // trim(path) // '/good_values" list')
3220 endif !}
3221 
3222 !
3223 ! Change to the new namelist, first saving the current list
3224 !
3225 
3226 save_current_list = fm_get_current_list()
3227 if (save_current_list .eq. ' ') then !{
3228  call mpp_error(fatal, trim(error_header) // ' Could not get the current list')
3229 endif !}
3230 
3231 if (.not. fm_change_list('/ocean_mod/namelists/' // trim(path_name))) then !{
3232  call mpp_error(fatal, trim(error_header) // ' Could not change to the namelist ' // trim(path_name))
3233 endif !}
3234 
3235 return
3236 
3237 end subroutine fm_util_start_namelist !}
3238 ! </SUBROUTINE> NAME="fm_util_start_namelist"
3239 
3240 
3241 !#######################################################################
3242 ! <SUBROUTINE NAME="fm_util_end_namelist">
3243 !
3244 ! <DESCRIPTION>
3245 ! Finish up processing a namelist
3246 ! </DESCRIPTION>
3247 !
3248 subroutine fm_util_end_namelist(path, name, caller, check) !{
3250 implicit none
3251 
3252 !
3253 ! arguments
3254 !
3255 
3256 character(len=*), intent(in) :: path
3257 character(len=*), intent(in) :: name
3258 character(len=*), intent(in), optional :: caller
3259 logical, intent(in), optional :: check
3260 
3261 !
3262 ! Local parameters
3263 !
3264 
3265 character(len=48), parameter :: sub_name = 'fm_util_end_namelist'
3266 
3267 !
3268 ! Local variables
3269 !
3270 
3271 character(len=fm_string_len), pointer, dimension(:) :: good_list => null()
3272 character(len=fm_path_name_len) :: path_name
3273 character(len=256) :: error_header
3274 character(len=256) :: warn_header
3275 character(len=256) :: note_header
3276 character(len=128) :: caller_str
3277 
3278 !
3279 ! set the caller string and headers
3280 !
3281 
3282 if (present(caller)) then !{
3283  caller_str = '[' // trim(caller) // ']'
3284 else !}{
3285  caller_str = fm_util_default_caller
3286 endif !}
3287 
3288 error_header = '==>Error from ' // trim(mod_name) // &
3289  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
3290 warn_header = '==>Warning from ' // trim(mod_name) // &
3291  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
3292 note_header = '==>Note from ' // trim(mod_name) // &
3293  '(' // trim(sub_name) // ')' // trim(caller_str) // ':'
3294 
3295 !
3296 ! check that a path is given (fatal if not)
3297 !
3298 
3299 if (name .eq. ' ') then !{
3300  call mpp_error(fatal, trim(error_header) // ' Empty name given')
3301 endif !}
3302 
3303 !
3304 ! Check that the path ane name match the preceding call to
3305 ! fm_util_start_namelist
3306 !
3307 
3308 if (path .ne. save_path) then !{
3309  call mpp_error(fatal, trim(error_header) // ' Path "' // trim(path) // '" does not match saved path "' // trim(save_path) // '"')
3310 elseif (name .ne. save_name) then !}{
3311  call mpp_error(fatal, trim(error_header) // ' Name "' // trim(name) // '" does not match saved name "' // trim(save_name) // '"')
3312 endif !}
3313 
3314 !
3315 ! Concatenate the path and name
3316 !
3317 
3318 if (path .eq. ' ') then !{
3319  path_name = name
3320 else !}{
3321  path_name = trim(path) // '/' // name
3322 endif !}
3323 save_path = ' '
3324 save_name = ' '
3325 
3326 !
3327 ! Check for any errors in the number of fields in this list
3328 !
3329 
3330 if (present(check)) then !{
3331  if (check) then !{
3332  if (caller_str .eq. ' ') then !{
3333  caller_str = trim(mod_name) // '(' // trim(sub_name) // ')'
3334  endif !}
3335  good_list => fm_util_get_string_array('/ocean_mod/GOOD/namelists/' // trim(path_name) // '/good_list', &
3336  caller = trim(mod_name) // '(' // trim(sub_name) // ')')
3337  if (associated(good_list)) then !{
3338  call fm_util_check_for_bad_fields('/ocean_mod/namelists/' // trim(path_name), good_list, caller = caller_str)
3339  deallocate(good_list)
3340  else !}{
3341  call mpp_error(fatal, trim(error_header) // ' Empty "' // trim(path_name) // '" list')
3342  endif !}
3343  endif !}
3344 endif !}
3345 
3346 !
3347 ! Change back to the saved list
3348 !
3349 
3350 if (save_current_list .ne. ' ') then !{
3351  if (.not. fm_change_list(save_current_list)) then !{
3352  call mpp_error(fatal, trim(error_header) // ' Could not change to the saved list: ' // trim(save_current_list))
3353  endif !}
3354 endif !}
3355 save_current_list = ' '
3356 
3357 !
3358 ! reset the default caller string
3359 !
3360 
3362 
3363 !
3364 ! reset the default no_overwrite string
3365 !
3366 
3368 
3369 !
3370 ! reset the default good_name_list string
3371 !
3372 
3374 
3375 return
3376 
3377 end subroutine fm_util_end_namelist !}
3378 ! </SUBROUTINE> NAME="fm_util_end_namelist"
3379 
3380 
3381 end module fm_util_mod !}
Definition: fms.F90:20
character(len=8) function, public fm_get_type(name)
integer function, public fm_get_index(name)
integer, parameter, public fm_path_name_len
integer function, public fm_util_get_index_string(name, string, caller)
Definition: fm_util.F90:637
logical function, public fm_util_get_logical(name, caller, index, default_value, scalar)
Definition: fm_util.F90:1348
subroutine, public fm_util_reset_no_overwrite
Definition: fm_util.F90:359
subroutine, public fm_util_set_value_string(name, value, caller, index, append, no_create, no_overwrite, good_name_list)
Definition: fm_util.F90:2902
subroutine, public fm_util_set_value_integer_array(name, value, length, caller, no_overwrite, good_name_list)
Definition: fm_util.F90:1711
character(len=128) save_default_caller
Definition: fm_util.F90:101
integer function, public fm_util_get_index_list(name, caller)
Definition: fm_util.F90:747
logical save_default_no_overwrite
Definition: fm_util.F90:105
subroutine, public fm_util_set_value_logical(name, value, caller, index, append, no_create, no_overwrite, good_name_list)
Definition: fm_util.F90:2546
integer function, public fm_new_list(name, create, keep)
character(len=fm_string_len) function, public fm_util_get_string(name, caller, index, default_value, scalar)
Definition: fm_util.F90:1592
subroutine, public fm_util_set_no_overwrite(no_overwrite)
Definition: fm_util.F90:317
integer, parameter, public fm_string_len
subroutine, public fm_util_set_value_real(name, value, caller, index, append, no_create, no_overwrite, good_name_list)
Definition: fm_util.F90:2724
logical function, public fm_change_list(name)
integer function, public fm_util_get_length(name, caller)
Definition: fm_util.F90:557
Definition: mpp.F90:39
subroutine, public fm_util_set_value_real_array(name, value, length, caller, no_overwrite, good_name_list)
Definition: fm_util.F90:2039
integer, parameter, public fm_type_name_len
subroutine, public fm_util_check_for_bad_fields(list, good_fields, caller)
Definition: fm_util.F90:393
subroutine, public fm_util_set_good_name_list(good_name_list)
Definition: fm_util.F90:240
real function, public fm_util_get_real(name, caller, index, default_value, scalar)
Definition: fm_util.F90:1467
subroutine, public fm_util_set_caller(caller)
Definition: fm_util.F90:159
subroutine, public fm_util_reset_caller
Definition: fm_util.F90:205
logical function, public fm_dump_list(name, recursive, unit)
logical function, public fm_exists(name)
subroutine, public fm_util_start_namelist(path, name, caller, no_overwrite, check)
Definition: fm_util.F90:3078
character(len=128), public fm_util_default_caller
Definition: fm_util.F90:89
logical default_no_overwrite
Definition: fm_util.F90:104
subroutine, public fm_util_reset_good_name_list
Definition: fm_util.F90:282
integer function, public fm_util_get_integer(name, caller, index, default_value, scalar)
Definition: fm_util.F90:1229
subroutine, public fm_util_set_value_string_array(name, value, length, caller, no_overwrite, good_name_list)
Definition: fm_util.F90:2203
subroutine, public fm_util_end_namelist(path, name, caller, check)
Definition: fm_util.F90:3249
character(len=fm_string_len) function, dimension(:), pointer, public fm_util_get_string_array(name, caller)
Definition: fm_util.F90:1131
character(len=fm_path_name_len) save_path
Definition: fm_util.F90:107
subroutine, public fm_util_set_value_logical_array(name, value, length, caller, no_overwrite, good_name_list)
Definition: fm_util.F90:1875
real function, dimension(:), pointer, public fm_util_get_real_array(name, caller)
Definition: fm_util.F90:1033
character(len=fm_path_name_len) save_current_list
Definition: fm_util.F90:106
character(len=48), parameter mod_name
Definition: fm_util.F90:95
character(len=128) default_good_name_list
Definition: fm_util.F90:102
character(len=fm_path_name_len) save_name
Definition: fm_util.F90:108
integer function, public fm_get_length(name)
subroutine, public fm_util_set_value_integer(name, value, caller, index, append, no_create, no_overwrite, good_name_list)
Definition: fm_util.F90:2368
integer function, dimension(:), pointer, public fm_util_get_integer_array(name, caller)
Definition: fm_util.F90:837
character(len=128) save_default_good_name_list
Definition: fm_util.F90:103
character(len=fm_path_name_len) function, public fm_get_current_list()
integer, parameter, public fm_field_name_len
logical function, dimension(:), pointer, public fm_util_get_logical_array(name, caller)
Definition: fm_util.F90:935