aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2002-07-15 17:16:10 +0000
committerschnetter <schnetter@b716e942-a2de-43ad-8f52-f3dfe468e4e7>2002-07-15 17:16:10 +0000
commit11ec06ae3f07d16d4da85e725941bb207b75ad11 (patch)
treed9204efc9069ad72ea40cdd5cc50d645db766cc6
parent42a98ca518377949e8cfc902343a4b2d97a25564 (diff)
Added file that I forgot, a Fortran interface to the hyperslabber.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinUtils/TGRtensor/trunk@3 b716e942-a2de-43ad-8f52-f3dfe468e4e7
-rw-r--r--src/hyperslab.c89
1 files changed, 89 insertions, 0 deletions
diff --git a/src/hyperslab.c b/src/hyperslab.c
new file mode 100644
index 0000000..8f03fd4
--- /dev/null
+++ b/src/hyperslab.c
@@ -0,0 +1,89 @@
+/* $Header$ */
+
+#include <assert.h>
+#include <stdlib.h>
+#include <string.h>
+
+#include "cctk.h"
+#include "cctk_Parameters.h"
+
+#include "Hyperslab.h"
+
+#if 0
+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*/]);
+#endif
+
+int Hyperslab_FillHyperslab (const cGH *cctkGH,
+ int target_proc,
+ int vindex,
+ int vtimelevel,
+ int hdim,
+ const int global_startpoint [/*vdim*/],
+ const int directions [/*vdim*/],
+ const int lengths [/*hdim*/],
+ const int downsample [/*hdim*/],
+ int nhdata,
+ void *hdata,
+ int hsize [/*hdim*/])
+{
+ void *myhdata;
+ int mynhdata;
+ int vtype;
+ int vtypesize;
+ int d;
+ int ierr;
+ myhdata = 0;
+ ierr = Hyperslab_GetHyperslab
+ (cctkGH, target_proc, vindex, vtimelevel, hdim,
+ global_startpoint, directions, lengths, downsample,
+ &myhdata, hsize);
+ if (myhdata) {
+ mynhdata = 1;
+ for (d=0; d<hdim; ++d) mynhdata *= hsize[d];
+ if (mynhdata > nhdata) {
+ CCTK_WARN (1, "Not enough space in hdata array");
+ mynhdata = nhdata;
+ }
+ assert (hdata);
+ vtype = CCTK_VarTypeI(vindex);
+ assert (vtype>=0);
+ vtypesize = CCTK_VarTypeSize(vtype);
+ assert (vtypesize>0);
+ memcpy (hdata, myhdata, mynhdata * vtypesize);
+ free (myhdata);
+ }
+ return ierr;
+}
+
+
+
+void CCTK_FCALL CCTK_FNAME(Hyperslab_FillHyperslab)
+ (int *ierr,
+ const cGH *cctkGH,
+ const int *target_proc,
+ const int *vindex,
+ const int *vtimelevel,
+ const int *hdim,
+ const int global_startpoint [/*vdim*/],
+ const int directions [/*vdim*/],
+ const int lengths [/*hdim*/],
+ const int downsample [/*hdim*/],
+ const int *nhdata,
+ void *hdata,
+ int hsize [/*hdim*/])
+{
+ *ierr = Hyperslab_FillHyperslab
+ (cctkGH, *target_proc, *vindex, *vtimelevel, *hdim,
+ global_startpoint, directions, lengths, downsample,
+ *nhdata, hdata, hsize);
+}