]> https://gitweb.dealii.org/ - dealii.git/commitdiff
pointer to triangulation is stored in CellDataStorage 8484/head
authorReza Rastak <rastak@stanford.edu>
Thu, 8 Aug 2019 22:18:54 +0000 (16:18 -0600)
committerReza Rastak <rastak@stanford.edu>
Fri, 9 Aug 2019 05:02:35 +0000 (23:02 -0600)
doc/news/changes/incompatibilities/20190808RezaRastak [new file with mode: 0644]
include/deal.II/base/quadrature_point_data.h

diff --git a/doc/news/changes/incompatibilities/20190808RezaRastak b/doc/news/changes/incompatibilities/20190808RezaRastak
new file mode 100644 (file)
index 0000000..fed6562
--- /dev/null
@@ -0,0 +1,5 @@
+Improved: CellDataStorage stores a reference to a single unique triangulation
+to ensure that cells from other triangulations cannot store any data on the
+current object.
+<br>
+(Reza Rastak, 2019/08/08)
index 758a0a698fb7491b7055516f0bd3ca74ce6bc8c3..f83a9b4c1adf57a3a3c74f49a9d2c65f6bad0f41 100644 (file)
@@ -91,6 +91,10 @@ public:
    * may reflect, for example, different constitutive models of continuum
    * mechanics in different parts of the domain.
    *
+   * @note The first time this method is called, it stores a SmartPointer to the
+   * Triangulation object that owns the cell. The future invocations of this
+   * method expects the cell to be from the same stored triangulation.
+   *
    * @pre The type @p T needs to either equal @p DataType, or be a class derived
    * from @p DataType. @p T needs to be default constructible.
    */
@@ -138,6 +142,9 @@ public:
    * cell-by-cell basis.
    *
    * @pre The type @p T needs to match the class provided to initialize() .
+   *
+   * @pre @p cell must be from the same Triangulation that is used to
+   * initialize() the cell data.
    */
   template <typename T = DataType>
   std::vector<std::shared_ptr<T>>
@@ -153,6 +160,9 @@ public:
    * cell-by-cell basis.
    *
    * @pre The type @p T needs to match the class provided to initialize() .
+   *
+   * @pre @p cell must be from the same Triangulation that is used to
+   * initialize() the cell data.
    */
   template <typename T = DataType>
   std::vector<std::shared_ptr<const T>>
@@ -160,31 +170,32 @@ public:
 
 private:
   /**
-   * Depending on @p CellIteratorType we need to define a triangulation type
+   * Number of dimensions
    */
-  using triangulation_type =
-    Triangulation<CellIteratorType::AccessorType::dimension,
-                  CellIteratorType::AccessorType::space_dimension>;
+  static constexpr unsigned int dimension =
+    CellIteratorType::AccessorType::dimension;
 
   /**
-   * The key is used as identified for a cell data structure.
-   * We need to use CellId as part of the key because it remains unique during
-   * adaptive refinement. We need to use a pointer to the triangulation in the
-   * key because we want to allow the data storage class to keep information on
-   * multiple meshes simultaneously.
+   * Number of space dimensions
    */
-  using key_type = std::pair<const triangulation_type *, CellId>;
+  static constexpr unsigned int space_dimension =
+    CellIteratorType::AccessorType::space_dimension;
 
   /**
-   * A map to store a vector of data on each cell.
+   * To ensure that all the cells in the CellDataStorage come from the same
+   * Triangulation, we need to store a reference to that Triangulation within
+   * the class.
    */
-  std::map<key_type, std::vector<std::shared_ptr<DataType>>> map;
+  SmartPointer<const Triangulation<dimension, space_dimension>,
+               CellDataStorage<CellIteratorType, DataType>>
+    tria;
 
   /**
-   * Computes the key from a cell iterator
+   * A map to store a vector of data on each cell.
+   * We need to use CellId as the key because it remains unique during
+   * adaptive refinement.
    */
-  static key_type
-  get_key(const CellIteratorType &cell);
+  std::map<CellId, std::vector<std::shared_ptr<DataType>>> map;
 
   /**
    * @addtogroup Exceptions
@@ -192,6 +203,13 @@ private:
   DeclExceptionMsg(
     ExcCellDataTypeMismatch,
     "Cell data is being retrieved with a type which is different than the type used to initialize it");
+
+  /**
+   * @addtogroup Exceptions
+   */
+  DeclExceptionMsg(
+    ExcTriangulationMismatch,
+    "The provided cell iterator does not belong to the triangulation that corresponds to the CellDataStorage object.");
 };
 
 
@@ -544,14 +562,6 @@ namespace parallel
 //                         CellDataStorage
 //--------------------------------------------------------------------
 
-template <typename CellIteratorType, typename DataType>
-typename CellDataStorage<CellIteratorType, DataType>::key_type
-CellDataStorage<CellIteratorType, DataType>::get_key(
-  const CellIteratorType &cell)
-{
-  return {&cell->get_triangulation(), cell->id()};
-}
-
 template <typename CellIteratorType, typename DataType>
 template <typename T>
 void
@@ -561,8 +571,13 @@ CellDataStorage<CellIteratorType, DataType>::initialize(
 {
   static_assert(std::is_base_of<DataType, T>::value,
                 "User's T class should be derived from user's DataType class");
+  // The first time this method is called, it has to initialize the reference
+  // to the triangulation object
+  if (!tria)
+    tria = &cell->get_triangulation();
+  Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
 
-  const auto key = get_key(cell);
+  const auto key = cell->id();
   if (map.find(key) == map.end())
     {
       map[key] = std::vector<std::shared_ptr<DataType>>(n_q_points);
@@ -596,10 +611,11 @@ template <typename CellIteratorType, typename DataType>
 bool
 CellDataStorage<CellIteratorType, DataType>::erase(const CellIteratorType &cell)
 {
-  const auto key = get_key(cell);
+  const auto key = cell->id();
   const auto it  = map.find(key);
   if (it == map.end())
     return false;
+  Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
   for (unsigned int i = 0; i < it->second.size(); i++)
     {
       Assert(
@@ -645,8 +661,9 @@ CellDataStorage<CellIteratorType, DataType>::get_data(
 {
   static_assert(std::is_base_of<DataType, T>::value,
                 "User's T class should be derived from user's DataType class");
+  Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
 
-  const auto it = map.find(get_key(cell));
+  const auto it = map.find(cell->id());
   Assert(it != map.end(), ExcMessage("Could not find data for the cell"));
 
   // It would be nice to have a specialized version of this function for
@@ -673,8 +690,9 @@ CellDataStorage<CellIteratorType, DataType>::get_data(
 {
   static_assert(std::is_base_of<DataType, T>::value,
                 "User's T class should be derived from user's DataType class");
+  Assert(&cell->get_triangulation() == tria, ExcTriangulationMismatch());
 
-  const auto it = map.find(get_key(cell));
+  const auto it = map.find(cell->id());
   Assert(it != map.end(), ExcMessage("Could not find QP data for the cell"));
 
   // Cast base class to the desired class. This has to be done irrespectively of

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.