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