The Elk Code
 
Loading...
Searching...
No Matches
rfmtctof.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 Lesser General Public
4! License. See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: rfmtctof
8! !INTERFACE:
9subroutine rfmtctof(rfmt)
10! !USES:
11use modmain
12use modomp
13! !INPUT/OUTPUT PARAMETERS:
14! rfmt : real muffin-tin function (in,real(npmtmax,natmtot))
15! !DESCRIPTION:
16! Converts a real muffin-tin function from a coarse to a fine radial mesh by
17! using cubic spline interpolation. See {\tt rfinterp} and {\tt spline}.
18!
19! !REVISION HISTORY:
20! Created October 2003 (JKD)
21!EOP
22!BOC
23implicit none
24! arguments
25real(8), intent(inout) :: rfmt(npmtmax,natmtot)
26! local variables
27integer is,ias,lm
28integer nr,nri,nro,iro
29integer nrc,nrci,nrco,irco
30integer i0,i1,nthd
31! automatic arrays
32real(8) rfmt1(npcmtmax),fi(nrcmtmax),fo(nrmtmax)
33if (lradstp == 1) return
34call holdthd(natmtot,nthd)
35!$OMP PARALLEL DO DEFAULT(SHARED) &
36!$OMP PRIVATE(rfmt1,fi,fo,is) &
37!$OMP PRIVATE(nr,nri,nro,iro) &
38!$OMP PRIVATE(nrc,nrci,nrco,irco) &
39!$OMP PRIVATE(lm,i0,i1) &
40!$OMP SCHEDULE(DYNAMIC) &
41!$OMP NUM_THREADS(nthd)
42do ias=1,natmtot
43 is=idxis(ias)
44 nr=nrmt(is)
45 nri=nrmti(is)
46 nro=nr-nri
47 iro=nri+1
48 nrc=nrcmt(is)
49 nrci=nrcmti(is)
50 nrco=nrc-nrci
51 irco=nrci+1
52! copy the input function
53 rfmt1(1:npcmt(is))=rfmt(1:npcmt(is),ias)
54! interpolate up to lmaxi over entire muffin-tin
55 do lm=1,lmmaxi
56 i1=lmmaxi*(nrci-1)+lm
57 fi(1:nrci)=rfmt1(lm:i1:lmmaxi)
58 i0=i1+lmmaxi
59 i1=lmmaxo*(nrc-irco)+i0
60 fi(irco:nrc)=rfmt1(i0:i1:lmmaxo)
61 call rfinterp(nrc,rcmt(:,is),wcrcmt(:,:,is),fi,nr,rlmt(:,1,is),fo)
62 i1=lmmaxi*(nri-1)+lm
63 rfmt(lm:i1:lmmaxi,ias)=fo(1:nri)
64 i0=i1+lmmaxi
65 i1=lmmaxo*(nr-iro)+i0
66 rfmt(i0:i1:lmmaxo,ias)=fo(iro:nr)
67 end do
68! interpolate up to lmaxo on outer part of muffin-tin
69 do lm=lmmaxi+1,lmmaxo
70 i0=lmmaxi*nrci+lm
71 i1=lmmaxo*(nrc-irco)+i0
72 fi(irco:nrc)=rfmt1(i0:i1:lmmaxo)
73 call rfinterp(nrco,rcmt(irco,is),wcrcmt(:,irco,is),fi(irco),nro, &
74 rsp(iro,is),fo(iro))
75 i0=lmmaxi*nri+lm
76 i1=lmmaxo*(nr-iro)+i0
77 rfmt(i0:i1:lmmaxo,ias)=fo(iro:nr)
78 end do
79end do
80!$OMP END PARALLEL DO
81call freethd(nthd)
82end subroutine
83!EOC
84
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
real(8), dimension(:,:), allocatable rcmt
Definition modmain.f90:177
integer lmmaxi
Definition modmain.f90:207
integer, dimension(maxspecies) nrcmt
Definition modmain.f90:173
real(8), dimension(:,:), allocatable rsp
Definition modmain.f90:135
integer, dimension(maxatoms *maxspecies) idxis
Definition modmain.f90:44
integer lradstp
Definition modmain.f90:171
real(8), dimension(:,:,:), allocatable wcrcmt
Definition modmain.f90:193
integer, dimension(maxspecies) npcmt
Definition modmain.f90:214
integer lmmaxo
Definition modmain.f90:203
integer, dimension(maxspecies) nrcmti
Definition modmain.f90:211
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
subroutine holdthd(nloop, nthd)
Definition modomp.f90:78
subroutine freethd(nthd)
Definition modomp.f90:106
subroutine rfinterp(ni, xi, wci, fi, no, xo, fo)
Definition rfinterp.f90:10
subroutine rfmtctof(rfmt)
Definition rfmtctof.f90:10