Reference documentation for deal.II version 9.0.0
boundary_lib.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2014 - 2018 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 #include <deal.II/opencascade/boundary_lib.h>
18 
19 #ifdef DEAL_II_WITH_OPENCASCADE
20 
21 
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>
32 #include <TopoDS.hxx>
33 
34 #include <Standard_Version.hxx>
35 #if (OCC_VERSION_MAJOR < 7)
36 # include <Handle_Adaptor3d_HCurve.hxx>
37 #endif
38 
39 
40 DEAL_II_NAMESPACE_OPEN
41 
42 
43 namespace OpenCASCADE
44 {
45 
46 
47  namespace
48  {
54  Handle_Adaptor3d_HCurve curve_adaptor(const TopoDS_Shape &shape)
55  {
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)));
63 
64  Assert(false, ExcInternalError());
65  return Handle(BRepAdaptor_HCurve)(new BRepAdaptor_HCurve());
66  }
67 
68 
69 
70 // Helper internal functions.
71  double shape_length(const TopoDS_Shape &sh)
72  {
73  Handle_Adaptor3d_HCurve adapt = curve_adaptor(sh);
74  return GCPnts_AbscissaPoint::Length(adapt->GetCurve());
75  }
76  }
77 
78  /*============================== NormalProjectionBoundary ==============================*/
79  template <int dim, int spacedim>
81  const double tolerance) :
82  sh(sh),
83  tolerance(tolerance)
84  {
85  Assert(spacedim == 3, ExcNotImplemented());
86  }
87 
88 
89 
90  template<int dim, int spacedim>
91  std::unique_ptr<Manifold<dim, spacedim> >
93  {
94  return std::unique_ptr<Manifold<dim,spacedim> >(new NormalProjectionBoundary(sh, tolerance));
95  }
96 
97 
98 
99  template <int dim, int spacedim>
101  project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
102  const Point<spacedim> &candidate) const
103  {
104  (void)surrounding_points;
105 #ifdef DEBUG
106  for (unsigned int i=0; i<surrounding_points.size(); ++i)
107  Assert(closest_point(sh, surrounding_points[i], tolerance)
108  .distance(surrounding_points[i]) <
109  std::max(tolerance*surrounding_points[i].norm(), tolerance),
110  ExcPointNotOnManifold<spacedim>(surrounding_points[i]));
111 #endif
112  return closest_point(sh, candidate,tolerance);
113  }
114 
115 
116  /*============================== DirectionalProjectionBoundary ==============================*/
117  template <int dim, int spacedim>
119  const Tensor<1,spacedim> &direction,
120  const double tolerance) :
121  sh(sh),
122  direction(direction),
123  tolerance(tolerance)
124  {
125  Assert(spacedim == 3, ExcNotImplemented());
126  }
127 
128 
129 
130  template<int dim, int spacedim>
131  std::unique_ptr<Manifold<dim, spacedim> >
133  {
134  return std::unique_ptr<Manifold<dim,spacedim> >
135  (new DirectionalProjectionBoundary(sh, direction, tolerance));
136  }
137 
138 
139 
140  template <int dim, int spacedim>
142  project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
143  const Point<spacedim> &candidate) const
144  {
145  (void)surrounding_points;
146 #ifdef DEBUG
147  for (unsigned int i=0; i<surrounding_points.size(); ++i)
148  Assert(closest_point(sh, surrounding_points[i],tolerance)
149  .distance(surrounding_points[i]) <
150  std::max(tolerance*surrounding_points[i].norm(), tolerance),
151  ExcPointNotOnManifold<spacedim>(surrounding_points[i]));
152 #endif
153  return line_intersection(sh, candidate, direction, tolerance);
154  }
155 
156 
157 
158  /*============================== NormalToMeshProjectionBoundary ==============================*/
159  template <int dim, int spacedim>
161  const double tolerance) :
162  sh(sh),
163  tolerance(tolerance)
164  {
165  Assert(spacedim == 3, ExcNotImplemented());
166  Assert(std::get<0>(count_elements(sh)) > 0,
167  ExcMessage("NormalToMeshProjectionBoundary needs a shape containing faces to operate."));
168  }
169 
170  template<int dim, int spacedim>
171  std::unique_ptr<Manifold<dim, spacedim> >
173  {
174  return std::unique_ptr<Manifold<dim, spacedim> >
176  }
177 
178 
179  template <int dim, int spacedim>
181  project_to_manifold (const ArrayView<const Point<spacedim>> &surrounding_points,
182  const Point<spacedim> &candidate) const
183  {
184  TopoDS_Shape out_shape;
185  Tensor<1,3> average_normal;
186 #ifdef DEBUG
187  for (unsigned int i=0; i<surrounding_points.size(); ++i)
188  {
189  Assert(closest_point(sh, surrounding_points[i], tolerance)
190  .distance(surrounding_points[i]) <
191  std::max(tolerance*surrounding_points[i].norm(), tolerance),
192  ExcPointNotOnManifold<spacedim>(surrounding_points[i]));
193  }
194 #endif
195 
196  switch (surrounding_points.size())
197  {
198  case 2:
199  {
200  for (unsigned int i=0; i<surrounding_points.size(); ++i)
201  {
202  std::tuple<Point<3>, Tensor<1,3>, double, double>
203  p_and_diff_forms =
205  surrounding_points[i],
206  tolerance);
207  average_normal += std::get<1>(p_and_diff_forms);
208  }
209 
210  average_normal/=2.0;
211 
212  Assert(average_normal.norm() > 1e-4,
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."));
214 
215  Tensor<1,3> T = surrounding_points[0]-surrounding_points[1];
216  T /= T.norm();
217  average_normal = average_normal-(average_normal*T)*T;
218  average_normal /= average_normal.norm();
219  break;
220  }
221  case 4:
222  {
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]};
226  Tensor<1,3> n1(n1_coords);
227  n1 = n1/n1.norm();
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]};
231  Tensor<1,3> n2(n2_coords);
232  n2 = n2/n2.norm();
233 
234  average_normal = (n1+n2)/2.0;
235 
236  Assert(average_normal.norm() > tolerance,
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."));
238 
239  average_normal /= average_normal.norm();
240  break;
241  }
242  case 8:
243  {
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]};
247  Tensor<1,3> n1(n1_coords);
248  n1 = n1/n1.norm();
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]};
252  Tensor<1,3> n2(n2_coords);
253  n2 = n2/n2.norm();
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]};
257  Tensor<1,3> n3(n3_coords);
258  n3 = n3/n3.norm();
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]};
262  Tensor<1,3> n4(n4_coords);
263  n4 = n4/n4.norm();
264 
265  average_normal = (n1+n2+n3+n4)/4.0;
266 
267  Assert(average_normal.norm() > tolerance,
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."));
269 
270  average_normal /= average_normal.norm();
271  break;
272  }
273  default:
274  {
275  AssertThrow(false, ExcNotImplemented());
276  break;
277  }
278  }
279 
280  return line_intersection(sh, candidate, average_normal, tolerance);
281  }
282 
283 
284  /*============================== ArclengthProjectionLineManifold ==============================*/
285  template <int dim, int spacedim>
287  const double tolerance):
288 
289  ChartManifold<dim,spacedim,1>(sh.Closed() ?
290  Point<1>(shape_length(sh)) :
291  Point<1>()),
292  sh(sh),
293  curve(curve_adaptor(sh)),
294  tolerance(tolerance),
295  length(shape_length(sh))
296  {
297  Assert(spacedim >= 2, ExcImpossibleInDimSpacedim(dim, spacedim));
298  }
299 
300 
301 
302  template<int dim, int spacedim>
303  std::unique_ptr<Manifold<dim, spacedim> >
305  {
306  return std::unique_ptr<Manifold<dim, spacedim> >
307  (new ArclengthProjectionLineManifold(sh,tolerance));
308  }
309 
310 
311 
312  template <int dim, int spacedim>
313  Point<1>
315  {
316  double t (0.0);
317  ShapeAnalysis_Curve curve_analysis;
318  gp_Pnt proj;
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));
321  (void)dist; // Silence compiler warning in Release mode.
322  return Point<1>(GCPnts_AbscissaPoint::Length(curve->GetCurve(),curve->GetCurve().FirstParameter(),t));
323  }
324 
325 
326 
327  template <int dim, int spacedim>
330  {
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);
334  }
335 
336  template <int dim, int spacedim>
338  NURBSPatchManifold( const TopoDS_Face &face,
339  const double tolerance)
340  :
341  face(face),
342  tolerance(tolerance)
343  {}
344 
345 
346 
347  template<int dim, int spacedim>
348  std::unique_ptr<Manifold<dim, spacedim> >
350  {
351  return std::unique_ptr<Manifold<dim, spacedim> >
352  (new NURBSPatchManifold<dim,spacedim>(face,tolerance));
353  }
354 
355 
356 
357  template <int dim, int spacedim> Point<2>
359  pull_back(const Point<spacedim> &space_point) const
360  {
361  Handle(Geom_Surface) SurfToProj = BRep_Tool::Surface(face);
362 
363  ShapeAnalysis_Surface projector(SurfToProj);
364  gp_Pnt2d proj_params = projector.ValueOfUV(point(space_point), tolerance);
365 
366  double u = proj_params.X();
367  double v = proj_params.Y();
368 
369  return Point<2>(u,v);
370  }
371 
372  template <int dim, int spacedim>
375  push_forward(const Point<2> &chart_point) const
376  {
377  return ::OpenCASCADE::push_forward<spacedim>(face, chart_point[0], chart_point[1]);
378  }
379 
380  template <int dim, int spacedim>
383  push_forward_gradient(const Point<2> &chart_point) const
384  {
386  Handle(Geom_Surface) surf = BRep_Tool::Surface(face);
387 
388  gp_Pnt q;
389  gp_Vec Du, Dv;
390  surf->D1(chart_point[0],chart_point[1], q, Du, Dv);
391 
392  DX[0][0] = Du.X();
393  DX[1][0] = Du.Y();
394  if (spacedim>2)
395  DX[2][0] = Du.Z();
396  else
397  Assert(std::abs(Du.Z()) < tolerance,
398  ExcMessage("Expecting derivative along Z to be zero! Bailing out."));
399  DX[0][1] = Dv.X();
400  DX[1][1] = Dv.Y();
401  if (spacedim>2)
402  DX[2][1] = Dv.Z();
403  else
404  Assert(std::abs(Dv.Z()) < tolerance,
405  ExcMessage("Expecting derivative along Z to be zero! Bailing out."));
406  return DX;
407  }
408 
409  template <int dim, int spacedim>
410  std::tuple<double, double, double, double>
413  {
414  Standard_Real umin, umax, vmin, vmax;
415  BRepTools::UVBounds(face, umin, umax, vmin, vmax);
416  return std::make_tuple(umin, umax, vmin, vmax);
417  }
418 
419 // Explicit instantiations
420 #include "boundary_lib.inst"
421 
422 } // end namespace OpenCASCADE
423 
424 DEAL_II_NAMESPACE_CLOSE
425 
426 #endif
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)
Definition: utilities.cc:89
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
Definition: boundary_lib.cc:92
Point< spacedim > point(const gp_Pnt &p, const double &tolerance=1e-10)
Definition: utilities.cc:183
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)
Definition: utilities.cc:418
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
numbers::NumberTraits< Number >::real_type norm() const
Definition: tensor.h:1273
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)
Definition: exceptions.h:1142
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)
Definition: utilities.cc:676
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)
Definition: boundary_lib.cc:80
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)
Definition: utilities.cc:686
static ::ExceptionBase & ExcInternalError()
DirectionalProjectionBoundary(const TopoDS_Shape &sh, const Tensor< 1, spacedim > &direction, const double tolerance=1e-7)