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.
9 subroutine vnlrho(tsh
,wfmt1
,wfmt2
,wfir1
,wfir2
,zrhomt
,zrhoir
)
12 ! !INPUT/OUTPUT PARAMETERS:
13 ! tsh : .true. if the muffin-tin density is to be in spherical harmonics
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 ! wfir1 : interstitial wavefunction 1 (in,complex(ngrtot))
20 ! wfir2 : interstitial wavefunction 2 (in,complex(ngrtot))
21 ! zrhomt : muffin-tin charge density in spherical harmonics/coordinates
22 ! (out,complex(lmmaxvr,nrcmtmax,natmtot))
23 ! zrhoir : interstitial charge density (out,complex(ngrtot))
25 ! Calculates the complex overlap charge density from two input wavefunctions:
26 ! $$ \rho({\bf r})\equiv\Psi_1^{\dag}({\bf r})\cdot\Psi_2({\bf r}). $$
27 ! Note that the muffin-tin wavefunctions are provided in spherical coordinates
28 ! and the returned density is either in terms of spherical harmonic
29 ! coefficients or spherical coordinates when {\tt tsh} is {\tt .true.} or
30 ! {\tt .false.}, respectively. See also the routine {\tt vnlrhomt}.
33 ! Created November 2004 (Sharma)
38 logical, intent(in
) :: tsh
39 complex(8), intent(in
) :: wfmt1(lmmaxvr
,nrcmtmax
,natmtot
,nspinor
)
40 complex(8), intent(in
) :: wfmt2(lmmaxvr
,nrcmtmax
,natmtot
,nspinor
)
41 complex(8), intent(in
) :: wfir1(ngrtot
,nspinor
)
42 complex(8), intent(in
) :: wfir2(ngrtot
,nspinor
)
43 complex(8), intent(out
) :: zrhomt(lmmaxvr
,nrcmtmax
,natmtot
)
44 complex(8), intent(out
) :: zrhoir(ngrtot
)
46 integer is
,ia
,ias
,nrc
,ir
48 complex(8), allocatable
:: zfmt(:,:)
49 if (spinpol
) allocate(zfmt(lmmaxvr
,nrcmtmax
))
55 call vnlrhomt(tsh
,is
,wfmt1(:,:,ias
,1),wfmt2(:,:,ias
,1),zrhomt(:,:,ias
))
58 call vnlrhomt(tsh
,is
,wfmt1(:,:,ias
,2),wfmt2(:,:,ias
,2),zfmt
)
59 zrhomt(:,1:nrc
,ias
)=zrhomt(:,1:nrc
,ias
)+zfmt(:,1:nrc
)
67 zrhoir(ir
)=conjg(wfir1(ir
,1))*wfir2(ir
,1)+conjg(wfir1(ir
,2))*wfir2(ir
,2)
72 zrhoir(ir
)=conjg(wfir1(ir
,1))*wfir2(ir
,1)
75 if (spinpol
) deallocate(zfmt
)