The Elk Code
 
Loading...
Searching...
No Matches
gndstulr.f90
Go to the documentation of this file.
1
2! Copyright (C) 2017 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 gndstulr
7use modmain
8use modulr
9use modmpi
10use modomp
11implicit none
12! local variables
13logical twrite
14integer ik0,ir,lp,nthd
15integer nmix,nwork,n
16real(8) dv
17character(64) str
18! allocatable arrays
19integer(omp_lock_kind), allocatable :: lock(:)
20real(8), allocatable :: work(:)
21complex(8), allocatable :: evecu(:,:)
22if (xctype(1) < 0) then
23 write(*,*)
24 write(*,'("Error(gndstulr): ultra long-range does not work with OEP")')
25 write(*,*)
26 stop
27end if
28if (spincore) then
29 write(*,*)
30 write(*,'("Error(gndstulr): ultra long-range does not work with &
31 &spin-polarised cores")')
32 write(*,*)
33 stop
34end if
35! no k-point reduction
37reducek=0
38! initialise global variables
39call init0
40call init1
41! write the κ-points to file
42call writekpa
43! write the k+κ-points to file
44call writekpts
45! read the regular Kohn-Sham potential from file
46call readstate
47! generate the first- and second-variational eigenvectors and eigenvalues for
48! the k+κ-point set
49call genvsig
50call gencore
51call linengy
52call genapwlofr
53call gensocfr
54call genevfsv
55call occupy
56! initialise the ultra long-range variables
57call initulr
58if (task == 700) then
59! initialise the long-range Kohn-Sham potential and magnetic field
60 call potuinit
61else
62! read in the potential and density from STATE_ULR.OUT
63 call readstulr
64end if
65! initialise the external Coulomb potential
66call vclqinit
67! initialise the external magnetic field if required
68if (spinpol) call bfcqinit
69! size of mixing vector (complex array)
70nmix=2*size(vsbsq)
71! determine the size of the mixer work array
72nwork=-1
73call mixerifc(mixtype,nmix,vsbsq,dv,nwork,work)
74allocate(work(nwork))
75! initialise the mixer
76iscl=0
77call mixerifc(mixtype,nmix,vsbsq,dv,nwork,work)
78! initialise the OpenMP locks
79allocate(lock(nqpt))
80do ir=1,nqpt
81 call omp_init_lock(lock(ir))
82end do
83! set last self-consistent loop flag
84tlast=.false.
85! begin the self-consistent loop
86if (mp_mpi) then
87! open ULR_INFO.OUT file
88 open(60,file='ULR_INFO.OUT',form='FORMATTED')
89! open RMSDVS.OUT
90 open(65,file='RMSDVS.OUT',form='FORMATTED')
91 call writeinfou(60)
92 call writebox(60,"Self-consistent loop started")
93end if
94do iscl=1,maxscl
95 if (mp_mpi) then
96 write(str,'("Loop number : ",I0)') iscl
97 call writebox(60,trim(str))
98 end if
99 if (iscl >= maxscl) then
100 if (mp_mpi) then
101 write(60,*)
102 write(60,'("Reached self-consistent loops maximum")')
103 end if
104 write(*,*)
105 write(*,'("Warning(gndstulr): failed to reach self-consistency after ",I4,&
106 &" loops")') iscl
107 tlast=.true.
108 end if
109! reset the OpenMP thread variables
110 call omp_reset
111! apply required local operations to the potential and magnetic field
112 call vblocalu
113! zero the density and magnetisation
114 rhormt(:,:,:)=0.d0
115 rhorir(:,:)=0.d0
116 if (spinpol) then
117 magrmt(:,:,:,:)=0.d0
118 magrir(:,:,:)=0.d0
119 end if
120! loop over original k-points
121 call holdthd(nkpt0/np_mpi,nthd)
122!$OMP PARALLEL DEFAULT(SHARED) &
123!$OMP PRIVATE(evecu) &
124!$OMP NUM_THREADS(nthd)
125 allocate(evecu(nstulr,nstulr))
126!$OMP DO SCHEDULE(DYNAMIC)
127 do ik0=1,nkpt0
128! distribute among MPI processes
129 if (mod(ik0-1,np_mpi) /= lp_mpi) cycle
130! solve the ultra long-range eigenvalue equation
131 call eveqnulr(ik0,evecu)
132! write the eigenvalues to file
133 call putevalu(ik0)
134! add to the density, magnetisation and current
135 call rhomaguk(ik0,lock,evecu)
136 end do
137!$OMP END DO
138 deallocate(evecu)
139!$OMP END PARALLEL
140 call freethd(nthd)
141 if (np_mpi > 1) then
142! broadcast eigenvalue array to every process
143 do ik0=1,nkpt0
144 lp=mod(ik0-1,np_mpi)
145 call mpi_bcast(evalu(:,ik0),nstulr,mpi_double_precision,lp,mpicom,ierror)
146 end do
147! add densities from each process and redistribute
149 call mpi_allreduce(mpi_in_place,rhormt,n,mpi_double_precision,mpi_sum, &
151 n=ngtc*nqpt
152 call mpi_allreduce(mpi_in_place,rhorir,n,mpi_double_precision,mpi_sum, &
154 if (spinpol) then
156 call mpi_allreduce(mpi_in_place,magrmt,n,mpi_double_precision,mpi_sum, &
158 n=ngtc*ndmag*nqpt
159 call mpi_allreduce(mpi_in_place,magrir,n,mpi_double_precision,mpi_sum, &
161 end if
162 end if
163! find the occupation numbers and Fermi energy
164 call occupyulr
165! synchronise MPI processes
166 call mpi_barrier(mpicom,ierror)
167! add the core density
168 call rhocoreu
169! perform partial Fourier transform to Q-space
170 call rhomagq
171! determine the muffin-tin and interstitial charges and moments
172 call chargeu
173 call momentu
174! compute the ultra long-range Kohn-Sham potential
175 call potksu
176! mix the old potential and field with the new
177 call mixerifc(mixtype,nmix,vsbsq,dv,nwork,work)
178! multiply the RMS change in potential by the number of Q-points
179 dv=dv*dble(nfqrz)
180! calculate and add the fixed spin moment effective field (after mixing)
181 call bfieldfsm
182 call addbfsmu
183! compute the energy components
184 call energyulr
185 if (mp_mpi) then
186! write eigenvalues to file
187 call writeevalu
188! output energy components
189 call writeengyu(60)
190! output charges
191 call writechg(60)
192! write muffin-tin charges for each R-vector
193 call writechgrmt
194 if (spinpol) then
195! output moments
196 call writemom(60)
197! write muffin-tin, interstitial and total moments for each R-vector
198 call writemomru
199! write the moments for each Q-vector
200 call writemomqu
201 end if
202! output effective fields for fixed spin moment calculations
203 if (fsmtype /= 0) call writefsm(60)
204! check for existence of the WRITE file
205 call checkwrite(twrite)
206! check self-consistent loop number modulo nwrite
207 if (nwrite >= 1) then
208 if (mod(iscl,nwrite) == 0) twrite=.true.
209 end if
210! write STATE_ULR.OUT file if required
211 if (twrite) then
212 call writestulr
213 write(60,*)
214 write(60,'("Wrote STATE_ULR.OUT")')
215 end if
216 end if
217! exit self-consistent loop if required
218 if (tlast) goto 10
219! check for convergence
220 if (iscl >= 2) then
221 if (mp_mpi) then
222 write(60,*)
223 write(60,'("RMS change in Kohn-Sham potential (target) : ",G18.10," (",&
224 &G18.10,")")') dv,epspot
225 flush(60)
226 write(65,'(G18.10)') dv
227 flush(65)
228 end if
229 if (dv < epspot) then
230 if (mp_mpi) then
231 write(60,*)
232 write(60,'("Convergence targets achieved")')
233 end if
234 tlast=.true.
235 end if
236 end if
237! check for STOP file
238 call checkstop
239 if (tstop) tlast=.true.
240! broadcast tlast from master process to all other processes
241 call mpi_bcast(tlast,1,mpi_logical,0,mpicom,ierror)
242! reset the OpenMP thread variables
243 call omp_reset
244end do
24510 continue
246if (mp_mpi) then
247! output timing information
248 write(60,*)
249 write(60,'("Timings (CPU seconds) :")')
250 write(60,'(" Hamiltonian matrix set up",T40,": ",F12.2)') timemat
251 write(60,'(" eigenvalue equation",T40,": ",F12.2)') timesv
252 write(60,'(" charge density calculation",T40,": ",F12.2)') timerho
253 write(60,'(" potential calculation",T40,": ",F12.2)') timepot
254 call writebox(60,"Self-consistent loop stopped")
255 if (maxscl > 1) then
256 call writestulr
257 write(60,*)
258 write(60,'("Wrote STATE_ULR.OUT")')
259 end if
260! close the ULR_INFO.OUT file
261 close(60)
262! close the RMSDVS.OUT file
263 close(65)
264end if
265! destroy the OpenMP locks
266do ir=1,nqpt
267 call omp_destroy_lock(lock(ir))
268end do
269deallocate(lock,work)
270! restore original parameters
272! synchronise MPI processes
273call mpi_barrier(mpicom,ierror)
274end subroutine
275
subroutine addbfsmu
Definition addbfsmu.f90:7
subroutine bfcqinit
Definition bfcqinit.f90:7
subroutine bfieldfsm
Definition bfieldfsm.f90:10
subroutine chargeu
Definition chargeu.f90:7
subroutine checkstop
Definition checkstop.f90:7
subroutine checkwrite(twrite)
Definition checkwrite.f90:7
subroutine energyulr
Definition energyulr.f90:7
subroutine eveqnulr(ik0, evecu)
Definition eveqnulr.f90:7
subroutine genapwlofr
Definition genapwlofr.f90:7
subroutine gencore
Definition gencore.f90:10
subroutine genevfsv
Definition genevfsv.f90:7
subroutine gensocfr
Definition gensocfr.f90:7
subroutine genvsig
Definition genvsig.f90:10
subroutine gndstulr
Definition gndstulr.f90:7
subroutine init0
Definition init0.f90:10
subroutine init1
Definition init1.f90:10
subroutine initulr
Definition initulr.f90:7
subroutine linengy
Definition linengy.f90:10
subroutine mixerifc(mtype, n, v, dv, nwork, work)
Definition mixerifc.f90:7
subroutine momentu
Definition momentu.f90:7
real(8) epspot
Definition modmain.f90:1058
integer nwrite
Definition modmain.f90:1056
integer nfqrz
Definition modmain.f90:539
real(8) timepot
Definition modmain.f90:1222
logical spinpol
Definition modmain.f90:228
logical spincore
Definition modmain.f90:940
integer fsmtype
Definition modmain.f90:251
integer, dimension(3) xctype
Definition modmain.f90:588
real(8) timerho
Definition modmain.f90:1220
integer mixtype
Definition modmain.f90:695
real(8) timemat
Definition modmain.f90:1214
integer ngtc
Definition modmain.f90:392
integer nqpt
Definition modmain.f90:525
logical tstop
Definition modmain.f90:1052
integer natmtot
Definition modmain.f90:40
integer npcmtmax
Definition modmain.f90:216
logical tlast
Definition modmain.f90:1050
integer task
Definition modmain.f90:1298
integer maxscl
Definition modmain.f90:1046
integer reducek
Definition modmain.f90:455
real(8) timesv
Definition modmain.f90:1218
integer iscl
Definition modmain.f90:1048
integer ndmag
Definition modmain.f90:238
integer reducek0
Definition modmain.f90:455
integer lp_mpi
Definition modmpi.f90:15
integer ierror
Definition modmpi.f90:19
integer mpicom
Definition modmpi.f90:11
integer np_mpi
Definition modmpi.f90:13
logical mp_mpi
Definition modmpi.f90:17
subroutine omp_reset
Definition modomp.f90:71
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
integer nstulr
Definition modulr.f90:94
real(8), dimension(:,:,:), allocatable rhormt
Definition modulr.f90:52
real(8), dimension(:,:,:,:), allocatable magrmt
Definition modulr.f90:53
real(8), dimension(:,:), allocatable evalu
Definition modulr.f90:96
complex(8), dimension(:), allocatable, target vsbsq
Definition modulr.f90:82
real(8), dimension(:,:,:), allocatable magrir
Definition modulr.f90:53
real(8), dimension(:,:), allocatable rhorir
Definition modulr.f90:52
integer nkpt0
Definition modulr.f90:18
subroutine occupy
Definition occupy.f90:10
subroutine occupyulr
Definition occupyulr.f90:7
subroutine potksu
Definition potksu.f90:7
subroutine potuinit
Definition potuinit.f90:7
subroutine putevalu(ik0)
Definition putevalu.f90:7
subroutine readstate
Definition readstate.f90:10
subroutine readstulr
Definition readstulr.f90:7
subroutine rhocoreu
Definition rhocoreu.f90:7
subroutine rhomagq
Definition rhomagq.f90:7
subroutine rhomaguk(ik0, lock, evecu)
Definition rhomaguk.f90:7
subroutine vblocalu
Definition vblocalu.f90:7
subroutine vclqinit
Definition vclqinit.f90:7
subroutine writebox(fnum, str)
Definition writebox.f90:7
subroutine writechg(fnum)
Definition writechg.f90:7
subroutine writechgrmt
subroutine writeengyu(fnum)
Definition writeengyu.f90:7
subroutine writeevalu
Definition writeevalu.f90:7
subroutine writefsm(fnum)
Definition writefsm.f90:7
subroutine writeinfou(fnum)
Definition writeinfou.f90:7
subroutine writekpa
Definition writekpa.f90:7
subroutine writekpts
Definition writekpts.f90:10
subroutine writemom(fnum)
Definition writemom.f90:7
subroutine writemomqu
Definition writemomqu.f90:7
subroutine writemomru
Definition writemomru.f90:7
subroutine writestulr
Definition writestulr.f90:7