File: rfmtctof.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (100 lines) | stat: -rw-r--r-- 2,281 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

!BOP
! !ROUTINE: rfmtctof
! !INTERFACE:
subroutine rfmtctof(rfmt)
! !USES:
use modmain
use modomp
! !INPUT/OUTPUT PARAMETERS:
!   rfmt : real muffin-tin function (in,real(npmtmax,natmtot))
! !DESCRIPTION:
!   Converts a real muffin-tin function from a coarse to a fine radial mesh by
!   using cubic spline interpolation. See {\tt rfinterp} and {\tt spline}.
!
! !REVISION HISTORY:
!   Created October 2003 (JKD)
!EOP
!BOC
implicit none
! arguments
real(8), intent(inout) :: rfmt(npmtmax,natmtot)
! local variables
integer is,ias,lm,nthd
integer nr,nri,nro
integer iro,ir,npi,i
integer nrc,nrci,nrco
integer irco,irc,npci
! allocatable arrays
real(8), allocatable :: fi(:),fo(:),rfmt1(:)
if (lradstp.eq.1) return
call omp_hold(natmtot,nthd)
!$OMP PARALLEL DEFAULT(SHARED) &
!$OMP PRIVATE(fi,fo,rfmt1,is,nr,nri,nro) &
!$OMP PRIVATE(iro,npi,nrc,nrci,nrco) &
!$OMP PRIVATE(irco,npci,lm,i,irc,ir) &
!$OMP NUM_THREADS(nthd)
!$OMP DO
do ias=1,natmtot
  allocate(fi(nrcmtmax),fo(nrmtmax),rfmt1(npmtmax))
  is=idxis(ias)
  nr=nrmt(is)
  nri=nrmti(is)
  nro=nr-nri
  iro=nri+1
  npi=npmti(is)
  nrc=nrcmt(is)
  nrci=nrcmti(is)
  nrco=nrc-nrci
  irco=nrci+1
  npci=npcmti(is)
! interpolate up to lmaxi over entire muffin-tin
  do lm=1,lmmaxi
    i=lm
    do irc=1,nrci
      fi(irc)=rfmt(i,ias)
      i=i+lmmaxi
    end do
    do irc=irco,nrc
      fi(irc)=rfmt(i,ias)
      i=i+lmmaxo
    end do
    call rfinterp(nrc,rcmt(:,is),fi,nr,rsp(:,is),fo)
    i=lm
    do ir=1,nri
      rfmt1(i)=fo(ir)
      i=i+lmmaxi
    end do
    do ir=iro,nr
      rfmt1(i)=fo(ir)
      i=i+lmmaxo
    end do
  end do
! interpolate up to lmaxo on outer part of muffin-tin
  do lm=lmmaxi+1,lmmaxo
    i=npci+lm
    do irc=irco,nrc
      fi(irc)=rfmt(i,ias)
      i=i+lmmaxo
    end do
    call rfinterp(nrco,rcmt(irco,is),fi(irco),nro,rsp(iro,is),fo(iro))
    i=npi+lm
    do ir=iro,nr
      rfmt1(i)=fo(ir)
      i=i+lmmaxo
    end do
  end do
  call dcopy(npmt(is),rfmt1,1,rfmt(:,ias),1)
  deallocate(fi,fo,rfmt1)
end do
!$OMP END DO
!$OMP END PARALLEL
call omp_free(nthd)
return
end subroutine
!EOC