From 31fb6006228454beb2e01deddddbbe41225c4a38 Mon Sep 17 00:00:00 2001 From: schnetter Date: Mon, 9 May 2005 20:57:58 +0000 Subject: Disable Ricci self check by default. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinUtils/TGRtensor/trunk@26 b716e942-a2de-43ad-8f52-f3dfe468e4e7 --- src/ricci.F90 | 10 +++++++++- src/ricci4.F90 | 50 ++++++++++++++++++++++++++++++++------------------ 2 files changed, 41 insertions(+), 19 deletions(-) diff --git a/src/ricci.F90 b/src/ricci.F90 index 83b144a..6575994 100644 --- a/src/ricci.F90 +++ b/src/ricci.F90 @@ -2,6 +2,8 @@ #include "cctk.h" +#undef DEBUG + module ricci implicit none private @@ -82,7 +84,9 @@ contains subroutine calc_ricci (gamma, dgamma, ri) CCTK_REAL, intent(in) :: gamma(3,3,3), dgamma(3,3,3,3) CCTK_REAL, intent(out) :: ri(3,3) - CCTK_REAL :: sum, cnt +#ifdef DEBUG + CCTK_REAL :: nrm, cnt +#endif integer :: i,j,k,l ! R_ij = Gamma^k_ij,k - Gamma^k_ik,j ! + Gamma^k_lk Gamma^l_ij - Gamma^k_lj Gamma^l_ki @@ -98,6 +102,7 @@ contains end do end do end do +#ifdef DEBUG ! check symmetries sum = 0 cnt = 0 @@ -108,12 +113,15 @@ contains end do end do if (sum > 1.0e-12 * cnt) call CCTK_WARN (0, "Ricci tensor is not symmetric") +#endif end subroutine calc_ricci subroutine calc_ricciderivs (gamma, dgamma, ddgamma, dri) CCTK_REAL, intent(in) :: gamma(3,3,3), dgamma(3,3,3,3), ddgamma(3,3,3,3,3) CCTK_REAL, intent(out) :: dri(3,3,3) +#if 0 CCTK_REAL :: nrm, cnt +#endif integer :: i,j,k,l,m ! R_ij = Gamma^k_ij,k - Gamma^k_ik,j ! + Gamma^k_lk Gamma^l_ij - Gamma^k_lj Gamma^l_ki diff --git a/src/ricci4.F90 b/src/ricci4.F90 index ada7733..ac3ad3a 100644 --- a/src/ricci4.F90 +++ b/src/ricci4.F90 @@ -2,6 +2,8 @@ #include "cctk.h" +#undef DEBUG + module ricci4 implicit none private @@ -56,7 +58,9 @@ contains subroutine calc_4ricci (gamma, dgamma, ri) CCTK_REAL, intent(in) :: gamma(4,4,4), dgamma(4,4,4,4) CCTK_REAL, intent(out) :: ri(4,4) - CCTK_REAL :: sum, cnt +#ifdef DEBUG + CCTK_REAL :: nrm, cnt +#endif integer :: i,j,k,l ! R_ij = Gamma^k_ij,k - Gamma^k_ik,j ! + Gamma^k_lk Gamma^l_ij - Gamma^k_lj Gamma^l_ki @@ -72,23 +76,27 @@ contains end do end do end do +#ifdef DEBUG ! check symmetries - sum = 0 + nrm = 0 cnt = 0 do i=1,4 do j=1,4 - sum = sum + (ri(i,j) - ri(j,i))**2 + nrm = nrm + (ri(i,j) - ri(j,i))**2 cnt = cnt + ri(i,j)**2 end do end do - if (sum > 1.0e-12 * cnt) call CCTK_WARN (0, "4-Ricci tensor is not symmetric") + if (nrm > 1.0e-12 * cnt) call CCTK_WARN (0, "4-Ricci tensor is not symmetric") +#endif end subroutine calc_4ricci subroutine calc_4riemann (gg, gamma, dgamma, rm) CCTK_REAL, intent(in) :: gg(4,4), gamma(4,4,4), dgamma(4,4,4,4) CCTK_REAL, intent(out) :: rm(4,4,4,4) CCTK_REAL :: rmu(4,4,4,4) - CCTK_REAL :: sum, cnt +#ifdef DEBUG + CCTK_REAL :: nrm, cnt +#endif integer :: i,j,k,l,m ! R^i_jkl = Gamma^i_jl,k - Gamma^i_jk,l ! + Gamma^m_jl Gamma^i_mk - Gamma^m_jk Gamma^i_ml @@ -119,30 +127,34 @@ contains end do end do end do +#ifdef DEBUG ! check symmetries - sum = 0 + nrm = 0 cnt = 0 do i=1,4 do j=1,4 do k=1,4 do l=1,4 - sum = sum + (rm(i,j,k,l) + rm(i,j,l,k))**2 - sum = sum + (rm(i,j,k,l) - rm(k,l,i,j))**2 - sum = sum + (rm(i,j,k,l) + rm(j,i,k,l))**2 - sum = sum + (rm(i,j,k,l) + rm(j,k,i,l) + rm(k,i,j,l))**2 + nrm = nrm + (rm(i,j,k,l) + rm(i,j,l,k))**2 + nrm = nrm + (rm(i,j,k,l) - rm(k,l,i,j))**2 + nrm = nrm + (rm(i,j,k,l) + rm(j,i,k,l))**2 + nrm = nrm + (rm(i,j,k,l) + rm(j,k,i,l) + rm(k,i,j,l))**2 cnt = cnt + rm(i,j,k,l)**2 end do end do end do end do - if (sum > 1.0e-12 * cnt) call CCTK_WARN (0, "4-Riemann tensor is not symmetric") + if (nrm > 1.0e-12 * cnt) call CCTK_WARN (0, "4-Riemann tensor is not symmetric") +#endif end subroutine calc_4riemann subroutine calc_4weyl (gg, rm, ri, rsc, we) CCTK_REAL, parameter :: one=1, half=one/2, third=one/3 CCTK_REAL, intent(in) :: gg(4,4), rm(4,4,4,4), ri(4,4), rsc CCTK_REAL, intent(out) :: we(4,4,4,4) - CCTK_REAL :: sum, cnt +#ifdef DEBUG + CCTK_REAL :: nrm, cnt +#endif integer :: i,j,k,l ! R_ijkl = C_ijkl + 2/(n-2) (g_i[k R_l]j - g_j[k R_l]i) ! - 2/(n-1)(n-2) R g_i[k g_l]j @@ -159,23 +171,25 @@ contains end do end do end do +#ifdef DEBUG ! check symmetries - sum = 0 + nrm = 0 cnt = 0 do i=1,4 do j=1,4 do k=1,4 do l=1,4 - sum = sum + (we(i,j,k,l) + we(i,j,l,k))**2 - sum = sum + (we(i,j,k,l) - we(k,l,i,j))**2 - sum = sum + (we(i,j,k,l) + we(j,i,k,l))**2 - sum = sum + (we(i,j,k,l) + we(j,k,i,l) + we(k,i,j,l))**2 + nrm = nrm + (we(i,j,k,l) + we(i,j,l,k))**2 + nrm = nrm + (we(i,j,k,l) - we(k,l,i,j))**2 + nrm = nrm + (we(i,j,k,l) + we(j,i,k,l))**2 + nrm = nrm + (we(i,j,k,l) + we(j,k,i,l) + we(k,i,j,l))**2 cnt = cnt + we(i,j,k,l)**2 end do end do end do end do - if (sum > 1.0e-12 * cnt) call CCTK_WARN (0, "4-Weyl tensor is not symmetric") + if (nrm > 1.0e-12 * cnt) call CCTK_WARN (0, "4-Weyl tensor is not symmetric") +#endif end subroutine calc_4weyl end module ricci4 -- cgit v1.2.3