Reference documentation for deal.II version 8.5.1
boundary_lib.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2014 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii__occ_boundary_lib_h
18 #define dealii__occ_boundary_lib_h
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_OPENCASCADE
23 
24 #include <deal.II/opencascade/utilities.h>
25 #include <deal.II/grid/tria_boundary.h>
26 #include <deal.II/grid/manifold.h>
27 
28 // opencascade needs "HAVE_CONFIG_H" to be exported...
29 #define HAVE_CONFIG_H
30 #include <Adaptor3d_Curve.hxx>
31 #include <Adaptor3d_HCurve.hxx>
32 #include <BRepAdaptor_Curve.hxx>
33 #undef HAVE_CONFIG_H
34 
35 DEAL_II_NAMESPACE_OPEN
36 
42 namespace OpenCASCADE
43 {
65  template <int dim, int spacedim>
66  class NormalProjectionBoundary : public Boundary<dim,spacedim>
67  {
68  public:
69 
77  NormalProjectionBoundary(const TopoDS_Shape &sh,
78  const double tolerance=1e-7);
79 
89  virtual Point<spacedim>
90  project_to_manifold (const std::vector<Point<spacedim> > &surrounding_points,
91  const Point<spacedim> &candidate) const;
92 
93 
94  private:
101  const TopoDS_Shape sh;
102 
106  const double tolerance;
107  };
108 
130  template <int dim, int spacedim>
131  class DirectionalProjectionBoundary : public Boundary<dim,spacedim>
132  {
133  public:
138  DirectionalProjectionBoundary(const TopoDS_Shape &sh,
140  const double tolerance=1e-7);
141 
151  virtual Point<spacedim>
152  project_to_manifold (const std::vector<Point<spacedim> > &surrounding_points,
153  const Point<spacedim> &candidate) const;
154 
155  private:
162  const TopoDS_Shape sh;
163 
168 
172  const double tolerance;
173  };
174 
175 
220  template <int dim, int spacedim>
221  class NormalToMeshProjectionBoundary : public Boundary<dim,spacedim>
222  {
223  public:
229  NormalToMeshProjectionBoundary(const TopoDS_Shape &sh,
230  const double tolerance=1e-7);
231 
238  virtual Point<spacedim>
239  project_to_manifold (const std::vector<Point<spacedim> > &surrounding_points,
240  const Point<spacedim> &candidate) const;
241 
242  private:
249  const TopoDS_Shape sh;
250 
255 
259  const double tolerance;
260  };
261 
281  template <int dim, int spacedim>
282  class ArclengthProjectionLineManifold : public ChartManifold<dim,spacedim,1>
283  {
284  public:
288  ArclengthProjectionLineManifold(const TopoDS_Shape &sh,
289  const double tolerance=1e-7);
290 
296  virtual Point<1>
297  pull_back(const Point<spacedim> &space_point) const;
298 
302  virtual Point<spacedim>
303  push_forward(const Point<1> &chart_point) const;
304 
305  private:
310  Handle_Adaptor3d_HCurve curve;
311 
315  const double tolerance;
316 
321  const double length;
322  };
323 
331  template <int dim, int spacedim>
332  class NURBSPatchManifold : public ChartManifold<dim, spacedim, 2>
333  {
334  public:
340  NURBSPatchManifold(const TopoDS_Face &face, const double tolerance = 1e-7);
341 
346  virtual Point<2>
347  pull_back(const Point<spacedim> &space_point) const;
348 
353  virtual Point<spacedim>
354  push_forward(const Point<2> &chart_point) const;
355 
368  virtual
370  push_forward_gradient(const Point<2> &chart_point) const;
371 
372  private:
377  std_cxx11::tuple<double, double, double, double>
378  get_uv_bounds() const;
379 
383  TopoDS_Face face;
384 
389  double tolerance;
390  };
391 
392 }
393 
396 DEAL_II_NAMESPACE_CLOSE
397 
398 
399 #endif // DEAL_II_WITH_OPENCASCADE
400 
401 /*------------------------------ occ_boundary_lib.h ------------------------------*/
402 #endif
403 /*------------------------------ occ_boundary_lib.h ------------------------------*/
virtual Point< 2 > pull_back(const Point< spacedim > &space_point) const
virtual Point< spacedim > push_forward(const Point< 2 > &chart_point) const
virtual DerivativeForm< 1, 2, spacedim > push_forward_gradient(const Point< 2 > &chart_point) const
Definition: tria.h:52
NormalToMeshProjectionBoundary(const TopoDS_Shape &sh, const double tolerance=1e-7)
virtual Point< spacedim > project_to_manifold(const std::vector< Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const
virtual Point< spacedim > project_to_manifold(const std::vector< Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const
std_cxx11::tuple< double, double, double, double > get_uv_bounds() const
virtual Point< spacedim > push_forward(const Point< 1 > &chart_point) const
ArclengthProjectionLineManifold(const TopoDS_Shape &sh, const double tolerance=1e-7)
virtual Point< spacedim > project_to_manifold(const std::vector< Point< spacedim > > &surrounding_points, const Point< spacedim > &candidate) const
Definition: boundary_lib.cc:93
NormalProjectionBoundary(const TopoDS_Shape &sh, const double tolerance=1e-7)
Definition: boundary_lib.cc:82
virtual Point< 1 > pull_back(const Point< spacedim > &space_point) const
DirectionalProjectionBoundary(const TopoDS_Shape &sh, const Tensor< 1, spacedim > &direction, const double tolerance=1e-7)