gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
OCCFace.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #include "GmshConfig.h"
7 #include "GmshMessage.h"
8 #include "GModel.h"
9 #include "GModelIO_OCC.h"
10 #include "GEdgeLoop.h"
11 #include "OCCEdge.h"
12 #include "OCCFace.h"
13 #include "Numeric.h"
14 #include "Context.h"
15 #include "robustPredicates.h"
16 
17 #if defined(HAVE_OCC)
18 
19 #include <BOPTools_AlgoTools.hxx>
20 #include <BOPTools_AlgoTools2D.hxx>
21 #include <BRepBndLib.hxx>
22 #include <BRepClass_FaceClassifier.hxx>
23 #include <BRepLProp_SLProps.hxx>
24 #include <BRepTools.hxx>
25 #include <BRep_Builder.hxx>
26 #include <Bnd_Box.hxx>
27 #include <Geom_BSplineSurface.hxx>
28 #include <Geom_BezierSurface.hxx>
29 #include <Geom_ConicalSurface.hxx>
30 #include <Geom_CylindricalSurface.hxx>
31 #include <Geom_Plane.hxx>
32 #include <Geom_SphericalSurface.hxx>
33 #include <Geom_SurfaceOfRevolution.hxx>
34 #include <Geom_ToroidalSurface.hxx>
35 #include <IntTools_Context.hxx>
36 #include <ShapeAnalysis.hxx>
37 #include <ShapeFix_Wire.hxx>
38 #include <Standard_Version.hxx>
39 #include <TopExp_Explorer.hxx>
40 #include <TopoDS.hxx>
41 #include <gp_Pln.hxx>
42 #include <gp_Sphere.hxx>
43 
44 OCCFace::OCCFace(GModel *m, TopoDS_Face s, int num)
45  : GFace(m, num), _s(s), _sf(s, Standard_True), _radius(-1)
46 {
47  _setup();
48 
49  if(CTX::instance()->debugSurface > 0 &&
50  tag() == CTX::instance()->debugSurface)
51  writeBREP("debugSurface.brep");
52 }
53 
54 void OCCFace::_setup()
55 {
56  edgeLoops.clear();
57  l_edges.clear();
58  l_dirs.clear();
59 
60  TopExp_Explorer exp2, exp3;
61  for(exp2.Init(_s.Oriented(TopAbs_FORWARD), TopAbs_WIRE); exp2.More();
62  exp2.Next()) {
63  TopoDS_Wire wire = TopoDS::Wire(exp2.Current());
64 
65  // it is crucial that the edges in the wire are ordered and oriented
66  // correctly (for non-periodic surfaces GEdgeLoop would correct it; but for
67  // periodic surfaces the location of degenerate edges linking 2 sides of the
68  // parametric space is crucial) - so always make sure to reorder the edges
69  ShapeFix_Wire sfw(wire, _s, CTX::instance()->geom.tolerance);
70  sfw.FixReorder();
71  wire = sfw.Wire();
72 
73  Msg::Debug("OCC surface %d - new wire", tag());
74  GEdgeLoop el;
75 
76  for(exp3.Init(wire, TopAbs_EDGE); exp3.More(); exp3.Next()) {
77  TopoDS_Edge edge = TopoDS::Edge(exp3.Current());
78  GEdge *e = nullptr;
79  if(model()->getOCCInternals())
80  e = model()->getOCCInternals()->getEdgeForOCCShape(model(), edge);
81  if(!e) { Msg::Error("Unknown curve in surface %d", tag()); }
82  else if(edge.Orientation() == TopAbs_INTERNAL &&
84  Msg::Debug("Adding embedded curve %d in surface %d", e->tag(), tag());
85  embedded_edges.push_back(e);
86  /*
87  if(e->getBeginVertex()) embedded_vertices.insert(e->getBeginVertex());
88  if(e->getEndVertex()) embedded_vertices.insert(e->getEndVertex());
89  */
90  }
91  else {
92  int ori = edge.Orientation() ? -1 : 1;
93  Msg::Debug("Curve %d (%d --> %d) ori %d", e->tag(),
94  e->getBeginVertex() ? e->getBeginVertex()->tag() : -1,
95  e->getEndVertex() ? e->getEndVertex()->tag() : -1, ori);
96  el.add(ori, e);
97  }
98  }
99 
100  if(!el.check()) {
101  // should not happen, since ShapeFix_Wire has been called before
102  Msg::Warning("Recomputing incorrect OpenCASCADE wire in surface %d", tag());
103  std::vector<GEdge*> edges;
104  el.getEdges(edges);
105  el.recompute(edges);
106  }
107 
108  for(GEdgeLoop::citer it = el.begin(); it != el.end(); ++it) {
109  l_edges.push_back(it->getEdge());
110  l_dirs.push_back(it->getSign());
111  if(el.count() == 2) {
112  it->getEdge()->meshAttributes.minimumMeshSegments =
113  std::max(it->getEdge()->meshAttributes.minimumMeshSegments, 2);
114  }
115  if(el.count() == 1) {
116  it->getEdge()->meshAttributes.minimumMeshSegments =
117  std::max(it->getEdge()->meshAttributes.minimumMeshSegments, 3);
118  }
119  }
120  edgeLoops.push_back(el);
121  }
122 
123  for(exp2.Init(_s.Oriented(TopAbs_FORWARD), TopAbs_VERTEX, TopAbs_EDGE);
124  exp2.More(); exp2.Next()) {
125  TopoDS_Vertex vertex = TopoDS::Vertex(exp2.Current());
126  GVertex *v = nullptr;
127  if(model()->getOCCInternals())
128  v = model()->getOCCInternals()->getVertexForOCCShape(model(), vertex);
129  if(!v) { Msg::Error("Unknown point in surface %d", tag()); }
130  else if(vertex.Orientation() == TopAbs_INTERNAL &&
132  Msg::Debug("Adding embedded point %d in surface %d", v->tag(), tag());
133  embedded_vertices.insert(v);
134  }
135  }
136 
137  BRepAdaptor_Surface surface(_s);
138  _periodic[0] = surface.IsUPeriodic();
139  _periodic[1] = surface.IsVPeriodic();
140  if(_periodic[0]) _period[0] = surface.UPeriod();
141  if(_periodic[1]) _period[1] = surface.VPeriod();
142 
143  ShapeAnalysis::GetFaceUVBounds(_s, _umin, _umax, _vmin, _vmax);
144  Msg::Debug("OCC surface %d with %d parameter bounds (%g,%g)(%g,%g)", tag(),
145  l_edges.size(), _umin, _umax, _vmin, _vmax);
146 
147  _occface = BRep_Tool::Surface(_s);
148 
149  // init projector, with little tolerance to converge on the borders of the
150  // surface
151  double umin = _umin;
152  double vmin = _vmin;
153  double umax = _umax;
154  double vmax = _vmax;
155  if(!_periodic[0]) {
156  const double du = _umax - _umin;
157  const double utol = std::max(fabs(du) * 1e-8, 1e-12);
158  umin -= utol;
159  umax += utol;
160  }
161  if(!_periodic[1]) {
162  const double dv = _vmax - _vmin;
163  const double vtol = std::max(fabs(dv) * 1e-8, 1e-12);
164  vmin -= vtol;
165  vmax += vtol;
166  }
167  _projector.Init(_occface, umin, umax, vmin, vmax);
168 
169  if(OCCFace::geomType() == GEntity::Sphere) {
170  BRepAdaptor_Surface surface(_s);
171  gp_Sphere sphere = surface.Sphere();
172  _radius = sphere.Radius();
173  gp_Pnt loc = sphere.Location();
174  _center = SPoint3(loc.X(), loc.Y(), loc.Z());
175  }
176 
177  // Only store references to this new face in edges at the end of the
178  // constructor, to avoid accessing it too early (e.g. when drawing an edge)
179  for(std::size_t i = 0; i < l_edges.size(); i++) {
180  GEdge *e = l_edges[i];
181  e->addFace(this);
182  OCCEdge *occe = dynamic_cast<OCCEdge *>(e);
183  if(occe && !e->is3D()) occe->setTrimmed(this);
184  }
185  for(std::size_t i = 0; i < embedded_edges.size(); i++) {
186  GEdge *e = embedded_edges[i];
187  // should not addFace(), as the edge is not on the boundary
188  //e->addFace(this);
189  OCCEdge *occe = dynamic_cast<OCCEdge *>(e);
190  if(occe && !e->is3D()) occe->setTrimmed(this);
191  }
192 }
193 
194 SBoundingBox3d OCCFace::bounds(bool fast)
195 {
196  if(CTX::instance()->geom.occBoundsUseSTL) {
197  buildSTLTriangulation();
198  // BRepBndLib can use the STL mesh if available, but unfortunately it
199  // enlarges the box with the mesh deflection tolerance and the shape
200  // tolerance, which makes it hard to get the expected minimal box in simple
201  // cases (e.g. for plane surfaces), and always leads to boxes that are too
202  // large; so we simply compute the box from the STL vertices. The downside
203  // of this approach is that the bbox might be *smaller* than the actual box
204  // for curved shapes, but this is preferable for us as boxes are mostly used
205  // to find/identify entities
206  if(stl_vertices_xyz.size()) {
207  SBoundingBox3d bbox;
208  for(std::size_t i = 0; i < stl_vertices_xyz.size(); i++)
209  bbox += stl_vertices_xyz[i];
210  return bbox;
211  }
212  }
213 
214  Bnd_Box b;
215  try {
216  BRepBndLib::Add(_s, b);
217  } catch(Standard_Failure &err) {
218  Msg::Error("OpenCASCADE exception %s", err.GetMessageString());
219  return SBoundingBox3d();
220  }
221  double xmin, ymin, zmin, xmax, ymax, zmax;
222  b.Get(xmin, ymin, zmin, xmax, ymax, zmax);
223  SBoundingBox3d bbox(xmin, ymin, zmin, xmax, ymax, zmax);
224  return bbox;
225 }
226 
227 Range<double> OCCFace::parBounds(int i) const
228 {
229  if(i == 0) return Range<double>(_umin, _umax);
230  return Range<double>(_vmin, _vmax);
231 }
232 
233 SVector3 OCCFace::normal(const SPoint2 &param) const
234 {
235  gp_Pnt pnt;
236  gp_Vec du, dv;
237 
238  _occface->D1(param.x(), param.y(), pnt, du, dv);
239 
240  SVector3 t1(du.X(), du.Y(), du.Z());
241  SVector3 t2(dv.X(), dv.Y(), dv.Z());
242  SVector3 n(crossprod(t1, t2));
243  n.normalize();
244  if(_s.Orientation() == TopAbs_REVERSED) return n * (-1.);
245  return n;
246 }
247 
248 Pair<SVector3, SVector3> OCCFace::firstDer(const SPoint2 &param) const
249 {
250  gp_Pnt pnt;
251  gp_Vec du, dv;
252  _occface->D1(param.x(), param.y(), pnt, du, dv);
253  return Pair<SVector3, SVector3>(SVector3(du.X(), du.Y(), du.Z()),
254  SVector3(dv.X(), dv.Y(), dv.Z()));
255 }
256 
257 void OCCFace::secondDer(const SPoint2 &param, SVector3 &dudu, SVector3 &dvdv,
258  SVector3 &dudv) const
259 {
260  gp_Pnt pnt;
261  gp_Vec du, dv, duu, dvv, duv;
262  _occface->D2(param.x(), param.y(), pnt, du, dv, duu, dvv, duv);
263 
264  dudu = SVector3(duu.X(), duu.Y(), duu.Z());
265  dvdv = SVector3(dvv.X(), dvv.Y(), dvv.Z());
266  dudv = SVector3(duv.X(), duv.Y(), duv.Z());
267 }
268 
269 GPoint OCCFace::point(double par1, double par2) const
270 {
271  double pp[2] = {par1, par2};
272  gp_Pnt val = _occface->Value(par1, par2);
273  return GPoint(val.X(), val.Y(), val.Z(), this, pp);
274 }
275 
276 bool OCCFace::_project(const double p[3], double uv[2], double xyz[3]) const
277 {
278  gp_Pnt pnt(p[0], p[1], p[2]);
279  _projector.Perform(pnt);
280  if(!_projector.NbPoints()) {
281  Msg::Debug("Projection of point (%g, %g, %g) on surface %d failed", p[0],
282  p[1], p[2], tag());
283  return false;
284  }
285  _projector.LowerDistanceParameters(uv[0], uv[1]);
286 
287  if(uv[0] < _umin || uv[0] > _umax || uv[1] < _vmin || uv[1] > _vmax)
288  Msg::Debug("Point projection is out of surface parameter bounds");
289 
290  if(xyz) {
291  pnt = _projector.NearestPoint();
292  xyz[0] = pnt.X();
293  xyz[1] = pnt.Y();
294  xyz[2] = pnt.Z();
295  }
296  return true;
297 }
298 
300  const double initialGuess[2]) const
301 {
302 #if defined(HAVE_ALGLIB)
303  // less robust but can be much faster
304  if(CTX::instance()->geom.occUseGenericClosestPoint)
305  return GFace::closestPoint(qp, initialGuess);
306 #endif
307  double uv[2], xyz[3];
308  if(_project(qp.data(), uv, xyz))
309  return GPoint(xyz[0], xyz[1], xyz[2], this, uv);
310  else
311  return GFace::closestPoint(qp, initialGuess);
312 }
313 
314 SPoint2 OCCFace::parFromPoint(const SPoint3 &qp, bool onSurface,
315  bool convTestXYZ) const
316 {
317  // less robust but can be much faster
318  if(CTX::instance()->geom.occUseGenericClosestPoint)
319  return GFace::parFromPoint(qp, onSurface, convTestXYZ);
320  double uv[2];
321  if(_project(qp.data(), uv, nullptr))
322  return SPoint2(uv[0], uv[1]);
323  else // fallback: force convergence test in XYZ coordinates
324  return GFace::parFromPoint(qp, true, true);
325 }
326 
327 GEntity::GeomType OCCFace::geomType() const
328 {
329  if(_occface->DynamicType() == STANDARD_TYPE(Geom_Plane))
330  return Plane;
331  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_ToroidalSurface))
332  return Torus;
333  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_BezierSurface))
334  return BezierSurface;
335  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_CylindricalSurface))
336  return Cylinder;
337  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_ConicalSurface))
338  return Cone;
339  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_SurfaceOfRevolution))
340  return SurfaceOfRevolution;
341  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_SphericalSurface))
342  return Sphere;
343  else if(_occface->DynamicType() == STANDARD_TYPE(Geom_BSplineSurface))
344  return BSplineSurface;
345  return Unknown;
346 }
347 
348 double OCCFace::curvatureMax(const SPoint2 &param) const
349 {
350  const double eps = 1.e-12;
351  BRepLProp_SLProps prop(_sf, 2, eps);
352  prop.SetParameters(param.x(), param.y());
353 
354  if(!prop.IsCurvatureDefined()) { return eps; }
355 
356  return std::max(fabs(prop.MinCurvature()), fabs(prop.MaxCurvature()));
357 }
358 
359 double OCCFace::curvatures(const SPoint2 &param, SVector3 &dirMax,
360  SVector3 &dirMin, double &curvMax,
361  double &curvMin) const
362 {
363  const double eps = 1.e-12;
364  BRepLProp_SLProps prop(_sf, 2, eps);
365  prop.SetParameters(param.x(), param.y());
366 
367  if(!prop.IsCurvatureDefined()) { return -1.; }
368 
369  curvMax = prop.MaxCurvature();
370  curvMin = prop.MinCurvature();
371 
372  gp_Dir dMax = gp_Dir();
373  gp_Dir dMin = gp_Dir();
374  prop.CurvatureDirections(dMax, dMin);
375 
376  dirMax[0] = dMax.X();
377  dirMax[1] = dMax.Y();
378  dirMax[2] = dMax.Z();
379  dirMin[0] = dMin.X();
380  dirMin[1] = dMin.Y();
381  dirMin[2] = dMin.Z();
382 
383  return curvMax;
384 }
385 
386 bool OCCFace::containsPoint(const SPoint3 &pt) const
387 {
388  const Standard_Real tolerance = BRep_Tool::Tolerance(_s);
389  BRepClass_FaceClassifier faceClassifier;
390  faceClassifier.Perform(_s, gp_Pnt{pt.x(), pt.y(), pt.z()}, tolerance);
391  const TopAbs_State state = faceClassifier.State();
392  return (state == TopAbs_IN || state == TopAbs_ON);
393 }
394 
395 bool OCCFace::containsParam(const SPoint2 &pt)
396 {
397  const Standard_Real tolerance = BRep_Tool::Tolerance(_s);
398  BRepClass_FaceClassifier faceClassifier;
399  faceClassifier.Perform(_s, gp_Pnt2d{pt.x(), pt.y()}, tolerance);
400  const TopAbs_State state = faceClassifier.State();
401  return (state == TopAbs_IN || state == TopAbs_ON);
402 }
403 
404 bool OCCFace::buildSTLTriangulation(bool force)
405 {
406  if(stl_triangles.size() && !force) return true;
407  stl_vertices_uv.clear();
408  stl_vertices_xyz.clear();
409  stl_triangles.clear();
410  if(!model()->getOCCInternals()->makeFaceSTL(
411  _s, stl_vertices_uv, stl_vertices_xyz, stl_normals, stl_triangles)) {
412  Msg::Info("OpenCASCADE triangulation of surface %d failed", tag());
413  // try the default algorithm in GFace
414  return GFace::buildSTLTriangulation(force);
415  }
416 
417  // compute the triangulation of the edges which are the boundaries of this
418  // face
419  std::vector<GEdge *> const &e = edges();
420  for(auto it = e.begin(); it != e.end(); it++) {
421  if((*it)->stl_vertices_xyz.size() == 0) {
422  const TopoDS_Edge *c = (TopoDS_Edge *)(*it)->getNativePtr();
423  model()->getOCCInternals()->makeEdgeSTLFromFace(
424  *c, _s, &((*it)->stl_vertices_xyz));
425  }
426  }
427 
428  return true;
429 }
430 
431 bool OCCFace::isSphere(double &radius, SPoint3 &center) const
432 {
433  switch(geomType()) {
434  case GEntity::Sphere:
435  radius = _radius;
436  center = _center;
437  return true;
438  default: return false;
439  }
440 }
441 
442 void OCCFace::writeBREP(const char *filename)
443 {
444  BRep_Builder b;
445  TopoDS_Compound c;
446  b.MakeCompound(c);
447  b.Add(c, _s);
448  BRepTools::Write(c, filename);
449 }
450 
451 #endif
crossprod
SVector3 crossprod(const SVector3 &a, const SVector3 &b)
Definition: SVector3.h:150
GFace
Definition: GFace.h:33
Msg::Info
static void Info(const char *fmt,...)
Definition: GmshMessage.cpp:587
SPoint2
Definition: SPoint2.h:12
Msg::Debug
static void Debug(const char *fmt,...)
Definition: GmshMessage.cpp:752
robustPredicates.h
c
static double c(int i, int j, fullMatrix< double > &CA, const std::vector< SPoint3 > &P, const std::vector< SPoint3 > &Q)
Definition: discreteFrechetDistance.cpp:15
Msg::Warning
static void Warning(const char *fmt,...)
Definition: GmshMessage.cpp:543
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
OCCFace.h
SPoint3
Definition: SPoint3.h:14
GEdgeLoop::getEdges
void getEdges(std::vector< GEdge * > &edges) const
Definition: GEdgeLoop.cpp:100
SVector3
Definition: SVector3.h:16
GEdgeLoop::count
int count(GEdge *) const
Definition: GEdgeLoop.cpp:86
GEdgeLoop.h
GEdgeLoop::begin
iter begin()
Definition: GEdgeLoop.h:48
GModelIO_OCC.h
GEdgeLoop::citer
std::list< GEdgeSigned >::const_iterator citer
Definition: GEdgeLoop.h:42
GmshMessage.h
edges
static int edges[6][2]
Definition: meshGRegionLocalMeshMod.cpp:23
GPoint
Definition: GPoint.h:13
GFace::closestPoint
virtual GPoint closestPoint(const SPoint3 &queryPoint, const double initialGuess[2]) const
Definition: GFace.cpp:1323
debugSurface
int debugSurface
Definition: meshGFace.cpp:3150
GEdge::addFace
void addFace(GFace *f)
Definition: GEdge.cpp:200
GEdgeLoop::recompute
void recompute(const std::vector< GEdge * > &wire)
Definition: GEdgeLoop.cpp:137
CTX::instance
static CTX * instance()
Definition: Context.cpp:122
SPoint3::x
double x(void) const
Definition: SPoint3.h:125
GEdgeLoop
Definition: GEdgeLoop.h:36
OCCEdge.h
Range
Definition: Range.h:10
SPoint3::data
const double * data() const
Definition: SPoint3.h:76
GVertex
Definition: GVertex.h:23
tolerance
#define tolerance
Definition: curvature.cpp:11
Numeric.h
GModel
Definition: GModel.h:44
GEdge::getBeginVertex
virtual GVertex * getBeginVertex() const
Definition: GEdge.h:63
GFace::parFromPoint
virtual SPoint2 parFromPoint(const SPoint3 &, bool onSurface=true, bool convTestXYZ=false) const
Definition: GFace.cpp:1254
SPoint3::y
double y(void) const
Definition: SPoint3.h:127
GEdge::is3D
virtual bool is3D() const
Definition: GEdge.h:167
GEntity::tag
int tag() const
Definition: GEntity.h:280
CTX::geom
contextGeometryOptions geom
Definition: Context.h:311
GEntity::GeomType
GeomType
Definition: GEntity.h:88
SPoint2::x
double x(void) const
Definition: SPoint2.h:86
Context.h
Pair
Definition: Pair.h:10
GEdge
Definition: GEdge.h:26
SPoint3::z
double z(void) const
Definition: SPoint3.h:129
GFace::buildSTLTriangulation
virtual bool buildSTLTriangulation(bool force=false)
Definition: GFace.cpp:1585
contextGeometryOptions::occAutoEmbed
int occAutoEmbed
Definition: Context.h:100
GModel.h
GEdge::getEndVertex
virtual GVertex * getEndVertex() const
Definition: GEdge.h:64
GEdgeLoop::end
iter end()
Definition: GEdgeLoop.h:49
GEntity::Sphere
@ Sphere
Definition: GEntity.h:108
GEdgeLoop::check
bool check()
Definition: GEdgeLoop.cpp:176
SBoundingBox3d
Definition: SBoundingBox3d.h:21
closestPoint
double closestPoint(const std::vector< SPoint3 > &P, const SPoint3 &p, SPoint3 &result)
Definition: hausdorffDistance.cpp:50
GEdgeLoop::add
void add(int ori, GEdge *ge)
Definition: GEdgeLoop.h:46
SPoint2::y
double y(void) const
Definition: SPoint2.h:88