TABLE OF CONTENTS
ABINIT/Psolver_hartree [ Functions ]
NAME
Psolver_hartree
FUNCTION
Given rho(r), compute Hartree potential considering the system as an isolated one. This potential is obtained from the convolution of 1/r and rho(r), treated in Fourier space. This method is a wrapper around Psolver() developped for BigDFT. It does not compute the xc energy nor potential. See psolver_rhohxc() to do it. WARNING : the XC energy and potential computation capability has been for spin-polarized case, as everything is done as if nspden=1
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR). 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)>=all input variables in this dataset mpi_enreg=MPI-parallelisation information. rhor(nfft,nspden)=electron density in real space in electrons/bohr**3
OUTPUT
enhartr=returned Hartree energy (hartree). vhartr(nfft)=Hartree potential. NOTE In PSolver, with nspden == 2, rhor(:,1) = density up and rhor(:,2) = density down. But in ABINIT (dtset%usewvl != 1) rhor(:,1) = total density and rhor(:,2) = density up . In ABINIT (dtset%usewvl != 1), the same convention is used as in PSolver.
PARENTS
mklocl_realspace,nres2vres
CHILDREN
h_potential,psolver_kernel,wrtout
SOURCE
46 #if defined HAVE_CONFIG_H 47 #include "config.h" 48 #endif 49 50 #include "abi_common.h" 51 52 53 subroutine psolver_hartree(enhartr, hgrid, icoulomb, me, mpi_comm, nfft, ngfft, nproc, & 54 & nscforder, nspden, rhor, vhartr, usewvl) 55 56 use defs_basis 57 use m_errors 58 use m_profiling_abi 59 #if defined HAVE_BIGDFT 60 use BigDFT_API, only : coulomb_operator 61 use poisson_solver, only : H_potential 62 #endif 63 64 !This section has been created automatically by the script Abilint (TD). 65 !Do not modify the following lines by hand. 66 #undef ABI_FUNC 67 #define ABI_FUNC 'psolver_hartree' 68 use interfaces_14_hidewrite 69 use interfaces_62_poisson, except_this_one => psolver_hartree 70 !End of the abilint section 71 72 implicit none 73 74 !Arguments ------------------------------------ 75 !scalars 76 integer, intent(in) :: nfft, nspden, icoulomb, usewvl, mpi_comm, me, nproc, nscforder 77 real(dp), intent(out) :: enhartr 78 !arrays 79 integer, intent(in) :: ngfft(3) 80 real(dp),intent(in) :: hgrid(3) 81 real(dp),intent(in) :: rhor(nfft,nspden) 82 real(dp),intent(out) :: vhartr(nfft) 83 84 !Local variables------------------------------- 85 #if defined HAVE_BIGDFT 86 !scalars 87 character(len=500) :: message 88 character(len = 1) :: datacode, bndcode 89 !arrays 90 real(dp), dimension(1) :: pot_ion_dummy 91 type(coulomb_operator):: kernel 92 #endif 93 94 ! ********************************************************************* 95 96 #if defined HAVE_BIGDFT 97 98 if (icoulomb == 0) then 99 ! The kernel is built with 'P'eriodic boundary counditions. 100 bndcode = 'P' 101 else if (icoulomb == 1) then 102 ! The kernel is built with 'F'ree boundary counditions. 103 bndcode = 'F' 104 else if (icoulomb == 2) then 105 ! The kernel is built with 'S'urface boundary counditions. 106 bndcode = 'S' 107 end if 108 109 if(nspden > 2 .and. usewvl/=0 )then 110 write(message, '(a,a,a,i0)' )& 111 & 'Only non-spin-polarised or collinear spin is allowed for wavelets,',ch10,& 112 & 'while the argument nspden = ', nspden 113 MSG_BUG(message) 114 end if 115 116 !We do the computation. 117 write(message, "(A,A,A,3I6)") "Psolver_hartree(): compute potential (Vhartree)...", ch10, & 118 & " | dimension:", ngfft(1:3) 119 call wrtout(std_out, message,'COLL') 120 121 if (usewvl == 0) then 122 vhartr(:) = rhor(:, 1) 123 124 datacode = 'G' 125 ! This may not work with MPI in the planewave code... 126 else 127 if(nspden==1)vhartr(:) = rhor(:, 1) 128 if(nspden==2)vhartr(:) = rhor(:, 1) + rhor(:, 2) 129 ! The data are 'D'istributed in the wavelet case or 'G'lobal otherwise. 130 if (nproc > 1) then 131 datacode = 'D' 132 else 133 datacode = 'G' 134 end if 135 end if 136 137 !We get the kernel. 138 call psolver_kernel( hgrid, 2, icoulomb, me, kernel, mpi_comm, ngfft, nproc, nscforder) 139 140 141 !We attack PSolver with the total density contained in vhartr. 142 !This is also valid for spin-polarized (collinear and non-collinear) 143 !systems. Thus we enter nspden (last arg of PSolver) as being 1. 144 !Warning : enxc and evxc are meaningless. 145 ! call psolver(bndcode, datacode, me, nproc, ngfft(1), ngfft(2), ngfft(3),& 146 !& 0, hgrid(1), hgrid(2), hgrid(3), vhartr, kernel%co%kernel, pot_ion_dummy, & 147 !& enhartr, enxc, evxc, 0.d0, .false., 1) 148 149 call H_potential(datacode,kernel,vhartr,pot_ion_dummy,& 150 & enhartr,zero,.false.) 151 152 153 #else 154 BIGDFT_NOTENABLED_ERROR() 155 if (.false.) write(std_out,*) nfft,nspden,icoulomb,usewvl,mpi_comm,me,nproc,nscforder,enhartr,& 156 & ngfft(1),hgrid(1),rhor(1,1),vhartr(1) 157 #endif 158 159 end subroutine psolver_hartree