TABLE OF CONTENTS


ABINIT/wvl_wfsinp_disk [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_wfsinp_disk

FUNCTION

 This method allocates and initialises wavefunctions with values from disk.
 See wvl_wfsinp_scratch() or wvl_wfsinp_reformat() from other initialisation
 routines.
 
 When initialised from scratch or from disk, wvl%wfs%[h]psi comes unallocated
 and will be allocated inside this routine.
 When initialised from memory (reformating), wvl%wfs%[h]psi will be reallocated.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DC)
 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

  dtset <type(dataset_type)>=input variables.
  hdr0 <type(hdr_type)>=the header of wf, den and pot files (read from restart)
  hdr <type(hdr_type)>=the header of wf, den and pot files
  mpi_enreg=informations about MPI parallelization
  option=1 for reading a file following ABINIT format, -1 for a BigDFT format.
  rprimd(3,3)=dimensional primitive translations in real space (bohr)
  wff <type(wffile_type)>= structure with informations on wf file.
  xred(3,natom)=reduced dimensionless atomic coordinates.

OUTPUT

SIDE EFFECTS

  wfs <type(wvl_projector_type)>=wavefunctions informations for wavelets.

PARENTS

      inwffil

CHILDREN

      first_orthon,wrtout,wvl_occ_abi2big,wvl_read

SOURCE

 44 #if defined HAVE_CONFIG_H
 45 #include "config.h"
 46 #endif
 47 
 48 #include "abi_common.h"
 49 
 50 
 51 subroutine wvl_wfsinp_disk(dtset, hdr0, hdr, mpi_enreg, occ, option, rprimd, wff, wfs, wvl, xred)
 52 
 53  use defs_basis
 54  use defs_abitypes
 55  use defs_wvltypes
 56  use m_wffile
 57  use m_profiling_abi
 58  use m_errors
 59  use m_xmpi
 60 
 61  use m_abi2big,          only : wvl_occ_abi2big
 62 
 63 #if defined HAVE_BIGDFT
 64  use BigDFT_API, only : first_orthon,sumrho,communicate_density,plot_density
 65  use dynamic_memory
 66 #endif
 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 'wvl_wfsinp_disk'
 72  use interfaces_14_hidewrite
 73  use interfaces_62_wvl_wfs
 74 !End of the abilint section
 75 
 76   implicit none
 77 
 78 !Arguments -------------------------------
 79   !scalars
 80   integer, intent(in)                       :: option
 81   type(dataset_type), intent(in)            :: dtset
 82   type(hdr_type), intent(in)                :: hdr0
 83   type(hdr_type), intent(in)                :: hdr
 84   type(MPI_type), intent(in)                :: mpi_enreg
 85   type(wffile_type), intent(in)             :: wff
 86   type(wvl_wf_type), intent(inout)          :: wfs
 87   type(wvl_internal_type), intent(inout)       :: wvl
 88   !type(wvl_denspot_type), intent(inout)       :: wvl_den
 89   !arrays
 90   real(dp), intent(inout) :: occ(dtset%mband*dtset%nkpt*dtset%nsppol)
 91   real(dp), intent(in)                      :: rprimd(3, 3)
 92   real(dp), intent(in)                      :: xred(3, dtset%natom)
 93 
 94 !Local variables-------------------------------
 95 #if defined HAVE_BIGDFT
 96  logical :: wvlbigdft=.false.
 97  integer :: comm,me,nproc
 98  character(len = 500)  :: message
 99  !real(dp), allocatable :: xcart(:,:)
