TABLE OF CONTENTS


ABINIT/aprxdr [ Functions ]

[ Top ] [ Functions ]

NAME

 aprxdr

FUNCTION

 Compute the approximative derivatives of the energy at different
 points along the line search, thanks to a finite-difference formula.
 This formula is the projection along the line search of the
 Eq.(11) in PRB54, 4383 (1996).

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (XG)
 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

 cplex: if 1, real space functions on FFT grid are REAL, if 2, COMPLEX
 choice= if==3, compute dedv_new, dedv_old, and dedv_mix,
 if/=3, compute only dedv_new and dedv_old.
 i_vresid and i_rhor, see the next lines.
 f_fftgr(nfft,nspden,n_fftgr)=different functions defined on the fft grid :
 The last residual potential is in f_fftgr(:,:,i_vresid(1)).
 The old  residual potential is in f_fftgr(:,:,i_vresid(2)).
 The previous old residual potential is in f_fftgr(:,:,i_vresid(3)).
 (needed only when choice==3)
 The old  density is in f_fftgr(:,:,i_rhor2).
 f_atm(3,natom,n_fftgr)=different functions defined for each atom :
 The last HF force is in f_atm(:,:,i_vresid(1)).
 The old  HF force is in f_fftgr(:,:,i_vresid(2)).
 The previous old HF force is in f_fftgr(:,:,i_vresid(3)).
 (needed only when choice==3)
 The old atomic positions are in f_atm(:,:,i_rhor2)
 moved_atm_inside: if==1, the atoms are allowed to move.
 mpicomm=the mpi communicator used for the summation
 mpi_summarize=set it to .true. if parallelisation is done over FFT
 natom=number of atoms in unit cell
 nfft=(effective) number of FFT grid points (for this processor)
 nfftot=total number of FFT grid points
 nspden=number of spin-density components
 rhor(nfft,nspden)=actual density
 xred(3,natom)=reduced atomic coordinates

OUTPUT

 dedv_mix=approximate derivative from previous old residual
 dedv_new=approximate derivative from new residual
 dedv_old=approximate derivative from old residual (output only when choice==3)

NOTES

 Should be OpenMP parallelized

PARENTS

      scfcge

CHILDREN

      dotprodm_vn

SOURCE

 62 #if defined HAVE_CONFIG_H
 63 #include "config.h"
 64 #endif
 65 
 66 #include "abi_common.h"
 67 
 68 subroutine aprxdr(cplex,choice,dedv_mix,dedv_new,dedv_old,&
 69 &  f_atm,f_fftgr,i_rhor2,i_vresid,moved_atm_inside,&
 70 &  mpicomm,mpi_summarize,natom,nfft,nfftot,nspden,n_fftgr,rhor,ucvol,xred)
 71 
 72  use m_profiling_abi
 73 
 74  use defs_basis
 75 
 76 !This section has been created automatically by the script Abilint (TD).
 77 !Do not modify the following lines by hand.
 78 #undef ABI_FUNC
 79 #define ABI_FUNC 'aprxdr'
 80  use interfaces_56_mixing, except_this_one => aprxdr
 81 !End of the abilint section
 82 
 83  implicit none
 84 
 85 !Arguments ------------------------------------
 86 !scalars
 87  integer,intent(in) :: choice,cplex,i_rhor2,moved_atm_inside,n_fftgr,natom,nfft
 88  integer,intent(in) :: mpicomm,nfftot,nspden
 89  logical, intent(in) :: mpi_summarize
 90  real(dp),intent(in) :: ucvol
 91  real(dp),intent(out) :: dedv_mix,dedv_new,dedv_old
 92 !arrays
 93  integer,intent(in) :: i_vresid(3)
 94  real(dp),intent(in) :: f_atm(3,natom,n_fftgr)
 95  real(dp),intent(in) :: f_fftgr(cplex*nfft,nspden,n_fftgr)
 96  real(dp),intent(in) :: rhor(cplex*nfft,nspden),xred(3,natom)
 97 
 98 !Local variables-------------------------------
 99 !scalars
100  integer :: iatom,idir
101 !arrays
102  real(dp) :: dedv_temp(1)
103  real(dp),allocatable :: ddens(:,:,:)
104 
105 ! *************************************************************************
106 
107  ABI_ALLOCATE(ddens,(cplex*nfft,nspden,1))
108 
109 !Compute approximative derivative of the energy
110 !with respect to change of potential
111 
112  ddens(:,:,1)=rhor(:,:)-f_fftgr(:,:,i_rhor2)
113 
114 !call dotprod_vn(cplex,1,ddens,dedv_old,nfft,nfftot,nspden,1,vresid,ucvol)
115 !Dot product ddens(:,:,1) f_fftgr(:,:,i_vresid(2))
116  call dotprodm_vn(cplex,1,ddens,dedv_temp,1,i_vresid(2),mpicomm,mpi_summarize,1,1,1,&
117 & nfft,nfftot,n_fftgr,nspden,f_fftgr,ucvol)
118  dedv_old = dedv_temp(1)
119 
120 !Dot product ddens(:,:,1) f_fftgr(:,:,i_vresid(1))
121  call dotprodm_vn(cplex,1,ddens,dedv_temp,1,i_vresid(1),mpicomm,mpi_summarize,1,1,1,&
122 & nfft,nfftot,n_fftgr,nspden,f_fftgr,ucvol)
123  dedv_new= dedv_temp(1)
124 
125  if(choice==3)then
126 !  Dot product ddens(:,:,1) f_fftgr(:,:,i_vresid(3))
127    call dotprodm_vn(cplex,1,ddens,dedv_temp,1,i_vresid(3),mpicomm,mpi_summarize,1,1,1,&
128 &   nfft,nfftot,n_fftgr,nspden,f_fftgr,ucvol)
129    dedv_mix = dedv_temp(1)
130  end if
131 
132  ABI_DEALLOCATE(ddens)
133 
134 !-------------------------------------------------------
135 
136 !Now, take care of eventual atomic displacements
137 
138  if(moved_atm_inside==1)then
139    do idir=1,3
140      do iatom=1,natom
141        dedv_new=dedv_new+&
142 &       f_atm(idir,iatom,i_vresid(1))*(xred(idir,iatom)-f_atm(idir,iatom,i_rhor2))
143        dedv_old=dedv_old+&
144 &       f_atm(idir,iatom,i_vresid(2))*(xred(idir,iatom)-f_atm(idir,iatom,i_rhor2))
145        if(choice==3) dedv_mix=dedv_mix+&
146 &       f_atm(idir,iatom,i_vresid(3))*(xred(idir,iatom)-f_atm(idir,iatom,i_rhor2))
147      end do
148    end do
149  end if
150 
151 end subroutine aprxdr