iexciting-0.9.224
[exciting.git] / src / vnlrhomt.f90
blob9baf041d2016ee892440a23885f79079906907a7
2 ! Copyright (C) 2002-2005 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.
6 !BOP
7 ! !ROUTINE: vnlrhomt
8 ! !INTERFACE:
9 subroutine vnlrhomt(tsh,is,wfmt1,wfmt2,zrhomt)
10 ! !USES:
11 use modmain
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! tsh : .true. if the density is to be in spherical harmonics (in,logical)
14 ! is : species number (in,integer)
15 ! wfmt1 : muffin-tin part of wavefunction 1 in spherical coordinates
16 ! (in,complex(lmmaxvr,nrcmtmax,natmtot,nspinor))
17 ! wfmt2 : muffin-tin part of wavefunction 2 in spherical coordinates
18 ! (in,complex(lmmaxvr,nrcmtmax,natmtot,nspinor))
19 ! zrhomt : muffin-tin charge density in spherical harmonics/coordinates
20 ! (out,complex(lmmaxvr,nrcmtmax))
21 ! !DESCRIPTION:
22 ! Calculates the complex overlap density in a single muffin-tin from two input
23 ! wavefunctions expressed in spherical coordinates. If {\tt tsh} is
24 ! {\tt .true.} then the output density is converted to a spherical harmonic
25 ! expansion. See routine {\tt vnlrho}.
27 ! !REVISION HISTORY:
28 ! Created November 2004 (Sharma)
29 !EOP
30 !BOC
31 implicit none
32 ! arguments
33 logical, intent(in) :: tsh
34 integer, intent(in) :: is
35 complex(8), intent(in) :: wfmt1(lmmaxvr,nrcmtmax)
36 complex(8), intent(in) :: wfmt2(lmmaxvr,nrcmtmax)
37 complex(8), intent(out) :: zrhomt(lmmaxvr,nrcmtmax)
38 ! local variables
39 integer irc
40 ! allocatable arrays
41 complex(8), allocatable :: zfmt(:,:)
42 if (tsh) then
43 ! output density in spherical harmonics
44 allocate(zfmt(lmmaxvr,nrcmtmax))
45 do irc=1,nrcmt(is)
46 zfmt(:,irc)=conjg(wfmt1(:,irc))*wfmt2(:,irc)
47 end do
48 call zgemm('N','N',lmmaxvr,nrcmt(is),lmmaxvr,zone,zfshtvr,lmmaxvr,zfmt, &
49 lmmaxvr,zzero,zrhomt,lmmaxvr)
50 deallocate(zfmt)
51 else
52 ! output density in spherical coordinates
53 do irc=1,nrcmt(is)
54 zrhomt(:,irc)=conjg(wfmt1(:,irc))*wfmt2(:,irc)
55 end do
56 end if
57 return
58 end subroutine
59 !EOC