aboutsummaryrefslogtreecommitdiff
path: root/Carpet
diff options
context:
space:
mode:
authorErik Schnetter <schnetter@gmail.com>2012-02-25 11:01:40 -0500
committerBarry Wardell <barry.wardell@gmail.com>2012-09-11 18:23:06 +0100
commit487280c1b638d05f4bbc2bc41c81d1fdb424cd4d (patch)
treed4ea7557f7cd42f74ce6f7a92b906e6401c69c97 /Carpet
parentae44c5e6d86a65cce382e313e2a10b214e9645da (diff)
CarpetLib: Change API to obtain pointer to grid function data
Change the API to obtain a pointer to grid function data: - Use a function "typed_data_pointer" instead of overloading the () operator (because this looks nicer) - Don't use a virtual function (because this isn't needed) - Update all uses
Diffstat (limited to 'Carpet')
-rw-r--r--Carpet/Carpet/src/Cycle.cc36
-rw-r--r--Carpet/Carpet/src/Storage.cc48
-rw-r--r--Carpet/Carpet/src/helpers.cc2
-rw-r--r--Carpet/Carpet/src/modes.cc7
-rw-r--r--Carpet/CarpetIOASCII/src/ioascii.cc12
-rw-r--r--Carpet/CarpetIOHDF5/src/Output.cc8
-rw-r--r--Carpet/CarpetIOHDF5/src/OutputSlice.cc2
-rw-r--r--Carpet/CarpetInterp/src/interp.cc7
-rw-r--r--Carpet/CarpetInterp2/src/fasterp.cc10
-rw-r--r--Carpet/CarpetLib/src/gf.cc4
-rw-r--r--Carpet/CarpetLib/src/gf.hh31
-rw-r--r--Carpet/CarpetLib/src/ggf.hh31
-rw-r--r--Carpet/CarpetReduce/src/reduce.cc10
-rw-r--r--Carpet/CarpetRegrid/src/automatic.cc2
-rw-r--r--Carpet/CarpetSlab/src/GetHyperslab.cc4
-rw-r--r--Carpet/CarpetSlab/src/slab.cc4
16 files changed, 136 insertions, 82 deletions
diff --git a/Carpet/Carpet/src/Cycle.cc b/Carpet/Carpet/src/Cycle.cc
index 64b7ff25b..449a5a7b3 100644
--- a/Carpet/Carpet/src/Cycle.cc
+++ b/Carpet/Carpet/src/Cycle.cc
@@ -82,11 +82,13 @@ namespace Carpet {
{
int const varindex = firstvarindex + var;
for (int tl=0; tl<numtimelevels; ++tl) {
- cctkGH->data[varindex][tl]
- = (tl < groupdata.AT(group).info.activetimelevels
- ? ((*arrdata.AT(group).AT(0).data.AT(var))
- (tl, 0, 0, 0)->storage())
- : NULL);
+ if (tl < groupdata.AT(group).info.activetimelevels) {
+ ggf *const ff = arrdata.AT(group).AT(0).data.AT(var);
+ void *const ptr = ff->data_pointer(tl, 0, 0, 0)->storage();
+ cctkGH->data[varindex][tl] = ptr;
+ } else {
+ cctkGH->data[varindex][tl] = NULL;
+ }
}
}
}
@@ -148,11 +150,13 @@ namespace Carpet {
{
int const varindex = firstvarindex + var;
for (int tl=0; tl<numtimelevels; ++tl) {
- cctkGH->data[varindex][tl]
- = (tl < groupdata.AT(group).info.activetimelevels
- ? ((*arrdata.AT(group).AT(0).data.AT(var))
- (tl, 0, 0, 0)->storage())
- : NULL);
+ if (tl < groupdata.AT(group).info.activetimelevels) {
+ ggf *const ff = arrdata.AT(group).AT(0).data.AT(var);
+ void *const ptr = ff->data_pointer(tl, 0, 0, 0)->storage();
+ cctkGH->data[varindex][tl] = ptr;
+ } else {
+ cctkGH->data[varindex][tl] = NULL;
+ }
}
}
}
@@ -205,11 +209,13 @@ namespace Carpet {
{
int const varindex = firstvarindex + var;
for (int tl=0; tl<numtimelevels; ++tl) {
- cctkGH->data[varindex][tl]
- = (tl < groupdata.AT(group).info.activetimelevels
- ? ((*arrdata.AT(group).AT(0).data.AT(var))
- (tl, 0, 0, 0)->storage())
- : NULL);
+ if (tl < groupdata.AT(group).info.activetimelevels) {
+ ggf *const ff = arrdata.AT(group).AT(0).data.AT(var);
+ void *const ptr = ff->data_pointer(tl, 0, 0, 0)->storage();
+ cctkGH->data[varindex][tl] = ptr;
+ } else {
+ cctkGH->data[varindex][tl] = NULL;
+ }
}
}
}
diff --git a/Carpet/Carpet/src/Storage.cc b/Carpet/Carpet/src/Storage.cc
index 5503936c8..2ebb4fc81 100644
--- a/Carpet/Carpet/src/Storage.cc
+++ b/Carpet/Carpet/src/Storage.cc
@@ -21,7 +21,7 @@ namespace Carpet {
static int
- GroupStorageCrease (const cGH* cgh, int n_groups, const int* groups,
+ GroupStorageCrease (const cGH* cctkGH, int n_groups, const int* groups,
const int* tls, int* status,
const bool inc);
@@ -33,13 +33,13 @@ namespace Carpet {
int
- GroupStorageCrease (const cGH* cgh, int n_groups, const int* groups,
+ GroupStorageCrease (const cGH* cctkGH, int n_groups, const int* groups,
const int* tls, int* status,
const bool inc)
{
DECLARE_CCTK_PARAMETERS;
- assert (cgh);
+ assert (cctkGH);
assert (n_groups >= 0);
assert (groups);
assert (tls);
@@ -184,11 +184,13 @@ namespace Carpet {
if (gp.grouptype != CCTK_GF) {
assert (rl==0 and m==0);
for (int tl=0; tl<gp.numtimelevels; ++tl) {
- cgh->data[varindex][tl]
- = (tl < groupdata.AT(group).activetimelevels.AT(ml).AT(rl)
- ? ((*arrdata.AT(group).AT(m).data.AT(var))
- (tl, 0, 0, 0)->storage())
- : NULL);
+ if (tl < groupdata.AT(group).activetimelevels.AT(ml).AT(rl)) {
+ ggf *const ff = arrdata.AT(group).AT(0).data.AT(var);
+ void *const ptr = ff->data_pointer(tl, 0, 0, 0)->storage();
+ cctkGH->data[varindex][tl] = ptr;
+ } else {
+ cctkGH->data[varindex][tl] = NULL;
+ }
}
} // if grouptype != GF
@@ -198,7 +200,7 @@ namespace Carpet {
} // if really change the number of active time levels
// Complain if there are not enough active time levels
- GroupStorageCheck (cgh, group, ml, rl);
+ GroupStorageCheck (cctkGH, group, ml, rl);
// Record current number of time levels
// Note: This adds the time levels of all refinement levels
@@ -218,7 +220,7 @@ namespace Carpet {
int
- GroupStorageIncrease (const cGH* cgh, int n_groups, const int* groups,
+ GroupStorageIncrease (const cGH* cctkGH, int n_groups, const int* groups,
const int* tls, int* status)
{
DECLARE_CCTK_PARAMETERS
@@ -227,13 +229,13 @@ namespace Carpet {
Checkpoint ("GroupStorageIncrease");
}
return
- GroupStorageCrease (cgh, n_groups, groups, tls, status, true);
+ GroupStorageCrease (cctkGH, n_groups, groups, tls, status, true);
}
int
- GroupStorageDecrease (const cGH* cgh, int n_groups, const int* groups,
+ GroupStorageDecrease (const cGH* cctkGH, int n_groups, const int* groups,
const int* tls, int* status)
{
DECLARE_CCTK_PARAMETERS
@@ -242,19 +244,19 @@ namespace Carpet {
Checkpoint ("GroupStorageDecrease");
}
return
- GroupStorageCrease (cgh, n_groups, groups, tls, status, false);
+ GroupStorageCrease (cctkGH, n_groups, groups, tls, status, false);
}
int
- EnableGroupStorage (const cGH* cgh, const char* groupname)
+ EnableGroupStorage (const cGH* cctkGH, const char* groupname)
{
const int group = CCTK_GroupIndex(groupname);
assert (group>=0 and group<CCTK_NumGroups());
const int tls = CCTK_MaxTimeLevelsGI(group);
int status;
- GroupStorageIncrease (cgh, 1, &group, &tls, &status);
+ GroupStorageIncrease (cctkGH, 1, &group, &tls, &status);
// Return whether storage was allocated previously
return status;
}
@@ -262,13 +264,13 @@ namespace Carpet {
int
- DisableGroupStorage (const cGH* cgh, const char* groupname)
+ DisableGroupStorage (const cGH* cctkGH, const char* groupname)
{
const int group = CCTK_GroupIndex(groupname);
assert (group>=0 and group<CCTK_NumGroups());
const int tls = 0;
int status;
- GroupStorageDecrease (cgh, 1, &group, &tls, &status);
+ GroupStorageDecrease (cctkGH, 1, &group, &tls, &status);
// Return whether storage was allocated previously
return status;
}
@@ -276,7 +278,7 @@ namespace Carpet {
int
- QueryGroupStorageB (const cGH* cgh, int group, const char* groupname)
+ QueryGroupStorageB (const cGH* cctkGH, int group, const char* groupname)
{
DECLARE_CCTK_PARAMETERS;
if (groupname) {
@@ -301,7 +303,7 @@ namespace Carpet {
const int*
- ArrayGroupSizeB (const cGH* cgh, int dir, int group, const char* groupname)
+ ArrayGroupSizeB (const cGH* cctkGH, int dir, int group, const char* groupname)
{
static const int zero = 0;
static const int error = 0;
@@ -323,7 +325,7 @@ namespace Carpet {
const int gpdim = groupdata.AT(group).info.dim;
assert (dir>=0 and dir<gpdim);
- if (CCTK_QueryGroupStorageI(cgh, group)) {
+ if (CCTK_QueryGroupStorageI(cctkGH, group)) {
return &groupdata.AT(group).info.lsh[dir];
@@ -337,7 +339,7 @@ namespace Carpet {
- int GroupDynamicData (const cGH* cgh, int group, cGroupDynamicData* data)
+ int GroupDynamicData (const cGH* cctkGH, int group, cGroupDynamicData* data)
{
// Return values:
// 0 for success
@@ -345,10 +347,10 @@ namespace Carpet {
// -3 if given GH pointer is invalid
// (-77 if group has zero variables)
// -78 if group does not exist
- if (not cgh) return -3;
+ if (not cctkGH) return -3;
if (not (group>=0 and group<CCTK_NumGroups())) return -78;
if (not data) return -1;
- assert (cgh);
+ assert (cctkGH);
assert (group>=0 and group<CCTK_NumGroups());
assert (data);
*data = groupdata.AT(group).info;
diff --git a/Carpet/Carpet/src/helpers.cc b/Carpet/Carpet/src/helpers.cc
index ccf5f71cc..760558368 100644
--- a/Carpet/Carpet/src/helpers.cc
+++ b/Carpet/Carpet/src/helpers.cc
@@ -124,7 +124,7 @@ namespace Carpet {
int const lc = arrdata.AT(groupindex).AT(m).hh->get_local_component(rl, c);
assert (lc >= 0 and
lc < arrdata.AT(groupindex).AT(m).hh->local_components(rl));
- gdata * const data = (*ff) (tl, rl, lc, mglevel);
+ gdata * const data = ff->data_pointer (tl, rl, lc, mglevel);
return data->storage();
} else {
return NULL;
diff --git a/Carpet/Carpet/src/modes.cc b/Carpet/Carpet/src/modes.cc
index ec26ada6a..b5c9a9fb6 100644
--- a/Carpet/Carpet/src/modes.cc
+++ b/Carpet/Carpet/src/modes.cc
@@ -197,7 +197,7 @@ namespace Carpet {
for (int tl=0; tl<max_tl; ++tl) {
if (ff and tl<active_tl) {
int const lc = 0;
- gdata * const data = (*ff) (tl, rl, lc, ml);
+ gdata * const data = ff->data_pointer (tl, rl, lc, ml);
assert (data);
cctkGH->data[firstvar+var][tl] = data->storage();
} else {
@@ -622,7 +622,7 @@ namespace Carpet {
int available_tl;
int tl_offset;
if (active_tl == 0) {
- // gropu has no storage
+ // group has no storage
available_tl = active_tl;
tl_offset = 0;
} else if (do_allow_past_timelevels) {
@@ -650,7 +650,8 @@ namespace Carpet {
for (int tl=0; tl<max_tl; ++tl) {
if (ff and tl<available_tl) {
gdata * const data =
- (*ff) (tl_offset+tl, reflevel, local_component, mglevel);
+ ff->data_pointer
+ (tl_offset+tl, reflevel, local_component, mglevel);
assert (data);
cctkGH->data[firstvar+var][tl] = data->storage();
} else {
diff --git a/Carpet/CarpetIOASCII/src/ioascii.cc b/Carpet/CarpetIOASCII/src/ioascii.cc
index fcd217d9a..cb29249df 100644
--- a/Carpet/CarpetIOASCII/src/ioascii.cc
+++ b/Carpet/CarpetIOASCII/src/ioascii.cc
@@ -605,7 +605,7 @@ namespace CarpetIOASCII {
if (dist::rank() == proc) {
const ggf* const ff =
arrdata.at(group).at(m).data.at(n + n_min);
- datas.at(n) = (*ff) (tl, rl, lc, ml);
+ datas.at(n) = ff->data_pointer (tl, rl, lc, ml);
} else {
datas.at(n) = NULL;
}
@@ -1375,7 +1375,7 @@ namespace CarpetIOASCII {
const vect<int,outdim> str = gfext.stride()[dirs];
const bbox<int,outdim> ext(lo,up,str);
- // Check whether the output origin is contained in the extent of
+ // check whether the output origin is contained in the extent of
// the data that should be output
ivect org1(org);
for (int d=0; d<outdim; ++d) org1[dirs[d]] = ext.lower()[d];
@@ -1428,10 +1428,10 @@ namespace CarpetIOASCII {
const gdata* gfdata = gfdatas.at(n);
os << (n==0 ? "\t" : " ");
switch (specific_cactus_type(vartype)) {
-#define TYPECASE(N,T) \
- case N: \
- os << (*(const data<T>*)gfdata)[index]; \
- break;
+#define TYPECASE(N,T) \
+ case N: \
+ os << (*(const data<T>*)gfdata)[index]; \
+ break;
#include "typecase.hh"
#undef TYPECASE
default:
diff --git a/Carpet/CarpetIOHDF5/src/Output.cc b/Carpet/CarpetIOHDF5/src/Output.cc
index d77324860..f8c580335 100644
--- a/Carpet/CarpetIOHDF5/src/Output.cc
+++ b/Carpet/CarpetIOHDF5/src/Output.cc
@@ -189,8 +189,8 @@ int WriteVarUnchunked (const cGH* const cctkGH,
const ggf* ff = arrdata.at(gindex).at(Carpet::map).data.at(var);
const gdata* const data =
local_component >= 0
- ? (*ff) (request->timelevel, refinementlevel,
- local_component, mglevel)
+ ? ff->data_pointer (request->timelevel, refinementlevel,
+ local_component, mglevel)
: NULL;
#if 0
// TODO: This does not work; data may be NULL
@@ -367,8 +367,8 @@ int WriteVarChunkedSequential (const cGH* const cctkGH,
const ggf* ff = arrdata.at(gindex).at(Carpet::map).data.at(var);
const gdata* const data =
local_component >= 0
- ? (*ff) (request->timelevel, refinementlevel,
- local_component, mglevel)
+ ? ff->data_pointer (request->timelevel, refinementlevel,
+ local_component, mglevel)
: NULL;
#if 0
// TODO: This does not work; data may be NULL
diff --git a/Carpet/CarpetIOHDF5/src/OutputSlice.cc b/Carpet/CarpetIOHDF5/src/OutputSlice.cc
index f0d570ced..9a691be5f 100644
--- a/Carpet/CarpetIOHDF5/src/OutputSlice.cc
+++ b/Carpet/CarpetIOHDF5/src/OutputSlice.cc
@@ -624,7 +624,7 @@ namespace CarpetIOHDF5 {
if (dist::rank() == proc) {
const ggf* const ff =
arrdata.at(group).at(m).data.at(n + n_min);
- datas.at(n) = (*ff) (tl, rl, lc, ml);
+ datas.at(n) = ff->data_pointer (tl, rl, lc, ml);
} else {
datas.at(n) = NULL;
}
diff --git a/Carpet/CarpetInterp/src/interp.cc b/Carpet/CarpetInterp/src/interp.cc
index 2a7a62fdb..864dd85d7 100644
--- a/Carpet/CarpetInterp/src/interp.cc
+++ b/Carpet/CarpetInterp/src/interp.cc
@@ -1706,9 +1706,10 @@ namespace CarpetInterp {
int const gi = CCTK_GroupIndexFromVarI (vi);
int const vi0 = CCTK_FirstVarIndexI (gi);
- input_arrays[n] =
- (*arrdata.AT(gi).AT(m).data.AT(vi-vi0))(my_tl, rl, lc, mglevel)->
- storage();
+ ggf const *const ff = arrdata.AT(gi).AT(m).data.AT(vi-vi0);
+ void const *const ptr =
+ ff->data_pointer(my_tl, rl, lc, mglevel)->storage();
+ input_arrays[n] = ptr;
}
} // for input arrays
diff --git a/Carpet/CarpetInterp2/src/fasterp.cc b/Carpet/CarpetInterp2/src/fasterp.cc
index 9b8fd956b..ba46eec69 100644
--- a/Carpet/CarpetInterp2/src/fasterp.cc
+++ b/Carpet/CarpetInterp2/src/fasterp.cc
@@ -508,7 +508,7 @@ namespace CarpetInterp2 {
fasterp_llocs_t local_locations (npoints);
if (CCTK_IsFunctionAliased ("MultiPatch_GlobalToLocal")) {
- // This is a muulti-patch simulation: convert global to local
+ // This is a multi-patch simulation: convert global to local
// coordinates
CCTK_REAL const * coords[dim];
@@ -1176,10 +1176,10 @@ namespace CarpetInterp2 {
assert (gi >= 0);
int const vi = varinds.AT(v) - CCTK_FirstVarIndexI (gi);
assert (vi >= 0);
- varptrs.AT(v) =
- (CCTK_REAL const *)
- (* Carpet::arrdata.AT(gi).AT(m).data.AT(vi))
- (tl, rl, lc, Carpet::mglevel)->storage();
+ ggf const *const ff = Carpet::arrdata.AT(gi).AT(m).data.AT(vi);
+ void const *const ptr =
+ ff->data_pointer(tl, rl, lc, Carpet::mglevel)->storage();
+ varptrs.AT(v) = (CCTK_REAL const *)ptr;
assert (varptrs.AT(v));
}
diff --git a/Carpet/CarpetLib/src/gf.cc b/Carpet/CarpetLib/src/gf.cc
index c9b003efc..8a27f9eef 100644
--- a/Carpet/CarpetLib/src/gf.cc
+++ b/Carpet/CarpetLib/src/gf.cc
@@ -16,7 +16,7 @@ gf<T>::gf (const int varindex_, const operator_type transport_operator_,
th& t_, dh& d_,
const int prolongation_order_time_,
const int vectorlength_, const int vectorindex_,
- gf* const vectorleader_)
+ ggf* const vectorleader_)
: ggf(varindex_, transport_operator_,
t_, d_, prolongation_order_time_,
vectorlength_, vectorindex_, vectorleader_)
@@ -39,6 +39,7 @@ gf<T>::~gf ()
+#if 0
// Access to the data
template<typename T>
const data<T>* gf<T>::operator() (int tl, int rl, int lc, int ml) const
@@ -65,6 +66,7 @@ data<T>* gf<T>::operator() (int tl, int rl, int lc, int ml)
assert (tl>=0 and tl<timelevels(ml, rl));
return (data<T>*)storage.AT(ml).AT(rl).AT(lc).AT(tl);
}
+#endif
diff --git a/Carpet/CarpetLib/src/gf.hh b/Carpet/CarpetLib/src/gf.hh
index a610ea46f..bf58ac7ff 100644
--- a/Carpet/CarpetLib/src/gf.hh
+++ b/Carpet/CarpetLib/src/gf.hh
@@ -37,7 +37,7 @@ public:
th& t, dh& d,
const int prolongation_order_time,
const int vectorlength, const int vectorindex,
- gf* const vectorleader);
+ ggf* const vectorleader);
// Destructors
virtual ~gf ();
@@ -50,7 +50,7 @@ public:
{
data<T>* const vl =
this->vectorleader
- ? (data<T>*)(*this->vectorleader)(tl,rl,lc,ml)
+ ? ((gf<T>*)this->vectorleader)->typed_data_pointer(tl,rl,lc,ml)
: NULL;
return new data<T>(this->varindex,
h.refcent, this->transport_operator,
@@ -58,8 +58,7 @@ public:
vl);
}
- virtual gdata*
- new_typed_data () const
+ virtual gdata* new_typed_data () const
{
return new data<T>(this->varindex,
h.refcent, this->transport_operator,
@@ -70,9 +69,29 @@ public:
// Access to the data
- virtual const data<T>* operator() (int tl, int rl, int lc, int ml) const;
+#if 0
+ virtual const data<T>* operator() (int tl, int rl, int lc, int ml) const
+ CCTK_MEMBER_ATTRIBUTE_PURE;
+ virtual data<T>* operator() (int tl, int rl, int lc, int ml)
+ CCTK_MEMBER_ATTRIBUTE_PURE;
+#endif
- virtual data<T>* operator() (int tl, int rl, int lc, int ml);
+ data<T> const* typed_data_pointer (int tl, int rl, int lc, int ml) const
+ {
+ assert (rl>=0 and rl<h.reflevels());
+ assert (lc>=0 and lc<h.local_components(rl));
+ assert (ml>=0 and ml<h.mglevels());
+ assert (tl>=0 and tl<timelevels(ml, rl));
+ return (data<T> const*)storage.AT(ml).AT(rl).AT(lc).AT(tl);
+ }
+ data<T>* typed_data_pointer (int tl, int rl, int lc, int ml)
+ {
+ assert (rl>=0 and rl<h.reflevels());
+ assert (lc>=0 and lc<h.local_components(rl));
+ assert (ml>=0 and ml<h.mglevels());
+ assert (tl>=0 and tl<timelevels(ml, rl));
+ return (data<T>*)storage.AT(ml).AT(rl).AT(lc).AT(tl);
+ }
diff --git a/Carpet/CarpetLib/src/ggf.hh b/Carpet/CarpetLib/src/ggf.hh
index 56878611c..e6b32713f 100644
--- a/Carpet/CarpetLib/src/ggf.hh
+++ b/Carpet/CarpetLib/src/ggf.hh
@@ -119,7 +119,7 @@ public:
// mg_prolongate, ref_restrict, or ref_prolongate.
// "Updating" means here that the boundaries have to be
- // synchronised. They don't need to be prolongated.
+ // synchronised. They don't need to be prolongated.
// Synchronise the boundaries of a component
void sync_all (comm_state& state, int tl, int rl, int ml);
@@ -153,9 +153,7 @@ public:
// Helpers
virtual gdata* typed_data (int tl, int rl, int lc, int ml) const = 0;
-
- virtual gdata*
- new_typed_data () const = 0;
+ virtual gdata* new_typed_data () const = 0;
@@ -199,8 +197,29 @@ protected:
public:
// Access to the data
- virtual const gdata* operator() (int tl, int rl, int lc, int ml) const CCTK_ATTRIBUTE_PURE = 0;
- virtual gdata* operator() (int tl, int rl, int lc, int ml) CCTK_ATTRIBUTE_PURE = 0;
+ const gdata* data_pointer (int const tl, int const rl, int const lc, int const ml) const
+ {
+ assert (rl>=0 and rl<h.reflevels());
+ assert (lc>=0 and lc<h.local_components(rl));
+ assert (ml>=0 and ml<h.mglevels());
+ assert (tl>=0 and tl<timelevels(ml, rl));
+ return storage.AT(ml).AT(rl).AT(lc).AT(tl);
+ }
+ gdata* data_pointer (int const tl, int const rl, int const lc, int const ml)
+ {
+ assert (rl>=0 and rl<h.reflevels());
+ assert (lc>=0 and lc<h.local_components(rl));
+ assert (ml>=0 and ml<h.mglevels());
+ assert (tl>=0 and tl<timelevels(ml, rl));
+ return storage.AT(ml).AT(rl).AT(lc).AT(tl);
+ }
+
+#if 0
+ virtual const gdata* operator() (int tl, int rl, int lc, int ml) const
+ CCTK_MEMBER_ATTRIBUTE_PURE = 0;
+ virtual gdata* operator() (int tl, int rl, int lc, int ml)
+ CCTK_MEMBER_ATTRIBUTE_PURE = 0;
+#endif
diff --git a/Carpet/CarpetReduce/src/reduce.cc b/Carpet/CarpetReduce/src/reduce.cc
index d18de8ce0..a849f100b 100644
--- a/Carpet/CarpetReduce/src/reduce.cc
+++ b/Carpet/CarpetReduce/src/reduce.cc
@@ -1464,9 +1464,13 @@ namespace CarpetReduce {
int const vi = invars[n];
int const gi = CCTK_GroupIndexFromVarI (vi);
int const vi0 = CCTK_FirstVarIndexI (gi);
- myinarrays.AT(tl).AT(n)
- = ((*arrdata.AT(gi).AT(Carpet::map).data.AT(vi-vi0))
- (tl, reflevel, local_component, mglevel)->storage());
+ ggf const *const ff =
+ arrdata.AT(gi).AT(Carpet::map).data.AT(vi-vi0);
+ void const *const ptr =
+ ff->
+ data_pointer(tl, reflevel, local_component, mglevel)->
+ storage();
+ myinarrays.AT(tl).AT(n) = ptr;
#endif
assert (myinarrays.AT(tl).AT(n));
}
diff --git a/Carpet/CarpetRegrid/src/automatic.cc b/Carpet/CarpetRegrid/src/automatic.cc
index 7248b732a..a524edae4 100644
--- a/Carpet/CarpetRegrid/src/automatic.cc
+++ b/Carpet/CarpetRegrid/src/automatic.cc
@@ -104,7 +104,7 @@ namespace CarpetRegrid {
assert (! region.extent.empty());
const data<CCTK_REAL>& errordata =
- *dynamic_cast<const data<CCTK_REAL>*> (errorgf(tl,rl,c,ml));
+ *errorgf.typed_data_pointer(tl,rl,c,ml);
Automatic_Recursive (cctkGH, hh, errordata, regl, region, reffact);
}
diff --git a/Carpet/CarpetSlab/src/GetHyperslab.cc b/Carpet/CarpetSlab/src/GetHyperslab.cc
index 27daeec1c..cd8cb6917 100644
--- a/Carpet/CarpetSlab/src/GetHyperslab.cc
+++ b/Carpet/CarpetSlab/src/GetHyperslab.cc
@@ -142,7 +142,7 @@ namespace CarpetSlab {
// Get sample data
const gdata* mydata;
- mydata = (*myff)(tl, rl, 0, 0);
+ mydata = myff->data_pointer(tl, rl, 0, 0);
// Stride of data in memory
const vect<int,dim> str = mydata->extent().stride();
@@ -183,7 +183,7 @@ namespace CarpetSlab {
BEGIN_LOCAL_COMPONENT_LOOP (cgh, gp.grouptype) {
// Get data object
- mydata = (*myff)(tl, rl, component, mglevel);
+ mydata = myff->data_pointer(tl, rl, component, mglevel);
// Calculate overlapping extents
bboxset<int,dim> const myextents =
diff --git a/Carpet/CarpetSlab/src/slab.cc b/Carpet/CarpetSlab/src/slab.cc
index 9a7ac6c4a..1825602c3 100644
--- a/Carpet/CarpetSlab/src/slab.cc
+++ b/Carpet/CarpetSlab/src/slab.cc
@@ -175,7 +175,7 @@ namespace CarpetSlab {
// Get sample data
const gdata* mydata;
- mydata = (*myff)(tl, rl, 0, 0);
+ mydata = myff->data_pointer(tl, rl, 0, 0);
// Stride of data in memory
const vect<int,dim> str = mydata->extent().stride();
@@ -216,7 +216,7 @@ namespace CarpetSlab {
BEGIN_COMPONENT_LOOP (cgh, gp.grouptype) {
// Get data object
- mydata = (*myff)(tl, rl, component, mglevel);
+ mydata = myff->data_pointer(tl, rl, component, mglevel);
// Calculate overlapping extents
const bboxset<int,dim> myextents =