From: Peter Munch Date: Wed, 8 Dec 2021 15:27:34 +0000 (+0100) Subject: Add Triangulation::all_reference_cells_are_simplex() and ::is_mixed_mesh() X-Git-Tag: v9.4.0-rc1~763^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=321eb1c4d8935f1e3e93e22be930cae7a27deb30;p=dealii.git Add Triangulation::all_reference_cells_are_simplex() and ::is_mixed_mesh() --- diff --git a/include/deal.II/grid/tria.h b/include/deal.II/grid/tria.h index 265b70d0ae..90d29f3379 100644 --- a/include/deal.II/grid/tria.h +++ b/include/deal.II/grid/tria.h @@ -3416,6 +3416,21 @@ public: bool all_reference_cells_are_hyper_cube() const; + /** + * Indicate if the triangulation only consists of simplex-like cells, i.e., + * lines, triangles, or tetrahedra. + */ + bool + all_reference_cells_are_simplex() const; + + /** + * Indicate if the triangulation consists of different cell types (mix of + * simplices, hypercubes, ...) or different face types, as in the case + * of pyramids or wedges.. + */ + bool + is_mixed_mesh() const; + #ifdef DOXYGEN /** * Write and read the data of this object from a stream for the purpose diff --git a/source/grid/tria.cc b/source/grid/tria.cc index 57fce4ad70..d279c78bb6 100644 --- a/source/grid/tria.cc +++ b/source/grid/tria.cc @@ -14655,7 +14655,36 @@ Triangulation::all_reference_cells_are_hyper_cube() const "cells used by this triangulation if the " "triangulation doesn't yet have any cells in it.")); return (this->reference_cells.size() == 1 && - this->reference_cells[0] == ReferenceCells::get_hypercube()); + this->reference_cells[0].is_hyper_cube()); +} + + + +template +bool +Triangulation::all_reference_cells_are_simplex() const +{ + Assert(this->reference_cells.size() > 0, + ExcMessage("You can't ask about the kinds of reference " + "cells used by this triangulation if the " + "triangulation doesn't yet have any cells in it.")); + return (this->reference_cells.size() == 1 && + this->reference_cells[0].is_simplex()); +} + + + +template +bool +Triangulation::is_mixed_mesh() const +{ + Assert(this->reference_cells.size() > 0, + ExcMessage("You can't ask about the kinds of reference " + "cells used by this triangulation if the " + "triangulation doesn't yet have any cells in it.")); + return reference_cells.size() > 1 || + ((reference_cells[0].is_hyper_cube() == false) && + (reference_cells[0].is_simplex() == false)); }