TABLE OF CONTENTS
ABINIT/wvl_wfsinp_scratch [ Functions ]
NAME
wvl_wfsinp_scratch
FUNCTION
This method allocates and initialises wavefunctions with values from input guess. See wvl_wfsinp_disk() 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 ireadwf=1 for reading from file, 0 otherwise. 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
wvl <type(wvl_data)>=wavefunctions & projectors informations for wavelets.
PARENTS
inwffil
CHILDREN
input_wf_diag,mklocl_wavelets,wrtout,wvl_occ_abi2big,wvl_occopt_abi2big xred2xcart
SOURCE
46 #if defined HAVE_CONFIG_H 47 #include "config.h" 48 #endif 49 50 #include "abi_common.h" 51 52 53 subroutine wvl_wfsinp_scratch(dtset, mpi_enreg, occ, rprimd, wvl, xred) 54 55 use defs_basis 56 use defs_abitypes 57 use defs_wvltypes 58 use m_wffile 59 use m_profiling_abi 60 use m_errors 61 use m_ab7_kpoints 62 use m_ab7_symmetry 63 use m_xmpi 64 65 use m_geometry, only : xred2xcart 66 use m_abi2big, only : wvl_occ_abi2big,wvl_occopt_abi2big 67 68 #if defined HAVE_BIGDFT 69 use BigDFT_API, only : createIonicPotential, input_wf_diag, gaussian_basis, & 70 & input_variables, calculate_rhocore, deallocate_Lzd_except_Glr, INPUT_IG_OFF,& 71 & SMEARING_DIST_ERF, PSPCODE_PAW 72 #endif 73 74 !This section has been created automatically by the script Abilint (TD). 75 !Do not modify the following lines by hand. 76 #undef ABI_FUNC 77 #define ABI_FUNC 'wvl_wfsinp_scratch' 78 use interfaces_14_hidewrite 79 use interfaces_67_common 80 !End of the abilint section 81 82 implicit none 83 84 !Arguments ------------------------------- 85 !scalars 86 type(dataset_type), intent(in) :: dtset 87 type(MPI_type), intent(inout) :: mpi_enreg 88 type(wvl_data), intent(inout) :: wvl 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 character(len = 500) :: message 97 integer :: comm,me,nproc 98 integer :: iscf_local 99 integer :: nvirt 100 integer :: ii,shift_vpsp,size_vpsp 101 logical :: onlywf=.false. ! find the wavefunctions and return 102 logical :: wvlbigdft=.false. 103 real(dp), allocatable :: xcart(:,:) 104 real(dp), allocatable :: rhor(:,:) 105 real(dp), pointer :: vpsp(:) 106 real(dp):: elecfield(3) 107 type(gaussian_basis) :: Gvirt 108 type(input_variables) :: in ! To be removed, waiting for BigDFT upgrade 109 #endif 110 111 ! ********************************************************************* 112 113 #if defined HAVE_BIGDFT 114 115 elecfield=zero 116 117 !If usewvl: wvlbigdft indicates that the BigDFT workflow will be followed 118 wvlbigdft=(dtset%usewvl==1.and.dtset%wvl_bigdft_comp==1) 119 120 write(message, '(a,a)' ) ch10,& 121 & ' wvl_wfsinp_scratch: wavefunction initialisation.' 122 call wrtout(std_out,message,'COLL') 123 124 comm=mpi_enreg%comm_wvl 125 me=xmpi_comm_rank(comm) 126 nproc=xmpi_comm_size(comm) 127 !Store xcart for each atom 128 ABI_ALLOCATE(xcart,(3, dtset%natom)) 129 call xred2xcart(dtset%natom, rprimd, xcart, xred) 130 131 !We allocate temporary arrays for rho and vpsp. 132 !allocate ionic potential 133 if (wvl%den%denspot%dpbox%n3pi > 0) then 134 size_vpsp=wvl%descr%Glr%d%n1i*wvl%descr%Glr%d%n2i*wvl%den%denspot%dpbox%n3pi 135 shift_vpsp=wvl%den%denspot%dpbox%ndims(1)*wvl%den%denspot%dpbox%ndims(2) & 136 & *wvl%den%denspot%dpbox%nscatterarr(me,4) 137 ABI_ALLOCATE(vpsp,(size_vpsp+shift_vpsp)) 138 else 139 ABI_ALLOCATE(vpsp,(1)) 140 end if 141 142 if(.not. onlywf) then 143 ABI_ALLOCATE(rhor,(dtset%nfft,dtset%nspden)) 144 call mklocl_wavelets(elecfield, xcart, mpi_enreg, dtset%natom, & 145 & dtset%nfft, dtset%nspden, 1, rprimd, vpsp, & 146 & wvl%den, wvl%descr, xcart) 147 ABI_DEALLOCATE(rhor) 148 end if 149 150 ! IMPORTANT: onlywf=.true. does not work yet, do not change this: 151 ! if(.not. wvlbigdft) onlywf=.true. !do not apply Hamiltonian inside input_wf_diag. 152 if(dtset%usepaw==1) wvl%descr%atoms%npspcode(:)=wvl%descr%npspcode_paw_init_guess(:) 153 iscf_local=dtset%iscf 154 if(.not. wvlbigdft) iscf_local=0 !important to have good occ values 155 156 !This routine allocates psi, hpsi and psit. 157 nvirt = 0 158 in%linear = INPUT_IG_OFF 159 in%nspin = dtset%nsppol 160 in%exctxpar = wvl%descr%exctxpar 161 in%itrpmax = dtset%nnsclo 162 in%iscf = iscf_local !dtset%iscf 163 in%Tel = dtset%tsmear 164 ! if (dtset%iscf == 0) in%Tel = zero 165 if (iscf_local == 0) in%Tel = zero 166 in%SIC = wvl%wfs%ks%SIC 167 in%orthpar = wvl%wfs%ks%orthpar 168 169 in%occopt=dtset%occopt 170 if(dtset%occopt>2 .and. dtset%occopt<7)then 171 call wvl_occopt_abi2big(in%occopt,in%occopt,1) 172 else 173 ! This will be used only for the initial wavefunctions: 174 in%occopt=SMEARING_DIST_ERF 175 end if 176 177 !Note: check if all required in "in" is passed. 178 !remove argument "in" from input_wf_diag 179 call input_wf_diag(me, nproc, & 180 & wvl%descr%atoms, wvl%den%denspot,& 181 & wvl%wfs%ks%orbs, nvirt, wvl%wfs%ks%comms, & 182 & wvl%wfs%ks%lzd, wvl%e%energs, xcart, & 183 & wvl%projectors%nlpsp, & 184 & dtset%ixc, wvl%wfs%ks%psi, wvl%wfs%ks%hpsi, wvl%wfs%ks%psit, & 185 & Gvirt, in%nspin, wvl%wfs%GPU, in, onlywf, wvl%projectors%G, wvl%descr%paw) 186 187 !This provisory: wvl%descr%paw could be passed as optional 188 !to input_wf_diag to allocate spsi inside this routine 189 if(dtset%usepaw==1) then 190 wvl%descr%atoms%npspcode(:)=PSPCODE_PAW 191 ABI_ALLOCATE(wvl%descr%paw%spsi,(max(wvl%wfs%ks%orbs%npsidim_orbs,wvl%wfs%ks%orbs%npsidim_comp))) 192 do ii=1,size(wvl%wfs%ks%psi) 193 wvl%descr%paw%spsi(ii)=wvl%wfs%ks%psi(ii) 194 wvl%wfs%ks%hpsi(ii)=wvl%wfs%ks%psi(ii) 195 end do 196 end if 197 198 if(wvlbigdft ) then 199 ! Copy occupations from BigDFT objects to ABINIT 200 call wvl_occ_abi2big(dtset%mband,dtset%nkpt,dtset%nsppol,occ,2,wvl%wfs) 201 else 202 ! Copy occupations from ABINIT to BigDFT objects 203 call wvl_occ_abi2big(dtset%mband,dtset%nkpt,dtset%nsppol,occ,1,wvl%wfs) 204 end if 205 206 write(message, '(a)' ) & 207 & ' | wavefunctions have been calculated.' 208 call wrtout(std_out,message,'COLL') 209 210 ABI_DEALLOCATE(xcart) 211 ABI_DEALLOCATE(vpsp) 212 213 #else 214 BIGDFT_NOTENABLED_ERROR() 215 if (.false.) write(std_out,*) dtset%nstep,mpi_enreg%me,wvl%wfs%ks,occ(1),rprimd(1,1),xred(1,1) 216 #endif 217 218 end subroutine wvl_wfsinp_scratch