The Elk Code
initulr.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2018 T. Mueller, 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 initulr
7 use modmain
8 use modulr
9 use modomp
10 implicit none
11 ! local variables
12 integer ik0,ik,ist,jst
13 integer iq,jq,ifq,ig
14 integer n,i,i1,i2,i3,nthd
15 ! allocatable arrays
16 integer, allocatable :: idx(:)
17 real(8), allocatable :: jlgqr(:,:)
18 ! combined target array for long-range density and magnetisation
19 if (allocated(rhmgr)) deallocate(rhmgr)
21 if (spinpol) n=n*(1+ndmag)
22 allocate(rhmgr(n))
23 ! associate pointer arrays with target
24 rhormt(1:npcmtmax,1:natmtot,1:nqpt) => rhmgr(1:)
25 i=size(rhormt)+1
26 rhorir(1:ngtc,1:nqpt) => rhmgr(i:)
27 if (spinpol) then
28  i=i+size(rhorir)
29  magrmt(1:npcmtmax,1:natmtot,1:ndmag,1:nqpt) => rhmgr(i:)
30  i=i+size(magrmt)
31  magrir(1:ngtc,1:ndmag,1:nqpt) => rhmgr(i:)
32 end if
33 if (allocated(rhoqmt)) deallocate(rhoqmt)
34 allocate(rhoqmt(npcmtmax,natmtot,nfqrz))
35 if (allocated(rhoqir)) deallocate(rhoqir)
36 allocate(rhoqir(ngtc,nfqrz))
37 if (allocated(chgmtru)) deallocate(chgmtru)
38 allocate(chgmtru(natmtot,nqpt))
39 if (allocated(magqmt)) deallocate(magqmt)
40 if (allocated(magqir)) deallocate(magqir)
41 if (allocated(mommtru)) deallocate(mommtru)
42 if (allocated(momirru)) deallocate(momirru)
43 if (allocated(momtotru)) deallocate(momtotru)
44 if (spinpol) then
45  allocate(magqmt(npcmtmax,natmtot,ndmag,nfqrz))
46  allocate(magqir(ngtc,ndmag,nfqrz))
47  allocate(mommtru(ndmag,natmtot,nqpt))
48  allocate(momirru(ndmag,nqpt))
49  allocate(momtotru(ndmag,nqpt))
50 end if
51 ! allocate Q-dependent potential and magnetic field arrays
52 if (allocated(vclq)) deallocate(vclq)
53 allocate(vclq(nfqrz))
54 if (allocated(bfcq)) deallocate(bfcq)
55 if (allocated(bfcmtq)) deallocate(bfcmtq)
56 if (allocated(bdipq)) deallocate(bdipq)
57 if (spinpol) then
58  allocate(bfcq(ndmag,nfqrz))
59  allocate(bfcmtq(natmtot,ndmag,nfqrz))
60  if (tbdipu) allocate(bdipq(ndmag,nfqrz))
61 end if
62 ! combined target array for Kohn-Sham potential and magnetic field
63 if (allocated(vsbsq)) deallocate(vsbsq)
65 if (spinpol) n=n*(1+ndmag)
66 allocate(vsbsq(n))
67 ! zero the array
68 vsbsq(1:n)=0.d0
69 ! associate pointer arrays with target
70 vsqmt(1:npcmtmax,1:natmtot,1:nfqrz) => vsbsq(1:)
71 i=size(vsqmt)+1
72 vsqir(1:ngtot,1:nfqrz) => vsbsq(i:)
73 if (spinpol) then
74  i=i+size(vsqir)
75  bsqmt(1:npcmtmax,1:natmtot,1:ndmag,1:nfqrz) => vsbsq(i:)
76  i=i+size(bsqmt)
77  bsqir(1:ngtot,1:ndmag,1:nfqrz) => vsbsq(i:)
78 end if
79 ! generate the Coulomb Green's function in Q-space with small Q cut-off
80 if (allocated(gclq)) deallocate(gclq)
81 allocate(gclq(nqpt))
82 call gengclqu
83 ! G+Q-vector arrays
84 if (allocated(vgqc)) deallocate(vgqc)
85 allocate(vgqc(3,ngvec,nfqrz))
86 if (allocated(gqc)) deallocate(gqc)
87 allocate(gqc(ngvec,nfqrz))
88 if (allocated(ylmgq)) deallocate(ylmgq)
89 allocate(ylmgq(lmmaxo,ngvec,nfqrz))
90 if (allocated(sfacgq)) deallocate(sfacgq)
91 allocate(sfacgq(ngvec,natmtot,nfqrz))
92 if (allocated(gclgq)) deallocate(gclgq)
93 allocate(gclgq(ngvec,nfqrz))
94 if (allocated(jlgqrmt)) deallocate(jlgqrmt)
95 allocate(jlgqrmt(0:lnpsd,ngvec,nspecies,nfqrz))
96 if (allocated(expqmt)) deallocate(expqmt)
97 allocate(expqmt(npcmtmax,natmtot,nqpt))
98 call holdthd(nfqrz,nthd)
99 !$OMP PARALLEL DEFAULT(SHARED) &
100 !$OMP PRIVATE(jlgqr,iq,ig) &
101 !$OMP PRIVATE(i1,i2,i3,jq) &
102 !$OMP NUM_THREADS(nthd)
103 allocate(jlgqr(njcmax,nspecies))
104 !$OMP DO SCHEDULE(DYNAMIC)
105 do ifq=1,nfqrz
106  iq=iqrzf(ifq)
107  do ig=1,ngvec
108 ! determine the G+Q-vectors
109  vgqc(1:3,ig,ifq)=vgc(1:3,ig)+vqc(1:3,iq)
110 ! G+Q-vector length
111  gqc(ig,ifq)=sqrt(vgqc(1,ig,ifq)**2+vgqc(2,ig,ifq)**2+vgqc(3,ig,ifq)**2)
112 ! spherical harmonics for G+Q-vectors
113  call genylmv(.true.,lmaxo,vgqc(:,ig,ifq),ylmgq(:,ig,ifq))
114  end do
115 ! generate the spherical Bessel functions j_l(|G+Q|r)
116  call genjlgpr(1,gqc(1,ifq),jlgqr)
117 ! structure factors for G+Q-vectors
118  call gensfacgp(ngvec,vgqc(:,:,ifq),ngvec,sfacgq(:,:,ifq))
119 ! generate the Coulomb Green's function in G+Q-space
120  call gengclgq(.true.,iq,ngvec,gqc(:,ifq),gclgq(:,ifq))
121 ! compute the spherical Bessel functions j_l(|G+Q|R_mt)
122  call genjlgprmt(lnpsd,ngvec,gqc(:,ifq),ngvec,jlgqrmt(:,:,:,ifq))
123 ! generate phase factor functions exp(iQ⋅r) in each muffin-tin
124  call genexpmt(1,jlgqr,ylmgq(:,:,ifq),ngvec,sfacgq(:,:,ifq),expqmt(:,:,iq))
125 ! store the phase factor function for -Q
126  i1=-ivq(1,iq); i2=-ivq(2,iq); i3=-ivq(3,iq)
127  if ((i1 >= intq(1,1)).and.(i1 <= intq(2,1)).and. &
128  (i2 >= intq(1,2)).and.(i2 <= intq(2,2)).and. &
129  (i3 >= intq(1,3)).and.(i3 <= intq(2,3)).and.(ifq > 1)) then
130  jq=ivqiq(i1,i2,i3)
131  expqmt(:,:,jq)=conjg(expqmt(:,:,iq))
132  end if
133 end do
134 !$OMP END DO
135 deallocate(jlgqr)
136 !$OMP END PARALLEL
137 call freethd(nthd)
138 ! number of long-range states
140 ! allocate eigenvalue array
141 if (allocated(evalu)) deallocate(evalu)
142 allocate(evalu(nstulr,nkpt0))
143 ! allocate the occupation number array
144 if (allocated(occulr)) deallocate(occulr)
145 allocate(occulr(nstulr,nkpt0))
146 ! initialise the occupation numbers if required
147 if (any(task == [700,701,720,725])) then
148  allocate(idx(nstulr))
149  do ik0=1,nkpt0
150  ik=(ik0-1)*nkpa+1
151  call sortidx(nstulr,occsv(1,ik),idx)
152  do ist=1,nstulr
153  i=idx(nstulr-ist+1)-1
154  ik=(ik0-1)*nkpa+i/nstsv+1
155  jst=mod(i,nstsv)+1
156  occulr(ist,ik0)=occsv(jst,ik)
157  end do
158  end do
159  deallocate(idx)
160 end if
161 ! zero the timing variables
162 timemat=0.d0
163 timesv=0.d0
164 timerho=0.d0
165 timepot=0.d0
166 end subroutine
167 
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition: modulr.f90:84
complex(8), dimension(:,:,:), allocatable sfacgq
Definition: modulr.f90:44
real(8), dimension(:,:,:,:), allocatable jlgqrmt
Definition: modulr.f90:40
integer npcmtmax
Definition: modmain.f90:216
pure subroutine gensfacgp(ngp, vgpc, ld, sfacgp)
Definition: gensfacgp.f90:10
integer task
Definition: modmain.f90:1286
integer ngtot
Definition: modmain.f90:390
integer lmmaxo
Definition: modmain.f90:203
integer ngtc
Definition: modmain.f90:392
logical spinpol
Definition: modmain.f90:228
integer nqpt
Definition: modmain.f90:525
complex(8), dimension(:,:), pointer, contiguous vsqir
Definition: modulr.f90:84
integer nstulr
Definition: modulr.f90:95
integer ndmag
Definition: modmain.f90:238
integer, dimension(:,:), allocatable ivq
Definition: modmain.f90:529
complex(8), dimension(:,:), allocatable bdipq
Definition: modulr.f90:81
Definition: modomp.f90:6
integer nkpt0
Definition: modulr.f90:18
real(8) timemat
Definition: modmain.f90:1204
integer, dimension(:), allocatable iqrzf
Definition: modmain.f90:543
integer lmaxo
Definition: modmain.f90:201
complex(8), dimension(:,:,:,:), allocatable magqmt
Definition: modulr.f90:61
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition: modulr.f90:85
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
real(8) timerho
Definition: modmain.f90:1210
complex(8), dimension(:,:,:), allocatable magqir
Definition: modulr.f90:61
subroutine gengclqu
Definition: gengclqu.f90:7
real(8), dimension(:,:), allocatable chgmtru
Definition: modulr.f90:56
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
real(8), dimension(:,:), allocatable momirru
Definition: modulr.f90:58
real(8) timepot
Definition: modmain.f90:1212
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition: gengclgq.f90:7
integer nstsv
Definition: modmain.f90:880
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
complex(8), dimension(:), allocatable vclq
Definition: modulr.f90:66
real(8), dimension(:,:,:,:), pointer, contiguous magrmt
Definition: modulr.f90:54
real(8), dimension(:,:), allocatable occulr
Definition: modulr.f90:99
integer ngvec
Definition: modmain.f90:396
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:896
real(8), dimension(:), allocatable, target rhmgr
Definition: modulr.f90:52
complex(8), dimension(:,:,:), allocatable ylmgq
Definition: modulr.f90:42
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
real(8), dimension(:,:,:), allocatable vgqc
Definition: modulr.f90:34
integer, dimension(:,:,:), allocatable ivqiq
Definition: modmain.f90:531
real(8), dimension(:,:,:), pointer, contiguous rhormt
Definition: modulr.f90:53
real(8), dimension(:,:), pointer, contiguous rhorir
Definition: modulr.f90:53
subroutine genexpmt(ngp, jlgpr, ylmgp, ld, sfacgp, expmt)
Definition: genexpmt.f90:7
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition: modulr.f90:85
real(8), dimension(:,:), allocatable gqc
Definition: modulr.f90:36
integer nfqrz
Definition: modmain.f90:539
integer, dimension(2, 3) intq
Definition: modmain.f90:517
complex(8), dimension(:,:,:), allocatable bfcmtq
Definition: modulr.f90:73
subroutine genjlgpr(ngp, gpc, jlgpr)
Definition: genjlgpr.f90:7
real(8), dimension(:), allocatable gclq
Definition: modmain.f90:553
complex(8), dimension(:,:,:), allocatable rhoqmt
Definition: modulr.f90:60
integer lnpsd
Definition: modmain.f90:626
real(8), dimension(:,:,:), pointer, contiguous magrir
Definition: modulr.f90:54
real(8), dimension(:,:,:), allocatable mommtru
Definition: modulr.f90:58
complex(8), dimension(:,:), allocatable bfcq
Definition: modulr.f90:71
integer nspecies
Definition: modmain.f90:34
subroutine freethd(nthd)
Definition: modomp.f90:112
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
subroutine genjlgprmt(lmax, ngp, gpc, ld, jlgprmt)
Definition: genjlgprmt.f90:10
integer njcmax
Definition: modmain.f90:1160
real(8) timesv
Definition: modmain.f90:1208
integer natmtot
Definition: modmain.f90:40
real(8), dimension(:,:), allocatable gclgq
Definition: modulr.f90:38
Definition: modulr.f90:6
integer nkpa
Definition: modulr.f90:24
complex(8), dimension(:), allocatable, target vsbsq
Definition: modulr.f90:83
real(8), dimension(:,:), allocatable evalu
Definition: modulr.f90:97
complex(8), dimension(:,:,:), allocatable expqmt
Definition: modulr.f90:46
real(8), dimension(:,:), allocatable momtotru
Definition: modulr.f90:58
subroutine initulr
Definition: initulr.f90:7
logical tbdipu
Definition: modulr.f90:79
complex(8), dimension(:,:), allocatable rhoqir
Definition: modulr.f90:60