aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authortradke <tradke@10716dce-81a3-4424-a2c8-48026a0d3035>2001-12-03 22:10:05 +0000
committertradke <tradke@10716dce-81a3-4424-a2c8-48026a0d3035>2001-12-03 22:10:05 +0000
commit9b7e0a9b724e4f9bd973c0b3320293eff7064227 (patch)
treefb7a748728630b27eca0a6a611c986814cde7eb1 /src
parent77d3ed2c0a6d387e06a798b3902ad72cef649dfb (diff)
Some small fixes for the old hyperslab API.
This is meant as preparing for the new API which should be committed soon now. git-svn-id: http://svn.cactuscode.org/arrangements/CactusPUGH/PUGHSlab/trunk@62 10716dce-81a3-4424-a2c8-48026a0d3035
Diffstat (limited to 'src')
-rw-r--r--src/DatatypeConversion.c65
-rw-r--r--src/Hyperslab.c339
-rw-r--r--src/NewHyperslab.c19
-rw-r--r--src/NewPUGHSlab.h19
-rw-r--r--src/PUGHSlab.h65
5 files changed, 300 insertions, 207 deletions
diff --git a/src/DatatypeConversion.c b/src/DatatypeConversion.c
index e80be03..478da9b 100644
--- a/src/DatatypeConversion.c
+++ b/src/DatatypeConversion.c
@@ -9,7 +9,8 @@
@@*/
#include "cctk.h"
-#include "NewPUGHSlab.h"
+#include "PUGHSlab.h"
+#include "PUGHSlabi.h"
/* the rcs ID and its dummy function to use it */
@@ -20,22 +21,26 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_DatatypeConversion_c)
/* macro to generate a predefined conversion function
along with its prototype */
#define CONVERSION_FUNCTION(src_type, dst_type, dst_scalartype, conversion) \
- static void Convert_##src_type##_to_##dst_type (void *_dst, \
- void *_src, \
- unsigned int nelems); \
+ static void Convert_##src_type##_to_##dst_type (const void *src, \
+ void *dst, \
+ CCTK_INT nelems, \
+ CCTK_INT src_stride, \
+ CCTK_INT dst_stride); \
\
- static void Convert_##src_type##_to_##dst_type (void *_dst, \
- void *_src, \
- unsigned int nelems) \
+ static void Convert_##src_type##_to_##dst_type (const void *src, \
+ void *dst, \
+ CCTK_INT nelems, \
+ CCTK_INT src_stride, \
+ CCTK_INT dst_stride) \
{ \
- dst_type *dst = (dst_type *) _dst; \
- src_type *src = (src_type *) _src; \
+ dst_type *_dst = (dst_type *) dst; \
+ const src_type *_src = (const src_type *) src; \
\
\
while (nelems--) \
{ \
- conversion (*src, *dst, dst_scalartype); \
- src++; dst++; \
+ conversion (*_src, *_dst, dst_scalartype); \
+ _src += src_stride; _dst += dst_stride; \
} \
}
@@ -85,12 +90,11 @@ CONVERSION_FUNCTION (CCTK_COMPLEX32, CCTK_COMPLEX8, CCTK_REAL4, CONVERT)
/* prototypes of routines defined in this source file */
-PUGHSlab_conversion_fn PUGHSlab_GetDatatypeConversionFn (int src, int dst);
-static int PUGHSlab_PrecisionVarType (int type);
+static int PUGHSlabi_PrecisionVarType (int type);
/*@@
- @routine PUGHSlab_GetDatatypeConversionFn
+ @routine PUGHSlabi_GetDatatypeConversionFn
@date Fri 23 Nov 2000
@author Thomas Radke
@desc
@@ -99,12 +103,12 @@ static int PUGHSlab_PrecisionVarType (int type);
@calls PUGHSlab_REAL8_to_REAL4
- @var src
+ @var src_type
@vdesc CCTK datatype of the source
@vtype int
@vio in
@endvar
- @var dst
+ @var dst_type
@vdesc CCTK datatype of the destination (hyperslab datatype)
@vtype int
@vio in
@@ -116,14 +120,15 @@ static int PUGHSlab_PrecisionVarType (int type);
or NULL if no predefined conversion function was found
@endreturndesc
@@*/
-PUGHSlab_conversion_fn PUGHSlab_GetDatatypeConversionFn (int src, int dst)
+t_hslabConversionFn PUGHSlabi_GetDatatypeConversionFn (int src_type,
+ int dst_type)
{
- PUGHSlab_conversion_fn retval;
+ t_hslabConversionFn retval;
/* get the precision datatype for generic CCTK_TYPE variable types */
- src = PUGHSlab_PrecisionVarType (src);
- dst = PUGHSlab_PrecisionVarType (dst);
+ src_type = PUGHSlabi_PrecisionVarType (src_type);
+ dst_type = PUGHSlabi_PrecisionVarType (dst_type);
/* get the appropriate conversion routine */
if (0)
@@ -131,55 +136,55 @@ PUGHSlab_conversion_fn PUGHSlab_GetDatatypeConversionFn (int src, int dst)
/* this is just because of the following #ifdef'ed else branches */
}
#if defined(CCTK_INT2) && defined(CCTK_INT4)
- else if (src == CCTK_VARIABLE_INT4 && dst == CCTK_VARIABLE_INT2)
+ else if (src_type == CCTK_VARIABLE_INT4 && dst_type == CCTK_VARIABLE_INT2)
{
retval = Convert_CCTK_INT4_to_CCTK_INT2;
}
#endif
#if defined(CCTK_INT2) && defined(CCTK_INT8)
- else if (src == CCTK_VARIABLE_INT8 && dst == CCTK_VARIABLE_INT2)
+ else if (src_type == CCTK_VARIABLE_INT8 && dst_type == CCTK_VARIABLE_INT2)
{
retval = Convert_CCTK_INT8_to_CCTK_INT2;
}
#endif
#if defined(CCTK_INT4) && defined(CCTK_INT8)
- else if (src == CCTK_VARIABLE_INT8 && dst == CCTK_VARIABLE_INT4)
+ else if (src_type == CCTK_VARIABLE_INT8 && dst_type == CCTK_VARIABLE_INT4)
{
retval = Convert_CCTK_INT8_to_CCTK_INT4;
}
#endif
#if defined(CCTK_REAL4) && defined(CCTK_REAL8)
- else if (src == CCTK_VARIABLE_REAL8 && dst == CCTK_VARIABLE_REAL4)
+ else if (src_type == CCTK_VARIABLE_REAL8 && dst_type == CCTK_VARIABLE_REAL4)
{
retval = Convert_CCTK_REAL8_to_CCTK_REAL4;
}
#endif
#if defined(CCTK_REAL4) && defined(CCTK_REAL16)
- else if (src == CCTK_VARIABLE_REAL16 && dst == CCTK_VARIABLE_REAL4)
+ else if (src_type == CCTK_VARIABLE_REAL16 && dst_type == CCTK_VARIABLE_REAL4)
{
retval = Convert_CCTK_REAL16_to_CCTK_REAL4;
}
#endif
#if defined(CCTK_REAL8) && defined(CCTK_REAL16)
- else if (src == CCTK_VARIABLE_REAL16 && dst == CCTK_VARIABLE_REAL8)
+ else if (src_type == CCTK_VARIABLE_REAL16 && dst_type == CCTK_VARIABLE_REAL8)
{
retval = Convert_CCTK_REAL16_to_CCTK_REAL8;
}
#endif
#if defined(CCTK_REAL4) && defined(CCTK_REAL8)
- else if (src == CCTK_VARIABLE_COMPLEX16 && dst == CCTK_VARIABLE_COMPLEX8)
+ else if (src_type == CCTK_VARIABLE_COMPLEX16 && dst_type == CCTK_VARIABLE_COMPLEX8)
{
retval = Convert_CCTK_COMPLEX16_to_CCTK_COMPLEX8;
}
#endif
#if defined(CCTK_REAL8) && defined(CCTK_REAL16)
- else if (src == CCTK_VARIABLE_COMPLEX32 && dst == CCTK_VARIABLE_COMPLEX8)
+ else if (src_type == CCTK_VARIABLE_COMPLEX32 && dst_type == CCTK_VARIABLE_COMPLEX8)
{
retval = Convert_CCTK_COMPLEX32_to_CCTK_COMPLEX8;
}
#endif
#if defined(CCTK_REAL8) && defined(CCTK_REAL16)
- else if (src == CCTK_VARIABLE_COMPLEX32 && dst == CCTK_VARIABLE_COMPLEX16)
+ else if (src_type == CCTK_VARIABLE_COMPLEX32 && dst_type == CCTK_VARIABLE_COMPLEX16)
{
retval = Convert_CCTK_COMPLEX32_to_CCTK_COMPLEX16;
}
@@ -197,7 +202,7 @@ PUGHSlab_conversion_fn PUGHSlab_GetDatatypeConversionFn (int src, int dst)
/* local functions */
/**************************************************************************/
/* get the precision datatype for generic CCTK_TYPE variable types */
-static int PUGHSlab_PrecisionVarType (int type)
+static int PUGHSlabi_PrecisionVarType (int type)
{
if (type == CCTK_VARIABLE_INT)
{
diff --git a/src/Hyperslab.c b/src/Hyperslab.c
index 6342fcb..3a8603a 100644
--- a/src/Hyperslab.c
+++ b/src/Hyperslab.c
@@ -1,4 +1,14 @@
-#include <stdio.h>
+ /*@@
+ @file Hyperslab.c
+ @date Thu 23 Nov 2000
+ @author Thomas Radke
+ @desc
+ Routines to extract hyperslabs from CCTK array variables
+ @enddesc
+ @version $Id$
+ @@*/
+
+
#include <stdlib.h>
#include <string.h>
@@ -36,20 +46,20 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_Hyperslab_c)
/* macro for copying all relevant data points
- (ranging from startpoint[] to endpoint[])
+ (ranging from origin[] to endpoint[])
into the typed hyperslab data buffer */
#define PICKUP_HYPERSLAB_DATA(cctk_type, vdim, vdata, hdata, point, \
- startpoint, endpoint, downsample,points_per_dim)\
+ origin, endpoint, downsample,points_per_dim)\
{ \
int i, dim, idx, dim0_elements; \
cctk_type *typed_vdata = (cctk_type *) vdata; \
cctk_type *typed_hdata = (cctk_type *) hdata; \
\
\
- /* set point to local startpoint */ \
- memcpy (point, startpoint, vdim * sizeof (point[0])); \
+ /* set point to local origin */ \
+ memcpy (point, origin, vdim * sizeof (point[0])); \
\
- dim0_elements = endpoint[0] - startpoint[0]; \
+ dim0_elements = endpoint[0] - origin[0]; \
\
/* do the nested loops starting with the innermost */ \
dim = 1; \
@@ -77,7 +87,7 @@ CCTK_FILEVERSION(CactusPUGH_PUGHSlab_Hyperslab_c)
/* reset innermost loopers */ \
for (dim--; dim > 0; dim--) \
{ \
- point[dim] = startpoint[dim]; \
+ point[dim] = origin[dim]; \
} \
dim = 1; \
} \
@@ -134,12 +144,22 @@ int Hyperslab_CollectData1D (const cGH *GH, int vindex, int vtimelvl,
/* routine to check parameters passed to the Hyperslab routines */
static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
int hdim,
- const int global_startpoint[/*vdim*/],
+ const int global_origin[/*vdim*/],
const int directions[/*vdim*/],
const int extents[/*hdim*/],
const int downsample_[/*hdim*/],
void **hdata,
int hsize[/*hdim*/]);
+#ifdef CCTK_MPI
+static void SortIntoUnchunked (int hdim,
+ const int hsize[],
+ int totals_global,
+ const CCTK_INT sizes_global[],
+ void **hdata,
+ void *chunked_hdata,
+ int vtypesize,
+ int nprocs);
+#endif /* CCTK_MPI */
/*@@
@@ -181,7 +201,7 @@ static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
@vtype int
@vio in
@endvar
- @var global_startpoint
+ @var global_origin
@vdesc global coordinates of the hyperslab origin
@vtype const int[dimensions of vindex]
@vio in
@@ -230,9 +250,11 @@ static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
@vio out
@endvar
@@*/
-int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
+int Hyperslab_GetLocalHyperslab (const cGH *GH,
+ int vindex,
+ int vtimelvl,
int hdim,
- const int global_startpoint[/*vdim*/],
+ const int global_origin[/*vdim*/],
const int directions[/*vdim*/],
const int extents[/*hdim*/],
const int downsample_[/*hdim*/],
@@ -242,10 +264,10 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
int hoffset_global[/*hdim*/])
{
int *point; /* looper over hyperslab dimensions */
- int *startpoint, /* hyperslab's local start and endpoint */
+ int *origin, /* hyperslab's local start and endpoint */
*endpoint; /* within the variable's grid dimensions */
int *downsample; /* the downsample_[] vector extended to vdim */
- int *my_global_startpoint, /* hyperslab's global start and endpoint */
+ int *my_global_origin, /* hyperslab's global start and endpoint */
*global_endpoint;
int *points_per_dim; /* points per subvolume */
int *do_dir; /* directions in which to span the hyperslab */
@@ -262,7 +284,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
/* do some plausibility checks */
errormsg = checkParameters (GH, vindex, vtimelvl, hdim,
- global_startpoint, directions, extents,
+ global_origin, directions, extents,
downsample_, hdata, hsize);
/* immediately return in case of errors */
@@ -280,9 +302,9 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
/* allocate the temporary arrays */
point = (int *) malloc (8 * vinfo.dim * sizeof (int));
- startpoint = point + 1*vinfo.dim;
+ origin = point + 1*vinfo.dim;
endpoint = point + 2*vinfo.dim;
- my_global_startpoint = point + 3*vinfo.dim;
+ my_global_origin = point + 3*vinfo.dim;
global_endpoint = point + 4*vinfo.dim;
downsample = point + 5*vinfo.dim;
do_dir = point + 6*vinfo.dim;
@@ -322,7 +344,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
if (do_dir[vdim])
{
global_endpoint[vdim] = extents[hdim] > 0 ?
- MIN (global_startpoint[vdim] + extents[hdim],
+ MIN (global_origin[vdim] + extents[hdim],
GA->extras->nsize[vdim]) :
GA->extras->nsize[vdim];
downsample[vdim] = downsample_[hdim];
@@ -330,7 +352,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
}
else
{
- global_endpoint[vdim] = global_startpoint[vdim] + 1;
+ global_endpoint[vdim] = global_origin[vdim] + 1;
downsample[vdim] = 1;
}
}
@@ -338,37 +360,37 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
/* get the local processor ID */
myproc = CCTK_MyProc (GH);
- /* compute my global startpoint from the global ranges */
+ /* compute my global origin from the global ranges */
for (vdim = 0; vdim < vinfo.dim; vdim++)
{
stagger_index = CCTK_StaggerDirIndex (vdim, vinfo.stagtype);
- if (global_startpoint[vdim] < MY_GLOBAL_EP (GA->extras, myproc,
+ if (global_origin[vdim] < MY_GLOBAL_EP (GA->extras, myproc,
stagger_index, vdim))
{
- if (global_startpoint[vdim] < MY_GLOBAL_SP (GA->extras, myproc,
+ if (global_origin[vdim] < MY_GLOBAL_SP (GA->extras, myproc,
stagger_index, vdim))
{
int npoints;
npoints = (MY_GLOBAL_SP (GA->extras, myproc, stagger_index, vdim)
- - global_startpoint[vdim]) / downsample[vdim];
+ - global_origin[vdim]) / downsample[vdim];
if ((MY_GLOBAL_SP (GA->extras, myproc, stagger_index, vdim)
- - global_startpoint[vdim]) % downsample[vdim])
+ - global_origin[vdim]) % downsample[vdim])
{
npoints++;
}
- my_global_startpoint[vdim] = global_startpoint[vdim] +
- npoints*downsample[vdim];
+ my_global_origin[vdim] = global_origin[vdim] +
+ npoints*downsample[vdim];
}
else
{
- my_global_startpoint[vdim] = global_startpoint[vdim];
+ my_global_origin[vdim] = global_origin[vdim];
}
}
else
{
- my_global_startpoint[vdim] = -1;
+ my_global_origin[vdim] = -1;
}
}
@@ -378,16 +400,16 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
{
stagger_index = CCTK_StaggerDirIndex (vdim, vinfo.stagtype);
- if (my_global_startpoint[vdim] >= 0 &&
- my_global_startpoint[vdim] < MY_GLOBAL_EP (GA->extras, myproc,
- stagger_index, vdim))
+ if (my_global_origin[vdim] >= 0 &&
+ my_global_origin[vdim] < MY_GLOBAL_EP (GA->extras, myproc,
+ stagger_index, vdim))
{
- startpoint[vdim] = my_global_startpoint[vdim] -
- GA->extras->lb[myproc][vdim];
+ origin[vdim] = my_global_origin[vdim] -
+ GA->extras->lb[myproc][vdim];
}
else
{
- startpoint[vdim] = -1;
+ origin[vdim] = -1;
}
if (global_endpoint[vdim] > MY_GLOBAL_SP (GA->extras, myproc,
@@ -403,21 +425,21 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#ifdef DEBUG
printf ("direction %d: local ranges[%d, %d)\n",
- vdim, startpoint[vdim], endpoint[vdim]);
+ vdim, origin[vdim], endpoint[vdim]);
#endif
- if (endpoint[vdim] < 0 || startpoint[vdim] < 0)
+ if (endpoint[vdim] < 0 || origin[vdim] < 0)
{
totals = 0;
- endpoint[vdim] = startpoint[vdim];
+ endpoint[vdim] = origin[vdim];
}
if (do_dir[vdim])
{
/* compute the global and local size in each hyperslab dimension */
- hsize_global[hdim] = (global_endpoint[vdim] - global_startpoint[vdim]) /
+ hsize_global[hdim] = (global_endpoint[vdim] - global_origin[vdim]) /
downsample[vdim];
- if ((global_endpoint[vdim] - global_startpoint[vdim]) %
+ if ((global_endpoint[vdim] - global_origin[vdim]) %
downsample[vdim])
{
hsize_global[hdim]++;
@@ -426,8 +448,8 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
{
hsize_global[hdim] -= 2 * GA->extras->nghostzones[vdim];
}
- hsize[hdim] = (endpoint[vdim] - startpoint[vdim]) / downsample[vdim];
- if ((endpoint[vdim] - startpoint[vdim]) % downsample[vdim])
+ hsize[hdim] = (endpoint[vdim] - origin[vdim]) / downsample[vdim];
+ if ((endpoint[vdim] - origin[vdim]) % downsample[vdim])
{
hsize[hdim]++;
}
@@ -441,7 +463,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#endif
/* nested loop over vinfo.dim dimensions */
- /* NOTE: the following code assumes startpoint[vdim] < endpoint[vdim] */
+ /* NOTE: the following code assumes origin[vdim] < endpoint[vdim] */
if (totals > 0)
{
void *vdata = CCTK_VarDataPtrI (GH, vtimelvl, vindex);
@@ -454,8 +476,8 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
{
if (do_dir[vdim])
{
- hoffset_global[hdim] = (my_global_startpoint[vdim] -
- global_startpoint[vdim]) / downsample[vdim];
+ hoffset_global[hdim] = (my_global_origin[vdim] -
+ global_origin[vdim]) / downsample[vdim];
if (GA->connectivity->perme[vdim])
{
hoffset_global[hdim] -= GA->extras->nghostzones[vdim];
@@ -487,32 +509,32 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
{
case CCTK_VARIABLE_CHAR:
PICKUP_HYPERSLAB_DATA (CCTK_BYTE, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
case CCTK_VARIABLE_INT:
PICKUP_HYPERSLAB_DATA (CCTK_INT, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
case CCTK_VARIABLE_REAL:
PICKUP_HYPERSLAB_DATA (CCTK_REAL, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
case CCTK_VARIABLE_COMPLEX:
PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, 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,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
#endif
@@ -520,7 +542,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#ifdef CCTK_INT4
case CCTK_VARIABLE_INT4:
PICKUP_HYPERSLAB_DATA (CCTK_INT4, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
#endif
@@ -528,7 +550,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#ifdef CCTK_INT8
case CCTK_VARIABLE_INT8:
PICKUP_HYPERSLAB_DATA (CCTK_INT8, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
#endif
@@ -536,13 +558,13 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#ifdef CCTK_REAL4
case CCTK_VARIABLE_REAL4:
PICKUP_HYPERSLAB_DATA (CCTK_REAL4, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
case CCTK_VARIABLE_COMPLEX8:
PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX8, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
#endif
@@ -550,13 +572,13 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#ifdef CCTK_REAL8
case CCTK_VARIABLE_REAL8:
PICKUP_HYPERSLAB_DATA (CCTK_REAL8, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
case CCTK_VARIABLE_COMPLEX16:
PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX16, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
#endif
@@ -564,13 +586,13 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
#ifdef CCTK_REAL16
case CCTK_VARIABLE_REAL16:
PICKUP_HYPERSLAB_DATA (CCTK_REAL16, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
case CCTK_VARIABLE_COMPLEX32:
PICKUP_HYPERSLAB_DATA (CCTK_COMPLEX32, vinfo.dim, vdata, *hdata,
- point, startpoint, endpoint, downsample,
+ point, origin, endpoint, downsample,
points_per_dim);
break;
#endif
@@ -629,7 +651,7 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
@vtype int
@vio in
@endvar
- @var global_startpoint
+ @var global_origin
@vdesc global coordinates of the hyperslab origin
@vtype const int[dimensions of vindex]
@vio in
@@ -667,9 +689,12 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
@vio out
@endvar
@@*/
-int Hyperslab_GetHyperslab (const cGH *GH, int target_proc, int vindex, int vtimelvl,
+int Hyperslab_GetHyperslab (const cGH *GH,
+ int target_proc,
+ int vindex,
+ int vtimelvl,
int hdim,
- const int global_startpoint[/*vdim*/],
+ const int global_origin[/*vdim*/],
const int directions[/*vdim*/],
const int extents[/*hdim*/],
const int downsample_[/*hdim*/],
@@ -701,7 +726,7 @@ int Hyperslab_GetHyperslab (const cGH *GH, int target_proc, int vindex, int vtim
/* do some plausibility checks */
errormsg = checkParameters (GH, vindex, vtimelvl, hdim,
- global_startpoint, directions, extents,
+ global_origin, directions, extents,
downsample_, hdata, hsize);
/* FIXME: hack for getting diagonals from 3D variables
@@ -727,7 +752,7 @@ int Hyperslab_GetHyperslab (const cGH *GH, int target_proc, int vindex, int vtim
{
length--;
}
- return (Hyperslab_CollectData1D (GH, vindex, vtimelvl, global_startpoint,
+ return (Hyperslab_CollectData1D (GH, vindex, vtimelvl, global_origin,
directions, downsample_[0], length,
hdata, hsize, target_proc));
}
@@ -782,7 +807,7 @@ int Hyperslab_GetHyperslab (const cGH *GH, int target_proc, int vindex, int vtim
/* get the processor-local hyperslab */
retval = Hyperslab_GetLocalHyperslab (GH, vindex, vtimelvl, hdim,
- global_startpoint, directions, extents,
+ global_origin, directions, extents,
downsample_, hdata_ptr, hsize_local,
hsize_global, hoffset_local);
@@ -901,124 +926,140 @@ int Hyperslab_GetHyperslab (const cGH *GH, int target_proc, int vindex, int vtim
The user wants it unchunked, so let's sort it... */
if (target_proc < 0 || target_proc == myproc)
{
- /* for hdim == 1 the hyperslab was already sorted by MPI_Gatherv()
- we just need to return that buffer */
- if (hdim == 1)
+ SortIntoUnchunked (hdim, hsize, totals_global, sizes_global, hdata, chunked_hdata, vtypesize, nprocs);
+
+ /* free allocated vectors for MPI communication */
+ free (displs);
+ }
+
+ /* free allocated arrays for communicating the hyperslab sizes of offsets */
+ free (sizes_local);
+ free (hsize_local);
+
+#endif /* CCTK_MPI */
+
+ return (retval);
+}
+
+
+/********************** local routines ************************************/
+#ifdef CCTK_MPI
+static void SortIntoUnchunked (int hdim,
+ const int hsize[],
+ int totals_global,
+ const CCTK_INT sizes_global[],
+ void **hdata,
+ void *chunked_hdata,
+ int vtypesize,
+ int nprocs)
+{
+ int i, j, proc;
+ int *point;
+ int *points_per_hdim;
+ char *copy_from, *copy_to;
+ const CCTK_INT *hsize_chunk, *hoffset_chunk;
+ int linear_hoffset;
+
+
+ /* for hdim == 1 the hyperslab was already sorted by MPI_Gatherv()
+ we just need to return that buffer */
+ if (hdim == 1)
+ {
+ *hdata = chunked_hdata;
+ }
+ else
+ {
+ /* allocate temporary vectors */
+ point = (int *) malloc (2 * hdim * sizeof (int));
+ points_per_hdim = point + hdim;
+
+ points_per_hdim[0] = 1;
+ for (i = 1; i < hdim; i++)
{
- *hdata = chunked_hdata;
+ points_per_hdim[i] = points_per_hdim[i-1] * hsize[i-1];
}
- else
- {
- int j;
- int *point;
- int *points_per_hdim;
- char *copy_from, *copy_to;
- CCTK_INT *hsize_chunk, *hoffset_chunk;
- int linear_hoffset;
+ /* allocate buffer for the returned hyperslab */
+ *hdata = malloc (totals_global * vtypesize);
- /* allocate temporary vectors */
- point = (int *) malloc (2 * hdim * sizeof (int));
- points_per_hdim = point + hdim;
+ /* use char pointers for easy incrementing when copying */
+ copy_from = (char *) chunked_hdata;
+ copy_to = (char *) *hdata;
- points_per_hdim[0] = 1;
- for (i = 1; i < hdim; i++)
+ /* now copy the chunks from each processor into the global hyperslab */
+ for (proc = 0; proc < nprocs; proc++)
+ {
+ /* skip processors which didn't contribute any data */
+ if (sizes_global[proc * (2*hdim + 1) + 2*hdim] <= 0)
{
- points_per_hdim[i] = points_per_hdim[i-1] * hsize[i-1];
+ continue;
}
- /* allocate buffer for the returned hyperslab */
- *hdata = malloc (totals_global * vtypesize);
+ hsize_chunk = sizes_global + proc * (2*hdim+1);
+ hoffset_chunk = hsize_chunk + hdim;
- /* use char pointers for easy incrementing when copying */
- copy_from = (char *) chunked_hdata;
- copy_to = (char *) *hdata;
+ /* set origin to zero */
+ memset (point, 0, hdim * sizeof (int));
- /* now copy the chunks from each processor into the global hyperslab */
- for (proc = 0; proc < nprocs; proc++)
+ i = 1;
+ while (1)
{
- /* skip processors which didn't contribute any data */
- if (sizes_global[proc * (2*hdim + 1) + 2*hdim] <= 0)
+ /* check for end of current loop */
+ if (point[i] >= hsize_chunk[i])
{
- continue;
- }
-
- hsize_chunk = sizes_global + proc * (2*hdim+1);
- hoffset_chunk = hsize_chunk + hdim;
-
- /* set startpoint to zero */
- memset (point, 0, hdim * sizeof (int));
-
- i = 1;
- while (1)
- {
- /* check for end of current loop */
- if (point[i] >= hsize_chunk[i])
+ /* increment outermost loopers */
+ for (i++; i < hdim; i++)
{
- /* increment outermost loopers */
- for (i++; i < hdim; i++)
- {
- if (++point[i] < hsize_chunk[i])
- {
- break;
- }
- }
-
- /* done if beyond outermost loop */
- if (i >= hdim)
+ if (++point[i] < hsize_chunk[i])
{
break;
}
+ }
- /* reset innermost loopers */
- for (i--; i > 0; i--)
- {
- point[i] = 0;
- }
- i = 1;
+ /* done if beyond outermost loop */
+ if (i >= hdim)
+ {
+ break;
}
- /* get the linear offset */
- linear_hoffset = hoffset_chunk[0];
- for (j = 1; j < hdim; j++)
+ /* reset innermost loopers */
+ for (i--; i > 0; i--)
{
- linear_hoffset += (hoffset_chunk[j] + point[j]) *
- points_per_hdim[j];
+ point[i] = 0;
}
- /* copy the line */
- memcpy (copy_to + linear_hoffset * vtypesize,
- copy_from, hsize_chunk[0] * vtypesize);
- copy_from += hsize_chunk[0] * vtypesize;
+ i = 1;
+ }
- /* increment current looper */
- point[i]++;
+ /* get the linear offset */
+ linear_hoffset = hoffset_chunk[0];
+ for (j = 1; j < hdim; j++)
+ {
+ linear_hoffset += (hoffset_chunk[j] + point[j]) *
+ points_per_hdim[j];
+ }
+ /* copy the line */
+ memcpy (copy_to + linear_hoffset * vtypesize,
+ copy_from, hsize_chunk[0] * vtypesize);
+ copy_from += hsize_chunk[0] * vtypesize;
- } /* end of nested loops over all dimensions */
- } /* end of loop over all processors */
+ /* increment current looper */
+ point[i]++;
- /* free allocated temporary vectors and the chunk buffer */
- free (chunked_hdata);
- free (point);
- }
+ } /* end of nested loops over all dimensions */
+ } /* end of loop over all processors */
- /* free allocated vectors for MPI communication */
- free (displs);
+ /* free allocated temporary vectors and the chunk buffer */
+ free (chunked_hdata);
+ free (point);
}
-
- /* free allocated arrays for communicating the hyperslab sizes of offsets */
- free (sizes_local);
- free (hsize_local);
-
-#endif /* CCTK_MPI */
-
- return (retval);
}
+#endif /* CCTK_MPI */
/********************** local routines ************************************/
static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
int hdim,
- const int global_startpoint[/*vdim*/],
+ const int global_origin[/*vdim*/],
const int directions[/*vdim*/],
const int extents[/*hdim*/],
const int downsample_[/*hdim*/],
@@ -1041,7 +1082,7 @@ static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
}
/* check the passed pointers */
- if (! global_startpoint || ! directions || ! extents || ! downsample_ ||
+ if (! global_origin || ! directions || ! extents || ! downsample_ ||
! hdata || ! hsize)
{
return ("NULL pointer(s) passed as parameters");
diff --git a/src/NewHyperslab.c b/src/NewHyperslab.c
index c92ac13..37aeeb0 100644
--- a/src/NewHyperslab.c
+++ b/src/NewHyperslab.c
@@ -17,6 +17,7 @@
#include "CactusPUGH/PUGH/src/include/pugh.h"
#include "NewPUGHSlab.h"
+#include "PUGHSlabi.h"
/* the rcs ID and its dummy function to use it */
static const char *rcsid = "$Header$";
@@ -119,7 +120,7 @@ static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
local hyperslab data and output it in parallel.
@enddesc
- @calls PUGHSlab_GetDatatypeConversionFn
+ @calls PUGHSlabi_GetDatatypeConversionFn
@var GH
@vdesc Pointer to CCTK grid hierarchy
@@ -148,7 +149,7 @@ static const char *checkParameters (const cGH *GH, int vindex, int vtimelvl,
@endvar
@var conversion_fn
@vdesc pointer to a user-supplied data conversion function
- @vtype PUGHSlab_conversion_fn
+ @vtype t_hslabConversionFn
@vio in
@endvar
@var global_startpoint
@@ -206,7 +207,7 @@ int NewHyperslab_GetLocalHyperslab (const cGH *GH,
int vtimelvl,
int hdim,
int htype,
- PUGHSlab_conversion_fn conversion_fn,
+ t_hslabConversionFn conversion_fn,
const int global_startpoint[/* vdim */],
const int directions[/* hdim*vdim */],
const int extents[/* hdim */],
@@ -277,7 +278,7 @@ int NewHyperslab_GetLocalHyperslab (const cGH *GH,
{
if (conversion_fn == NULL)
{
- conversion_fn = PUGHSlab_GetDatatypeConversionFn (vinfo.vartype, htype);
+ conversion_fn = PUGHSlabi_GetDatatypeConversionFn (vinfo.vartype, htype);
if (! conversion_fn)
{
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
@@ -589,7 +590,7 @@ int NewHyperslab_GetLocalHyperslab (const cGH *GH,
{
if (conversion_fn)
{
- conversion_fn (typed_hdata, typed_vdata, dim0_points);
+ conversion_fn (typed_vdata, typed_hdata, dim0_points, 1, 1);
}
else
{
@@ -600,11 +601,9 @@ int NewHyperslab_GetLocalHyperslab (const cGH *GH,
{
if (conversion_fn)
{
- for (i = 0; i < dim0_hsize; i += hdata_size)
- {
- (*conversion_fn) (typed_hdata + i, typed_vdata, 1);
- typed_vdata += downsample[0];
- }
+ conversion_fn (typed_vdata, typed_hdata, dim0_points,
+ downsample[0], 1);
+ typed_vdata += downsample[0] * dim0_points;
}
else
{
diff --git a/src/NewPUGHSlab.h b/src/NewPUGHSlab.h
index 7fcffcb..3d76628 100644
--- a/src/NewPUGHSlab.h
+++ b/src/NewPUGHSlab.h
@@ -1,5 +1,5 @@
/*@@
- @header PUGHSlab.h
+ @header NewPUGHSlab.h
@date Sun 28 May 2000
@author Thomas Radke
@desc
@@ -10,19 +10,16 @@
@@*/
-#ifndef _PUGHSLAB_PUGHSLAB_H_
-#define _PUGHSLAB_PUGHSLAB_H_
+#ifndef _PUGHSLAB_NEWPUGHSLAB_H_
+#define _PUGHSLAB_NEWPUGHSLAB_H_
+
+#include "PUGHSlab.h"
#ifdef __cplusplus
extern "C"
{
#endif
-/* prototype of datatype conversion routines */
-typedef void (*PUGHSlab_conversion_fn) (void *hdata,
- void *vdata,
- unsigned int nelems);
-
/* function prototypes */
int NewHyperslab_GetLocalHyperslab (const cGH *GH,
@@ -30,7 +27,7 @@ int NewHyperslab_GetLocalHyperslab (const cGH *GH,
int vtimelvl,
int hdim,
int htype,
- PUGHSlab_conversion_fn copy_fn,
+ t_hslabConversionFn copy_fn,
const int global_startpoint[/* vdim */],
const int directions[/* hdim * vdim */],
const int lengths[/* hdim */],
@@ -57,11 +54,9 @@ int NewHyperslab_GetHyperslab (cGH *GH,
int hsize[/* hdim */]);
#endif
-PUGHSlab_conversion_fn PUGHSlab_GetDatatypeConversionFn (int vtype, int htype);
-
#ifdef __cplusplus
}
#endif
-#endif /* _PUGHSLAB_PUGHSLAB_H_ */
+#endif /* _PUGHSLAB_NEWPUGHSLAB_H_ */
diff --git a/src/PUGHSlab.h b/src/PUGHSlab.h
index c401602..c2ff990 100644
--- a/src/PUGHSlab.h
+++ b/src/PUGHSlab.h
@@ -1,17 +1,61 @@
/*@@
- @header Hyperslab.h
+ @header PUGHSlab.h
@date Sun 28 May 2000
@author Thomas Radke
@desc
- Function declarations of thorn Hyperslab
+ Global function declarations of thorn PUGHSlab
@enddesc
- @history
- @endhistory
+ @version $Header$
@@*/
+#ifndef _PUGHSLAB_PUGHSLAB_H_
+#define _PUGHSLAB_PUGHSLAB_H_
+
+#ifdef __cplusplus
+extern "C"
+{
+#endif
+
+/* prototype of a datatype conversion routine */
+typedef void (*t_hslabConversionFn) (const void *src,
+ void *dst,
+ CCTK_INT nelems,
+ CCTK_INT src_stride,
+ CCTK_INT dst_stride);
+
/* function prototypes */
-int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
+CCTK_INT Hyperslab_Get (const cGH *GH,
+ CCTK_INT mapping_handle,
+ const CCTK_INT vindex,
+ const CCTK_INT timelevel,
+ const CCTK_INT hdatatype,
+ void *hdata);
+CCTK_INT Hyperslab_GetList (const cGH *GH,
+ CCTK_INT mapping_handle,
+ CCTK_INT num_arrays,
+ const CCTK_INT *vindices /* num_arrays */,
+ const CCTK_INT *timelevels /* num_arrays */,
+ const CCTK_INT *hdatatypes /* num_arrays */,
+ void *const *hdata /* num_arrays */);
+CCTK_INT Hyperslab_DefineGlobalMappingByIndex (
+ const cGH *GH,
+ CCTK_INT vindex,
+ CCTK_INT hdim,
+ const CCTK_INT *direction /* vdim*hdim */,
+ const CCTK_INT *origin /* vdim */,
+ const CCTK_INT *extent /* hdim */,
+ const CCTK_INT *downsample /* hdim */,
+ CCTK_INT table_handle,
+ CCTK_INT target_proc,
+ t_hslabConversionFn conversion_fn,
+ CCTK_INT *hsize /* hdim */);
+CCTK_INT Hyperslab_FreeMapping (CCTK_INT mapping_handle);
+
+
+int Hyperslab_GetLocalHyperslab (const cGH *GH,
+ int vindex,
+ int vtimelvl,
int hdim,
const int global_startpoint [/*vdim*/],
const int directions [/*vdim*/],
@@ -20,10 +64,19 @@ int Hyperslab_GetLocalHyperslab (const cGH *GH, int vindex, int vtimelvl,
void **hdata,
int hsize [/*hdim*/], int ghsize [/*hdim*/],
int hoffset [/*hdim*/]);
-int Hyperslab_GetHyperslab (const cGH *GH, int target_proc, int vindex, int vtimelvl,
+int Hyperslab_GetHyperslab (const cGH *GH,
+ int target_proc,
+ int vindex,
+ int vtimelvl,
int hdim,
const int global_startpoint [/*vdim*/],
const int directions [/*vdim*/],
const int lengths [/*hdim*/],
const int downsample_ [/*hdim*/],
void **hdata, int hsize [/*hdim*/]);
+
+#ifdef __cplusplus
+}
+#endif
+
+#endif /* _PUGHSLAB_PUGHSLAB_H_ */