12#ifdef MFEM_USE_MOONOLITH
18#include "moonolith_build_quadrature.hpp"
20using namespace mfem::internal;
27template <
class Polytope>
class CutGeneric :
public Cut
30 using Point =
typename Polytope::Point;
31 static const int Dim = Point::n_dims;
32 using Quadrature_t = moonolith::Quadrature<double, Dim>;
33 using BuildQuadrature_t = moonolith::BuildQuadrature<Polytope>;
36 const int from_elem_idx,
41 void Describe()
const override;
44 inline Quadrature_t &GetQRule() {
return q_rule_; }
45 inline int GetOrder()
const {
return order_; }
47 inline void SetOrder(
const int order) { order_ = order; }
49 virtual void MakePolytope(
Mesh &mesh,
const int elem_idx,
50 Polytope &polygon) = 0;
52 bool IsValidPhysicalPoint(
const Vector &p_mfem)
const
56 for (
int d = 0; d < Dim; ++d)
61 moonolith::HPolytope<double, Dim> poly;
64 if (!poly.contains(
p, 1e-8)) {
return false; }
68 return poly.contains(
p, 1e-8);
81 assert(
p.x <= 1 + 1e-8);
82 assert(
p.y <= 1 + 1e-8);
83 assert(
p.z <= 1 + 1e-8);
87 ok = ok && (
p.x ==
p.x);
88 ok = ok && (
p.y ==
p.y);
89 ok = ok && (
p.z ==
p.z);
91 ok = ok && (
p.x >= -1e-8);
92 ok = ok && (
p.y >= -1e-8);
93 ok = ok && (
p.z >= -1e-8);
95 ok = ok && (
p.x <= 1 + 1e-8);
96 ok = ok && (
p.y <= 1 + 1e-8);
97 ok = ok && (
p.z <= 1 + 1e-8);
101 static void ConvertQRule(
const IntegrationRule &ir, Quadrature_t &result)
103 const int size = ir.
Size();
107 for (
int k = 0; k < size; ++k)
110 result.points[k][0] = qp.x;
111 result.points[k][1] = qp.y;
112 if constexpr (Dim == 3)
114 result.points[k][12] = qp.z;
117 result.weights[k] = qp.weight;
121 inline const Polytope & from()
const {
return from_; }
122 inline const Polytope & to()
const {
return to_; }
126 Quadrature_t q_rule_;
127 Quadrature_t physical_quadrature_;
128 BuildQuadrature_t builder_;
130 double intersection_measure_{0};
133class Cut2D :
public CutGeneric<moonolith::Polygon<double, 2>>
136 using Polygon_t = moonolith::Polygon<double, 2>;
137 using Quadrature_t = moonolith::Quadrature<double, 2>;
138 using BuildQuadrature_t = moonolith::BuildQuadrature<Polygon_t>;
140 void SetIntegrationOrder(
const int order)
override;
144 void MakePolytope(
Mesh &mesh,
const int elem_idx,
145 Polygon_t &polygon)
override;
151class Cut3D :
public CutGeneric<moonolith::Polyhedron<double>>
154 using Polyhedron_t = moonolith::Polyhedron<double>;
155 using Quadrature_t = moonolith::Quadrature<double, 3>;
156 using BuildQuadrature_t = moonolith::BuildQuadrature<Polyhedron_t>;
158 void SetIntegrationOrder(
const int order)
override;
162 void MakePolytope(
Mesh &mesh,
const int elem_idx,
163 Polyhedron_t &polyhedron)
override;
172 const Vector &physical_p,
const double &w,
179 assert(ref_p.
x >= -1e-8);
180 assert(ref_p.
y >= -1e-8);
181 assert(ref_p.
z >= -1e-8);
183 assert(ref_p.
x <= 1 + 1e-8);
184 assert(ref_p.
y <= 1 + 1e-8);
185 assert(ref_p.
z <= 1 + 1e-8);
187#ifdef MFEM_DEBUG_MOONOLITH
192 physical_test_p -= physical_p;
194 assert(physical_test_p.
Norml2() < 1e-10);
210 assert(
dim == 2 ||
dim == 3);
214template <
class Polytope>
216 const int from_elem_idx,
218 const int to_elem_idx,
223 MakePolytope(*from_space.
GetMesh(), from_elem_idx, from_);
224 MakePolytope(*to_space.
GetMesh(), to_elem_idx, to_);
226 if (!builder_.apply(q_rule_, from_, to_, physical_quadrature_))
232#ifdef MFEM_DEBUG_MOONOLITH
235 static int counter = 0;
237 moonolith::MatlabScripter script;
239 script.plot(from_,
"\'g-\'");
240 script.plot(from_,
"\'g.\'");
241 script.plot(to_,
"\'r-\'");
242 script.plot(to_,
"\'r.\'");
243 script.plot(physical_quadrature_.points,
"\'*b\'");
245 mfem::out <<
"measure(" << counter <<
"): " << moonolith::measure(
246 physical_quadrature_) <<
"\n";
248 script.save(
"out_" + std::to_string(counter++) +
".m");
256 const int n_qp = physical_quadrature_.n_points();
266 double from_measure = moonolith::measure(from_);
267 double to_measure = moonolith::measure(to_);
270 for (
int qp = 0; qp < n_qp; ++qp)
273 for (
int d = 0; d < Dim; ++d)
275 p(d) = physical_quadrature_.points[qp][d];
278 double w = physical_quadrature_.weights[qp];
279 intersection_measure_ += w;
281 assert(IsValidPhysicalPoint(
p));
284 from_quadrature[qp]);
290 assert(IsValidQPoint(from_quadrature[qp]));
291 assert(IsValidQPoint(to_quadrature[qp]));
297template <
class Polytope>
void CutGeneric<Polytope>::Describe()
const
299 mfem::out <<
"Cut measure " << intersection_measure_ <<
'\n';
302template class CutGeneric<::moonolith::Polygon<double, 2>>;
303template class CutGeneric<::moonolith::Polyhedron<double>>;
305void Cut2D::MakePolytope(
Mesh &mesh,
const int elem_idx, Polygon_t &polygon)
310 const int n_points = buffer_pts.
Width();
311 polygon.resize(n_points);
313 for (
int k = 0; k < n_points; ++k)
315 for (
int d = 0; d < 2; ++d)
317 polygon.points[k][d] = buffer_pts(d, k);
321 assert(polygon.check_convexity());
322 assert(::moonolith::measure(polygon) > 0.0);
329 const int size = ir.
Size();
330 auto &q_rule = this->GetQRule();
335 for (
int k = 0; k < size; ++k)
338 q_rule.points[k][0] = qp.x;
339 q_rule.points[k][1] = qp.y;
340 q_rule.weights[k] = qp.weight;
348void Cut2D::SetIntegrationOrder(
const int order)
350 if (this->GetOrder() != order)
354 SetQuadratureRule(ir);
358void Cut3D::SetIntegrationOrder(
const int order)
360 if (this->GetOrder() != order)
364 SetQuadratureRule(ir);
368void Cut3D::MakePolytope(
Mesh &mesh,
const int elem_idx,
369 Polyhedron_t &polyhedron)
377 const int e_type = e.
GetType();
383 const int n_faces = buffer_faces.
Size();
385 polyhedron.el_ptr.resize(n_faces + 1);
386 polyhedron.points.resize(buffer_vertices.
Size());
387 polyhedron.el_index.resize(MaxVertsXFace(e_type) * n_faces);
388 polyhedron.el_ptr[0] = 0;
390 for (
int i = 0; i < buffer_vertices.
Size(); ++i)
392 for (
int j = 0; j <
dim; ++j)
394 polyhedron.points[i][j] = buffer_pts(j, i);
399 for (
int i = 0; i < buffer_faces.
Size(); ++i)
402 const int eptr = polyhedron.el_ptr[i];
404 for (
int j = 0; j < f2v.
Size(); ++j)
406 const int v_offset = buffer_vertices.
Find(f2v[j]);
407 polyhedron.el_index[eptr + j] = v_offset;
410 polyhedron.el_ptr[i + 1] = polyhedron.el_ptr[i] + f2v.
Size();
413 polyhedron.fix_ordering();
419 const int size = ir.
Size();
421 auto &q_rule = this->GetQRule();
427 for (
int k = 0; k < size; ++k)
430 q_rule.points[k][0] = qp.x;
431 q_rule.points[k][1] = qp.y;
432 q_rule.points[k][2] = qp.z;
433 q_rule.weights[k] = qp.weight;
445 return std::make_shared<Cut2D>();
449 return std::make_shared<Cut3D>();
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
int Find(const T &el) const
Return the first index where 'el' is found; return -1 if not found.
All subclasses of Cut will implement intersection routines and quadrature point generation within the...
Data type dense matrix using column-major storage.
Abstract data type element.
virtual Type GetType() const =0
Returns element's type.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Mesh * GetMesh() const
Returns the mesh.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetOrder() const
Returns the order of the integration rule.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void GetElementVertices(int i, Array< int > &v) const
Returns the indices of the vertices of element i.
const Element * GetElement(int i) const
Return pointer to the i'th element object.
int Dimension() const
Dimension of the reference space used within the elements.
void GetElementFaces(int i, Array< int > &faces, Array< int > &ori) const
Return the indices and the orientations of all faces of element i.
void GetFaceVertices(int i, Array< int > &vert) const
Returns the indices of the vertices of face i.
void GetPointMatrix(int i, DenseMatrix &pointmat) const
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
real_t Norml2() const
Returns the l2 norm of the vector.
int Size() const
Returns the size of the vector.
std::shared_ptr< Cut > NewCut(const int dim)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void TransformToReference(ElementTransformation &Trans, int type, const Vector &physical_p, const double &w, IntegrationPoint &ref_p)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p(const Vector &x, real_t t)