The Elk Code
gendmftm.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 J. K. Dewhurst, S. Sharma, E. K. U. Gross and L. Nordstrom.
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 gendmftm
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
14 ! automatic arrays
15 real(8) wkpr(-lmmaxdm:lmmaxdm)
16 complex(8) dm(lmmaxdm,2,lmmaxdm,2)
17 ! allocate global array
18 if (allocated(dmftm)) deallocate(dmftm)
19 allocate(dmftm(lmmaxdm,2,lmmaxdm,2,natmtot))
20 ! zero the fixed tensor moment density matrices
21 dmftm(:,:,:,:,:)=0.d0
22 do i=1,ntmfix
23  is=itmfix(1,i)
24  if (is > nspecies) then
25  write(*,*)
26  write(*,'("Error(gendmftm): invalid species number : ",I8)') is
27  write(*,*)
28  stop
29  end if
30  ia=itmfix(2,i)
31  if (ia > natoms(is)) then
32  write(*,*)
33  write(*,'("Error(gendmftm): invalid atom number : ",I8)') ia
34  write(*,'(" for species ",I4)') is
35  write(*,*)
36  stop
37  end if
38  ias=idxas(ia,is)
39  l=itmfix(3,i)
40  if (l > lmaxdm) then
41  write(*,*)
42  write(*,'("Error(gendmftm): l > lmaxdm ",2I8)') l,lmaxdm
43  write(*,'(" for species ",I4," and atom ",I4)') is,ia
44  write(*,*)
45  stop
46  end if
47 ! generate the 3-index density matrix
48  k=itmfix(4,i)
49  p=itmfix(5,i)
50  r=itmfix(6,i)
51  t=itmfix(7,i)
52  if (abs(t) > lmmaxdm) then
53  write(*,*)
54  write(*,'("Error(gendmftm): invalid t : ",I8)') t
55  write(*,'(" for tensor moment entry ",I3)') i
56  write(*,*)
57  stop
58  end if
59 ! scale factor for conventional normalisation
60  t0=sqrt(dble((2*l+1)*2))
61  wkpr(:)=0.d0
62  wkpr(t)=wkprfix(i)/t0
63  call tm3todm(l,k,p,r,lmmaxdm,wkpr,dm)
64  dmftm(:,:,:,:,ias)=dmftm(:,:,:,:,ias)+dm(:,:,:,:)
65 end do
66 end subroutine
67 
integer, dimension(:,:), allocatable itmfix
Definition: moddftu.f90:76
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
complex(8), dimension(:,:,:,:,:), allocatable dmftm
Definition: moddftu.f90:80
real(8), dimension(:), allocatable wkprfix
Definition: moddftu.f90:78
subroutine gendmftm
Definition: gendmftm.f90:7
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer nspecies
Definition: modmain.f90:34
integer natmtot
Definition: modmain.f90:40
subroutine tm3todm(l, k, p, r, ld, wkpr, dm)
Definition: tm3todm.f90:10
integer, parameter lmaxdm
Definition: moddftu.f90:13
integer ntmfix
Definition: moddftu.f90:72