File: grad2rfmt.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 (77 lines) | stat: -rw-r--r-- 1,613 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

! Copyright (C) 2009 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
! This file is distributed under the terms of the GNU Lesser General Public
! License. See the file COPYING for license details.

subroutine grad2rfmt(nr,nri,r,rfmt,g2rfmt)
use modmain
implicit none
! arguments
integer, intent(in) :: nr,nri
real(8), intent(in) :: r(nr)
real(8), intent(in) :: rfmt(*)
real(8), intent(out) :: g2rfmt(*)
! local variables
integer nro,iro,ir
integer l,m,lm,npi,i
real(8) t1
! automatic arrays
real(8) ri(nr),ri2(nr),fr(nr),cf(3,nr)
! tabulate 1/r and 1/r^2
do ir=1,nr
  ri(ir)=1.d0/r(ir)
  ri2(ir)=ri(ir)**2
end do
nro=nr-nri
iro=nri+1
npi=lmmaxi*nri
lm=0
do l=0,lmaxi
  t1=-dble(l*(l+1))
  do m=-l,l
    lm=lm+1
! use a cubic spline to compute radial derivatives
    i=lm
    do ir=1,nri
      fr(ir)=rfmt(i)
      i=i+lmmaxi
    end do
    do ir=iro,nr
      fr(ir)=rfmt(i)
      i=i+lmmaxo
    end do
    call spline(nr,r,fr,cf)
! apply Laplacian
    i=lm
    do ir=1,nri
      g2rfmt(i)=2.d0*(ri(ir)*cf(1,ir)+cf(2,ir))+ri2(ir)*t1*rfmt(i)
      i=i+lmmaxi
    end do
    do ir=iro,nr
      g2rfmt(i)=2.d0*(ri(ir)*cf(1,ir)+cf(2,ir))+ri2(ir)*t1*rfmt(i)
      i=i+lmmaxo
    end do
  end do
end do
do l=lmaxi+1,lmaxo
  t1=-dble(l*(l+1))
  do m=-l,l
    lm=lm+1
    i=npi+lm
    do ir=iro,nr
      fr(ir)=rfmt(i)
      i=i+lmmaxo
    end do
    call spline(nro,r(iro),fr(iro),cf(1,iro))
    i=npi+lm
    do ir=iro,nr
      g2rfmt(i)=2.d0*(ri(ir)*cf(1,ir)+cf(2,ir))+ri2(ir)*t1*rfmt(i)
      i=i+lmmaxo
    end do
  end do
end do
! apply smoothing if required
call rfmtsm(msmooth,nr,nri,g2rfmt)
return
end subroutine