diff options
author | allen <allen@6a3ddf76-46e1-4315-99d9-bc56cac1ef84> | 1999-03-15 15:06:33 +0000 |
---|---|---|
committer | allen <allen@6a3ddf76-46e1-4315-99d9-bc56cac1ef84> | 1999-03-15 15:06:33 +0000 |
commit | fa4dcb58bbb9348ceaa571780d42e8033056af65 (patch) | |
tree | 41cac10802ea1fc25b9d4541f445326abd4dfa62 | |
parent | 00543fe0be420faddc2cacd95be540144ccd74e7 (diff) |
Fixing some small problems
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/IDAnalyticBH/trunk@5 6a3ddf76-46e1-4315-99d9-bc56cac1ef84
-rw-r--r-- | param.ccl | 4 | ||||
-rw-r--r-- | schedule.ccl | 7 | ||||
-rw-r--r-- | src/BrillLindquist.F | 6 | ||||
-rw-r--r-- | src/Misner.F | 6 | ||||
-rw-r--r-- | src/ParamChecker.F | 6 | ||||
-rw-r--r-- | src/Schwarzschild.F | 9 |
6 files changed, 25 insertions, 13 deletions
@@ -5,7 +5,7 @@ private: REAL mass "Mass of black hole" { -: :: "Not sure if it can be negative or no" +: :: "Not sure if it can be negative or not" } 2.0 @@ -103,7 +103,7 @@ REAL bl_M_4 "Mass of 4th BL hole" #{ #} -#EXTENDS KEYWORD model "" +#EXTENDS KEYWORD initial_data "" #{ # "schwarzschild" :: "One Schwarzshild black hole" #} diff --git a/schedule.ccl b/schedule.ccl index bde3600..6812c02 100644 --- a/schedule.ccl +++ b/schedule.ccl @@ -1,7 +1,8 @@ # Schedule definitions for thorn IDAnalyticBH # $Header$ -if (CCTK_Equals(model,"schwarzschild")) + +if (CCTK_Equals(initial_data,"schwarzschild")) { if (use_conformal) { @@ -15,7 +16,7 @@ if (CCTK_Equals(model,"schwarzschild")) } -if (CCTK_Equals(model,"bl_bh")) +if (CCTK_Equals(initial_data,"bl_bh")) { STORAGE: confac, confac_1derivs, confac_2derivs @@ -26,7 +27,7 @@ if (CCTK_Equals(model,"bl_bh")) } -if (CCTK_Equals(model,"misner_bh")) +if (CCTK_Equals(initial_data,"misner_bh")) { STORAGE: confac, confac_1derivs, confac_2derivs diff --git a/src/BrillLindquist.F b/src/BrillLindquist.F index 109ed5f..38db1d8 100644 --- a/src/BrillLindquist.F +++ b/src/BrillLindquist.F @@ -2,6 +2,8 @@ #include "declare_arguments.h" #include "declare_parameters.h" +#include "../../Einstein/src/Einstein.h" + subroutine BrillLindquist(CCTK_FARGUMENTS) implicit none @@ -17,9 +19,9 @@ c Must use conformal metric c ------------------------- - if (conformal_state == 0) then + if (conformal_state == NOCONFORMAL_METRIC) then call CCTK_Warn(1,"Changing to conformal metric in BrillLindquist") - conformal_state = 1 + conformal_state = CONFORMAL_METRIC end if c Put parameters into arrays for following calculations diff --git a/src/Misner.F b/src/Misner.F index fa65290..532c1ce 100644 --- a/src/Misner.F +++ b/src/Misner.F @@ -4,6 +4,8 @@ C This wrapper written by Carsten Gundlach. #include "declare_arguments.h" #include "declare_parameters.h" +#include "../../Einstein/src/Einstein.h" + subroutine Misner(CCTK_FARGUMENTS) implicit none @@ -22,9 +24,9 @@ C Finite differencing step. c Must use conformal metric c ------------------------- - if (conformal_state == 0) then + if (conformal_state == NOCONFORMAL_METRIC) then call CCTK_Warn(1,"Changing to conformal metric in Misner") - conformal_state = 1 + conformal_state = CONFORMAL_METRIC end if C Initialize some C global variables before calls to diff --git a/src/ParamChecker.F b/src/ParamChecker.F index 6cbd070..10c9caa 100644 --- a/src/ParamChecker.F +++ b/src/ParamChecker.F @@ -11,20 +11,20 @@ integer CCTK_Equals - if (CCTK_Equals(model,"schwarzschild")) then + if (CCTK_Equals(initial_data,"schwarzschild")==1) then write(*,*) 'One Black Hole: Throat at :',mass/2d0 if (use_conformal == 1) then write(*,*) "Uses conformal metric" else write(*,*) "Uses non-conformal metric" end if - else if (CCTK_Equals(model,"bl")) then + else if (CCTK_Equals(initial_data,"bl")==1) then write(*,*) "Setting up Brill Lindquist Solution" write(*,*) " Uses conformal metric" write(*,*) " Number of black holes (bl_nbh) :",bl_nbh write(*,*) " Black hole masses (bl_M_?) :",bl_M_1,bl_M_2, & bl_M_3,bl_M_4 - else if (CCTK_Equals(model,"misner")) then + else if (CCTK_Equals(initial_data,"misner")==1) then write(*,*) "Setting up Misner Solution" write(*,*) " Uses conformal metric" end if diff --git a/src/Schwarzschild.F b/src/Schwarzschild.F index 83b2cd9..d597e45 100644 --- a/src/Schwarzschild.F +++ b/src/Schwarzschild.F @@ -2,6 +2,9 @@ #include "declare_arguments.h" #include "declare_parameters.h" +c Need include file from Einstein +#include "../../packages/CactusEinstein/Einstein/src/Einstein.h" + subroutine Schwarzschild(CCTK_FARGUMENTS) implicit none @@ -24,7 +27,9 @@ three = 3.d0 c conformal metric flag - if (use_conformal==1) then + if (use_conformal == 1) then + + conformal_state = CONFORMAL_METRIC c Compute conformal factor psi = ( one + mass/two/r) @@ -51,6 +56,8 @@ c derivatives of psi / psi else + conformal_state = NOCONFORMAL_METRIC + gxx = ( one + mass/two/r)**4 gyy = gxx gzz = gxx |