XBeach
|
00001 module getkey_module 00002 ! This module makes parameters available through introspection/reflection 00003 ! You can ask for all available parameters using getkeys 00004 ! You can ask for the value of a specific parameter using getkey 00005 ! The naming is still a bit inconsistent so that might change in the future 00006 use typesandkinds, only: slen 00007 use mnemmodule, only: maxrank 00008 implicit none 00009 private 00010 public parameter, getkey, getkeys, getkey_indextype, getnkeys 00011 save 00012 ! analogue to the arraytype in mnemonic 00013 type parameter 00014 character type ! 'c', 'i' or 'r': integer or real*8 00015 integer rank ! 0,1,2,3,4 00016 character(slen) :: name ! 'v','ve', ..... 00017 character(slen) :: units ! m, following udunits convention 00018 character(slen) :: description 00019 character(slen), dimension(maxrank) :: dimensions ! the dimensions of the variable, for example (s%nx, s%ny) 00020 00021 character(slen), pointer :: c0 ! pointer to characters 00022 character(slen), dimension(:), pointer :: c1 ! pointer to array of characters 00023 00024 real*8, pointer :: r0 ! pointer to real8 scalar 00025 real*8, dimension(:), pointer :: r1 ! pointer to real8 (:) 00026 00027 integer, pointer :: i0 ! pointer to integer scalar 00028 integer, dimension(:), pointer :: i1 ! pointer to integer (:) 00029 end type parameter 00030 00031 contains 00032 ! Lookup the index and valuetype of a parameter, currently only works for rank 0 00033 subroutine getkey_indextype(par, key, index, type) 00034 use params, only: parameters 00035 implicit none 00036 type(parameters), intent(in) :: par 00037 character(len=*), intent(in) :: key 00038 integer, intent(out) :: index 00039 character, intent(out) :: type 00040 integer :: i 00041 include 'getkey.inc' 00042 index = -1 00043 type = '' 00044 do i=1,nintegerkeys 00045 if (integerkeys(i) == key) then 00046 index = i 00047 type = 'i' 00048 return 00049 end if 00050 end do 00051 do i=1,nrealkeys 00052 if (realkeys(i) == key) then 00053 index = i 00054 type = 'r' 00055 return 00056 end if 00057 end do 00058 do i=1,ncharacterkeys 00059 if (characterkeys(i) == key) then 00060 index = i 00061 type = 'c' 00062 return 00063 end if 00064 end do 00065 end subroutine getkey_indextype 00066 00067 integer function getkey(par, key, value) 00068 ! This subroutine returns a value (parameter), given a key (parameter name). 00069 use params, only: parameters 00070 use mnemmodule 00071 implicit none 00072 00073 type(parameters), intent(in) :: par 00074 character(len=*), intent(in) :: key 00075 type(parameter), intent(inout) :: value 00076 ! Ok this is ugly, I changed the targets to save, so the pointer is still available after the function ends 00077 ! But this means that if you call getkey twice, then your old value will now refer to the new value 00078 character(slen),target, save :: charvalue 00079 integer, target, save :: intvalue 00080 real*8, target, save :: realvalue 00081 integer :: index 00082 character :: type 00083 include 'getkey.inc' 00084 ! This sets index and type 00085 call getkey_indextype(par, key, index, type) 00086 getkey = -1 00087 if (index .eq. -1 ) return 00088 value%name = key 00089 value%rank = 0 00090 value%type = type 00091 if (type == 'c') then 00092 charvalue = charactervalues(index) 00093 value%c0 => charvalue 00094 elseif (type == 'i') then 00095 intvalue = integervalues(index) 00096 value%i0 => intvalue 00097 elseif (type == 'r') then 00098 realvalue = realvalues(index) 00099 value%r0 => realvalue 00100 end if 00101 getkey = 0 00102 end function getkey 00103 00104 subroutine getnkeys(par, n) 00105 use params, only: parameters 00106 implicit none 00107 type(parameters), intent(in) :: par 00108 integer, intent(out) :: n 00109 include 'getkey.inc' 00110 n = ncharacterkeys+nintegerkeys+nrealkeys 00111 end subroutine getnkeys 00112 00113 subroutine getkeys(par, keys) 00114 ! This subroutine gives a list of all available keys (parameter names) 00115 use params, only: parameters 00116 implicit none 00117 type(parameters), intent(in) :: par 00118 character(slen), dimension(:), allocatable, intent(out) :: keys 00119 00120 include 'getkey.inc' 00121 allocate(keys(ncharacterkeys+nintegerkeys+nrealkeys)) 00122 keys(1:ncharacterkeys) = characterkeys(1:ncharacterkeys) 00123 keys(ncharacterkeys+1:ncharacterkeys+nintegerkeys) = integerkeys(1:nintegerkeys) 00124 keys(ncharacterkeys+nintegerkeys+1:ncharacterkeys+nintegerkeys+nrealkeys) = realkeys(1:nrealkeys) 00125 00126 00127 end subroutine getkeys 00128 00129 00130 00131 00132 00133 end module getkey_module 00134