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