--- /dev/null
+New: Add method to get the coarse-grid cell from CellID.
+<br>
+(Peter Munch, 2019/08/15)
# define DEAL_II_DOF_INDEX_MPI_TYPE MPI_UNSIGNED
#endif
+#ifdef DEAL_II_WITH_64BIT_INDICES
+ /**
+ * The type used for coarse-cell ids.
+ */
+ using coarse_cell_id = unsigned long long int;
+#else
+ /**
+ * The type used for coarse-cell ids.
+ */
+ using coarse_cell_id = unsigned int;
+#endif
+
/**
* The type used to denote boundary indicators associated with every piece
* of the boundary and, in the case of meshes that describe manifolds in
const types::global_dof_index invalid_dof_index =
static_cast<types::global_dof_index>(-1);
+ /**
+ * An invalid value for coarse-cell ids.
+ */
+ const types::coarse_cell_id invalid_coarse_cell_id =
+ static_cast<types::coarse_cell_id>(-1);
+
/**
* Invalid material_id which we need in several places as a default value.
* We assume that all material_ids lie in the range [0,
* and the number of children of a cell in the current space dimension (i.e.,
* GeometryInfo<dim>::max_children_per_cell).
*/
- CellId(const unsigned int coarse_cell_id,
+ CellId(const types::coarse_cell_id coarse_cell_id,
const std::vector<std::uint8_t> &child_indices);
/**
* GeometryInfo<dim>::max_children_per_cell). The array
* @p child_indices must have at least @p n_child_indices valid entries.
*/
- CellId(const unsigned int coarse_cell_id,
- const unsigned int n_child_indices,
- const std::uint8_t *child_indices);
+ CellId(const types::coarse_cell_id coarse_cell_id,
+ const unsigned int n_child_indices,
+ const std::uint8_t * child_indices);
/**
* Construct a CellId object with a given binary representation that was
void
serialize(Archive &ar, const unsigned int version);
+ /**
+ * Return the id of the coarse cell.
+ */
+ types::coarse_cell_id
+ get_coarse_cell_id() const;
+
private:
/**
* The number of the coarse cell within whose tree the cell
* represented by the current object is located.
*/
- unsigned int coarse_cell_id;
+ types::coarse_cell_id coarse_cell_id;
/**
* The number of child indices stored in the child_indices array. This is
return true; // other.id is longer
}
+
+inline types::coarse_cell_id
+CellId::get_coarse_cell_id() const
+{
+ return coarse_cell_id;
+}
+
+
DEAL_II_NAMESPACE_CLOSE
#endif
std::pair<std::pair<cell_iterator, unsigned int>, std::bitset<3>>> &
get_periodic_face_map() const;
+ /**
+ * Translate the unique id of a coarse cell to its index.
+ *
+ * @note For serial and shared triangulation both id and index are the same.
+ * For distributed triangulations setting both might differ, since the
+ * id might correspond to a global id and the index to a local id.
+ *
+ * @param coarse_cell_id Unique id of the coarse cell.
+ * @return Index of the coarse cell within the current triangulation.
+ */
+ virtual unsigned int
+ coarse_cell_id_to_coarse_cell_index(
+ const types::coarse_cell_id coarse_cell_id) const;
+
+
+ /**
+ * Translate the index of coarse cell to its unique id.
+ *
+ * @note: See the note of the method
+ * translate_coarse_cell_id_to_coarse_cell_index().
+ *
+ * @param coarse_cell_index Index of the coarse cell.
+ * @return Id of the coarse cell.
+ */
+ virtual types::coarse_cell_id
+ coarse_cell_index_to_coarse_cell_id(
+ const unsigned int coarse_cell_index) const;
+
BOOST_SERIALIZATION_SPLIT_MEMBER()
}
+
+template <int dim, int spacedim>
+inline unsigned int
+Triangulation<dim, spacedim>::coarse_cell_id_to_coarse_cell_index(
+ const types::coarse_cell_id coarse_cell_id) const
+{
+ return coarse_cell_id;
+}
+
+
+
+template <int dim, int spacedim>
+inline types::coarse_cell_id
+Triangulation<dim, spacedim>::coarse_cell_index_to_coarse_cell_id(
+ const unsigned int coarse_cell_index) const
+{
+ return coarse_cell_index;
+}
+
+
+
/* -------------- declaration of explicit specializations ------------- */
template <>
CellId::CellId()
- : coarse_cell_id(numbers::invalid_unsigned_int)
+ : coarse_cell_id(numbers::invalid_coarse_cell_id)
, n_child_indices(numbers::invalid_unsigned_int)
{
// initialize the child indices to invalid values
-CellId::CellId(const unsigned int coarse_cell_id,
+CellId::CellId(const types::coarse_cell_id coarse_cell_id,
const std::vector<std::uint8_t> &id)
: coarse_cell_id(coarse_cell_id)
, n_child_indices(id.size())
-CellId::CellId(const unsigned int coarse_cell_id,
- const unsigned int n_child_indices,
- const std::uint8_t *id)
+CellId::CellId(const types::coarse_cell_id coarse_cell_id,
+ const unsigned int n_child_indices,
+ const std::uint8_t * id)
: coarse_cell_id(coarse_cell_id)
, n_child_indices(n_child_indices)
{