TABLE OF CONTENTS
ABINIT/find_getdtset [ Functions ]
NAME
find_getdtset
FUNCTION
Find the number of the dataset (iget) for a given value of a "get" variable (getvalue) of name getname, given the number of the current dataset (idtset). Also find the coefficients of mixing of the images of the old dataset, to initialize the new dataset images (use a linear interpolation)
COPYRIGHT
Copyright (C) 2009-2017 ABINIT group (XG) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
dtsets(0:ndtset_alloc)=<type datasets_type>contains all input variables getvalue=value of the get variable getname=name of the get variable idtset=number of the current dataset mxnimage=dimension of miximage ndtset_alloc=dimension of dtsets
OUTPUT
iget=number of the dataset from which the value must be get, 0 if the data should not be got from another dataset miximage(mxnimage,mxnimage)=coefficients of mixing of the images of the old dataset, to initialize the new dataset images
PARENTS
driver
CHILDREN
wrtout
SOURCE
39 #if defined HAVE_CONFIG_H 40 #include "config.h" 41 #endif 42 43 #include "abi_common.h" 44 45 46 subroutine find_getdtset(dtsets,getvalue,getname,idtset,iget,miximage,mxnimage,ndtset_alloc) 47 48 use defs_basis 49 use defs_abitypes 50 use m_profiling_abi 51 use m_errors 52 53 !This section has been created automatically by the script Abilint (TD). 54 !Do not modify the following lines by hand. 55 #undef ABI_FUNC 56 #define ABI_FUNC 'find_getdtset' 57 use interfaces_14_hidewrite 58 !End of the abilint section 59 60 implicit none 61 62 !Arguments ------------------------------------ 63 !scalars 64 integer, intent(in) :: getvalue,idtset,mxnimage,ndtset_alloc 65 integer, intent(out) :: iget 66 real(dp), intent(out) :: miximage(mxnimage,mxnimage) 67 character(len=*),intent(in) :: getname 68 type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc) 69 70 !Local variables------------------------------- 71 integer :: iimage 72 real(dp) :: newimage_get,ratio 73 character(len=500) :: message 74 75 ! ************************************************************************* 76 77 iget=0 78 if(getvalue>0 .or. (getvalue<0 .and. idtset+getvalue>0) )then 79 ! In case getvalue is a negative number (so must add to idtset) 80 if(getvalue<0 .and. idtset+getvalue>0) iget=idtset+getvalue 81 if(getvalue>0)then 82 do iget=1,idtset 83 if( dtsets(iget)%jdtset==getvalue )exit 84 end do 85 if(iget==idtset)then 86 ! The index of the dataset, from which the data ought to be taken, 87 ! does not correspond to a previous dataset. 88 write(message, '(a,i3,4a,i3,7a)' )& 89 & 'The component number',idtset,' of the input variable ',trim(getname),',',& 90 & ' equal to',getvalue,',',ch10,& 91 & 'does not correspond to an existing index.',ch10,& 92 & 'Action : correct ',trim(getname),' or jdtset in your input file.' 93 MSG_ERROR(message) 94 end if 95 end if 96 write(message, '(3a,i3,2a)' )& 97 & ' find_getdtset : ',trim(getname),'/=0, take data from output of dataset with index',dtsets(iget)%jdtset,'.',ch10 98 call wrtout(ab_out,message,'COLL') 99 call wrtout(std_out,message,'COLL') 100 end if 101 102 !For the time being, uses a simple interpolation when the images do not match. If only one image, take the first get image. 103 miximage(:,:)=zero 104 if(dtsets(idtset)%nimage==1)then 105 miximage(1,1)=one 106 else 107 do iimage=1,dtsets(idtset)%nimage 108 ratio=(iimage-one)/real(dtsets(idtset)%nimage-one) 109 newimage_get=one+ratio*(dtsets(iget)%nimage-one) 110 if(abs(newimage_get-nint(newimage_get))<tol8)then 111 miximage(iimage,nint(newimage_get))=one 112 else 113 miximage(iimage,floor(newimage_get))=one-(newimage_get-floor(newimage_get)) 114 miximage(iimage,ceiling(newimage_get))=one-miximage(iimage,floor(newimage_get)) 115 end if 116 end do 117 end if 118 119 end subroutine find_getdtset