aboutsummaryrefslogtreecommitdiff
path: root/CarpetAttic
diff options
context:
space:
mode:
authoreschnett <>2001-03-16 23:35:00 +0000
committereschnett <>2001-03-16 23:35:00 +0000
commit62eb5283a3076fef3035c35e4ff6203c6cbea975 (patch)
treeec0328eedd76240378d555b08f884cf8e5e19f36 /CarpetAttic
parent9106019dcedd29654b5d6b62ad1c2f23278e7762 (diff)
Fixed FlexIO output.
darcs-hash:20010316233532-f6438-e89c26eba7e91077b819ce9c64db0dc192a7d6db.gz
Diffstat (limited to 'CarpetAttic')
-rw-r--r--CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc23
1 files changed, 13 insertions, 10 deletions
diff --git a/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc b/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc
index e0494da20..9d2b6d6f7 100644
--- a/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc
+++ b/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc
@@ -28,7 +28,7 @@
#include "ioflexio.hh"
-static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc,v 1.2 2001/03/16 21:32:17 eschnett Exp $";
+static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIO/src/ioflexio.cc,v 1.3 2001/03/17 00:35:32 eschnett Exp $";
@@ -198,7 +198,9 @@ namespace CarpetIOFlexIO {
amrwriter = new AMRwriter(*writer);
// Set datatype
- assert (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL8);
+ assert (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL8
+ || (sizeof(CCTK_REAL) == sizeof(CCTK_REAL8)
+ &&CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL));
// TODO: Set datatype correctly
amrwriter->setType (IObase::Float64);
@@ -210,9 +212,9 @@ namespace CarpetIOFlexIO {
(cgh, &lower[d], &upper[d], d+1, 0, "cart3d");
assert (ierr==0);
origin[d] = lower[d];
- delta[d] = cgh->cctk_delta_space[d] / reflevelfact;
+ delta[d] = cgh->cctk_delta_space[d];
}
- timestep = cgh->cctk_delta_time / reflevelfact;
+ timestep = base_delta_time;
amrwriter->setTopLevelParameters
(dim, origin, delta, timestep, maxreflevels);
@@ -226,10 +228,10 @@ namespace CarpetIOFlexIO {
(interlevel_timerefinement, interlevel_spacerefinement);
// Set level
- amrwriter->setLevel(reflevel);
+ amrwriter->setLevel (reflevel);
// Set current time
- amrwriter->setTime(cgh->cctk_iteration);
+ amrwriter->setTime (cgh->cctk_iteration);
}
// Traverse all components on this refinement and multigrid
@@ -259,7 +261,7 @@ namespace CarpetIOFlexIO {
// Make temporary copy on processor 0
const bbox<int,dim> ext = data->extent();
- generic_data<dim>* tmp = data->make_typed ();
+ generic_data<dim>* const tmp = data->make_typed ();
tmp->allocate (ext, 0);
tmp->copy_from (data, ext);
@@ -267,10 +269,10 @@ namespace CarpetIOFlexIO {
if (CCTK_MyProc(cgh)==0) {
int origin[dim], dims[dim];
for (int d=0; d<dim; ++d) {
- origin[d] = (ext.lower() / ext.stride())[d];
- dims[d] = ((ext.upper() - ext.lower()) / ext.stride() + 1)[d];
+ origin[d] = (ext.lower() / ext.stride())[d];
+ dims[d] = (ext.shape() / ext.stride())[d];
}
- amrwriter->write (origin, dims, (void*)data->storage());
+ amrwriter->write (origin, dims, (void*)tmp->storage());
}
// Delete temporary copy
@@ -286,6 +288,7 @@ namespace CarpetIOFlexIO {
writer = 0;
}
+ break;
} // ARRAY or GROUP
default: