TABLE OF CONTENTS
ABINIT/wvl_hpsitopsi [ Functions ]
NAME
wvl_hpsitopsi
FUNCTION
Heart of the wavelet resolution, compute new wavefunctions mixed witf previous by computing the gradient of the wavefunctions knowing the external potential.
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. istep=id of the current iteration (first is 1). mpi_enreg=informations about MPI parallelization proj <type(wvl_projector_type)>=projectors informations for wavelets. vtrial(dtset%nfft)=external potential. xcart(3,natom)=cartesian atomic coordinates
OUTPUT
SIDE EFFECTS
energies <type(energies_type)>=storage for energies computed here : | e_kinetic(OUT)=kinetic energy part of total energy | e_localpsp(OUT)=local pseudopotential part of total energy | e_nonlocalpsp(OUT)=nonlocal pseudopotential part of total energy residm=max value for gradient in the minimisation process. rhor(dtset%nfft)=electron density in r space wfs <type(wvl_projector_type)>=wavefunctions informations for wavelets.
PARENTS
vtorho
CHILDREN
calculate_energy_and_gradient,hpsitopsi,pawcprj_alloc,wrtout
SOURCE
43 #if defined HAVE_CONFIG_H 44 #include "config.h" 45 #endif 46 47 #include "abi_common.h" 48 49 50 subroutine wvl_hpsitopsi(cprj,dtset,energies,istep,mcprj,mpi_enreg,residm,wvl,xcart) 51 52 use m_profiling_abi 53 54 use defs_basis 55 use defs_abitypes 56 use defs_wvltypes 57 use m_errors 58 use m_xmpi 59 60 use m_energies, only : energies_type 61 use m_pawcprj, only : pawcprj_type, pawcprj_alloc 62 #if defined HAVE_BIGDFT 63 use BigDFT_API, only : hpsitopsi, calculate_energy_and_gradient 64 #endif 65 66 !This section has been created automatically by the script Abilint (TD). 67 !Do not modify the following lines by hand. 68 #undef ABI_FUNC 69 #define ABI_FUNC 'wvl_hpsitopsi' 70 use interfaces_14_hidewrite 71 !End of the abilint section 72 73 implicit none 74 75 !Arguments ------------------------------- 76 type(dataset_type), intent(in) :: dtset 77 type(energies_type), intent(inout) :: energies 78 integer, intent(in) :: istep,mcprj 79 type(MPI_type), intent(in) :: mpi_enreg 80 real(dp), intent(inout) :: residm 81 type(wvl_data), intent(inout) :: wvl 82 !arrays 83 real(dp), intent(in) :: xcart(3, dtset%natom) 84 type(pawcprj_type),dimension(dtset%natom,mcprj),intent(inout)::cprj 85 86 !Local variables------------------------------- 87 #if defined HAVE_BIGDFT 88 integer :: iatom,icprj 89 character(len = 500) :: message 90 real(dp), save :: etotal_local 91 integer, save :: ids 92 real(dp) :: gnrm_zero 93 integer :: comm,me,nproc 94 integer :: nlmn(dtset%natom) 95 #endif 96 97 98 ! ********************************************************************* 99 100 DBG_ENTER("COLL") 101 102 #if defined HAVE_BIGDFT 103 104 if(wvl%wfs%ks%orthpar%methOrtho .ne. 0) then 105 write(message,'(2a)') ch10,& 106 & 'wvl_hpsitopsi: the only orthogonalization method supported for PAW+WVL is Cholesky' 107 MSG_ERROR(message) 108 end if 109 110 write(message, '(a,a)' ) ch10,& 111 & ' wvl_hpsitopsi: compute the new wavefunction from the trial potential.' 112 call wrtout(std_out,message,'COLL') 113 114 comm=mpi_enreg%comm_wvl 115 me=xmpi_comm_rank(comm) 116 nproc=xmpi_comm_size(comm) 117 118 !Initialisation of mixing parameter 119 if (istep == 1) then 120 etotal_local = real(1.d100, dp) 121 ids = dtset%nwfshist 122 end if 123 124 !WARNING! e_hartree is taken from the previous iteration as e_xc 125 !Update physical values 126 127 !Precondition, minimise (DIIS or steepest descent) and ortho. 128 !Compute also the norm of the gradient. 129 if(dtset%usepaw==1) then 130 call calculate_energy_and_gradient(istep, me, nproc, wvl%wfs%GPU, dtset%wvl_nprccg, & 131 & dtset%iscf, wvl%e%energs, wvl%wfs%ks, residm, gnrm_zero,wvl%descr%paw) 132 else 133 call calculate_energy_and_gradient(istep, me, nproc, wvl%wfs%GPU, dtset%wvl_nprccg, & 134 & dtset%iscf, wvl%e%energs, wvl%wfs%ks, residm, gnrm_zero) 135 end if 136 etotal_local = wvl%wfs%ks%diis%energy 137 138 if(dtset%usepaw==1) then 139 call hpsitopsi(me, nproc, istep, ids, wvl%wfs%ks,& 140 & wvl%descr%atoms,wvl%projectors%nlpsp,& 141 & wvl%descr%paw,xcart,energies%e_nonlocalpsp,wvl%projectors%G) 142 else 143 call hpsitopsi(me, nproc, istep, ids, wvl%wfs%ks,& 144 & wvl%descr%atoms,wvl%projectors%nlpsp) 145 end if 146 147 if(dtset%usepaw==1) then 148 ! PENDING : cprj should not be copied' 149 150 ! Cannot use pawcprj_copy because cprj and paw%cprj are not the same objects 151 ! Get nlmn from bigdft cprj, and allocate our copy 152 do iatom=1,dtset%natom 153 nlmn(iatom) = wvl%descr%paw%cprj(iatom,1)%nlmn 154 end do 155 call pawcprj_alloc(cprj,mcprj,nlmn) 156 157 do iatom=1,dtset%natom 158 do icprj=1,mcprj 159 cprj(iatom,icprj)%cp(:,:)= wvl%descr%paw%cprj(iatom,icprj)%cp(:,:) 160 end do 161 end do 162 163 164 end if 165 166 #else 167 BIGDFT_NOTENABLED_ERROR() 168 if (.false.) write(std_out,*) dtset%nstep,energies%e_ewald,istep,mcprj,mpi_enreg%nproc,residm,& 169 & wvl%wfs%ks,xcart(1,1),cprj(1,1)%nlmn 170 #endif 171 172 DBG_EXIT("COLL") 173 174 end subroutine wvl_hpsitopsi