TABLE OF CONTENTS
ABINIT/macroin [ Functions ]
NAME
macroin
FUNCTION
Treat "macro" input variables, that can : - initialize several other input variables for one given dataset - initialize several other input variables for a set of datasets. Note that the treatment of these different types of macro input variables is different. Documentation of such input variables is very important, including the proper echo, in the output file, of what such input variables have done. Important information : all the "macro" input variables should be properly identifiable to be so, and it is proposed to make them start with the string "macro".
COPYRIGHT
Copyright (C) 2009-2018 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
ndtset_alloc=number of datasets, corrected for allocation of at least one data set. ecut_tmp(3,2,10)= possible ecut values as read in psp files
OUTPUT
dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables, some of which are given a value here. The dataset with number 0 should NOT be modified in the present routine.
PARENTS
m_ab7_invars_f90
CHILDREN
intagm
SOURCE
42 #if defined HAVE_CONFIG_H 43 #include "config.h" 44 #endif 45 46 #include "abi_common.h" 47 48 49 subroutine macroin(dtsets,ecut_tmp,lenstr,ndtset_alloc,string) 50 51 use defs_basis 52 use defs_abitypes 53 use m_profiling_abi 54 use m_xmpi 55 use m_errors 56 57 use m_parser, only : intagm 58 59 !This section has been created automatically by the script Abilint (TD). 60 !Do not modify the following lines by hand. 61 #undef ABI_FUNC 62 #define ABI_FUNC 'macroin' 63 !End of the abilint section 64 65 implicit none 66 67 !Arguments ------------------------------------ 68 !scalars 69 integer,intent(in) :: ndtset_alloc,lenstr 70 character(len=*),intent(inout) :: string 71 !arrays 72 real(dp),intent(in) :: ecut_tmp(3,2,10) 73 type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc) !vz_i ziontypat 74 75 !Local variables ------------------------------- 76 !scalars 77 integer :: idtset,iatom,jdtset,marr,tread 78 !!arrays 79 integer,allocatable :: intarr(:) 80 real(dp) :: ecutmax(3),ecutdgmax(3) 81 real(dp),allocatable :: dprarr(:) 82 character(len=500) :: message 83 !****************************************************************** 84 85 do idtset=1,ndtset_alloc 86 jdtset=dtsets(idtset)%jdtset 87 if (dtsets(idtset)%macro_uj>0) then 88 dtsets(idtset)%irdwfk = 1 ! preconverged wave function compulsory 89 ! dtsets(idtset)%nline = maxval((/ int(dtsets(idtset)%natom/2) , 6 /)) ! using default value: \DeltaU< 1% 90 ! dtsets(idtset)%nnsclo = 4 ! using default value: \DeltaU< 1% 91 dtsets(idtset)%tolvrs = 10d-8 ! convergence on the potential; 10d-8^= 10d-5 on occupation 92 dtsets(idtset)%diemix = 0.45_dp ! fastest convergence: dn= E^(-istep * 0.229 ) 93 dtsets(idtset)%dmatpuopt= 3 ! normalization of the occupation operator 94 ! dtsets(idtset)%nstep = 255 ! expected convergence after 10 \pm 3, 30 as in default normally suficient 95 ! dtsets(idtset)%iscf = 17 ! mixing on potential, 17: default for PAW 96 end if ! macro_uj 97 98 !Read parameters 99 marr=dtsets(idtset)%npsp;if (dtsets(idtset)%npsp<3) marr=3 100 marr=max(marr,dtsets(idtset)%nimage) 101 ABI_ALLOCATE(intarr,(marr)) 102 ABI_ALLOCATE(dprarr,(marr)) 103 104 call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),"accuracy",tread,'INT') 105 106 ecutmax=-one 107 ecutdgmax=-one 108 do iatom=1,dtsets(idtset)%natom 109 ecutmax(:)=max(ecutmax(:),ecut_tmp(:,1,dtsets(idtset)%typat(iatom))) 110 ecutdgmax(:)=max(ecutdgmax(:),ecut_tmp(:,2,dtsets(idtset)%typat(iatom))) 111 end do 112 113 if(tread==1) then 114 dtsets(idtset)%accuracy=intarr(1) 115 if (dtsets(idtset)%accuracy==1) then 116 if (ecutmax(1)>zero) dtsets(idtset)%ecut=ecutmax(1) 117 if (ecutdgmax(1)>zero.and.dtsets(idtset)%usepaw==1) dtsets(idtset)%pawecutdg=ecutdgmax(1) 118 dtsets(idtset)%boxcutmin=1.5_dp 119 if (dtsets(idtset)%usepaw==1) then 120 dtsets(idtset)%bxctmindg=1.5_dp 121 dtsets(idtset)%pawxcdev=1 122 dtsets(idtset)%pawmixdg=0 123 dtsets(idtset)%pawovlp=10 124 dtsets(idtset)%pawnhatxc=0 125 dtsets(idtset)%mqgriddg=0 126 end if 127 dtsets(idtset)%mqgrid=0 128 dtsets(idtset)%tolimg=5.0d-5 129 dtsets(idtset)%tolvrs=tol3 130 dtsets(idtset)%tolmxf=1.0d-3 131 dtsets(idtset)%toldff=zero 132 dtsets(idtset)%optforces=1 133 dtsets(idtset)%timopt=0 134 dtsets(idtset)%npulayit=4 135 dtsets(idtset)%nstep=30 136 dtsets(idtset)%prteig=0 137 dtsets(idtset)%prtden=0 138 else if (dtsets(idtset)%accuracy==2) then 139 if (ecutmax(2)>zero) dtsets(idtset)%ecut=ecutmax(2) 140 if (ecutdgmax(2)>zero.and.dtsets(idtset)%usepaw==1) dtsets(idtset)%pawecutdg=ecutdgmax(2) 141 dtsets(idtset)%boxcutmin=1.8_dp 142 if (dtsets(idtset)%usepaw==1) then 143 dtsets(idtset)%bxctmindg=1.8_dp 144 dtsets(idtset)%pawxcdev=1 145 dtsets(idtset)%pawmixdg=0 146 dtsets(idtset)%pawovlp=7 147 dtsets(idtset)%pawnhatxc=1 148 dtsets(idtset)%mqgriddg=0 149 end if 150 dtsets(idtset)%mqgrid=0 151 dtsets(idtset)%tolimg=5.0d-5 152 dtsets(idtset)%tolvrs=tol5 153 dtsets(idtset)%tolmxf=5.0d-4 154 dtsets(idtset)%toldff=zero 155 dtsets(idtset)%optforces=1 156 dtsets(idtset)%timopt=0 157 dtsets(idtset)%npulayit=7 158 dtsets(idtset)%nstep=30 159 dtsets(idtset)%prteig=0 160 dtsets(idtset)%prtden=0 161 else if (dtsets(idtset)%accuracy==3) then 162 if (ecutmax(2)>zero) dtsets(idtset)%ecut=ecutmax(2) 163 if (ecutdgmax(2)>zero.and.dtsets(idtset)%usepaw==1) dtsets(idtset)%pawecutdg=ecutdgmax(2) 164 dtsets(idtset)%boxcutmin=1.8_dp 165 if (dtsets(idtset)%usepaw==1) then 166 dtsets(idtset)%bxctmindg=1.8_dp 167 dtsets(idtset)%pawxcdev=1 168 dtsets(idtset)%pawmixdg=0 169 dtsets(idtset)%pawovlp=7 170 dtsets(idtset)%pawnhatxc=1 171 dtsets(idtset)%mqgriddg=0 172 end if 173 dtsets(idtset)%mqgrid=0 174 dtsets(idtset)%tolimg=5.0d-5 175 dtsets(idtset)%tolvrs=tol7 176 dtsets(idtset)%tolmxf=1.0d-4 177 dtsets(idtset)%toldff=zero 178 dtsets(idtset)%optforces=2 179 dtsets(idtset)%timopt=1 180 if(xmpi_paral==1) dtsets(idtset)%timopt = 0 181 dtsets(idtset)%npulayit=7 182 dtsets(idtset)%nstep=30 183 dtsets(idtset)%prteig=1 184 dtsets(idtset)%prtden=1 185 else if (dtsets(idtset)%accuracy==4) then 186 if (ecutmax(3)>zero) dtsets(idtset)%ecut=ecutmax(3) 187 if (ecutdgmax(3)>zero.and.dtsets(idtset)%usepaw==1) dtsets(idtset)%pawecutdg=ecutdgmax(3) 188 dtsets(idtset)%boxcutmin=two 189 if (dtsets(idtset)%usepaw==1) then 190 dtsets(idtset)%bxctmindg=two 191 dtsets(idtset)%pawxcdev=1 192 dtsets(idtset)%pawmixdg=0 193 dtsets(idtset)%pawovlp=5 194 dtsets(idtset)%pawnhatxc=1 195 dtsets(idtset)%mqgriddg=0 196 end if 197 dtsets(idtset)%mqgrid=0 198 dtsets(idtset)%tolimg=5.0d-5 199 dtsets(idtset)%tolvrs=tol9 200 dtsets(idtset)%tolmxf=5.0d-5 201 dtsets(idtset)%toldff=zero 202 dtsets(idtset)%optforces=2 203 dtsets(idtset)%timopt=1 204 if(xmpi_paral==1) dtsets(idtset)%timopt = 0 205 dtsets(idtset)%npulayit=7 206 dtsets(idtset)%nstep=30 207 dtsets(idtset)%prteig=1 208 dtsets(idtset)%prtden=1 209 else if (dtsets(idtset)%accuracy==5) then 210 if (ecutmax(2)>zero) dtsets(idtset)%ecut=ecutmax(2) 211 if (ecutdgmax(2)>zero.and.dtsets(idtset)%usepaw==1) dtsets(idtset)%pawecutdg=ecutdgmax(2) 212 dtsets(idtset)%boxcutmin=two 213 if (dtsets(idtset)%usepaw==1) then 214 dtsets(idtset)%bxctmindg=two 215 dtsets(idtset)%pawxcdev=2 216 dtsets(idtset)%pawmixdg=1 217 dtsets(idtset)%pawovlp=5 218 dtsets(idtset)%pawnhatxc=1 219 dtsets(idtset)%mqgriddg=0 220 end if 221 dtsets(idtset)%mqgrid=0 222 dtsets(idtset)%tolimg=5.0d-5 223 dtsets(idtset)%tolvrs=tol10 224 dtsets(idtset)%tolmxf=1.0d-6 225 dtsets(idtset)%toldff=zero 226 dtsets(idtset)%optforces=2 227 dtsets(idtset)%timopt=1 228 if(xmpi_paral==1) dtsets(idtset)%timopt = 0 229 dtsets(idtset)%npulayit=15 230 dtsets(idtset)%nstep=50 231 dtsets(idtset)%prteig=1 232 dtsets(idtset)%prtden=1 233 else if (dtsets(idtset)%accuracy==6) then 234 if (ecutmax(3)>zero) dtsets(idtset)%ecut=ecutmax(3) 235 if (ecutdgmax(3)>zero.and.dtsets(idtset)%usepaw==1) dtsets(idtset)%pawecutdg=ecutdgmax(3) 236 dtsets(idtset)%boxcutmin=two 237 if (dtsets(idtset)%usepaw==1) then 238 dtsets(idtset)%bxctmindg=two 239 dtsets(idtset)%pawxcdev=2 240 dtsets(idtset)%pawmixdg=1 241 dtsets(idtset)%pawovlp=5 242 dtsets(idtset)%pawnhatxc=1 243 dtsets(idtset)%mqgriddg=0 244 end if 245 dtsets(idtset)%mqgrid=0 246 dtsets(idtset)%tolimg=5.0d-5 247 dtsets(idtset)%tolvrs=tol12 248 dtsets(idtset)%tolmxf=1.0d-6 249 dtsets(idtset)%toldff=zero 250 dtsets(idtset)%optforces=2 251 dtsets(idtset)%timopt=1 252 if(xmpi_paral==1) dtsets(idtset)%timopt = 0 253 dtsets(idtset)%npulayit=15 254 dtsets(idtset)%nstep=50 255 dtsets(idtset)%prteig=1 256 dtsets(idtset)%prtden=1 257 elseif(dtsets(idtset)%accuracy>6)then 258 write(message, '(a,a,a)' )& 259 & 'accuracy >6 is forbiden !',ch10,& 260 & 'Action : check your input data file.' 261 MSG_ERROR(message) 262 end if 263 else 264 if (ecutmax(3)>zero) dtsets(idtset)%ecut=ecutmax(3) 265 end if 266 ABI_DEALLOCATE(intarr) 267 ABI_DEALLOCATE(dprarr) 268 end do 269 270 end subroutine macroin