TABLE OF CONTENTS


ABINIT/m_fold2block [ Modules ]

[ Top ] [ Modules ]

NAME

  m_fold2block

FUNCTION

  This module contains basic tools to operate on vectors expressed in reduced coordinates.

COPYRIGHT

 Copyright (C) 2008-2018 ABINIT group (AB)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 MODULE m_fold2block
24 
25  use defs_basis
26  use m_abicore
27  use m_errors
28 
29  implicit none
30 
31  public
32  !private
33 
34 CONTAINS  !===========================================================

m_fold2block/getargs [ Functions ]

[ Top ] [ m_fold2block ] [ Functions ]

NAME

 getargs

FUNCTION

 Part of fold2Bloch that interprets the command line
 arguments and checks for their corectness.

INPUTS

 fname: Name of input file
 folds: Number of folds in x, y, and z directions.

OUTPUT

 folds: Number of folds in x, y, and z directions.
 Returns the folds array filled with each element.

PARENTS

      fold2Bloch

CHILDREN

SOURCE

273 subroutine getargs(folds, fname)
274 
275 
276 !This section has been created automatically by the script Abilint (TD).
277 !Do not modify the following lines by hand.
278 #undef ABI_FUNC
279 #define ABI_FUNC 'getargs'
280 !End of the abilint section
281 
282  implicit none
283 
284 !Arguments ------------------------------------
285 !scalars
286  character(len=fnlen), intent(inout) :: fname
287 !arrays
288  integer, intent(inout) :: folds(3)
289 
290 !Local variables-------------------------------
291 !scalars
292  integer :: num_args, argcount, ii, ios
293  character(len=fnlen) :: argfolds
294  logical :: dir
295 !arrays
296  character(len=fnlen), allocatable :: args(:)
297 
298 ! *************************************************************************
299 
300  !get number of args
301  num_args=command_argument_count()
302  ABI_ALLOCATE(args,(num_args))
303 
304  !Check for errors in arguments
305  if (num_args>2) then
306    write(std_out,*) "     Too many arguments."
307    write(std_out,*) "     Usage: $fold2Bloch file_WFK x:y:z (folds)"
308    MSG_ERROR("Aborting now")
309  elseif (num_args<2) then
310    if (num_args==1) then
311      call get_command_argument(1,args(1))
312      if ((args(1)=="-h").or.(args(1)=="--help")) then
313        write(std_out,*) "     Usage: fold2Bloch file_WFK x:y:z (folds)"
314        write(std_out,*) "     file_WFK is the WFK file name (ex. Ex_WFK)"
315        write(std_out,*) "     x:y:z integers, greater than 0 that represent a multiplicity"
316        write(std_out,*) "     in the corresponding directions used when constructing the supercell."
317        MSG_ERROR("Aborting now")
318      else
319        write(std_out,*) "     Not all arguments are present."
320        write(std_out,*) "     Make sure that file name and number of folds are indicated."
321        write(std_out,*) "     Usage: $fold2Bloch file_WFK x:y:z (folds)"
322        MSG_ERROR("Aborting now")
323      end if
324    else
325      write(std_out,*) "     Not all arguments are present."
326      write(std_out,*) "     Make sure that file name and number of folds are indicated."
327      write(std_out,*) "     Usage: $fold2Bloch file_WFK x:y:z (folds)"
328      MSG_ERROR("Aborting now")
329    end if
330  else
331    do argcount=1, num_args
332      call get_command_argument(argcount,args(argcount))
333      if (argcount==1) then
334        read(args(argcount), "(a)") fname !read first argument
335      else
336        read(args(argcount), "(a)") argfolds !read second argument
337      end if
338    end do
339  end if
340 
341  ! Does intput file exist?
342  inquire(file=fname, exist=dir)
343  if (.not.(dir)) then
344    write(std_out,*) "     Case file not found: ", trim(fname)
345    write(std_out,*) "     Usage: $fold2Bloch file_WFK x:y:z (folds)"
346    MSG_ERROR("Aborting now")
347  end if
348 
349  !Was the number of folds entered in correct format?
350  ii=0
351  ii=INDEX(argfolds, ':') !look for first ":" to differentiate between axis
352  if (ii==0) then
353    write(std_out,*) "     Unknown number of folds. See below or type:"
354    write(std_out,*) "     fold2Bloch <-h> or fold2Bloch <--help> for more information."
355    write(std_out,*) '     Usage: $fold2Bloch file_WFK x:y:z (folds)'
356    MSG_ERROR("Aborting now")
357  end if
358  read (argfolds(1:ii-1), *, iostat=ios) folds(1) !read X folds
359  if ((ios/=0).or.(folds(1)<=0)) then
360    MSG_ERROR('Number of folds has to be a positive integer greater than 0')
361  end if
362  argfolds=argfolds(ii+1:) !Start argfolds from the first ":"
363  ii=0
364  ii=INDEX(argfolds, ':') !look for second ":"
365  if (ii==0) then
366    write(std_out,*) '     Unknown number of folds. See below or type:'
367    write(std_out,*) "     fold2Bloch <-h> or fold2Bloch <--help> for more information."
368    write(std_out,*) '     Usage: $fold2Bloch file_WFK x:y:z (folds)'
369    MSG_ERROR("Aborting now")
370  end if
371  read (argfolds(1:ii-1),*, iostat=ios) folds(2) !read Y folds
372  if ((ios/=0).or.(folds(2)<=0)) then
373    MSG_ERROR('Number of folds has to be a positive integer greater than 0')
374  end if
375  read(argfolds(ii+1:),*, iostat=ios) Folds(3) !read Z folds
376  if ((ios/=0).or.(folds(3)<=0)) then
377    MSG_ERROR('Number of folds has to be a positive integer greater than 0')
378  end if
379 
380  ABI_DEALLOCATE(args)
381 
382 end subroutine getargs

