File: genzrm.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (52 lines) | stat: -rw-r--r-- 1,181 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52

! Copyright (C) 2016 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

subroutine genzrm(n,wf11,wf12,wf21,wf22,zrho,ld,zmag)
use modmain
implicit none
! arguments
integer, intent(in) :: n
complex(8), intent(in) :: wf11(n),wf12(n),wf21(n),wf22(n)
complex(8), intent(out) :: zrho(n)
integer, intent(in) :: ld
complex(8), intent(out) :: zmag(ld,ndmag)
! local variables
integer i
complex(8) z11,z12,z21,z22,z1,z2
if (ncmag) then
! non-collinear case
  do i=1,n
    z11=wf11(i)
    z12=wf12(i)
    z21=wf21(i)
    z22=wf22(i)
! up-dn spin density
    z1=conjg(z11)*z22
! dn-up spin density
    z2=conjg(z12)*z21
! x-component: up-dn + dn-up
    zmag(i,1)=z1+z2
! y-component: i*(dn-up - up-dn)
    z1=z2-z1
    zmag(i,2)=cmplx(-aimag(z1),dble(z1),8)
    z1=conjg(z11)*z21
    z2=conjg(z12)*z22
! z-component: up-up - dn-dn
    zmag(i,3)=z1-z2
! density: up-up + dn-dn
    zrho(i)=z1+z2
  end do
else
! collinear case
  do i=1,n
    z1=conjg(wf11(i))*wf21(i)
    z2=conjg(wf12(i))*wf22(i)
    zmag(i,1)=z1-z2
    zrho(i)=z1+z2
  end do
end if
return
end subroutine