aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2005-10-21 14:04:55 +0000
committerschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2005-10-21 14:04:55 +0000
commita23a211054dde7ae7bc968a4b63ecef756c34587 (patch)
tree8a1de316bf73f68d445a81812fe61f8fc5049e25
parentb5796e9102b6a39b53a9cf0e4bde9b855974a0d3 (diff)
Correct more errors in the third derivatives
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinUtils/TGRtensor/trunk@35 b716e942-a2de-43ad-8f52-f3dfe468e4e7
-rw-r--r--src/derivs.F908
1 files changed, 4 insertions, 4 deletions
diff --git a/src/derivs.F90 b/src/derivs.F90
index 2e509d2..95e0fe5 100644
--- a/src/derivs.F90
+++ b/src/derivs.F90
@@ -47,18 +47,18 @@ contains
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(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(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,1,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(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,1,2)
f(2,3,2) = f(3,2,2)
@@ -67,7 +67,7 @@ contains
f(:,2,3) = f(:,3,2)
f(1,3,3) = f(3,1,3)
f(2,3,3) = f(3,2,3)
- 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)
+ 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