aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorRoland Haas <rhaas@tapir.caltech.edu>2013-09-04 14:45:56 -0700
committerRoland Haas <rhaas@tapir.caltech.edu>2013-09-27 10:41:38 -0700
commit0e435fca01c0ce3356565c661afb0ff9d64880a1 (patch)
tree78aa48df07d6def528a5b4368b96691e9de57722
parent38923afb0681b02fae2171cf54671cbe587ff748 (diff)
CarpetLib: add some code for electric fence like functionality
-rw-r--r--Carpet/CarpetLib/param.ccl9
-rw-r--r--Carpet/CarpetLib/src/data.cc13
-rw-r--r--Carpet/CarpetLib/src/data.hh1
-rw-r--r--Carpet/CarpetLib/src/gdata.hh2
-rw-r--r--Carpet/CarpetLib/src/mem.cc55
-rw-r--r--Carpet/CarpetLib/src/mem.hh3
6 files changed, 80 insertions, 3 deletions
diff --git a/Carpet/CarpetLib/param.ccl b/Carpet/CarpetLib/param.ccl
index 541ec2f49..5f1d23b82 100644
--- a/Carpet/CarpetLib/param.ccl
+++ b/Carpet/CarpetLib/param.ccl
@@ -76,6 +76,15 @@ BOOLEAN poison_new_memory "Try to catch uninitialised data by setting newly allo
{
} "no"
+BOOLEAN electric_fence "Surround each allocated memory block by canaries to check for out-of-bounds accesses" STEERABLE=recover
+{
+} "no"
+
+CCTK_INT fence_width "number of guard cells to use" STEERABLE=recover
+{
+ 1:* :: "any number of cells"
+} 1
+
RESTRICTED:
CCTK_INT poison_value "Integer value (0..255) used to poison new timelevels (with memset)" STEERABLE=always
diff --git a/Carpet/CarpetLib/src/data.cc b/Carpet/CarpetLib/src/data.cc
index 0d1b29d9e..9ed352df2 100644
--- a/Carpet/CarpetLib/src/data.cc
+++ b/Carpet/CarpetLib/src/data.cc
@@ -561,6 +561,19 @@ size_t data<T>::allocsize (const ibbox & extent_, const int proc_) const
return vectorlength * prod (pad_shape (extent_)) * sizeof (T);
}
+template<typename T>
+bool data<T>::check_fence (const int upperlower) const
+{
+ assert ((vectorleader != NULL) xor (vectorindex == 0));
+ bool retval = true;
+ // since vectors share _memory we only check the first/last vector element
+ if ((vectorindex == 0 and upperlower == 0) or
+ (vectorindex == vectorlength-1 and upperlower == 1)) {
+ retval = _memory->is_fence_intact(upperlower);
+ }
+ return retval;
+}
+
// Data manipulators
diff --git a/Carpet/CarpetLib/src/data.hh b/Carpet/CarpetLib/src/data.hh
index e64899ae0..bcce28e8f 100644
--- a/Carpet/CarpetLib/src/data.hh
+++ b/Carpet/CarpetLib/src/data.hh
@@ -65,6 +65,7 @@ public:
void* const memptr = NULL, size_t const memsize = 0);
virtual void free ();
virtual size_t allocsize (const ibbox& extent, const int proc) const;
+ virtual bool check_fence (const int upperlower) const;
public:
diff --git a/Carpet/CarpetLib/src/gdata.hh b/Carpet/CarpetLib/src/gdata.hh
index 5b5decba6..c000a4d74 100644
--- a/Carpet/CarpetLib/src/gdata.hh
+++ b/Carpet/CarpetLib/src/gdata.hh
@@ -105,6 +105,8 @@ public:
0;
virtual void free () = 0;
virtual size_t allocsize (const ibbox& extent, const int proc) const = 0;
+ // true if fence is intact
+ virtual bool check_fence (const int upperlower) const = 0;
// Accessors
bool has_storage () const {
diff --git a/Carpet/CarpetLib/src/mem.cc b/Carpet/CarpetLib/src/mem.cc
index ddbe0c73b..15435c726 100644
--- a/Carpet/CarpetLib/src/mem.cc
+++ b/Carpet/CarpetLib/src/mem.cc
@@ -95,6 +95,7 @@ mem (size_t const vectorlength, size_t const nelems,
#else
size_t const vector_size = 1;
#endif
+ size_t const canary = electric_fence ? 2*fence_width : 0;
size_t const final_padding = vector_size - 1;
size_t const max_cache_linesize = get_max_cache_linesize();
size_t const alignment =
@@ -103,7 +104,8 @@ mem (size_t const vectorlength, size_t const nelems,
// Safety check
assert(alignment <= 1024);
- const size_t nbytes = (vectorlength * nelems + final_padding) * sizeof (T);
+ const size_t nbytes = (vectorlength * nelems + canary + final_padding) *
+ sizeof (T);
if (max_allowed_memory_MB > 0 and
(total_allocated_bytes + nbytes > MEGA * max_allowed_memory_MB))
{
@@ -133,7 +135,7 @@ mem (size_t const vectorlength, size_t const nelems,
}
storage_base_ = (T*)ptr;
- storage_ = (T*)align_up(size_t(ptr), alignment);
+ storage_ = (T*)align_up(size_t(storage_base_ + canary/2), alignment);
assert(size_t(storage_) % alignment == 0);
owns_storage_ = true;
@@ -142,6 +144,14 @@ mem (size_t const vectorlength, size_t const nelems,
if (poison_new_memory) {
memset (storage_, poison_value, nbytes);
}
+ if (electric_fence) {
+ // put poison just before and just after payload region
+ // FIXME: this will not work with alignment for vectorizing. Not sure how
+ // to support that as well as protect grid scalars.
+ memset (storage_ - fence_width, poison_value, fence_width*sizeof(T));
+ memset (storage_ + vectorlength_ * nelems_, poison_value,
+ fence_width*sizeof(T));
+ }
} else {
assert (memsize >= vectorlength * nelems * sizeof (T));
@@ -158,8 +168,13 @@ template<typename T>
mem<T>::
~mem ()
{
+ DECLARE_CCTK_PARAMETERS;
assert (not has_clients());
if (owns_storage_) {
+ // do we really want this automated check in the destructor?
+ // what if we are already terminating to to a failed fence check?
+ if (electric_fence)
+ assert(is_fence_intact(0) && is_fence_intact(1) );
free(storage_base_);
#if VECTORISE
@@ -167,10 +182,11 @@ mem<T>::
#else
size_t const vector_size = 1;
#endif
+ size_t const canary = electric_fence ? 2*fence_width : 0;
size_t const final_padding = vector_size - 1;
const size_t nbytes =
- (vectorlength_ * nelems_ + final_padding) * sizeof (T);
+ (vectorlength_ * nelems_ + canary + final_padding) * sizeof (T);
total_allocated_bytes -= nbytes;
}
-- total_allocated_objects;
@@ -179,6 +195,38 @@ mem<T>::
template<typename T>
+bool mem<T>::
+is_fence_intact (const int upperlower) const
+{
+ DECLARE_CCTK_PARAMETERS;
+ bool retval = true;
+
+ if (owns_storage_) {
+ assert(storage_ and storage_base_);
+ if (electric_fence) {
+ T worm;
+ memset (&worm, poison_value, sizeof (T));
+ if (upperlower) {
+ for(int i=0; i<fence_width; ++i) {
+ retval = retval &&
+ (memcmp (&worm, storage_ + vectorlength_ * nelems_ + i,
+ sizeof (T)) == 0);
+ }
+ } else {
+ for(int i=0; i<fence_width; ++i) {
+ retval = retval &&
+ (memcmp (&worm, storage_ - 1 - i, sizeof (T)) == 0);
+ }
+ }
+ }
+ }
+
+ return retval;
+}
+
+
+
+template<typename T>
void mem<T>::
register_client (size_t const vectorindex)
{
@@ -248,6 +296,7 @@ mempool::
}
}
+// TODO: add electric fence
void *
mempool::
alloc (size_t nbytes)
diff --git a/Carpet/CarpetLib/src/mem.hh b/Carpet/CarpetLib/src/mem.hh
index 807daff70..5ae98ff75 100644
--- a/Carpet/CarpetLib/src/mem.hh
+++ b/Carpet/CarpetLib/src/mem.hh
@@ -60,6 +60,9 @@ public:
assert (clients_.AT(vectorindex));
return & storage_ [vectorindex * nelems_];
}
+
+ // return true if fence is intact, ie poison value is still there
+ bool is_fence_intact (const int upperlower) const;
void register_client (size_t vectorindex);
void unregister_client (size_t vectorindex);