The Elk Code
vmatmtftm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 L. Nordstrom, 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 
6 subroutine vmatmtftm
7 use modmain
8 use moddftu
9 implicit none
10 ! local variables
11 integer is,ia,ias,i
12 integer l,k,p,r,t
13 real(8) t0,t1
14 ! automatic arrays
15 real(8) wkpr(-lmmaxdm:lmmaxdm)
16 complex(8) dm(lmmaxdm,2,lmmaxdm,2)
17 ! loop over fixed tensor moment entries
18 do i=1,ntmfix
19  is=itmfix(1,i)
20  ia=itmfix(2,i)
21  ias=idxas(ia,is)
22  l=itmfix(3,i)
23  k=itmfix(4,i)
24  p=itmfix(5,i)
25  r=itmfix(6,i)
26  t=itmfix(7,i)
27 ! decompose density matrix in 3-index tensor moment components
28  call dmtotm3(l,k,p,r,lmmaxdm,dmatmt(:,:,:,:,ias),wkpr)
29 ! scale factor for conventional normalisation
30  t0=sqrt(dble((2*l+1)*2))
31 ! take difference between current and target moment
32  t1=wkpr(t)-wkprfix(i)/t0
33  wkpr(:)=0.d0
34  wkpr(t)=tauftm*t1
35 ! compute new density matrix
36  call tm3todm(l,k,p,r,lmmaxdm,wkpr,dm)
37 ! add to global FTM potential matrix
38  vmftm(:,:,:,:,ias)=vmftm(:,:,:,:,ias)+dm(:,:,:,:)
39 end do
40 ! add to muffin-tin potential matrix (fix by Leon Kerber)
41 do i=1,ntmfix
42  is=itmfix(1,i)
43  ia=itmfix(2,i)
44  ias=idxas(ia,is)
45  vmatmt(:,:,:,:,ias)=vmatmt(:,:,:,:,ias)+vmftm(:,:,:,:,ias)
46 end do
47 end subroutine
48 
integer, dimension(:,:), allocatable itmfix
Definition: moddftu.f90:76
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
real(8), dimension(:), allocatable wkprfix
Definition: moddftu.f90:78
complex(8), dimension(:,:,:,:,:), allocatable dmatmt
Definition: moddftu.f90:16
real(8) tauftm
Definition: moddftu.f90:84
complex(8), dimension(:,:,:,:,:), allocatable vmatmt
Definition: moddftu.f90:20
subroutine dmtotm3(l, k, p, r, ld, dm, wkpr)
Definition: dmtotm3.f90:10
complex(8), dimension(:,:,:,:,:), allocatable vmftm
Definition: moddftu.f90:82
subroutine vmatmtftm
Definition: vmatmtftm.f90:7
subroutine tm3todm(l, k, p, r, ld, wkpr, dm)
Definition: tm3todm.f90:10
integer ntmfix
Definition: moddftu.f90:72