Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
dictionary.F90
Go to the documentation of this file.
2
3 use iso_c_binding, only : c_int, c_float, c_double, c_bool
5 use exceptions
6 use logfiles
7 use fstring
9 use exceptions
10 implicit none
11
12 private
13
14 type, public :: dict_entry_t
15
16 character (len=:), allocatable :: key
17 character (len=:), allocatable :: secondary_key
18
19 type (fstring_list_t) :: sl
20 type (dict_entry_t), pointer :: previous => null()
21 type (dict_entry_t), pointer :: next => null()
22
23 contains
24
25 procedure :: add_key_sub
26 generic :: add_key => add_key_sub
27
28 procedure :: add_string_sub
29 procedure :: add_integer_sub
30 procedure :: add_float_sub
31 procedure :: add_double_sub
32 procedure :: add_logical_sub
33 generic :: add_value => add_string_sub, &
38 end type dict_entry_t
39
40
41 type, public :: dict_t
42
43 type (dict_entry_t), pointer :: first => null()
44 type (dict_entry_t), pointer :: last => null()
45 type (dict_entry_t), pointer :: current => null()
46 integer (c_int) :: count = 0
47
48 contains
49
50 procedure, private :: get_entry_by_key_fn
51 procedure, private :: get_entry_by_index_fn
52 generic :: get_entry => get_entry_by_key_fn, &
54
55 procedure, private :: get_next_entry_by_key_fn
56 procedure, private :: get_next_entry_fn
57 generic :: get_next_entry => get_next_entry_by_key_fn, &
59
60 procedure, private :: add_entry_to_dict_sub
61 generic :: add_entry => add_entry_to_dict_sub
62
63 procedure, private :: key_name_already_in_use_fn
64 generic :: key_already_in_use => key_name_already_in_use_fn
65
66 procedure :: find_dict_entry_fn
67 generic :: find_dict_entry => find_dict_entry_fn
68
69 procedure, private :: delete_entry_by_key_sub
70 generic :: delete_entry => delete_entry_by_key_sub
71
73 generic :: print_all => print_all_dictionary_entries_sub
74
75 procedure, private :: grep_dictionary_key_names_fn
76 generic :: grep_keys => grep_dictionary_key_names_fn
77
78 procedure, public :: get_value => get_value_as_string_sub
79
80 procedure, private :: get_values_as_int_sub
81 procedure, private :: get_values_as_float_sub
82 procedure, private :: get_values_as_logical_sub
83 procedure, private :: get_values_as_string_list_sub
88 generic :: get_values => get_values_as_int_sub, &
96
97 end type dict_t
98
99 ! CF = "Control File"; this dictionary will be populated elsewhere with all of the
100 ! directives found in the SWB control file.
101 type (dict_t), public :: cf_dict
102 type (dict_entry_t), public, pointer :: cf_entry
103
104contains
105
106 !!
107 !! This section contains methods bound to the DICT_ENTRY_T class
108 !!
109
110 subroutine add_key_sub(this, sKey)
111
112 class(dict_entry_t) :: this
113 character (len=*), intent(in) :: sKey
114
115 this%key = trim(skey)
116
117 end subroutine add_key_sub
118
119!--------------------------------------------------------------------------------------------------
120
121 subroutine add_string_sub(this, sValue)
122
123 class(dict_entry_t) :: this
124 character (len=*), intent(in) :: sValue
125
126 call this%sl%append(svalue)
127
128 end subroutine add_string_sub
129
130!--------------------------------------------------------------------------------------------------
131
132 subroutine add_float_sub(this, fValue)
133
134 class(dict_entry_t) :: this
135 real (c_float), intent(in) :: fValue
136
137 call this%sl%append( ascharacter( fvalue ) )
138
139 end subroutine add_float_sub
140
141!--------------------------------------------------------------------------------------------------
142
143 subroutine add_integer_sub(this, iValue)
144
145 class(dict_entry_t) :: this
146 integer (c_int), intent(in) :: iValue
147
148 call this%sl%append( ascharacter( ivalue ) )
149
150 end subroutine add_integer_sub
151
152!--------------------------------------------------------------------------------------------------
153
154 subroutine add_double_sub(this, dValue)
155
156 class(dict_entry_t) :: this
157 real (c_double), intent(in) :: dValue
158
159 call this%sl%append( ascharacter( dvalue ) )
160
161 end subroutine add_double_sub
162
163!--------------------------------------------------------------------------------------------------
164
165 subroutine add_logical_sub(this, lValue)
166
167 class(dict_entry_t) :: this
168 logical (c_bool), intent(in) :: lValue
169
170 call this%sl%append( ascharacter( lvalue ) )
171
172 end subroutine add_logical_sub
173
174!--------------------------------------------------------------------------------------------------
175
176 !!
177 !! The methods in the section below are bound to the DICT_T type, which is essentially a
178 !! linked list of DICT_ENTRY_T objects
179 !!
180
181 function get_entry_by_key_fn(this, sKey) result( pDict )
182
183 class(dict_t) :: this
184 character (len=*), intent(in) :: skey
185 type (dict_entry_t), pointer :: pdict
186
187 this%current => this%first
188
189 do while ( associated( this%current ) )
190
191 if ( this%current%key .strapprox. skey ) exit
192
193 this%current => this%current%next
194
195 enddo
196
197 if (.not. associated( this%current ) ) &
198 call warn( smessage="Failed to find a dictionary entry with a key value of "//dquote(skey), &
199 iloglevel=log_debug, &
200 lecho=false )
201
202 pdict => this%current
203
204 end function get_entry_by_key_fn
205
206!--------------------------------------------------------------------------------------------------
207
208 function get_entry_by_index_fn(this, iIndex) result( pDict )
209
210 class(dict_t) :: this
211 integer (c_int), intent(in) :: iindex
212 type (dict_entry_t), pointer :: pdict
213
214 ! [ LOCALS ]
215 integer (c_int) :: icurrentindex
216
217 icurrentindex = 1
218
219 this%current => this%first
220
221 do
222
223 if ( .not. associated( this%current ) ) exit
224
225 if ( icurrentindex == iindex ) exit
226
227 this%current => this%current%next
228
229 icurrentindex = icurrentindex + 1
230
231 enddo
232
233 if (.not. associated( this%current ) ) &
234 call warn( smessage="Failed to find a dictionary entry with a index value of " &
235 //ascharacter(iindex), &
236 iloglevel=log_debug, &
237 lecho=false )
238
239 pdict => this%current
240
241 end function get_entry_by_index_fn
242
243!--------------------------------------------------------------------------------------------------
244
245 function get_next_entry_by_key_fn(this, sKey) result( pDict )
246
247 class(dict_t) :: this
248 character (len=*), intent(in) :: skey
249 type (dict_entry_t), pointer :: pdict
250
251 ! if "current" location is not null, it will point to the location of the
252 ! last key value found. move forward by one before examining the next key...
253 if (associated( this%current) ) pdict => this%current%next
254
255 do while ( associated( pdict ) )
256
257 if ( pdict%key .strequal. skey ) exit
258
259 pdict => pdict%next
260 this%current => pdict
261
262 enddo
263
264 if (.not. associated( pdict ) ) &
265 call warn( smessage="Failed to find another dictionary entry with a key value of "//dquote(skey), &
266 iloglevel=log_debug, &
267 lecho=false )
268
269 end function get_next_entry_by_key_fn
270
271!--------------------------------------------------------------------------------------------------
272
273 function get_next_entry_fn(this) result( pDict )
274
275 class(dict_t) :: this
276 type (dict_entry_t), pointer :: pdict
277
278 pdict => null()
279
280 ! if "current" location is not null, it will point to the location of the
281 ! last key value found. move forward by one before examining the next key...
282 if (associated( this%current) ) then
283 pdict => this%current%next
284 this%current => this%current%next
285 endif
286
287 if (.not. associated( pdict ) ) &
288 call warn( smessage="Reached end of dictionary.", &
289 iloglevel=log_debug, &
290 lecho=false )
291
292 end function get_next_entry_fn
293
294!--------------------------------------------------------------------------------------------------
295
296 function grep_dictionary_key_names_fn(this, sKey) result( slString )
297
298 class(dict_t) :: this
299 character (len=*), intent(in) :: skey
300 type (fstring_list_t) :: slstring
301
302 ! [ LOCALS ]
303 integer (c_int) :: iindex
304
305 this%current => this%first
306
307 do while ( associated(this%current ) )
308
309 if ( this%current%key .containssimilar. skey ) &
310 call slstring%append(this%current%key)
311
312 this%current => this%current%next
313
314 enddo
315
316 if ( slstring%get(1) == '<NA>' ) &
317 call warn(smessage="Failed to find a dictionary entry associated with a key value of " &
318 //dquote(skey)//".", smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false )
319
321
322!--------------------------------------------------------------------------------------------------
323
324function key_name_already_in_use_fn(this, sKey) result( in_use )
325
326 class(dict_t) :: this
327 character (len=*), intent(in) :: skey
328 logical (c_bool) :: in_use
329
330 ! [ LOCALS ]
331 integer (c_int) :: iindex
332 integer (c_int) :: icount
333
334 this%current => this%first
335 icount = 0
336
337 do while ( associated(this%current ) )
338
339 !iIndex = index(string=this%current%key, substring=asUppercase(sKey) )
340
341 !if ( iIndex > 0 ) iCount = iCount + 1
342
343 if ( this%current%key .strapprox. skey ) icount = icount + 1
344
345 this%current => this%current%next
346
347 enddo
348
349 if ( icount == 0 ) then
350
351 in_use = false
352
353 else
354
355 in_use = true
356
357 endif
358
360
361
362!--------------------------------------------------------------------------------------------------
363
364function find_dict_entry_fn(this, sSearchKey) result( pDict )
365
366 class(dict_t) :: this
367 character (len=*), intent(in) :: ssearchkey
368 type (dict_entry_t), pointer :: pdict
369
370 ! [ LOCALS ]
371 integer (c_int) :: iindex
372
373 pdict => null()
374 this%current => this%first
375
376 do while ( associated(this%current ) )
377
378 !iIndex = index(string=this%current%key, substring=asUppercase(sSearchKey) )
379
380 !if ( iIndex > 0 ) then
381
382 if (this%current%key .strapprox. ssearchkey ) then
383
384 pdict => this%current
385 exit
386
387 endif
388
389 this%current => this%current%next
390
391 enddo
392
393end function find_dict_entry_fn
394
395!--------------------------------------------------------------------------------------------------
396
397 subroutine add_entry_to_dict_sub(this, dict_entry)
398
399 class(dict_t) :: this
400 type (DICT_ENTRY_T), pointer :: dict_entry
401
402 type (DICT_ENTRY_T), pointer :: temp_dict_entry
403 integer (c_int) :: iIndex
404
405 temp_dict_entry => null()
406
407 if ( associated(dict_entry) ) then
408
409 temp_dict_entry => this%find_dict_entry( dict_entry%key )
410
411 ! the dictionary already contains a dictionary entry with this
412 ! key value; append the values to the existing entry
413 if ( associated( temp_dict_entry ) ) then
414
415 do iindex=1, dict_entry%sl%count
416
417 call temp_dict_entry%sl%append( dict_entry%sl%get( iindex ) )
418
419 enddo
420
421 temp_dict_entry => null()
422
423 else
424
425 this%count = this%count + 1
426
427 if ( associated( this%last) ) then
428
429 ! dictionary has at least one entry
430
431 dict_entry%previous => this%last
432 dict_entry%next => null()
433 this%last%next => dict_entry
434 this%last => dict_entry
435 this%current => dict_entry
436
437 else ! this is the first dictionary entry
438
439 this%first => dict_entry
440 this%last => dict_entry
441 this%current => dict_entry
442 dict_entry%previous => null()
443 dict_entry%next => null()
444
445 endif
446
447 endif
448
449 else
450
451 call warn( smessage="Internal programming error: dictionary entry is null", &
452 smodule=__file__, iline=__line__ )
453
454 endif
455
456 end subroutine add_entry_to_dict_sub
457
458!--------------------------------------------------------------------------------------------------
459
460 subroutine delete_entry_by_key_sub(this, sKey)
461
462 class(dict_t) :: this
463 character (len=*), intent(in) :: sKey
464
465 ! [ LOCALS ]
466 type (DICT_ENTRY_T), pointer :: pTarget
467 type (DICT_ENTRY_T), pointer :: pTemp
468
469 ptarget => this%get_entry(skey)
470
471 if ( associated( ptarget ) ) then
472
473 ! first remove object from linked list
474 ptemp => ptarget%previous
475
476 if (associated( ptemp) ) then
477 ! if pTemp is unassociated, it means we are at the head of the list
478 ptemp%next => ptarget%next
479 ptarget%next%previous => ptemp
480
481 else
482
483 ptemp => ptarget%next
484 this%first => ptemp
485
486 endif
487
488 ! set "current" pointer to entry that was just before the now-deleted entry
489 this%current => ptemp
490
491 call ptarget%sl%clear()
492
493 endif
494
495 end subroutine delete_entry_by_key_sub
496
497!--------------------------------------------------------------------------------------------------
498
499 subroutine get_values_as_int_sub(this, sKey, iValues, is_fatal )
500
501 class(dict_t) :: this
502 character (len=*) :: sKey
503 integer (c_int), allocatable, intent(out) :: iValues(:)
504 logical (c_bool), optional :: is_fatal
505
506 ! [ LOCALS ]
507 type (DICT_ENTRY_T), pointer :: pTarget
508 integer (c_int) :: iStat
509 logical (c_bool) :: is_fatal_l
510 logical (c_bool) :: empty_entries
511
512 if ( present( is_fatal ) ) then
513 is_fatal_l = is_fatal
514 else
515 is_fatal_l = false
516 endif
517
518 ptarget => this%get_entry(skey)
519
520 if ( associated( ptarget ) ) then
521
522 if ( is_fatal_l ) then
523 empty_entries = ptarget%sl%empty_entries_present()
524 if( empty_entries ) call warn(smessage="There are missing values associated" &
525 //" with the key value of "//dquote(skey), &
526 iloglevel=log_all, lecho=true, lfatal=true )
527 endif
528
529 ivalues = ptarget%sl%get_integer()
530
531 else
532
533 allocate(ivalues(1), stat=istat)
534 call assert(istat == 0, "Failed to allocate memory to iValues array", &
535 __file__, __line__)
536
537 call warn(smessage="Failed to find a dictionary entry associated with key value of "//dquote(skey), &
538 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
539
540 ivalues = itinyval
541
542 endif
543
544 end subroutine get_values_as_int_sub
545
546!--------------------------------------------------------------------------------------------------
547
548 subroutine get_values_as_logical_sub(this, sKey, lValues, is_fatal)
549
550 class(dict_t) :: this
551 character (len=*) :: sKey
552 logical (c_bool), allocatable, intent(out) :: lValues(:)
553 logical (c_bool), optional :: is_fatal
554
555 ! [ LOCALS ]
556 type (DICT_ENTRY_T), pointer :: pTarget
557 integer (c_int) :: iStat
558 logical (c_bool) :: is_fatal_l
559 logical (c_bool) :: empty_entries
560
561 if ( present( is_fatal ) ) then
562 is_fatal_l = is_fatal
563 else
564 is_fatal_l = false
565 endif
566
567 ptarget => this%get_entry(skey)
568
569 if ( associated( ptarget ) ) then
570
571 if ( is_fatal_l ) then
572 empty_entries = ptarget%sl%empty_entries_present()
573 if( empty_entries ) call warn(smessage="There are missing values associated" &
574 //" with the key value of "//dquote(skey), &
575 iloglevel=log_all, lecho=true, lfatal=true )
576 endif
577
578 lvalues = ptarget%sl%get_logical()
579
580 else
581
582 allocate(lvalues(1), stat=istat)
583 call assert(istat == 0, "Failed to allocate memory to lValues array", &
584 __file__, __line__)
585
586 call warn(smessage="Failed to find a dictionary entry associated with key value of "//dquote(skey), &
587 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
588
589 lvalues = false
590
591 endif
592
593 end subroutine get_values_as_logical_sub
594
595!--------------------------------------------------------------------------------------------------
596
597 !> Search through keys for a match; return logical values.
598 !!
599 !! THis routine allows for multiple header values to be supplied
600 !! in the search for the appropriate column.
601 !! @param[in] this Object of DICT_T class.
602 !! @param[in] slKeys String list containing one or more possible key values to
603 !! search for.
604 !! @param[out] iValues Integer vector of values associated with one of the provided keys.
605
606 subroutine get_values_as_logical_given_list_of_keys_sub(this, slKeys, lValues, is_fatal)
607
608 class(dict_t) :: this
609 type (FSTRING_LIST_T) :: slKeys
610 logical (c_bool), allocatable, intent(out) :: lValues(:)
611 logical (c_bool), optional :: is_fatal
612
613 ! [ LOCALS ]
614 type (DICT_ENTRY_T), pointer :: pTarget
615 integer (c_int) :: iStat
616 integer (c_int) :: iCount
617 character (len=:), allocatable :: sText
618 logical (c_bool) :: is_fatal_l
619 logical (c_bool) :: empty_entries
620
621 if ( present( is_fatal ) ) then
622 is_fatal_l = is_fatal
623 else
624 is_fatal_l = false
625 endif
626
627 icount = 0
628
629 do while ( icount < slkeys%count )
630
631 icount = icount + 1
632
633 stext = slkeys%get( icount)
634
635 ptarget => this%get_entry( stext )
636
637 if ( associated( ptarget ) ) exit
638
639 enddo
640
641 if ( associated( ptarget ) ) then
642
643 if ( is_fatal_l ) then
644 empty_entries = ptarget%sl%empty_entries_present()
645 if( empty_entries ) call warn(smessage="There are missing values associated" &
646 //" with the key values of "//slkeys%list_all(), &
647 iloglevel=log_all, lecho=true, lfatal=true )
648 endif
649
650 lvalues = ptarget%sl%get_logical()
651
652 else
653
654 allocate(lvalues(1), stat=istat)
655 call assert(istat == 0, "Failed to allocate memory to lValues array", &
656 __file__, __line__)
657
658 call warn(smessage="Failed to find a dictionary entry associated with key value(s) of: "//dquote(slkeys%list_all()), &
659 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
660
661 lvalues = false
662
663 endif
664
666
667!--------------------------------------------------------------------------------------------------
668
669 subroutine get_values_as_string_list_given_list_of_keys_sub(this, slKeys, slString, is_fatal )
670
671 class(dict_t) :: this
672 type (FSTRING_LIST_T) :: slKeys
673 type ( FSTRING_LIST_T ), intent(out) :: slString
674 logical (c_bool), optional :: is_fatal
675
676 ! [ LOCALS ]
677 type (DICT_ENTRY_T), pointer :: pTarget
678 integer (c_int) :: iStat
679 integer (c_int) :: iCount
680 character (len=:), allocatable :: sText
681 logical (c_bool) :: is_fatal_l
682 logical (c_bool) :: empty_entries
683
684 if ( present( is_fatal ) ) then
685 is_fatal_l = is_fatal
686 else
687 is_fatal_l = false
688 endif
689
690 icount = 0
691
692 do while ( icount < slkeys%count )
693
694 icount = icount + 1
695
696 stext = slkeys%get( icount)
697
698 ptarget => this%get_entry( stext )
699
700 if ( associated( ptarget ) ) exit
701
702 enddo
703
704 if ( associated( ptarget ) ) then
705
706 if ( is_fatal_l ) then
707 empty_entries = ptarget%sl%empty_entries_present()
708 if( empty_entries ) call warn(smessage="There are missing values associated" &
709 //" with the key values of "//slkeys%list_all(), &
710 iloglevel=log_all, lecho=true, lfatal=true )
711 endif
712
713 slstring = ptarget%sl
714
715 else
716
717 call slstring%append("<NA>")
718 call warn(smessage="Failed to find a dictionary entry associated with key value(s) of: "//dquote(slkeys%list_all()), &
719 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
720
721 endif
722
724
725!--------------------------------------------------------------------------------------------------
726
727 subroutine get_value_as_string_sub(this, sText, sKey, iIndex, is_fatal )
728
729 class(dict_t) :: this
730 character(len=:), allocatable, intent(out) :: sText
731 character (len=*), intent(in), optional :: sKey
732 integer (c_int), intent(in), optional :: iIndex
733 logical (c_bool), optional :: is_fatal
734
735 ! [ LOCALS ]
736 type (DICT_ENTRY_T), pointer :: pTarget
737 integer (c_int) :: iStat
738 logical (c_bool) :: is_fatal_l
739 logical (c_bool) :: empty_entries
740
741 if ( present( is_fatal ) ) then
742 is_fatal_l = is_fatal
743 else
744 is_fatal_l = false
745 endif
746
747 if ( present( skey ) ) this%current => this%get_entry(skey)
748
749 if ( present( iindex ) ) this%current => this%get_entry( iindex )
750
751 if ( associated( this%current ) ) then
752
753 stext = this%current%sl%get( 1, this%current%sl%count )
754
755 else
756
757 stext= ""
758
759 endif
760
761 end subroutine get_value_as_string_sub
762
763!--------------------------------------------------------------------------------------------------
764
765 subroutine get_values_as_string_list_sub(this, sKey, slString, is_fatal)
766
767 class(dict_t) :: this
768 character (len=*), intent(in) :: sKey
769 type (FSTRING_LIST_T), intent(out) :: slString
770 logical (c_bool), optional :: is_fatal
771
772 ! [ LOCALS ]
773 type (DICT_ENTRY_T), pointer :: pTarget
774 integer (c_int) :: iStat
775 logical (c_bool) :: is_fatal_l
776 logical (c_bool) :: empty_entries
777
778 ptarget => this%get_entry(skey)
779
780 if ( present(is_fatal) ) then
781 is_fatal_l = is_fatal
782 else
783 is_fatal_l = false
784 endif
785
786 if ( associated( ptarget ) ) then
787
788 if ( is_fatal_l ) then
789 empty_entries = ptarget%sl%empty_entries_present()
790 if( empty_entries ) call warn(smessage="There are missing values associated" &
791 //" with the key value of "//dquote(skey), &
792 iloglevel=log_all, lecho=true, lfatal=true )
793 endif
794
795 slstring = ptarget%sl
796
797 else
798
799 call slstring%append("<NA>")
800 call warn(smessage="Failed to find a dictionary entry associated with key value of "//dquote(skey), &
801 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
802
803 endif
804
805
806 end subroutine get_values_as_string_list_sub
807
808!--------------------------------------------------------------------------------------------------
809
810 !> Search through keys for a match; return integer values.
811 !!
812 !! THis routine allows for multiple header values to be supplied
813 !! in the search for the appropriate column.
814 !! @param[in] this Object of DICT_T class.
815 !! @param[in] slKeys String list containing one or more possible key values to
816 !! search for.
817 !! @param[out] iValues Integer vector of values associated with one of the provided keys.
818
819 subroutine get_values_as_int_given_list_of_keys_sub(this, slKeys, iValues, is_fatal)
820
821 class(dict_t) :: this
822 type (FSTRING_LIST_T) :: slKeys
823 integer (c_int), allocatable, intent(out) :: iValues(:)
824 logical (c_bool), optional :: is_fatal
825
826 ! [ LOCALS ]
827 type (DICT_ENTRY_T), pointer :: pTarget
828 integer (c_int) :: iStat
829 integer (c_int) :: iCount
830 character (len=256) :: sText
831 logical (c_bool) :: is_fatal_l
832 logical (c_bool) :: empty_entries
833
834 if ( present( is_fatal ) ) then
835 is_fatal_l = is_fatal
836 else
837 is_fatal_l = false
838 endif
839
840 icount = 0
841
842 ptarget => null()
843
844 do while ( icount < slkeys%count )
845
846 icount = icount + 1
847
848 stext = slkeys%get( icount)
849
850 ptarget => this%get_entry( stext )
851
852 if ( associated( ptarget ) ) exit
853
854 enddo
855
856 if ( associated( ptarget ) ) then
857
858 if ( is_fatal_l ) then
859 empty_entries = ptarget%sl%empty_entries_present()
860 if( empty_entries ) call warn(smessage="There are missing values associated" &
861 //" with the key values of "//slkeys%list_all(), &
862 iloglevel=log_all, lecho=true, lfatal=true )
863 endif
864
865 ivalues = ptarget%sl%get_integer()
866
867 else
868
869 allocate(ivalues(1), stat=istat)
870 call assert(istat == 0, "Failed to allocate memory to iValues array", &
871 __file__, __line__)
872
873 call warn(smessage="Failed to find a dictionary entry associated with key value(s) of: "//dquote(slkeys%list_all()), &
874 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
875
876 ivalues = itinyval
877
878 endif
879
881
882!--------------------------------------------------------------------------------------------------
883
884 !> Search through keys for a match; return float values.
885 !!
886 !! THis routine allows for multiple header values to be supplied
887 !! in the search for the appropriate column.
888 !! @param[in] this Object of DICT_T class.
889 !! @param[in] slKeys String list containing one or more possible key values to
890 !! search for.
891 !! @param[out] fValues Float vector of values associated with one of the provided keys.
892
893 subroutine get_values_as_float_given_list_of_keys_sub(this, slKeys, fValues, is_fatal)
894
895 class(dict_t) :: this
896 type (FSTRING_LIST_T) :: slKeys
897 real (c_float), allocatable, intent(out) :: fValues(:)
898 logical (c_bool), optional :: is_fatal
899
900 ! [ LOCALS ]
901 type (DICT_ENTRY_T), pointer :: pTarget
902 integer (c_int) :: iStat
903 integer (c_int) :: iCount
904 character (len=:), allocatable :: sText
905 logical (c_bool) :: is_fatal_l
906 logical (c_bool) :: empty_entries
907
908 if ( present( is_fatal ) ) then
909 is_fatal_l = is_fatal
910 else
911 is_fatal_l = false
912 endif
913
914 icount = 0
915 ptarget => null()
916
917 do while ( icount < slkeys%count )
918
919 icount = icount + 1
920
921 stext = slkeys%get( icount)
922
923 ptarget => this%get_entry( stext )
924
925 if ( associated( ptarget ) ) exit
926
927 enddo
928
929 if ( associated( ptarget ) ) then
930
931 if ( is_fatal_l ) then
932 empty_entries = ptarget%sl%empty_entries_present()
933 if( empty_entries ) call warn(smessage="There are missing values associated" &
934 //" with the key values of "//slkeys%list_all(), &
935 iloglevel=log_all, lecho=true, lfatal=true )
936 endif
937
938 fvalues = ptarget%sl%get_float()
939
940 else
941
942 allocate(fvalues(1), stat=istat)
943 call assert(istat == 0, "Failed to allocate memory to fValues array", &
944 __file__, __line__)
945
946 call warn(smessage="Failed to find a dictionary entry associated with key value(s) of: "//dquote(slkeys%list_all()), &
947 smodule=__file__, iline=__line__, iloglevel=log_debug, lecho=false)
948
949
950 fvalues = ftinyval
951
952 endif
953
954
956
957!--------------------------------------------------------------------------------------------------
958
959 subroutine get_values_as_float_sub(this, sKey, fValues, is_fatal)
960
961 class(dict_t) :: this
962 character (len=*), intent(in) :: sKey
963 real (c_float), allocatable, intent(out) :: fValues(:)
964 logical (c_bool), optional :: is_fatal
965
966 ! [ LOCALS ]
967 type (DICT_ENTRY_T), pointer :: pTarget
968 integer (c_int) :: iStat
969 logical (c_bool) :: is_fatal_l
970 logical (c_bool) :: empty_entries
971
972 if ( present( is_fatal ) ) then
973 is_fatal_l = is_fatal
974 else
975 is_fatal_l = false
976 endif
977
978 ptarget => this%get_entry(skey)
979
980 if ( associated( ptarget ) ) then
981
982 if ( is_fatal_l ) then
983 empty_entries = ptarget%sl%empty_entries_present()
984 if( empty_entries ) call warn(smessage="There are missing values associated" &
985 //" with the key value of "//dquote(skey), &
986 iloglevel=log_all, lecho=true, lfatal=true )
987 endif
988
989 fvalues = ptarget%sl%get_float()
990
991 else
992
993 allocate(fvalues(1), stat=istat)
994 call assert(istat == 0, "Failed to allocate memory to iValues array", &
995 __file__, __line__)
996
997 call warn(smessage="Failed to find a dictionary entry associated with key value of "//dquote(skey), &
998 smodule=__file__, iline=__line__, iloglevel=log_debug , lecho=false)
999
1000 fvalues = ftinyval
1001
1002 endif
1003
1004 end subroutine get_values_as_float_sub
1005
1006!--------------------------------------------------------------------------------------------------
1007
1008 subroutine print_all_dictionary_entries_sub(this, iLogLevel, sDescription, lEcho )
1009
1010 class(dict_t) :: this
1011 integer (c_int), intent(in), optional :: iLogLevel
1012 character (len=*), intent(in), optional :: sDescription
1013 logical (c_bool), intent(in), optional :: lEcho
1014
1015 ! [ LOCALS ]
1016 type (DICT_ENTRY_T), pointer :: current
1017 character (len=512) :: sTempBuf
1018 character (len=:), allocatable :: sDescription_
1019 integer (c_int) :: iCount
1020 integer (c_int) :: iIndex
1021 integer (c_int) :: iLogLevel_l
1022 logical (c_bool) :: lEcho_l
1023
1024 if ( present( iloglevel ) ) then
1025 iloglevel_l = iloglevel
1026 else
1027 iloglevel_l = logs%iLogLevel
1028 end if
1029
1030 if ( present( lecho ) ) then
1031 lecho_l = lecho
1032 else
1033 lecho_l = false
1034 end if
1035
1036 if ( present( sdescription ) ) then
1037 sdescription_ = trim(sdescription)
1038 else
1039 sdescription_ = "dictionary data structure"
1040 endif
1041
1042 current => this%first
1043 icount = 0
1044
1045 call logs%write( "### Summary of all items stored in "//trim(sdescription_), &
1046 iloglevel=iloglevel_l, lecho=lecho_l )
1047
1048 do while ( associated( current ) )
1049
1050 icount = icount + 1
1051 stempbuf = current%key
1052
1053 call logs%write( ascharacter(icount)//") KEY: "//dquote(stempbuf), &
1054 iloglevel=iloglevel_l, lecho=lecho_l, itab=2, ilinesbefore=1, ilinesafter=1 )
1055
1056 select case ( iloglevel_l )
1057
1058 case ( log_general )
1059
1060 call current%sl%print( lu=logs%iUnitNum( log_general ) )
1061
1062 case ( log_debug )
1063
1064 call current%sl%print( lu=logs%iUnitNum( log_debug ) )
1065
1066 case ( log_all )
1067
1068 call current%sl%print( lu=logs%iUnitNum( log_general ) )
1069 call current%sl%print( lu=logs%iUnitNum( log_debug ) )
1070
1071 case default
1072
1073 end select
1074
1075 if ( lecho_l ) call current%sl%print()
1076
1077 current => current%next
1078
1079 enddo
1080
1082
1083end module dictionary
This module contains physical constants and convenience functions aimed at performing unit conversion...
logical(c_bool), parameter, public true
logical(c_bool), parameter, public false
integer(c_int), parameter, public itinyval
subroutine get_values_as_string_list_sub(this, skey, slstring, is_fatal)
type(dict_entry_t) function, pointer get_entry_by_key_fn(this, skey)
type(dict_entry_t) function, pointer get_next_entry_by_key_fn(this, skey)
subroutine add_key_sub(this, skey)
subroutine get_values_as_float_given_list_of_keys_sub(this, slkeys, fvalues, is_fatal)
Search through keys for a match; return float values.
type(dict_entry_t) function, pointer get_entry_by_index_fn(this, iindex)
subroutine add_integer_sub(this, ivalue)
subroutine get_values_as_logical_given_list_of_keys_sub(this, slkeys, lvalues, is_fatal)
Search through keys for a match; return logical values.
type(dict_entry_t), pointer, public cf_entry
subroutine print_all_dictionary_entries_sub(this, iloglevel, sdescription, lecho)
subroutine get_values_as_int_given_list_of_keys_sub(this, slkeys, ivalues, is_fatal)
Search through keys for a match; return integer values.
type(dict_entry_t) function, pointer get_next_entry_fn(this)
type(dict_t), public cf_dict
logical(c_bool) function key_name_already_in_use_fn(this, skey)
subroutine delete_entry_by_key_sub(this, skey)
subroutine get_values_as_int_sub(this, skey, ivalues, is_fatal)
subroutine add_double_sub(this, dvalue)
subroutine get_value_as_string_sub(this, stext, skey, iindex, is_fatal)
type(dict_entry_t) function, pointer find_dict_entry_fn(this, ssearchkey)
subroutine get_values_as_string_list_given_list_of_keys_sub(this, slkeys, slstring, is_fatal)
type(fstring_list_t) function grep_dictionary_key_names_fn(this, skey)
subroutine add_entry_to_dict_sub(this, dict_entry)
subroutine add_float_sub(this, fvalue)
subroutine get_values_as_float_sub(this, skey, fvalues, is_fatal)
subroutine add_logical_sub(this, lvalue)
subroutine add_string_sub(this, svalue)
subroutine get_values_as_logical_sub(this, skey, lvalues, is_fatal)
subroutine, public warn(smessage, smodule, iline, shints, lfatal, iloglevel, lecho)
type(logfile_t), public logs
Definition logfiles.F90:62
@ log_general
Definition logfiles.F90:22