TABLE OF CONTENTS


ABINIT/initwf [ Functions ]

[ Top ] [ Functions ]

NAME

 initwf

FUNCTION

 Initialization of wavefunctions.
 If formeig==1, and partially filled case, I am not sure that the eig_k are initialized properly ...
 formeig option (format of the eigenvalues and eigenvector) :
   0 => ground-state format (initialisation of
        eigenvectors with random numbers, vector of eigenvalues)
   1 => respfn format (initialisation of
        eigenvectors with 0 s, hermitian matrix of eigenvalues)

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR)
 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

 formeig=see above
 headform=header format (might be needed to read the block of wfs)
 icg=shift to be given to the location of the data in the array cg
 ikpt= number of the k point of which the wf is initialised
 spin=spin index
 mcg=dimension of the cg array
 mpi_enreg=information about MPI parallelization
 nband_k=number of bands at this particular k point
 nkpt=number of k points
 npw=number of plane waves
 nspinor=number of spinorial components of the wavefunctions (on current proc)
 wff1=structure info for file containing wavefunctions (when needed)

OUTPUT

 cg(2,mcg)=complex wf array
  if ground state format (formeig=0):
    eig_k(nband_k)=list of eigenvalues (input or init to large number), hartree
  if respfn format (formeig=1):
    eig_k(2*nband_k*nband_k)= matrix of eigenvalues (input or init to large number), hartree

SIDE EFFECTS

 Input/output:
 occ_k(nband_k)=list of occupations (input or left to their initial value)
 ikptsp_old=number of the previous spin-k point, or 0 if first call of present file

PARENTS

      wfsinp

CHILDREN

      rwwf,timab,wffreadskipk,wfk_close,wfk_open_read,wfk_read_band_block
      wrtout

SOURCE

 58 #if defined HAVE_CONFIG_H
 59 #include "config.h"
 60 #endif
 61 
 62 #include "abi_common.h"
 63 
 64 subroutine initwf(cg,eig_k,formeig,headform,icg,ikpt,ikptsp_old,&
 65 &  spin,mcg,mpi_enreg,nband_k,nkpt,npw,nspinor,occ_k,wff1)
 66 
 67  use defs_basis
 68  use defs_abitypes
 69  use m_profiling_abi
 70  use m_errors
 71  use m_xmpi
 72  use m_wfk
 73  use m_wffile
 74 
 75  use m_io_tools,     only : get_unit
 76 
 77 !This section has been created automatically by the script Abilint (TD).
 78 !Do not modify the following lines by hand.
 79 #undef ABI_FUNC
 80 #define ABI_FUNC 'initwf'
 81  use interfaces_14_hidewrite
 82  use interfaces_18_timing
 83  use interfaces_56_io_mpi
 84  use interfaces_62_iowfdenpot, except_this_one => initwf
 85 !End of the abilint section
 86 
 87  implicit none
 88 
 89 !Arguments ------------------------------------
 90 !scalars
 91  integer,intent(in) :: formeig,headform,icg,ikpt,spin,mcg,nband_k,nkpt,npw,nspinor
 92  integer,intent(inout) :: ikptsp_old
 93  type(MPI_type),intent(in) :: mpi_enreg
 94  type(wffile_type),intent(inout) :: wff1
 95 !arrays
 96  real(dp),intent(inout) :: occ_k(nband_k)
 97  real(dp),intent(inout) :: cg(2,mcg),eig_k((2*nband_k)**formeig*nband_k) !vz_i
 98 
 99 !Local variables-------------------------------
