1 // Copyright (C) 2007-2014 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 // File : StdMeshers_Cartesian_3D.cxx
25 #include "StdMeshers_Cartesian_3D.hxx"
27 #include "SMDS_MeshNode.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Comment.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_CartesianParameters3D.hxx"
36 #include <utilities.h>
37 #include <Utils_ExceptHandlers.hxx>
38 #include <Basics_OCCTVersion.hxx>
40 #include <GEOMUtils.hxx>
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepBuilderAPI_Copy.hxx>
46 #include <BRepBuilderAPI_MakeFace.hxx>
47 #include <BRepTools.hxx>
48 #include <BRep_Builder.hxx>
49 #include <BRep_Tool.hxx>
50 #include <Bnd_B3d.hxx>
51 #include <Bnd_Box.hxx>
53 #include <GCPnts_UniformDeflection.hxx>
54 #include <Geom2d_BSplineCurve.hxx>
55 #include <Geom2d_BezierCurve.hxx>
56 #include <Geom2d_TrimmedCurve.hxx>
57 #include <GeomAPI_ProjectPointOnSurf.hxx>
58 #include <GeomLib.hxx>
59 #include <Geom_BSplineCurve.hxx>
60 #include <Geom_BSplineSurface.hxx>
61 #include <Geom_BezierCurve.hxx>
62 #include <Geom_BezierSurface.hxx>
63 #include <Geom_RectangularTrimmedSurface.hxx>
64 #include <Geom_TrimmedCurve.hxx>
65 #include <IntAna_IntConicQuad.hxx>
66 #include <IntAna_IntLinTorus.hxx>
67 #include <IntAna_Quadric.hxx>
68 #include <IntCurveSurface_TransitionOnCurve.hxx>
69 #include <IntCurvesFace_Intersector.hxx>
70 #include <Poly_Triangulation.hxx>
71 #include <Precision.hxx>
73 #include <TopExp_Explorer.hxx>
74 #include <TopLoc_Location.hxx>
75 #include <TopTools_MapOfShape.hxx>
77 #include <TopoDS_Compound.hxx>
78 #include <TopoDS_Face.hxx>
79 #include <TopoDS_TShape.hxx>
80 #include <gp_Cone.hxx>
81 #include <gp_Cylinder.hxx>
84 #include <gp_Pnt2d.hxx>
85 #include <gp_Sphere.hxx>
86 #include <gp_Torus.hxx>
92 #include <tbb/parallel_for.h>
93 //#include <tbb/enumerable_thread_specific.h>
102 #if OCC_VERSION_LARGE <= 0x06050300
103 // workaround is required only for OCCT6.5.3 and older (see OCC22809)
104 #define ELLIPSOLID_WORKAROUND
107 #ifdef ELLIPSOLID_WORKAROUND
108 #include <BRepIntCurveSurface_Inter.hxx>
109 #include <BRepTopAdaptor_TopolTool.hxx>
110 #include <BRepAdaptor_HSurface.hxx>
113 //=============================================================================
117 //=============================================================================
119 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
120 :SMESH_3D_Algo(hypId, studyId, gen)
122 _name = "Cartesian_3D";
123 _shapeType = (1 << TopAbs_SOLID); // 1 bit /shape type
124 _compatibleHypothesis.push_back("CartesianParameters3D");
126 _onlyUnaryInput = false; // to mesh all SOLIDs at once
127 _requireDiscreteBoundary = false; // 2D mesh not needed
128 _supportSubmeshes = false; // do not use any existing mesh
131 //=============================================================================
133 * Check presence of a hypothesis
135 //=============================================================================
137 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh& aMesh,
138 const TopoDS_Shape& aShape,
139 Hypothesis_Status& aStatus)
141 aStatus = SMESH_Hypothesis::HYP_MISSING;
143 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
144 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
145 if ( h == hyps.end())
150 for ( ; h != hyps.end(); ++h )
152 if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
154 aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
159 return aStatus == HYP_OK;
166 //=============================================================================
167 // Definitions of internal utils
168 // --------------------------------------------------------------------------
170 Trans_TANGENT = IntCurveSurface_Tangent,
171 Trans_IN = IntCurveSurface_In,
172 Trans_OUT = IntCurveSurface_Out,
175 // --------------------------------------------------------------------------
177 * \brief Common data of any intersection between a Grid and a shape
179 struct B_IntersectPoint
181 mutable const SMDS_MeshNode* _node;
182 mutable vector< TGeomID > _faceIDs;
184 B_IntersectPoint(): _node(NULL) {}
185 void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
186 int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
187 bool IsOnFace( int faceID ) const;
188 virtual ~B_IntersectPoint() {}
190 // --------------------------------------------------------------------------
192 * \brief Data of intersection between a GridLine and a TopoDS_Face
194 struct F_IntersectPoint : public B_IntersectPoint
197 mutable Transition _transition;
198 mutable size_t _indexOnLine;
200 bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
202 // --------------------------------------------------------------------------
204 * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
206 struct E_IntersectPoint : public B_IntersectPoint
212 // --------------------------------------------------------------------------
214 * \brief A line of the grid and its intersections with 2D geometry
219 double _length; // line length
220 multiset< F_IntersectPoint > _intPoints;
222 void RemoveExcessIntPoints( const double tol );
223 bool GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut );
225 // --------------------------------------------------------------------------
227 * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
231 gp_XYZ _uNorm, _vNorm, _zNorm;
232 vector< gp_XYZ > _origins; // origin points of all planes in one direction
233 vector< double > _zProjs; // projections of origins to _zNorm
235 // --------------------------------------------------------------------------
237 * \brief Iterator on the parallel grid lines of one direction
243 size_t _iVar1, _iVar2, _iConst;
244 string _name1, _name2, _nameConst;
246 LineIndexer( size_t sz1, size_t sz2, size_t sz3,
247 size_t iv1, size_t iv2, size_t iConst,
248 const string& nv1, const string& nv2, const string& nConst )
250 _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
251 _curInd[0] = _curInd[1] = _curInd[2] = 0;
252 _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
253 _name1 = nv1; _name2 = nv2; _nameConst = nConst;
256 size_t I() const { return _curInd[0]; }
257 size_t J() const { return _curInd[1]; }
258 size_t K() const { return _curInd[2]; }
259 void SetIJK( size_t i, size_t j, size_t k )
261 _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
265 if ( ++_curInd[_iVar1] == _size[_iVar1] )
266 _curInd[_iVar1] = 0, ++_curInd[_iVar2];
268 bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
269 size_t LineIndex () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
270 size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
271 size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
272 size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
273 void SetIndexOnLine (size_t i) { _curInd[ _iConst ] = i; }
274 size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
276 // --------------------------------------------------------------------------
278 * \brief Container of GridLine's
282 vector< double > _coords[3]; // coordinates of grid nodes
283 gp_XYZ _axes [3]; // axis directions
284 vector< GridLine > _lines [3]; // in 3 directions
285 double _tol, _minCellSize;
287 gp_Mat _invB; // inverted basis of _axes
288 //bool _isOrthogonalAxes;
290 vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
291 vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
293 list< E_IntersectPoint > _edgeIntP; // intersections with EDGEs
294 TopTools_IndexedMapOfShape _shapes;
296 SMESH_MesherHelper* _helper;
298 size_t CellIndex( size_t i, size_t j, size_t k ) const
300 return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
302 size_t NodeIndex( size_t i, size_t j, size_t k ) const
304 return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
306 size_t NodeIndexDX() const { return 1; }
307 size_t NodeIndexDY() const { return _coords[0].size(); }
308 size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
310 LineIndexer GetLineIndexer(size_t iDir) const;
312 void SetCoordinates(const vector<double>& xCoords,
313 const vector<double>& yCoords,
314 const vector<double>& zCoords,
315 const double* axesDirs,
316 const Bnd_Box& bndBox );
317 void ComputeUVW(const gp_XYZ& p, double uvw[3]);
318 void ComputeNodes(SMESH_MesherHelper& helper);
320 #ifdef ELLIPSOLID_WORKAROUND
321 // --------------------------------------------------------------------------
323 * \brief struct temporary replacing IntCurvesFace_Intersector until
324 * OCCT bug 0022809 is fixed
325 * http://tracker.dev.opencascade.org/view.php?id=22809
327 struct TMP_IntCurvesFace_Intersector
329 BRepAdaptor_Surface _surf;
331 BRepIntCurveSurface_Inter _intcs;
332 vector<IntCurveSurface_IntersectionPoint> _points;
333 BRepTopAdaptor_TopolTool _clsf;
335 TMP_IntCurvesFace_Intersector(const TopoDS_Face& face, const double tol)
336 :_surf( face ), _tol( tol ), _clsf( new BRepAdaptor_HSurface(_surf) ) {}
337 Bnd_Box Bounding() const { Bnd_Box b; BRepBndLib::Add (_surf.Face(), b); return b; }
338 void Perform( const gp_Lin& line, const double w0, const double w1 )
341 for ( _intcs.Init( _surf.Face(), line, _tol ); _intcs.More(); _intcs.Next() )
342 if ( w0 <= _intcs.W() && _intcs.W() <= w1 )
343 _points.push_back( _intcs.Point() );
345 bool IsDone() const { return true; }
346 int NbPnt() const { return _points.size(); }
347 IntCurveSurface_TransitionOnCurve Transition( const int i ) const { return _points[ i-1 ].Transition(); }
348 double WParameter( const int i ) const { return _points[ i-1 ].W(); }
349 TopAbs_State ClassifyUVPoint(const gp_Pnt2d& p) { return _clsf.Classify( p, _tol ); }
351 #define __IntCurvesFace_Intersector TMP_IntCurvesFace_Intersector
353 #define __IntCurvesFace_Intersector IntCurvesFace_Intersector
355 // --------------------------------------------------------------------------
357 * \brief Intersector of TopoDS_Face with all GridLine's
359 struct FaceGridIntersector
365 __IntCurvesFace_Intersector* _surfaceInt;
366 vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
368 FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
370 bool IsInGrid(const Bnd_Box& gridBox);
372 void StoreIntersections()
374 for ( size_t i = 0; i < _intersections.size(); ++i )
376 multiset< F_IntersectPoint >::iterator ip =
377 _intersections[i].first->_intPoints.insert( _intersections[i].second );
378 ip->_faceIDs.reserve( 1 );
379 ip->_faceIDs.push_back( _faceID );
382 const Bnd_Box& GetFaceBndBox()
384 GetCurveFaceIntersector();
387 __IntCurvesFace_Intersector* GetCurveFaceIntersector()
391 _surfaceInt = new __IntCurvesFace_Intersector( _face, Precision::PConfusion() );
392 _bndBox = _surfaceInt->Bounding();
393 if ( _bndBox.IsVoid() )
394 BRepBndLib::Add (_face, _bndBox);
398 bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
400 // --------------------------------------------------------------------------
402 * \brief Intersector of a surface with a GridLine
404 struct FaceLineIntersector
407 double _u, _v, _w; // params on the face and the line
408 Transition _transition; // transition of at intersection (see IntCurveSurface.cdl)
409 Transition _transIn, _transOut; // IN and OUT transitions depending of face orientation
412 gp_Cylinder _cylinder;
416 __IntCurvesFace_Intersector* _surfaceInt;
418 vector< F_IntersectPoint > _intPoints;
420 void IntersectWithPlane (const GridLine& gridLine);
421 void IntersectWithCylinder(const GridLine& gridLine);
422 void IntersectWithCone (const GridLine& gridLine);
423 void IntersectWithSphere (const GridLine& gridLine);
424 void IntersectWithTorus (const GridLine& gridLine);
425 void IntersectWithSurface (const GridLine& gridLine);
427 bool UVIsOnFace() const;
428 void addIntPoint(const bool toClassify=true);
429 bool isParamOnLineOK( const double linLength )
431 return -_tol < _w && _w < linLength + _tol;
433 FaceLineIntersector():_surfaceInt(0) {}
434 ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
436 // --------------------------------------------------------------------------
438 * \brief Class representing topology of the hexahedron and creating a mesh
439 * volume basing on analysis of hexahedron intersection with geometry
443 // --------------------------------------------------------------------------------
446 // --------------------------------------------------------------------------------
447 struct _Node //!< node either at a hexahedron corner or at intersection
449 const SMDS_MeshNode* _node; // mesh node at hexahedron corner
450 const B_IntersectPoint* _intPoint;
451 const _Face* _usedInFace;
453 _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
454 :_node(n), _intPoint(ip), _usedInFace(0) {}
455 const SMDS_MeshNode* Node() const
456 { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
457 //const F_IntersectPoint* FaceIntPnt() const
458 //{ return static_cast< const F_IntersectPoint* >( _intPoint ); }
459 const E_IntersectPoint* EdgeIntPnt() const
460 { return static_cast< const E_IntersectPoint* >( _intPoint ); }
461 bool IsUsedInFace( const _Face* polygon = 0 )
463 return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
465 void Add( const E_IntersectPoint* ip )
470 else if ( !_intPoint->_node ) {
471 ip->Add( _intPoint->_faceIDs );
475 _intPoint->Add( ip->_faceIDs );
478 int IsLinked( const B_IntersectPoint* other,
479 int avoidFace=-1 ) const // returns id of a common face
481 return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
483 bool IsOnFace( int faceID ) const // returns true if faceID is found
485 return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
489 if ( const SMDS_MeshNode* n = Node() )
490 return SMESH_TNodeXYZ( n );
491 if ( const E_IntersectPoint* eip =
492 dynamic_cast< const E_IntersectPoint* >( _intPoint ))
494 return gp_Pnt( 1e100, 0, 0 );
497 // --------------------------------------------------------------------------------
498 struct _Link // link connecting two _Node's
501 _Face* _faces[2]; // polygons sharing a link
502 vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
503 vector< _Node* > _fIntNodes; // _Node's at _fIntPoints
504 vector< _Link > _splits;
505 _Link() { _faces[0] = 0; }
507 // --------------------------------------------------------------------------------
512 _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
513 void Reverse() { _reverse = !_reverse; }
514 int NbResultLinks() const { return _link->_splits.size(); }
515 _OrientedLink ResultLink(int i) const
517 return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
519 _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
520 _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
521 operator bool() const { return _link; }
522 vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
524 vector< TGeomID > faces;
525 const B_IntersectPoint *ip0, *ip1;
526 if (( ip0 = _link->_nodes[0]->_intPoint ) &&
527 ( ip1 = _link->_nodes[1]->_intPoint ))
529 for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
530 if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
531 !usedIDs.count( ip0->_faceIDs[i] ) )
532 faces.push_back( ip0->_faceIDs[i] );
536 bool HasEdgeNodes() const
538 return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
539 dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
543 return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
545 void AddFace( _Face* f )
547 if ( _link->_faces[0] )
549 _link->_faces[1] = f;
553 _link->_faces[0] = f;
554 _link->_faces[1] = 0;
557 void RemoveFace( _Face* f )
559 if ( !_link->_faces[0] ) return;
561 if ( _link->_faces[1] == f )
563 _link->_faces[1] = 0;
565 else if ( _link->_faces[0] == f )
568 if ( _link->_faces[1] )
570 _link->_faces[0] = _link->_faces[1];
571 _link->_faces[1] = 0;
576 // --------------------------------------------------------------------------------
579 vector< _OrientedLink > _links; // links on GridLine's
580 vector< _Link > _polyLinks; // links added to close a polygonal face
581 vector< _Node* > _eIntNodes; // nodes at intersection with EDGEs
582 bool isPolyLink( const _OrientedLink& ol )
584 return _polyLinks.empty() ? false :
585 ( &_polyLinks[0] <= ol._link && ol._link <= &_polyLinks.back() );
588 // --------------------------------------------------------------------------------
589 struct _volumeDef // holder of nodes of a volume mesh element
591 vector< _Node* > _nodes;
592 vector< int > _quantities;
593 typedef boost::shared_ptr<_volumeDef> Ptr;
594 void set( const vector< _Node* >& nodes,
595 const vector< int >& quant = vector< int >() )
596 { _nodes = nodes; _quantities = quant; }
599 // topology of a hexahedron
602 _Link _hexLinks [12];
605 // faces resulted from hexahedron intersection
606 vector< _Face > _polygons;
608 // intresections with EDGEs
609 vector< const E_IntersectPoint* > _eIntPoints;
611 // additional nodes created at intersection points
612 vector< _Node > _intNodes;
614 // nodes inside the hexahedron (at VERTEXes)
615 vector< _Node* > _vIntNodes;
617 // computed volume elements
618 //vector< _volumeDef::Ptr > _volumeDefs;
619 _volumeDef _volumeDefs;
622 double _sizeThreshold, _sideLength[3];
623 int _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
624 int _origNodeInd; // index of _hexNodes[0] node within the _grid
628 Hexahedron(const double sizeThreshold, Grid* grid);
629 int MakeElements(SMESH_MesherHelper& helper,
630 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
631 void ComputeElements();
632 void Init() { init( _i, _j, _k ); }
635 Hexahedron(const Hexahedron& other );
636 void init( size_t i, size_t j, size_t k );
637 void init( size_t i );
638 void addEdges(SMESH_MesherHelper& helper,
639 vector< Hexahedron* >& intersectedHex,
640 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
641 gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
642 double proj, BRepAdaptor_Curve& curve,
643 const gp_XYZ& axis, const gp_XYZ& origin );
644 int getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
645 bool addIntersection( const E_IntersectPoint& ip,
646 vector< Hexahedron* >& hexes,
647 int ijk[], int dIJK[] );
648 bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
649 bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
650 int addElements(SMESH_MesherHelper& helper);
651 bool is1stNodeOut( _Link& link ) const;
652 bool isInHole() const;
653 bool checkPolyhedronSize() const;
658 bool debugDumpLink( _Link* link );
659 _Node* FindEqualNode( vector< _Node* >& nodes,
660 const E_IntersectPoint* ip,
663 for ( size_t i = 0; i < nodes.size(); ++i )
664 if ( nodes[i]->EdgeIntPnt() == ip ||
665 nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
672 // --------------------------------------------------------------------------
674 * \brief Hexahedron computing volumes in one thread
676 struct ParallelHexahedron
678 vector< Hexahedron* >& _hexVec;
680 ParallelHexahedron( vector< Hexahedron* >& hv, vector<int>& ind): _hexVec(hv), _index(ind) {}
681 void operator() ( const tbb::blocked_range<size_t>& r ) const
683 for ( size_t i = r.begin(); i != r.end(); ++i )
684 if ( Hexahedron* hex = _hexVec[ _index[i]] )
685 hex->ComputeElements();
688 // --------------------------------------------------------------------------
690 * \brief Structure intersecting certain nb of faces with GridLine's in one thread
692 struct ParallelIntersector
694 vector< FaceGridIntersector >& _faceVec;
695 ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
696 void operator() ( const tbb::blocked_range<size_t>& r ) const
698 for ( size_t i = r.begin(); i != r.end(); ++i )
699 _faceVec[i].Intersect();
704 //=============================================================================
705 // Implementation of internal utils
706 //=============================================================================
708 * \brief adjust \a i to have \a val between values[i] and values[i+1]
710 inline void locateValue( int & i, double val, const vector<double>& values,
711 int& di, double tol )
713 //val += values[0]; // input \a val is measured from 0.
714 if ( i > values.size()-2 )
717 while ( i+2 < values.size() && val > values[ i+1 ])
719 while ( i > 0 && val < values[ i ])
722 if ( i > 0 && val - values[ i ] < tol )
724 else if ( i+2 < values.size() && values[ i+1 ] - val < tol )
729 //=============================================================================
731 * Remove coincident intersection points
733 void GridLine::RemoveExcessIntPoints( const double tol )
735 if ( _intPoints.size() < 2 ) return;
737 set< Transition > tranSet;
738 multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
739 while ( ip2 != _intPoints.end() )
743 while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
745 tranSet.insert( ip1->_transition );
746 tranSet.insert( ip2->_transition );
747 ip2->Add( ip1->_faceIDs );
748 _intPoints.erase( ip1 );
751 if ( tranSet.size() > 1 ) // points with different transition coincide
753 bool isIN = tranSet.count( Trans_IN );
754 bool isOUT = tranSet.count( Trans_OUT );
756 (*ip1)._transition = Trans_TANGENT;
758 (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
762 //================================================================================
764 * Return "is OUT" state for nodes before the given intersection point
766 bool GridLine::GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut )
768 if ( ip->_transition == Trans_IN )
770 if ( ip->_transition == Trans_OUT )
772 if ( ip->_transition == Trans_APEX )
774 // singularity point (apex of a cone)
775 if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
777 multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
778 if ( ipAft == _intPoints.end() )
781 if ( ipBef->_transition != ipAft->_transition )
782 return ( ipBef->_transition == Trans_OUT );
783 return ( ipBef->_transition != Trans_OUT );
785 // _transition == Trans_TANGENT
788 //================================================================================
792 void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
793 const SMDS_MeshNode* n) const
795 if ( _faceIDs.empty() )
798 for ( size_t i = 0; i < fIDs.size(); ++i )
800 vector< TGeomID >::iterator it =
801 std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
802 if ( it == _faceIDs.end() )
803 _faceIDs.push_back( fIDs[i] );
808 //================================================================================
810 * Returns index of a common face if any, else zero
812 int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
815 for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
816 if ( avoidFace != other->_faceIDs[i] &&
817 IsOnFace ( other->_faceIDs[i] ))
818 return other->_faceIDs[i];
821 //================================================================================
823 * Returns \c true if \a faceID in in this->_faceIDs
825 bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
827 vector< TGeomID >::const_iterator it =
828 std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
829 return ( it != _faceIDs.end() );
831 //================================================================================
833 * Return an iterator on GridLine's in a given direction
835 LineIndexer Grid::GetLineIndexer(size_t iDir) const
837 const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
838 const string s [] = { "X", "Y", "Z" };
839 LineIndexer li( _coords[0].size(), _coords[1].size(), _coords[2].size(),
840 indices[iDir*3], indices[iDir*3+1], indices[iDir*3+2],
841 s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
844 //=============================================================================
846 * Creates GridLine's of the grid
848 void Grid::SetCoordinates(const vector<double>& xCoords,
849 const vector<double>& yCoords,
850 const vector<double>& zCoords,
851 const double* axesDirs,
852 const Bnd_Box& shapeBox)
854 _coords[0] = xCoords;
855 _coords[1] = yCoords;
856 _coords[2] = zCoords;
858 _axes[0].SetCoord( axesDirs[0],
861 _axes[1].SetCoord( axesDirs[3],
864 _axes[2].SetCoord( axesDirs[6],
867 _axes[0].Normalize();
868 _axes[1].Normalize();
869 _axes[2].Normalize();
871 _invB.SetCols( _axes[0], _axes[1], _axes[2] );
874 // _isOrthogonalAxes = ( Abs( _axes[0] * _axes[1] ) < 1e-20 &&
875 // Abs( _axes[1] * _axes[2] ) < 1e-20 &&
876 // Abs( _axes[2] * _axes[0] ) < 1e-20 );
879 _minCellSize = Precision::Infinite();
880 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
882 for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
884 double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
885 if ( cellLen < _minCellSize )
886 _minCellSize = cellLen;
889 if ( _minCellSize < Precision::Confusion() )
890 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
891 SMESH_Comment("Too small cell size: ") << _minCellSize );
892 _tol = _minCellSize / 1000.;
894 // attune grid extremities to shape bounding box
896 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
897 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
898 double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
899 &_coords[0].back(), &_coords[1].back(), &_coords[2].back() };
900 for ( int i = 0; i < 6; ++i )
901 if ( fabs( sP[i] - *cP[i] ) < _tol )
902 *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
904 for ( int iDir = 0; iDir < 3; ++iDir )
906 if ( _coords[iDir][0] - sP[iDir] > _tol )
908 _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
909 _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
911 if ( sP[iDir+3] - _coords[iDir].back() > _tol )
913 _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
914 _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
917 _tol = _minCellSize / 1000.;
919 _origin = ( _coords[0][0] * _axes[0] +
920 _coords[1][0] * _axes[1] +
921 _coords[2][0] * _axes[2] );
924 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
926 LineIndexer li = GetLineIndexer( iDir );
927 _lines[iDir].resize( li.NbLines() );
928 double len = _coords[ iDir ].back() - _coords[iDir].front();
929 for ( ; li.More(); ++li )
931 GridLine& gl = _lines[iDir][ li.LineIndex() ];
932 gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
933 _coords[1][li.J()] * _axes[1] +
934 _coords[2][li.K()] * _axes[2] );
935 gl._line.SetDirection( _axes[ iDir ]);
940 //================================================================================
942 * Computes coordinates of a point in the grid CS
944 void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
946 // gp_XYZ p = P - _origin;
947 // UVW[ 0 ] = p.X() * _invB( 1, 1 ) + p.Y() * _invB( 1, 2 ) + p.Z() * _invB( 1, 3 );
948 // UVW[ 1 ] = p.X() * _invB( 2, 1 ) + p.Y() * _invB( 2, 2 ) + p.Z() * _invB( 2, 3 );
949 // UVW[ 2 ] = p.X() * _invB( 3, 1 ) + p.Y() * _invB( 3, 2 ) + p.Z() * _invB( 3, 3 );
950 // UVW[ 0 ] += _coords[0][0];
951 // UVW[ 1 ] += _coords[1][0];
952 // UVW[ 2 ] += _coords[2][0];
953 gp_XYZ p = P * _invB;
954 p.Coord( UVW[0], UVW[1], UVW[2] );
956 //================================================================================
960 void Grid::ComputeNodes(SMESH_MesherHelper& helper)
962 // state of each node of the grid relative to the geometry
963 const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
964 vector< bool > isNodeOut( nbGridNodes, false );
965 _nodes.resize( nbGridNodes, 0 );
966 _gridIntP.resize( nbGridNodes, NULL );
968 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
970 LineIndexer li = GetLineIndexer( iDir );
972 // find out a shift of node index while walking along a GridLine in this direction
973 li.SetIndexOnLine( 0 );
974 size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
975 li.SetIndexOnLine( 1 );
976 const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
978 const vector<double> & coords = _coords[ iDir ];
979 for ( ; li.More(); ++li ) // loop on lines in iDir
981 li.SetIndexOnLine( 0 );
982 nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
984 GridLine& line = _lines[ iDir ][ li.LineIndex() ];
985 const gp_XYZ lineLoc = line._line.Location().XYZ();
986 const gp_XYZ lineDir = line._line.Direction().XYZ();
987 line.RemoveExcessIntPoints( _tol );
988 multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
989 multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
992 const double* nodeCoord = & coords[0];
993 const double* coord0 = nodeCoord;
994 const double* coordEnd = coord0 + coords.size();
995 double nodeParam = 0;
996 for ( ; ip != intPnts.end(); ++ip )
998 // set OUT state or just skip IN nodes before ip
999 if ( nodeParam < ip->_paramOnLine - _tol )
1001 isOut = line.GetIsOutBefore( ip, isOut );
1003 while ( nodeParam < ip->_paramOnLine - _tol )
1006 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
1007 if ( ++nodeCoord < coordEnd )
1008 nodeParam = *nodeCoord - *coord0;
1012 if ( nodeCoord == coordEnd ) break;
1014 // create a mesh node on a GridLine at ip if it does not coincide with a grid node
1015 if ( nodeParam > ip->_paramOnLine + _tol )
1017 // li.SetIndexOnLine( 0 );
1018 // double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1019 // xyz[ li._iConst ] += ip->_paramOnLine;
1020 gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
1021 ip->_node = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1022 ip->_indexOnLine = nodeCoord-coord0-1;
1024 // create a mesh node at ip concident with a grid node
1027 int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1028 if ( !_nodes[ nodeIndex ] )
1030 //li.SetIndexOnLine( nodeCoord-coord0 );
1031 //double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1032 gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1033 _nodes [ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1034 _gridIntP[ nodeIndex ] = & * ip;
1036 if ( _gridIntP[ nodeIndex ] )
1037 _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1039 _gridIntP[ nodeIndex ] = & * ip;
1040 // ip->_node = _nodes[ nodeIndex ]; -- to differ from ip on links
1041 ip->_indexOnLine = nodeCoord-coord0;
1042 if ( ++nodeCoord < coordEnd )
1043 nodeParam = *nodeCoord - *coord0;
1046 // set OUT state to nodes after the last ip
1047 for ( ; nodeCoord < coordEnd; ++nodeCoord )
1048 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
1052 // Create mesh nodes at !OUT nodes of the grid
1054 for ( size_t z = 0; z < _coords[2].size(); ++z )
1055 for ( size_t y = 0; y < _coords[1].size(); ++y )
1056 for ( size_t x = 0; x < _coords[0].size(); ++x )
1058 size_t nodeIndex = NodeIndex( x, y, z );
1059 if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
1061 //_nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
1062 gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1063 _coords[1][y] * _axes[1] +
1064 _coords[2][z] * _axes[2] );
1065 _nodes[ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1070 // check validity of transitions
1071 const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1072 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1074 LineIndexer li = GetLineIndexer( iDir );
1075 for ( ; li.More(); ++li )
1077 multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1078 if ( intPnts.empty() ) continue;
1079 if ( intPnts.size() == 1 )
1081 if ( intPnts.begin()->_transition != Trans_TANGENT &&
1082 intPnts.begin()->_transition != Trans_APEX )
1083 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1084 SMESH_Comment("Wrong SOLE transition of GridLine (")
1085 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1086 << ") along " << li._nameConst
1087 << ": " << trName[ intPnts.begin()->_transition] );
1091 if ( intPnts.begin()->_transition == Trans_OUT )
1092 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1093 SMESH_Comment("Wrong START transition of GridLine (")
1094 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1095 << ") along " << li._nameConst
1096 << ": " << trName[ intPnts.begin()->_transition ]);
1097 if ( intPnts.rbegin()->_transition == Trans_IN )
1098 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1099 SMESH_Comment("Wrong END transition of GridLine (")
1100 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1101 << ") along " << li._nameConst
1102 << ": " << trName[ intPnts.rbegin()->_transition ]);
1109 //=============================================================================
1111 * Checks if the face is encosed by the grid
1113 bool FaceGridIntersector::IsInGrid(const Bnd_Box& gridBox)
1115 // double x0,y0,z0, x1,y1,z1;
1116 // const Bnd_Box& faceBox = GetFaceBndBox();
1117 // faceBox.Get(x0,y0,z0, x1,y1,z1);
1119 // if ( !gridBox.IsOut( gp_Pnt( x0,y0,z0 )) &&
1120 // !gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
1123 // double X0,Y0,Z0, X1,Y1,Z1;
1124 // gridBox.Get(X0,Y0,Z0, X1,Y1,Z1);
1125 // double faceP[6] = { x0,y0,z0, x1,y1,z1 };
1126 // double gridP[6] = { X0,Y0,Z0, X1,Y1,Z1 };
1127 // gp_Dir axes[3] = { gp::DX(), gp::DY(), gp::DZ() };
1128 // for ( int iDir = 0; iDir < 6; ++iDir )
1130 // if ( iDir < 3 && gridP[ iDir ] <= faceP[ iDir ] ) continue;
1131 // if ( iDir >= 3 && gridP[ iDir ] >= faceP[ iDir ] ) continue;
1133 // // check if the face intersects a side of a gridBox
1135 // gp_Pnt p = iDir < 3 ? gp_Pnt( X0,Y0,Z0 ) : gp_Pnt( X1,Y1,Z1 );
1136 // gp_Ax1 norm( p, axes[ iDir % 3 ] );
1137 // if ( iDir < 3 ) norm.Reverse();
1139 // gp_XYZ O = norm.Location().XYZ(), N = norm.Direction().XYZ();
1141 // TopLoc_Location loc = _face.Location();
1142 // Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(_face,loc);
1143 // if ( !aPoly.IsNull() )
1145 // if ( !loc.IsIdentity() )
1147 // norm.Transform( loc.Transformation().Inverted() );
1148 // O = norm.Location().XYZ(), N = norm.Direction().XYZ();
1150 // const double deflection = aPoly->Deflection();
1152 // const TColgp_Array1OfPnt& nodes = aPoly->Nodes();
1153 // for ( int i = nodes.Lower(); i <= nodes.Upper(); ++i )
1154 // if (( nodes( i ).XYZ() - O ) * N > _grid->_tol + deflection )
1159 // BRepAdaptor_Surface surf( _face );
1160 // double u0, u1, v0, v1, du, dv, u, v;
1161 // BRepTools::UVBounds( _face, u0, u1, v0, v1);
1162 // if ( surf.GetType() == GeomAbs_Plane ) {
1163 // du = u1 - u0, dv = v1 - v0;
1166 // du = surf.UResolution( _grid->_minCellSize / 10. );
1167 // dv = surf.VResolution( _grid->_minCellSize / 10. );
1169 // for ( u = u0, v = v0; u <= u1 && v <= v1; u += du, v += dv )
1171 // gp_Pnt p = surf.Value( u, v );
1172 // if (( p.XYZ() - O ) * N > _grid->_tol )
1174 // TopAbs_State state = GetCurveFaceIntersector()->ClassifyUVPoint(gp_Pnt2d( u, v ));
1175 // if ( state == TopAbs_IN || state == TopAbs_ON )
1183 //=============================================================================
1185 * Intersects TopoDS_Face with all GridLine's
1187 void FaceGridIntersector::Intersect()
1189 FaceLineIntersector intersector;
1190 intersector._surfaceInt = GetCurveFaceIntersector();
1191 intersector._tol = _grid->_tol;
1192 intersector._transOut = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1193 intersector._transIn = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1195 typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1196 PIntFun interFunction;
1198 BRepAdaptor_Surface surf( _face );
1199 switch ( surf.GetType() ) {
1201 intersector._plane = surf.Plane();
1202 interFunction = &FaceLineIntersector::IntersectWithPlane;
1204 case GeomAbs_Cylinder:
1205 intersector._cylinder = surf.Cylinder();
1206 interFunction = &FaceLineIntersector::IntersectWithCylinder;
1209 intersector._cone = surf.Cone();
1210 interFunction = &FaceLineIntersector::IntersectWithCone;
1212 case GeomAbs_Sphere:
1213 intersector._sphere = surf.Sphere();
1214 interFunction = &FaceLineIntersector::IntersectWithSphere;
1217 intersector._torus = surf.Torus();
1218 interFunction = &FaceLineIntersector::IntersectWithTorus;
1221 interFunction = &FaceLineIntersector::IntersectWithSurface;
1224 _intersections.clear();
1225 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1227 if ( surf.GetType() == GeomAbs_Plane )
1229 // check if all lines in this direction are parallel to a plane
1230 if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1231 Precision::Angular()))
1233 // find out a transition, that is the same for all lines of a direction
1234 gp_Dir plnNorm = intersector._plane.Axis().Direction();
1235 gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1236 intersector._transition =
1237 ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1239 if ( surf.GetType() == GeomAbs_Cylinder )
1241 // check if all lines in this direction are parallel to a cylinder
1242 if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1243 Precision::Angular()))
1247 // intersect the grid lines with the face
1248 for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1250 GridLine& gridLine = _grid->_lines[iDir][iL];
1251 if ( _bndBox.IsOut( gridLine._line )) continue;
1253 intersector._intPoints.clear();
1254 (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1255 for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1256 _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1260 //================================================================================
1262 * Return true if (_u,_v) is on the face
1264 bool FaceLineIntersector::UVIsOnFace() const
1266 TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1267 return ( state == TopAbs_IN || state == TopAbs_ON );
1269 //================================================================================
1271 * Store an intersection if it is IN or ON the face
1273 void FaceLineIntersector::addIntPoint(const bool toClassify)
1275 if ( !toClassify || UVIsOnFace() )
1278 p._paramOnLine = _w;
1279 p._transition = _transition;
1280 _intPoints.push_back( p );
1283 //================================================================================
1285 * Intersect a line with a plane
1287 void FaceLineIntersector::IntersectWithPlane (const GridLine& gridLine)
1289 IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1290 _w = linPlane.ParamOnConic(1);
1291 if ( isParamOnLineOK( gridLine._length ))
1293 ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1297 //================================================================================
1299 * Intersect a line with a cylinder
1301 void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1303 IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1304 if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1306 _w = linCylinder.ParamOnConic(1);
1307 if ( linCylinder.NbPoints() == 1 )
1308 _transition = Trans_TANGENT;
1310 _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1311 if ( isParamOnLineOK( gridLine._length ))
1313 ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1316 if ( linCylinder.NbPoints() > 1 )
1318 _w = linCylinder.ParamOnConic(2);
1319 if ( isParamOnLineOK( gridLine._length ))
1321 ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1322 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1328 //================================================================================
1330 * Intersect a line with a cone
1332 void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1334 IntAna_IntConicQuad linCone(gridLine._line,_cone);
1335 if ( !linCone.IsDone() ) return;
1337 gp_Vec du, dv, norm;
1338 for ( int i = 1; i <= linCone.NbPoints(); ++i )
1340 _w = linCone.ParamOnConic( i );
1341 if ( !isParamOnLineOK( gridLine._length )) continue;
1342 ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1345 ElSLib::D1( _u, _v, _cone, P, du, dv );
1347 double normSize2 = norm.SquareMagnitude();
1348 if ( normSize2 > Precision::Angular() * Precision::Angular() )
1350 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1351 cos /= sqrt( normSize2 );
1352 if ( cos < -Precision::Angular() )
1353 _transition = _transIn;
1354 else if ( cos > Precision::Angular() )
1355 _transition = _transOut;
1357 _transition = Trans_TANGENT;
1361 _transition = Trans_APEX;
1363 addIntPoint( /*toClassify=*/false);
1367 //================================================================================
1369 * Intersect a line with a sphere
1371 void FaceLineIntersector::IntersectWithSphere (const GridLine& gridLine)
1373 IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1374 if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1376 _w = linSphere.ParamOnConic(1);
1377 if ( linSphere.NbPoints() == 1 )
1378 _transition = Trans_TANGENT;
1380 _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1381 if ( isParamOnLineOK( gridLine._length ))
1383 ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1386 if ( linSphere.NbPoints() > 1 )
1388 _w = linSphere.ParamOnConic(2);
1389 if ( isParamOnLineOK( gridLine._length ))
1391 ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1392 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1398 //================================================================================
1400 * Intersect a line with a torus
1402 void FaceLineIntersector::IntersectWithTorus (const GridLine& gridLine)
1404 IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1405 if ( !linTorus.IsDone()) return;
1407 gp_Vec du, dv, norm;
1408 for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1410 _w = linTorus.ParamOnLine( i );
1411 if ( !isParamOnLineOK( gridLine._length )) continue;
1412 linTorus.ParamOnTorus( i, _u,_v );
1415 ElSLib::D1( _u, _v, _torus, P, du, dv );
1417 double normSize = norm.Magnitude();
1418 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1420 if ( cos < -Precision::Angular() )
1421 _transition = _transIn;
1422 else if ( cos > Precision::Angular() )
1423 _transition = _transOut;
1425 _transition = Trans_TANGENT;
1426 addIntPoint( /*toClassify=*/false);
1430 //================================================================================
1432 * Intersect a line with a non-analytical surface
1434 void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1436 _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1437 if ( !_surfaceInt->IsDone() ) return;
1438 for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1440 _transition = Transition( _surfaceInt->Transition( i ) );
1441 _w = _surfaceInt->WParameter( i );
1442 addIntPoint(/*toClassify=*/false);
1445 //================================================================================
1447 * check if its face can be safely intersected in a thread
1449 bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1454 TopLoc_Location loc;
1455 Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1456 Handle(Geom_RectangularTrimmedSurface) ts =
1457 Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1458 while( !ts.IsNull() ) {
1459 surf = ts->BasisSurface();
1460 ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1462 if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1463 surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1464 if ( !noSafeTShapes.insert((const Standard_Transient*) _face.TShape() ).second )
1468 TopExp_Explorer exp( _face, TopAbs_EDGE );
1469 for ( ; exp.More(); exp.Next() )
1471 bool edgeIsSafe = true;
1472 const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1475 Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1478 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1479 while( !tc.IsNull() ) {
1480 c = tc->BasisCurve();
1481 tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1483 if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1484 c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1491 Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1494 Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1495 while( !tc.IsNull() ) {
1496 c2 = tc->BasisCurve();
1497 tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1499 if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1500 c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1504 if ( !edgeIsSafe && !noSafeTShapes.insert((const Standard_Transient*) e.TShape() ).second )
1509 //================================================================================
1511 * \brief Creates topology of the hexahedron
1513 Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1514 : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbFaceIntNodes(0)
1516 _polygons.reserve(100); // to avoid reallocation;
1518 //set nodes shift within grid->_nodes from the node 000
1519 size_t dx = _grid->NodeIndexDX();
1520 size_t dy = _grid->NodeIndexDY();
1521 size_t dz = _grid->NodeIndexDZ();
1523 size_t i100 = i000 + dx;
1524 size_t i010 = i000 + dy;
1525 size_t i110 = i010 + dx;
1526 size_t i001 = i000 + dz;
1527 size_t i101 = i100 + dz;
1528 size_t i011 = i010 + dz;
1529 size_t i111 = i110 + dz;
1530 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1531 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1532 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1533 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1534 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1535 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1536 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1537 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1539 vector< int > idVec;
1540 // set nodes to links
1541 for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1543 SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1544 _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1545 link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1546 link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1549 // set links to faces
1550 int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1551 for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1553 SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1554 _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1555 bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1556 faceID == SMESH_Block::ID_Fx1z ||
1557 faceID == SMESH_Block::ID_F0yz );
1558 quad._links.resize(4);
1559 vector<_OrientedLink>::iterator frwLinkIt = quad._links.begin();
1560 vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1561 for ( int i = 0; i < 4; ++i )
1563 bool revLink = revFace;
1564 if ( i > 1 ) // reverse links u1 and v0
1566 _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1567 link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1572 //================================================================================
1574 * \brief Copy constructor
1576 Hexahedron::Hexahedron( const Hexahedron& other )
1577 :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbFaceIntNodes(0)
1579 _polygons.reserve(100); // to avoid reallocation;
1581 for ( int i = 0; i < 8; ++i )
1582 _nodeShift[i] = other._nodeShift[i];
1584 for ( int i = 0; i < 12; ++i )
1586 const _Link& srcLink = other._hexLinks[ i ];
1587 _Link& tgtLink = this->_hexLinks[ i ];
1588 tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1589 tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1592 for ( int i = 0; i < 6; ++i )
1594 const _Face& srcQuad = other._hexQuads[ i ];
1595 _Face& tgtQuad = this->_hexQuads[ i ];
1596 tgtQuad._links.resize(4);
1597 for ( int j = 0; j < 4; ++j )
1599 const _OrientedLink& srcLink = srcQuad._links[ j ];
1600 _OrientedLink& tgtLink = tgtQuad._links[ j ];
1601 tgtLink._reverse = srcLink._reverse;
1602 tgtLink._link = _hexLinks + ( srcLink._link - other._hexLinks );
1607 //================================================================================
1609 * \brief Initializes its data by given grid cell
1611 void Hexahedron::init( size_t i, size_t j, size_t k )
1613 _i = i; _j = j; _k = k;
1614 // set nodes of grid to nodes of the hexahedron and
1615 // count nodes at hexahedron corners located IN and ON geometry
1616 _nbCornerNodes = _nbBndNodes = 0;
1617 _origNodeInd = _grid->NodeIndex( i,j,k );
1618 for ( int iN = 0; iN < 8; ++iN )
1620 _hexNodes[iN]._node = _grid->_nodes [ _origNodeInd + _nodeShift[iN] ];
1621 _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
1622 _nbCornerNodes += bool( _hexNodes[iN]._node );
1623 _nbBndNodes += bool( _hexNodes[iN]._intPoint );
1626 _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1627 _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1628 _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1633 if ( _nbFaceIntNodes + _eIntPoints.size() > 0 &&
1634 _nbFaceIntNodes + _nbCornerNodes + _eIntPoints.size() > 3)
1636 _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
1639 // create sub-links (_splits) by splitting links with _fIntPoints
1640 for ( int iLink = 0; iLink < 12; ++iLink )
1642 _Link& link = _hexLinks[ iLink ];
1643 link._fIntNodes.resize( link._fIntPoints.size() );
1644 for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
1646 _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
1647 link._fIntNodes[ i ] = & _intNodes.back();
1650 link._splits.clear();
1651 split._nodes[ 0 ] = link._nodes[0];
1652 bool isOut = ( ! link._nodes[0]->Node() ); // is1stNodeOut( iLink );
1653 bool checkTransition;
1654 for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
1656 if ( link._fIntNodes[i]->Node() ) // intersection non-coinsident with a grid node
1658 if ( split._nodes[ 0 ]->Node() && !isOut )
1660 split._nodes[ 1 ] = link._fIntNodes[i];
1661 link._splits.push_back( split );
1663 split._nodes[ 0 ] = link._fIntNodes[i];
1664 checkTransition = true;
1666 else // FACE intersection coinsident with a grid node
1668 checkTransition = ( link._nodes[0]->Node() );
1670 if ( checkTransition )
1672 switch ( link._fIntPoints[i]->_transition ) {
1673 case Trans_OUT: isOut = true; break;
1674 case Trans_IN : isOut = false; break;
1676 if ( !link._fIntNodes[i]->Node() && i == 0 )
1677 isOut = is1stNodeOut( link );
1679 ; // isOut remains the same
1683 if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
1685 split._nodes[ 1 ] = link._nodes[1];
1686 link._splits.push_back( split );
1690 // Create _Node's at intersections with EDGEs.
1692 const double tol2 = _grid->_tol * _grid->_tol;
1693 int facets[3], nbFacets, subEntity;
1695 for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
1697 nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
1698 _Node* equalNode = 0;
1699 switch( nbFacets ) {
1700 case 1: // in a _Face
1702 _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
1703 equalNode = FindEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1705 equalNode->Add( _eIntPoints[ iP ] );
1708 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1709 quad._eIntNodes.push_back( & _intNodes.back() );
1713 case 2: // on a _Link
1715 _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
1716 if ( link._splits.size() > 0 )
1718 equalNode = FindEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
1720 equalNode->Add( _eIntPoints[ iP ] );
1724 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1725 for ( int iF = 0; iF < 2; ++iF )
1727 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1728 equalNode = FindEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1730 equalNode->Add( _eIntPoints[ iP ] );
1733 quad._eIntNodes.push_back( & _intNodes.back() );
1739 case 3: // at a corner
1741 _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
1742 if ( node.Node() > 0 )
1744 if ( node._intPoint )
1745 node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
1749 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1750 for ( int iF = 0; iF < 3; ++iF )
1752 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1753 equalNode = FindEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1755 equalNode->Add( _eIntPoints[ iP ] );
1758 quad._eIntNodes.push_back( & _intNodes.back() );
1764 } // switch( nbFacets )
1766 if ( nbFacets == 0 ||
1767 _grid->_shapes( _eIntPoints[ iP ]->_shapeID ).ShapeType() == TopAbs_VERTEX )
1769 equalNode = FindEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
1771 equalNode->Add( _eIntPoints[ iP ] );
1774 if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
1775 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1776 _vIntNodes.push_back( & _intNodes.back() );
1779 } // loop on _eIntPoints
1781 else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
1784 // create sub-links (_splits) of whole links
1785 for ( int iLink = 0; iLink < 12; ++iLink )
1787 _Link& link = _hexLinks[ iLink ];
1788 link._splits.clear();
1789 if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
1791 split._nodes[ 0 ] = link._nodes[0];
1792 split._nodes[ 1 ] = link._nodes[1];
1793 link._splits.push_back( split );
1799 //================================================================================
1801 * \brief Initializes its data by given grid cell (countered from zero)
1803 void Hexahedron::init( size_t iCell )
1805 size_t iNbCell = _grid->_coords[0].size() - 1;
1806 size_t jNbCell = _grid->_coords[1].size() - 1;
1807 _i = iCell % iNbCell;
1808 _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1809 _k = iCell / iNbCell / jNbCell;
1813 //================================================================================
1815 * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1817 void Hexahedron::ComputeElements()
1821 int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
1822 if ( _nbCornerNodes + nbIntersections < 4 )
1825 if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
1829 _polygons.reserve( 20 );
1831 // Create polygons from quadrangles
1832 // --------------------------------
1835 vector< _OrientedLink > splits;
1836 vector<_Node*> chainNodes, usedEdgeNodes;
1837 _Face* coplanarPolyg;
1839 bool hasEdgeIntersections = !_eIntPoints.empty();
1841 for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1843 _Face& quad = _hexQuads[ iF ] ;
1845 _polygons.resize( _polygons.size() + 1 );
1846 _Face* polygon = &_polygons.back();
1847 polygon->_polyLinks.reserve( 20 );
1850 for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1851 for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
1852 splits.push_back( quad._links[ iE ].ResultLink( iS ));
1854 // add splits of links to a polygon and add _polyLinks to make
1855 // polygon's boundary closed
1857 int nbSplits = splits.size();
1858 if ( nbSplits < 2 && quad._eIntNodes.empty() )
1862 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1863 if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
1864 quad._eIntNodes[ iP ]->_usedInFace = 0;
1866 int nbUsedEdgeNodes = 0;
1868 while ( nbSplits > 0 )
1871 while ( !splits[ iS ] )
1874 if ( !polygon->_links.empty() )
1876 _polygons.resize( _polygons.size() + 1 );
1877 polygon = &_polygons.back();
1878 polygon->_polyLinks.reserve( 20 );
1880 polygon->_links.push_back( splits[ iS ] );
1881 splits[ iS++ ]._link = 0;
1884 _Node* nFirst = polygon->_links.back().FirstNode();
1885 _Node *n1,*n2 = polygon->_links.back().LastNode();
1886 for ( ; nFirst != n2 && iS < splits.size(); ++iS )
1888 _OrientedLink& split = splits[ iS ];
1889 if ( !split ) continue;
1891 n1 = split.FirstNode();
1894 // try to connect to intersections with EDGEs
1895 if ( quad._eIntNodes.size() > nbUsedEdgeNodes &&
1896 findChain( n2, n1, quad, chainNodes ))
1898 for ( size_t i = 1; i < chainNodes.size(); ++i )
1900 polyLink._nodes[0] = chainNodes[i-1];
1901 polyLink._nodes[1] = chainNodes[i];
1902 polygon->_polyLinks.push_back( polyLink );
1903 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1904 nbUsedEdgeNodes += ( polyLink._nodes[1]->IsUsedInFace( polygon ));
1906 if ( chainNodes.back() != n1 )
1908 n2 = chainNodes.back();
1913 // try to connect to a split ending on the same FACE
1916 _OrientedLink foundSplit;
1917 for ( int i = iS; i < splits.size() && !foundSplit; ++i )
1918 if (( foundSplit = splits[ i ]) &&
1919 ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
1921 polyLink._nodes[0] = n2;
1922 polyLink._nodes[1] = foundSplit.FirstNode();
1923 polygon->_polyLinks.push_back( polyLink );
1924 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1929 foundSplit._link = 0;
1933 n2 = foundSplit.FirstNode();
1938 if ( n2->IsLinked( nFirst->_intPoint ))
1940 polyLink._nodes[0] = n2;
1941 polyLink._nodes[1] = n1;
1942 polygon->_polyLinks.push_back( polyLink );
1943 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1947 polygon->_links.push_back( split );
1950 n2 = polygon->_links.back().LastNode();
1954 if ( nFirst != n2 ) // close a polygon
1956 if ( !findChain( n2, nFirst, quad, chainNodes ))
1958 if ( !closePolygon( polygon, chainNodes ))
1959 chainNodes.push_back( nFirst );
1961 for ( size_t i = 1; i < chainNodes.size(); ++i )
1963 polyLink._nodes[0] = chainNodes[i-1];
1964 polyLink._nodes[1] = chainNodes[i];
1965 polygon->_polyLinks.push_back( polyLink );
1966 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1967 nbUsedEdgeNodes += bool( polyLink._nodes[1]->IsUsedInFace( polygon ));
1971 if ( polygon->_links.size() < 3 && nbSplits > 0 )
1973 polygon->_polyLinks.clear();
1974 polygon->_links.clear();
1976 } // while ( nbSplits > 0 )
1978 // if ( quad._eIntNodes.size() > nbUsedEdgeNodes )
1980 // // make _vIntNodes from not used _eIntNodes
1981 // const double tol = 0.05 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
1982 // for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1984 // if ( quad._eIntNodes[ iP ]->IsUsedInFace() ) continue;
1985 // _Node* equalNode =
1986 // FindEqualNode( _vIntNodes, quad._eIntNodes[ iP ].EdgeIntPnt(), tol*tol );
1988 // equalNode->Add( quad._eIntNodes[ iP ].EdgeIntPnt() );
1990 // _vIntNodes.push_back( quad._eIntNodes[ iP ]);
1994 if ( polygon->_links.size() < 3 )
1996 _polygons.pop_back();
1997 //usedEdgeNodes.resize( usedEdgeNodes.size() - nbUsedEdgeNodes );
1999 } // loop on 6 hexahedron sides
2001 // Create polygons closing holes in a polyhedron
2002 // ----------------------------------------------
2004 // clear _usedInFace
2005 for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2006 _intNodes[ iN ]._usedInFace = 0;
2008 // add polygons to their links and mark used nodes
2009 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2011 _Face& polygon = _polygons[ iP ];
2012 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2014 polygon._links[ iL ].AddFace( &polygon );
2015 polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
2019 vector< _OrientedLink* > freeLinks;
2020 freeLinks.reserve(20);
2021 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2023 _Face& polygon = _polygons[ iP ];
2024 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2025 if ( polygon._links[ iL ].NbFaces() < 2 )
2027 freeLinks.push_back( & polygon._links[ iL ]);
2028 freeLinks.back()->FirstNode()->IsUsedInFace() == true;
2031 int nbFreeLinks = freeLinks.size();
2032 if ( nbFreeLinks > 0 && nbFreeLinks < 3 ) return;
2034 // put not used intersection nodes to _vIntNodes
2035 int nbVertexNodes = 0; // nb not used vertex nodes
2037 for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2038 nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
2040 const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
2041 for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
2043 if ( _intNodes[ iN ].IsUsedInFace() ) continue;
2044 if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
2046 FindEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
2047 if ( !equalNode /*|| equalNode->IsUsedInFace()*/ )
2049 _vIntNodes.push_back( &_intNodes[ iN ]);
2055 set<TGeomID> usedFaceIDs;
2056 TGeomID curFace = 0;
2057 const size_t nbQuadPolygons = _polygons.size();
2059 // create polygons by making closed chains of free links
2060 size_t iPolygon = _polygons.size();
2061 while ( nbFreeLinks > 0 )
2063 if ( iPolygon == _polygons.size() )
2064 _polygons.resize( _polygons.size() + 1 );
2065 _Face& polygon = _polygons[ iPolygon ];
2066 polygon._polyLinks.reserve( 20 );
2067 polygon._links.reserve( 20 );
2069 _OrientedLink* curLink = 0;
2071 if (( !hasEdgeIntersections ) ||
2072 ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
2074 // get a remaining link to start from
2075 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2076 if (( curLink = freeLinks[ iL ] ))
2077 freeLinks[ iL ] = 0;
2078 polygon._links.push_back( *curLink );
2082 // find all links connected to curLink
2083 curNode = curLink->FirstNode();
2085 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2086 if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
2088 curLink = freeLinks[ iL ];
2089 freeLinks[ iL ] = 0;
2091 polygon._links.push_back( *curLink );
2093 } while ( curLink );
2095 else // there are intersections with EDGEs
2097 // get a remaining link to start from, one lying on minimal
2100 vector< pair< TGeomID, int > > facesOfLink[3];
2101 pair< TGeomID, int > faceOfLink( -1, -1 );
2102 vector< TGeomID > faces;
2103 for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2104 if ( freeLinks[ iL ] )
2106 faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2107 if ( faces.size() == 1 )
2109 faceOfLink = make_pair( faces[0], iL );
2110 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2112 facesOfLink[0].push_back( faceOfLink );
2114 else if ( facesOfLink[0].empty() )
2116 faceOfLink = make_pair(( faces.empty() ? -1 : faces[0]), iL );
2117 facesOfLink[ 1 + faces.empty() ].push_back( faceOfLink );
2120 for ( int i = 0; faceOfLink.second < 0 && i < 3; ++i )
2121 if ( !facesOfLink[i].empty() )
2122 faceOfLink = facesOfLink[i][0];
2124 if ( faceOfLink.first < 0 ) // all faces used
2126 for ( size_t i = 0; i < facesOfLink[2].size() && faceOfLink.first < 1; ++i )
2128 curLink = freeLinks[ facesOfLink[2][i].second ];
2129 faceOfLink.first = curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2131 usedFaceIDs.clear();
2133 curFace = faceOfLink.first;
2134 curLink = freeLinks[ faceOfLink.second ];
2135 freeLinks[ faceOfLink.second ] = 0;
2137 usedFaceIDs.insert( curFace );
2138 polygon._links.push_back( *curLink );
2141 // find all links bounding a FACE of curLink
2144 // go forward from curLink
2145 curNode = curLink->LastNode();
2147 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2148 if ( freeLinks[ iL ] &&
2149 freeLinks[ iL ]->FirstNode() == curNode &&
2150 freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2152 curLink = freeLinks[ iL ];
2153 freeLinks[ iL ] = 0;
2154 polygon._links.push_back( *curLink );
2157 } while ( curLink );
2159 std::reverse( polygon._links.begin(), polygon._links.end() );
2161 curLink = & polygon._links.back();
2164 // go backward from curLink
2165 curNode = curLink->FirstNode();
2167 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2168 if ( freeLinks[ iL ] &&
2169 freeLinks[ iL ]->LastNode() == curNode &&
2170 freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2172 curLink = freeLinks[ iL ];
2173 freeLinks[ iL ] = 0;
2174 polygon._links.push_back( *curLink );
2177 } while ( curLink );
2179 curNode = polygon._links.back().FirstNode();
2181 if ( polygon._links[0].LastNode() != curNode )
2183 if ( nbVertexNodes > 0 )
2185 // add links with _vIntNodes if not already used
2186 for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2187 if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2188 _vIntNodes[ iN ]->IsOnFace( curFace ))
2190 _vIntNodes[ iN ]->_usedInFace = &polygon;
2192 polyLink._nodes[0] = _vIntNodes[ iN ];
2193 polyLink._nodes[1] = curNode;
2194 polygon._polyLinks.push_back( polyLink );
2195 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
2196 freeLinks.push_back( &polygon._links.back() );
2198 curNode = _vIntNodes[ iN ];
2199 // TODO: to reorder _vIntNodes within polygon, if there are several ones
2202 // if ( polygon._links.size() > 1 )
2204 polyLink._nodes[0] = polygon._links[0].LastNode();
2205 polyLink._nodes[1] = curNode;
2206 polygon._polyLinks.push_back( polyLink );
2207 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
2208 freeLinks.push_back( &polygon._links.back() );
2212 } // if there are intersections with EDGEs
2214 if ( polygon._links.size() < 2 ||
2215 polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
2216 return; // closed polygon not found -> invalid polyhedron
2218 if ( polygon._links.size() == 2 )
2220 if ( freeLinks.back() == &polygon._links.back() )
2222 freeLinks.pop_back();
2225 if ( polygon._links.front().NbFaces() > 0 )
2226 polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
2227 if ( polygon._links.back().NbFaces() > 0 )
2228 polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
2230 _polygons.pop_back();
2232 else // polygon._links.size() >= 2
2234 // add polygon to its links
2235 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2237 polygon._links[ iL ].AddFace( &polygon );
2238 polygon._links[ iL ].Reverse();
2240 if ( hasEdgeIntersections && iPolygon == _polygons.size() - 1 )
2242 // check that a polygon does not lie in the plane of another polygon
2244 for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
2246 if ( polygon._links[ iL ].NbFaces() < 2 )
2247 continue; // it's a just added free link
2248 // look for a polygon made on a hexa side and sharing
2249 // two or more haxa links
2251 coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
2252 for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
2253 if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
2254 !coplanarPolyg->isPolyLink( polygon._links[ iL2 ]) &&
2255 coplanarPolyg < & _polygons[ nbQuadPolygons ])
2257 if ( iL2 == polygon._links.size() )
2260 if ( 0 /*coplanarPolyg*/ ) // coplanar polygon found
2262 freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
2263 nbFreeLinks -= polygon._polyLinks.size();
2265 // fill freeLinks with links not shared by coplanarPolyg and polygon
2266 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2267 if ( polygon._links[ iL ]._link->_faces[1] &&
2268 polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
2270 _Face* p = polygon._links[ iL ]._link->_faces[0];
2271 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2272 if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
2274 freeLinks.push_back( & p->_links[ iL2 ] );
2276 freeLinks.back()->RemoveFace( &polygon );
2280 for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
2281 if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
2282 coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
2284 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
2285 if ( p == coplanarPolyg )
2286 p = coplanarPolyg->_links[ iL ]._link->_faces[1];
2287 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2288 if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
2290 freeLinks.push_back( & p->_links[ iL2 ] );
2292 freeLinks.back()->RemoveFace( coplanarPolyg );
2296 // set coplanarPolyg to be re-created next
2297 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2298 if ( coplanarPolyg == & _polygons[ iP ] )
2301 _polygons[ iPolygon ]._links.clear();
2302 _polygons[ iPolygon ]._polyLinks.clear();
2305 if ( freeLinks.back() == &polygon._links.back() )
2307 freeLinks.pop_back();
2310 _polygons.pop_back();
2311 usedFaceIDs.erase( curFace );
2313 } // if ( coplanarPolyg )
2314 } // if ( hasEdgeIntersections )
2316 iPolygon = _polygons.size();
2318 } // end of case ( polygon._links.size() > 2 )
2319 } // while ( nbFreeLinks > 0 )
2321 if ( ! checkPolyhedronSize() )
2326 // create a classic cell if possible
2327 const int nbNodes = _nbCornerNodes + nbIntersections;
2328 bool isClassicElem = false;
2329 if ( nbNodes == 8 && _polygons.size() == 6 ) isClassicElem = addHexa();
2330 else if ( nbNodes == 4 && _polygons.size() == 4 ) isClassicElem = addTetra();
2331 else if ( nbNodes == 6 && _polygons.size() == 5 ) isClassicElem = addPenta();
2332 else if ( nbNodes == 5 && _polygons.size() == 5 ) isClassicElem = addPyra ();
2333 if ( !isClassicElem )
2335 _volumeDefs._nodes.clear();
2336 _volumeDefs._quantities.clear();
2338 for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2340 const size_t nbLinks = _polygons[ iF ]._links.size();
2341 _volumeDefs._quantities.push_back( nbLinks );
2342 for ( size_t iL = 0; iL < nbLinks; ++iL )
2343 _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
2347 //================================================================================
2349 * \brief Create elements in the mesh
2351 int Hexahedron::MakeElements(SMESH_MesherHelper& helper,
2352 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2354 SMESHDS_Mesh* mesh = helper.GetMeshDS();
2356 size_t nbCells[3] = { _grid->_coords[0].size() - 1,
2357 _grid->_coords[1].size() - 1,
2358 _grid->_coords[2].size() - 1 };
2359 const size_t nbGridCells = nbCells[0] * nbCells[1] * nbCells[2];
2360 vector< Hexahedron* > intersectedHex( nbGridCells, 0 );
2363 // set intersection nodes from GridLine's to links of intersectedHex
2364 int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
2365 for ( int iDir = 0; iDir < 3; ++iDir )
2367 int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
2368 dInd[1][ iDirOther[iDir][0] ] = -1;
2369 dInd[2][ iDirOther[iDir][1] ] = -1;
2370 dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
2371 // loop on GridLine's parallel to iDir
2372 LineIndexer lineInd = _grid->GetLineIndexer( iDir );
2373 for ( ; lineInd.More(); ++lineInd )
2375 GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
2376 multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2377 for ( ; ip != line._intPoints.end(); ++ip )
2379 // if ( !ip->_node ) continue; // intersection at a grid node
2380 lineInd.SetIndexOnLine( ip->_indexOnLine );
2381 for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
2383 i = int(lineInd.I()) + dInd[iL][0];
2384 j = int(lineInd.J()) + dInd[iL][1];
2385 k = int(lineInd.K()) + dInd[iL][2];
2386 if ( i < 0 || i >= nbCells[0] ||
2387 j < 0 || j >= nbCells[1] ||
2388 k < 0 || k >= nbCells[2] ) continue;
2390 const size_t hexIndex = _grid->CellIndex( i,j,k );
2391 Hexahedron *& hex = intersectedHex[ hexIndex ];
2394 hex = new Hexahedron( *this );
2400 const int iLink = iL + iDir * 4;
2401 hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
2402 //hex->_hexLinks[iLink]._fIntNodes.push_back( _Node( 0, &(*ip) ));
2403 hex->_nbFaceIntNodes += bool( ip->_node );
2409 // implement geom edges into the mesh
2410 addEdges( helper, intersectedHex, edge2faceIDsMap );
2412 // add not split hexadrons to the mesh
2414 vector<int> intHexInd( nbIntHex );
2416 for ( size_t i = 0; i < intersectedHex.size(); ++i )
2418 Hexahedron * & hex = intersectedHex[ i ];
2421 intHexInd[ nbIntHex++ ] = i;
2422 if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
2423 continue; // treat intersected hex later
2424 this->init( hex->_i, hex->_j, hex->_k );
2430 if (( _nbCornerNodes == 8 ) &&
2431 ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
2433 // order of _hexNodes is defined by enum SMESH_Block::TShapeID
2434 SMDS_MeshElement* el =
2435 mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
2436 _hexNodes[3].Node(), _hexNodes[1].Node(),
2437 _hexNodes[4].Node(), _hexNodes[6].Node(),
2438 _hexNodes[7].Node(), _hexNodes[5].Node() );
2439 mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
2444 intersectedHex[ i ] = 0;
2448 else if ( _nbCornerNodes > 3 && !hex )
2450 // all intersection of hex with geometry are at grid nodes
2451 hex = new Hexahedron( *this );
2456 intHexInd.push_back(0);
2457 intHexInd[ nbIntHex++ ] = i;
2461 // add elements resulted from hexadron intersection
2463 intHexInd.resize( nbIntHex );
2464 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, nbIntHex ),
2465 ParallelHexahedron( intersectedHex, intHexInd ),
2466 tbb::simple_partitioner()); // ComputeElements() is called here
2467 for ( size_t i = 0; i < intHexInd.size(); ++i )
2468 if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
2469 nbAdded += hex->addElements( helper );
2471 for ( size_t i = 0; i < intHexInd.size(); ++i )
2472 if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
2474 hex->ComputeElements();
2475 nbAdded += hex->addElements( helper );
2479 for ( size_t i = 0; i < intersectedHex.size(); ++i )
2480 if ( intersectedHex[ i ] )
2481 delete intersectedHex[ i ];
2486 //================================================================================
2488 * \brief Implements geom edges into the mesh
2490 void Hexahedron::addEdges(SMESH_MesherHelper& helper,
2491 vector< Hexahedron* >& hexes,
2492 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2494 if ( edge2faceIDsMap.empty() ) return;
2496 // Prepare planes for intersecting with EDGEs
2499 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
2501 GridPlanes& planes = pln[ iDirZ ];
2502 int iDirX = ( iDirZ + 1 ) % 3;
2503 int iDirY = ( iDirZ + 2 ) % 3;
2504 // planes._uNorm = ( _grid->_axes[ iDirY ] ^ _grid->_axes[ iDirZ ] ).Normalized();
2505 // planes._vNorm = ( _grid->_axes[ iDirZ ] ^ _grid->_axes[ iDirX ] ).Normalized();
2506 planes._zNorm = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
2507 planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
2508 planes._zProjs [0] = 0;
2509 const double zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2510 const vector< double > & u = _grid->_coords[ iDirZ ];
2511 for ( int i = 1; i < planes._zProjs.size(); ++i )
2513 planes._zProjs [i] = zFactor * ( u[i] - u[0] );
2517 const double deflection = _grid->_minCellSize / 20.;
2518 const double tol = _grid->_tol;
2519 E_IntersectPoint ip;
2521 // Intersect EDGEs with the planes
2522 map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
2523 for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
2525 const TGeomID edgeID = e2fIt->first;
2526 const TopoDS_Edge & E = TopoDS::Edge( _grid->_shapes( edgeID ));
2527 BRepAdaptor_Curve curve( E );
2528 TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
2529 TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
2531 ip._faceIDs = e2fIt->second;
2532 ip._shapeID = edgeID;
2534 // discretize the EGDE
2535 GCPnts_UniformDeflection discret( curve, deflection, true );
2536 if ( !discret.IsDone() || discret.NbPoints() < 2 )
2539 // perform intersection
2540 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
2542 GridPlanes& planes = pln[ iDirZ ];
2543 int iDirX = ( iDirZ + 1 ) % 3;
2544 int iDirY = ( iDirZ + 2 ) % 3;
2545 double xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
2546 double yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
2547 double zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
2548 //double zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2549 int dIJK[3], d000[3] = { 0,0,0 };
2550 double o[3] = { _grid->_coords[0][0],
2551 _grid->_coords[1][0],
2552 _grid->_coords[2][0] };
2554 // locate the 1st point of a segment within the grid
2555 gp_XYZ p1 = discret.Value( 1 ).XYZ();
2556 double u1 = discret.Parameter( 1 );
2557 double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
2559 _grid->ComputeUVW( p1, ip._uvw );
2560 int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
2561 int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
2562 int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
2563 locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
2564 locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
2565 locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
2567 int ijk[3]; // grid index where a segment intersect a plane
2572 // add the 1st vertex point to a hexahedron
2576 ip._shapeID = _grid->_shapes.Add( v1 );
2577 _grid->_edgeIntP.push_back( ip );
2578 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2579 _grid->_edgeIntP.pop_back();
2580 ip._shapeID = edgeID;
2582 for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
2584 // locate the 2nd point of a segment within the grid
2585 gp_XYZ p2 = discret.Value( iP ).XYZ();
2586 double u2 = discret.Parameter( iP );
2587 double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
2589 if ( Abs( zProj2 - zProj1 ) <= std::numeric_limits<double>::min() )
2591 locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
2593 // treat intersections with planes between 2 end points of a segment
2594 int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
2595 int iZ = iZ1 + ( iZ1 < iZ2 );
2596 for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
2598 ip._point = findIntPoint( u1, zProj1, u2, zProj2,
2599 planes._zProjs[ iZ ],
2600 curve, planes._zNorm, _grid->_origin );
2601 _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
2602 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2603 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2606 // add ip to hex "above" the plane
2607 _grid->_edgeIntP.push_back( ip );
2609 bool added = addIntersection(_grid->_edgeIntP.back(), hexes, ijk, dIJK);
2611 // add ip to hex "below" the plane
2612 ijk[ iDirZ ] = iZ-1;
2613 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, dIJK ) &&
2615 _grid->_edgeIntP.pop_back();
2622 // add the 2nd vertex point to a hexahedron
2625 ip._shapeID = _grid->_shapes.Add( v2 );
2627 _grid->ComputeUVW( p1, ip._uvw );
2628 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2629 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2631 _grid->_edgeIntP.push_back( ip );
2632 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2633 _grid->_edgeIntP.pop_back();
2634 ip._shapeID = edgeID;
2636 } // loop on 3 grid directions
2639 // Create nodes at found intersections
2640 // const E_IntersectPoint* eip;
2641 // for ( size_t i = 0; i < hexes.size(); ++i )
2643 // Hexahedron* h = hexes[i];
2644 // if ( !h ) continue;
2645 // for ( int iF = 0; iF < 6; ++iF )
2647 // _Face& quad = h->_hexQuads[ iF ];
2648 // for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2649 // if ( !quad._eIntNodes[ iP ]._node )
2650 // if (( eip = quad._eIntNodes[ iP ].EdgeIntPnt() ))
2651 // quad._eIntNodes[ iP ]._intPoint->_node = helper.AddNode( eip->_point.X(),
2653 // eip->_point.Z() );
2655 // for ( size_t iP = 0; iP < hexes[i]->_vIntNodes.size(); ++iP )
2656 // if (( eip = h->_vIntNodes[ iP ].EdgeIntPnt() ))
2657 // h->_vIntNodes[ iP ]._intPoint->_node = helper.AddNode( eip->_point.X(),
2659 // eip->_point.Z() );
2663 //================================================================================
2665 * \brief Finds intersection of a curve with a plane
2666 * \param [in] u1 - parameter of one curve point
2667 * \param [in] proj1 - projection of the curve point to the plane normal
2668 * \param [in] u2 - parameter of another curve point
2669 * \param [in] proj2 - projection of the other curve point to the plane normal
2670 * \param [in] proj - projection of a point where the curve intersects the plane
2671 * \param [in] curve - the curve
2672 * \param [in] axis - the plane normal
2673 * \param [in] origin - the plane origin
2674 * \return gp_Pnt - the found intersection point
2676 gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
2677 double u2, double proj2,
2679 BRepAdaptor_Curve& curve,
2681 const gp_XYZ& origin)
2683 double r = (( proj - proj1 ) / ( proj2 - proj1 ));
2684 double u = u1 * ( 1 - r ) + u2 * r;
2685 gp_Pnt p = curve.Value( u );
2686 double newProj = axis * ( p.XYZ() - origin );
2687 if ( Abs( proj - newProj ) > _grid->_tol / 10. )
2690 return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
2692 return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
2697 //================================================================================
2699 * \brief Returns indices of a hexahedron sub-entities holding a point
2700 * \param [in] ip - intersection point
2701 * \param [out] facets - 0-3 facets holding a point
2702 * \param [out] sub - index of a vertex or an edge holding a point
2703 * \return int - number of facets holding a point
2705 int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
2707 enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
2709 int vertex = 0, egdeMask = 0;
2711 if ( Abs( _grid->_coords[0][ _i ] - ip->_uvw[0] ) < _grid->_tol ) {
2712 facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
2715 else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
2716 facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
2720 if ( Abs( _grid->_coords[1][ _j ] - ip->_uvw[1] ) < _grid->_tol ) {
2721 facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
2724 else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
2725 facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
2729 if ( Abs( _grid->_coords[2][ _k ] - ip->_uvw[2] ) < _grid->_tol ) {
2730 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
2733 else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
2734 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
2741 case 0: sub = 0; break;
2742 case 1: sub = facets[0]; break;
2744 const int edge [3][8] = {
2745 { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
2746 SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
2747 { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
2748 SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
2749 { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
2750 SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
2752 switch ( egdeMask ) {
2753 case X | Y: sub = edge[ 0 ][ vertex ]; break;
2754 case X | Z: sub = edge[ 1 ][ vertex ]; break;
2755 default: sub = edge[ 2 ][ vertex ];
2761 sub = vertex + SMESH_Block::ID_FirstV;
2766 //================================================================================
2768 * \brief Adds intersection with an EDGE
2770 bool Hexahedron::addIntersection( const E_IntersectPoint& ip,
2771 vector< Hexahedron* >& hexes,
2772 int ijk[], int dIJK[] )
2776 size_t hexIndex[4] = {
2777 _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
2778 dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
2779 dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
2780 dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
2782 for ( int i = 0; i < 4; ++i )
2784 if ( /*0 <= hexIndex[i] &&*/ hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
2786 Hexahedron* h = hexes[ hexIndex[i] ];
2787 // check if ip is really inside the hex
2789 if (( _grid->_coords[0][ h->_i ] - _grid->_tol > ip._uvw[0] ) ||
2790 ( _grid->_coords[0][ h->_i+1 ] + _grid->_tol < ip._uvw[0] ) ||
2791 ( _grid->_coords[1][ h->_j ] - _grid->_tol > ip._uvw[1] ) ||
2792 ( _grid->_coords[1][ h->_j+1 ] + _grid->_tol < ip._uvw[1] ) ||
2793 ( _grid->_coords[2][ h->_k ] - _grid->_tol > ip._uvw[2] ) ||
2794 ( _grid->_coords[2][ h->_k+1 ] + _grid->_tol < ip._uvw[2] ))
2795 throw SALOME_Exception("ip outside a hex");
2797 h->_eIntPoints.push_back( & ip );
2803 //================================================================================
2805 * \brief Finds nodes at a path from one node to another via intersections with EDGEs
2807 bool Hexahedron::findChain( _Node* n1,
2810 vector<_Node*>& chn )
2813 chn.push_back( n1 );
2814 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2815 if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2816 n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
2817 n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2819 chn.push_back( quad._eIntNodes[ iP ]);
2820 chn.push_back( n2 );
2821 quad._eIntNodes[ iP ]->_usedInFace = &quad;
2828 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2829 if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2830 chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2832 chn.push_back( quad._eIntNodes[ iP ]);
2833 found = quad._eIntNodes[ iP ]->_usedInFace = &quad;
2836 } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
2838 if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
2839 chn.push_back( n2 );
2841 return chn.size() > 1;
2843 //================================================================================
2845 * \brief Try to heal a polygon whose ends are not connected
2847 bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
2849 int i = -1, nbLinks = polygon->_links.size();
2852 vector< _OrientedLink > newLinks;
2853 // find a node lying on the same FACE as the last one
2854 _Node* node = polygon->_links.back().LastNode();
2855 int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
2856 for ( i = nbLinks - 2; i >= 0; --i )
2857 if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
2861 for ( ; i < nbLinks; ++i )
2862 newLinks.push_back( polygon->_links[i] );
2866 // find a node lying on the same FACE as the first one
2867 node = polygon->_links[0].FirstNode();
2868 avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
2869 for ( i = 1; i < nbLinks; ++i )
2870 if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
2873 for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
2874 newLinks.push_back( polygon->_links[i] );
2876 if ( newLinks.size() > 1 )
2878 polygon->_links.swap( newLinks );
2880 chainNodes.push_back( polygon->_links.back().LastNode() );
2881 chainNodes.push_back( polygon->_links[0].FirstNode() );
2886 //================================================================================
2888 * \brief Checks transition at the 1st node of a link
2890 bool Hexahedron::is1stNodeOut( _Link& link /*int iLink*/ ) const
2892 // new version is for the case: tangent transition at the 1st node
2894 if ( link._fIntNodes.size() > 1 )
2896 // check transition at the next intersection
2897 switch ( link._fIntPoints[1]->_transition ) {
2898 case Trans_OUT: return false;
2899 case Trans_IN : return true;
2900 default: ; // tangent transition
2903 gp_Pnt p1 = link._nodes[0]->Point();
2904 gp_Pnt p2 = link._nodes[1]->Point();
2905 gp_Pnt testPnt = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
2907 TGeomID faceID = link._fIntPoints[0]->_faceIDs[0];
2908 const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( faceID ));
2909 TopLoc_Location loc;
2910 GeomAPI_ProjectPointOnSurf& proj =
2911 _grid->_helper->GetProjector( face, loc, 0.1*_grid->_tol );
2912 testPnt.Transform( loc );
2913 proj.Perform( testPnt );
2914 if ( proj.IsDone() &&
2915 proj.NbPoints() > 0 &&
2916 proj.LowerDistance() > _grid->_tol )
2918 Quantity_Parameter u,v;
2919 proj.LowerDistanceParameters( u,v );
2921 if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
2926 if ( face.Orientation() == TopAbs_REVERSED )
2928 gp_Vec v( proj.NearestPoint(), testPnt );
2929 return v * normal > 0;
2932 // if ( !_hexLinks[ iLink ]._nodes[0]->Node() ) // no node
2934 // if ( !_hexLinks[ iLink ]._nodes[0]->_intPoint ) // no intersection with geometry
2936 // switch ( _hexLinks[ iLink ]._nodes[0]->FaceIntPnt()->_transition ) {
2937 // case Trans_OUT: return true;
2938 // case Trans_IN : return false;
2939 // default: ; // tangent transition
2942 // // ijk of a GridLine corresponding to the link
2943 // int iDir = iLink / 4;
2944 // int indSub = iLink % 4;
2945 // LineIndexer li = _grid->GetLineIndexer( iDir );
2946 // li.SetIJK( _i,_j,_k );
2947 // size_t lineIndex[4] = { li.LineIndex (),
2948 // li.LineIndex10(),
2949 // li.LineIndex01(),
2950 // li.LineIndex11() };
2951 // GridLine& line = _grid->_lines[ iDir ][ lineIndex[ indSub ]];
2953 // // analyze transition of previous ip
2954 // bool isOut = true;
2955 // multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2956 // for ( ; ip != line._intPoints.end(); ++ip )
2958 // if ( &(*ip) == _hexLinks[ iLink ]._nodes[0]->_intPoint )
2960 // switch ( ip->_transition ) {
2961 // case Trans_OUT: isOut = true;
2962 // case Trans_IN : isOut = false;
2967 // if ( ip == line._intPoints.end() )
2968 // cout << "BUG: Wrong GridLine. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl;
2972 //================================================================================
2974 * \brief Adds computed elements to the mesh
2976 int Hexahedron::addElements(SMESH_MesherHelper& helper)
2979 // add elements resulted from hexahedron intersection
2980 //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
2982 vector< const SMDS_MeshNode* > nodes( _volumeDefs._nodes.size() );
2983 for ( size_t iN = 0; iN < nodes.size(); ++iN )
2984 if ( !( nodes[iN] = _volumeDefs._nodes[iN]->Node() ))
2986 if ( const E_IntersectPoint* eip = _volumeDefs._nodes[iN]->EdgeIntPnt() )
2987 nodes[iN] = _volumeDefs._nodes[iN]->_intPoint->_node =
2988 helper.AddNode( eip->_point.X(),
2992 throw SALOME_Exception("Bug: no node at intersection point");
2995 if ( !_volumeDefs._quantities.empty() )
2997 helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
3001 switch ( nodes.size() )
3003 case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
3004 nodes[4],nodes[5],nodes[6],nodes[7] );
3006 case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
3008 case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
3011 helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
3015 nbAdded += int ( _volumeDefs._nodes.size() > 0 );
3020 //================================================================================
3022 * \brief Return true if the element is in a hole
3024 bool Hexahedron::isInHole() const
3026 if ( !_vIntNodes.empty() )
3029 const int ijk[3] = { _i, _j, _k };
3030 F_IntersectPoint curIntPnt;
3032 // consider a cell to be in a hole if all links in any direction
3033 // comes OUT of geometry
3034 for ( int iDir = 0; iDir < 3; ++iDir )
3036 const vector<double>& coords = _grid->_coords[ iDir ];
3037 LineIndexer li = _grid->GetLineIndexer( iDir );
3038 li.SetIJK( _i,_j,_k );
3039 size_t lineIndex[4] = { li.LineIndex (),
3043 bool allLinksOut = true, hasLinks = false;
3044 for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
3046 const _Link& link = _hexLinks[ iL + 4*iDir ];
3047 // check transition of the first node of a link
3048 const F_IntersectPoint* firstIntPnt = 0;
3049 if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
3051 curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0];
3052 const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
3053 multiset< F_IntersectPoint >::const_iterator ip =
3054 line._intPoints.upper_bound( curIntPnt );
3056 firstIntPnt = &(*ip);
3058 else if ( !link._fIntPoints.empty() )
3060 firstIntPnt = link._fIntPoints[0];
3066 allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
3069 if ( hasLinks && allLinksOut )
3075 //================================================================================
3077 * \brief Return true if a polyhedron passes _sizeThreshold criterion
3079 bool Hexahedron::checkPolyhedronSize() const
3082 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3084 const _Face& polygon = _polygons[iP];
3085 gp_XYZ area (0,0,0);
3086 gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
3087 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3089 gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
3093 volume += p1 * area;
3097 double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
3099 return volume > initVolume / _sizeThreshold;
3101 //================================================================================
3103 * \brief Tries to create a hexahedron
3105 bool Hexahedron::addHexa()
3107 if ( _polygons[0]._links.size() != 4 ||
3108 _polygons[1]._links.size() != 4 ||
3109 _polygons[2]._links.size() != 4 ||
3110 _polygons[3]._links.size() != 4 ||
3111 _polygons[4]._links.size() != 4 ||
3112 _polygons[5]._links.size() != 4 )
3116 for ( int iL = 0; iL < 4; ++iL )
3119 nodes[iL] = _polygons[0]._links[iL].FirstNode();
3122 // find a top node above the base node
3123 _Link* link = _polygons[0]._links[iL]._link;
3124 //ASSERT( link->_faces.size() > 1 );
3125 if ( !link->_faces[0] || !link->_faces[1] )
3126 return debugDumpLink( link );
3127 // a quadrangle sharing <link> with _polygons[0]
3128 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
3129 for ( int i = 0; i < 4; ++i )
3130 if ( quad->_links[i]._link == link )
3132 // 1st node of a link opposite to <link> in <quad>
3133 nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
3139 _volumeDefs.set( vector< _Node* >( nodes, nodes+8 ));
3143 //================================================================================
3145 * \brief Tries to create a tetrahedron
3147 bool Hexahedron::addTetra()
3150 nodes[0] = _polygons[0]._links[0].FirstNode();
3151 nodes[1] = _polygons[0]._links[1].FirstNode();
3152 nodes[2] = _polygons[0]._links[2].FirstNode();
3154 _Link* link = _polygons[0]._links[0]._link;
3155 //ASSERT( link->_faces.size() > 1 );
3156 if ( !link->_faces[0] || !link->_faces[1] )
3157 return debugDumpLink( link );
3159 // a triangle sharing <link> with _polygons[0]
3160 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
3161 for ( int i = 0; i < 3; ++i )
3162 if ( tria->_links[i]._link == link )
3164 nodes[3] = tria->_links[(i+1)%3].LastNode();
3165 _volumeDefs.set( vector< _Node* >( nodes, nodes+4 ));
3171 //================================================================================
3173 * \brief Tries to create a pentahedron
3175 bool Hexahedron::addPenta()
3177 // find a base triangular face
3179 for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
3180 if ( _polygons[ iF ]._links.size() == 3 )
3182 if ( iTri < 0 ) return false;
3187 for ( int iL = 0; iL < 3; ++iL )
3190 nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
3193 // find a top node above the base node
3194 _Link* link = _polygons[ iTri ]._links[iL]._link;
3195 //ASSERT( link->_faces.size() > 1 );
3196 if ( !link->_faces[0] || !link->_faces[1] )
3197 return debugDumpLink( link );
3198 // a quadrangle sharing <link> with a base triangle
3199 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
3200 if ( quad->_links.size() != 4 ) return false;
3201 for ( int i = 0; i < 4; ++i )
3202 if ( quad->_links[i]._link == link )
3204 // 1st node of a link opposite to <link> in <quad>
3205 nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
3211 _volumeDefs.set( vector< _Node* >( nodes, nodes+6 ));
3213 return ( nbN == 6 );
3215 //================================================================================
3217 * \brief Tries to create a pyramid
3219 bool Hexahedron::addPyra()
3221 // find a base quadrangle
3223 for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
3224 if ( _polygons[ iF ]._links.size() == 4 )
3226 if ( iQuad < 0 ) return false;
3230 nodes[0] = _polygons[iQuad]._links[0].FirstNode();
3231 nodes[1] = _polygons[iQuad]._links[1].FirstNode();
3232 nodes[2] = _polygons[iQuad]._links[2].FirstNode();
3233 nodes[3] = _polygons[iQuad]._links[3].FirstNode();
3235 _Link* link = _polygons[iQuad]._links[0]._link;
3236 //ASSERT( link->_faces.size() > 1 );
3237 if ( !link->_faces[0] || !link->_faces[1] )
3238 return debugDumpLink( link );
3240 // a triangle sharing <link> with a base quadrangle
3241 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
3242 if ( tria->_links.size() != 3 ) return false;
3243 for ( int i = 0; i < 3; ++i )
3244 if ( tria->_links[i]._link == link )
3246 nodes[4] = tria->_links[(i+1)%3].LastNode();
3247 _volumeDefs.set( vector< _Node* >( nodes, nodes+5 ));
3253 //================================================================================
3255 * \brief Dump a link and return \c false
3257 bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
3260 gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
3261 cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
3262 << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
3263 << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
3268 //================================================================================
3270 * \brief computes exact bounding box with axes parallel to given ones
3272 //================================================================================
3274 void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
3275 const double* axesDirs,
3279 TopoDS_Compound allFacesComp;
3280 b.MakeCompound( allFacesComp );
3281 for ( size_t iF = 0; iF < faceVec.size(); ++iF )
3282 b.Add( allFacesComp, faceVec[ iF ] );
3284 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
3285 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
3287 for ( int i = 0; i < 6; ++i )
3288 farDist = Max( farDist, 10 * sP[i] );
3290 gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
3291 gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
3292 gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
3293 axis[0].Normalize();
3294 axis[1].Normalize();
3295 axis[2].Normalize();
3297 gp_Mat basis( axis[0], axis[1], axis[2] );
3298 gp_Mat bi = basis.Inverted();
3301 for ( int iDir = 0; iDir < 3; ++iDir )
3303 gp_XYZ axis0 = axis[ iDir ];
3304 gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
3305 gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
3306 for ( int isMax = 0; isMax < 2; ++isMax )
3308 double shift = isMax ? farDist : -farDist;
3309 gp_XYZ orig = shift * axis0;
3310 gp_XYZ norm = axis1 ^ axis2;
3311 gp_Pln pln( orig, norm );
3312 norm = pln.Axis().Direction().XYZ();
3313 BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
3315 gp_Pnt& pAxis = isMax ? pMax : pMin;
3316 gp_Pnt pPlane, pFaces;
3317 double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
3322 for ( int i = 0; i < 2; ++i ) {
3323 corner.SetCoord( 1, sP[ i*3 ]);
3324 for ( int j = 0; j < 2; ++j ) {
3325 corner.SetCoord( 2, sP[ i*3 + 1 ]);
3326 for ( int k = 0; k < 2; ++k )
3328 corner.SetCoord( 3, sP[ i*3 + 2 ]);
3334 corner = isMax ? bb.CornerMax() : bb.CornerMin();
3335 pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
3339 gp_XYZ pf = pFaces.XYZ() * bi;
3340 pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
3346 shapeBox.Add( pMin );
3347 shapeBox.Add( pMax );
3354 //=============================================================================
3356 * \brief Generates 3D structured Cartesian mesh in the internal part of
3357 * solid shapes and polyhedral volumes near the shape boundary.
3358 * \param theMesh - mesh to fill in
3359 * \param theShape - a compound of all SOLIDs to mesh
3360 * \retval bool - true in case of success
3362 //=============================================================================
3364 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh & theMesh,
3365 const TopoDS_Shape & theShape)
3367 // The algorithm generates the mesh in following steps:
3369 // 1) Intersection of grid lines with the geometry boundary.
3370 // This step allows to find out if a given node of the initial grid is
3371 // inside or outside the geometry.
3373 // 2) For each cell of the grid, check how many of it's nodes are outside
3374 // of the geometry boundary. Depending on a result of this check
3375 // - skip a cell, if all it's nodes are outside
3376 // - skip a cell, if it is too small according to the size threshold
3377 // - add a hexahedron in the mesh, if all nodes are inside
3378 // - add a polyhedron in the mesh, if some nodes are inside and some outside
3380 _computeCanceled = false;
3382 SMESH_MesherHelper helper( theMesh );
3387 grid._helper = &helper;
3389 vector< TopoDS_Shape > faceVec;
3391 TopTools_MapOfShape faceMap;
3392 TopExp_Explorer fExp;
3393 for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
3394 if ( !faceMap.Add( fExp.Current() ))
3395 faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
3397 for ( fExp.ReInit(); fExp.More(); fExp.Next() )
3398 if ( faceMap.Contains( fExp.Current() ))
3399 faceVec.push_back( fExp.Current() );
3401 vector<FaceGridIntersector> facesItersectors( faceVec.size() );
3402 map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
3403 TopExp_Explorer eExp;
3405 for ( int i = 0; i < faceVec.size(); ++i )
3407 facesItersectors[i]._face = TopoDS::Face ( faceVec[i] );
3408 facesItersectors[i]._faceID = grid._shapes.Add( faceVec[i] );
3409 facesItersectors[i]._grid = &grid;
3410 shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
3412 if ( _hyp->GetToAddEdges() )
3414 helper.SetSubShape( faceVec[i] );
3415 for ( eExp.Init( faceVec[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
3417 const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
3418 if ( !SMESH_Algo::isDegenerated( edge ) &&
3419 !helper.IsRealSeam( edge ))
3420 edge2faceIDsMap[ grid._shapes.Add( edge )].push_back( facesItersectors[i]._faceID );
3425 getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
3427 vector<double> xCoords, yCoords, zCoords;
3428 _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
3430 grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
3432 if ( _computeCanceled ) return false;
3435 { // copy partner faces and curves of not thread-safe types
3436 set< const Standard_Transient* > tshapes;
3437 BRepBuilderAPI_Copy copier;
3438 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3440 if ( !facesItersectors[i].IsThreadSafe(tshapes) )
3442 copier.Perform( facesItersectors[i]._face );
3443 facesItersectors[i]._face = TopoDS::Face( copier );
3447 // Intersection of grid lines with the geometry boundary.
3448 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
3449 ParallelIntersector( facesItersectors ),
3450 tbb::simple_partitioner());
3452 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3453 facesItersectors[i].Intersect();
3456 // put interesection points onto the GridLine's; this is done after intersection
3457 // to avoid contention of facesItersectors for writing into the same GridLine
3458 // in case of parallel work of facesItersectors
3459 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3460 facesItersectors[i].StoreIntersections();
3462 TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
3463 helper.SetSubShape( solidExp.Current() );
3464 helper.SetElementsOnShape( true );
3466 if ( _computeCanceled ) return false;
3468 // create nodes on the geometry
3469 grid.ComputeNodes(helper);
3471 if ( _computeCanceled ) return false;
3473 // create volume elements
3474 Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
3475 int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
3477 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
3480 // make all SOLIDs computed
3481 if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
3483 SMDS_ElemIteratorPtr volIt = sm1->GetElements();
3484 for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
3486 const SMDS_MeshElement* vol = volIt->next();
3487 sm1->RemoveElement( vol, /*isElemDeleted=*/false );
3488 meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
3491 // make other sub-shapes computed
3492 setSubmeshesComputed( theMesh, theShape );
3495 // remove free nodes
3496 if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
3498 TIDSortedNodeSet nodesToRemove;
3499 // get intersection nodes
3500 for ( int iDir = 0; iDir < 3; ++iDir )
3502 vector< GridLine >& lines = grid._lines[ iDir ];
3503 for ( size_t i = 0; i < lines.size(); ++i )
3505 multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
3506 for ( ; ip != lines[i]._intPoints.end(); ++ip )
3507 if ( ip->_node && ip->_node->NbInverseElements() == 0 )
3508 nodesToRemove.insert( nodesToRemove.end(), ip->_node );
3512 for ( size_t i = 0; i < grid._nodes.size(); ++i )
3513 if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
3514 nodesToRemove.insert( nodesToRemove.end(), grid._nodes[i] );
3517 TIDSortedNodeSet::iterator n = nodesToRemove.begin();
3518 for ( ; n != nodesToRemove.end(); ++n )
3519 meshDS->RemoveFreeNode( *n, smDS, /*fromGroups=*/false );
3525 // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
3526 catch ( SMESH_ComputeError& e)
3528 return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
3533 //=============================================================================
3537 //=============================================================================
3539 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh & theMesh,
3540 const TopoDS_Shape & theShape,
3541 MapShapeNbElems& theResMap)
3544 // std::vector<int> aResVec(SMDSEntity_Last);
3545 // for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
3546 // if(IsQuadratic) {
3547 // aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3548 // int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
3549 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
3552 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
3553 // aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3555 // SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
3556 // aResMap.insert(std::make_pair(sm,aResVec));
3561 //=============================================================================
3565 * \brief Event listener setting/unsetting _alwaysComputed flag to
3566 * submeshes of inferior levels to prevent their computing
3568 struct _EventListener : public SMESH_subMeshEventListener
3572 _EventListener(const string& algoName):
3573 SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
3576 // --------------------------------------------------------------------------------
3577 // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
3579 static void setAlwaysComputed( const bool isComputed,
3580 SMESH_subMesh* subMeshOfSolid)
3582 SMESH_subMeshIteratorPtr smIt =
3583 subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3584 while ( smIt->more() )
3586 SMESH_subMesh* sm = smIt->next();
3587 sm->SetIsAlwaysComputed( isComputed );
3589 subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
3592 // --------------------------------------------------------------------------------
3593 // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
3595 virtual void ProcessEvent(const int event,
3596 const int eventType,
3597 SMESH_subMesh* subMeshOfSolid,
3598 SMESH_subMeshEventListenerData* data,
3599 const SMESH_Hypothesis* hyp = 0)
3601 if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
3603 setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
3608 SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
3609 if ( !algo3D || _algoName != algo3D->GetName() )
3610 setAlwaysComputed( false, subMeshOfSolid );
3614 // --------------------------------------------------------------------------------
3615 // set the event listener
3617 static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
3619 subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
3624 }; // struct _EventListener
3628 //================================================================================
3630 * \brief Sets event listener to submeshes if necessary
3631 * \param subMesh - submesh where algo is set
3632 * This method is called when a submesh gets HYP_OK algo_state.
3633 * After being set, event listener is notified on each event of a submesh.
3635 //================================================================================
3637 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
3639 _EventListener::SetOn( subMesh, GetName() );
3642 //================================================================================
3644 * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
3646 //================================================================================
3648 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh& theMesh,
3649 const TopoDS_Shape& theShape)
3651 for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
3652 _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));