TABLE OF CONTENTS


ABINIT/inqpt [ Functions ]

[ Top ] [ Functions ]

NAME

 inqpt

FUNCTION

 Initialize the q point
 for one particular dataset, characterized by jdtset.

COPYRIGHT

 Copyright (C) 2011-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

  chksymbreak= if 1, will check whether the k point grid is symmetric, and stop if not.
  iout=unit number for echoed output
  jdtset=number of the dataset looked for
  lenstr=actual length of the string
  msym=default maximal number of symmetries
  natom=number of atoms
  rprimd(3,3)=dimensional real space primitive translations (bohr)
  spinat(3,1:natom)=spin-magnetization of the atoms
  string*(*)=character string containing all the input data.
             Initialized previously in instrng.
  typat(natom)=type for each atom
  vacuum(3)=for each direction, 0 if no vacuum, 1 if vacuum
  xred(3,natom,nimage) =reduced coordinates of atoms

OUTPUT

  qptn(3)=reduced coordinates of eventual q point (normalisation is already included)
  kptrlatt(3,3)=q-point lattice specification (if kptopt/=0)
  wtqc=weigth of the eventual current q point

SIDE EFFECTS

PARENTS

      invars1

CHILDREN

      getkgrid,intagm,metric,symfind,symlatt

SOURCE

 48 #if defined HAVE_CONFIG_H
 49 #include "config.h"
 50 #endif
 51 
 52 #include "abi_common.h"
 53 
 54 subroutine inqpt(chksymbreak,iout,jdtset,lenstr,msym,natom,qptn,wtqc,rprimd,spinat,string,typat,vacuum,xred,qptrlatt)
 55 
 56  use defs_basis
 57  use m_errors
 58  use m_profiling_abi
 59 
 60  use m_geometry,     only : metric
 61  use m_parser,  only : intagm
 62 
 63 !This section has been created automatically by the script Abilint (TD).
 64 !Do not modify the following lines by hand.
 65 #undef ABI_FUNC
 66 #define ABI_FUNC 'inqpt'
 67  use interfaces_41_geometry
 68  use interfaces_56_recipspace
 69 !End of the abilint section
 70 
 71  implicit none
 72 
 73 !Arguments ------------------------------------
 74 !scalars
 75  integer,intent(in)   :: chksymbreak,iout,jdtset,lenstr,msym,natom
 76  real(dp),intent(inout) :: wtqc
 77  character(len=*),intent(in) :: string
 78 !arrays
 79  integer,intent(in) :: typat(natom),vacuum(3)
 80  real(dp),intent(out) :: qptn(3)
 81  integer,intent(inout) :: qptrlatt(3,3) !vz_i
 82  real(dp),intent(in) :: rprimd(3,3)
 83  real(dp),intent(in) :: spinat(3,natom)
 84  real(dp),intent(in) :: xred(3,natom)
 85 
 86 !Local variables-------------------------------
 87 !scalars
 88  integer :: ii,iqpt,iscf_fake,marr,nptsym,nqpt,nqpt_computed,nshiftq,nsym_new,qptopt
 89  integer :: tread,tread_qptrlatt,tread_ngqpt,use_inversion
 90  real(dp) :: qptnrm,qptrlen,tolsym,ucvol
 91  character(len=500) :: message
 92 !arrays
 93  integer :: bravais(11)
 94  integer :: ngqpt(3)
 95  integer, allocatable :: symafm_new(:)
 96  integer, allocatable :: ptsymrel(:,:,:),symrel_new(:,:,:)
 97  integer,allocatable :: intarr(:)
 98  real(dp) :: gmet(3,3),gprimd(3,3),qpt(3),rmet(3,3),shiftq(3,210)
 99  real(dp),allocatable :: qpts(:,:),tnons_new(:,:),wtq(:)
