diff options
author | schnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7> | 2006-08-16 20:02:03 +0000 |
---|---|---|
committer | schnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7> | 2006-08-16 20:02:03 +0000 |
commit | 91059f0111b815b4a11ee3712c453353bd66a929 (patch) | |
tree | d40d56990a89bc787b552597a306805326607cc2 | |
parent | 7daf9b0034c36ef4627a2568ba191555f2c4ac77 (diff) |
Add routines to calculate the time derivatives of the metric from the
extrinsic curvature and the Einstein equations.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinUtils/TGRtensor/trunk@40 b716e942-a2de-43ad-8f52-f3dfe468e4e7
-rw-r--r-- | src/adm_metric.F90 | 75 |
1 files changed, 75 insertions, 0 deletions
diff --git a/src/adm_metric.F90 b/src/adm_metric.F90 index 3df077e..9e284cd 100644 --- a/src/adm_metric.F90 +++ b/src/adm_metric.F90 @@ -14,6 +14,12 @@ module adm_metric public calc_3metricdot public calc_extcurv + + public calc_3metricdot_simple + public calc_3metricderivdot_simple + public calc_extcurvdot_simple + public calc_3metricdot2_simple + contains subroutine calc_4metric (gg, alfa, beta, g4) @@ -273,5 +279,74 @@ contains end do end subroutine calc_extcurv + + + + subroutine calc_3metricdot_simple (kk, gg_dot) + CCTK_REAL, intent(in) :: kk(3,3) + CCTK_REAL, intent(out) :: gg_dot(3,3) + integer :: i,j + + ! d/dt g_ij = -2 K_ij + + do i=1,3 + do j=1,3 + gg_dot(i,j) = -2 * kk(i,j) + end do + end do + end subroutine calc_3metricdot_simple + + subroutine calc_3metricderivdot_simple (dkk, dgg_dot) + CCTK_REAL, intent(in) :: dkk(3,3,3) + CCTK_REAL, intent(out) :: dgg_dot(3,3,3) + integer :: i,j,k + ! d/dt g_ij,k = -2 K_ij,k + + do i=1,3 + do j=1,3 + do k=1,3 + dgg_dot(i,j,k) = -2 * dkk(i,j,k) + end do + end do + 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) + integer :: i,j + + ! d/dt g_ij = -2 K_ij + ! d^2/dt^2 g_ij = -2 d/dt K_ij + + do i=1,3 + do j=1,3 + gg_dot2(i,j) = -2 * kk_dot(i,j) + end do + end do + end subroutine calc_3metricdot2_simple + end module adm_metric |