20#ifdef DEAL_II_WITH_OPENCASCADE
26# include <IGESControl_Controller.hxx>
27# include <IGESControl_Reader.hxx>
28# include <IGESControl_Writer.hxx>
29# include <STEPControl_Controller.hxx>
30# include <STEPControl_Reader.hxx>
31# include <STEPControl_Writer.hxx>
32# include <TopExp_Explorer.hxx>
34# include <TopoDS_Edge.hxx>
35# include <TopoDS_Face.hxx>
36# include <TopoDS_Shape.hxx>
40# if DEAL_II_OPENCASCADE_VERSION_GTE(7, 0, 0)
41# include <Standard_Transient.hxx>
43# include <Handle_Standard_Transient.hxx>
46# include <BRepAdaptor_Curve.hxx>
47# if DEAL_II_OPENCASCADE_VERSION_GTE(7, 6, 0)
48# include <BRepAlgoAPI_Section.hxx>
50# include <BRepAdaptor_HCompCurve.hxx>
51# include <BRepAdaptor_HCurve.hxx>
52# include <BRepAlgo_Section.hxx>
54# include <BRepAdaptor_Surface.hxx>
55# include <BRepBndLib.hxx>
56# include <BRepBuilderAPI_MakeEdge.hxx>
57# include <BRepBuilderAPI_Sewing.hxx>
58# include <BRepBuilderAPI_Transform.hxx>
59# include <BRepMesh_IncrementalMesh.hxx>
60# include <BRepTools.hxx>
61# include <BRep_Builder.hxx>
62# include <GCPnts_AbscissaPoint.hxx>
63# include <GeomAPI_Interpolate.hxx>
64# include <GeomAPI_ProjectPointOnCurve.hxx>
65# include <GeomAPI_ProjectPointOnSurf.hxx>
66# include <GeomConvert_CompCurveToBSplineCurve.hxx>
67# include <GeomLProp_SLProps.hxx>
68# include <Geom_BoundedCurve.hxx>
69# include <Geom_Plane.hxx>
70# include <IntCurvesFace_ShapeIntersector.hxx>
71# include <Poly_Triangulation.hxx>
72# include <ShapeAnalysis_Surface.hxx>
73# include <StlAPI_Reader.hxx>
74# include <StlAPI_Writer.hxx>
75# include <TColStd_HSequenceOfTransient.hxx>
76# include <TColStd_SequenceOfTransient.hxx>
77# include <TColgp_HArray1OfPnt.hxx>
78# include <TopLoc_Location.hxx>
91 std::tuple<unsigned int, unsigned int, unsigned int>
95 unsigned int n_faces = 0, n_edges = 0, n_vertices = 0;
96 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next(), ++n_faces)
98 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next(), ++n_edges)
100 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next(), ++n_vertices)
102 return std::tuple<unsigned int, unsigned int, unsigned int>(n_faces,
109 std::vector<TopoDS_Face> & faces,
110 std::vector<TopoDS_Edge> & edges,
111 std::vector<TopoDS_Vertex> &
vertices)
118 for (exp.Init(shape, TopAbs_FACE); exp.More(); exp.Next())
120 faces.push_back(TopoDS::Face(exp.Current()));
122 for (exp.Init(shape, TopAbs_EDGE); exp.More(); exp.Next())
124 edges.push_back(TopoDS::Edge(exp.Current()));
126 for (exp.Init(shape, TopAbs_VERTEX); exp.More(); exp.Next())
128 vertices.push_back(TopoDS::Vertex(exp.Current()));
135 std::vector<TopoDS_Compound> & compounds,
136 std::vector<TopoDS_CompSolid> &compsolids,
137 std::vector<TopoDS_Solid> & solids,
138 std::vector<TopoDS_Shell> & shells,
139 std::vector<TopoDS_Wire> & wires)
142 compsolids.resize(0);
148 for (exp.Init(shape, TopAbs_COMPOUND); exp.More(); exp.Next())
150 compounds.push_back(TopoDS::Compound(exp.Current()));
152 for (exp.Init(shape, TopAbs_COMPSOLID); exp.More(); exp.Next())
154 compsolids.push_back(TopoDS::CompSolid(exp.Current()));
156 for (exp.Init(shape, TopAbs_SOLID); exp.More(); exp.Next())
158 solids.push_back(TopoDS::Solid(exp.Current()));
160 for (exp.Init(shape, TopAbs_SHELL); exp.More(); exp.Next())
162 shells.push_back(TopoDS::Shell(exp.Current()));
164 for (exp.Init(shape, TopAbs_WIRE); exp.More(); exp.Next())
166 wires.push_back(TopoDS::Wire(exp.Current()));
170 template <
int spacedim>
177 return gp_Pnt(p[0], 0, 0);
179 return gp_Pnt(p[0], p[1], 0);
181 return gp_Pnt(p[0], p[1], p[2]);
187 template <
int spacedim>
189 point(
const gp_Pnt &p,
const double tolerance)
197 "Cannot convert OpenCASCADE point to 1d if p.Y() != 0."));
200 "Cannot convert OpenCASCADE point to 1d if p.Z() != 0."));
205 "Cannot convert OpenCASCADE point to 2d if p.Z() != 0."));
219 const double tolerance)
221 const double rel_tol =
223 if (direction.
norm() > 0.0)
224 return (p1 * direction < p2 * direction - rel_tol);
226 for (
int d = dim; d >= 0; --d)
227 if (p1[d] < p2[d] - rel_tol)
229 else if (p2[d] < p1[d] - rel_tol)
239 read_IGES(
const std::string &filename,
const double scale_factor)
241 IGESControl_Reader reader;
242 IFSelect_ReturnStatus stat;
243 stat = reader.ReadFile(filename.c_str());
246 Standard_Boolean failsonly = Standard_False;
247 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
248 reader.PrintCheckLoad(failsonly, mode);
250 Standard_Integer nRoots = reader.TransferRoots();
259 scale.SetScale(Origin, scale_factor);
261 TopoDS_Shape sh = reader.OneShape();
262 BRepBuilderAPI_Transform trans(sh, scale);
264 return trans.Shape();
268 write_IGES(
const TopoDS_Shape &shape,
const std::string &filename)
270 IGESControl_Controller::Init();
271 IGESControl_Writer ICW(
"MM", 0);
272 Standard_Boolean ok = ICW.AddShape(shape);
275 Standard_Boolean OK = ICW.Write(filename.c_str());
283 StlAPI_Reader reader;
285 reader.Read(shape, filename.c_str());
292 const std::string & filename,
293 const double deflection,
294 const bool sew_different_faces,
295 const double sewer_tolerance,
296 const bool is_relative,
297 const double angular_deflection,
298 const bool in_parallel)
301 std::vector<TopoDS_Vertex>
vertices;
302 std::vector<TopoDS_Edge> edges;
303 std::vector<TopoDS_Face> faces;
305 const bool mesh_is_present =
306 std::none_of(faces.begin(), faces.end(), [&Loc](
const TopoDS_Face &face) {
307 Handle(Poly_Triangulation) theTriangulation =
308 BRep_Tool::Triangulation(face, Loc);
309 return theTriangulation.IsNull();
311 TopoDS_Shape shape_to_be_written = shape;
312 if (!mesh_is_present)
314 if (sew_different_faces)
316 BRepBuilderAPI_Sewing sewer(sewer_tolerance);
317 sewer.Add(shape_to_be_written);
319 shape_to_be_written = sewer.SewedShape();
322 shape_to_be_written = shape;
326 BRepMesh_IncrementalMesh mesh_im(shape_to_be_written,
333 StlAPI_Writer writer;
335# if DEAL_II_OPENCASCADE_VERSION_GTE(6, 9, 0)
337 const auto error = writer.Write(shape_to_be_written, filename.c_str());
340# if !DEAL_II_OPENCASCADE_VERSION_GTE(7, 2, 0)
351 writer.Write(shape_to_be_written, filename.c_str());
356 read_STEP(
const std::string &filename,
const double scale_factor)
358 STEPControl_Reader reader;
359 IFSelect_ReturnStatus stat;
360 stat = reader.ReadFile(filename.c_str());
363 Standard_Boolean failsonly = Standard_False;
364 IFSelect_PrintCount mode = IFSelect_ItemsByEntity;
365 reader.PrintCheckLoad(failsonly, mode);
367 Standard_Integer nRoots = reader.TransferRoots();
376 scale.SetScale(Origin, scale_factor);
378 TopoDS_Shape sh = reader.OneShape();
379 BRepBuilderAPI_Transform trans(sh, scale);
381 return trans.Shape();
385 write_STEP(
const TopoDS_Shape &shape,
const std::string &filename)
387 STEPControl_Controller::Init();
388 STEPControl_Writer SCW;
389 IFSelect_ReturnStatus status;
390 status = SCW.Transfer(shape, STEPControl_AsIs);
392 ExcMessage(
"Failed to add shape to STEP controller."));
394 status = SCW.Write(filename.c_str());
397 ExcMessage(
"Failed to write translated shape to STEP file."));
403 double tolerance = 0.0;
405 std::vector<TopoDS_Face> faces;
406 std::vector<TopoDS_Edge> edges;
407 std::vector<TopoDS_Vertex>
vertices;
412 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(vertex));
414 for (
const auto &edge : edges)
415 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(edge));
417 for (
const auto &face : faces)
418 tolerance = std::fmax(tolerance, BRep_Tool::Tolerance(face));
432 Handle(Geom_Plane) plane =
new Geom_Plane(c_x, c_y, c_z, c);
433# if DEAL_II_OPENCASCADE_VERSION_GTE(7, 6, 0)
434 BRepAlgoAPI_Section section(in_shape, plane);
436 BRepAlgo_Section section(in_shape, plane);
438 TopoDS_Shape edges = section.Shape();
443 join_edges(
const TopoDS_Shape &in_shape,
const double tolerance)
445 TopoDS_Edge out_shape;
446 const TopoDS_Shape & edges = in_shape;
447 std::vector<Handle_Geom_BoundedCurve> intersections;
451 gp_Pnt PIn(0.0, 0.0, 0.0);
452 gp_Pnt PFin(0.0, 0.0, 0.0);
453 gp_Pnt PMid(0.0, 0.0, 0.0);
454 TopExp_Explorer edgeExplorer(edges, TopAbs_EDGE);
456 while (edgeExplorer.More())
458 edge = TopoDS::Edge(edgeExplorer.Current());
459 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, L, First, Last);
460 intersections.push_back(Handle(Geom_BoundedCurve)::DownCast(curve));
466 unsigned int numIntersEdges = intersections.size();
469 GeomConvert_CompCurveToBSplineCurve convert_bspline(intersections[0]);
471 bool check =
false, one_added =
true, one_failed =
true;
472 std::vector<bool> added(numIntersEdges,
false);
474 while (one_added ==
true)
478 for (
unsigned int i = 1; i < numIntersEdges; ++i)
479 if (added[i] ==
false)
481 Handle(Geom_Curve) curve = intersections[i];
482 Handle(Geom_BoundedCurve) bcurve =
483 Handle(Geom_BoundedCurve)::DownCast(curve);
484 check = convert_bspline.Add(bcurve, tolerance,
false,
true, 0);
489 Handle(Geom_BoundedCurve) bcurve =
490 Handle(Geom_BoundedCurve)::DownCast(curve);
492 convert_bspline.Add(bcurve, tolerance,
false,
true, 0);
494 one_failed = one_failed || (check ==
false);
495 one_added = one_added || (check ==
true);
500 Assert(one_failed ==
false,
501 ExcMessage(
"Joining some of the Edges failed."));
503 Handle(Geom_Curve) bspline = convert_bspline.BSplineCurve();
505 out_shape = BRepBuilderAPI_MakeEdge(bspline);
514 const double tolerance)
518 gp_Pnt P0 =
point(origin);
521 dim > 1 ? direction[1] : 0,
522 dim > 2 ? direction[2] : 0));
526 gp_Pnt Pproj(0.0, 0.0, 0.0);
530 IntCurvesFace_ShapeIntersector Inters;
531 Inters.Load(in_shape, tolerance);
535 Inters.Perform(line, -RealLast(), +RealLast());
538 double minDistance = 1e7;
540 for (
int i = 0; i < Inters.NbPnt(); ++i)
542 const double distance =
point(origin).Distance(Inters.Pnt(i + 1));
545 if (distance < minDistance)
547 minDistance = distance;
548 result = point<dim>(Inters.Pnt(i + 1));
560 const double tolerance)
562 unsigned int n_vertices = curve_points.size();
564 if (direction * direction > 0)
566 std::sort(curve_points.begin(),
569 return OpenCASCADE::point_compare(p1,
577 Handle(TColgp_HArray1OfPnt)
vertices =
578 new TColgp_HArray1OfPnt(1, n_vertices);
579 for (
unsigned int vertex = 0; vertex < n_vertices; ++vertex)
585 GeomAPI_Interpolate bspline_generator(
vertices, closed, tolerance);
586 bspline_generator.Perform();
587 Assert((bspline_generator.IsDone()),
588 ExcMessage(
"Interpolated bspline generation failed"));
590 Handle(Geom_BSplineCurve) bspline = bspline_generator.Curve();
591 TopoDS_Edge out_shape = BRepBuilderAPI_MakeEdge(bspline);
592 out_shape.Closed(closed);
598 template <
int spacedim>
599 std::vector<TopoDS_Edge>
607 std::map<unsigned int, std::pair<unsigned int, unsigned int>> vert_to_faces;
608 std::map<unsigned int, std::pair<unsigned int, unsigned int>> face_to_verts;
609 std::map<unsigned int, bool> visited_faces;
610 std::map<unsigned int, Point<spacedim>> vert_to_point;
612 unsigned int face_index;
614 for (
const auto &cell :
triangulation.active_cell_iterators())
616 if (cell->face(f)->at_boundary())
619 face_index = cell->face(f)->index();
620 const unsigned int v0 = cell->face(f)->vertex_index(0);
621 const unsigned int v1 = cell->face(f)->vertex_index(1);
622 face_to_verts[face_index].first =
v0;
623 face_to_verts[face_index].second =
v1;
624 visited_faces[face_index] =
false;
629 f, 0,
true,
false,
false)];
631 f, 1,
true,
false,
false)];
634 if (vert_to_faces.find(
v0) == vert_to_faces.end())
636 vert_to_faces[
v0].first = face_index;
640 vert_to_faces[
v0].second = face_index;
642 if (vert_to_faces.find(
v1) == vert_to_faces.end())
644 vert_to_faces[
v1].first = face_index;
648 vert_to_faces[
v1].second = face_index;
654 std::vector<TopoDS_Edge> interpolation_curves;
655 bool finished = (face_to_verts.size() == 0);
656 face_index = finished ? 0 : face_to_verts.begin()->first;
658 while (finished ==
false)
660 const unsigned int start_point_index = face_to_verts[face_index].first;
661 unsigned int point_index = start_point_index;
664 std::vector<Point<spacedim>> pointlist;
667 visited_faces[face_index] =
true;
668 auto current_point = vert_to_point[point_index];
669 pointlist.push_back(current_point);
672 if (face_to_verts[face_index].
first != point_index)
673 point_index = face_to_verts[face_index].first;
675 point_index = face_to_verts[face_index].second;
678 if (vert_to_faces[point_index].
first != face_index)
679 face_index = vert_to_faces[point_index].first;
681 face_index = vert_to_faces[point_index].second;
683 while (point_index != start_point_index);
685 interpolation_curves.push_back(
689 for (
const auto &f : visited_faces)
690 if (f.second ==
false)
692 face_index = f.first;
697 return interpolation_curves;
702 std::tuple<Point<dim>, TopoDS_Shape, double,
double>
705 const double tolerance)
708 gp_Pnt Pproj =
point(origin);
710 double minDistance = 1e7;
711 gp_Pnt tmp_proj(0.0, 0.0, 0.0);
713# ifdef DEAL_II_HAVE_CXX17
714 [[maybe_unused]]
unsigned int counter = 0;
716 unsigned int counter = 0;
719 unsigned int face_counter = 0;
721 TopoDS_Shape out_shape;
725 for (exp.Init(in_shape, TopAbs_FACE); exp.More(); exp.Next())
727 TopoDS_Face face = TopoDS::Face(exp.Current());
731 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
733 ShapeAnalysis_Surface projector(SurfToProj);
734 gp_Pnt2d proj_params = projector.ValueOfUV(
point(origin), tolerance);
736 SurfToProj->D0(proj_params.X(), proj_params.Y(), tmp_proj);
738 double distance = point<dim>(tmp_proj).distance(origin);
739 if (distance < minDistance)
741 minDistance = distance;
755 if (face_counter == 0)
756 for (exp.Init(in_shape, TopAbs_EDGE); exp.More(); exp.Next())
758 TopoDS_Edge edge = TopoDS::Edge(exp.Current());
759 if (!BRep_Tool::Degenerated(edge))
767 Handle(Geom_Curve) CurveToProj =
768 BRep_Tool::Curve(edge, L, First, Last);
770 GeomAPI_ProjectPointOnCurve Proj(
point(origin), CurveToProj);
771 unsigned int num_proj_points = Proj.NbPoints();
772 if ((num_proj_points > 0) && (Proj.LowerDistance() < minDistance))
774 minDistance = Proj.LowerDistance();
775 Pproj = Proj.NearestPoint();
777 u = Proj.LowerDistanceParameter();
784 return std::tuple<Point<dim>, TopoDS_Shape, double,
double>(
785 point<dim>(Pproj), out_shape, u, v);
793 const double tolerance)
795 std::tuple<Point<dim>, TopoDS_Shape, double,
double> ref =
797 return std::get<0>(ref);
803 const double tolerance)
806 std::tuple<Point<3>, TopoDS_Shape, double,
double> shape_and_params =
809 TopoDS_Shape &out_shape = std::get<1>(shape_and_params);
810 double & u = std::get<2>(shape_and_params);
811 double & v = std::get<3>(shape_and_params);
815 std::tuple<unsigned int, unsigned int, unsigned int>
numbers =
822 "Could not find normal: the shape containing the closest point has 0 faces."));
826 "Could not find normal: the shape containing the closest point has more than 1 face."));
830 exp.Init(out_shape, TopAbs_FACE);
831 TopoDS_Face face = TopoDS::Face(exp.Current());
837 push_forward(
const TopoDS_Shape &in_shape,
const double u,
const double v)
839 switch (in_shape.ShapeType())
843 BRepAdaptor_Surface surf(TopoDS::Face(in_shape));
844 return point<dim>(surf.Value(u, v));
848 BRepAdaptor_Curve curve(TopoDS::Edge(in_shape));
849 return point<dim>(curve.Value(u));
863 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
864 GeomLProp_SLProps props(SurfToProj, u, v, 1, 1e-7);
865 gp_Pnt Value = props.Value();
867 gp_Dir Normal = props.Normal();
868 Assert(props.IsCurvatureDefined(),
869 ExcMessage(
"Curvature is not well defined!"));
870 Standard_Real Min_Curvature = props.MinCurvature();
871 Standard_Real Max_Curvature = props.MaxCurvature();
878 if (face.Orientation() == TopAbs_REVERSED)
885 return std::tuple<Point<3>,
Tensor<1, 3>, double,
double>(point<3>(Value),
893 template <
int spacedim>
898 BRepAdaptor_Surface surf(face);
899 const double u0 = surf.FirstUParameter();
900 const double u1 = surf.LastUParameter();
901 const double v0 = surf.FirstVParameter();
902 const double v1 = surf.LastVParameter();
904 std::vector<CellData<2>> cells;
905 std::vector<Point<spacedim>>
vertices;
908 vertices.push_back(point<spacedim>(surf.Value(u0,
v0)));
909 vertices.push_back(point<spacedim>(surf.Value(u1,
v0)));
910 vertices.push_back(point<spacedim>(surf.Value(u0,
v1)));
911 vertices.push_back(point<spacedim>(surf.Value(u1,
v1)));
914 for (
unsigned int i = 0; i < 4; ++i)
917 cells.push_back(cell);
921# include "utilities.inst"
Abstract base class for mapping classes.
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
numbers::NumberTraits< Number >::real_type norm() const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcUnsupportedShape()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
TopoDS_Edge interpolation_curve(std::vector< Point< dim > > &curve_points, const Tensor< 1, dim > &direction=Tensor< 1, dim >(), const bool closed=false, const double tolerance=1e-7)
std::tuple< Point< 3 >, Tensor< 1, 3 >, double, double > push_forward_and_differential_forms(const TopoDS_Face &face, const double u, const double v, const double tolerance=1e-7)
void extract_compound_shapes(const TopoDS_Shape &shape, std::vector< TopoDS_Compound > &compounds, std::vector< TopoDS_CompSolid > &compsolids, std::vector< TopoDS_Solid > &solids, std::vector< TopoDS_Shell > &shells, std::vector< TopoDS_Wire > &wires)
Point< dim > push_forward(const TopoDS_Shape &in_shape, const double u, const double v)
std::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
void write_STEP(const TopoDS_Shape &shape, const std::string &filename)
std::tuple< Point< 3 >, Tensor< 1, 3 >, double, double > closest_point_and_differential_forms(const TopoDS_Shape &in_shape, const Point< 3 > &origin, const double tolerance=1e-7)
void extract_geometrical_shapes(const TopoDS_Shape &shape, std::vector< TopoDS_Face > &faces, std::vector< TopoDS_Edge > &edges, std::vector< TopoDS_Vertex > &vertices)
TopoDS_Edge join_edges(const TopoDS_Shape &in_shape, const double tolerance=1e-7)
TopoDS_Shape read_STEP(const std::string &filename, const double scale_factor=1e-3)
bool point_compare(const Point< dim > &p1, const Point< dim > &p2, const Tensor< 1, dim > &direction=Tensor< 1, dim >(), const double tolerance=1e-10)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Point< dim > closest_point(const TopoDS_Shape &in_shape, const Point< dim > &origin, const double tolerance=1e-7)
std::tuple< Point< dim >, TopoDS_Shape, double, double > project_point_and_pull_back(const TopoDS_Shape &in_shape, const Point< dim > &origin, const double tolerance=1e-7)
Point< dim > line_intersection(const TopoDS_Shape &in_shape, const Point< dim > &origin, const Tensor< 1, dim > &direction, const double tolerance=1e-7)
void write_STL(const TopoDS_Shape &shape, const std::string &filename, const double deflection, const bool sew_different_faces=false, const double sewer_tolerance=1e-6, const bool is_relative=false, const double angular_deflection=0.5, const bool in_parallel=false)
void create_triangulation(const TopoDS_Face &face, Triangulation< 2, spacedim > &tria)
TopoDS_Shape read_STL(const std::string &filename)
TopoDS_Shape intersect_plane(const TopoDS_Shape &in_shape, const double c_x, const double c_y, const double c_z, const double c, const double tolerance=1e-7)
std::vector< TopoDS_Edge > create_curves_from_triangulation_boundary(const Triangulation< 2, spacedim > &triangulation, const Mapping< 2, spacedim > &mapping=StaticMappingQ1< 2, spacedim >::mapping)
double get_shape_tolerance(const TopoDS_Shape &shape)
void write_IGES(const TopoDS_Shape &shape, const std::string &filename)
TopoDS_Shape read_IGES(const std::string &filename, const double scale_factor=1e-3)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
std::vector< unsigned int > vertices
static unsigned int face_to_cell_vertices(const unsigned int face, const unsigned int vertex, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
const ::Triangulation< dim, spacedim > & tria