16 #include "moonolith_build_quadrature.hpp"
18 using namespace mfem::internal;
25 template <
class Polytope>
class CutGeneric :
public Cut
28 using Point =
typename Polytope::Point;
29 static const int Dim = Point::n_dims;
30 using Quadrature_t = moonolith::Quadrature<double, Dim>;
31 using BuildQuadrature_t = moonolith::BuildQuadrature<Polytope>;
33 bool BuildQuadrature(
const FiniteElementSpace &from_space,
34 const int from_elem_idx,
35 const FiniteElementSpace &to_space,
36 const int to_elem_idx, IntegrationRule &from_quadrature,
37 IntegrationRule &to_quadrature)
override;
39 void Describe()
const override;
42 inline Quadrature_t &GetQRule() {
return q_rule_; }
43 inline int GetOrder()
const {
return order_; }
45 inline void SetOrder(
const int order) { order_ = order; }
47 virtual void MakePolytope(Mesh &mesh,
const int elem_idx,
48 Polytope &polygon) = 0;
50 bool IsValidPhysicalPoint(
const Vector &p_mfem)
const
54 for(
int d = 0; d < Dim; ++d) {
58 moonolith::HPolytope<double, Dim> poly;
61 if(!poly.contains(p, 1e-8))
return false;
65 return poly.contains(p, 1e-8);
68 bool IsValidQPoint(
const IntegrationPoint &p)
const
78 assert(p.x <= 1 + 1e-8);
79 assert(p.y <= 1 + 1e-8);
80 assert(p.z <= 1 + 1e-8);
84 ok = ok && (p.x == p.x);
85 ok = ok && (p.y == p.y);
86 ok = ok && (p.z == p.z);
88 ok = ok && (p.x >= -1e-8);
89 ok = ok && (p.y >= -1e-8);
90 ok = ok && (p.z >= -1e-8);
92 ok = ok && (p.x <= 1 + 1e-8);
93 ok = ok && (p.y <= 1 + 1e-8);
94 ok = ok && (p.z <= 1 + 1e-8);
98 static void ConvertQRule(
const IntegrationRule &ir, Quadrature_t &result)
100 const int size = ir.Size();
104 for (
int k = 0; k < size; ++k)
107 result.points[k][0] = qp.x;
108 result.points[k][1] = qp.y;
109 if constexpr (Dim == 3) {
110 result.points[k][12] = qp.z;
113 result.weights[k] = qp.weight;
117 inline const Polytope & from()
const {
return from_; }
118 inline const Polytope & to()
const {
return to_; }
122 Quadrature_t q_rule_;
123 Quadrature_t physical_quadrature_;
124 BuildQuadrature_t builder_;
126 double intersection_measure_{0};
129 class Cut2D :
public CutGeneric<moonolith::Polygon<double, 2>>
132 using Polygon_t = moonolith::Polygon<double, 2>;
133 using Quadrature_t = moonolith::Quadrature<double, 2>;
134 using BuildQuadrature_t = moonolith::BuildQuadrature<Polygon_t>;
136 void SetIntegrationOrder(
const int order)
override;
139 void SetQuadratureRule(
const IntegrationRule &ir)
override;
140 void MakePolytope(Mesh &mesh,
const int elem_idx,
141 Polygon_t &polygon)
override;
144 DenseMatrix buffer_pts;
147 class Cut3D :
public CutGeneric<moonolith::Polyhedron<double>>
150 using Polyhedron_t = moonolith::Polyhedron<double>;
151 using Quadrature_t = moonolith::Quadrature<double, 3>;
152 using BuildQuadrature_t = moonolith::BuildQuadrature<Polyhedron_t>;
154 void SetIntegrationOrder(
const int order)
override;
157 void SetQuadratureRule(
const IntegrationRule &ir)
override;
158 void MakePolytope(Mesh &mesh,
const int elem_idx,
159 Polyhedron_t &polyhedron)
override;
162 DenseMatrix buffer_pts;
163 Array<int> buffer_vertices;
164 Array<int> buffer_faces, buffer_cor;
168 const Vector &physical_p,
const double &w,
175 assert(ref_p.
x >= -1e-8);
176 assert(ref_p.
y >= -1e-8);
177 assert(ref_p.
z >= -1e-8);
179 assert(ref_p.
x <= 1 + 1e-8);
180 assert(ref_p.
y <= 1 + 1e-8);
181 assert(ref_p.
z <= 1 + 1e-8);
183 #ifdef MFEM_DEBUG_MOONOLITH
188 physical_test_p -= physical_p;
190 assert(physical_test_p.
Norml2() < 1e-10);
206 assert(dim == 2 || dim == 3);
210 template <
class Polytope>
211 bool CutGeneric<Polytope>::BuildQuadrature(
const FiniteElementSpace &from_space,
212 const int from_elem_idx,
213 const FiniteElementSpace &to_space,
214 const int to_elem_idx,
215 IntegrationRule &from_quadrature,
216 IntegrationRule &to_quadrature)
219 MakePolytope(*from_space.GetMesh(), from_elem_idx, from_);
220 MakePolytope(*to_space.GetMesh(), to_elem_idx, to_);
222 if (!builder_.apply(q_rule_, from_, to_, physical_quadrature_))
228 #ifdef MFEM_DEBUG_MOONOLITH
231 static int counter = 0;
233 moonolith::MatlabScripter script;
235 script.plot(from_,
"\'g-\'");
236 script.plot(from_,
"\'g.\'");
237 script.plot(to_,
"\'r-\'");
238 script.plot(to_,
"\'r.\'");
239 script.plot(physical_quadrature_.points,
"\'*b\'");
241 std::cout <<
"measure(" << counter <<
"): " << moonolith::measure(physical_quadrature_) <<
"\n";
243 script.save(
"out_" + std::to_string(counter++) +
".m");
246 #endif // MFEM_DEBUG_MOONOLITH
248 int from_type = from_space.GetFE(from_elem_idx)->GetGeomType();
249 int to_type = to_space.GetFE(to_elem_idx)->GetGeomType();
251 const int n_qp = physical_quadrature_.n_points();
253 from_quadrature.SetSize(n_qp);
254 to_quadrature.SetSize(n_qp);
256 ElementTransformation &from_trans =
257 *from_space.GetElementTransformation(from_elem_idx);
258 ElementTransformation &to_trans =
259 *to_space.GetElementTransformation(to_elem_idx);
261 double from_measure = moonolith::measure(from_);
262 double to_measure = moonolith::measure(to_);
265 for (
int qp = 0; qp < n_qp; ++qp)
268 for (
int d = 0; d < Dim; ++d)
270 p(d) = physical_quadrature_.points[qp][d];
273 double w = physical_quadrature_.weights[qp];
274 intersection_measure_ += w;
276 assert(IsValidPhysicalPoint(p));
279 from_quadrature[qp]);
285 assert(IsValidQPoint(from_quadrature[qp]));
286 assert(IsValidQPoint(to_quadrature[qp]));
292 template <
class Polytope>
void CutGeneric<Polytope>::Describe()
const
294 mfem::out <<
"Cut measure " << intersection_measure_ <<
'\n';
297 template class CutGeneric<::moonolith::Polygon<double, 2>>;
298 template class CutGeneric<::moonolith::Polyhedron<double>>;
300 void Cut2D::MakePolytope(Mesh &mesh,
const int elem_idx, Polygon_t &polygon)
303 mesh.GetPointMatrix(elem_idx, buffer_pts);
305 const int n_points = buffer_pts.Width();
306 polygon.resize(n_points);
308 for (
int k = 0; k < n_points; ++k)
310 for (
int d = 0; d < 2; ++d)
312 polygon.points[k][d] = buffer_pts(d, k);
316 assert(polygon.check_convexity());
317 assert(::moonolith::measure(polygon) > 0.0);
322 void Cut2D::SetQuadratureRule(
const IntegrationRule &ir)
324 const int size = ir.Size();
325 auto &q_rule = this->GetQRule();
330 for (
int k = 0; k < size; ++k)
333 q_rule.points[k][0] = qp.x;
334 q_rule.points[k][1] = qp.y;
335 q_rule.weights[k] = qp.weight;
339 this->SetOrder(ir.GetOrder());
343 void Cut2D::SetIntegrationOrder(
const int order)
345 if (this->GetOrder() != order)
348 assert(ir.GetOrder() >= order);
349 SetQuadratureRule(ir);
353 void Cut3D::SetIntegrationOrder(
const int order)
355 if (this->GetOrder() != order)
358 assert(ir.GetOrder() >= order);
359 SetQuadratureRule(ir);
363 void Cut3D::MakePolytope(Mesh &mesh,
const int elem_idx,
364 Polyhedron_t &polyhedron)
367 const int dim = mesh.Dimension();
369 assert(mesh.GetElement(elem_idx));
371 const Element &e = *mesh.GetElement(elem_idx);
372 const int e_type = e.GetType();
374 mesh.GetElementFaces(elem_idx, buffer_faces, buffer_cor);
375 mesh.GetPointMatrix(elem_idx, buffer_pts);
376 mesh.GetElementVertices(elem_idx, buffer_vertices);
378 const int n_faces = buffer_faces.Size();
380 polyhedron.el_ptr.resize(n_faces + 1);
381 polyhedron.points.resize(buffer_vertices.Size());
382 polyhedron.el_index.resize(MaxVertsXFace(e_type) * n_faces);
383 polyhedron.el_ptr[0] = 0;
385 for (
int i = 0; i < buffer_vertices.Size(); ++i)
387 for (
int j = 0; j <
dim; ++j)
389 polyhedron.points[i][j] = buffer_pts(j, i);
394 for (
int i = 0; i < buffer_faces.Size(); ++i)
396 mesh.GetFaceVertices(buffer_faces[i], f2v);
397 const int eptr = polyhedron.el_ptr[i];
399 for (
int j = 0; j < f2v.Size(); ++j)
401 const int v_offset = buffer_vertices.Find(f2v[j]);
402 polyhedron.el_index[eptr + j] = v_offset;
405 polyhedron.el_ptr[i + 1] = polyhedron.el_ptr[i] + f2v.Size();
408 polyhedron.fix_ordering();
412 void Cut3D::SetQuadratureRule(
const IntegrationRule &ir)
414 const int size = ir.Size();
416 auto &q_rule = this->GetQRule();
422 for (
int k = 0; k < size; ++k)
425 q_rule.points[k][0] = qp.x;
426 q_rule.points[k][1] = qp.y;
427 q_rule.points[k][2] = qp.z;
428 q_rule.weights[k] = qp.weight;
432 this->SetOrder(ir.GetOrder());
440 return std::make_shared<Cut2D>();
444 return std::make_shared<Cut3D>();
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
double Norml2() const
Returns the l2 norm of the vector.
int Size() const
Returns the size of the vector.
double p(const Vector &x, double t)
std::shared_ptr< Cut > NewCut(const int dim)
Class for integration point with weight.
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)