aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authortradke <tradke@10716dce-81a3-4424-a2c8-48026a0d3035>2001-06-27 16:53:19 +0000
committertradke <tradke@10716dce-81a3-4424-a2c8-48026a0d3035>2001-06-27 16:53:19 +0000
commit28da073c156187bd544660b0fb6bb68de8d236c3 (patch)
tree17179563a9ab1f87ce4d78999445256411c0a614
parent8a318e7efda20f397db2769f6dc92263dfc65b75 (diff)
Hyperslabs can be done now for every possible CCTK datatype.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusPUGH/PUGHSlab/trunk@58 10716dce-81a3-4424-a2c8-48026a0d3035
-rw-r--r--src/Hyperslab.c95
-rw-r--r--src/NewHyperslab.c15
2 files changed, 77 insertions, 33 deletions
diff --git a/src/Hyperslab.c b/src/Hyperslab.c
index 5f44e78..919e384 100644
--- a/src/Hyperslab.c
+++ b/src/Hyperslab.c
@@ -41,7 +41,7 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_Hyperslab_c)
#define PICKUP_HYPERSLAB_DATA(cctk_type, vdim, vdata, hdata, point, \
startpoint, endpoint, downsample,points_per_dim)\
{ \
- int i, dim, vindex, dim0_elements; \
+ int i, dim, idx, dim0_elements; \
cctk_type *typed_vdata = (cctk_type *) vdata; \
cctk_type *typed_hdata = (cctk_type *) hdata; \
\
@@ -83,10 +83,10 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_Hyperslab_c)
} \
\
/* get the linear offset */ \
- vindex = point[0]; \
+ idx = point[0]; \
for (i = 1; i < vdim; i++) \
{ \
- vindex += point[i] * points_per_dim[i]; \
+ idx += point[i] * points_per_dim[i]; \
} \
\
/* copy the data in lowest dimension */ \
@@ -94,7 +94,7 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_Hyperslab_c)
otherwise copy element-wise in a loop */ \
if (downsample[0] == 1) \
{ \
- memcpy (typed_hdata, typed_vdata + vindex, \
+ memcpy (typed_hdata, typed_vdata + idx, \
dim0_elements * sizeof (*typed_hdata)); \
typed_hdata += dim0_elements; \
} \
@@ -102,7 +102,7 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_Hyperslab_c)
{ \
for (i = 0; i < dim0_elements; i += downsample[0]) \
{ \
- *typed_hdata++ = typed_vdata[vindex + i]; \
+ *typed_hdata++ = typed_vdata[idx + i]; \
} \
} \
\
@@ -491,21 +491,91 @@ int Hyperslab_GetLocalHyperslab (cGH *GH, int vindex, int vtimelvl,
point, startpoint, endpoint, downsample,
points_per_dim);
break;
+
case CCTK_VARIABLE_INT:
PICKUP_HYPERSLAB_DATA (CCTK_INT, vinfo.dim, vdata, *hdata,
point, startpoint, endpoint, downsample,
points_per_dim);
break;
+
case CCTK_VARIABLE_REAL:
PICKUP_HYPERSLAB_DATA (CCTK_REAL, vinfo.dim, vdata, *hdata,
point, startpoint, endpoint, downsample,
points_per_dim);
break;
+
case CCTK_VARIABLE_COMPLEX:
PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX, vinfo.dim, vdata, *hdata,
point, startpoint, endpoint, downsample,
points_per_dim);
break;
+
+#ifdef CCTK_INT2
+ case CCTK_VARIABLE_INT2:
+ PICKUP_HYPERSLAB_DATA (CCTK_INT2, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+#endif
+
+#ifdef CCTK_INT4
+ case CCTK_VARIABLE_INT4:
+ PICKUP_HYPERSLAB_DATA (CCTK_INT4, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+#endif
+
+#ifdef CCTK_INT8
+ case CCTK_VARIABLE_INT8:
+ PICKUP_HYPERSLAB_DATA (CCTK_INT8, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+#endif
+
+#ifdef CCTK_REAL4
+ case CCTK_VARIABLE_REAL4:
+ PICKUP_HYPERSLAB_DATA (CCTK_REAL4, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+
+ case CCTK_VARIABLE_COMPLEX8:
+ PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX8, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+#endif
+
+#ifdef CCTK_REAL8
+ case CCTK_VARIABLE_REAL8:
+ PICKUP_HYPERSLAB_DATA (CCTK_REAL8, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+
+ case CCTK_VARIABLE_COMPLEX16:
+ PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX16, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+#endif
+
+#ifdef CCTK_REAL16
+ case CCTK_VARIABLE_REAL16:
+ PICKUP_HYPERSLAB_DATA (CCTK_REAL16, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+
+ case CCTK_VARIABLE_COMPLEX32:
+ PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX32, vinfo.dim, vdata, *hdata,
+ point, startpoint, endpoint, downsample,
+ points_per_dim);
+ break;
+#endif
+
default:
CCTK_WARN (1, "Unsupported variable type");
retval = -1;
@@ -807,20 +877,7 @@ int Hyperslab_GetHyperslab (cGH *GH, int target_proc, int vindex, int vtimelvl,
}
/* detect the MPI datatype to use */
- mpi_vtype = 0;
- switch (vinfo.vartype)
- {
- case CCTK_VARIABLE_CHAR:
- mpi_vtype = PUGH_MPI_CHAR; break;
- case CCTK_VARIABLE_INT:
- mpi_vtype = PUGH_MPI_INT; break;
- case CCTK_VARIABLE_REAL:
- mpi_vtype = PUGH_MPI_REAL; break;
- case CCTK_VARIABLE_COMPLEX:
- mpi_vtype = PUGH_pGH (GH)->PUGH_mpi_complex; break;
- default:
- CCTK_WARN (1, "Unsupported variable type"); break;
- }
+ mpi_vtype = PUGH_MPIDataType (vinfo.vartype);;
/* collect the hyperslab chunks from all processors */
if (target_proc < 0)
diff --git a/src/NewHyperslab.c b/src/NewHyperslab.c
index 050d745..b0aa21b 100644
--- a/src/NewHyperslab.c
+++ b/src/NewHyperslab.c
@@ -956,20 +956,7 @@ int NewHyperslab_GetHyperslab (cGH *GH,
}
/* detect the MPI datatype to use */
- mpi_vtype = 0;
- switch (htype)
- {
- case CCTK_VARIABLE_CHAR:
- mpi_vtype = PUGH_MPI_CHAR; break;
- case CCTK_VARIABLE_INT:
- mpi_vtype = PUGH_MPI_INT; break;
- case CCTK_VARIABLE_REAL:
- mpi_vtype = PUGH_MPI_REAL; break;
- case CCTK_VARIABLE_COMPLEX:
- mpi_vtype = PUGH_pGH (GH)->PUGH_mpi_complex; break;
- default:
- CCTK_WARN (1, "Unsupported variable type"); break;
- }
+ mpi_vtype = PUGH_MPIDataType (htype);
/* collect the hyperslab chunks from all processors */
if (target_proc < 0)