The Elk Code
Loading...
Searching...
No Matches
genwkpr0.f90
Go to the documentation of this file.
1
2
! Copyright (C) 2022 Leon Kerber, J. K. Dewhurst and S. Sharma.
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
genwkpr0
7
use
modmain
8
use
modtddft
9
use
moddftu
10
implicit none
11
integer
is,ia,ias,idu
12
integer
n,l,p,k,r,i
13
! determine the total number of tensor moments
14
n=0
15
do
idu=1,
ndftu
16
is=
isldu
(1,idu)
17
l=
isldu
(2,idu)
18
do
ia=1,
natoms
(is)
19
do
k=0,2*l
20
do
p=0,1
21
do
r=abs(k-p),k+p
22
n=n+1
23
end do
24
end do
25
end do
26
end do
27
end do
28
! allocate the t=0 tensor moment global array
29
if
(
allocated
(
wkpr0
))
deallocate
(
wkpr0
)
30
allocate
(
wkpr0
(-
lmmaxdm
:
lmmaxdm
,n))
31
i=0
32
do
idu=1,
ndftu
33
is=
isldu
(1,idu)
34
l=
isldu
(2,idu)
35
do
ia=1,
natoms
(is)
36
ias=
idxas
(ia,is)
37
do
k=0,2*l
38
do
p=0,1
39
do
r=abs(k-p),k+p
40
i=i+1
41
! decompose density matrix in 3-index tensor moment components
42
call
dmtotm3
(l,k,p,r,
lmmaxdm
,
dmatmt
(:,:,:,:,ias),
wkpr0
(:,i))
43
end do
44
end do
45
end do
46
end do
47
end do
48
end subroutine
49
dmtotm3
subroutine dmtotm3(l, k, p, r, ld, dm, wkpr)
Definition
dmtotm3.f90:10
genwkpr0
subroutine genwkpr0
Definition
genwkpr0.f90:7
moddftu
Definition
moddftu.f90:6
moddftu::ndftu
integer ndftu
Definition
moddftu.f90:38
moddftu::lmmaxdm
integer, parameter lmmaxdm
Definition
moddftu.f90:14
moddftu::wkpr0
real(8), dimension(:,:), allocatable wkpr0
Definition
moddftu.f90:86
moddftu::isldu
integer, dimension(2, maxdftu) isldu
Definition
moddftu.f90:40
moddftu::dmatmt
complex(8), dimension(:,:,:,:,:), allocatable dmatmt
Definition
moddftu.f90:16
modmain
Definition
modmain.f90:6
modmain::natoms
integer, dimension(maxspecies) natoms
Definition
modmain.f90:36
modmain::idxas
integer, dimension(maxatoms, maxspecies) idxas
Definition
modmain.f90:42
modtddft
Definition
modtddft.f90:6
genwkpr0.f90
Generated by
1.9.8