100 !scalars
101  integer,parameter :: nkpt_max=50
102  integer :: ikpt0,nband_disk,tim_rwwf 
103  character(len=500) :: msg
104 !arrays
105  integer,allocatable :: kg_dum(:,:)
106  real(dp) :: tsec(2)
107 #if 0
108  integer :: iomode,comm,funt,ierr
109  type(wfk_t) :: Wfk
110 #endif
111 
112 ! *************************************************************************
113 
114 !DEBUG
115 !write(std_out,*)' initwf : enter, ikptsp_old,ikpt,spin,nkpt= ',ikptsp_old,ikpt,spin,nkpt
116 !stop
117 !ENDDEBUG
118 
119 #if 0
120  MSG_WARNING("Entering new IO section")
121 !call WffClose(wff1,ierr)
122  comm   = MPI_enreg%comm_cell
123  iomode = iomode_from_fname(wff1%fname)
124  call wfk_open_read(Wfk,wff1%fname,formeig,iomode,get_unit(),comm)
125  call wfk_read_band_block(Wfk,(/1,nband_k/),ikpt,spin,xmpio_at,cg_k=cg(1:,icg:),eig_k=eig_k,occ_k=occ_k)
126  call wfk_close(Wfk)
127 !call clsopn(wff1)
128  RETURN
129 #endif
130 
131  call timab(770,1,tsec)
132  call timab(771,1,tsec)
133 
134  ABI_ALLOCATE(kg_dum,(3,0))
135 
136 !Skip wavefunctions for k-points not treated by this proc.
137 !(from ikptsp_old+1 to ikpt+(spin-1)*nkpt-1)
138  if (ikptsp_old<ikpt+(spin-1)*nkpt-1) then
139    do ikpt0=ikptsp_old+1,ikpt+(spin-1)*nkpt-1
140      call WffReadSkipK(formeig,headform,ikpt0,spin,mpi_enreg,wff1)
141    end do
142  end if
143 
144 !DEBUG
145 !write(std_out,*)' initwf : before rwwf'
146 !write(std_out,*)' formeig,icg,ikpt,spin=',formeig,icg,ikpt,spin
147 !write(std_out,*)' nband_k,nband_disk,npw,nspinor=',nband_k,nband_disk,npw,nspinor
148 !write(std_out,*)' unwff1=',unwff1
149 !stop
150 !ENDDEBUG
151 
152  if(mpi_enreg%paralbd==0)tim_rwwf=2
153  if(mpi_enreg%paralbd==1)tim_rwwf=20
154 
155  call timab(771,2,tsec)
156 
157  call rwwf(cg,eig_k,formeig,headform,icg,ikpt,spin,kg_dum,nband_k,mcg,mpi_enreg,nband_k,nband_disk,&
158 & npw,nspinor,occ_k,1,0,tim_rwwf,wff1)
159 
160  call timab(772,1,tsec)
161 
162  if (ikpt<=nkpt_max) then
163    write(msg,'(3(a,i0))')' initwf: disk file gives npw= ',npw,' nband= ',nband_disk,' for kpt number= ',ikpt
164    call wrtout(std_out,msg,'PERS')
165  else if (ikpt==nkpt_max+1) then
166    call wrtout(std_out,' initwf: the number of similar message is sufficient... stop printing them','PERS')
167  end if
168 
169 !Check the number of bands on disk file against desired number. These are not required to agree)
170  if (nband_disk/=nband_k) then
171    write(msg,'(2(a,i0),3a,i0,3a)')&
172 &   'For kpt number ',ikpt,' disk file has ',nband_disk,' bands',ch10,&
173 &   'but input file gave nband= ',nband_k,'.',ch10,&
174 &   'This is not fatal. Bands are skipped or filled with random numbers.'
175    MSG_COMMENT(msg)
176  end if
177 
178  if (ikpt<=nkpt_max) then
179    write(msg,'(a,i0,a)')' initwf: ',nband_disk,' bands have been initialized from disk'
180    call wrtout(std_out,msg,'PERS')
181  end if
182 
183  ikptsp_old=ikpt+(spin-1)*nkpt
184 
185  ABI_DEALLOCATE(kg_dum)
186 
187  call timab(772,2,tsec)
188  call timab(770,2,tsec)
189 
190 end subroutine initwf