100 #if defined DEBUG_MODE
101  !integer, parameter :: ndebug = 5  !5 will not work for wavelets compiling with debug=naughty
102  integer,parameter :: ndebug = 0
103 #else
104  integer, parameter :: ndebug = 0
105 #endif
106 #endif
107 
108 ! *********************************************************************
109 
110 #if defined HAVE_BIGDFT
111 
112  write(message, '(a,a)' ) ch10,' wvl_wfsinp_disk: wavefunction initialisation.'
113  call wrtout(std_out,message,'COLL')
114 
115 !If usewvl: wvlbigdft indicates that the BigDFT workflow will be followed
116  wvlbigdft=(dtset%usewvl==1.and.dtset%wvl_bigdft_comp==1)
117 
118  comm=mpi_enreg%comm_wvl
119  me=xmpi_comm_rank(comm)
120  nproc=xmpi_comm_size(comm)
121 !We allocate psi.
122 !ABI_ALLOCATE(wfs%ks%psi,( max(wfs%ks%orbs%npsidim_comp,wfs%ks%orbs%npsidim_orbs)+ndebug) )
123  wfs%ks%psi=f_malloc_ptr(max(wfs%ks%orbs%npsidim_comp,wfs%ks%orbs%npsidim_orbs)+ndebug,id='psi')
124 
125  write(message, '(a,a,a,a,I0)' ) ch10, &
126 & ' wvl_wfsinp_disk: allocate wavefunctions,', ch10, &
127 & '  size of the compressed array per proc: ', &
128 & product(shape(wfs%ks%psi))
129  call wrtout(std_out,message,'COLL')
130 
131  call wvl_read(dtset, hdr0, hdr, mpi_enreg, option, rprimd, wff, wfs, wvl, xred)
132 
133 !We orthogonalise,only for NC.
134  if(wvl%paw%usepaw==0 .and. wvlbigdft) then
135    call first_orthon(me, nproc, wfs%ks%orbs, wfs%ks%lzd, wfs%ks%comms, &
136 &   wfs%ks%psi, wfs%ks%hpsi, wfs%ks%psit, wfs%ks%orthpar,wvl%paw)
137  else
138 !  ABI_ALLOCATE(wfs%ks%hpsi,(max(wfs%ks%orbs%npsidim_orbs,wfs%ks%orbs%npsidim_comp)))
139    wfs%ks%hpsi=f_malloc_ptr(max(wfs%ks%orbs%npsidim_orbs,wfs%ks%orbs%npsidim_comp),id='hpsi')
140    if(wvl%paw%usepaw==1) then
141      ABI_ALLOCATE(wvl%paw%spsi,(max(wfs%ks%orbs%npsidim_orbs,wfs%ks%orbs%npsidim_comp)))
142    end if
143 
144 !  Set orbs%eval=-0.5.
145 !  This will be done in LDiagHam
146 !  For the moment we skip this, since hpsi is not yet calculated
147 !  and it an input argument in LDiagHam.
148    wfs%ks%orbs%eval(:)=-0.5d0
149 
150 !  Copy occupations from BigDFT objects to ABINIT
151    call wvl_occ_abi2big(dtset%mband,dtset%nkpt,dtset%nsppol,occ,2,wfs)
152 
153 
154  end if
155 
156 #else
157  BIGDFT_NOTENABLED_ERROR()
158  if (.false.) write(std_out,*) option,dtset%nstep,hdr0%ecut,hdr%ecut,mpi_enreg%me,wff%me,&
159 & wfs%ks,wvl%h(1),occ(1),rprimd(1,1),xred(1,1)
160 #endif
161 
162 !for testing
163 !!Plot the density
164 !call sumrho(wvl_den%denspot%dpbox,wfs%ks%orbs,&
165 !& wfs%GPU,& wvl%atoms%sym,&
166 !& wvl_den%denspot%rhod,wfs%psi,wvl_den%denspot%rho_psi)
167 !call communicate_density(wvl_den%denspot%dpbox,wfs%ks%orbs%nspin,&
168 !& wvl_den%denspot%rhod,wvl_den%denspot%dpcom%nscatterarr,&
169 !& wvl_den%denspot%rho_psi,wvl_den%denspot%rhov)
170 !call plot_density('electronic_density',&
171 !& me,nproc,wfs%Lzd%Glr%d%n1,wfs%Lzd%Glr%d%n2,wfs%Lzd%Glr%d%n3,&
172 !& wfs%Lzd%Glr%d%n1i,wfs%Lzd%Glr%d%n2i,wfs%Lzd%Glr%d%n3i,&
173 !& wvl_den%denspot%dpcom%nscatterarr(me,2),  &
174 !& wfs%orbs%nspin,&
175 !& wvl_den%denspot%hgrids(1),wvl_den%denspot%hgrids(2),wvl_den%denspot%hgrids(3),&
176 !& wvl%atoms,xcart,wvl_den%denspot%dpcom%ngatherarr,&
177 !& wvl_den%denspot%rhov(1+wvl_den%denspot%dpcom%nscatterarr(me,4)*wfs%Lzd%Glr%d%n1i*wfs%Lzd%Glr%d%n2i))
178 !ABI_DEALLOCATE(xcart)
179 !end of debug
180 
181 
182 end subroutine wvl_wfsinp_disk