aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorgoodale <goodale@5c0f84ea-6048-4d6e-bfa4-55cd5f2e0dd7>2002-04-26 16:22:11 +0000
committergoodale <goodale@5c0f84ea-6048-4d6e-bfa4-55cd5f2e0dd7>2002-04-26 16:22:11 +0000
commit302ba3da79031ba347157c4c28ca82a9cf3ed4e6 (patch)
treef6abb1779185b2a6912235be8e5a4c6846158bf3
parentcb8aa9cd754854014d4a44a77d514fc30db3e010 (diff)
Converted to new ADMBase, StaticConformal stuff.
Tom git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/IDLinearWaves/trunk@68 5c0f84ea-6048-4d6e-bfa4-55cd5f2e0dd7
-rw-r--r--interface.ccl2
-rw-r--r--param.ccl8
-rw-r--r--schedule.ccl12
-rw-r--r--src/ParamCheck.c62
-rw-r--r--src/make.code.defn2
-rw-r--r--src/planewaves.F7741
-rw-r--r--src/teukwaves.F7740
-rw-r--r--test/test_pw_ADM_leap.par23
-rw-r--r--test/test_pw_ADM_sl.par25
9 files changed, 160 insertions, 55 deletions
diff --git a/interface.ccl b/interface.ccl
index 189dea2..319510e 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -1,4 +1,4 @@
# Interface definition for thorn IDLinearWaves
implements: idlinearwaves
-inherits: einstein
+inherits: ADMBase StaticConformal Grid
diff --git a/param.ccl b/param.ccl
index 65b7bb6..ab334be 100644
--- a/param.ccl
+++ b/param.ccl
@@ -1,6 +1,6 @@
# Parameter definitions for thorn IDLinearWaves
-shares:einstein
+shares: ADMBase
EXTENDS KEYWORD initial_data
{
@@ -8,9 +8,11 @@ EXTENDS KEYWORD initial_data
"planewaves" :: "linear waves initial data- plane waves"
}
-USES BOOLEAN use_conformal
+USES KEYWORD metric_type
-USES BOOLEAN use_conformal_derivs
+shares: StaticConformal
+
+USES KEYWORD conformal_storage
private:
diff --git a/schedule.ccl b/schedule.ccl
index 5c536b2..d29b43d 100644
--- a/schedule.ccl
+++ b/schedule.ccl
@@ -1,7 +1,14 @@
# Schedule definitions for thorn IDLinearWaves
+
+
if (CCTK_Equals(initial_data,"planewaves"))
{
+ schedule IDLinearWaves_ParamChecker at CCTK_PARAMCHECK
+ {
+ LANG: C
+ } "Check that the metric_type is recognised"
+
schedule planewaves at CCTK_INITIAL
{
LANG: Fortran
@@ -10,6 +17,11 @@ if (CCTK_Equals(initial_data,"planewaves"))
if (CCTK_Equals(initial_data,"teukwaves"))
{
+ schedule IDLinearWaves_ParamChecker at CCTK_PARAMCHECK
+ {
+ LANG: C
+ } "Check that the metric_type is recognised"
+
schedule teukwaves at CCTK_INITIAL
{
LANG: Fortran
diff --git a/src/ParamCheck.c b/src/ParamCheck.c
new file mode 100644
index 0000000..c67605f
--- /dev/null
+++ b/src/ParamCheck.c
@@ -0,0 +1,62 @@
+ /*@@
+ @file ParamCheck.c
+ @date Fri Apr 26 18:03:09 2002
+ @author Tom Goodale
+ @desc
+ Check the parameters for IDLinearwaves
+ @enddesc
+ @version $Header$
+ @@*/
+
+#include "cctk.h"
+
+#include "cctk_Arguments.h"
+#include "cctk_Parameters.h"
+
+static const char *rcsid = "$Header$";
+
+CCTK_FILEVERSION(CactusEinstein_IDLinearWaves_ParamCheck_c)
+
+/********************************************************************
+ ********************* Local Data Types ***********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Local Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ***************** Scheduled Routine Prototypes *********************
+ ********************************************************************/
+
+void IDLinearWaves_ParamChecker(CCTK_ARGUMENTS);
+
+/********************************************************************
+ ********************* Other Routine Prototypes *********************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* Local Data *****************************
+ ********************************************************************/
+
+/********************************************************************
+ ********************* External Routines **********************
+ ********************************************************************/
+
+void IDLinearWaves_ParamChecker(CCTK_ARGUMENTS)
+{
+ DECLARE_CCTK_ARGUMENTS
+ DECLARE_CCTK_PARAMETERS
+
+ /* Do we know how to deal with this type of metric ? */
+ if(! CCTK_EQUALS(metric_type, "physical") &&
+ ! CCTK_EQUALS(metric_type, "static conformal"))
+ {
+ CCTK_PARAMWARN("Unknown ADMBase::metric_type - known types are \"physical\" and \"static conformal\"");
+ }
+}
+
+/********************************************************************
+ ********************* Local Routines *************************
+ ********************************************************************/
+
diff --git a/src/make.code.defn b/src/make.code.defn
index a152d03..9969001 100644
--- a/src/make.code.defn
+++ b/src/make.code.defn
@@ -2,7 +2,7 @@
# $Header$
# Source files in this directory
-SRCS = \
+SRCS = ParamCheck.c\
planewaves.F77\
teukwaves.F77
diff --git a/src/planewaves.F77 b/src/planewaves.F77
index 79c1d4e..880ff58 100644
--- a/src/planewaves.F77
+++ b/src/planewaves.F77
@@ -37,8 +37,6 @@
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
-c Using macro definitions from Einstein
-#include "CactusEinstein/Einstein/src/Einstein.h"
subroutine planewaves(CCTK_ARGUMENTS)
@@ -229,28 +227,41 @@ c loop over sh ends here:
enddo
c initialize the conformal factor
- if (use_conformal == 1) then
- 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;
+
+ else if(CCTK_EQUALS(conformal_storage,"factor+derivs+2nd derivs")) then
+ conformal_state = 3;
+ end if
+
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
psi(i,j,k) = 1d0
- psix(i,j,k) = 0d0
- psiy(i,j,k) = 0d0
- psiz(i,j,k) = 0d0
- psixy(i,j,k) = 0d0
- psixz(i,j,k) = 0d0
- psiyz(i,j,k) = 0d0
- psixx(i,j,k) = 0d0
- psiyy(i,j,k) = 0d0
- psizz(i,j,k) = 0d0
+ if(conformal_state .gt. 1) then
+ psix(i,j,k) = 0d0
+ psiy(i,j,k) = 0d0
+ psiz(i,j,k) = 0d0
+ endif
+ if(conformal_state .gt. 2) then
+ psixy(i,j,k) = 0d0
+ psixz(i,j,k) = 0d0
+ psiyz(i,j,k) = 0d0
+ psixx(i,j,k) = 0d0
+ psiyy(i,j,k) = 0d0
+ psizz(i,j,k) = 0d0
+ endif
end do
end do
end do
- else
- conformal_state = NOCONFORMAL_METRIC
end if
diff --git a/src/teukwaves.F77 b/src/teukwaves.F77
index 9e02eb0..dd31347 100644
--- a/src/teukwaves.F77
+++ b/src/teukwaves.F77
@@ -569,30 +569,42 @@ c time symmetry
enddo
c initialize the conformal factor
- if (use_conformal == 1) then
- 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;
+
+ else if(CCTK_EQUALS(conformal_storage,"factor+derivs+2nd derivs")) then
+ conformal_state = 3;
+ end if
+
do k=1,cctk_lsh(3)
do j=1,cctk_lsh(2)
do i=1,cctk_lsh(1)
psi(i,j,k) = 1d0
- psix(i,j,k) = 0d0
- psiy(i,j,k) = 0d0
- psiz(i,j,k) = 0d0
- psixy(i,j,k) = 0d0
- psixz(i,j,k) = 0d0
- psiyz(i,j,k) = 0d0
- psixx(i,j,k) = 0d0
- psiyy(i,j,k) = 0d0
- psizz(i,j,k) = 0d0
+ if(conformal_state .gt. 1) then
+ psix(i,j,k) = 0d0
+ psiy(i,j,k) = 0d0
+ psiz(i,j,k) = 0d0
+ endif
+ if(conformal_state .gt. 2) then
+ psixy(i,j,k) = 0d0
+ psixz(i,j,k) = 0d0
+ psiyz(i,j,k) = 0d0
+ psixx(i,j,k) = 0d0
+ psiyy(i,j,k) = 0d0
+ psizz(i,j,k) = 0d0
+ endif
end do
end do
end do
- else
- conformal_state = NOCONFORMAL_METRIC
end if
-
return
end
diff --git a/test/test_pw_ADM_leap.par b/test/test_pw_ADM_leap.par
index fda8c39..7a5db46 100644
--- a/test/test_pw_ADM_leap.par
+++ b/test/test_pw_ADM_leap.par
@@ -3,7 +3,7 @@
#########################################################
# Required thorns
-ActiveThorns = "Boundary time pugh pughslab pughreduce cartgrid3d einstein ADM IDLinearWaves ioascii ioutil iobasic"
+ActiveThorns = "Boundary time pugh pughslab pughreduce cartgrid3d ADMBase ADM ADMMacros StaticConformal CoordGauge IDLinearWaves ioascii ioutil iobasic"
# GENERAL
@@ -22,12 +22,15 @@ cactus::cctk_itlast = 100
# Einstein
-einstein::initial_data = "planewaves"
-einstein::evolution_system = "ADM"
+admbase::metric_type = "static conformal"
+admbase::initial_data = "planewaves"
+admbase::evolution_method = "ADM"
adm::method = "leapfrog"
adm::bound = "flat"
-einstein::slicing = "geodesic"
+
+admbase::lapse_evolution_method = "coordgauge"
+coordgauge::lapse_list = "geodesic"
# IDLinearwaves
idlinearwaves::amplitude = 0.001
@@ -45,15 +48,15 @@ IO::out_fileinfo = "none"
IO::parfile_write = "no"
IOBasic::outScalar_every = 10
-IOBasic::outScalar_vars = "einstein::gzz einstein::kzz einstein::alp"
+IOBasic::outScalar_vars = "admbase::gzz admbase::kzz admbase::alp"
IOASCII::out1D_every = 10
-IOASCII::out1D_vars = "einstein::gxx einstein::gxy einstein::gxz
- einstein::gyy einstein::gyz einstein::gzz
- einstein::kxx einstein::kxy einstein::kxz
- einstein::kyy einstein::kyz einstein::kzz"
+IOASCII::out1D_vars = "admbase::gxx admbase::gxy admbase::gxz
+ admbase::gyy admbase::gyz admbase::gzz
+ admbase::kxx admbase::kxy admbase::kxz
+ admbase::kyy admbase::kyz admbase::kzz"
IOBasic::outInfo_every = 5
-IOBasic::outInfo_vars = "einstein::gzz"
+IOBasic::outInfo_vars = "admbase::gzz"
##################################################################
diff --git a/test/test_pw_ADM_sl.par b/test/test_pw_ADM_sl.par
index 31f4963..394787f 100644
--- a/test/test_pw_ADM_sl.par
+++ b/test/test_pw_ADM_sl.par
@@ -3,7 +3,7 @@
#########################################################
# Required thorns
-ActiveThorns = "boundary time pugh pughslab pughreduce cartgrid3d einstein ADM IDLinearWaves ioascii ioutil iobasic"
+ActiveThorns = "boundary time pugh pughslab pughreduce cartgrid3d ADMBase StaticConformal ADM ADMMacros ADMAnalysis CoordGauge IDLinearWaves ioascii ioutil iobasic"
# GENERAL
@@ -22,12 +22,15 @@ cactus::cctk_itlast = 100
# Einstein
-einstein::initial_data = "planewaves"
-einstein::evolution_system = "ADM"
+admbase::metric_type = "static conformal"
+admbase::initial_data = "planewaves"
+admbase::evolution_method = "ADM"
adm::method = "stagleap"
adm::bound = "flat"
-einstein::slicing = "geodesic"
+
+admbase::lapse_evolution_method = "coordgauge"
+coordgauge::lapse_list = "geodesic"
# IDLinearwaves
idlinearwaves::amplitude = 0.001
@@ -45,15 +48,15 @@ IO::out_fileinfo = "none"
IO::parfile_write = "no"
IOBasic::outScalar_every = 10
-IOBasic::outScalar_vars = "einstein::gxx einstein::kxx
- einstein::grr einstein::alp"
+IOBasic::outScalar_vars = "admbase::gxx admbase::kxx
+ admanalysis::grr admbase::alp"
IOASCII::out1D_every = 10
-IOASCII::out1D_vars = "einstein::gxx einstein::gxy einstein::gxz
- einstein::gyy einstein::gyz einstein::gzz
- einstein::kxx einstein::kxy einstein::kxz
- einstein::kyy einstein::kyz einstein::kzz"
+IOASCII::out1D_vars = "admbase::gxx admbase::gxy admbase::gxz
+ admbase::gyy admbase::gyz admbase::gzz
+ admbase::kxx admbase::kxy admbase::kxz
+ admbase::kyy admbase::kyz admbase::kzz"
IOBasic::outInfo_every = 5
-IOBasic::outInfo_vars = "einstein::gxx"
+IOBasic::outInfo_vars = "admbase::gxx"
##################################################################