100  real(dp),allocatable :: dprarr(:)
101 
102 ! *************************************************************************
103 
104 !Compute the maximum size of arrays intarr and dprarr (nshiftq is 210 at maximum)
105  marr=630
106  ABI_ALLOCATE(intarr,(marr))
107  ABI_ALLOCATE(dprarr,(marr))
108 
109 !Find the method to generate the k points
110  qptopt=0
111  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'qptopt',tread,'INT')
112  if(tread==1)qptopt=intarr(1)
113 
114  if(qptopt==0)then
115 
116 !  Read qpt and qptnrm
117    qpt=zero
118    call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'qpt',tread,'DPR')
119    if(tread==1) qpt(1:3)=dprarr(1:3)
120 
121    qptnrm=1.0_dp
122    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'qptnrm',tread,'DPR')
123 
124    if(tread==1) qptnrm=dprarr(1)
125    if(qptnrm<tol10)then
126      write(message, '(5a)' )&
127 &     'The input variable qptnrm is lower than 1.0d-10,',ch10,&
128 &     'while it must be a positive, non-zero number.   ',ch10,&
129 &     'Action: correct the qptnrm in the input file.'
130      MSG_ERROR(message)
131    end if
132 
133    qptn(:)=qpt(:)/qptnrm
134 
135 !  DBSP: one could want ot define wtq in order to reproduce what is obtained
136 !  with ngqpt but without having to do initialize the qgrid (extremly slow in
137 !  case of large grid > 50x50x50
138    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'wtq',tread,'DPR')
139    if(tread==1) wtqc=dprarr(1)
140 
141  else if(qptopt>=1 .and. qptopt<=4)then
142 
143    ngqpt(:)=0
144    call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'ngqpt',tread_ngqpt,'INT')
145 
146    if(tread_ngqpt==1)then
147      ngqpt(1:3)=intarr(1:3)
148      do ii=1,3
149        if(ngqpt(ii)<1)then
150          write(message,  '(a,i0,a,a,a,i0,a,a,a)' ) &
151 &         'The input variable ngqpt(',ii,') must be strictly positive,',ch10,&
152 &         'while it is found to be',ngqpt(ii),'.',ch10,&
153 &         'Action: change it in your input file, or change qptopt.'
154          MSG_ERROR(message)
155        end if
156      end do
157    end if
158 
159    call intagm(dprarr,intarr,jdtset,marr,9,string(1:lenstr),'qptrlatt',tread_qptrlatt,'INT')
160    if(tread_qptrlatt==1) qptrlatt(:,:)=reshape(intarr(1:9), (/3,3/) )
161 
162    if(tread_ngqpt==1 .and. tread_qptrlatt==1)then
163      write(message,  '(a,a,a,a,a)' ) &
164 &     'The input variables ngqpt and qptrlatt cannot both ',ch10,&
165 &     'be defined in the input file.',ch10,&
166 &     'Action: change one of ngqpt or qptrlatt in your input file.'
167      MSG_ERROR(message)
168    else if(tread_ngqpt==1)then
169      qptrlatt(:,:)=0
170      qptrlatt(1,1)=ngqpt(1)
171      qptrlatt(2,2)=ngqpt(2)
172      qptrlatt(3,3)=ngqpt(3)
173    end if
174 
175    nshiftq=1
176    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nshiftq',tread,'INT')
177    if(tread==1)nshiftq=intarr(1)
178 
179    if(nshiftq<1 .or. nshiftq>210 )then
180      write(message,  '(3a,i4,a,a,a)' )&
181 &     'The only allowed values of nshiftq are between 1 and 210,',ch10,&
182 &     'while it is found to be',nshiftq,'.',ch10,&
183 &     'Action: change the value of nshiftq in your input file, or change qptopt.'
184      MSG_ERROR(message)
185    end if
186 
187    shiftq=half
188    call intagm(dprarr,intarr,jdtset,marr,3*nshiftq,string(1:lenstr),'shiftq',tread,'DPR')
189 
190    if(tread==1)then
191      shiftq(:,1:nshiftq)=reshape( dprarr(1:3*nshiftq), (/3,nshiftq/) )
192    else
193      if(nshiftq/=1)then
194        write(message,  '(3a,i4,a,a)' )&
195 &       'When nshiftq is not equal to 1, shiftq must be defined in the input file.',ch10,&
196 &       'However, shiftq is not defined, while nshiftq=',nshiftq,ch10,&
197 &       'Action: change the value of nshiftq in your input file, or define shiftq.'
198        MSG_ERROR(message)
199      end if
200    end if
201 
202 !  Re-generate symmetry operations from the lattice and atomic coordinates
203 !  This is a fundamental difference with respect to the k point generation.
204    tolsym=tol8
205    ABI_ALLOCATE(ptsymrel,(3,3,msym))
206    ABI_ALLOCATE(symafm_new,(msym))
207    ABI_ALLOCATE(symrel_new,(3,3,msym))
208    ABI_ALLOCATE(tnons_new,(3,msym))
209    call symlatt(bravais,msym,nptsym,ptsymrel,rprimd,tolsym)
210    use_inversion=1
211    call metric(gmet,gprimd,-1,rmet,rprimd,ucvol)
212    call symfind(0,(/zero,zero,zero/),gprimd,0,msym,natom,0,nptsym,nsym_new,0,&
213 &   ptsymrel,spinat,symafm_new,symrel_new,tnons_new,tolsym,typat,use_inversion,xred)
214 
215 !  Prepare to compute the q-point grid in the ZB or IZB
216    iscf_fake=0 ! Do not need the weights
217 
218 !  Compute the maximum number of q points
219    nqpt=0
220    ABI_ALLOCATE(qpts,(3,nqpt))
221    ABI_ALLOCATE(wtq,(nqpt))
222    call getkgrid(chksymbreak,0,iscf_fake,qpts,qptopt,qptrlatt,qptrlen,&
223 &   msym,nqpt,nqpt_computed,nshiftq,nsym_new,rprimd,&
224 &   shiftq,symafm_new,symrel_new,vacuum,wtq)
225    nqpt=nqpt_computed
226    ABI_DEALLOCATE(qpts)
227    ABI_DEALLOCATE(wtq)
228 
229 !  Find the index of the q point within the set of q points that will be generated
230    iqpt=0
231    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'iqpt',tread,'INT')
232    if(tread==1)iqpt=intarr(1)
233 
234 !  Checks that iqpt is among the computed q points
235    if(iqpt<0)then
236      write(message, '(a,i5,3a)' )&
237 &     'The input variable iqpt,',iqpt,' is negative, while it should be 0 or positive.',ch10,&
238 &     'Action: correct iqpt in the input file.'
239      MSG_ERROR(message)
240    end if
241 
242    if(iqpt>nqpt_computed)then
243      write(message, '(a,i5,3a,i5,7a)' )&
244 &     'The input variable iqpt,',iqpt,' is bigger than the computed number of q-points in the grid,',ch10,&
245 &     'which is ',nqpt,'.',ch10,&
246 &     'The latter has been computed from the input variables qptrlatt, ngqpt, nshiftq,',ch10,&
247 &     'shiftq, as well as qptopt, the symmetries of the lattice, and spinat.',ch10,&
248 &     'Action: correct iqpt in the input file, or correct the computed q-point grid.'
249      MSG_ERROR(message)
250    end if
251 
252 !  Compute the q-point grid in the BZ or the IBZ
253    ABI_ALLOCATE(qpts,(3,nqpt))
254    ABI_ALLOCATE(wtq,(nqpt))
255 
256    call getkgrid(chksymbreak,iout,iscf_fake,qpts,qptopt,qptrlatt,qptrlen,&
257 &   msym,nqpt,nqpt_computed,nshiftq,nsym_new,rprimd,&
258 &   shiftq,symafm_new,symrel_new,vacuum,wtq)
259 
260 !  Transfer to qptn, and deallocate
261    qptn(:)=zero
262    if(iqpt/=0)then
263      qptn(:)=qpts(:,iqpt)
264      wtqc = wtq(iqpt)
265    end if
266 
267    ABI_DEALLOCATE(ptsymrel)
268    ABI_DEALLOCATE(symafm_new)
269    ABI_DEALLOCATE(symrel_new)
270    ABI_DEALLOCATE(tnons_new)
271    ABI_DEALLOCATE(qpts)
272    ABI_DEALLOCATE(wtq)
273 
274  else
275 
276    write(message,  '(a,a,a,i4,a,a,a)' ) &
277 &   'The only values of qptopt allowed are smaller than 4.',ch10,&
278 &   'The input value of qptopt is',qptopt,'.',ch10,&
279 &   'Action: change qptopt in your input file.'
280    MSG_ERROR(message)
281  end if
282 
283  ABI_DEALLOCATE(intarr)
284  ABI_DEALLOCATE(dprarr)
285 
286 end subroutine inqpt