TABLE OF CONTENTS


ABINIT/vn_nl_rec [ Functions ]

[ Top ] [ Functions ]

NAME

 vn_nl_rec

FUNCTION

 this routine computes the contribution to the vector vn, during
 recursion, due to the non-local psp.

COPYRIGHT

 Copyright (C) 2008-2018 ABINIT group ( ).
 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

  vn(:,:,:)=the vector on the real-space grid.
  inf_ucvol=volume of infinitesimal cell
  natom=number of atoms
  typat(natom)=the type of psps associated to the atoms
  ngfftrec(3)=first 3 components of ngfftrec (truncated box, if different from ngfft) for the real-space grid
  nlrec<type(nlpsprec_type)> in recursion_type containing information concerning psp 
  projec(ngfftrec(1),ngfftrec(2),ngfftrec(3),lmnmax,natom) is the  vector, on the ngfftrec grid containing 
  the non-lacal projector $Y_{lm}(r-R_A)f_{lk}(r-R_A)

OUTPUT

 vn_nl(:,:,:)=the non_local contribution to vn

PARENTS

      recursion,recursion_nl

CHILDREN

      timab

NOTES

SOURCE

 41 #if defined HAVE_CONFIG_H
 42 #include "config.h"
 43 #endif
 44 
 45 #include "abi_common.h"
 46 
 47 
 48 subroutine vn_nl_rec(vn,natom,typat,ngfftrec,&
 49   &                    inf_ucvol,nlrec,projec)
 50 
 51  use m_profiling_abi
 52  
 53  use defs_basis
 54  use defs_rectypes
 55  use m_per_cond,only             :per_cond
 56  use m_linalg_interfaces
 57 
 58 !This section has been created automatically by the script Abilint (TD).
 59 !Do not modify the following lines by hand.
 60 #undef ABI_FUNC
 61 #define ABI_FUNC 'vn_nl_rec'
 62  use interfaces_18_timing
 63 !End of the abilint section
 64 
 65  implicit none
 66  
 67 !Arguments -------------------------------
 68 !scalars
 69  integer,intent(in) :: natom
 70  real(dp),intent(in) :: inf_ucvol
 71  type(nlpsprec_type),intent(in) :: nlrec
 72 !arrays
 73  integer,intent(in) :: ngfftrec(3),typat(natom)
 74  real(dp),intent(in) :: projec(0:,0:,0:,1:,1:)
 75  real(dp),intent(inout):: vn(0:ngfftrec(1)*ngfftrec(2)*ngfftrec(3)-1)
 76 !Local variables-------------------------------
 77 !scalars
 78  integer :: iatom,nfftrec
 79  integer :: jlmn,il,in,jn
 80  integer :: ipsp,ilmn
 81  integer :: npsp,lmnmax
 82  real(dp):: vn_nl_loc 
 83 !arrays
 84  real(dp):: vn_nl(0:ngfftrec(1)-1,0:ngfftrec(2)-1,0:ngfftrec(3)-1)
 85  real(dp):: vtempo(0:ngfftrec(1)-1,0:ngfftrec(2)-1,0:ngfftrec(3)-1)
 86  real(dp):: tsec(2)
 87 ! *************************************************************************
 88 
 89  call timab(615,1,tsec)
 90 !--Initialisation
 91  
 92  vn_nl = zero
 93  npsp = nlrec%npsp
 94  lmnmax = nlrec%lmnmax
 95  nfftrec = product(ngfftrec)
 96  vtempo(:,:,:) = reshape(source=vn,shape=ngfftrec(:3))
 97 
 98 !--Sum_iatom \int dr1 E(r-r_a,r1-r_a)vn(r1) *infucvol
 99  do iatom=1,natom !--Loop on atoms
100    ipsp = typat(natom)
101 
102 !  --If psp(typat(iatom)) is local then cycle
103    if(all(nlrec%pspinfo(:,ipsp)==0))  cycle 
104 
105    
106 !  write(std_out,*)'lmnmax',nlrec%lmnmax,lmnmax
107 
108    do ilmn = 1, lmnmax
109      do jlmn = 1,lmnmax
110        if(nlrec%indlmn(4,ilmn,ipsp)==nlrec%indlmn(4,jlmn,ipsp)) then
111          il = 1+nlrec%indlmn(1,jlmn,ipsp)
112          in = nlrec%indlmn(3,ilmn,ipsp)
113          jn = nlrec%indlmn(3,jlmn,ipsp)
114          vn_nl_loc = ddot(nfftrec,projec(:,:,:,jlmn,iatom),1,vtempo,1)
115          vn_nl = vn_nl+projec(:,:,:,ilmn,iatom)*vn_nl_loc*nlrec%mat_exp_psp_nl(in,jn,il,ipsp)
116        end if
117      end do
118    end do
119  end do !--End loop on atoms
120  vtempo = vtempo + vn_nl*inf_ucvol
121 
122  vn = reshape(source=vtempo,shape=(/nfftrec/))
123 
124  call timab(615,2,tsec)
125 end subroutine vn_nl_rec