aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorjthorn <jthorn@df1f8a13-aa1d-4dd4-9681-27ded5b42416>2002-05-14 10:46:39 +0000
committerjthorn <jthorn@df1f8a13-aa1d-4dd4-9681-27ded5b42416>2002-05-14 10:46:39 +0000
commit2365beb7250709bd29ec6c7b10548a8798e3e83e (patch)
treeb4f61eb1409284f58f148767a04abb356b0931a8 /src
parent0a060daa1fbbb675a8fee5b6ff909931102c57a1 (diff)
rework pointer casting for storing results of complex interpolations
-- now we use the same scheme as when fetching complex gridfns, namely we pointer-cast our void * pointer to an array of 2-element arrays of real values, then subscript this [pt][part] where pt = the point index and part = 0 for real part and 1 for imaginary part. git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/LocalInterp/trunk@32 df1f8a13-aa1d-4dd4-9681-27ded5b42416
Diffstat (limited to 'src')
-rw-r--r--src/GeneralizedPolynomial-Uniform/template.c33
1 files changed, 13 insertions, 20 deletions
diff --git a/src/GeneralizedPolynomial-Uniform/template.c b/src/GeneralizedPolynomial-Uniform/template.c
index 0e6f004..0943c6b 100644
--- a/src/GeneralizedPolynomial-Uniform/template.c
+++ b/src/GeneralizedPolynomial-Uniform/template.c
@@ -514,7 +514,7 @@ CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
/*
* save origin/delta variables, precompute 1/delta factors
- * ... in theory we could only compute those factors we're going to use,
+ * ... in theory we could compute only those factors we're going to use,
* but it's not worth the trouble, so we just compute them all
*/
{
@@ -552,7 +552,6 @@ CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
/*
* interpolate at each point
*/
- {
int pt;
for (pt = 0 ; pt < N_interp_points ; ++pt)
{
@@ -604,7 +603,6 @@ int pt;
switch (interp_coords_type_code)
{
-
case CCTK_VARIABLE_REAL:
{
const CCTK_REAL *const interp_coords_ptr_real
@@ -1191,20 +1189,18 @@ case CCTK_VARIABLE_REAL16:
case CCTK_VARIABLE_COMPLEX:
{
- CCTK_COMPLEX *const output_array_ptr_complex
- = (CCTK_COMPLEX *) output_arrays[out];
- ((CCTK_REAL *) & output_array_ptr_complex[pt]) [part]
- = result;
+ CCTK_REAL (*const output_array_ptr_complex)[COMPLEX_N_PARTS]
+ = (CCTK_REAL (*)[COMPLEX_N_PARTS]) output_arrays[out];
+ output_array_ptr_complex[pt][part] = result;
break;
}
#ifdef HAVE_CCTK_COMPLEX8
case CCTK_VARIABLE_COMPLEX8:
{
- CCTK_COMPLEX8 *const output_array_ptr_complex8
- = (CCTK_COMPLEX8 *) output_arrays[out];
- ((CCTK_REAL4 *) & output_array_ptr_complex8[pt]) [part]
- = result;
+ CCTK_REAL4 (*const output_array_ptr_complex8)[COMPLEX_N_PARTS]
+ = (CCTK_REAL4 (*)[COMPLEX_N_PARTS]) output_arrays[out];
+ output_array_ptr_complex8[pt][part] = result;
break;
}
#endif /* HAVE_CCTK_COMPLEX8 */
@@ -1212,10 +1208,9 @@ case CCTK_VARIABLE_COMPLEX8:
#ifdef HAVE_CCTK_COMPLEX16
case CCTK_VARIABLE_COMPLEX16:
{
- CCTK_COMPLEX16 *const output_array_ptr_complex16
- = (CCTK_COMPLEX16 *) output_arrays[out];
- ((CCTK_REAL8 *) & output_array_ptr_complex16[pt]) [part]
- = result;
+ CCTK_REAL8 (*const output_array_ptr_complex16)[COMPLEX_N_PARTS]
+ = (CCTK_REAL8 (*)[COMPLEX_N_PARTS]) output_arrays[out];
+ output_array_ptr_complex16[pt][part] = result;
break;
}
#endif /* HAVE_CCTK_COMPLEX16 */
@@ -1223,10 +1218,9 @@ case CCTK_VARIABLE_COMPLEX16:
#ifdef HAVE_CCTK_COMPLEX32
case CCTK_VARIABLE_COMPLEX32:
{
- CCTK_COMPLEX32 *const output_array_ptr_complex32
- = (CCTK_COMPLEX32 *) output_arrays[out];
- ((CCTK_REAL16 *) & output_array_ptr_complex32[pt]) [part]
- = result;
+ CCTK_REAL16 (*const output_array_ptr_complex32)[COMPLEX_N_PARTS]
+ = (CCTK_REAL16 (*)[COMPLEX_N_PARTS]) output_arrays[out];
+ output_array_ptr_complex32[pt][part] = result;
break;
}
#endif /* HAVE_CCTK_COMPLEX32 */
@@ -1254,5 +1248,4 @@ default:
return 0; /*** NORMAL RETURN ***/
}
}
- }
}