aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorguzman <guzman@b6f3ac56-194f-0410-8878-cdf6079d7f1b>2002-05-03 14:10:24 +0000
committerguzman <guzman@b6f3ac56-194f-0410-8878-cdf6079d7f1b>2002-05-03 14:10:24 +0000
commit4d10c7d660916c9acecf050533b003f3fbc492d0 (patch)
tree288d78fd0daab700194d7ea398132447d7ff7920
parent6675c2aa8ec85069fe4b1f7d71f994e073317a47 (diff)
Converted version for Einstein 2
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/IDAxiOddBrillBH/trunk@40 b6f3ac56-194f-0410-8878-cdf6079d7f1b
-rw-r--r--src/IDAxiOddBrillBH.F31
1 files changed, 27 insertions, 4 deletions
diff --git a/src/IDAxiOddBrillBH.F b/src/IDAxiOddBrillBH.F
index c49bf48..02865e3 100644
--- a/src/IDAxiOddBrillBH.F
+++ b/src/IDAxiOddBrillBH.F
@@ -25,9 +25,6 @@ c
c @endhistory
c@@*/
-c Need include file from Einstein
-#include "CactusEinstein/Einstein/src/Einstein.h"
-
subroutine IDAxiOddBrillBH(CCTK_ARGUMENTS)
implicit none
@@ -62,6 +59,7 @@ c Perhaps this and others should go into cctk.h
integer :: nx,ny,nz
integer :: i,j,k,it,ier,nquads,ntries
integer :: npoints,handle,ierror
+ integer make_conformal_derivs
pi = 4.0d0*atan(1.0d0)
c DO NOT use integer*4
@@ -85,7 +83,32 @@ c --------------
write(*,*) 'Bowen & York Momenta: :',byJ
write(*,*) 'Brill Wave: sin^n theta :',n
- conformal_state = CONFORMAL_METRIC
+c Check if we should create and store conformal factor stuff
+
+ if (CCTK_EQUALS(metric_type, "static conformal")) then
+
+ conformal_state = 1
+
+ if(CCTK_EQUALS(conformal_storage,"factor+derivs")) then
+
+ conformal_state = 2
+ make_conformal_derivs = 1
+
+ else if (CCTK_EQUALS(conformal_storage,"factor+derivs+2nd derivs")) then
+
+ conformal_state = 3
+ make_conformal_derivs = 1
+
+ endif
+
+ else
+
+ make_conformal_derivs = 0
+
+ endif
+
+c conformal_state = CONFORMAL_METRIC
+
c Sovle on this sized cartesian grid
c 3D grid size NE x NT x NP
c Add 2 zones for eta coordinate and 4 for theta