TABLE OF CONTENTS


ABINIT/wvl_psitohpsi [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_psitohpsi

FUNCTION

 Compute new trial potential and calculate the hamiltionian application into hpsi.

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, MT)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/dpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  mpi_enreg=information about MPI parallelization

OUTPUT

  vxc(nfft,nspden)=exchange-correlation potential (hartree)
  vtrial(nfft,nspden)=new potential

NOTES

PARENTS

      afterscfloop,rhotov,setvtr,vtorho

CHILDREN

      psitohpsi,total_energies,wrtout,wvl_vtrial_abi2big,wvl_vxc_abi2big

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 
 40 subroutine wvl_psitohpsi(alphamix,eexctX, eexcu, ehart, ekin_sum, epot_sum, eproj_sum, eSIC_DC, &
 41      & itrp, iter, iscf, me, natom, nfft, nproc, nspden, rpnrm, scf, &
 42      & vexcu, wvl, wvlbigdft, xcart, xcstr,vtrial,vxc)
 43 
 44  use defs_basis
 45  use defs_wvltypes
 46  use m_profiling_abi
 47  use m_errors
 48  use m_abi2big, only: wvl_vxc_abi2big,wvl_vtrial_abi2big
 49 #if defined HAVE_BIGDFT
 50  use BigDFT_API, only: psitohpsi, KS_POTENTIAL, total_energies
 51 #endif
 52 
 53 !This section has been created automatically by the script Abilint (TD).
 54 !Do not modify the following lines by hand.
 55 #undef ABI_FUNC
 56 #define ABI_FUNC 'wvl_psitohpsi'
 57  use interfaces_14_hidewrite
 58 !End of the abilint section
 59 
 60  implicit none
 61 
 62 !Arguments-------------------------------
 63 !scalars
 64  integer, intent(in) :: me, nproc, itrp, iter, iscf, natom, nfft, nspden
 65  real(dp), intent(in) :: alphamix
 66  real(dp), intent(out) :: rpnrm
 67  logical, intent(in) :: scf
 68  logical, intent(in) :: wvlbigdft
 69  type(wvl_data), intent(inout) :: wvl
 70  real(dp), intent(inout) :: eexctX,eSIC_DC,ehart,eexcu,vexcu, ekin_sum, epot_sum, eproj_sum
 71  real(dp), dimension(6), intent(out) :: xcstr
 72  real(dp), intent(inout) :: xcart(3, natom)
 73 !arrays
 74  real(dp),intent(out), optional :: vxc(nfft,nspden)
 75  real(dp),intent(out), optional :: vtrial(nfft,nspden)
 76 
 77 !Local variables-------------------------------
 78 !scalars
 79 #if defined HAVE_BIGDFT
 80  character(len=500) :: message
 81  integer :: linflag = 0
 82  character(len=3), parameter :: unblock_comms = "OFF"
 83 #endif
 84 
 85 ! *************************************************************************
 86 
 87  DBG_ENTER("COLL")
 88 
 89 #if defined HAVE_BIGDFT
 90 
 91  if(wvl%descr%atoms%npspcode(1)==7) then
 92    call psitohpsi(me,nproc,wvl%descr%atoms,scf,wvl%den%denspot, &
 93 &   itrp, iter, iscf, alphamix,&
 94 &   wvl%projectors%nlpsp,xcart,linflag,unblock_comms, &
 95 &   wvl%wfs%GPU,wvl%wfs%ks,wvl%e%energs,rpnrm,xcstr,&
 96 &   wvl%projectors%G,wvl%descr%paw)
 97  else
 98    call psitohpsi(me,nproc,wvl%descr%atoms,scf,wvl%den%denspot, &
 99 &   itrp, iter, iscf, alphamix,&
100 &   wvl%projectors%nlpsp,xcart,linflag,unblock_comms, &
101 &   wvl%wfs%GPU,wvl%wfs%ks,wvl%e%energs,rpnrm,xcstr)
102  end if
103 
104  if(scf) then
105    ehart     = wvl%e%energs%eh
106    eexcu     = wvl%e%energs%exc
107    vexcu     = wvl%e%energs%evxc
108  end if
109  eexctX    = wvl%e%energs%eexctX
110  eSIC_DC   = wvl%e%energs%evsic
111  ekin_sum  = wvl%e%energs%ekin
112  eproj_sum = wvl%e%energs%eproj
113  epot_sum  = wvl%e%energs%epot
114 
115 !Correct local potential, since in BigDFT 
116 !this variable contains more terms
117 !Do the following only if sumpion==.true. in psolver_rhohxc.
118 !For the moment it is set to false.
119 
120  epot_sum=epot_sum-real(2,dp)*wvl%e%energs%eh
121  epot_sum=epot_sum-wvl%e%energs%evxc
122 
123  if(wvlbigdft) then
124    call total_energies(wvl%e%energs, iter, me)
125  end if
126 
127 !Note: if evxc is not rested here,
128 !we have to rest this from etotal in prtene, afterscfcv and etotfor.
129 !check ABINIT-6.15.1.
130 
131  if(scf) then
132    if (present(vxc)) then
133      write(message, '(a,a,a,a)' ) ch10, ' wvl_psitohpsi : but why are you copying vxc :..o('
134      call wrtout(std_out,message,'COLL')
135      call wvl_vxc_abi2big(2,vxc,wvl%den)
136    end if
137    if (wvl%den%denspot%rhov_is == KS_POTENTIAL .and. present(vtrial)) then
138      write(message, '(a,a,a,a)' ) ch10, ' wvl_psitohpsi : but why are you copying vtrial :..o('
139      call wrtout(std_out,message,'COLL')
140      call wvl_vtrial_abi2big(2,vtrial,wvl%den)
141    end if
142  end if
143 
144 #else
145  BIGDFT_NOTENABLED_ERROR()
146  if (.false.) write(std_out,*) me,nproc,itrp,iter,iscf,natom,nfft,nspden,alphamix,rpnrm,scf,&
147 & wvlbigdft,wvl%wfs%ks,eexctX,eSIC_DC,ehart,eexcu,vexcu,ekin_sum,&
148 & epot_sum,eproj_sum,xcstr(1),xcart(1,1),vxc(1,1),vtrial(1,1)
149 #endif
150 
151  DBG_EXIT("COLL")
152 
153 end subroutine wvl_psitohpsi