TABLE OF CONTENTS
ABINIT/simple_j_dia [ Functions ]
NAME
simple_j_dia
FUNCTION
simple test current for H atoms
COPYRIGHT
Copyright (C) 2005-2018 ABINIT group (JJ,MT) 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
NOTES
PARENTS
CHILDREN
acrossb
SOURCE
28 #if defined HAVE_CONFIG_H 29 #include "config.h" 30 #endif 31 32 #include "abi_common.h" 33 34 subroutine simple_j_dia(jdia,natom,nfft,pawfgrtab) 35 36 37 use defs_basis 38 use m_profiling_abi 39 use m_errors 40 41 use m_geometry, only : acrossb 42 use m_pawfgrtab, only : pawfgrtab_type 43 44 !This section has been created automatically by the script Abilint (TD). 45 !Do not modify the following lines by hand. 46 #undef ABI_FUNC 47 #define ABI_FUNC 'simple_j_dia' 48 !End of the abilint section 49 50 implicit none 51 52 !Arguments ------------------------------------ 53 !scalars 54 integer,intent(in) :: natom,nfft 55 !arrays 56 real(dp),intent(out) :: jdia(3,3,nfft) 57 type(pawfgrtab_type),intent(in) :: pawfgrtab(natom) 58 59 !Local variables------------------------------- 60 !scalars 61 integer :: iatom,idir,ifgd,ifftsph 62 real(dp) :: Bvecsum,nrm,nrmmax,psi2,scale 63 !arrays 64 real(dp) :: axb(3),B(3,3),rvec(3) 65 real(dp),allocatable :: Bvec(:,:) 66 67 ! ************************************************************************ 68 69 DBG_ENTER("COLL") 70 71 !make sure current field is set to zero everywhere 72 jdia(:,:,:) = zero 73 74 !define the directions of the external B field 75 B = zero 76 B(1,1) = one; B(2,2) = one; B(3,3) = one; 77 78 !loop over atoms in cell 79 scale = half 80 nrm = zero 81 nrmmax = zero 82 Bvecsum = zero 83 do iatom = 1, natom 84 ABI_ALLOCATE(Bvec,(3,pawfgrtab(iatom)%nfgd)) 85 do ifgd=1, pawfgrtab(iatom)%nfgd 86 ifftsph = pawfgrtab(iatom)%ifftsph(ifgd) 87 rvec(:) = pawfgrtab(iatom)%rfgd(:,ifgd) 88 nrm = sqrt(dot_product(rvec,rvec)) 89 if (nrm > nrmmax) nrmmax = nrm 90 psi2=exp(-two*nrm/scale)/(pi*scale*scale*scale) 91 do idir=1, 3 92 call acrossb(B(idir,:),rvec,axb) 93 jdia(idir,:,ifftsph) = -psi2*axb(:)/(two*Sp_Lt) 94 end do ! end loop over idir for jdia 95 if (nrm > zero) then 96 call acrossb(rvec,jdia(3,:,ifftsph),axb) 97 Bvec(:,ifgd) = axb(:)/(Sp_Lt*nrm*nrm*nrm) 98 end if 99 Bvecsum = Bvecsum + Bvec(3,ifgd) 100 end do ! end loop over nfgd points in sphere 101 write(std_out,'(i8,f12.8,f12.8)')pawfgrtab(iatom)%nfgd,Bvecsum,Bvecsum*(four_pi*third*nrm**3)/pawfgrtab(iatom)%nfgd 102 ABI_DEALLOCATE(Bvec) 103 end do ! end loop over atoms in cell 104 105 DBG_EXIT("COLL") 106 107 end subroutine simple_j_dia