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,ifq,ig,n,i,nthd
14 ! allocatable arrays
15 integer, allocatable :: idx(:)
16 ! combined target array for long-range density and magnetisation
17 if (allocated(rhmgr)) deallocate(rhmgr)
19 if (spinpol) n=n*(1+ndmag)
20 allocate(rhmgr(n))
21 ! associate pointer arrays with target
22 rhormt(1:npcmtmax,1:natmtot,1:nqpt) => rhmgr(1:)
23 i=size(rhormt)+1
24 rhorir(1:ngtc,1:nqpt) => rhmgr(i:)
25 if (spinpol) then
26  i=i+size(rhorir)
27  magrmt(1:npcmtmax,1:natmtot,1:ndmag,1:nqpt) => rhmgr(i:)
28  i=i+size(magrmt)
29  magrir(1:ngtc,1:ndmag,1:nqpt) => rhmgr(i:)
30 end if
31 if (allocated(rhoqmt)) deallocate(rhoqmt)
32 allocate(rhoqmt(npcmtmax,natmtot,nfqrz))
33 if (allocated(rhoqir)) deallocate(rhoqir)
34 allocate(rhoqir(ngtc,nfqrz))
35 if (allocated(chgmtru)) deallocate(chgmtru)
36 allocate(chgmtru(natmtot,nqpt))
37 if (allocated(magqmt)) deallocate(magqmt)
38 if (allocated(magqir)) deallocate(magqir)
39 if (allocated(mommtru)) deallocate(mommtru)
40 if (allocated(momirru)) deallocate(momirru)
41 if (allocated(momtotru)) deallocate(momtotru)
42 if (spinpol) then
43  allocate(magqmt(npcmtmax,natmtot,ndmag,nfqrz))
44  allocate(magqir(ngtc,ndmag,nfqrz))
45  allocate(mommtru(ndmag,natmtot,nqpt))
46  allocate(momirru(ndmag,nqpt))
47  allocate(momtotru(ndmag,nqpt))
48 end if
49 ! allocate Q-dependent potential and magnetic field arrays
50 if (allocated(vclq)) deallocate(vclq)
51 allocate(vclq(nfqrz))
52 if (allocated(bfcq)) deallocate(bfcq)
53 if (allocated(bfcmtq)) deallocate(bfcmtq)
54 if (allocated(bdipq)) deallocate(bdipq)
55 if (spinpol) then
56  allocate(bfcq(ndmag,nfqrz))
57  allocate(bfcmtq(natmtot,ndmag,nfqrz))
58  if (tbdipu) allocate(bdipq(ndmag,nfqrz))
59 end if
60 ! combined target array for Kohn-Sham potential and magnetic field
61 if (allocated(vsbsq)) deallocate(vsbsq)
63 if (spinpol) n=n*(1+ndmag)
64 allocate(vsbsq(n))
65 ! zero the array
66 vsbsq(1:n)=0.d0
67 ! associate pointer arrays with target
68 vsqmt(1:npcmtmax,1:natmtot,1:nfqrz) => vsbsq(1:)
69 i=size(vsqmt)+1
70 vsqir(1:ngtot,1:nfqrz) => vsbsq(i:)
71 if (spinpol) then
72  i=i+size(vsqir)
73  bsqmt(1:npcmtmax,1:natmtot,1:ndmag,1:nfqrz) => vsbsq(i:)
74  i=i+size(bsqmt)
75  bsqir(1:ngtot,1:ndmag,1:nfqrz) => vsbsq(i:)
76 end if
77 ! generate the Coulomb Green's function in Q-space with small Q cut-off
78 if (allocated(gclq)) deallocate(gclq)
79 allocate(gclq(nqpt))
80 call gengclqu
81 ! G+Q-vector arrays
82 if (allocated(vgqc)) deallocate(vgqc)
83 allocate(vgqc(3,ngvec,nfqrz))
84 if (allocated(gqc)) deallocate(gqc)
85 allocate(gqc(ngvec,nfqrz))
86 if (allocated(ylmgq)) deallocate(ylmgq)
87 allocate(ylmgq(lmmaxo,ngvec,nfqrz))
88 if (allocated(sfacgq)) deallocate(sfacgq)
89 allocate(sfacgq(ngvec,natmtot,nfqrz))
90 if (allocated(gclgq)) deallocate(gclgq)
91 allocate(gclgq(ngvec,nfqrz))
92 if (allocated(jlgqrmt)) deallocate(jlgqrmt)
93 allocate(jlgqrmt(0:lnpsd,ngvec,nspecies,nfqrz))
94 call holdthd(nfqrz,nthd)
95 !$OMP PARALLEL DO DEFAULT(SHARED) &
96 !$OMP PRIVATE(iq,ig) &
97 !$OMP SCHEDULE(DYNAMIC) &
98 !$OMP NUM_THREADS(nthd)
99 do ifq=1,nfqrz
100  iq=iqrzf(ifq)
101  do ig=1,ngvec
102 ! determine the G+Q-vectors
103  vgqc(1:3,ig,ifq)=vgc(1:3,ig)+vqc(1:3,iq)
104 ! G+Q-vector length
105  gqc(ig,ifq)=sqrt(vgqc(1,ig,ifq)**2+vgqc(2,ig,ifq)**2+vgqc(3,ig,ifq)**2)
106 ! spherical harmonics for G+Q-vectors
107  call genylmv(.true.,lmaxo,vgqc(:,ig,ifq),ylmgq(:,ig,ifq))
108  end do
109 ! structure factors for G+Q-vectors
110  call gensfacgp(ngvec,vgqc(:,:,ifq),ngvec,sfacgq(:,:,ifq))
111 ! generate the Coulomb Green's function in G+Q-space
112  call gengclgq(.true.,iq,ngvec,gqc(:,ifq),gclgq(:,ifq))
113 ! compute the spherical Bessel functions j_l(|G+Q|R_mt)
114  call genjlgprmt(lnpsd,ngvec,gqc(:,ifq),ngvec,jlgqrmt(:,:,:,ifq))
115 end do
116 !$OMP END PARALLEL DO
117 call freethd(nthd)
118 ! number of long-range states
120 ! allocate eigenvalue array
121 if (allocated(evalu)) deallocate(evalu)
122 allocate(evalu(nstulr,nkpt0))
123 ! allocate the occupation number array
124 if (allocated(occulr)) deallocate(occulr)
125 allocate(occulr(nstulr,nkpt0))
126 ! initialise the occupation numbers if required
127 if (any(task == [700,701,720,725])) then
128  allocate(idx(nstulr))
129  do ik0=1,nkpt0
130  ik=(ik0-1)*nkpa+1
131  call sortidx(nstulr,occsv(1,ik),idx)
132  do ist=1,nstulr
133  i=idx(nstulr-ist+1)-1
134  ik=(ik0-1)*nkpa+i/nstsv+1
135  jst=mod(i,nstsv)+1
136  occulr(ist,ik0)=occsv(jst,ik)
137  end do
138  end do
139  deallocate(idx)
140 end if
141 ! zero the timing variables
142 timemat=0.d0
143 timesv=0.d0
144 timerho=0.d0
145 timepot=0.d0
146 end subroutine
147 
complex(8), dimension(:,:,:), pointer, contiguous vsqmt
Definition: modulr.f90:82
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:1297
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:82
integer nstulr
Definition: modulr.f90:93
integer ndmag
Definition: modmain.f90:238
complex(8), dimension(:,:), allocatable bdipq
Definition: modulr.f90:79
Definition: modomp.f90:6
integer nkpt0
Definition: modulr.f90:18
real(8) timemat
Definition: modmain.f90:1215
integer, dimension(:), allocatable iqrzf
Definition: modmain.f90:543
integer lmaxo
Definition: modmain.f90:201
complex(8), dimension(:,:,:,:), allocatable magqmt
Definition: modulr.f90:59
complex(8), dimension(:,:,:), pointer, contiguous bsqir
Definition: modulr.f90:83
pure subroutine genylmv(t4pil, lmax, v, ylm)
Definition: genylmv.f90:10
real(8) timerho
Definition: modmain.f90:1221
complex(8), dimension(:,:,:), allocatable magqir
Definition: modulr.f90:59
subroutine gengclqu
Definition: gengclqu.f90:7
real(8), dimension(:,:), allocatable chgmtru
Definition: modulr.f90:54
real(8), dimension(:,:), allocatable vgc
Definition: modmain.f90:420
real(8), dimension(:,:), allocatable momirru
Definition: modulr.f90:56
real(8) timepot
Definition: modmain.f90:1223
pure subroutine gengclgq(treg, iq, ngq, gqc, gclgq)
Definition: gengclgq.f90:7
integer nstsv
Definition: modmain.f90:887
real(8), dimension(:,:), allocatable vqc
Definition: modmain.f90:547
complex(8), dimension(:), allocatable vclq
Definition: modulr.f90:64
real(8), dimension(:,:,:,:), pointer, contiguous magrmt
Definition: modulr.f90:52
real(8), dimension(:,:), allocatable occulr
Definition: modulr.f90:97
integer ngvec
Definition: modmain.f90:396
real(8), dimension(:,:), allocatable occsv
Definition: modmain.f90:903
real(8), dimension(:), allocatable, target rhmgr
Definition: modulr.f90:50
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
real(8), dimension(:,:,:), pointer, contiguous rhormt
Definition: modulr.f90:51
real(8), dimension(:,:), pointer, contiguous rhorir
Definition: modulr.f90:51
complex(8), dimension(:,:,:,:), pointer, contiguous bsqmt
Definition: modulr.f90:83
real(8), dimension(:,:), allocatable gqc
Definition: modulr.f90:36
integer nfqrz
Definition: modmain.f90:539
complex(8), dimension(:,:,:), allocatable bfcmtq
Definition: modulr.f90:71
real(8), dimension(:), allocatable gclq
Definition: modmain.f90:553
complex(8), dimension(:,:,:), allocatable rhoqmt
Definition: modulr.f90:58
integer lnpsd
Definition: modmain.f90:628
real(8), dimension(:,:,:), pointer, contiguous magrir
Definition: modulr.f90:52
real(8), dimension(:,:,:), allocatable mommtru
Definition: modulr.f90:56
complex(8), dimension(:,:), allocatable bfcq
Definition: modulr.f90:69
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
real(8) timesv
Definition: modmain.f90:1219
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:81
real(8), dimension(:,:), allocatable evalu
Definition: modulr.f90:95
real(8), dimension(:,:), allocatable momtotru
Definition: modulr.f90:56
subroutine initulr
Definition: initulr.f90:7
logical tbdipu
Definition: modulr.f90:77
complex(8), dimension(:,:), allocatable rhoqir
Definition: modulr.f90:58