The Elk Code
rvfcross.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2007 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 !BOP
7 ! !ROUTINE: rvfcross
8 ! !INTERFACE:
9 subroutine rvfcross(rvfmt1,rvfir1,rvfmt2,rvfir2,rvfmt3,rvfir3)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! rvfmt1 : first input muffin-tin field (in,real(npmtmax,natmtot,3))
14 ! rvfir1 : first input interstitial field (in,real(ngtot,3))
15 ! rvfmt2 : second input muffin-tin field (in,real(npmtmax,natmtot,3))
16 ! rvfir2 : second input interstitial field (in,real(ngtot,3))
17 ! rvfmt3 : output muffin-tin field (out,real(npmtmax,natmtot,3))
18 ! rvfir3 : output interstitial field (out,real(ngtot,3))
19 ! !DESCRIPTION:
20 ! Given two real vector fields, ${\bf f}_1$ and ${\bf f}_2$, defined over the
21 ! entire unit cell, this routine computes the local cross product
22 ! $$ {\bf f}_3({\bf r})\equiv{\bf f}_1({\bf r})\times{\bf f}_2({\bf r}). $$
23 !
24 ! !REVISION HISTORY:
25 ! Created February 2007 (JKD)
26 !EOP
27 !BOC
28 implicit none
29 ! arguments
30 real(8), intent(in) :: rvfmt1(npmtmax,natmtot,3),rvfir1(ngtot,3)
31 real(8), intent(in) :: rvfmt2(npmtmax,natmtot,3),rvfir2(ngtot,3)
32 real(8), intent(out) :: rvfmt3(npmtmax,natmtot,3),rvfir3(ngtot,3)
33 ! local variables
34 integer is,ias,nr,nri,ir,i
35 real(8) v1(3),v2(3),v3(3)
36 ! allocatable arrays
37 real(8), allocatable :: rvfmt4(:,:),rvfmt5(:,:)
38 !---------------------------!
39 ! muffin-tin region !
40 !---------------------------!
41 allocate(rvfmt4(npmtmax,3),rvfmt5(npmtmax,3))
42 do ias=1,natmtot
43  is=idxis(ias)
44  nr=nrmt(is)
45  nri=nrmti(is)
46  do i=1,3
47  call rbsht(nr,nri,rvfmt1(:,ias,i),rvfmt4(:,i))
48  call rbsht(nr,nri,rvfmt2(:,ias,i),rvfmt5(:,i))
49  end do
50  do i=1,npmt(is)
51  v1(:)=rvfmt4(i,:)
52  v2(:)=rvfmt5(i,:)
53  call r3cross(v1,v2,v3)
54  rvfmt4(i,:)=v3(:)
55  end do
56  do i=1,3
57  call rfsht(nr,nri,rvfmt4(:,i),rvfmt3(:,ias,i))
58  end do
59 end do
60 deallocate(rvfmt4,rvfmt5)
61 !-----------------------------!
62 ! interstitial region !
63 !-----------------------------!
64 do ir=1,ngtot
65  v1(:)=rvfir1(ir,:)
66  v2(:)=rvfir2(ir,:)
67  call r3cross(v1,v2,v3)
68  rvfir3(ir,:)=v3(:)
69 end do
70 end subroutine
71 !EOC
72 
integer, dimension(maxspecies) npmt
Definition: modmain.f90:213
subroutine rvfcross(rvfmt1, rvfir1, rvfmt2, rvfir2, rvfmt3, rvfir3)
Definition: rvfcross.f90:10
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition: rbsht.f90:7
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
subroutine rfsht(nr, nri, rfmt1, rfmt2)
Definition: rfsht.f90:7
pure subroutine r3cross(x, y, z)
Definition: r3cross.f90:10
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150