* 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.
*/
* 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>>
* 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>>
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
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.");
};
// 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
{
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);
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(
{
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
{
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