aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@b2a53a04-0f4f-0410-87ed-f9f25ced00cf>2004-10-08 19:07:34 +0000
committerschnetter <schnetter@b2a53a04-0f4f-0410-87ed-f9f25ced00cf>2004-10-08 19:07:34 +0000
commit9bfb372b21b2a6a3fcdc34c4fc988a22c2bf101a (patch)
tree9c574555d06c027892ea06ad51caca9ebe95d7b9
parent56edb0f45bea22c4405509665b921476d6461dda (diff)
Use boolean instead of integer for parameter use_sources.
Fix segmentation fault when there is no matter. Add prototype for function relax. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinInitialData/TwoPunctures/trunk@24 b2a53a04-0f4f-0410-87ed-f9f25ced00cf
-rw-r--r--param.ccl7
-rw-r--r--src/FuncAndJacobian.c7
-rw-r--r--src/Newton.c3
-rw-r--r--src/ParamCheck.c10
4 files changed, 14 insertions, 13 deletions
diff --git a/param.ccl b/param.ccl
index 82463f6..076c4ff 100644
--- a/param.ccl
+++ b/param.ccl
@@ -103,8 +103,7 @@ REAL par_S_minus[3] "spin of the m- puncture"
} 0.0
-INT par_use_sources "Use sources?"
-{
- 0:1 :: "0 for no (default), 1 for yes"
-} 0
+BOOLEAN use_sources "Use sources?"
+{
+} "no"
diff --git a/src/FuncAndJacobian.c b/src/FuncAndJacobian.c
index 4f0d13c..e105883 100644
--- a/src/FuncAndJacobian.c
+++ b/src/FuncAndJacobian.c
@@ -192,7 +192,8 @@ void
F_of_v (CCTK_POINTER_TO_CONST cctkGH,
int nvar, int n1, int n2, int n3, derivs v, double *F,
derivs u)
-{ // Calculates the left hand sides of the non-linear equations F_m(v_n)=0
+{
+ // Calculates the left hand sides of the non-linear equations F_m(v_n)=0
// and the function u (u.d0[]) as well as its derivatives
// (u.d1[], u.d2[], u.d3[], u.d11[], u.d12[], u.d13[], u.d22[], u.d23[], u.d33[])
// at interior points and at the boundaries "+/-"
@@ -206,7 +207,7 @@ F_of_v (CCTK_POINTER_TO_CONST cctkGH,
allocate_derivs (&U, nvar);
sources=calloc(n1*n2*n3, sizeof(CCTK_REAL));
- if (par_use_sources)
+ if (use_sources)
{
CCTK_REAL *s_x, *s_y, *s_z;
CCTK_INT i3D;
@@ -259,7 +260,7 @@ F_of_v (CCTK_POINTER_TO_CONST cctkGH,
for (i = 0; i < n1; i++)
for (j = 0; j < n2; j++)
for (k = 0; k < n3; k++)
- sources[i*n1*n2 + j*n2 + k]=0.0;
+ sources[i*n2*n3 + j*n3 + k]=0.0;
Derivatives_AB3 (nvar, n1, n2, n3, v);
diff --git a/src/Newton.c b/src/Newton.c
index 0dbb627..c1af7f3 100644
--- a/src/Newton.c
+++ b/src/Newton.c
@@ -9,6 +9,9 @@
#include "TP_utilities.h"
#include "TwoPunctures.h"
+static void relax (double *dv, int nvar, int n1, int n2, int n3, double *rhs,
+ int *ncols, int **cols, double **JFD);
+
// -----------------------------------------------------------------------------------
void
resid (double *res, int ntotal, double *dv, double *rhs,
diff --git a/src/ParamCheck.c b/src/ParamCheck.c
index 9ec622e..9a13add 100644
--- a/src/ParamCheck.c
+++ b/src/ParamCheck.c
@@ -18,15 +18,13 @@ TwoPunctures_ParamCheck (CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
- if (par_use_sources)
+ if (use_sources)
{
CCTK_INFO("Solving for BH-NS");
- if (CCTK_IsFunctionAliased ("Set_Rho_ADM"))
- CCTK_INFO("Aliased Functions found");
- else
+ if (! CCTK_IsFunctionAliased ("Set_Rho_ADM"))
CCTK_WARN(0, "I found no (aliased) function for matter sources, but "
- "was said to use matter.\n");
+ "was said to use matter.");
}
else
- CCTK_INFO("not using sources (only BHs)");
+ CCTK_INFO("Solving for BHs");
}