The Elk Code
 
Loading...
Searching...
No Matches
gencrm.f90
Go to the documentation of this file.
1
2! Copyright (C) 2016 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6pure subroutine gencrm(n,wf11,wf12,wf21,wf22,crho,ld,cmag)
7use modmain
8implicit none
9! arguments
10integer, intent(in) :: n
11complex(4), intent(in) :: wf11(n),wf12(n),wf21(n),wf22(n)
12complex(4), intent(out) :: crho(n)
13integer, intent(in) :: ld
14complex(4), intent(out) :: cmag(ld,ndmag)
15! local variables
16integer i
17complex(4) c11,c12,c21,c22,c1,c2
18if (ncmag) then
19! non-collinear case
20!$OMP SIMD PRIVATE(c11,c12,c21,c22,c1,c2) SIMDLEN(8)
21 do i=1,n
22 c11=wf11(i); c12=wf12(i)
23 c21=wf21(i); c22=wf22(i)
24! up-dn spin density
25 c1=conjg(c11)*c22
26! dn-up spin density
27 c2=conjg(c12)*c21
28! x-component: up-dn + dn-up
29 cmag(i,1)=c1+c2
30! y-component: i*(dn-up - up-dn)
31 cmag(i,2)=ci*(c2-c1)
32 c1=conjg(c11)*c21
33 c2=conjg(c12)*c22
34! z-component: up-up - dn-dn
35 cmag(i,3)=c1-c2
36! density: up-up + dn-dn
37 crho(i)=c1+c2
38 end do
39else
40! collinear case
41!$OMP SIMD PRIVATE(c1,c2) SIMDLEN(8)
42 do i=1,n
43 c1=conjg(wf11(i))*wf21(i)
44 c2=conjg(wf12(i))*wf22(i)
45 cmag(i,1)=c1-c2
46 crho(i)=c1+c2
47 end do
48end if
49end subroutine
50
pure subroutine gencrm(n, wf11, wf12, wf21, wf22, crho, ld, cmag)
Definition gencrm.f90:7
logical ncmag
Definition modmain.f90:240
complex(4), parameter ci
Definition modmain.f90:1237
integer ndmag
Definition modmain.f90:238