m_fold2block/NewK [ Functions ]

[ Top ] [ m_fold2block ] [ Functions ]

NAME

 NewK

FUNCTION

 Part of fold2Bloch that determines the unfolded
 positions of each K point.

INPUTS

 FX, FY, FZ= Number of folds in each dimention
 XX, YY, ZZ= Original K point coordinates

OUTPUT

 nkval= Array of determined coordinates of unfolded K point locaations.
 Depends on the number of folds the WFK file was structured with.

PARENTS

      fold2Bloch

CHILDREN

SOURCE

178  subroutine NewK(XX, YY, ZZ, FX, FY, FZ, NKVal)
179 
180 
181 !This section has been created automatically by the script Abilint (TD).
182 !Do not modify the following lines by hand.
183 #undef ABI_FUNC
184 #define ABI_FUNC 'NewK'
185 !End of the abilint section
186 
187  implicit none
188 
189 !Arguments ------------------------------------
190 !scalars
191  integer, intent(in) :: FX, FY, FZ
192  real(dp), intent(in) :: XX, YY, ZZ
193  real(dp), intent(inout) :: NKVal (3,FX*FY*FZ)
194 !arrays
195 
196 !Local variables-------------------------------
197 !scalars
198  real(dp) :: field_x, field_y, field_z
199  real(dp) :: TX, TY, TZ
200  integer :: loop, ii, jj, kk, size
201 !arrays
202 
203 ! *************************************************************************
204 
205  size=FX*FY*FZ
206 
207  !Brillouin zone range
208  field_x=0.5*FX
209  field_y=0.5*FY
210  field_z=0.5*FZ
211  loop=1
212 
213  !Determine new K point states
214  do ii=0, (FX-1)
215    if ((XX+ii)>(field_x)) then !if falls outside of field, loop it around
216      TX=XX-(field_x*2)+ii
217    else
218      TX=XX+ii
219    end if
220    do jj=0, (FY-1)
221      if ((YY+jj)>(field_y)) then
222        TY=YY-(field_y*2)+jj
223      else
224        TY=YY+jj
225      end if
226      do kk=0,(FZ-1)
227        if ((ZZ+kk)>(field_z)) then
228          TZ=ZZ-(field_z*2)+kk
229        else
230          TZ=ZZ+kk
231        end if
232        NKVal(1,loop)=TX !Assign to an output array
233        NKVal(2,loop)=TY
234        NKVal(3,loop)=TZ
235        loop=loop+1
236      end do
237    end do
238  end do
239 
240  !reduce the values to fit in Brillouin zone
241  do ii=1, FX*FY*FZ
242    nkval(1,ii)=nkval(1,ii)/FX
243    nkval(2,ii)=nkval(2,ii)/FY
244    nkval(3,ii)=nkval(3,ii)/FZ
245  end do
246 
247  END SUBROUTINE NewK

m_fold2block/progress [ Functions ]

[ Top ] [ m_fold2block ] [ Functions ]

NAME

 progress

FUNCTION

 Part of fold2Bloch that unfolds the wavefunction and calculates
 the weight of each band.

INPUTS

 ikpt: K point index
 nkpt: Total number of K points in the function
 kpt: Current K point being procesed.

OUTPUT

 Write to screen K point being processed and percent complete.

PARENTS

      fold2Bloch

CHILDREN

SOURCE

