17 #include <deal.II/opencascade/boundary_lib.h> 19 #ifdef DEAL_II_WITH_OPENCASCADE 22 #include <GCPnts_AbscissaPoint.hxx> 23 #include <BRepAdaptor_Curve.hxx> 24 #include <BRepAdaptor_CompCurve.hxx> 25 #include <BRepAdaptor_HCurve.hxx> 26 #include <BRepAdaptor_HCompCurve.hxx> 27 #include <GCPnts_AbscissaPoint.hxx> 28 #include <ShapeAnalysis_Curve.hxx> 29 #include <BRep_Tool.hxx> 30 #include <BRepTools.hxx> 31 #include <ShapeAnalysis_Surface.hxx> 34 #include <Standard_Version.hxx> 35 #if (OCC_VERSION_MAJOR < 7) 36 # include <Handle_Adaptor3d_HCurve.hxx> 40 DEAL_II_NAMESPACE_OPEN
54 Handle_Adaptor3d_HCurve curve_adaptor(
const TopoDS_Shape &shape)
56 Assert( (shape.ShapeType() == TopAbs_WIRE) ||
57 (shape.ShapeType() == TopAbs_EDGE),
59 if (shape.ShapeType() == TopAbs_WIRE)
60 return Handle(BRepAdaptor_HCompCurve)(
new BRepAdaptor_HCompCurve(TopoDS::Wire(shape)));
61 else if (shape.ShapeType() == TopAbs_EDGE)
62 return Handle(BRepAdaptor_HCurve)(
new BRepAdaptor_HCurve(TopoDS::Edge(shape)));
65 return Handle(BRepAdaptor_HCurve)(
new BRepAdaptor_HCurve());
71 double shape_length(
const TopoDS_Shape &sh)
73 Handle_Adaptor3d_HCurve adapt = curve_adaptor(sh);
74 return GCPnts_AbscissaPoint::Length(adapt->GetCurve());
79 template <
int dim,
int spacedim>
81 const double tolerance) :
90 template<
int dim,
int spacedim>
91 std::unique_ptr<Manifold<dim, spacedim> >
99 template <
int dim,
int spacedim>
104 (void)surrounding_points;
106 for (
unsigned int i=0; i<surrounding_points.size(); ++i)
109 std::max(tolerance*surrounding_points[i].norm(), tolerance),
110 ExcPointNotOnManifold<spacedim>(surrounding_points[i]));
117 template <
int dim,
int spacedim>
120 const double tolerance) :
122 direction(direction),
130 template<
int dim,
int spacedim>
131 std::unique_ptr<Manifold<dim, spacedim> >
134 return std::unique_ptr<Manifold<dim,spacedim> >
140 template <
int dim,
int spacedim>
145 (void)surrounding_points;
147 for (
unsigned int i=0; i<surrounding_points.size(); ++i)
150 std::max(tolerance*surrounding_points[i].norm(), tolerance),
151 ExcPointNotOnManifold<spacedim>(surrounding_points[i]));
159 template <
int dim,
int spacedim>
161 const double tolerance) :
167 ExcMessage(
"NormalToMeshProjectionBoundary needs a shape containing faces to operate."));
170 template<
int dim,
int spacedim>
171 std::unique_ptr<Manifold<dim, spacedim> >
174 return std::unique_ptr<Manifold<dim, spacedim> >
179 template <
int dim,
int spacedim>
184 TopoDS_Shape out_shape;
187 for (
unsigned int i=0; i<surrounding_points.size(); ++i)
190 .distance(surrounding_points[i]) <
191 std::max(tolerance*surrounding_points[i].norm(), tolerance),
192 ExcPointNotOnManifold<spacedim>(surrounding_points[i]));
196 switch (surrounding_points.size())
200 for (
unsigned int i=0; i<surrounding_points.size(); ++i)
205 surrounding_points[i],
207 average_normal += std::get<1>(p_and_diff_forms);
213 ExcMessage(
"Failed to refine cell: the average of the surface normals at the surrounding edge turns out to be a null vector, making the projection direction undetermined."));
215 Tensor<1,3> T = surrounding_points[0]-surrounding_points[1];
217 average_normal = average_normal-(average_normal*T)*T;
218 average_normal /= average_normal.
norm();
223 Tensor<1,3> u = surrounding_points[1]-surrounding_points[0];
224 Tensor<1,3> v = surrounding_points[2]-surrounding_points[0];
225 const double n1_coords[3] = {u[1] *v[2]-u[2] *v[1],u[2] *v[0]-u[0] *v[2],u[0] *v[1]-u[1] *v[0]};
228 u = surrounding_points[2]-surrounding_points[3];
229 v = surrounding_points[1]-surrounding_points[3];
230 const double n2_coords[3] = {u[1] *v[2]-u[2] *v[1],u[2] *v[0]-u[0] *v[2],u[0] *v[1]-u[1] *v[0]};
234 average_normal = (n1+n2)/2.0;
237 ExcMessage(
"Failed to refine cell: the normal estimated via the surrounding points turns out to be a null vector, making the projection direction undetermined."));
239 average_normal /= average_normal.
norm();
244 Tensor<1,3> u = surrounding_points[1]-surrounding_points[0];
245 Tensor<1,3> v = surrounding_points[2]-surrounding_points[0];
246 const double n1_coords[3] = {u[1] *v[2]-u[2] *v[1],u[2] *v[0]-u[0] *v[2],u[0] *v[1]-u[1] *v[0]};
249 u = surrounding_points[2]-surrounding_points[3];
250 v = surrounding_points[1]-surrounding_points[3];
251 const double n2_coords[3] = {u[1] *v[2]-u[2] *v[1],u[2] *v[0]-u[0] *v[2],u[0] *v[1]-u[1] *v[0]};
254 u = surrounding_points[4]-surrounding_points[7];
255 v = surrounding_points[6]-surrounding_points[7];
256 const double n3_coords[3] = {u[1] *v[2]-u[2] *v[1],u[2] *v[0]-u[0] *v[2],u[0] *v[1]-u[1] *v[0]};
259 u = surrounding_points[6]-surrounding_points[7];
260 v = surrounding_points[5]-surrounding_points[7];
261 const double n4_coords[3] = {u[1] *v[2]-u[2] *v[1],u[2] *v[0]-u[0] *v[2],u[0] *v[1]-u[1] *v[0]};
265 average_normal = (n1+n2+n3+n4)/4.0;
268 ExcMessage(
"Failed to refine cell: the normal estimated via the surrounding points turns out to be a null vector, making the projection direction undetermined."));
270 average_normal /= average_normal.
norm();
285 template <
int dim,
int spacedim>
287 const double tolerance):
290 Point<1>(shape_length(sh)) :
293 curve(curve_adaptor(sh)),
294 tolerance(tolerance),
295 length(shape_length(sh))
302 template<
int dim,
int spacedim>
303 std::unique_ptr<Manifold<dim, spacedim> >
306 return std::unique_ptr<Manifold<dim, spacedim> >
312 template <
int dim,
int spacedim>
317 ShapeAnalysis_Curve curve_analysis;
319 const double dist = curve_analysis.Project(curve->GetCurve(),
point(space_point), tolerance, proj, t,
true);
320 Assert(dist < tolerance*length, ExcPointNotOnManifold<spacedim>(space_point));
322 return Point<1>(GCPnts_AbscissaPoint::Length(curve->GetCurve(),curve->GetCurve().FirstParameter(),t));
327 template <
int dim,
int spacedim>
331 GCPnts_AbscissaPoint AP(curve->GetCurve(), chart_point[0], curve->GetCurve().FirstParameter());
332 gp_Pnt P = curve->GetCurve().Value(AP.Parameter());
333 return point<spacedim>(P);
336 template <
int dim,
int spacedim>
339 const double tolerance)
347 template<
int dim,
int spacedim>
348 std::unique_ptr<Manifold<dim, spacedim> >
351 return std::unique_ptr<Manifold<dim, spacedim> >
357 template <
int dim,
int spacedim>
Point<2> 361 Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
363 ShapeAnalysis_Surface projector(SurfToProj);
364 gp_Pnt2d proj_params = projector.ValueOfUV(
point(space_point), tolerance);
366 double u = proj_params.X();
367 double v = proj_params.Y();
372 template <
int dim,
int spacedim>
377 return ::OpenCASCADE::push_forward<spacedim>(face, chart_point[0], chart_point[1]);
380 template <
int dim,
int spacedim>
386 Handle(Geom_Surface) surf = BRep_Tool::Surface(face);
390 surf->D1(chart_point[0],chart_point[1], q, Du, Dv);
397 Assert(std::abs(Du.Z()) < tolerance,
398 ExcMessage(
"Expecting derivative along Z to be zero! Bailing out."));
404 Assert(std::abs(Dv.Z()) < tolerance,
405 ExcMessage(
"Expecting derivative along Z to be zero! Bailing out."));
409 template <
int dim,
int spacedim>
410 std::tuple<double, double, double, double>
414 Standard_Real umin, umax, vmin, vmax;
415 BRepTools::UVBounds(face, umin, umax, vmin, vmax);
416 return std::make_tuple(umin, umax, vmin, vmax);
420 #include "boundary_lib.inst" 424 DEAL_II_NAMESPACE_CLOSE
virtual Point< 2 > pull_back(const Point< spacedim > &space_point) const override
std::tuple< unsigned int, unsigned int, unsigned int > count_elements(const TopoDS_Shape &shape)
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
static ::ExceptionBase & ExcImpossibleInDimSpacedim(int arg1, int arg2)
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Point< spacedim > point(const gp_Pnt &p, const double &tolerance=1e-10)
virtual Point< spacedim > push_forward(const Point< 1 > &chart_point) const override
Point< dim > line_intersection(const TopoDS_Shape &in_shape, const Point< dim > &origin, const Tensor< 1, dim > &direction, const double tolerance=1e-7)
#define AssertThrow(cond, exc)
numbers::NumberTraits< Number >::real_type norm() const
virtual Point< spacedim > push_forward(const Point< 2 > &chart_point) const override
static ::ExceptionBase & ExcMessage(std::string arg1)
NormalToMeshProjectionBoundary(const TopoDS_Shape &sh, const double tolerance=1e-7)
#define Assert(cond, exc)
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
virtual std::unique_ptr< Manifold< dim, spacedim > > clone() const override
Point< dim > closest_point(const TopoDS_Shape &in_shape, const Point< dim > &origin, const double tolerance=1e-7)
virtual Point< spacedim > project_to_manifold(const ArrayView< const Point< spacedim >> &surrounding_points, const Point< spacedim > &candidate) const override
numbers::NumberTraits< Number >::real_type distance(const Point< dim, Number > &p) const
static ::ExceptionBase & ExcNotImplemented()
ArclengthProjectionLineManifold(const TopoDS_Shape &sh, const double tolerance=1e-7)
virtual Point< 1 > pull_back(const Point< spacedim > &space_point) const override
static ::ExceptionBase & ExcUnsupportedShape()
std::tuple< double, double, double, double > get_uv_bounds() const
virtual DerivativeForm< 1, 2, spacedim > push_forward_gradient(const Point< 2 > &chart_point) const override
NormalProjectionBoundary(const TopoDS_Shape &sh, const double tolerance=1e-7)
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)
static ::ExceptionBase & ExcInternalError()
DirectionalProjectionBoundary(const TopoDS_Shape &sh, const Tensor< 1, spacedim > &direction, const double tolerance=1e-7)