TABLE OF CONTENTS
ABINIT/metric_so [ Functions ]
NAME
metric_so
FUNCTION
Computes Pauli matrices and antisymmetric tensor needed for spin-orbit.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, GZ, 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 . For the initials of contributors, see ~abinit/doc/developers/contributors.txt.
INPUTS
gprimd(3,3)=dimensional primitive translations for reciprocal space (bohr**-1)
OUTPUT
amet(2,3,3,2,2)=the antisymmetric tensor A(Re/Im,y,y'',s,s'') pauli(2,2,2,3)=Pauli matrixes
NOTES
(the preprocessing based on CPP does not allow isolated quotes, so that they have been doubled in the following latex equations) $2 Pauli(s,s'',1) & = & 0 & 1 \nonumber & & 1 & 0 \nonumber 2 Pauli(s,s'',2) & = & 0 & -i \nonumber & & i & 0 \nonumber 2 Pauli(s,s'',3) & = & 1 & 0 \nonumber & & 0 & -1 \end{eqnarray} }} $Amet(y,y'',s,s'') & = & -i Pauli(s,s'',n) E(n,m,m'') Gprimd(m,y) Gprimd(m'',y'') \end{eqnarray} }} E(n,m,m''): full antisymetric tensor s,s'': spin indices (1..2) y,y'',m,m'',n: metric indices (1..3) a,b: strain indices (1..3) Amet and Pauli are complex
PARENTS
nonlop_pl
CHILDREN
SOURCE
53 #if defined HAVE_CONFIG_H 54 #include "config.h" 55 #endif 56 57 #include "abi_common.h" 58 59 60 subroutine metric_so(amet,gprimd,pauli) 61 62 use defs_basis 63 use m_profiling_abi 64 65 !This section has been created automatically by the script Abilint (TD). 66 !Do not modify the following lines by hand. 67 #undef ABI_FUNC 68 #define ABI_FUNC 'metric_so' 69 !End of the abilint section 70 71 implicit none 72 73 !Arguments ------------------------------------ 74 !arrays 75 real(dp),intent(in) :: gprimd(3,3) 76 real(dp),intent(out) :: amet(2,3,3,2,2),pauli(2,2,2,3) 77 78 !Local variables------------------------------- 79 !scalars 80 integer :: iy1,iy2,m1,m2,n 81 !arrays 82 real(dp) :: buffer1(3,3,2,2) !,buffer2(3,3,3,3,2,2) 83 84 ! ********************************************************************** 85 86 !Fill in Pauli matrices and make them spin matrices: 87 !Pauli(Re/Im,up,down,coord): 88 pauli(:,:,:,:)=0.d0 89 pauli(1,1,2,1)= 1.d0;pauli(1,2,1,1)= 1.d0 90 pauli(2,1,2,2)=-1.d0;pauli(2,2,1,2)= 1.d0 91 pauli(1,1,1,3)= 1.d0;pauli(1,2,2,3)=-1.d0 92 pauli(:,:,:,:)= 0.5d0*pauli(:,:,:,:) 93 94 !Construct the antisymmetric tensor: 95 amet(:,:,:,:,:)=0.d0 96 do iy2=1,3 97 do iy1=1,3 98 do n=1,3 99 m1=mod(n ,3)+1 ! n,m1,m2 is an even permutation 100 m2=mod(m1,3)+1 101 amet(1:2,iy1,iy2,1:2,1:2) = amet(:,iy1,iy2,:,:) & 102 & + pauli(:,:,:,n) & 103 & *(gprimd(m1,iy1)*gprimd(m2,iy2) & 104 & -gprimd(m2,iy1)*gprimd(m1,iy2)) 105 end do 106 end do 107 end do 108 !amet= -i amet 109 buffer1(:,:,:,:)=-amet(1,:,:,:,:) 110 amet(1,:,:,:,:) = amet(2,:,:,:,:) 111 amet(2,:,:,:,:) =buffer1(:,:,:,:) 112 113 !DEBUG 114 !!Eventually construct the gradients of Amet wrt strains: 115 !! DAmet(y,y'',a,b,s,s'')= d[Amet(y,y'',s,s'')]/d[strain(a,b)] 116 !! -i Pauli(s,s'',n)* 117 !! ( E(n,a,m'')*Gprimd(b,y )*Gprimd(m'',y'') 118 !! +E(n,m,a )*Gprimd(b,y'')*Gprimd(m ,y ) ) 119 !damet(:,:,:,:,:,:,:)=0.d0 120 !do ib=1,3 121 !do ia=1,3 122 !m1=mod(ia,3)+1 ! ia,m1,m2 is an even permutation 123 !m2=mod(m1,3)+1 124 !do iy2=1,3 125 !do iy1=1,3 126 !damet(:,iy1,iy2,ia,ib,:,:) = damet(:,iy1,iy2,ia,ib,:,:) & 127 !& + (pauli(:,:,:,m2)*gprimd(m1,iy2) & 128 !-pauli(:,:,:,m1)*gprimd(m2,iy2))*gprimd(ib,iy1) & 129 !& + (pauli(:,:,:,m1)*gprimd(m2,iy1) & 130 !-pauli(:,:,:,m2)*gprimd(m1,iy1))*gprimd(ib,iy2) 131 !end do 132 !end do 133 !end do 134 !end do 135 !! damet= i damet 136 !buffer2(:,:,:,:,:,:)= damet(1,:,:,:,:,:,:) 137 !damet(1,:,:,:,:,:,:)= -damet(2,:,:,:,:,:,:) 138 !damet(2,:,:,:,:,:,:)=buffer2(:,:,:,:,:,:) 139 !! Symetrize damet(:,:,:,a,b,:,:) 140 !damet(:,:,:,1,2,:,:)=0.5d0*(damet(:,:,:,1,2,:,:)+damet(:,:,:,2,1,:,:)) 141 !damet(:,:,:,1,3,:,:)=0.5d0*(damet(:,:,:,1,3,:,:)+damet(:,:,:,3,1,:,:)) 142 !damet(:,:,:,2,3,:,:)=0.5d0*(damet(:,:,:,2,3,:,:)+damet(:,:,:,3,2,:,:)) 143 !damet(:,:,:,2,1,:,:)=damet(:,:,:,1,2,:,:) 144 !damet(:,:,:,3,1,:,:)=damet(:,:,:,1,3,:,:) 145 !damet(:,:,:,3,2,:,:)=damet(:,:,:,2,3,:,:) 146 !ENDDEBUG 147 148 end subroutine metric_so