408  subroutine progress(ikpt, nkpt, kpt)
409 
410 
411 !This section has been created automatically by the script Abilint (TD).
412 !Do not modify the following lines by hand.
413 #undef ABI_FUNC
414 #define ABI_FUNC 'progress'
415 !End of the abilint section
416 
417  implicit none
418 
419 !Arguments ------------------------------------
420 !scalars
421  integer, intent(in) :: ikpt, nkpt
422 !arrays
423  real(dp), intent(in) :: kpt(3)
424 
425 ! *************************************************************************
426 
427  write(std_out, '(a,i10, a1)', advance='no') ''//achar(27)//'[97m',100*ikpt/nkpt,'%'//achar(27)//'[0m'
428  write(std_out, '(a25,3f12.6,a)') ''//achar(27)//'[32m Processing K point: ', kpt,''//achar(27)//'[0m'
429 
430 end subroutine progress

m_fold2block/SortC [ Functions ]

[ Top ] [ m_fold2block ] [ Functions ]

NAME

 SortC

FUNCTION

 Part of fold2Bloch that unfolds the wavefunction and calculates
 the weight of each band.

INPUTS

 FX, FY, FZ= Number of folds in each dimention
 vector= Array of energy coefficient position vectors
 coefc= Array of energy coefficients of in a certain band
 NV= Number of vectors/coefficients

OUTPUT

 wegihts= Calculated weights of a band in an unfolded state.
 Depends on the number of folds the WFK file was structured with.

PARENTS

      fold2Bloch

CHILDREN

SOURCE

 65  subroutine SortC(FX, FY, FZ, Vector, CoefC, NV, Weights)
 66 
 67 
 68 !This section has been created automatically by the script Abilint (TD).
 69 !Do not modify the following lines by hand.
 70 #undef ABI_FUNC
 71 #define ABI_FUNC 'SortC'
 72 !End of the abilint section
 73 
 74  implicit none
 75 
 76 !Arguments ------------------------------------
 77 !scalars
 78  integer, intent(in) :: FX, FY, FZ, NV
 79 !arrays
 80  integer, intent(in) :: Vector(3,NV)
 81  real(dp), intent(in) :: CoefC(2,NV)
 82  real(dp), intent(inout) ::  Weights(FX*FY*FZ)
 83 
 84 !Local variables-------------------------------
 85 !scalars
 86  real :: sumtot
 87  integer :: remainder_x, remainder_y, remainder_z, jj, kk, ll, pp, el
 88 !arrays
 89  real(dp), allocatable :: Sums(:)
 90  real(dp), allocatable :: coefsqr(:),TGroupC(:,:,:,:)
 91  integer, allocatable :: counter(:,:,:)
 92 
 93 ! *************************************************************************
 94 
 95  ABI_ALLOCATE(TGroupC,(FX, FY, FZ, NV))
 96  ABI_ALLOCATE(Sums,(FX*FY*FZ))
 97  ABI_ALLOCATE(counter,(FX,FY,FZ))
 98  ABI_ALLOCATE(coefsqr,(NV))
 99 
100  !Convert to sum of squares
101  do jj=1, NV
102    coefsqr(jj)=coefc(1,jj)*coefc(1,jj)+coefc(2,jj)*coefc(2,jj)
103  end do
104 
105  !Initiates the counter and TGroupC elements at 0
106  do jj=1,FX
107    do kk=1,FY
108      do ll=1,FZ
109        counter(jj,kk,ll)=0
110        do pp=1,NV
111          TGroupC(jj,kk,ll,pp)=0.0
112        end do
113      end do
114    end do
115  end do
116 
117  !Sorts the energy coeeficients
118  do jj=1, (NV)
119    remainder_x=MODULO(Vector(1,jj), FX)
120    remainder_y=MODULO(Vector(2,jj), FY)
121    remainder_z=MODULO(Vector(3,jj), FZ)
122    counter(remainder_x+1, remainder_y+1, remainder_z+1)=counter(remainder_x+1, remainder_y+1, remainder_z+1)+1
123    TGroupC(remainder_x+1, remainder_y+1, remainder_z+1, counter(remainder_x+1, remainder_y+1, remainder_z+1))=Coefsqr(jj)
124  end do
125 
126  !Sums the squares  of all coefficients per group
127  el=1
128  do jj=1, FX
129    do kk=1, FY
130      do ll=1, FZ
131        if (counter(jj, kk, ll)>0) then
132          Sums(el)=SUM(TGroupC(jj, kk, ll,1:counter(jj, kk, ll)))
133          el=el+1
134        else
135          Sums(el)=0.0
136          el=el+1
137        end if
138      end do
139    end do
140  end do
141 
142  !Assign weights to an array
143  sumtot=SUM(Sums)
144  do jj=1, (FX*FY*FZ)
145    Weights(jj)=Sums(jj)/sumtot
146  end do
147  ABI_DEALLOCATE(TGroupC)
148  ABI_DEALLOCATE(Sums)
149  ABI_DEALLOCATE(counter)
150  ABI_DEALLOCATE(coefsqr)
151 
152  end subroutine SortC