aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authoreschnett <eschnett@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2011-01-21 03:14:15 +0000
committereschnett <eschnett@578cdeb0-5ea1-4b81-8215-5a3b8777ee0b>2011-01-21 03:14:15 +0000
commit61a33fbee9c1b55997fe7f2e0b8c9c25f4bca9ba (patch)
treedcaa6ae57505eccbf0049355bc0ba570e865346d
parentb530d68aa224e23e1ae0b0629faa680678d99767 (diff)
Correct signatures of scheduled functions
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/MoL/trunk@144 578cdeb0-5ea1-4b81-8215-5a3b8777ee0b
-rw-r--r--src/Counter.c14
-rw-r--r--src/IndexArrays.c4
-rw-r--r--src/ParamCheck.c7
-rw-r--r--src/RHSNaNCheck.c6
-rw-r--r--src/RKCoefficients.c7
-rw-r--r--src/SetTime.c22
6 files changed, 19 insertions, 41 deletions
diff --git a/src/Counter.c b/src/Counter.c
index 486b024..44515de 100644
--- a/src/Counter.c
+++ b/src/Counter.c
@@ -30,9 +30,9 @@ CCTK_FILEVERSION(CactusBase_MoL_Counter_c);
***************** Scheduled Routine Prototypes *********************
********************************************************************/
-int MoL_SetCounter(CCTK_ARGUMENTS);
+void MoL_SetCounter(CCTK_ARGUMENTS);
-int MoL_DecrementCounter(CCTK_ARGUMENTS);
+void MoL_DecrementCounter(CCTK_ARGUMENTS);
/********************************************************************
********************* Other Routine Prototypes *********************
@@ -61,7 +61,7 @@ int MoL_DecrementCounter(CCTK_ARGUMENTS);
@@*/
-int MoL_SetCounter(CCTK_ARGUMENTS)
+void MoL_SetCounter(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -89,10 +89,6 @@ int MoL_SetCounter(CCTK_ARGUMENTS)
}
}
/* #endif */
-
-
-
- return 0;
}
/*@@
@@ -110,7 +106,7 @@ int MoL_SetCounter(CCTK_ARGUMENTS)
@@*/
-int MoL_DecrementCounter(CCTK_ARGUMENTS)
+void MoL_DecrementCounter(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -137,8 +133,6 @@ int MoL_DecrementCounter(CCTK_ARGUMENTS)
}
}
/* #endif */
-
- return 0;
}
/********************************************************************
diff --git a/src/IndexArrays.c b/src/IndexArrays.c
index 3c7c1e0..d3ceccb 100644
--- a/src/IndexArrays.c
+++ b/src/IndexArrays.c
@@ -36,7 +36,7 @@ CCTK_FILEVERSION(CactusBase_MoL_IndexArrays_c);
void MoL_SetupIndexArrays(CCTK_ARGUMENTS);
-void MoL_FreeIndexArrays(void);
+void MoL_FreeIndexArrays(CCTK_ARGUMENTS);
/********************************************************************
********************* Other Routine Prototypes *********************
@@ -366,7 +366,7 @@ void MoL_SetupIndexArrays(CCTK_ARGUMENTS)
@@*/
-void MoL_FreeIndexArrays(void)
+void MoL_FreeIndexArrays(CCTK_ARGUMENTS)
{
if (EvolvedVariableIndex)
diff --git a/src/ParamCheck.c b/src/ParamCheck.c
index 301a252..0257eb2 100644
--- a/src/ParamCheck.c
+++ b/src/ParamCheck.c
@@ -31,7 +31,7 @@ CCTK_FILEVERSION(CactusBase_MoL_ParamCheck_c);
***************** Scheduled Routine Prototypes *********************
********************************************************************/
-int MoL_ParamCheck(CCTK_ARGUMENTS);
+void MoL_ParamCheck(CCTK_ARGUMENTS);
/********************************************************************
********************* Other Routine Prototypes *********************
@@ -60,7 +60,7 @@ int MoL_ParamCheck(CCTK_ARGUMENTS);
@@*/
-int MoL_ParamCheck(CCTK_ARGUMENTS)
+void MoL_ParamCheck(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
@@ -171,9 +171,6 @@ int MoL_ParamCheck(CCTK_ARGUMENTS)
CCTK_PARAMWARN("Adaptive time step sizes are only possible with the RK45, RK45CK, RK65, and RK87 solvers");
}
}
-
- return 0;
-
}
/********************************************************************
diff --git a/src/RHSNaNCheck.c b/src/RHSNaNCheck.c
index d149d4f..ebdfca7 100644
--- a/src/RHSNaNCheck.c
+++ b/src/RHSNaNCheck.c
@@ -32,7 +32,7 @@ CCTK_FILEVERSION(CactusBase_MoL_RHSNaNCheck_c);
***************** Scheduled Routine Prototypes *********************
********************************************************************/
-int MoL_NaNCheck(CCTK_ARGUMENTS);
+void MoL_NaNCheck(CCTK_ARGUMENTS);
/********************************************************************
********************* Other Routine Prototypes *********************
@@ -61,7 +61,7 @@ int MoL_NaNCheck(CCTK_ARGUMENTS);
@@*/
-int MoL_NaNCheck(CCTK_ARGUMENTS)
+void MoL_NaNCheck(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -92,8 +92,6 @@ int MoL_NaNCheck(CCTK_ARGUMENTS)
"NaNs were found on iteration %d inside MoL",
(int)(MoL_Intermediate_Steps - *MoL_Intermediate_Step + 1));
}
-
- return 0;
}
/********************************************************************
diff --git a/src/RKCoefficients.c b/src/RKCoefficients.c
index 57ac38a..36299c7 100644
--- a/src/RKCoefficients.c
+++ b/src/RKCoefficients.c
@@ -33,7 +33,7 @@ CCTK_FILEVERSION(CactusBase_MoL_RKCoefficients_c);
***************** Scheduled Routine Prototypes *********************
********************************************************************/
-int MoL_SetupRKCoefficients(CCTK_ARGUMENTS);
+void MoL_SetupRKCoefficients(CCTK_ARGUMENTS);
/********************************************************************
********************* Other Routine Prototypes *********************
@@ -64,7 +64,7 @@ int MoL_SetupRKCoefficients(CCTK_ARGUMENTS);
@@*/
-int MoL_SetupRKCoefficients(CCTK_ARGUMENTS)
+void MoL_SetupRKCoefficients(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
@@ -245,9 +245,6 @@ int MoL_SetupRKCoefficients(CCTK_ARGUMENTS)
CCTK_WARN(0, "RKCoefficients does not recognize the value "
"of Generic_Type");
}
-
- return 0;
-
}
/********************************************************************
diff --git a/src/SetTime.c b/src/SetTime.c
index 5866d0d..891bb19 100644
--- a/src/SetTime.c
+++ b/src/SetTime.c
@@ -33,11 +33,11 @@ CCTK_FILEVERSION(CactusBase_MoL_SetTime_c);
***************** Scheduled Routine Prototypes *********************
********************************************************************/
-int MoL_SetTime(CCTK_ARGUMENTS);
+void MoL_SetTime(CCTK_ARGUMENTS);
-int MoL_ResetTime(CCTK_ARGUMENTS);
+void MoL_ResetTime(CCTK_ARGUMENTS);
-int MoL_ResetDeltaTime(CCTK_ARGUMENTS);
+void MoL_ResetDeltaTime(CCTK_ARGUMENTS);
/********************************************************************
********************* Other Routine Prototypes *********************
@@ -116,15 +116,13 @@ static const CCTK_REAL alpha_array65[8] = {
@@*/
-int MoL_SetInitialTime(CCTK_ARGUMENTS)
+void MoL_SetInitialTime(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
*Original_Time = cctkGH->cctk_time;
*Original_Delta_Time = cctkGH->cctk_delta_time;
-
- return 0;
}
/*@@
@@ -142,7 +140,7 @@ int MoL_SetInitialTime(CCTK_ARGUMENTS)
@@*/
-int MoL_SetTime(CCTK_ARGUMENTS)
+void MoL_SetTime(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -168,8 +166,6 @@ int MoL_SetTime(CCTK_ARGUMENTS)
cctkGH->cctk_delta_time = beta*(*Original_Delta_Time);
}
-
- return 0;
}
/*@@
@@ -189,7 +185,7 @@ int MoL_SetTime(CCTK_ARGUMENTS)
@@*/
-int MoL_ResetTime(CCTK_ARGUMENTS)
+void MoL_ResetTime(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -335,8 +331,6 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
free(previous_times);
previous_times = NULL;
-
- return 0;
}
/*@@
@@ -356,7 +350,7 @@ int MoL_ResetTime(CCTK_ARGUMENTS)
@@*/
-int MoL_ResetDeltaTime(CCTK_ARGUMENTS)
+void MoL_ResetDeltaTime(CCTK_ARGUMENTS)
{
DECLARE_CCTK_ARGUMENTS;
@@ -445,8 +439,6 @@ int MoL_ResetDeltaTime(CCTK_ARGUMENTS)
cctkGH->cctk_delta_time/cctkGH->cctk_timefac);
fflush(stdout);
#endif
-
- return 0;
}
/********************************************************************