The Elk Code
zfmtctof.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2013 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU Lesser General Public
4 ! License. See the file COPYING for license details.
5 
6 subroutine zfmtctof(zfmt)
7 use modmain
8 use modomp
9 implicit none
10 ! arguments
11 complex(8), intent(inout) :: zfmt(npmtmax,natmtot)
12 ! local variables
13 integer is,ias,lm
14 integer nr,nri,nro,iro
15 integer nrc,nrci,nrco,irco
16 integer i0,i1,nthd
17 ! automatic arrays
18 real(8) fi1(nrcmtmax),fi2(nrcmtmax)
19 real(8) fo1(nrmtmax),fo2(nrmtmax)
20 complex(8) zfmt1(npcmtmax)
21 if (lradstp == 1) return
22 call holdthd(natmtot,nthd)
23 !$OMP PARALLEL DO DEFAULT(SHARED) &
24 !$OMP PRIVATE(zfmt1,fi1,fi2,fo1,fo2,is) &
25 !$OMP PRIVATE(nr,nri,nro,iro) &
26 !$OMP PRIVATE(nrc,nrci,nrco,irco) &
27 !$OMP PRIVATE(lm,i0,i1) &
28 !$OMP SCHEDULE(DYNAMIC) &
29 !$OMP NUM_THREADS(nthd)
30 do ias=1,natmtot
31  is=idxis(ias)
32  nr=nrmt(is)
33  nri=nrmti(is)
34  nro=nr-nri
35  iro=nri+1
36  nrc=nrcmt(is)
37  nrci=nrcmti(is)
38  nrco=nrc-nrci
39  irco=nrci+1
40 ! copy the input function
41  zfmt1(1:npcmt(is))=zfmt(1:npcmt(is),ias)
42 ! interpolate up to lmaxi over entire muffin-tin
43  do lm=1,lmmaxi
44  i1=lmmaxi*(nrci-1)+lm
45  fi1(1:nrci)=dble(zfmt1(lm:i1:lmmaxi))
46  fi2(1:nrci)=aimag(zfmt1(lm:i1:lmmaxi))
47  i0=i1+lmmaxi
48  i1=lmmaxo*(nrc-irco)+i0
49  fi1(irco:nrc)=dble(zfmt1(i0:i1:lmmaxo))
50  fi2(irco:nrc)=aimag(zfmt1(i0:i1:lmmaxo))
51  call rfinterp(nrc,rcmt(:,is),wcrcmt(:,:,is),fi1,nr,rlmt(:,1,is),fo1)
52  call rfinterp(nrc,rcmt(:,is),wcrcmt(:,:,is),fi2,nr,rlmt(:,1,is),fo2)
53  i1=lmmaxi*(nri-1)+lm
54  zfmt(lm:i1:lmmaxi,ias)=cmplx(fo1(1:nri),fo2(1:nri),8)
55  i0=i1+lmmaxi
56  i1=lmmaxo*(nr-iro)+i0
57  zfmt(i0:i1:lmmaxo,ias)=cmplx(fo1(iro:nr),fo2(iro:nr),8)
58  end do
59 ! interpolate up to lmaxo on outer part of muffin-tin
60  do lm=lmmaxi+1,lmmaxo
61  i0=lmmaxi*nrci+lm
62  i1=lmmaxo*(nrc-irco)+i0
63  fi1(irco:nrc)=dble(zfmt1(i0:i1:lmmaxo))
64  fi2(irco:nrc)=aimag(zfmt1(i0:i1:lmmaxo))
65  call rfinterp(nrco,rcmt(irco,is),wcrcmt(:,irco,is),fi1(irco),nro, &
66  rsp(iro,is),fo1(iro))
67  call rfinterp(nrco,rcmt(irco,is),wcrcmt(:,irco,is),fi2(irco),nro, &
68  rsp(iro,is),fo2(iro))
69  i0=lmmaxi*nri+lm
70  i1=lmmaxo*(nr-iro)+i0
71  zfmt(i0:i1:lmmaxo,ias)=cmplx(fo1(iro:nr),fo2(iro:nr),8)
72  end do
73 end do
74 !$OMP END PARALLEL DO
75 call freethd(nthd)
76 end subroutine
77 
real(8), dimension(:,:), allocatable rcmt
Definition: modmain.f90:177
subroutine zfmtctof(zfmt)
Definition: zfmtctof.f90:7
integer, dimension(maxspecies) npcmt
Definition: modmain.f90:214
integer lmmaxo
Definition: modmain.f90:203
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
Definition: modomp.f90:6
subroutine rfinterp(ni, xi, wci, fi, no, xo, fo)
Definition: rfinterp.f90:10
integer lradstp
Definition: modmain.f90:171
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
real(8), dimension(:,:,:), allocatable wcrcmt
Definition: modmain.f90:193
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
integer, dimension(maxspecies) nrcmt
Definition: modmain.f90:173
integer, dimension(maxspecies) nrcmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150