TABLE OF CONTENTS


ABINIT/wvl_tail_corrections [ Functions ]

[ Top ] [ Functions ]

NAME

 wvl_tail_corrections

FUNCTION

 Perform a minimization on the wavefunctions (especially the treatment
 of the kinetic operator) with exponentialy decreasing functions on
 boundaries.

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

      afterscfloop

CHILDREN

      calculatetailcorrection,dcopy,wrtout,xmpi_allgatherv

SOURCE

 33 #if defined HAVE_CONFIG_H
 34 #include "config.h"
 35 #endif
 36 
 37 #include "abi_common.h"
 38 
 39 subroutine wvl_tail_corrections(dtset, energies, etotal, mpi_enreg, psps, wvl, xcart)
 40 
 41  use defs_basis
 42  use defs_datatypes
 43  use defs_abitypes
 44  use defs_wvltypes
 45  use m_xmpi
 46  use m_errors
 47  use m_profiling_abi
 48 
 49  use m_energies, only : energies_type
 50 #if defined HAVE_BIGDFT
 51   use BigDFT_API, only: CalculateTailCorrection
 52 #endif
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'wvl_tail_corrections'
 58  use interfaces_14_hidewrite
 59 !End of the abilint section
 60 
 61  implicit none
 62 
 63 !Arguments ------------------------------------
 64 !scalars
 65  real(dp),intent(out) :: etotal
 66  type(MPI_type),intent(in) :: mpi_enreg
 67  type(dataset_type),intent(in) :: dtset
 68  type(energies_type),intent(inout) :: energies
 69  type(pseudopotential_type),intent(in) :: psps
 70  type(wvl_data),intent(inout) :: wvl
 71 !arrays
 72  real(dp),intent(in) :: xcart(3,dtset%natom)
 73 
 74 !Local variables-------------------------------
 75 #if defined HAVE_BIGDFT
 76 !scalars
 77  integer :: ierr,me,nbuf,nproc,nsize,spaceComm
 78  real(dp) :: ekin_sum,epot_sum,eproj_sum
 79  logical :: parallel
 80  character(len=500) :: message
 81 !arrays
 82  integer :: ntails(3)
 83  real(dp) :: atails(3)
 84 #endif
 85 
 86 ! *************************************************************************
 87 
 88 #if defined HAVE_BIGDFT
 89 
 90  spaceComm=mpi_enreg%comm_wvl
 91  me=xmpi_comm_rank(spaceComm)
 92  nproc=xmpi_comm_size(spaceComm)
 93  parallel = (nproc > 1)
 94 
 95 !Write a message with the total energy before tail corrections.
 96  etotal = energies%e_kinetic + energies%e_hartree + energies%e_xc + &
 97 & energies%e_localpsp + energies%e_corepsp + energies%e_fock+&
 98 & energies%e_entropy + energies%e_elecfield + energies%e_magfield+&
 99 & energies%e_ewald + energies%e_chempot + energies%e_vdw_dftd
