1 // Copyright (C) 2007-2016 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 "SMESHDS_Mesh.hxx"
29 #include "SMESH_Block.hxx"
30 #include "SMESH_Comment.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_MesherHelper.hxx"
33 #include "SMESH_subMesh.hxx"
34 #include "SMESH_subMeshEventListener.hxx"
35 #include "StdMeshers_CartesianParameters3D.hxx"
37 #include <utilities.h>
38 #include <Utils_ExceptHandlers.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 <BRepTopAdaptor_FClass2d.hxx>
49 #include <BRep_Builder.hxx>
50 #include <BRep_Tool.hxx>
51 #include <Bnd_B3d.hxx>
52 #include <Bnd_Box.hxx>
54 #include <GCPnts_UniformDeflection.hxx>
55 #include <Geom2d_BSplineCurve.hxx>
56 #include <Geom2d_BezierCurve.hxx>
57 #include <Geom2d_TrimmedCurve.hxx>
58 #include <GeomAPI_ProjectPointOnSurf.hxx>
59 #include <GeomLib.hxx>
60 #include <Geom_BSplineCurve.hxx>
61 #include <Geom_BSplineSurface.hxx>
62 #include <Geom_BezierCurve.hxx>
63 #include <Geom_BezierSurface.hxx>
64 #include <Geom_RectangularTrimmedSurface.hxx>
65 #include <Geom_TrimmedCurve.hxx>
66 #include <IntAna_IntConicQuad.hxx>
67 #include <IntAna_IntLinTorus.hxx>
68 #include <IntAna_Quadric.hxx>
69 #include <IntCurveSurface_TransitionOnCurve.hxx>
70 #include <IntCurvesFace_Intersector.hxx>
71 #include <Poly_Triangulation.hxx>
72 #include <Precision.hxx>
74 #include <TopExp_Explorer.hxx>
75 #include <TopLoc_Location.hxx>
76 #include <TopTools_IndexedMapOfShape.hxx>
77 #include <TopTools_MapOfShape.hxx>
79 #include <TopoDS_Compound.hxx>
80 #include <TopoDS_Face.hxx>
81 #include <TopoDS_TShape.hxx>
82 #include <gp_Cone.hxx>
83 #include <gp_Cylinder.hxx>
86 #include <gp_Pnt2d.hxx>
87 #include <gp_Sphere.hxx>
88 #include <gp_Torus.hxx>
94 #include <tbb/parallel_for.h>
95 //#include <tbb/enumerable_thread_specific.h>
104 //=============================================================================
108 //=============================================================================
110 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
111 :SMESH_3D_Algo(hypId, studyId, gen)
113 _name = "Cartesian_3D";
114 _shapeType = (1 << TopAbs_SOLID); // 1 bit /shape type
115 _compatibleHypothesis.push_back("CartesianParameters3D");
117 _onlyUnaryInput = false; // to mesh all SOLIDs at once
118 _requireDiscreteBoundary = false; // 2D mesh not needed
119 _supportSubmeshes = false; // do not use any existing mesh
122 //=============================================================================
124 * Check presence of a hypothesis
126 //=============================================================================
128 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh& aMesh,
129 const TopoDS_Shape& aShape,
130 Hypothesis_Status& aStatus)
132 aStatus = SMESH_Hypothesis::HYP_MISSING;
134 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
135 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
136 if ( h == hyps.end())
141 for ( ; h != hyps.end(); ++h )
143 if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
145 aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
150 return aStatus == HYP_OK;
157 //=============================================================================
158 // Definitions of internal utils
159 // --------------------------------------------------------------------------
161 Trans_TANGENT = IntCurveSurface_Tangent,
162 Trans_IN = IntCurveSurface_In,
163 Trans_OUT = IntCurveSurface_Out,
166 // --------------------------------------------------------------------------
168 * \brief Common data of any intersection between a Grid and a shape
170 struct B_IntersectPoint
172 mutable const SMDS_MeshNode* _node;
173 mutable vector< TGeomID > _faceIDs;
175 B_IntersectPoint(): _node(NULL) {}
176 void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
177 int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
178 bool IsOnFace( int faceID ) const;
179 virtual ~B_IntersectPoint() {}
181 // --------------------------------------------------------------------------
183 * \brief Data of intersection between a GridLine and a TopoDS_Face
185 struct F_IntersectPoint : public B_IntersectPoint
188 mutable Transition _transition;
189 mutable size_t _indexOnLine;
191 bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
193 // --------------------------------------------------------------------------
195 * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
197 struct E_IntersectPoint : public B_IntersectPoint
203 // --------------------------------------------------------------------------
205 * \brief A line of the grid and its intersections with 2D geometry
210 double _length; // line length
211 multiset< F_IntersectPoint > _intPoints;
213 void RemoveExcessIntPoints( const double tol );
214 bool GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut );
216 // --------------------------------------------------------------------------
218 * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
223 vector< gp_XYZ > _origins; // origin points of all planes in one direction
224 vector< double > _zProjs; // projections of origins to _zNorm
226 // --------------------------------------------------------------------------
228 * \brief Iterator on the parallel grid lines of one direction
234 size_t _iVar1, _iVar2, _iConst;
235 string _name1, _name2, _nameConst;
237 LineIndexer( size_t sz1, size_t sz2, size_t sz3,
238 size_t iv1, size_t iv2, size_t iConst,
239 const string& nv1, const string& nv2, const string& nConst )
241 _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
242 _curInd[0] = _curInd[1] = _curInd[2] = 0;
243 _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
244 _name1 = nv1; _name2 = nv2; _nameConst = nConst;
247 size_t I() const { return _curInd[0]; }
248 size_t J() const { return _curInd[1]; }
249 size_t K() const { return _curInd[2]; }
250 void SetIJK( size_t i, size_t j, size_t k )
252 _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
256 if ( ++_curInd[_iVar1] == _size[_iVar1] )
257 _curInd[_iVar1] = 0, ++_curInd[_iVar2];
259 bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
260 size_t LineIndex () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
261 size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
262 size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
263 size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
264 void SetIndexOnLine (size_t i) { _curInd[ _iConst ] = i; }
265 size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
267 // --------------------------------------------------------------------------
269 * \brief Container of GridLine's
273 vector< double > _coords[3]; // coordinates of grid nodes
274 gp_XYZ _axes [3]; // axis directions
275 vector< GridLine > _lines [3]; // in 3 directions
276 double _tol, _minCellSize;
278 gp_Mat _invB; // inverted basis of _axes
280 vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
281 vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
283 list< E_IntersectPoint > _edgeIntP; // intersections with EDGEs
284 TopTools_IndexedMapOfShape _shapes;
286 SMESH_MesherHelper* _helper;
288 size_t CellIndex( size_t i, size_t j, size_t k ) const
290 return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
292 size_t NodeIndex( size_t i, size_t j, size_t k ) const
294 return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
296 size_t NodeIndexDX() const { return 1; }
297 size_t NodeIndexDY() const { return _coords[0].size(); }
298 size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
300 LineIndexer GetLineIndexer(size_t iDir) const;
302 void SetCoordinates(const vector<double>& xCoords,
303 const vector<double>& yCoords,
304 const vector<double>& zCoords,
305 const double* axesDirs,
306 const Bnd_Box& bndBox );
307 void ComputeUVW(const gp_XYZ& p, double uvw[3]);
308 void ComputeNodes(SMESH_MesherHelper& helper);
310 // --------------------------------------------------------------------------
312 * \brief Intersector of TopoDS_Face with all GridLine's
314 struct FaceGridIntersector
320 IntCurvesFace_Intersector* _surfaceInt;
321 vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
323 FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
326 void StoreIntersections()
328 for ( size_t i = 0; i < _intersections.size(); ++i )
330 multiset< F_IntersectPoint >::iterator ip =
331 _intersections[i].first->_intPoints.insert( _intersections[i].second );
332 ip->_faceIDs.reserve( 1 );
333 ip->_faceIDs.push_back( _faceID );
336 const Bnd_Box& GetFaceBndBox()
338 GetCurveFaceIntersector();
341 IntCurvesFace_Intersector* GetCurveFaceIntersector()
345 _surfaceInt = new IntCurvesFace_Intersector( _face, Precision::PConfusion() );
346 _bndBox = _surfaceInt->Bounding();
347 if ( _bndBox.IsVoid() )
348 BRepBndLib::Add (_face, _bndBox);
352 bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
354 // --------------------------------------------------------------------------
356 * \brief Intersector of a surface with a GridLine
358 struct FaceLineIntersector
361 double _u, _v, _w; // params on the face and the line
362 Transition _transition; // transition of at intersection (see IntCurveSurface.cdl)
363 Transition _transIn, _transOut; // IN and OUT transitions depending of face orientation
366 gp_Cylinder _cylinder;
370 IntCurvesFace_Intersector* _surfaceInt;
372 vector< F_IntersectPoint > _intPoints;
374 void IntersectWithPlane (const GridLine& gridLine);
375 void IntersectWithCylinder(const GridLine& gridLine);
376 void IntersectWithCone (const GridLine& gridLine);
377 void IntersectWithSphere (const GridLine& gridLine);
378 void IntersectWithTorus (const GridLine& gridLine);
379 void IntersectWithSurface (const GridLine& gridLine);
381 bool UVIsOnFace() const;
382 void addIntPoint(const bool toClassify=true);
383 bool isParamOnLineOK( const double linLength )
385 return -_tol < _w && _w < linLength + _tol;
387 FaceLineIntersector():_surfaceInt(0) {}
388 ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
390 // --------------------------------------------------------------------------
392 * \brief Class representing topology of the hexahedron and creating a mesh
393 * volume basing on analysis of hexahedron intersection with geometry
397 // --------------------------------------------------------------------------------
400 // --------------------------------------------------------------------------------
401 struct _Node //!< node either at a hexahedron corner or at intersection
403 const SMDS_MeshNode* _node; // mesh node at hexahedron corner
404 const B_IntersectPoint* _intPoint;
405 const _Face* _usedInFace;
407 _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
408 :_node(n), _intPoint(ip), _usedInFace(0) {}
409 const SMDS_MeshNode* Node() const
410 { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
411 const E_IntersectPoint* EdgeIntPnt() const
412 { return static_cast< const E_IntersectPoint* >( _intPoint ); }
413 bool IsUsedInFace( const _Face* polygon = 0 )
415 return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
417 void Add( const E_IntersectPoint* ip )
422 else if ( !_intPoint->_node ) {
423 ip->Add( _intPoint->_faceIDs );
427 _intPoint->Add( ip->_faceIDs );
430 TGeomID IsLinked( const B_IntersectPoint* other,
431 TGeomID avoidFace=-1 ) const // returns id of a common face
433 return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
435 bool IsOnFace( TGeomID faceID ) const // returns true if faceID is found
437 return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
441 if ( const SMDS_MeshNode* n = Node() )
442 return SMESH_TNodeXYZ( n );
443 if ( const E_IntersectPoint* eip =
444 dynamic_cast< const E_IntersectPoint* >( _intPoint ))
446 return gp_Pnt( 1e100, 0, 0 );
448 TGeomID ShapeID() const
450 if ( const E_IntersectPoint* eip = dynamic_cast< const E_IntersectPoint* >( _intPoint ))
451 return eip->_shapeID;
455 // --------------------------------------------------------------------------------
456 struct _Link // link connecting two _Node's
459 _Face* _faces[2]; // polygons sharing a link
460 vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
461 vector< _Node* > _fIntNodes; // _Node's at _fIntPoints
462 vector< _Link > _splits;
463 _Link() { _faces[0] = 0; }
465 // --------------------------------------------------------------------------------
470 _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
471 void Reverse() { _reverse = !_reverse; }
472 int NbResultLinks() const { return _link->_splits.size(); }
473 _OrientedLink ResultLink(int i) const
475 return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
477 _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
478 _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
479 operator bool() const { return _link; }
480 vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
482 vector< TGeomID > faces;
483 const B_IntersectPoint *ip0, *ip1;
484 if (( ip0 = _link->_nodes[0]->_intPoint ) &&
485 ( ip1 = _link->_nodes[1]->_intPoint ))
487 for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
488 if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
489 !usedIDs.count( ip0->_faceIDs[i] ) )
490 faces.push_back( ip0->_faceIDs[i] );
494 bool HasEdgeNodes() const
496 return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
497 dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
501 return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
503 void AddFace( _Face* f )
505 if ( _link->_faces[0] )
507 _link->_faces[1] = f;
511 _link->_faces[0] = f;
512 _link->_faces[1] = 0;
515 void RemoveFace( _Face* f )
517 if ( !_link->_faces[0] ) return;
519 if ( _link->_faces[1] == f )
521 _link->_faces[1] = 0;
523 else if ( _link->_faces[0] == f )
525 _link->_faces[0] = 0;
526 if ( _link->_faces[1] )
528 _link->_faces[0] = _link->_faces[1];
529 _link->_faces[1] = 0;
534 // --------------------------------------------------------------------------------
537 vector< _OrientedLink > _links; // links on GridLine's
538 vector< _Link > _polyLinks; // links added to close a polygonal face
539 vector< _Node* > _eIntNodes; // nodes at intersection with EDGEs
540 bool IsPolyLink( const _OrientedLink& ol )
542 return _polyLinks.empty() ? false :
543 ( &_polyLinks[0] <= ol._link && ol._link <= &_polyLinks.back() );
545 void AddPolyLink(_Node* n0, _Node* n1, _Face* faceToFindEqual=0)
547 if ( faceToFindEqual && faceToFindEqual != this ) {
548 for ( size_t iL = 0; iL < faceToFindEqual->_polyLinks.size(); ++iL )
549 if ( faceToFindEqual->_polyLinks[iL]._nodes[0] == n1 &&
550 faceToFindEqual->_polyLinks[iL]._nodes[1] == n0 )
553 ( _OrientedLink( & faceToFindEqual->_polyLinks[iL], /*reverse=*/true ));
560 _polyLinks.push_back( l );
561 _links.push_back( _OrientedLink( &_polyLinks.back() ));
564 // --------------------------------------------------------------------------------
565 struct _volumeDef // holder of nodes of a volume mesh element
567 vector< _Node* > _nodes;
568 vector< int > _quantities;
569 typedef boost::shared_ptr<_volumeDef> Ptr;
570 void set( const vector< _Node* >& nodes,
571 const vector< int >& quant = vector< int >() )
572 { _nodes = nodes; _quantities = quant; }
573 void set( _Node** nodes, int nb )
574 { _nodes.assign( nodes, nodes + nb ); }
577 // topology of a hexahedron
580 _Link _hexLinks [12];
583 // faces resulted from hexahedron intersection
584 vector< _Face > _polygons;
586 // intresections with EDGEs
587 vector< const E_IntersectPoint* > _eIntPoints;
589 // additional nodes created at intersection points
590 vector< _Node > _intNodes;
592 // nodes inside the hexahedron (at VERTEXes)
593 vector< _Node* > _vIntNodes;
595 // computed volume elements
596 //vector< _volumeDef::Ptr > _volumeDefs;
597 _volumeDef _volumeDefs;
600 double _sizeThreshold, _sideLength[3];
601 int _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
602 int _origNodeInd; // index of _hexNodes[0] node within the _grid
606 Hexahedron(const double sizeThreshold, Grid* grid);
607 int MakeElements(SMESH_MesherHelper& helper,
608 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
609 void ComputeElements();
610 void Init() { init( _i, _j, _k ); }
613 Hexahedron(const Hexahedron& other );
614 void init( size_t i, size_t j, size_t k );
615 void init( size_t i );
616 void addEdges(SMESH_MesherHelper& helper,
617 vector< Hexahedron* >& intersectedHex,
618 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
619 gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
620 double proj, BRepAdaptor_Curve& curve,
621 const gp_XYZ& axis, const gp_XYZ& origin );
622 int getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
623 bool addIntersection( const E_IntersectPoint& ip,
624 vector< Hexahedron* >& hexes,
625 int ijk[], int dIJK[] );
626 bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
627 bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
628 bool findChainOnEdge( const vector< _OrientedLink >& splits,
629 const _OrientedLink& prevSplit,
630 const _OrientedLink& avoidSplit,
633 vector<_Node*>& chn);
634 int addElements(SMESH_MesherHelper& helper);
635 bool isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper ) const;
636 void sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID face);
637 bool isInHole() const;
638 bool checkPolyhedronSize() const;
643 bool debugDumpLink( _Link* link );
644 _Node* findEqualNode( vector< _Node* >& nodes,
645 const E_IntersectPoint* ip,
648 for ( size_t i = 0; i < nodes.size(); ++i )
649 if ( nodes[i]->EdgeIntPnt() == ip ||
650 nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
654 bool isImplementEdges() const { return !_grid->_edgeIntP.empty(); }
655 bool isOutParam(const double uvw[3]) const;
659 // --------------------------------------------------------------------------
661 * \brief Hexahedron computing volumes in one thread
663 struct ParallelHexahedron
665 vector< Hexahedron* >& _hexVec;
666 ParallelHexahedron( vector< Hexahedron* >& hv ): _hexVec(hv) {}
667 void operator() ( const tbb::blocked_range<size_t>& r ) const
669 for ( size_t i = r.begin(); i != r.end(); ++i )
670 if ( Hexahedron* hex = _hexVec[ i ] )
671 hex->ComputeElements();
674 // --------------------------------------------------------------------------
676 * \brief Structure intersecting certain nb of faces with GridLine's in one thread
678 struct ParallelIntersector
680 vector< FaceGridIntersector >& _faceVec;
681 ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
682 void operator() ( const tbb::blocked_range<size_t>& r ) const
684 for ( size_t i = r.begin(); i != r.end(); ++i )
685 _faceVec[i].Intersect();
690 //=============================================================================
691 // Implementation of internal utils
692 //=============================================================================
694 * \brief adjust \a i to have \a val between values[i] and values[i+1]
696 inline void locateValue( int & i, double val, const vector<double>& values,
697 int& di, double tol )
699 //val += values[0]; // input \a val is measured from 0.
700 if ( i > (int) values.size()-2 )
703 while ( i+2 < (int) values.size() && val > values[ i+1 ])
705 while ( i > 0 && val < values[ i ])
708 if ( i > 0 && val - values[ i ] < tol )
710 else if ( i+2 < (int) values.size() && values[ i+1 ] - val < tol )
715 //=============================================================================
717 * Remove coincident intersection points
719 void GridLine::RemoveExcessIntPoints( const double tol )
721 if ( _intPoints.size() < 2 ) return;
723 set< Transition > tranSet;
724 multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
725 while ( ip2 != _intPoints.end() )
729 while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
731 tranSet.insert( ip1->_transition );
732 tranSet.insert( ip2->_transition );
733 ip2->Add( ip1->_faceIDs );
734 _intPoints.erase( ip1 );
737 if ( tranSet.size() > 1 ) // points with different transition coincide
739 bool isIN = tranSet.count( Trans_IN );
740 bool isOUT = tranSet.count( Trans_OUT );
742 (*ip1)._transition = Trans_TANGENT;
744 (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
748 //================================================================================
750 * Return "is OUT" state for nodes before the given intersection point
752 bool GridLine::GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut )
754 if ( ip->_transition == Trans_IN )
756 if ( ip->_transition == Trans_OUT )
758 if ( ip->_transition == Trans_APEX )
760 // singularity point (apex of a cone)
761 if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
763 multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
764 if ( ipAft == _intPoints.end() )
767 if ( ipBef->_transition != ipAft->_transition )
768 return ( ipBef->_transition == Trans_OUT );
769 return ( ipBef->_transition != Trans_OUT );
771 // _transition == Trans_TANGENT
774 //================================================================================
778 void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
779 const SMDS_MeshNode* n) const
781 if ( _faceIDs.empty() )
784 for ( size_t i = 0; i < fIDs.size(); ++i )
786 vector< TGeomID >::iterator it =
787 std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
788 if ( it == _faceIDs.end() )
789 _faceIDs.push_back( fIDs[i] );
794 //================================================================================
796 * Returns index of a common face if any, else zero
798 int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
801 for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
802 if ( avoidFace != other->_faceIDs[i] &&
803 IsOnFace ( other->_faceIDs[i] ))
804 return other->_faceIDs[i];
807 //================================================================================
809 * Returns \c true if \a faceID in in this->_faceIDs
811 bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
813 vector< TGeomID >::const_iterator it =
814 std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
815 return ( it != _faceIDs.end() );
817 //================================================================================
819 * Return an iterator on GridLine's in a given direction
821 LineIndexer Grid::GetLineIndexer(size_t iDir) const
823 const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
824 const string s [] = { "X", "Y", "Z" };
825 LineIndexer li( _coords[0].size(), _coords[1].size(), _coords[2].size(),
826 indices[iDir*3], indices[iDir*3+1], indices[iDir*3+2],
827 s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
830 //=============================================================================
832 * Creates GridLine's of the grid
834 void Grid::SetCoordinates(const vector<double>& xCoords,
835 const vector<double>& yCoords,
836 const vector<double>& zCoords,
837 const double* axesDirs,
838 const Bnd_Box& shapeBox)
840 _coords[0] = xCoords;
841 _coords[1] = yCoords;
842 _coords[2] = zCoords;
844 _axes[0].SetCoord( axesDirs[0],
847 _axes[1].SetCoord( axesDirs[3],
850 _axes[2].SetCoord( axesDirs[6],
853 _axes[0].Normalize();
854 _axes[1].Normalize();
855 _axes[2].Normalize();
857 _invB.SetCols( _axes[0], _axes[1], _axes[2] );
861 _minCellSize = Precision::Infinite();
862 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
864 for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
866 double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
867 if ( cellLen < _minCellSize )
868 _minCellSize = cellLen;
871 if ( _minCellSize < Precision::Confusion() )
872 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
873 SMESH_Comment("Too small cell size: ") << _minCellSize );
874 _tol = _minCellSize / 1000.;
876 // attune grid extremities to shape bounding box
878 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
879 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
880 double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
881 &_coords[0].back(), &_coords[1].back(), &_coords[2].back() };
882 for ( int i = 0; i < 6; ++i )
883 if ( fabs( sP[i] - *cP[i] ) < _tol )
884 *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
886 for ( int iDir = 0; iDir < 3; ++iDir )
888 if ( _coords[iDir][0] - sP[iDir] > _tol )
890 _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
891 _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
893 if ( sP[iDir+3] - _coords[iDir].back() > _tol )
895 _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
896 _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
899 _tol = _minCellSize / 1000.;
901 _origin = ( _coords[0][0] * _axes[0] +
902 _coords[1][0] * _axes[1] +
903 _coords[2][0] * _axes[2] );
906 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
908 LineIndexer li = GetLineIndexer( iDir );
909 _lines[iDir].resize( li.NbLines() );
910 double len = _coords[ iDir ].back() - _coords[iDir].front();
911 for ( ; li.More(); ++li )
913 GridLine& gl = _lines[iDir][ li.LineIndex() ];
914 gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
915 _coords[1][li.J()] * _axes[1] +
916 _coords[2][li.K()] * _axes[2] );
917 gl._line.SetDirection( _axes[ iDir ]);
922 //================================================================================
924 * Computes coordinates of a point in the grid CS
926 void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
928 gp_XYZ p = P * _invB;
929 p.Coord( UVW[0], UVW[1], UVW[2] );
931 //================================================================================
935 void Grid::ComputeNodes(SMESH_MesherHelper& helper)
937 // state of each node of the grid relative to the geometry
938 const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
939 vector< bool > isNodeOut( nbGridNodes, false );
940 _nodes.resize( nbGridNodes, 0 );
941 _gridIntP.resize( nbGridNodes, NULL );
943 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
945 LineIndexer li = GetLineIndexer( iDir );
947 // find out a shift of node index while walking along a GridLine in this direction
948 li.SetIndexOnLine( 0 );
949 size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
950 li.SetIndexOnLine( 1 );
951 const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
953 const vector<double> & coords = _coords[ iDir ];
954 for ( ; li.More(); ++li ) // loop on lines in iDir
956 li.SetIndexOnLine( 0 );
957 nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
959 GridLine& line = _lines[ iDir ][ li.LineIndex() ];
960 const gp_XYZ lineLoc = line._line.Location().XYZ();
961 const gp_XYZ lineDir = line._line.Direction().XYZ();
962 line.RemoveExcessIntPoints( _tol );
963 multiset< F_IntersectPoint >& intPnts = line._intPoints;
964 multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
967 const double* nodeCoord = & coords[0];
968 const double* coord0 = nodeCoord;
969 const double* coordEnd = coord0 + coords.size();
970 double nodeParam = 0;
971 for ( ; ip != intPnts.end(); ++ip )
973 // set OUT state or just skip IN nodes before ip
974 if ( nodeParam < ip->_paramOnLine - _tol )
976 isOut = line.GetIsOutBefore( ip, isOut );
978 while ( nodeParam < ip->_paramOnLine - _tol )
981 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
982 if ( ++nodeCoord < coordEnd )
983 nodeParam = *nodeCoord - *coord0;
987 if ( nodeCoord == coordEnd ) break;
989 // create a mesh node on a GridLine at ip if it does not coincide with a grid node
990 if ( nodeParam > ip->_paramOnLine + _tol )
992 // li.SetIndexOnLine( 0 );
993 // double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
994 // xyz[ li._iConst ] += ip->_paramOnLine;
995 gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
996 ip->_node = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
997 ip->_indexOnLine = nodeCoord-coord0-1;
999 // create a mesh node at ip concident with a grid node
1002 int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1003 if ( !_nodes[ nodeIndex ] )
1005 //li.SetIndexOnLine( nodeCoord-coord0 );
1006 //double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1007 gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1008 _nodes [ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1009 _gridIntP[ nodeIndex ] = & * ip;
1011 if ( _gridIntP[ nodeIndex ] )
1012 _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1014 _gridIntP[ nodeIndex ] = & * ip;
1015 // ip->_node = _nodes[ nodeIndex ]; -- to differ from ip on links
1016 ip->_indexOnLine = nodeCoord-coord0;
1017 if ( ++nodeCoord < coordEnd )
1018 nodeParam = *nodeCoord - *coord0;
1021 // set OUT state to nodes after the last ip
1022 for ( ; nodeCoord < coordEnd; ++nodeCoord )
1023 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
1027 // Create mesh nodes at !OUT nodes of the grid
1029 for ( size_t z = 0; z < _coords[2].size(); ++z )
1030 for ( size_t y = 0; y < _coords[1].size(); ++y )
1031 for ( size_t x = 0; x < _coords[0].size(); ++x )
1033 size_t nodeIndex = NodeIndex( x, y, z );
1034 if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
1036 //_nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
1037 gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1038 _coords[1][y] * _axes[1] +
1039 _coords[2][z] * _axes[2] );
1040 _nodes[ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1045 // check validity of transitions
1046 const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1047 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1049 LineIndexer li = GetLineIndexer( iDir );
1050 for ( ; li.More(); ++li )
1052 multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1053 if ( intPnts.empty() ) continue;
1054 if ( intPnts.size() == 1 )
1056 if ( intPnts.begin()->_transition != Trans_TANGENT &&
1057 intPnts.begin()->_transition != Trans_APEX )
1058 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1059 SMESH_Comment("Wrong SOLE transition of GridLine (")
1060 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1061 << ") along " << li._nameConst
1062 << ": " << trName[ intPnts.begin()->_transition] );
1066 if ( intPnts.begin()->_transition == Trans_OUT )
1067 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1068 SMESH_Comment("Wrong START transition of GridLine (")
1069 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1070 << ") along " << li._nameConst
1071 << ": " << trName[ intPnts.begin()->_transition ]);
1072 if ( intPnts.rbegin()->_transition == Trans_IN )
1073 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1074 SMESH_Comment("Wrong END transition of GridLine (")
1075 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1076 << ") along " << li._nameConst
1077 << ": " << trName[ intPnts.rbegin()->_transition ]);
1084 //=============================================================================
1086 * Intersects TopoDS_Face with all GridLine's
1088 void FaceGridIntersector::Intersect()
1090 FaceLineIntersector intersector;
1091 intersector._surfaceInt = GetCurveFaceIntersector();
1092 intersector._tol = _grid->_tol;
1093 intersector._transOut = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1094 intersector._transIn = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1096 typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1097 PIntFun interFunction;
1099 bool isDirect = true;
1100 BRepAdaptor_Surface surf( _face );
1101 switch ( surf.GetType() ) {
1103 intersector._plane = surf.Plane();
1104 interFunction = &FaceLineIntersector::IntersectWithPlane;
1105 isDirect = intersector._plane.Direct();
1107 case GeomAbs_Cylinder:
1108 intersector._cylinder = surf.Cylinder();
1109 interFunction = &FaceLineIntersector::IntersectWithCylinder;
1110 isDirect = intersector._cylinder.Direct();
1113 intersector._cone = surf.Cone();
1114 interFunction = &FaceLineIntersector::IntersectWithCone;
1115 //isDirect = intersector._cone.Direct();
1117 case GeomAbs_Sphere:
1118 intersector._sphere = surf.Sphere();
1119 interFunction = &FaceLineIntersector::IntersectWithSphere;
1120 isDirect = intersector._sphere.Direct();
1123 intersector._torus = surf.Torus();
1124 interFunction = &FaceLineIntersector::IntersectWithTorus;
1125 //isDirect = intersector._torus.Direct();
1128 interFunction = &FaceLineIntersector::IntersectWithSurface;
1131 std::swap( intersector._transOut, intersector._transIn );
1133 _intersections.clear();
1134 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1136 if ( surf.GetType() == GeomAbs_Plane )
1138 // check if all lines in this direction are parallel to a plane
1139 if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1140 Precision::Angular()))
1142 // find out a transition, that is the same for all lines of a direction
1143 gp_Dir plnNorm = intersector._plane.Axis().Direction();
1144 gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1145 intersector._transition =
1146 ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1148 if ( surf.GetType() == GeomAbs_Cylinder )
1150 // check if all lines in this direction are parallel to a cylinder
1151 if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1152 Precision::Angular()))
1156 // intersect the grid lines with the face
1157 for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1159 GridLine& gridLine = _grid->_lines[iDir][iL];
1160 if ( _bndBox.IsOut( gridLine._line )) continue;
1162 intersector._intPoints.clear();
1163 (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1164 for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1165 _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1169 //================================================================================
1171 * Return true if (_u,_v) is on the face
1173 bool FaceLineIntersector::UVIsOnFace() const
1175 TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1176 return ( state == TopAbs_IN || state == TopAbs_ON );
1178 //================================================================================
1180 * Store an intersection if it is IN or ON the face
1182 void FaceLineIntersector::addIntPoint(const bool toClassify)
1184 if ( !toClassify || UVIsOnFace() )
1187 p._paramOnLine = _w;
1188 p._transition = _transition;
1189 _intPoints.push_back( p );
1192 //================================================================================
1194 * Intersect a line with a plane
1196 void FaceLineIntersector::IntersectWithPlane(const GridLine& gridLine)
1198 IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1199 _w = linPlane.ParamOnConic(1);
1200 if ( isParamOnLineOK( gridLine._length ))
1202 ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1206 //================================================================================
1208 * Intersect a line with a cylinder
1210 void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1212 IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1213 if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1215 _w = linCylinder.ParamOnConic(1);
1216 if ( linCylinder.NbPoints() == 1 )
1217 _transition = Trans_TANGENT;
1219 _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1220 if ( isParamOnLineOK( gridLine._length ))
1222 ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1225 if ( linCylinder.NbPoints() > 1 )
1227 _w = linCylinder.ParamOnConic(2);
1228 if ( isParamOnLineOK( gridLine._length ))
1230 ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1231 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1237 //================================================================================
1239 * Intersect a line with a cone
1241 void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1243 IntAna_IntConicQuad linCone(gridLine._line,_cone);
1244 if ( !linCone.IsDone() ) return;
1246 gp_Vec du, dv, norm;
1247 for ( int i = 1; i <= linCone.NbPoints(); ++i )
1249 _w = linCone.ParamOnConic( i );
1250 if ( !isParamOnLineOK( gridLine._length )) continue;
1251 ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1254 ElSLib::D1( _u, _v, _cone, P, du, dv );
1256 double normSize2 = norm.SquareMagnitude();
1257 if ( normSize2 > Precision::Angular() * Precision::Angular() )
1259 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1260 cos /= sqrt( normSize2 );
1261 if ( cos < -Precision::Angular() )
1262 _transition = _transIn;
1263 else if ( cos > Precision::Angular() )
1264 _transition = _transOut;
1266 _transition = Trans_TANGENT;
1270 _transition = Trans_APEX;
1272 addIntPoint( /*toClassify=*/false);
1276 //================================================================================
1278 * Intersect a line with a sphere
1280 void FaceLineIntersector::IntersectWithSphere (const GridLine& gridLine)
1282 IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1283 if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1285 _w = linSphere.ParamOnConic(1);
1286 if ( linSphere.NbPoints() == 1 )
1287 _transition = Trans_TANGENT;
1289 _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1290 if ( isParamOnLineOK( gridLine._length ))
1292 ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1295 if ( linSphere.NbPoints() > 1 )
1297 _w = linSphere.ParamOnConic(2);
1298 if ( isParamOnLineOK( gridLine._length ))
1300 ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1301 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1307 //================================================================================
1309 * Intersect a line with a torus
1311 void FaceLineIntersector::IntersectWithTorus (const GridLine& gridLine)
1313 IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1314 if ( !linTorus.IsDone()) return;
1316 gp_Vec du, dv, norm;
1317 for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1319 _w = linTorus.ParamOnLine( i );
1320 if ( !isParamOnLineOK( gridLine._length )) continue;
1321 linTorus.ParamOnTorus( i, _u,_v );
1324 ElSLib::D1( _u, _v, _torus, P, du, dv );
1326 double normSize = norm.Magnitude();
1327 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1329 if ( cos < -Precision::Angular() )
1330 _transition = _transIn;
1331 else if ( cos > Precision::Angular() )
1332 _transition = _transOut;
1334 _transition = Trans_TANGENT;
1335 addIntPoint( /*toClassify=*/false);
1339 //================================================================================
1341 * Intersect a line with a non-analytical surface
1343 void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1345 _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1346 if ( !_surfaceInt->IsDone() ) return;
1347 for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1349 _transition = Transition( _surfaceInt->Transition( i ) );
1350 _w = _surfaceInt->WParameter( i );
1351 addIntPoint(/*toClassify=*/false);
1354 //================================================================================
1356 * check if its face can be safely intersected in a thread
1358 bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1363 TopLoc_Location loc;
1364 Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1365 Handle(Geom_RectangularTrimmedSurface) ts =
1366 Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1367 while( !ts.IsNull() ) {
1368 surf = ts->BasisSurface();
1369 ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1371 if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1372 surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1373 if ( !noSafeTShapes.insert( _face.TShape().get() ).second )
1377 TopExp_Explorer exp( _face, TopAbs_EDGE );
1378 for ( ; exp.More(); exp.Next() )
1380 bool edgeIsSafe = true;
1381 const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1384 Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1387 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1388 while( !tc.IsNull() ) {
1389 c = tc->BasisCurve();
1390 tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1392 if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1393 c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1400 Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1403 Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1404 while( !tc.IsNull() ) {
1405 c2 = tc->BasisCurve();
1406 tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1408 if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1409 c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1413 if ( !edgeIsSafe && !noSafeTShapes.insert( e.TShape().get() ).second )
1418 //================================================================================
1420 * \brief Creates topology of the hexahedron
1422 Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1423 : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbFaceIntNodes(0)
1425 _polygons.reserve(100); // to avoid reallocation;
1427 //set nodes shift within grid->_nodes from the node 000
1428 size_t dx = _grid->NodeIndexDX();
1429 size_t dy = _grid->NodeIndexDY();
1430 size_t dz = _grid->NodeIndexDZ();
1432 size_t i100 = i000 + dx;
1433 size_t i010 = i000 + dy;
1434 size_t i110 = i010 + dx;
1435 size_t i001 = i000 + dz;
1436 size_t i101 = i100 + dz;
1437 size_t i011 = i010 + dz;
1438 size_t i111 = i110 + dz;
1439 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1440 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1441 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1442 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1443 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1444 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1445 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1446 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1448 vector< int > idVec;
1449 // set nodes to links
1450 for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1452 SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1453 _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1454 link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1455 link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1458 // set links to faces
1459 int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1460 for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1462 SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1463 _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1464 bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1465 faceID == SMESH_Block::ID_Fx1z ||
1466 faceID == SMESH_Block::ID_F0yz );
1467 quad._links.resize(4);
1468 vector<_OrientedLink>::iterator frwLinkIt = quad._links.begin();
1469 vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1470 for ( int i = 0; i < 4; ++i )
1472 bool revLink = revFace;
1473 if ( i > 1 ) // reverse links u1 and v0
1475 _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1476 link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1481 //================================================================================
1483 * \brief Copy constructor
1485 Hexahedron::Hexahedron( const Hexahedron& other )
1486 :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbFaceIntNodes(0)
1488 _polygons.reserve(100); // to avoid reallocation;
1490 for ( int i = 0; i < 8; ++i )
1491 _nodeShift[i] = other._nodeShift[i];
1493 for ( int i = 0; i < 12; ++i )
1495 const _Link& srcLink = other._hexLinks[ i ];
1496 _Link& tgtLink = this->_hexLinks[ i ];
1497 tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1498 tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1501 for ( int i = 0; i < 6; ++i )
1503 const _Face& srcQuad = other._hexQuads[ i ];
1504 _Face& tgtQuad = this->_hexQuads[ i ];
1505 tgtQuad._links.resize(4);
1506 for ( int j = 0; j < 4; ++j )
1508 const _OrientedLink& srcLink = srcQuad._links[ j ];
1509 _OrientedLink& tgtLink = tgtQuad._links[ j ];
1510 tgtLink._reverse = srcLink._reverse;
1511 tgtLink._link = _hexLinks + ( srcLink._link - other._hexLinks );
1516 //================================================================================
1518 * \brief Initializes its data by given grid cell
1520 void Hexahedron::init( size_t i, size_t j, size_t k )
1522 _i = i; _j = j; _k = k;
1523 // set nodes of grid to nodes of the hexahedron and
1524 // count nodes at hexahedron corners located IN and ON geometry
1525 _nbCornerNodes = _nbBndNodes = 0;
1526 _origNodeInd = _grid->NodeIndex( i,j,k );
1527 for ( int iN = 0; iN < 8; ++iN )
1529 _hexNodes[iN]._node = _grid->_nodes [ _origNodeInd + _nodeShift[iN] ];
1530 _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
1531 _nbCornerNodes += bool( _hexNodes[iN]._node );
1532 _nbBndNodes += bool( _hexNodes[iN]._intPoint );
1534 _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1535 _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1536 _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1541 if ( _nbFaceIntNodes + _eIntPoints.size() > 0 &&
1542 _nbFaceIntNodes + _nbCornerNodes + _eIntPoints.size() > 3)
1544 _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
1546 // this method can be called in parallel, so use own helper
1547 SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
1549 // create sub-links (_splits) by splitting links with _fIntPoints
1551 for ( int iLink = 0; iLink < 12; ++iLink )
1553 _Link& link = _hexLinks[ iLink ];
1554 link._fIntNodes.resize( link._fIntPoints.size() );
1555 for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
1557 _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
1558 link._fIntNodes[ i ] = & _intNodes.back();
1561 link._splits.clear();
1562 split._nodes[ 0 ] = link._nodes[0];
1563 bool isOut = ( ! link._nodes[0]->Node() );
1564 bool checkTransition;
1565 for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
1567 const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
1568 if ( !isGridNode ) // intersection non-coincident with a grid node
1570 if ( split._nodes[ 0 ]->Node() && !isOut )
1572 split._nodes[ 1 ] = link._fIntNodes[i];
1573 link._splits.push_back( split );
1575 split._nodes[ 0 ] = link._fIntNodes[i];
1576 checkTransition = true;
1578 else // FACE intersection coincident with a grid node (at link ends)
1580 checkTransition = ( i == 0 && link._nodes[0]->Node() );
1582 if ( checkTransition )
1584 if ( link._fIntPoints[i]->_faceIDs.size() > 1 || _eIntPoints.size() > 0 )
1585 isOut = isOutPoint( link, i, helper );
1587 switch ( link._fIntPoints[i]->_transition ) {
1588 case Trans_OUT: isOut = true; break;
1589 case Trans_IN : isOut = false; break;
1591 isOut = isOutPoint( link, i, helper );
1595 if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
1597 split._nodes[ 1 ] = link._nodes[1];
1598 link._splits.push_back( split );
1602 // Create _Node's at intersections with EDGEs.
1604 const double tol2 = _grid->_tol * _grid->_tol;
1605 int facets[3], nbFacets, subEntity;
1607 for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
1609 nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
1610 _Node* equalNode = 0;
1611 switch( nbFacets ) {
1612 case 1: // in a _Face
1614 _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
1615 equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1617 equalNode->Add( _eIntPoints[ iP ] );
1620 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1621 quad._eIntNodes.push_back( & _intNodes.back() );
1625 case 2: // on a _Link
1627 _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
1628 if ( link._splits.size() > 0 )
1630 equalNode = findEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
1632 equalNode->Add( _eIntPoints[ iP ] );
1636 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1637 for ( int iF = 0; iF < 2; ++iF )
1639 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1640 equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1642 equalNode->Add( _eIntPoints[ iP ] );
1645 quad._eIntNodes.push_back( & _intNodes.back() );
1651 case 3: // at a corner
1653 _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
1654 if ( node.Node() > 0 )
1656 if ( node._intPoint )
1657 node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
1661 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1662 for ( int iF = 0; iF < 3; ++iF )
1664 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1665 equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1667 equalNode->Add( _eIntPoints[ iP ] );
1670 quad._eIntNodes.push_back( & _intNodes.back() );
1676 } // switch( nbFacets )
1678 if ( nbFacets == 0 ||
1679 _grid->_shapes( _eIntPoints[ iP ]->_shapeID ).ShapeType() == TopAbs_VERTEX )
1681 equalNode = findEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
1683 equalNode->Add( _eIntPoints[ iP ] );
1685 else if ( nbFacets == 0 ) {
1686 if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
1687 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1688 _vIntNodes.push_back( & _intNodes.back() );
1691 } // loop on _eIntPoints
1693 else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
1696 // create sub-links (_splits) of whole links
1697 for ( int iLink = 0; iLink < 12; ++iLink )
1699 _Link& link = _hexLinks[ iLink ];
1700 link._splits.clear();
1701 if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
1703 split._nodes[ 0 ] = link._nodes[0];
1704 split._nodes[ 1 ] = link._nodes[1];
1705 link._splits.push_back( split );
1711 //================================================================================
1713 * \brief Initializes its data by given grid cell (countered from zero)
1715 void Hexahedron::init( size_t iCell )
1717 size_t iNbCell = _grid->_coords[0].size() - 1;
1718 size_t jNbCell = _grid->_coords[1].size() - 1;
1719 _i = iCell % iNbCell;
1720 _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1721 _k = iCell / iNbCell / jNbCell;
1725 //================================================================================
1727 * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1729 void Hexahedron::ComputeElements()
1733 int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
1734 if ( _nbCornerNodes + nbIntersections < 4 )
1737 if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
1741 _polygons.reserve( 20 );
1743 // Create polygons from quadrangles
1744 // --------------------------------
1746 vector< _OrientedLink > splits;
1747 vector<_Node*> chainNodes;
1748 _Face* coplanarPolyg;
1750 bool hasEdgeIntersections = !_eIntPoints.empty();
1752 for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1754 _Face& quad = _hexQuads[ iF ] ;
1756 _polygons.resize( _polygons.size() + 1 );
1757 _Face* polygon = &_polygons.back();
1758 polygon->_polyLinks.reserve( 20 );
1761 for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1762 for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
1763 splits.push_back( quad._links[ iE ].ResultLink( iS ));
1765 // add splits of links to a polygon and add _polyLinks to make
1766 // polygon's boundary closed
1768 int nbSplits = splits.size();
1769 if (( nbSplits == 1 ) &&
1770 ( quad._eIntNodes.empty() ||
1771 splits[0].FirstNode()->IsLinked( splits[0].LastNode()->_intPoint )))
1772 //( quad._eIntNodes.empty() || _nbCornerNodes + nbIntersections > 6 ))
1776 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1777 if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
1778 quad._eIntNodes[ iP ]->_usedInFace = 0;
1780 size_t nbUsedEdgeNodes = 0;
1781 _Face* prevPolyg = 0; // polygon previously created from this quad
1783 while ( nbSplits > 0 )
1786 while ( !splits[ iS ] )
1789 if ( !polygon->_links.empty() )
1791 _polygons.resize( _polygons.size() + 1 );
1792 polygon = &_polygons.back();
1793 polygon->_polyLinks.reserve( 20 );
1795 polygon->_links.push_back( splits[ iS ] );
1796 splits[ iS++ ]._link = 0;
1799 _Node* nFirst = polygon->_links.back().FirstNode();
1800 _Node *n1,*n2 = polygon->_links.back().LastNode();
1801 for ( ; nFirst != n2 && iS < splits.size(); ++iS )
1803 _OrientedLink& split = splits[ iS ];
1804 if ( !split ) continue;
1806 n1 = split.FirstNode();
1809 n1->_intPoint->_faceIDs.size() > 1 )
1811 // n1 is at intersection with EDGE
1812 if ( findChainOnEdge( splits, polygon->_links.back(), split, iS, quad, chainNodes ))
1814 for ( size_t i = 1; i < chainNodes.size(); ++i )
1815 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
1816 prevPolyg = polygon;
1817 n2 = chainNodes.back();
1821 else if ( n1 != n2 )
1823 // try to connect to intersections with EDGEs
1824 if ( quad._eIntNodes.size() > nbUsedEdgeNodes &&
1825 findChain( n2, n1, quad, chainNodes ))
1827 for ( size_t i = 1; i < chainNodes.size(); ++i )
1829 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i] );
1830 nbUsedEdgeNodes += ( chainNodes[i]->IsUsedInFace( polygon ));
1832 if ( chainNodes.back() != n1 )
1834 n2 = chainNodes.back();
1839 // try to connect to a split ending on the same FACE
1842 _OrientedLink foundSplit;
1843 for ( size_t i = iS; i < splits.size() && !foundSplit; ++i )
1844 if (( foundSplit = splits[ i ]) &&
1845 ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
1851 foundSplit._link = 0;
1855 if ( n2 != foundSplit.FirstNode() )
1857 polygon->AddPolyLink( n2, foundSplit.FirstNode() );
1858 n2 = foundSplit.FirstNode();
1864 if ( n2->IsLinked( nFirst->_intPoint ))
1866 polygon->AddPolyLink( n2, n1, prevPolyg );
1869 } // if ( n1 != n2 )
1871 polygon->_links.push_back( split );
1874 n2 = polygon->_links.back().LastNode();
1878 if ( nFirst != n2 ) // close a polygon
1880 if ( !findChain( n2, nFirst, quad, chainNodes ))
1882 if ( !closePolygon( polygon, chainNodes ))
1883 if ( !isImplementEdges() )
1884 chainNodes.push_back( nFirst );
1886 for ( size_t i = 1; i < chainNodes.size(); ++i )
1888 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
1889 nbUsedEdgeNodes += bool( chainNodes[i]->IsUsedInFace( polygon ));
1893 if ( polygon->_links.size() < 3 && nbSplits > 0 )
1895 polygon->_polyLinks.clear();
1896 polygon->_links.clear();
1898 } // while ( nbSplits > 0 )
1900 if ( polygon->_links.size() < 3 )
1902 _polygons.pop_back();
1904 } // loop on 6 hexahedron sides
1906 // Create polygons closing holes in a polyhedron
1907 // ----------------------------------------------
1909 // clear _usedInFace
1910 for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
1911 _intNodes[ iN ]._usedInFace = 0;
1913 // add polygons to their links and mark used nodes
1914 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1916 _Face& polygon = _polygons[ iP ];
1917 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1919 polygon._links[ iL ].AddFace( &polygon );
1920 polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
1924 vector< _OrientedLink* > freeLinks;
1925 freeLinks.reserve(20);
1926 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1928 _Face& polygon = _polygons[ iP ];
1929 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1930 if ( polygon._links[ iL ].NbFaces() < 2 )
1931 freeLinks.push_back( & polygon._links[ iL ]);
1933 int nbFreeLinks = freeLinks.size();
1934 if ( nbFreeLinks == 1 ) return;
1936 // put not used intersection nodes to _vIntNodes
1937 int nbVertexNodes = 0; // nb not used vertex nodes
1939 for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
1940 nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
1942 const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
1943 for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
1945 if ( _intNodes[ iN ].IsUsedInFace() ) continue;
1946 if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
1948 findEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
1951 _vIntNodes.push_back( &_intNodes[ iN ]);
1957 set<TGeomID> usedFaceIDs;
1958 vector< TGeomID > faces;
1959 TGeomID curFace = 0;
1960 const size_t nbQuadPolygons = _polygons.size();
1961 E_IntersectPoint ipTmp;
1963 // create polygons by making closed chains of free links
1964 size_t iPolygon = _polygons.size();
1965 while ( nbFreeLinks > 0 )
1967 if ( iPolygon == _polygons.size() )
1969 _polygons.resize( _polygons.size() + 1 );
1970 _polygons[ iPolygon ]._polyLinks.reserve( 20 );
1971 _polygons[ iPolygon ]._links.reserve( 20 );
1973 _Face& polygon = _polygons[ iPolygon ];
1975 _OrientedLink* curLink = 0;
1977 if (( !hasEdgeIntersections ) ||
1978 ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
1980 // get a remaining link to start from
1981 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1982 if (( curLink = freeLinks[ iL ] ))
1983 freeLinks[ iL ] = 0;
1984 polygon._links.push_back( *curLink );
1988 // find all links connected to curLink
1989 curNode = curLink->FirstNode();
1991 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1992 if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
1994 curLink = freeLinks[ iL ];
1995 freeLinks[ iL ] = 0;
1997 polygon._links.push_back( *curLink );
1999 } while ( curLink );
2001 else // there are intersections with EDGEs
2003 // get a remaining link to start from, one lying on minimal nb of FACEs
2005 typedef pair< TGeomID, int > TFaceOfLink;
2006 TFaceOfLink faceOfLink( -1, -1 );
2007 TFaceOfLink facesOfLink[3] = { faceOfLink, faceOfLink, faceOfLink };
2008 for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2009 if ( freeLinks[ iL ] )
2011 faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2012 if ( faces.size() == 1 )
2014 faceOfLink = TFaceOfLink( faces[0], iL );
2015 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2017 facesOfLink[0] = faceOfLink;
2019 else if ( facesOfLink[0].first < 0 )
2021 faceOfLink = TFaceOfLink(( faces.empty() ? -1 : faces[0]), iL );
2022 facesOfLink[ 1 + faces.empty() ] = faceOfLink;
2025 for ( int i = 0; faceOfLink.first < 0 && i < 3; ++i )
2026 faceOfLink = facesOfLink[i];
2028 if ( faceOfLink.first < 0 ) // all faces used
2030 for ( size_t iL = 0; iL < freeLinks.size() && faceOfLink.first < 1; ++iL )
2031 if (( curLink = freeLinks[ iL ]))
2034 curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2035 faceOfLink.second = iL;
2037 usedFaceIDs.clear();
2039 curFace = faceOfLink.first;
2040 curLink = freeLinks[ faceOfLink.second ];
2041 freeLinks[ faceOfLink.second ] = 0;
2043 usedFaceIDs.insert( curFace );
2044 polygon._links.push_back( *curLink );
2047 // find all links lying on a curFace
2050 // go forward from curLink
2051 curNode = curLink->LastNode();
2053 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2054 if ( freeLinks[ iL ] &&
2055 freeLinks[ iL ]->FirstNode() == curNode &&
2056 freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2058 curLink = freeLinks[ iL ];
2059 freeLinks[ iL ] = 0;
2060 polygon._links.push_back( *curLink );
2063 } while ( curLink );
2065 std::reverse( polygon._links.begin(), polygon._links.end() );
2067 curLink = & polygon._links.back();
2070 // go backward from curLink
2071 curNode = curLink->FirstNode();
2073 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2074 if ( freeLinks[ iL ] &&
2075 freeLinks[ iL ]->LastNode() == curNode &&
2076 freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2078 curLink = freeLinks[ iL ];
2079 freeLinks[ iL ] = 0;
2080 polygon._links.push_back( *curLink );
2083 } while ( curLink );
2085 curNode = polygon._links.back().FirstNode();
2087 if ( polygon._links[0].LastNode() != curNode )
2089 if ( nbVertexNodes > 0 )
2091 // add links with _vIntNodes if not already used
2093 for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2094 if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2095 _vIntNodes[ iN ]->IsOnFace( curFace ))
2097 _vIntNodes[ iN ]->_usedInFace = &polygon;
2098 chainNodes.push_back( _vIntNodes[ iN ] );
2100 if ( chainNodes.size() > 1 )
2102 sortVertexNodes( chainNodes, curNode, curFace );
2104 for ( size_t i = 0; i < chainNodes.size(); ++i )
2106 polygon.AddPolyLink( chainNodes[ i ], curNode );
2107 curNode = chainNodes[ i ];
2108 freeLinks.push_back( &polygon._links.back() );
2111 nbVertexNodes -= chainNodes.size();
2113 // if ( polygon._links.size() > 1 )
2115 polygon.AddPolyLink( polygon._links[0].LastNode(), curNode );
2116 freeLinks.push_back( &polygon._links.back() );
2120 } // if there are intersections with EDGEs
2122 if ( polygon._links.size() < 2 ||
2123 polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
2124 return; // closed polygon not found -> invalid polyhedron
2126 if ( polygon._links.size() == 2 )
2128 if ( freeLinks.back() == &polygon._links.back() )
2130 freeLinks.pop_back();
2133 if ( polygon._links.front().NbFaces() > 0 )
2134 polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
2135 if ( polygon._links.back().NbFaces() > 0 )
2136 polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
2138 if ( iPolygon == _polygons.size()-1 )
2139 _polygons.pop_back();
2141 else // polygon._links.size() >= 2
2143 // add polygon to its links
2144 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2146 polygon._links[ iL ].AddFace( &polygon );
2147 polygon._links[ iL ].Reverse();
2149 if ( /*hasEdgeIntersections &&*/ iPolygon == _polygons.size() - 1 )
2151 // check that a polygon does not lie on a hexa side
2153 for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
2155 if ( polygon._links[ iL ].NbFaces() < 2 )
2156 continue; // it's a just added free link
2157 // look for a polygon made on a hexa side and sharing
2158 // two or more haxa links
2160 coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
2161 for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
2162 if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
2163 !coplanarPolyg->IsPolyLink( polygon._links[ iL ]) &&
2164 !coplanarPolyg->IsPolyLink( polygon._links[ iL2 ]) &&
2165 coplanarPolyg < & _polygons[ nbQuadPolygons ])
2167 if ( iL2 == polygon._links.size() )
2170 if ( coplanarPolyg ) // coplanar polygon found
2172 freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
2173 nbFreeLinks -= polygon._polyLinks.size();
2175 // an E_IntersectPoint used to mark nodes of coplanarPolyg
2176 // as lying on curFace while they are not at intersection with geometry
2177 ipTmp._faceIDs.resize(1);
2178 ipTmp._faceIDs[0] = curFace;
2180 // fill freeLinks with links not shared by coplanarPolyg and polygon
2181 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2182 if ( polygon._links[ iL ]._link->_faces[1] &&
2183 polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
2185 _Face* p = polygon._links[ iL ]._link->_faces[0];
2186 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2187 if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
2189 freeLinks.push_back( & p->_links[ iL2 ] );
2191 freeLinks.back()->RemoveFace( &polygon );
2195 for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
2196 if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
2197 coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
2199 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
2200 if ( p == coplanarPolyg )
2201 p = coplanarPolyg->_links[ iL ]._link->_faces[1];
2202 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2203 if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
2205 // set links of coplanarPolyg in place of used freeLinks
2206 // to re-create coplanarPolyg next
2208 for ( ; iL3 < freeLinks.size() && freeLinks[ iL3 ]; ++iL3 );
2209 if ( iL3 < freeLinks.size() )
2210 freeLinks[ iL3 ] = ( & p->_links[ iL2 ] );
2212 freeLinks.push_back( & p->_links[ iL2 ] );
2214 freeLinks[ iL3 ]->RemoveFace( coplanarPolyg );
2215 // mark nodes of coplanarPolyg as lying on curFace
2216 for ( int iN = 0; iN < 2; ++iN )
2218 _Node* n = freeLinks[ iL3 ]->_link->_nodes[ iN ];
2219 if ( n->_intPoint ) n->_intPoint->Add( ipTmp._faceIDs );
2220 else n->_intPoint = &ipTmp;
2225 // set coplanarPolyg to be re-created next
2226 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2227 if ( coplanarPolyg == & _polygons[ iP ] )
2230 _polygons[ iPolygon ]._links.clear();
2231 _polygons[ iPolygon ]._polyLinks.clear();
2234 _polygons.pop_back();
2235 usedFaceIDs.erase( curFace );
2237 } // if ( coplanarPolyg )
2238 } // if ( hasEdgeIntersections ) - search for coplanarPolyg
2240 iPolygon = _polygons.size();
2242 } // end of case ( polygon._links.size() > 2 )
2243 } // while ( nbFreeLinks > 0 )
2245 if ( ! checkPolyhedronSize() )
2250 for ( size_t i = 0; i < 8; ++i )
2251 if ( _hexNodes[ i ]._intPoint == &ipTmp )
2252 _hexNodes[ i ]._intPoint = 0;
2254 // create a classic cell if possible
2257 for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2258 nbPolygons += (_polygons[ iF ]._links.size() > 0 );
2260 //const int nbNodes = _nbCornerNodes + nbIntersections;
2262 for ( size_t i = 0; i < 8; ++i )
2263 nbNodes += _hexNodes[ i ].IsUsedInFace();
2264 for ( size_t i = 0; i < _intNodes.size(); ++i )
2265 nbNodes += _intNodes[ i ].IsUsedInFace();
2267 bool isClassicElem = false;
2268 if ( nbNodes == 8 && nbPolygons == 6 ) isClassicElem = addHexa();
2269 else if ( nbNodes == 4 && nbPolygons == 4 ) isClassicElem = addTetra();
2270 else if ( nbNodes == 6 && nbPolygons == 5 ) isClassicElem = addPenta();
2271 else if ( nbNodes == 5 && nbPolygons == 5 ) isClassicElem = addPyra ();
2272 if ( !isClassicElem )
2274 _volumeDefs._nodes.clear();
2275 _volumeDefs._quantities.clear();
2277 for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2279 const size_t nbLinks = _polygons[ iF ]._links.size();
2280 if ( nbLinks == 0 ) continue;
2281 _volumeDefs._quantities.push_back( nbLinks );
2282 for ( size_t iL = 0; iL < nbLinks; ++iL )
2283 _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
2287 //================================================================================
2289 * \brief Create elements in the mesh
2291 int Hexahedron::MakeElements(SMESH_MesherHelper& helper,
2292 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2294 SMESHDS_Mesh* mesh = helper.GetMeshDS();
2296 size_t nbCells[3] = { _grid->_coords[0].size() - 1,
2297 _grid->_coords[1].size() - 1,
2298 _grid->_coords[2].size() - 1 };
2299 const size_t nbGridCells = nbCells[0] * nbCells[1] * nbCells[2];
2300 vector< Hexahedron* > allHexa( nbGridCells, 0 );
2303 // set intersection nodes from GridLine's to links of allHexa
2304 int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
2305 for ( int iDir = 0; iDir < 3; ++iDir )
2307 int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
2308 dInd[1][ iDirOther[iDir][0] ] = -1;
2309 dInd[2][ iDirOther[iDir][1] ] = -1;
2310 dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
2311 // loop on GridLine's parallel to iDir
2312 LineIndexer lineInd = _grid->GetLineIndexer( iDir );
2313 for ( ; lineInd.More(); ++lineInd )
2315 GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
2316 multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2317 for ( ; ip != line._intPoints.end(); ++ip )
2319 // if ( !ip->_node ) continue; // intersection at a grid node
2320 lineInd.SetIndexOnLine( ip->_indexOnLine );
2321 for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
2323 i = int(lineInd.I()) + dInd[iL][0];
2324 j = int(lineInd.J()) + dInd[iL][1];
2325 k = int(lineInd.K()) + dInd[iL][2];
2326 if ( i < 0 || i >= (int) nbCells[0] ||
2327 j < 0 || j >= (int) nbCells[1] ||
2328 k < 0 || k >= (int) nbCells[2] ) continue;
2330 const size_t hexIndex = _grid->CellIndex( i,j,k );
2331 Hexahedron *& hex = allHexa[ hexIndex ];
2334 hex = new Hexahedron( *this );
2340 const int iLink = iL + iDir * 4;
2341 hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
2342 hex->_nbFaceIntNodes += bool( ip->_node );
2348 // implement geom edges into the mesh
2349 addEdges( helper, allHexa, edge2faceIDsMap );
2351 // add not split hexadrons to the mesh
2353 vector< Hexahedron* > intHexa( nbIntHex, (Hexahedron*) NULL );
2354 for ( size_t i = 0; i < allHexa.size(); ++i )
2356 Hexahedron * & hex = allHexa[ i ];
2359 intHexa.push_back( hex );
2360 if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
2361 continue; // treat intersected hex later
2362 this->init( hex->_i, hex->_j, hex->_k );
2368 if (( _nbCornerNodes == 8 ) &&
2369 ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
2371 // order of _hexNodes is defined by enum SMESH_Block::TShapeID
2372 SMDS_MeshElement* el =
2373 mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
2374 _hexNodes[3].Node(), _hexNodes[1].Node(),
2375 _hexNodes[4].Node(), _hexNodes[6].Node(),
2376 _hexNodes[7].Node(), _hexNodes[5].Node() );
2377 mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
2382 else if ( _nbCornerNodes > 3 && !hex )
2384 // all intersection of hex with geometry are at grid nodes
2385 hex = new Hexahedron( *this );
2389 intHexa.push_back( hex );
2393 // add elements resulted from hexadron intersection
2395 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, intHexa.size() ),
2396 ParallelHexahedron( intHexa ),
2397 tbb::simple_partitioner()); // ComputeElements() is called here
2398 for ( size_t i = 0; i < intHexa.size(); ++i )
2399 if ( Hexahedron * hex = intHexa[ i ] )
2400 nbAdded += hex->addElements( helper );
2402 for ( size_t i = 0; i < intHexa.size(); ++i )
2403 if ( Hexahedron * hex = intHexa[ i ] )
2405 hex->ComputeElements();
2406 nbAdded += hex->addElements( helper );
2410 for ( size_t i = 0; i < allHexa.size(); ++i )
2412 delete allHexa[ i ];
2417 //================================================================================
2419 * \brief Implements geom edges into the mesh
2421 void Hexahedron::addEdges(SMESH_MesherHelper& helper,
2422 vector< Hexahedron* >& hexes,
2423 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2425 if ( edge2faceIDsMap.empty() ) return;
2427 // Prepare planes for intersecting with EDGEs
2430 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
2432 GridPlanes& planes = pln[ iDirZ ];
2433 int iDirX = ( iDirZ + 1 ) % 3;
2434 int iDirY = ( iDirZ + 2 ) % 3;
2435 planes._zNorm = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
2436 planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
2437 planes._zProjs [0] = 0;
2438 const double zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2439 const vector< double > & u = _grid->_coords[ iDirZ ];
2440 for ( size_t i = 1; i < planes._zProjs.size(); ++i )
2442 planes._zProjs [i] = zFactor * ( u[i] - u[0] );
2446 const double deflection = _grid->_minCellSize / 20.;
2447 const double tol = _grid->_tol;
2448 E_IntersectPoint ip;
2450 // Intersect EDGEs with the planes
2451 map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
2452 for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
2454 const TGeomID edgeID = e2fIt->first;
2455 const TopoDS_Edge & E = TopoDS::Edge( _grid->_shapes( edgeID ));
2456 BRepAdaptor_Curve curve( E );
2457 TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
2458 TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
2460 ip._faceIDs = e2fIt->second;
2461 ip._shapeID = edgeID;
2463 // discretize the EDGE
2464 GCPnts_UniformDeflection discret( curve, deflection, true );
2465 if ( !discret.IsDone() || discret.NbPoints() < 2 )
2468 // perform intersection
2469 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
2471 GridPlanes& planes = pln[ iDirZ ];
2472 int iDirX = ( iDirZ + 1 ) % 3;
2473 int iDirY = ( iDirZ + 2 ) % 3;
2474 double xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
2475 double yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
2476 double zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
2477 int dIJK[3], d000[3] = { 0,0,0 };
2478 double o[3] = { _grid->_coords[0][0],
2479 _grid->_coords[1][0],
2480 _grid->_coords[2][0] };
2482 // locate the 1st point of a segment within the grid
2483 gp_XYZ p1 = discret.Value( 1 ).XYZ();
2484 double u1 = discret.Parameter( 1 );
2485 double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
2487 _grid->ComputeUVW( p1, ip._uvw );
2488 int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
2489 int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
2490 int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
2491 locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
2492 locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
2493 locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
2495 int ijk[3]; // grid index where a segment intersect a plane
2500 // add the 1st vertex point to a hexahedron
2504 ip._shapeID = _grid->_shapes.Add( v1 );
2505 _grid->_edgeIntP.push_back( ip );
2506 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2507 _grid->_edgeIntP.pop_back();
2508 ip._shapeID = edgeID;
2510 for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
2512 // locate the 2nd point of a segment within the grid
2513 gp_XYZ p2 = discret.Value( iP ).XYZ();
2514 double u2 = discret.Parameter( iP );
2515 double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
2517 if ( Abs( zProj2 - zProj1 ) > std::numeric_limits<double>::min() )
2519 locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
2521 // treat intersections with planes between 2 end points of a segment
2522 int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
2523 int iZ = iZ1 + ( iZ1 < iZ2 );
2524 for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
2526 ip._point = findIntPoint( u1, zProj1, u2, zProj2,
2527 planes._zProjs[ iZ ],
2528 curve, planes._zNorm, _grid->_origin );
2529 _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
2530 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2531 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2534 // add ip to hex "above" the plane
2535 _grid->_edgeIntP.push_back( ip );
2537 bool added = addIntersection(_grid->_edgeIntP.back(), hexes, ijk, dIJK);
2539 // add ip to hex "below" the plane
2540 ijk[ iDirZ ] = iZ-1;
2541 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, dIJK ) &&
2543 _grid->_edgeIntP.pop_back();
2551 // add the 2nd vertex point to a hexahedron
2554 ip._shapeID = _grid->_shapes.Add( v2 );
2556 _grid->ComputeUVW( p1, ip._uvw );
2557 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2558 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2560 _grid->_edgeIntP.push_back( ip );
2561 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2562 _grid->_edgeIntP.pop_back();
2563 ip._shapeID = edgeID;
2565 } // loop on 3 grid directions
2570 //================================================================================
2572 * \brief Finds intersection of a curve with a plane
2573 * \param [in] u1 - parameter of one curve point
2574 * \param [in] proj1 - projection of the curve point to the plane normal
2575 * \param [in] u2 - parameter of another curve point
2576 * \param [in] proj2 - projection of the other curve point to the plane normal
2577 * \param [in] proj - projection of a point where the curve intersects the plane
2578 * \param [in] curve - the curve
2579 * \param [in] axis - the plane normal
2580 * \param [in] origin - the plane origin
2581 * \return gp_Pnt - the found intersection point
2583 gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
2584 double u2, double proj2,
2586 BRepAdaptor_Curve& curve,
2588 const gp_XYZ& origin)
2590 double r = (( proj - proj1 ) / ( proj2 - proj1 ));
2591 double u = u1 * ( 1 - r ) + u2 * r;
2592 gp_Pnt p = curve.Value( u );
2593 double newProj = axis * ( p.XYZ() - origin );
2594 if ( Abs( proj - newProj ) > _grid->_tol / 10. )
2597 return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
2599 return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
2604 //================================================================================
2606 * \brief Returns indices of a hexahedron sub-entities holding a point
2607 * \param [in] ip - intersection point
2608 * \param [out] facets - 0-3 facets holding a point
2609 * \param [out] sub - index of a vertex or an edge holding a point
2610 * \return int - number of facets holding a point
2612 int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
2614 enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
2616 int vertex = 0, edgeMask = 0;
2618 if ( Abs( _grid->_coords[0][ _i ] - ip->_uvw[0] ) < _grid->_tol ) {
2619 facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
2622 else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
2623 facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
2627 if ( Abs( _grid->_coords[1][ _j ] - ip->_uvw[1] ) < _grid->_tol ) {
2628 facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
2631 else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
2632 facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
2636 if ( Abs( _grid->_coords[2][ _k ] - ip->_uvw[2] ) < _grid->_tol ) {
2637 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
2640 else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
2641 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
2648 case 0: sub = 0; break;
2649 case 1: sub = facets[0]; break;
2651 const int edge [3][8] = {
2652 { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
2653 SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
2654 { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
2655 SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
2656 { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
2657 SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
2659 switch ( edgeMask ) {
2660 case X | Y: sub = edge[ 0 ][ vertex ]; break;
2661 case X | Z: sub = edge[ 1 ][ vertex ]; break;
2662 default: sub = edge[ 2 ][ vertex ];
2668 sub = vertex + SMESH_Block::ID_FirstV;
2673 //================================================================================
2675 * \brief Adds intersection with an EDGE
2677 bool Hexahedron::addIntersection( const E_IntersectPoint& ip,
2678 vector< Hexahedron* >& hexes,
2679 int ijk[], int dIJK[] )
2683 size_t hexIndex[4] = {
2684 _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
2685 dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
2686 dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
2687 dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
2689 for ( int i = 0; i < 4; ++i )
2691 if ( /*0 <= hexIndex[i] &&*/ hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
2693 Hexahedron* h = hexes[ hexIndex[i] ];
2694 // check if ip is really inside the hex
2696 if ( h->isOutParam( ip._uvw ))
2697 throw SALOME_Exception("ip outside a hex");
2699 h->_eIntPoints.push_back( & ip );
2705 //================================================================================
2707 * \brief Finds nodes at a path from one node to another via intersections with EDGEs
2709 bool Hexahedron::findChain( _Node* n1,
2712 vector<_Node*>& chn )
2715 chn.push_back( n1 );
2716 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2717 if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2718 n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
2719 n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2721 chn.push_back( quad._eIntNodes[ iP ]);
2722 chn.push_back( n2 );
2723 quad._eIntNodes[ iP ]->_usedInFace = &quad;
2730 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2731 if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2732 chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2734 chn.push_back( quad._eIntNodes[ iP ]);
2735 found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
2738 } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
2740 if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
2741 chn.push_back( n2 );
2743 return chn.size() > 1;
2745 //================================================================================
2747 * \brief Try to heal a polygon whose ends are not connected
2749 bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
2751 int i = -1, nbLinks = polygon->_links.size();
2754 vector< _OrientedLink > newLinks;
2755 // find a node lying on the same FACE as the last one
2756 _Node* node = polygon->_links.back().LastNode();
2757 int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
2758 for ( i = nbLinks - 2; i >= 0; --i )
2759 if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
2763 for ( ; i < nbLinks; ++i )
2764 newLinks.push_back( polygon->_links[i] );
2768 // find a node lying on the same FACE as the first one
2769 node = polygon->_links[0].FirstNode();
2770 avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
2771 for ( i = 1; i < nbLinks; ++i )
2772 if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
2775 for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
2776 newLinks.push_back( polygon->_links[i] );
2778 if ( newLinks.size() > 1 )
2780 polygon->_links.swap( newLinks );
2782 chainNodes.push_back( polygon->_links.back().LastNode() );
2783 chainNodes.push_back( polygon->_links[0].FirstNode() );
2788 //================================================================================
2790 * \brief Finds nodes on the same EDGE as the first node of avoidSplit.
2792 * This function is for a case where an EDGE lies on a quad which lies on a FACE
2793 * so that a part of quad in ON and another part in IN
2795 bool Hexahedron::findChainOnEdge( const vector< _OrientedLink >& splits,
2796 const _OrientedLink& prevSplit,
2797 const _OrientedLink& avoidSplit,
2800 vector<_Node*>& chn )
2802 if ( !isImplementEdges() )
2805 _Node* pn1 = prevSplit.FirstNode();
2806 _Node* pn2 = prevSplit.LastNode();
2807 int avoidFace = pn1->IsLinked( pn2->_intPoint ); // FACE under the quad
2808 if ( avoidFace < 1 && pn1->_intPoint )
2811 _Node* n, *stopNode = avoidSplit.LastNode();
2814 if ( !quad._eIntNodes.empty() )
2816 chn.push_back( pn2 );
2821 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2822 if (( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad )) &&
2823 ( chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint, avoidFace )) &&
2824 ( !avoidFace || quad._eIntNodes[ iP ]->IsOnFace( avoidFace )))
2826 chn.push_back( quad._eIntNodes[ iP ]);
2827 found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
2835 for ( i = splits.size()-1; i >= 0; --i )
2840 n = splits[i].LastNode();
2841 if ( n == stopNode )
2844 ( n->IsLinked( pn2->_intPoint, avoidFace )) &&
2845 ( !avoidFace || n->IsOnFace( avoidFace )))
2848 n = splits[i].FirstNode();
2849 if ( n == stopNode )
2851 if (( n->IsLinked( pn2->_intPoint, avoidFace )) &&
2852 ( !avoidFace || n->IsOnFace( avoidFace )))
2856 if ( n && n != stopNode)
2859 chn.push_back( pn2 );
2866 //================================================================================
2868 * \brief Checks transition at the ginen intersection node of a link
2870 bool Hexahedron::isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper ) const
2874 const bool moreIntPoints = ( iP+1 < (int) link._fIntPoints.size() );
2877 _Node* n1 = link._fIntNodes[ iP ];
2879 n1 = link._nodes[0];
2880 _Node* n2 = moreIntPoints ? link._fIntNodes[ iP+1 ] : 0;
2881 if ( !n2 || !n2->Node() )
2882 n2 = link._nodes[1];
2886 // get all FACEs under n1 and n2
2887 set< TGeomID > faceIDs;
2888 if ( moreIntPoints ) faceIDs.insert( link._fIntPoints[iP+1]->_faceIDs.begin(),
2889 link._fIntPoints[iP+1]->_faceIDs.end() );
2890 if ( n2->_intPoint ) faceIDs.insert( n2->_intPoint->_faceIDs.begin(),
2891 n2->_intPoint->_faceIDs.end() );
2892 if ( faceIDs.empty() )
2893 return false; // n2 is inside
2894 if ( n1->_intPoint ) faceIDs.insert( n1->_intPoint->_faceIDs.begin(),
2895 n1->_intPoint->_faceIDs.end() );
2896 faceIDs.insert( link._fIntPoints[iP]->_faceIDs.begin(),
2897 link._fIntPoints[iP]->_faceIDs.end() );
2899 // get a point between 2 nodes
2900 gp_Pnt p1 = n1->Point();
2901 gp_Pnt p2 = n2->Point();
2902 gp_Pnt pOnLink = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
2904 TopLoc_Location loc;
2906 set< TGeomID >::iterator faceID = faceIDs.begin();
2907 for ( ; faceID != faceIDs.end(); ++faceID )
2909 // project pOnLink on a FACE
2910 if ( *faceID < 1 ) continue;
2911 const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( *faceID ));
2912 GeomAPI_ProjectPointOnSurf& proj =
2913 helper.GetProjector( face, loc, 0.1*_grid->_tol );
2914 gp_Pnt testPnt = pOnLink.Transformed( loc.Transformation().Inverted() );
2915 proj.Perform( testPnt );
2916 if ( proj.IsDone() && proj.NbPoints() > 0 )
2919 proj.LowerDistanceParameters( u,v );
2921 if ( proj.LowerDistance() <= 0.1 * _grid->_tol )
2927 // find isOut by normals
2929 if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
2934 if ( face.Orientation() == TopAbs_REVERSED )
2936 gp_Vec v( proj.NearestPoint(), testPnt );
2937 isOut = ( v * normal > 0 );
2942 // classify a projection
2943 if ( !n1->IsOnFace( *faceID ) || !n2->IsOnFace( *faceID ))
2945 BRepTopAdaptor_FClass2d cls( face, Precision::Confusion() );
2946 TopAbs_State state = cls.Perform( gp_Pnt2d( u,v ));
2947 if ( state == TopAbs_OUT )
2959 //================================================================================
2961 * \brief Sort nodes on a FACE
2963 void Hexahedron::sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID faceID)
2965 if ( nodes.size() > 20 ) return;
2967 // get shapes under nodes
2968 TGeomID nShapeIds[20], *nShapeIdsEnd = &nShapeIds[0] + nodes.size();
2969 for ( size_t i = 0; i < nodes.size(); ++i )
2970 if ( !( nShapeIds[i] = nodes[i]->ShapeID() ))
2973 // get shapes of the FACE
2974 const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( faceID ));
2975 list< TopoDS_Edge > edges;
2976 list< int > nbEdges;
2977 int nbW = SMESH_Block::GetOrderedEdges (face, edges, nbEdges);
2979 // select a WIRE - remove EDGEs of irrelevant WIREs from edges
2980 list< TopoDS_Edge >::iterator e = edges.begin(), eEnd = e;
2981 list< int >::iterator nE = nbEdges.begin();
2982 for ( ; nbW > 0; ++nE, --nbW )
2984 std::advance( eEnd, *nE );
2985 for ( ; e != eEnd; ++e )
2986 for ( int i = 0; i < 2; ++i )
2989 _grid->_shapes.FindIndex( *e ) :
2990 _grid->_shapes.FindIndex( SMESH_MesherHelper::IthVertex( 0, *e ));
2992 ( std::find( &nShapeIds[0], nShapeIdsEnd, id ) != nShapeIdsEnd ))
2994 edges.erase( eEnd, edges.end() ); // remove rest wires
2995 e = eEnd = edges.end();
3002 edges.erase( edges.begin(), eEnd ); // remove a current irrelevant wire
3005 // rotate edges to have the first one at least partially out of the hexa
3006 list< TopoDS_Edge >::iterator e = edges.begin(), eMidOut = edges.end();
3007 for ( ; e != edges.end(); ++e )
3009 if ( !_grid->_shapes.FindIndex( *e ))
3014 for ( int i = 0; i < 2 && !isOut; ++i )
3018 TopoDS_Vertex v = SMESH_MesherHelper::IthVertex( 0, *e );
3019 p = BRep_Tool::Pnt( v );
3021 else if ( eMidOut == edges.end() )
3023 TopLoc_Location loc;
3024 Handle(Geom_Curve) c = BRep_Tool::Curve( *e, loc, f, l);
3025 if ( c.IsNull() ) break;
3026 p = c->Value( 0.5 * ( f + l )).Transformed( loc );
3033 _grid->ComputeUVW( p.XYZ(), uvw );
3034 if ( isOutParam( uvw ))
3045 if ( e != edges.end() )
3046 edges.splice( edges.end(), edges, edges.begin(), e );
3047 else if ( eMidOut != edges.end() )
3048 edges.splice( edges.end(), edges, edges.begin(), eMidOut );
3050 // sort nodes according to the order of edges
3051 _Node* orderNodes [20];
3052 //TGeomID orderShapeIDs[20];
3054 TGeomID id, *pID = 0;
3055 for ( e = edges.begin(); e != edges.end(); ++e )
3057 if (( id = _grid->_shapes.FindIndex( SMESH_MesherHelper::IthVertex( 0, *e ))) &&
3058 (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
3060 //orderShapeIDs[ nbN ] = id;
3061 orderNodes [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
3064 if (( id = _grid->_shapes.FindIndex( *e )) &&
3065 (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
3067 //orderShapeIDs[ nbN ] = id;
3068 orderNodes [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
3072 if ( nbN != nodes.size() )
3075 bool reverse = ( orderNodes[0 ]->Point().SquareDistance( curNode->Point() ) >
3076 orderNodes[nbN-1]->Point().SquareDistance( curNode->Point() ));
3078 for ( size_t i = 0; i < nodes.size(); ++i )
3079 nodes[ i ] = orderNodes[ reverse ? nbN-1-i : i ];
3082 //================================================================================
3084 * \brief Adds computed elements to the mesh
3086 int Hexahedron::addElements(SMESH_MesherHelper& helper)
3089 // add elements resulted from hexahedron intersection
3090 //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
3092 vector< const SMDS_MeshNode* > nodes( _volumeDefs._nodes.size() );
3093 for ( size_t iN = 0; iN < nodes.size(); ++iN )
3094 if ( !( nodes[iN] = _volumeDefs._nodes[iN]->Node() ))
3096 if ( const E_IntersectPoint* eip = _volumeDefs._nodes[iN]->EdgeIntPnt() )
3097 nodes[iN] = _volumeDefs._nodes[iN]->_intPoint->_node =
3098 helper.AddNode( eip->_point.X(),
3102 throw SALOME_Exception("Bug: no node at intersection point");
3105 if ( !_volumeDefs._quantities.empty() )
3107 helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
3111 switch ( nodes.size() )
3113 case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
3114 nodes[4],nodes[5],nodes[6],nodes[7] );
3116 case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
3118 case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
3121 helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
3125 nbAdded += int ( _volumeDefs._nodes.size() > 0 );
3130 //================================================================================
3132 * \brief Return true if the element is in a hole
3134 bool Hexahedron::isInHole() const
3136 if ( !_vIntNodes.empty() )
3139 const size_t ijk[3] = { _i, _j, _k };
3140 F_IntersectPoint curIntPnt;
3142 // consider a cell to be in a hole if all links in any direction
3143 // comes OUT of geometry
3144 for ( int iDir = 0; iDir < 3; ++iDir )
3146 const vector<double>& coords = _grid->_coords[ iDir ];
3147 LineIndexer li = _grid->GetLineIndexer( iDir );
3148 li.SetIJK( _i,_j,_k );
3149 size_t lineIndex[4] = { li.LineIndex (),
3153 bool allLinksOut = true, hasLinks = false;
3154 for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
3156 const _Link& link = _hexLinks[ iL + 4*iDir ];
3157 // check transition of the first node of a link
3158 const F_IntersectPoint* firstIntPnt = 0;
3159 if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
3161 curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0] + _grid->_tol;
3162 const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
3163 multiset< F_IntersectPoint >::const_iterator ip =
3164 line._intPoints.upper_bound( curIntPnt );
3166 firstIntPnt = &(*ip);
3168 else if ( !link._fIntPoints.empty() )
3170 firstIntPnt = link._fIntPoints[0];
3176 allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
3179 if ( hasLinks && allLinksOut )
3185 //================================================================================
3187 * \brief Return true if a polyhedron passes _sizeThreshold criterion
3189 bool Hexahedron::checkPolyhedronSize() const
3192 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3194 const _Face& polygon = _polygons[iP];
3195 if ( polygon._links.empty() )
3197 gp_XYZ area (0,0,0);
3198 gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
3199 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3201 gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
3205 volume += p1 * area;
3209 double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
3211 return volume > initVolume / _sizeThreshold;
3213 //================================================================================
3215 * \brief Tries to create a hexahedron
3217 bool Hexahedron::addHexa()
3219 int nbQuad = 0, iQuad = -1;
3220 for ( size_t i = 0; i < _polygons.size(); ++i )
3222 if ( _polygons[i]._links.empty() )
3224 if ( _polygons[i]._links.size() != 4 )
3235 for ( int iL = 0; iL < 4; ++iL )
3238 nodes[iL] = _polygons[iQuad]._links[iL].FirstNode();
3241 // find a top node above the base node
3242 _Link* link = _polygons[iQuad]._links[iL]._link;
3243 if ( !link->_faces[0] || !link->_faces[1] )
3244 return debugDumpLink( link );
3245 // a quadrangle sharing <link> with _polygons[iQuad]
3246 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[iQuad] )];
3247 for ( int i = 0; i < 4; ++i )
3248 if ( quad->_links[i]._link == link )
3250 // 1st node of a link opposite to <link> in <quad>
3251 nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
3257 _volumeDefs.set( &nodes[0], 8 );
3261 //================================================================================
3263 * \brief Tries to create a tetrahedron
3265 bool Hexahedron::addTetra()
3268 for ( size_t i = 0; i < _polygons.size() && iTria < 0; ++i )
3269 if ( _polygons[i]._links.size() == 3 )
3275 nodes[0] = _polygons[iTria]._links[0].FirstNode();
3276 nodes[1] = _polygons[iTria]._links[1].FirstNode();
3277 nodes[2] = _polygons[iTria]._links[2].FirstNode();
3279 _Link* link = _polygons[iTria]._links[0]._link;
3280 if ( !link->_faces[0] || !link->_faces[1] )
3281 return debugDumpLink( link );
3283 // a triangle sharing <link> with _polygons[0]
3284 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[iTria] )];
3285 for ( int i = 0; i < 3; ++i )
3286 if ( tria->_links[i]._link == link )
3288 nodes[3] = tria->_links[(i+1)%3].LastNode();
3289 _volumeDefs.set( &nodes[0], 4 );
3295 //================================================================================
3297 * \brief Tries to create a pentahedron
3299 bool Hexahedron::addPenta()
3301 // find a base triangular face
3303 for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
3304 if ( _polygons[ iF ]._links.size() == 3 )
3306 if ( iTri < 0 ) return false;
3311 for ( int iL = 0; iL < 3; ++iL )
3314 nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
3317 // find a top node above the base node
3318 _Link* link = _polygons[ iTri ]._links[iL]._link;
3319 if ( !link->_faces[0] || !link->_faces[1] )
3320 return debugDumpLink( link );
3321 // a quadrangle sharing <link> with a base triangle
3322 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
3323 if ( quad->_links.size() != 4 ) return false;
3324 for ( int i = 0; i < 4; ++i )
3325 if ( quad->_links[i]._link == link )
3327 // 1st node of a link opposite to <link> in <quad>
3328 nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
3334 _volumeDefs.set( &nodes[0], 6 );
3336 return ( nbN == 6 );
3338 //================================================================================
3340 * \brief Tries to create a pyramid
3342 bool Hexahedron::addPyra()
3344 // find a base quadrangle
3346 for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
3347 if ( _polygons[ iF ]._links.size() == 4 )
3349 if ( iQuad < 0 ) return false;
3353 nodes[0] = _polygons[iQuad]._links[0].FirstNode();
3354 nodes[1] = _polygons[iQuad]._links[1].FirstNode();
3355 nodes[2] = _polygons[iQuad]._links[2].FirstNode();
3356 nodes[3] = _polygons[iQuad]._links[3].FirstNode();
3358 _Link* link = _polygons[iQuad]._links[0]._link;
3359 if ( !link->_faces[0] || !link->_faces[1] )
3360 return debugDumpLink( link );
3362 // a triangle sharing <link> with a base quadrangle
3363 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
3364 if ( tria->_links.size() != 3 ) return false;
3365 for ( int i = 0; i < 3; ++i )
3366 if ( tria->_links[i]._link == link )
3368 nodes[4] = tria->_links[(i+1)%3].LastNode();
3369 _volumeDefs.set( &nodes[0], 5 );
3375 //================================================================================
3377 * \brief Dump a link and return \c false
3379 bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
3382 gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
3383 cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
3384 << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
3385 << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
3389 //================================================================================
3391 * \brief Classify a point by grid parameters
3393 bool Hexahedron::isOutParam(const double uvw[3]) const
3395 return (( _grid->_coords[0][ _i ] - _grid->_tol > uvw[0] ) ||
3396 ( _grid->_coords[0][ _i+1 ] + _grid->_tol < uvw[0] ) ||
3397 ( _grid->_coords[1][ _j ] - _grid->_tol > uvw[1] ) ||
3398 ( _grid->_coords[1][ _j+1 ] + _grid->_tol < uvw[1] ) ||
3399 ( _grid->_coords[2][ _k ] - _grid->_tol > uvw[2] ) ||
3400 ( _grid->_coords[2][ _k+1 ] + _grid->_tol < uvw[2] ));
3403 //================================================================================
3405 * \brief computes exact bounding box with axes parallel to given ones
3407 //================================================================================
3409 void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
3410 const double* axesDirs,
3414 TopoDS_Compound allFacesComp;
3415 b.MakeCompound( allFacesComp );
3416 for ( size_t iF = 0; iF < faceVec.size(); ++iF )
3417 b.Add( allFacesComp, faceVec[ iF ] );
3419 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
3420 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
3422 for ( int i = 0; i < 6; ++i )
3423 farDist = Max( farDist, 10 * sP[i] );
3425 gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
3426 gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
3427 gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
3428 axis[0].Normalize();
3429 axis[1].Normalize();
3430 axis[2].Normalize();
3432 gp_Mat basis( axis[0], axis[1], axis[2] );
3433 gp_Mat bi = basis.Inverted();
3436 for ( int iDir = 0; iDir < 3; ++iDir )
3438 gp_XYZ axis0 = axis[ iDir ];
3439 gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
3440 gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
3441 for ( int isMax = 0; isMax < 2; ++isMax )
3443 double shift = isMax ? farDist : -farDist;
3444 gp_XYZ orig = shift * axis0;
3445 gp_XYZ norm = axis1 ^ axis2;
3446 gp_Pln pln( orig, norm );
3447 norm = pln.Axis().Direction().XYZ();
3448 BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
3450 gp_Pnt& pAxis = isMax ? pMax : pMin;
3451 gp_Pnt pPlane, pFaces;
3452 double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
3457 for ( int i = 0; i < 2; ++i ) {
3458 corner.SetCoord( 1, sP[ i*3 ]);
3459 for ( int j = 0; j < 2; ++j ) {
3460 corner.SetCoord( 2, sP[ i*3 + 1 ]);
3461 for ( int k = 0; k < 2; ++k )
3463 corner.SetCoord( 3, sP[ i*3 + 2 ]);
3469 corner = isMax ? bb.CornerMax() : bb.CornerMin();
3470 pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
3474 gp_XYZ pf = pFaces.XYZ() * bi;
3475 pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
3481 shapeBox.Add( pMin );
3482 shapeBox.Add( pMax );
3489 //=============================================================================
3491 * \brief Generates 3D structured Cartesian mesh in the internal part of
3492 * solid shapes and polyhedral volumes near the shape boundary.
3493 * \param theMesh - mesh to fill in
3494 * \param theShape - a compound of all SOLIDs to mesh
3495 * \retval bool - true in case of success
3497 //=============================================================================
3499 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh & theMesh,
3500 const TopoDS_Shape & theShape)
3502 // The algorithm generates the mesh in following steps:
3504 // 1) Intersection of grid lines with the geometry boundary.
3505 // This step allows to find out if a given node of the initial grid is
3506 // inside or outside the geometry.
3508 // 2) For each cell of the grid, check how many of it's nodes are outside
3509 // of the geometry boundary. Depending on a result of this check
3510 // - skip a cell, if all it's nodes are outside
3511 // - skip a cell, if it is too small according to the size threshold
3512 // - add a hexahedron in the mesh, if all nodes are inside
3513 // - add a polyhedron in the mesh, if some nodes are inside and some outside
3515 _computeCanceled = false;
3517 SMESH_MesherHelper helper( theMesh );
3522 grid._helper = &helper;
3524 vector< TopoDS_Shape > faceVec;
3526 TopTools_MapOfShape faceMap;
3527 TopExp_Explorer fExp;
3528 for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
3529 if ( !faceMap.Add( fExp.Current() ))
3530 faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
3532 for ( fExp.ReInit(); fExp.More(); fExp.Next() )
3533 if ( faceMap.Contains( fExp.Current() ))
3534 faceVec.push_back( fExp.Current() );
3536 vector<FaceGridIntersector> facesItersectors( faceVec.size() );
3537 map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
3538 TopExp_Explorer eExp;
3540 for ( size_t i = 0; i < faceVec.size(); ++i )
3542 facesItersectors[i]._face = TopoDS::Face ( faceVec[i] );
3543 facesItersectors[i]._faceID = grid._shapes.Add( faceVec[i] );
3544 facesItersectors[i]._grid = &grid;
3545 shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
3547 if ( _hyp->GetToAddEdges() )
3549 helper.SetSubShape( faceVec[i] );
3550 for ( eExp.Init( faceVec[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
3552 const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
3553 if ( !SMESH_Algo::isDegenerated( edge ) &&
3554 !helper.IsRealSeam( edge ))
3555 edge2faceIDsMap[ grid._shapes.Add( edge )].push_back( facesItersectors[i]._faceID );
3560 getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
3562 vector<double> xCoords, yCoords, zCoords;
3563 _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
3565 grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
3567 if ( _computeCanceled ) return false;
3570 { // copy partner faces and curves of not thread-safe types
3571 set< const Standard_Transient* > tshapes;
3572 BRepBuilderAPI_Copy copier;
3573 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3575 if ( !facesItersectors[i].IsThreadSafe(tshapes) )
3577 copier.Perform( facesItersectors[i]._face );
3578 facesItersectors[i]._face = TopoDS::Face( copier );
3582 // Intersection of grid lines with the geometry boundary.
3583 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
3584 ParallelIntersector( facesItersectors ),
3585 tbb::simple_partitioner());
3587 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3588 facesItersectors[i].Intersect();
3591 // put intersection points onto the GridLine's; this is done after intersection
3592 // to avoid contention of facesItersectors for writing into the same GridLine
3593 // in case of parallel work of facesItersectors
3594 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3595 facesItersectors[i].StoreIntersections();
3597 TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
3598 helper.SetSubShape( solidExp.Current() );
3599 helper.SetElementsOnShape( true );
3601 if ( _computeCanceled ) return false;
3603 // create nodes on the geometry
3604 grid.ComputeNodes(helper);
3606 if ( _computeCanceled ) return false;
3608 // create volume elements
3609 Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
3610 int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
3612 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
3615 // make all SOLIDs computed
3616 if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
3618 SMDS_ElemIteratorPtr volIt = sm1->GetElements();
3619 for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
3621 const SMDS_MeshElement* vol = volIt->next();
3622 sm1->RemoveElement( vol, /*isElemDeleted=*/false );
3623 meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
3626 // make other sub-shapes computed
3627 setSubmeshesComputed( theMesh, theShape );
3630 // remove free nodes
3631 if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
3633 TIDSortedNodeSet nodesToRemove;
3634 // get intersection nodes
3635 for ( int iDir = 0; iDir < 3; ++iDir )
3637 vector< GridLine >& lines = grid._lines[ iDir ];
3638 for ( size_t i = 0; i < lines.size(); ++i )
3640 multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
3641 for ( ; ip != lines[i]._intPoints.end(); ++ip )
3642 if ( ip->_node && ip->_node->NbInverseElements() == 0 )
3643 nodesToRemove.insert( nodesToRemove.end(), ip->_node );
3647 for ( size_t i = 0; i < grid._nodes.size(); ++i )
3648 if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
3649 nodesToRemove.insert( nodesToRemove.end(), grid._nodes[i] );
3652 TIDSortedNodeSet::iterator n = nodesToRemove.begin();
3653 for ( ; n != nodesToRemove.end(); ++n )
3654 meshDS->RemoveFreeNode( *n, smDS, /*fromGroups=*/false );
3660 // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
3661 catch ( SMESH_ComputeError& e)
3663 return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
3668 //=============================================================================
3672 //=============================================================================
3674 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh & theMesh,
3675 const TopoDS_Shape & theShape,
3676 MapShapeNbElems& theResMap)
3679 // std::vector<int> aResVec(SMDSEntity_Last);
3680 // for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
3681 // if(IsQuadratic) {
3682 // aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3683 // int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
3684 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
3687 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
3688 // aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3690 // SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
3691 // aResMap.insert(std::make_pair(sm,aResVec));
3696 //=============================================================================
3700 * \brief Event listener setting/unsetting _alwaysComputed flag to
3701 * submeshes of inferior levels to prevent their computing
3703 struct _EventListener : public SMESH_subMeshEventListener
3707 _EventListener(const string& algoName):
3708 SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
3711 // --------------------------------------------------------------------------------
3712 // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
3714 static void setAlwaysComputed( const bool isComputed,
3715 SMESH_subMesh* subMeshOfSolid)
3717 SMESH_subMeshIteratorPtr smIt =
3718 subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3719 while ( smIt->more() )
3721 SMESH_subMesh* sm = smIt->next();
3722 sm->SetIsAlwaysComputed( isComputed );
3724 subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
3727 // --------------------------------------------------------------------------------
3728 // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
3730 virtual void ProcessEvent(const int event,
3731 const int eventType,
3732 SMESH_subMesh* subMeshOfSolid,
3733 SMESH_subMeshEventListenerData* data,
3734 const SMESH_Hypothesis* hyp = 0)
3736 if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
3738 setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
3743 SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
3744 if ( !algo3D || _algoName != algo3D->GetName() )
3745 setAlwaysComputed( false, subMeshOfSolid );
3749 // --------------------------------------------------------------------------------
3750 // set the event listener
3752 static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
3754 subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
3759 }; // struct _EventListener
3763 //================================================================================
3765 * \brief Sets event listener to submeshes if necessary
3766 * \param subMesh - submesh where algo is set
3767 * This method is called when a submesh gets HYP_OK algo_state.
3768 * After being set, event listener is notified on each event of a submesh.
3770 //================================================================================
3772 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
3774 _EventListener::SetOn( subMesh, GetName() );
3777 //================================================================================
3779 * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
3781 //================================================================================
3783 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh& theMesh,
3784 const TopoDS_Shape& theShape)
3786 for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
3787 _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));