The Elk Code
genhafspt.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2023 E. Harris-Lee, S. Sharma and J. K. Dewhurst.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine genhafspt(evecsv,pmat,h)
7 use modmain
8 use modtddft
9 implicit none
10 ! arguments
11 complex(8), intent(in) :: evecsv(nstsv,nstsv),pmat(nstsv,nstsv,3)
12 complex(8), intent(inout) :: h(nstsv,nstsv)
13 ! local variables
14 integer jst,i,j,l
15 real(8) ca,v(3)
16 complex(8) zm(2,2)
17 ! allocatable arrays
18 complex(8), allocatable :: smat(:,:,:,:)
19 allocate(smat(nstsv,nstsv,2,2))
20 ! generate the spin operator matrix elements in the second-variational basis
21 call gensmatk(evecsv,smat)
22 ! coupling constant of the external spin-polarised A-field (-1/c)
23 ca=-1.d0/solsc
24 ! add the spin-current operator to the Hamiltonian
25 do l=1,3
26  v(1:3)=ca*afspt(l,1:3,itimes)
27  zm(1,1)=v(3)
28  zm(1,2)=cmplx(v(1),-v(2),8)
29  zm(2,1)=cmplx(v(1),v(2),8)
30  zm(2,2)=-v(3)
31  do j=1,2
32  do i=1,2
33  call zhemm('L','U',nstsv,nstsv,zm(i,j),pmat(:,:,l),nstsv,smat(:,:,i,j), &
34  nstsv,zone,h,nstsv)
35  end do
36  end do
37 end do
38 ! add effective magnetic field to Hamiltonian if required (E. Harris-Lee)
39 if (tbaspat) then
40  call r3mtv(afspt(:,:,itimes),afieldt(:,itimes),v)
41  v(1:3)=(ca**2)*v(1:3)
42  zm(1,1)=v(3)
43  zm(1,2)=cmplx(v(1),-v(2),8)
44  zm(2,1)=cmplx(v(1),v(2),8)
45  zm(2,2)=-v(3)
46  do j=1,2
47  do i=1,2
48  do jst=1,nstsv
49  h(1:jst,jst)=h(1:jst,jst)+zm(i,j)*smat(1:jst,jst,i,j)
50  end do
51  end do
52  end do
53 end if
54 deallocate(smat)
55 end subroutine
56 
pure subroutine r3mtv(a, x, y)
Definition: r3mtv.f90:10
subroutine genhafspt(evecsv, pmat, h)
Definition: genhafspt.f90:7
complex(8), parameter zone
Definition: modmain.f90:1240
subroutine gensmatk(evecsv, smat)
Definition: gensmatk.f90:7
logical tbaspat
Definition: modtddft.f90:65
real(8) solsc
Definition: modmain.f90:1253
real(8), dimension(:,:), allocatable afieldt
Definition: modtddft.f90:58
integer itimes
Definition: modtddft.f90:46
real(8), dimension(:,:,:), allocatable afspt
Definition: modtddft.f90:62