100  if (dtset%usepaw==0) etotal = etotal + energies%e_nonlocalpsp
101  if (dtset%usepaw/=0) etotal = etotal + energies%e_paw
102  write(message,'(a,2x,e19.12)') ' Total energy before tail correction', etotal
103  call wrtout(std_out, message, 'COLL')
104 
105 !Calculate kinetic energy correction due to boundary conditions
106  nbuf = nint(dtset%tl_radius / dtset%wvl_hgrid)
107  ntails = (/ wvl%descr%Glr%d%n1, wvl%descr%Glr%d%n2, wvl%descr%Glr%d%n3 /) + 2 * nbuf
108  atails = real(ntails, dp) * dtset%wvl_hgrid
109  write(message,'(a,a,i6,a,A,A,3F12.6,A,A,3I12,A)') ch10,&
110 & ' Tail requires ',nbuf,' additional grid points around cell.', ch10, &
111 & '  | new acell:', atails, ch10, &
112 & '  | new box size for wavelets:', ntails, ch10
113  call wrtout(std_out,message,'COLL')
114  call wrtout(ab_out,message,'COLL')
115 
116 
117 !Calculate energy correction due to finite size effects
118 !---reformat potential
119  nsize = wvl%descr%Glr%d%n1i * wvl%descr%Glr%d%n2i
120  ABI_ALLOCATE(wvl%den%denspot%pot_work, (nsize * wvl%descr%Glr%d%n3i * dtset%nsppol))
121 
122  if (parallel) then
123    call xmpi_allgatherv(wvl%den%denspot%rhov, &
124 &   nsize * wvl%den%denspot%dpbox%nscatterarr(me, 2), &
125 &   wvl%den%denspot%pot_work, nsize * wvl%den%denspot%dpbox%ngatherarr(:,2), &
126 &   nsize * wvl%den%denspot%dpbox%ngatherarr(:,3),spaceComm,ierr)
127  else
128    call dcopy(wvl%descr%Glr%d%n1i * wvl%descr%Glr%d%n2i * &
129 &   wvl%descr%Glr%d%n3i * dtset%nsppol,wvl%den%denspot%rhov,1,wvl%den%denspot%pot_work,1)
130  end if
131 
132  if(dtset%usepaw==1) then
133    call CalculateTailCorrection(me, nproc, wvl%descr%atoms, dtset%tl_radius, &
134 &   wvl%wfs%ks%orbs, wvl%wfs%ks%lzd%Glr, wvl%projectors%nlpsp, dtset%tl_nprccg, &
135 &   wvl%den%denspot%pot_work, dtset%wvl_hgrid, xcart, psps%gth_params%radii_cf, &
136 &   dtset%wvl_crmult, dtset%wvl_frmult, dtset%nsppol, &
137 &   wvl%wfs%ks%psi, .false., ekin_sum, epot_sum, eproj_sum,&
138 &   wvl%projectors%G,wvl%descr%paw)
139  else
140    call CalculateTailCorrection(me, nproc, wvl%descr%atoms, dtset%tl_radius, &
141 &   wvl%wfs%ks%orbs, wvl%wfs%ks%lzd%Glr, wvl%projectors%nlpsp, dtset%tl_nprccg, &
142 &   wvl%den%denspot%pot_work, dtset%wvl_hgrid, xcart, psps%gth_params%radii_cf, &
143 &   dtset%wvl_crmult, dtset%wvl_frmult, dtset%nsppol, &
144 &   wvl%wfs%ks%psi, .false., ekin_sum, epot_sum, eproj_sum)
145  end if
146 
147  ABI_DEALLOCATE(wvl%den%denspot%pot_work)
148 
149  energies%e_kinetic = ekin_sum
150  energies%e_localpsp = epot_sum - two * energies%e_hartree
151  energies%e_nonlocalpsp = eproj_sum
152  energies%e_corepsp = zero
153  energies%e_chempot = zero
154 #if defined HAVE_BIGDFT
155  energies%e_localpsp = energies%e_localpsp - wvl%e%energs%evxc
156 #endif
157 
158  write(message,'(a,3(1x,e18.11))') ' ekin_sum,epot_sum,eproj_sum',  &
159  ekin_sum,epot_sum,eproj_sum
160  call wrtout(std_out, message, 'COLL')
161  write(message,'(a,2(1x,e18.11))') ' ehart,eexcu', &
162 & energies%e_hartree,energies%e_xc
163  call wrtout(std_out, message, 'COLL')
164 
165  etotal = energies%e_kinetic + energies%e_hartree + energies%e_xc + &
166 & energies%e_localpsp + energies%e_corepsp + energies%e_fock+&
167 & energies%e_entropy + energies%e_elecfield + energies%e_magfield+&
168 & energies%e_ewald + energies%e_vdw_dftd
169  if (dtset%usepaw==0) etotal = etotal + energies%e_nonlocalpsp
170  if (dtset%usepaw/=0) etotal = etotal + energies%e_paw
171 
172  write(message,'(a,2x,e19.12)') ' Total energy with tail correction', etotal
173  call wrtout(std_out, message, 'COLL')
174 
175 !--- End if of tail calculation
176 
177 #else
178  BIGDFT_NOTENABLED_ERROR()
179  if (.false.) write(std_out,*) etotal,mpi_enreg%nproc,dtset%nstep,energies%e_ewald,psps%npsp,&
180 & wvl%wfs%ks,xcart(1,1)
181 #endif
182 
183 end subroutine wvl_tail_corrections