The Elk Code
 
Loading...
Searching...
No Matches
oepmain.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6subroutine oepmain
7use modmain
8use modmpi
9use modomp
10implicit none
11! local variables
12integer ik,idm,is,ias
13integer nrc,nrci,np,npc
14integer n,nthd,it
15real(8) resp,t1
16! automatic arrays
17real(8) rfmt2(npcmtmax)
18! allocatable arrays
19real(8), allocatable :: rfmt1(:,:),rfir(:)
20real(8), allocatable :: rvfmt(:,:,:),rvfir(:,:)
21complex(8), allocatable :: vclcv(:,:,:,:),vclvv(:,:,:)
22! external functions
23real(8), external :: rfint,rfinpc
24if (iscl < 1) return
25! calculate Coulomb matrix elements
26allocate(vclcv(ncrmax,natmtot,nstsv,nkpt),vclvv(nstsv,nstsv,nkpt))
27call oepvcl(vclcv,vclvv)
28! allocate local arrays
29allocate(rfmt1(npmtmax,natmtot),rfir(ngtot))
30if (spinpol) then
31 allocate(rvfmt(npmtmax,natmtot,ndmag),rvfir(ngtot,ndmag))
32end if
33! zero initial exchange potential and magnetic field
34vxmt(:,:)=0.d0
35vxir(:)=0.d0
36if (spinpol) then
37 bxmt(:,:,:)=0.d0
38 bxir(:,:)=0.d0
39end if
40!------------------------------!
41! start iteration loop !
42!------------------------------!
43do it=1,maxitoep
44 if (mp_mpi.and.(mod(it,10) == 0)) then
45 write(*,'("Info(oepmain): done ",I4," iterations of ",I4)') it,maxitoep
46 end if
47! zero the residuals
48 dvxmt(:,:)=0.d0
49 dvxir(:)=0.d0
50 if (spinpol) then
51 dbxmt(:,:,:)=0.d0
52 dbxir(:,:)=0.d0
53 end if
54! calculate the k-dependent residuals
55 call holdthd(nkpt/np_mpi,nthd)
56!$OMP PARALLEL DO DEFAULT(SHARED) &
57!$OMP SCHEDULE(DYNAMIC) &
58!$OMP NUM_THREADS(nthd)
59 do ik=1,nkpt
60! distribute among MPI processes
61 if (mod(ik-1,np_mpi) /= lp_mpi) cycle
62 call oepresk(ik,vclcv,vclvv)
63 end do
64!$OMP END PARALLEL DO
65 call freethd(nthd)
66! add residuals from each process and redistribute
67 if (np_mpi > 1) then
68 n=npcmtmax*natmtot
69 call mpi_allreduce(mpi_in_place,dvxmt,n,mpi_double_precision,mpi_sum, &
71 call mpi_allreduce(mpi_in_place,dvxir,ngtot,mpi_double_precision, &
72 mpi_sum,mpicom,ierror)
73 if (spinpol) then
74 n=n*ndmag
75 call mpi_allreduce(mpi_in_place,dbxmt,n,mpi_double_precision,mpi_sum, &
77 n=ngtot*ndmag
78 call mpi_allreduce(mpi_in_place,dbxir,n,mpi_double_precision,mpi_sum, &
80 end if
81 end if
82! convert muffin-tin residuals to spherical harmonics
83 call holdthd(natmtot,nthd)
84!$OMP PARALLEL DO DEFAULT(SHARED) &
85!$OMP PRIVATE(is,nrc,nrci,idm) &
86!$OMP SCHEDULE(DYNAMIC) &
87!$OMP NUM_THREADS(nthd)
88 do ias=1,natmtot
89 is=idxis(ias)
90 nrc=nrcmt(is)
91 nrci=nrcmti(is)
92 call rfsht(nrc,nrci,dvxmt(:,ias),rfmt1(:,ias))
93 do idm=1,ndmag
94 call rfsht(nrc,nrci,dbxmt(:,ias,idm),rvfmt(:,ias,idm))
95 end do
96 end do
97!$OMP END PARALLEL DO
98 call freethd(nthd)
99! symmetrise the residuals
101 rfmt1,dvxir)
102 if (spinpol) then
104 igrzf,npmtmax,rvfmt,ngtot,dbxir)
105 end if
106! magnitude of residuals
107 resoep=rfinpc(npmtmax,rfmt1,dvxir,rfmt1,dvxir)
108 do idm=1,ndmag
109 t1=rfinpc(npmtmax,rvfmt(:,:,idm),dbxir(:,idm),rvfmt(:,:,idm),dbxir(:,idm))
110 resoep=resoep+t1
111 end do
112 resoep=sqrt(resoep)/omega
113! update the step size
114 if (it >= 2) then
115! check residual against previous value
116 if (resoep < resp) then
118 else
120 end if
121 end if
122! store previous residual
123 resp=resoep
124! update exchange potential and magnetic field
125 call holdthd(natmtot,nthd)
126!$OMP PARALLEL DO DEFAULT(SHARED) &
127!$OMP PRIVATE(rfmt2,is,nrc,nrci,npc,idm) &
128!$OMP SCHEDULE(DYNAMIC) &
129!$OMP NUM_THREADS(nthd)
130 do ias=1,natmtot
131 is=idxis(ias)
132 nrc=nrcmt(is)
133 nrci=nrcmti(is)
134 npc=npcmt(is)
135! convert residual to spherical coordinates
136 call rbsht(nrc,nrci,rfmt1(:,ias),rfmt2)
137! subtract from exchange potential
138 vxmt(1:npc,ias)=vxmt(1:npc,ias)-tauoep*rfmt2(1:npc)
139! repeat for exchange magnetic field
140 do idm=1,ndmag
141 call rbsht(nrc,nrci,rvfmt(:,ias,idm),rfmt2)
142 bxmt(1:npc,ias,idm)=bxmt(1:npc,ias,idm)-tauoep*rfmt2(1:npc)
143 end do
144 end do
145!$OMP END PARALLEL DO
146 call freethd(nthd)
147 vxir(:)=vxir(:)-tauoep*dvxir(:)
148 do idm=1,ndmag
149 bxir(:,idm)=bxir(:,idm)-tauoep*dbxir(:,idm)
150 end do
151! end iteration loop
152end do
153! convert the exchange potential and field to spherical harmonics
154call holdthd(natmtot,nthd)
155!$OMP PARALLEL DO DEFAULT(SHARED) &
156!$OMP PRIVATE(is,nrc,nrci,idm) &
157!$OMP SCHEDULE(DYNAMIC) &
158!$OMP NUM_THREADS(nthd)
159do ias=1,natmtot
160 is=idxis(ias)
161 nrc=nrcmt(is)
162 nrci=nrcmti(is)
163 call rfsht(nrc,nrci,vxmt(:,ias),rfmt1(:,ias))
164 do idm=1,ndmag
165 call rfsht(nrc,nrci,bxmt(:,ias,idm),rvfmt(:,ias,idm))
166 end do
167end do
168!$OMP END PARALLEL DO
169call freethd(nthd)
170! convert potential and field from a coarse to a fine radial mesh
171call rfmtctof(rfmt1)
172do idm=1,ndmag
173 call rfmtctof(rvfmt(:,:,idm))
174end do
175! add to existing (density derived) correlation potential and field
176do ias=1,natmtot
177 is=idxis(ias)
178 np=npmt(is)
179 vxcmt(1:np,ias)=vxcmt(1:np,ias)+rfmt1(1:np,ias)
180 do idm=1,ndmag
181 bxcmt(1:np,ias,idm)=bxcmt(1:np,ias,idm)+rvfmt(1:np,ias,idm)
182 end do
183end do
184vxcir(:)=vxcir(:)+vxir(:)
185do idm=1,ndmag
186 bxcir(:,idm)=bxcir(:,idm)+bxir(:,idm)
187end do
188! symmetrise the exchange potential and field
190 vxcir)
191if (spinpol) then
194end if
195deallocate(rfmt1,rfir,vclcv,vclvv)
196if (spinpol) deallocate(rvfmt,rvfir)
197! set the constant part of the exchange potential equal to zero
198call rfint0(0.d0,vxcmt,vxcir)
199end subroutine
200
real(8) resoep
Definition modmain.f90:1147
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable bxcmt
Definition modmain.f90:636
integer ngtot
Definition modmain.f90:390
integer, dimension(3) ngridg
Definition modmain.f90:386
logical ncmag
Definition modmain.f90:240
integer nfgrz
Definition modmain.f90:412
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:,:), allocatable dbxir
Definition modmain.f90:1145
real(8), dimension(:,:), allocatable bxir
Definition modmain.f90:1143
logical spinpol
Definition modmain.f90:228
integer, dimension(:), allocatable igrzf
Definition modmain.f90:416
real(8) omega
Definition modmain.f90:20
integer ngvec
Definition modmain.f90:396
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
real(8), dimension(:,:), allocatable bxcir
Definition modmain.f90:636
real(8), dimension(:), allocatable vxcir
Definition modmain.f90:634
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer nkpt
Definition modmain.f90:461
real(8), dimension(:), allocatable vxir
Definition modmain.f90:1143
integer natmtot
Definition modmain.f90:40
real(8) tau0oep
Definition modmain.f90:1141
integer maxitoep
Definition modmain.f90:1139
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer, dimension(:), allocatable igfft
Definition modmain.f90:406
real(8) tauoep
Definition modmain.f90:1141
integer nstsv
Definition modmain.f90:886
integer ncrmax
Definition modmain.f90:1137
integer npmtmax
Definition modmain.f90:216
real(8), dimension(:,:,:), allocatable bxmt
Definition modmain.f90:1143
integer, dimension(maxspecies) npmt
Definition modmain.f90:213
real(8), dimension(:,:), allocatable dvxmt
Definition modmain.f90:1145
real(8), dimension(:,:), allocatable vxcmt
Definition modmain.f90:634
real(8), dimension(:,:), allocatable vxmt
Definition modmain.f90:1143
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
integer iscl
Definition modmain.f90:1048
real(8), dimension(:), allocatable dvxir
Definition modmain.f90:1145
integer ndmag
Definition modmain.f90:238
real(8), dimension(:,:,:), allocatable dbxmt
Definition modmain.f90:1145
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 holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine oepmain
Definition oepmain.f90:7
subroutine oepresk(ik, vclcv, vclvv)
Definition oepresk.f90:7
subroutine oepvcl(vclcv, vclvv)
Definition oepvcl.f90:7
subroutine rbsht(nr, nri, rfmt1, rfmt2)
Definition rbsht.f90:7
real(8) function rfinpc(ld, rfmt1, rfir1, rfmt2, rfir2)
Definition rfinpc.f90:7
subroutine rfint0(rf0, rfmt, rfir)
Definition rfint0.f90:7
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10
subroutine rfsht(nr, nri, rfmt1, rfmt2)
Definition rfsht.f90:7
subroutine symrf(nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld, rfmt, rfir)
Definition symrf.f90:11
subroutine symrvf(tspin, tnc, nrmt_, nrmti_, npmt_, ngridg_, ngtot_, ngvec_, nfgrz_, igfft_, igrzf_, ld1, rvfmt, ld2, rvfir)
Definition symrvf.f90:11