aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2005-10-19 15:15:55 +0000
committerschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2005-10-19 15:15:55 +0000
commit5d7f29b265659df0b093a44f629b6c3d1c8b3b81 (patch)
tree115d941c570dde72ad50e484d4329daa0b835512
parent90cc3406685472af510d6d810f154ccc6f8220fe (diff)
Correct third derivative operators again
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinUtils/TGRtensor/trunk@33 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 22351a8..2e509d2 100644
--- a/src/derivs.F90
+++ b/src/derivs.F90
@@ -53,20 +53,20 @@ contains
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(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(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(1,3,2) = f(3,1,2)
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(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)
end subroutine get_derivs3
#ifndef TGR_INCLUDED