TABLE OF CONTENTS
ABINIT/wvl_wfsinp_reformat [ Functions ]
NAME
wvl_wfsinp_reformat
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%ks%[h]psi comes unallocated and will be allocated inside this routine. When initialised from memory (reformating), wvl%wfs%ks%[h]psi will be reallocated. The projectors are also recomputed. The scalar arrays should be reallocated using dtset%nfft after a call to this routine.
COPYRIGHT
Copyright (C) 2005-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 .
INPUTS
OUTPUT
SIDE EFFECTS
NOTES
PARENTS
mover
CHILDREN
copy_old_wavefunctions,deallocate_wfd,first_orthon local_potential_dimensions,nullify_gaussian_basis,psolver_kernel reformatmywaves,wrtout,wvl_denspot_free,wvl_denspot_set wvl_projectors_free,wvl_projectors_set,wvl_setboxgeometry,wvl_setngfft wvl_wfs_free,wvl_wfs_lr_copy,wvl_wfs_set,xred2xcart
SOURCE
45 #if defined HAVE_CONFIG_H 46 #include "config.h" 47 #endif 48 49 #include "abi_common.h" 50 51 52 subroutine wvl_wfsinp_reformat(dtset, mpi_enreg, psps, rprimd, wvl, xred, xred_old) 53 54 use defs_basis 55 use defs_datatypes 56 use defs_abitypes 57 use defs_wvltypes 58 use m_profiling_abi 59 use m_errors 60 use m_xmpi 61 62 use m_geometry, only : xred2xcart 63 64 #if defined HAVE_BIGDFT 65 use BigDFT_API, only : copy_old_wavefunctions, reformatmywaves, first_orthon, & 66 & deallocate_wfd, wavefunctions_descriptors, deallocate_lr, & 67 & local_potential_dimensions, copy_coulomb_operator, & 68 & deallocate_coulomb_operator, nullify_gaussian_basis 69 use dynamic_memory 70 #endif 71 72 !This section has been created automatically by the script Abilint (TD). 73 !Do not modify the following lines by hand. 74 #undef ABI_FUNC 75 #define ABI_FUNC 'wvl_wfsinp_reformat' 76 use interfaces_14_hidewrite 77 use interfaces_43_wvl_wrappers 78 use interfaces_62_poisson 79 !End of the abilint section 80 81 implicit none 82 83 !Arguments ------------------------------------ 84 type(dataset_type), intent(inout) :: dtset 85 type(MPI_type), intent(inout) :: mpi_enreg 86 type(pseudopotential_type), intent(in) :: psps 87 type(wvl_data), intent(inout) :: wvl 88 real(dp), intent(inout) :: rprimd(3,3) 89 real(dp), intent(inout) :: xred_old(3, dtset%natom) 90 real(dp), intent(inout) :: xred(3, dtset%natom) 91 92 !Local variables------------------------------- 93 #if defined HAVE_BIGDFT 94 integer :: itypat 95 integer :: nSize_old(3) 96 real(dp) :: hgrid_old(3) 97 real(dp), allocatable :: xcart(:,:), xcart_old(:,:) 98 real(dp), pointer :: psi_old(:), eigen_old(:) 99 integer :: comm,me,nproc,icoulomb 100 type(coulomb_operator)::kernel 101 type(wavefunctions_descriptors) :: keys_old 102 character(len=500) :: message 103 #if defined DEBUG_MODE 104 !integer, parameter :: ndebug = 5 !5 will not work for wavelets compiling with debug=naughty 105 integer,parameter :: ndebug = 0 106 #else 107 integer, parameter :: ndebug = 0 108 #endif 109 #endif 110 111 ! ********************************************************************* 112 113 #if defined HAVE_BIGDFT 114 115 write(message, '(a,a)' ) ch10,& 116 & ' wvl_wfsinp_reformat: reformat the wavefunctions.' 117 call wrtout(std_out, message, 'COLL') 118 119 comm=mpi_enreg%comm_wvl 120 me=xmpi_comm_rank(comm) 121 nproc=xmpi_comm_size(comm) 122 !Convert input xred_old (reduced coordinates) to xcart_old (cartesian) 123 ABI_ALLOCATE(xcart_old,(3, dtset%natom)) 124 call xred2xcart(dtset%natom, rprimd, xcart_old, xred_old) 125 126 !Copy current to old. 127 ABI_ALLOCATE(eigen_old,(wvl%wfs%ks%orbs%norb)) 128 eigen_old = wvl%wfs%ks%orbs%eval 129 hgrid_old = wvl%descr%h 130 call copy_old_wavefunctions(nproc, wvl%wfs%ks%orbs, & 131 & wvl%descr%Glr%d%n1, wvl%descr%Glr%d%n2, wvl%descr%Glr%d%n3, & 132 & wvl%wfs%ks%lzd%Glr%wfd, wvl%wfs%ks%psi, nSize_old(1), nSize_old(2), nSize_old(3), & 133 & keys_old, psi_old) 134 !Patch because copy_old_wavefunctions() free wvl%wfs%ks%lzd%Glr%wfd but don't nullify it. 135 nullify(wvl%wfs%ks%lzd%glr%wfd%keyglob) 136 nullify(wvl%wfs%ks%lzd%glr%wfd%keygloc) 137 nullify(wvl%wfs%ks%lzd%glr%wfd%keyvloc) 138 nullify(wvl%wfs%ks%lzd%glr%wfd%keyvglob) 139 140 !We deallocate the previous projectors. 141 call wvl_projectors_free(wvl%projectors) 142 143 !Deallocate old wavefunctions 144 call wvl_wfs_free(wvl%wfs) 145 146 !Deallocate old denspot 147 call wvl_denspot_free(wvl%den) 148 149 !We change the box geometry. 150 call wvl_setBoxGeometry(dtset%prtvol, psps%gth_params%radii_cf, rprimd, xred, & 151 & wvl%descr, dtset%wvl_crmult, dtset%wvl_frmult) 152 call wvl_denspot_set(wvl%den, psps%gth_params, dtset%ixc, dtset%natom, dtset%nsppol, & 153 & rprimd, wvl%descr, dtset%wvl_crmult, dtset%wvl_frmult, comm, xred) 154 if (wvl%descr%atoms%astruct%geocode == "F") then 155 icoulomb = 1 156 else if (wvl%descr%atoms%astruct%geocode == "S") then 157 icoulomb = 2 158 else 159 icoulomb = 0 160 end if 161 !calculation of the Poisson kernel anticipated to reduce memory peak for small systems 162 call psolver_kernel( wvl%den%denspot%dpbox%hgrids, 2, icoulomb, me, kernel, & 163 & comm, wvl%den%denspot%dpbox%ndims, nproc, dtset%nscforder) 164 ! Shallow copy of the kernel (still owned by ABINIT). 165 wvl%den%denspot%pkernel = kernel 166 wvl%den%denspot%pkernelseq = kernel 167 !Associate the denspot distribution into mpi_enreg. 168 mpi_enreg%nscatterarr => wvl%den%denspot%dpbox%nscatterarr 169 mpi_enreg%ngatherarr => wvl%den%denspot%dpbox%ngatherarr 170 mpi_enreg%ngfft3_ionic = wvl%den%denspot%dpbox%n3pi 171 call wvl_setngfft(me, dtset%mgfft, dtset%nfft, & 172 & dtset%ngfft, nproc, wvl%den%denspot%dpbox%ndims(1), & 173 & wvl%den%denspot%dpbox%ndims(2), & 174 & wvl%den%denspot%dpbox%ndims(3),wvl%den%denspot%dpbox%n3d) 175 176 !We copy the geometry structure. 177 call wvl_wfs_lr_copy(wvl%wfs, wvl%descr) 178 !Reallocate them with new size. 179 call wvl_wfs_set(dtset%strprecon,dtset%spinmagntarget, dtset%kpt, me, dtset%natom, sum(dtset%nband), dtset%nkpt, & 180 & nproc, dtset%nspinor, dtset%nsppol, dtset%nwfshist, dtset%occ_orig, psps, rprimd, & 181 & wvl%wfs, dtset%wtk, wvl%descr, dtset%wvl_crmult, dtset%wvl_frmult, xred) 182 183 !Recopy old eval for precond. 184 wvl%wfs%ks%orbs%eval = eigen_old 185 ABI_DEALLOCATE(eigen_old) 186 187 !We allocate psi. 188 !ABI_ALLOCATE(wvl%wfs%ks%psi,( max(wvl%wfs%ks%orbs%npsidim_comp,wvl%wfs%ks%orbs%npsidim_orbs)+ndebug) ) 189 wvl%wfs%ks%psi=f_malloc_ptr(max(wvl%wfs%ks%orbs%npsidim_comp,wvl%wfs%ks%orbs%npsidim_orbs)+ndebug,id='psi') 190 write(message, '(a,a,a,a,I0)' ) ch10, & 191 & ' wvl_wfsinp_reformat: allocate wavefunctions,', ch10, & 192 & ' size of the compressed array per proc: ', & 193 & product(shape(wvl%wfs%ks%psi)) 194 call wrtout(std_out,message,'COLL') 195 196 !Convert input xred (reduced coordinates) to xcart (cartesian) 197 ABI_ALLOCATE(xcart,(3, dtset%natom)) 198 call xred2xcart(dtset%natom, rprimd, xcart, xred) 199 200 !We transfer the old wavefunctions to the new ones. 201 call reformatmywaves(me, wvl%wfs%ks%orbs, wvl%descr%atoms, & 202 & hgrid_old(1), hgrid_old(2), hgrid_old(3), nSize_old(1), nSize_old(2), & 203 & nSize_old(3), xcart_old, keys_old, psi_old, wvl%descr%h(1), wvl%descr%h(2), & 204 & wvl%descr%h(3), wvl%descr%Glr%d%n1, wvl%descr%Glr%d%n2, wvl%descr%Glr%d%n3, xcart, & 205 & wvl%wfs%ks%lzd%Glr%wfd, wvl%wfs%ks%psi) 206 ABI_DEALLOCATE(xcart) 207 ABI_DEALLOCATE(xcart_old) 208 209 !We free the old descriptors and arrays. 210 ABI_DEALLOCATE(psi_old) 211 call deallocate_wfd(keys_old) 212 213 call local_potential_dimensions(me,wvl%wfs%ks%lzd,wvl%wfs%ks%orbs,wvl%den%denspot%xc,& 214 & wvl%den%denspot%dpbox%ngatherarr(0,1)) 215 216 !it seems that the table "wvl%projectors%G" is no more used 217 !but it's not allocated -> fortran runtime error 218 #if defined HAVE_BIGDFT 219 ABI_DATATYPE_ALLOCATE(wvl%projectors%G,(dtset%ntypat)) 220 do itypat=1,dtset%ntypat 221 call nullify_gaussian_basis(wvl%projectors%G(itypat)) 222 end do 223 #endif 224 225 !Reallocate projectors for the new positions. 226 call wvl_projectors_set(me, dtset%natom, wvl%projectors, psps, rprimd, & 227 & wvl%wfs, wvl%descr, dtset%wvl_frmult, xred) 228 229 !Orthogonilise new wavefunctions. 230 call first_orthon(me, nproc, wvl%wfs%ks%orbs, wvl%wfs%ks%lzd, wvl%wfs%ks%comms, & 231 & wvl%wfs%ks%psi, wvl%wfs%ks%hpsi, wvl%wfs%ks%psit, wvl%wfs%ks%orthpar,wvl%descr%paw) 232 233 #else 234 BIGDFT_NOTENABLED_ERROR() 235 if (.false.) write(std_out,*) dtset%nstep,mpi_enreg%me,psps%npsp,wvl%wfs%ks,rprimd(1,1),& 236 & xred_old(1,1),xred(1,1) 237 #endif 238 239 end subroutine wvl_wfsinp_reformat