aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2010-04-07 22:59:56 +0000
committerschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2010-04-07 22:59:56 +0000
commit356d0361bd3d586d1761eae4ea81537cf95bd46c (patch)
tree5704045d9d54ea852afb51de1872d350ccd13e19
parent1839a7bfae273274a5738b7debdb5d2b80207b70 (diff)
Simplify code.
Take matter contribution into account when evaluating the Einstein equations. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinUtils/TGRtensor/trunk@51 b716e942-a2de-43ad-8f52-f3dfe468e4e7
-rw-r--r--src/adm_metric_simple.F9070
1 files changed, 37 insertions, 33 deletions
diff --git a/src/adm_metric_simple.F90 b/src/adm_metric_simple.F90
index 29d1b84..89ca057 100644
--- a/src/adm_metric_simple.F90
+++ b/src/adm_metric_simple.F90
@@ -5,19 +5,21 @@
module adm_metric_simple
! This module assumes that alpha=1 and beta^i=0
+ use constants
use tensor
implicit none
private
public calc_4metric_simple
public calc_4metricderivs_simple
public calc_4metricderivs2_simple
+
public calc_3metric_simple
public calc_3metricderivs_simple
+ public calc_extcurv_simple
+ public calc_extcurvdot_simple
public calc_3metricdot_simple
- public calc_extcurv_simple
public calc_3metricderivdot_simple
- public calc_extcurvdot_simple
public calc_3metricdot2_simple
contains
@@ -71,7 +73,7 @@ contains
CCTK_REAL, intent(out) :: ddg4(0:3,0:3,0:3,0:3)
CCTK_REAL :: d4gg(3,3,0:3)
CCTK_REAL :: dd4gg(3,3,0:3,0:3)
- integer :: i,j,k
+ integer :: i,j
! 4-metric
g4(0 ,0 ) = -1
@@ -118,7 +120,6 @@ contains
subroutine calc_3metric_simple (g4, gg)
CCTK_REAL, intent(in) :: g4(0:3,0:3)
CCTK_REAL, intent(out) :: gg(3,3)
- CCTK_REAL :: dtg, gu(3,3)
! ds^2 = -alpha^2 dt^2 + g_ij (dx^i + beta^i dt) (dx^j + beta^j dt)
gg = g4(1:3,1:3)
@@ -130,7 +131,7 @@ contains
CCTK_REAL, intent(out) :: dgg(3,3,3)
CCTK_REAL, intent(out) :: gg_dot(3,3)
CCTK_REAL :: d4gg(3,3,0:3)
- integer :: i,j
+ integer :: i
! ds^2 = -alpha^2 dt^2 + g_ij (dx^i + beta^i dt) (dx^j + beta^j dt)
gg = g4(1:3,1:3)
@@ -142,22 +143,47 @@ contains
dgg = d4gg(:,:,1:3)
end subroutine calc_3metricderivs_simple
-
-
- subroutine calc_extcurv_simple (gg, dgg, gg_dot, kk)
- CCTK_REAL, intent(in) :: gg(3,3), dgg(3,3,3), gg_dot(3,3)
+ subroutine calc_extcurv_simple (gg_dot, kk)
+ CCTK_REAL, intent(in) :: gg_dot(3,3)
CCTK_REAL, intent(out) :: kk(3,3)
integer :: i,j
- ! d/dt g_ij = -2 alpha K_ij + g_kj beta^k,i + g_ik beta^k,j + beta^k g_ij,k
+ ! d/dt g_ij = -2 K_ij
do i=1,3
do j=1,3
- kk(i,j) = - gg_dot(i,j) / 2
+ kk(i,j) = - gg_dot(i,j)
+ kk(i,j) = kk(i,j) / 2
end do
end do
+
end subroutine calc_extcurv_simple
+ subroutine calc_extcurvdot_simple (gg,gu,ri, kk, tt, kk_dot)
+ CCTK_REAL, intent(in) :: gg(3,3), gu(3,3)
+ CCTK_REAL, intent(in) :: ri(3,3)
+ CCTK_REAL, intent(in) :: kk(3,3)
+ CCTK_REAL, intent(in) :: tt(3,3)
+ CCTK_REAL, intent(out) :: kk_dot(3,3)
+ integer :: i,j,l,m
+
+ ! d/dt K_ij = R_ij - 2 K_il g^lm K_mj + g^lm K_lm K_ij
+ ! - 8 pi (T_ij - 1/2 g_ij T)
+
+ do i=1,3
+ do j=1,3
+ kk_dot(i,j) = ri(i,j) - 8*pi * tt(i,j)
+ do l=1,3
+ do m=1,3
+ kk_dot(i,j) = kk_dot(i,j) &
+ + gu(l,m) * (-2 * kk(i,l) * kk(m,j) + kk(l,m) * kk(i,j)) &
+ + 4*pi * gg(i,j) * gu(l,m) * tt(l,m)
+ end do
+ end do
+ end do
+ end do
+ end subroutine calc_extcurvdot_simple
+
subroutine calc_3metricdot_simple (kk, gg_dot)
@@ -190,28 +216,6 @@ contains
end do
end subroutine calc_3metricderivdot_simple
- subroutine calc_extcurvdot_simple (gu,ri, kk, kk_dot)
- CCTK_REAL, intent(in) :: gu(3,3)
- CCTK_REAL, intent(in) :: ri(3,3)
- CCTK_REAL, intent(in) :: kk(3,3)
- CCTK_REAL, intent(out) :: kk_dot(3,3)
- integer :: i,j,l,m
-
- ! d/dt K_ij = R_ij - 2 K_il g^lm K_mj + g^lm K_lm K_ij
-
- do i=1,3
- do j=1,3
- kk_dot(i,j) = ri(i,j)
- do l=1,3
- do m=1,3
- kk_dot(i,j) = kk_dot(i,j) + gu(l,m) &
- * (-2 * kk(i,l) * kk(m,j) + kk(l,m) * kk(i,j))
- end do
- end do
- end do
- end do
- end subroutine calc_extcurvdot_simple
-
subroutine calc_3metricdot2_simple (kk_dot, gg_dot2)
CCTK_REAL, intent(in) :: kk_dot(3,3)
CCTK_REAL, intent(out) :: gg_dot2(3,3)