aboutsummaryrefslogtreecommitdiff
path: root/src/derivs.F90
blob: 22351a8dea9f95818d47e9d7ff0545e6fce49580 (plain)
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
! $Header$

#ifndef TGR_INCLUDED

#include "cctk.h"
#include "cctk_Parameters.h"

module derivs
  implicit none
  private
  public get_derivs
  public get_derivs2
  public get_derivs3
contains
#endif
  subroutine get_derivs (a, f, pos, off, dx)
    CCTK_REAL, intent(in)  :: a(*)
    CCTK_REAL, intent(out) :: f(3)
    integer,   intent(in)  :: pos, off(3)
    CCTK_REAL, intent(in)  :: dx(3)
    integer :: i
    do i=1,3
       f(i) = (a(pos+off(i)) - a(pos-off(i))) / (2*dx(i))
    end do
  end subroutine get_derivs
  subroutine get_derivs2 (a, f, pos, off, dx)
    CCTK_REAL, intent(in)  :: a(*)
    CCTK_REAL, intent(out) :: f(3,3)
    integer,   intent(in)  :: pos, off(3)
    CCTK_REAL, intent(in)  :: dx(3)
    integer :: i
    do i=1,3
       f(i,i) = (a(pos+off(i)) - 2*a(pos) + a(pos-off(i))) / dx(i)**2
    end do
    f(1,2) = (  a(pos-off(1)-off(2)) - a(pos+off(1)-off(2)) &
         &    - a(pos-off(1)+off(2)) + a(pos+off(1)+off(2))) / (4*dx(1)*dx(2))
    f(2,1) = f(1,2)
    f(1,3) = (  a(pos-off(1)-off(3)) - a(pos+off(1)-off(3)) &
         &    - a(pos-off(1)+off(3)) + a(pos+off(1)+off(3))) / (4*dx(1)*dx(3))
    f(3,1) = f(1,3)
    f(2,3) = (  a(pos-off(2)-off(3)) - a(pos+off(2)-off(3)) &
         &    - a(pos-off(2)+off(3)) + a(pos+off(2)+off(3))) / (4*dx(2)*dx(3))
    f(3,2) = f(2,3)
  end subroutine get_derivs2
  subroutine get_derivs3 (a, f, pos, off, dx)
    CCTK_REAL, intent(in)  :: a(*)
    CCTK_REAL, intent(out) :: f(3,3,3)
    integer,   intent(in)  :: pos, off(3)
    CCTK_REAL, intent(in)  :: dx(3)
    f(1,1,1) = (- a(pos+2*off(1)) + 2*a(pos+off(1)) - 2*a(pos-off(1)) + a(pos-2*off(1))) / (2*dx(1)**3)
    f(2,1,1) = (a(pos+off(1)+off(2)) - 2*a(pos+off(2)) + a(pos-off(1)+off(2)) - a(pos+off(1)-off(2)) + 2*a(pos-off(2)) - a(pos-off(1)-off(2))) / (2*dx(1)*dx(1)*dx(2))
    f(3,1,1) = (a(pos+off(1)+off(3)) - 2*a(pos+off(3)) + a(pos-off(1)+off(3)) - a(pos+off(1)-off(3)) + 2*a(pos-off(3)) - a(pos-off(1)-off(3))) / (2*dx(1)*dx(1)*dx(3))
    f(1,2,1) = f(2,1,1)
    f(2,2,1) = (a(pos+off(1)+off(2)) - 2*a(pos+off(1)) + a(pos+off(1)-off(2)) - a(pos-off(1)+off(2)) + 2*a(pos-off(1)) - a(pos-off(1)-off(2))) / (2*dx(1)*dx(2)*dx(2))
    f(3,2,1) = (a(pos+off(1)+off(2)+off(3)) - a(pos-off(1)+off(2)+off(3)) - a(pos+off(1)-off(2)+off(3)) + a(pos-off(1)-off(2)+off(3)) - a(pos+off(1)+off(2)-off(3)) + a(pos-off(1)+off(2)-off(3)) + a(pos+off(1)-off(2)-off(3)) - a(pos-off(1)-off(2)-off(3))) * (8*dx(1)*dx(2)*dx(3))
    f(1,3,1) = f(3,3,1)
    f(2,3,1) = f(3,2,1)
    f(3,3,1) = (a(pos+off(1)+off(3)) - 2*a(pos+off(1)) + a(pos+off(1)-off(3)) - a(pos-off(1)+off(3)) + 2*a(pos-off(1)) - a(pos-off(1)-off(3))) / (2*dx(1)*dx(3)*dx(3))
    f(:,1,2) = f(:,2,1)
    f(1,2,2) = f(2,2,1)
    f(2,2,2) = (- a(pos+2*off(2)) + 2*a(pos+off(2)) - 2*a(pos-off(2)) + a(pos-2*off(2))) / (2*dx(2)**3)
    f(3,2,2) = (a(pos+off(2)+off(3)) - 2*a(pos+off(3)) + a(pos-off(2)+off(3)) - a(pos+off(2)-off(3)) + 2*a(pos-off(3)) - a(pos-off(2)-off(3))) / (2*dx(2)*dx(2)*dx(3))
    f(1,3,2) = f(3,2,1)
    f(2,3,2) = f(3,2,2)
    f(3,3,2) = (a(pos+off(2)+off(3)) - 2*a(pos+off(2)) + a(pos+off(2)-off(3)) - a(pos-off(2)+off(3)) + 2*a(pos-off(2)) - a(pos-off(2)-off(3))) / (2*dx(2)*dx(3)*dx(3))
    f(:,1,3) = f(:,3,1)
    f(:,2,3) = f(:,3,2)
    f(1,3,3) = f(3,3,1)
    f(2,3,3) = f(3,3,2)
    f(3,3,3) = (- a(pos+2*off(3)) + 2*a(pos+off(3)) - 2*a(pos-off(3)) + a(pos-2*off(3))) / (2*dx(3)**3)
  end subroutine get_derivs3
#ifndef TGR_INCLUDED
end module derivs
#endif