1 // Copyright (C) 2007-2015 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 // File : StdMeshers_Cartesian_3D.cxx
25 #include "StdMeshers_Cartesian_3D.hxx"
27 #include "SMDS_MeshNode.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Comment.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_CartesianParameters3D.hxx"
36 #include <utilities.h>
37 #include <Utils_ExceptHandlers.hxx>
38 #include <Basics_OCCTVersion.hxx>
40 #include <GEOMUtils.hxx>
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepBuilderAPI_Copy.hxx>
46 #include <BRepBuilderAPI_MakeFace.hxx>
47 #include <BRepTools.hxx>
48 #include <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_MapOfShape.hxx>
78 #include <TopoDS_Compound.hxx>
79 #include <TopoDS_Face.hxx>
80 #include <TopoDS_TShape.hxx>
81 #include <gp_Cone.hxx>
82 #include <gp_Cylinder.hxx>
85 #include <gp_Pnt2d.hxx>
86 #include <gp_Sphere.hxx>
87 #include <gp_Torus.hxx>
93 #include <tbb/parallel_for.h>
94 //#include <tbb/enumerable_thread_specific.h>
103 //=============================================================================
107 //=============================================================================
109 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
110 :SMESH_3D_Algo(hypId, studyId, gen)
112 _name = "Cartesian_3D";
113 _shapeType = (1 << TopAbs_SOLID); // 1 bit /shape type
114 _compatibleHypothesis.push_back("CartesianParameters3D");
116 _onlyUnaryInput = false; // to mesh all SOLIDs at once
117 _requireDiscreteBoundary = false; // 2D mesh not needed
118 _supportSubmeshes = false; // do not use any existing mesh
121 //=============================================================================
123 * Check presence of a hypothesis
125 //=============================================================================
127 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh& aMesh,
128 const TopoDS_Shape& aShape,
129 Hypothesis_Status& aStatus)
131 aStatus = SMESH_Hypothesis::HYP_MISSING;
133 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
134 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
135 if ( h == hyps.end())
140 for ( ; h != hyps.end(); ++h )
142 if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
144 aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
149 return aStatus == HYP_OK;
156 //=============================================================================
157 // Definitions of internal utils
158 // --------------------------------------------------------------------------
160 Trans_TANGENT = IntCurveSurface_Tangent,
161 Trans_IN = IntCurveSurface_In,
162 Trans_OUT = IntCurveSurface_Out,
165 // --------------------------------------------------------------------------
167 * \brief Common data of any intersection between a Grid and a shape
169 struct B_IntersectPoint
171 mutable const SMDS_MeshNode* _node;
172 mutable vector< TGeomID > _faceIDs;
174 B_IntersectPoint(): _node(NULL) {}
175 void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
176 int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
177 bool IsOnFace( int faceID ) const;
178 virtual ~B_IntersectPoint() {}
180 // --------------------------------------------------------------------------
182 * \brief Data of intersection between a GridLine and a TopoDS_Face
184 struct F_IntersectPoint : public B_IntersectPoint
187 mutable Transition _transition;
188 mutable size_t _indexOnLine;
190 bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
192 // --------------------------------------------------------------------------
194 * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
196 struct E_IntersectPoint : public B_IntersectPoint
202 // --------------------------------------------------------------------------
204 * \brief A line of the grid and its intersections with 2D geometry
209 double _length; // line length
210 multiset< F_IntersectPoint > _intPoints;
212 void RemoveExcessIntPoints( const double tol );
213 bool GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut );
215 // --------------------------------------------------------------------------
217 * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
222 vector< gp_XYZ > _origins; // origin points of all planes in one direction
223 vector< double > _zProjs; // projections of origins to _zNorm
225 // --------------------------------------------------------------------------
227 * \brief Iterator on the parallel grid lines of one direction
233 size_t _iVar1, _iVar2, _iConst;
234 string _name1, _name2, _nameConst;
236 LineIndexer( size_t sz1, size_t sz2, size_t sz3,
237 size_t iv1, size_t iv2, size_t iConst,
238 const string& nv1, const string& nv2, const string& nConst )
240 _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
241 _curInd[0] = _curInd[1] = _curInd[2] = 0;
242 _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
243 _name1 = nv1; _name2 = nv2; _nameConst = nConst;
246 size_t I() const { return _curInd[0]; }
247 size_t J() const { return _curInd[1]; }
248 size_t K() const { return _curInd[2]; }
249 void SetIJK( size_t i, size_t j, size_t k )
251 _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
255 if ( ++_curInd[_iVar1] == _size[_iVar1] )
256 _curInd[_iVar1] = 0, ++_curInd[_iVar2];
258 bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
259 size_t LineIndex () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
260 size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
261 size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
262 size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
263 void SetIndexOnLine (size_t i) { _curInd[ _iConst ] = i; }
264 size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
266 // --------------------------------------------------------------------------
268 * \brief Container of GridLine's
272 vector< double > _coords[3]; // coordinates of grid nodes
273 gp_XYZ _axes [3]; // axis directions
274 vector< GridLine > _lines [3]; // in 3 directions
275 double _tol, _minCellSize;
277 gp_Mat _invB; // inverted basis of _axes
279 vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
280 vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
282 list< E_IntersectPoint > _edgeIntP; // intersections with EDGEs
283 TopTools_IndexedMapOfShape _shapes;
285 SMESH_MesherHelper* _helper;
287 size_t CellIndex( size_t i, size_t j, size_t k ) const
289 return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
291 size_t NodeIndex( size_t i, size_t j, size_t k ) const
293 return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
295 size_t NodeIndexDX() const { return 1; }
296 size_t NodeIndexDY() const { return _coords[0].size(); }
297 size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
299 LineIndexer GetLineIndexer(size_t iDir) const;
301 void SetCoordinates(const vector<double>& xCoords,
302 const vector<double>& yCoords,
303 const vector<double>& zCoords,
304 const double* axesDirs,
305 const Bnd_Box& bndBox );
306 void ComputeUVW(const gp_XYZ& p, double uvw[3]);
307 void ComputeNodes(SMESH_MesherHelper& helper);
309 // --------------------------------------------------------------------------
311 * \brief Intersector of TopoDS_Face with all GridLine's
313 struct FaceGridIntersector
319 IntCurvesFace_Intersector* _surfaceInt;
320 vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
322 FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
325 void StoreIntersections()
327 for ( size_t i = 0; i < _intersections.size(); ++i )
329 multiset< F_IntersectPoint >::iterator ip =
330 _intersections[i].first->_intPoints.insert( _intersections[i].second );
331 ip->_faceIDs.reserve( 1 );
332 ip->_faceIDs.push_back( _faceID );
335 const Bnd_Box& GetFaceBndBox()
337 GetCurveFaceIntersector();
340 IntCurvesFace_Intersector* GetCurveFaceIntersector()
344 _surfaceInt = new IntCurvesFace_Intersector( _face, Precision::PConfusion() );
345 _bndBox = _surfaceInt->Bounding();
346 if ( _bndBox.IsVoid() )
347 BRepBndLib::Add (_face, _bndBox);
351 bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
353 // --------------------------------------------------------------------------
355 * \brief Intersector of a surface with a GridLine
357 struct FaceLineIntersector
360 double _u, _v, _w; // params on the face and the line
361 Transition _transition; // transition of at intersection (see IntCurveSurface.cdl)
362 Transition _transIn, _transOut; // IN and OUT transitions depending of face orientation
365 gp_Cylinder _cylinder;
369 IntCurvesFace_Intersector* _surfaceInt;
371 vector< F_IntersectPoint > _intPoints;
373 void IntersectWithPlane (const GridLine& gridLine);
374 void IntersectWithCylinder(const GridLine& gridLine);
375 void IntersectWithCone (const GridLine& gridLine);
376 void IntersectWithSphere (const GridLine& gridLine);
377 void IntersectWithTorus (const GridLine& gridLine);
378 void IntersectWithSurface (const GridLine& gridLine);
380 bool UVIsOnFace() const;
381 void addIntPoint(const bool toClassify=true);
382 bool isParamOnLineOK( const double linLength )
384 return -_tol < _w && _w < linLength + _tol;
386 FaceLineIntersector():_surfaceInt(0) {}
387 ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
389 // --------------------------------------------------------------------------
391 * \brief Class representing topology of the hexahedron and creating a mesh
392 * volume basing on analysis of hexahedron intersection with geometry
396 // --------------------------------------------------------------------------------
399 // --------------------------------------------------------------------------------
400 struct _Node //!< node either at a hexahedron corner or at intersection
402 const SMDS_MeshNode* _node; // mesh node at hexahedron corner
403 const B_IntersectPoint* _intPoint;
404 const _Face* _usedInFace;
406 _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
407 :_node(n), _intPoint(ip), _usedInFace(0) {}
408 const SMDS_MeshNode* Node() const
409 { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
410 const E_IntersectPoint* EdgeIntPnt() const
411 { return static_cast< const E_IntersectPoint* >( _intPoint ); }
412 bool IsUsedInFace( const _Face* polygon = 0 )
414 return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
416 void Add( const E_IntersectPoint* ip )
421 else if ( !_intPoint->_node ) {
422 ip->Add( _intPoint->_faceIDs );
426 _intPoint->Add( ip->_faceIDs );
429 TGeomID IsLinked( const B_IntersectPoint* other,
430 TGeomID avoidFace=-1 ) const // returns id of a common face
432 return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
434 bool IsOnFace( TGeomID faceID ) const // returns true if faceID is found
436 return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
440 if ( const SMDS_MeshNode* n = Node() )
441 return SMESH_TNodeXYZ( n );
442 if ( const E_IntersectPoint* eip =
443 dynamic_cast< const E_IntersectPoint* >( _intPoint ))
445 return gp_Pnt( 1e100, 0, 0 );
447 TGeomID ShapeID() const
449 if ( const E_IntersectPoint* eip = dynamic_cast< const E_IntersectPoint* >( _intPoint ))
450 return eip->_shapeID;
454 // --------------------------------------------------------------------------------
455 struct _Link // link connecting two _Node's
458 _Face* _faces[2]; // polygons sharing a link
459 vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
460 vector< _Node* > _fIntNodes; // _Node's at _fIntPoints
461 vector< _Link > _splits;
462 _Link() { _faces[0] = 0; }
464 // --------------------------------------------------------------------------------
469 _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
470 void Reverse() { _reverse = !_reverse; }
471 int NbResultLinks() const { return _link->_splits.size(); }
472 _OrientedLink ResultLink(int i) const
474 return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
476 _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
477 _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
478 operator bool() const { return _link; }
479 vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
481 vector< TGeomID > faces;
482 const B_IntersectPoint *ip0, *ip1;
483 if (( ip0 = _link->_nodes[0]->_intPoint ) &&
484 ( ip1 = _link->_nodes[1]->_intPoint ))
486 for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
487 if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
488 !usedIDs.count( ip0->_faceIDs[i] ) )
489 faces.push_back( ip0->_faceIDs[i] );
493 bool HasEdgeNodes() const
495 return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
496 dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
500 return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
502 void AddFace( _Face* f )
504 if ( _link->_faces[0] )
506 _link->_faces[1] = f;
510 _link->_faces[0] = f;
511 _link->_faces[1] = 0;
514 void RemoveFace( _Face* f )
516 if ( !_link->_faces[0] ) return;
518 if ( _link->_faces[1] == f )
520 _link->_faces[1] = 0;
522 else if ( _link->_faces[0] == f )
524 _link->_faces[0] = 0;
525 if ( _link->_faces[1] )
527 _link->_faces[0] = _link->_faces[1];
528 _link->_faces[1] = 0;
533 // --------------------------------------------------------------------------------
536 vector< _OrientedLink > _links; // links on GridLine's
537 vector< _Link > _polyLinks; // links added to close a polygonal face
538 vector< _Node* > _eIntNodes; // nodes at intersection with EDGEs
539 bool IsPolyLink( const _OrientedLink& ol )
541 return _polyLinks.empty() ? false :
542 ( &_polyLinks[0] <= ol._link && ol._link <= &_polyLinks.back() );
544 void AddPolyLink(_Node* n0, _Node* n1, _Face* faceToFindEqual=0)
546 if ( faceToFindEqual && faceToFindEqual != this ) {
547 for ( size_t iL = 0; iL < faceToFindEqual->_polyLinks.size(); ++iL )
548 if ( faceToFindEqual->_polyLinks[iL]._nodes[0] == n1 &&
549 faceToFindEqual->_polyLinks[iL]._nodes[1] == n0 )
552 ( _OrientedLink( & faceToFindEqual->_polyLinks[iL], /*reverse=*/true ));
559 _polyLinks.push_back( l );
560 _links.push_back( _OrientedLink( &_polyLinks.back() ));
563 // --------------------------------------------------------------------------------
564 struct _volumeDef // holder of nodes of a volume mesh element
566 vector< _Node* > _nodes;
567 vector< int > _quantities;
568 typedef boost::shared_ptr<_volumeDef> Ptr;
569 void set( const vector< _Node* >& nodes,
570 const vector< int >& quant = vector< int >() )
571 { _nodes = nodes; _quantities = quant; }
572 void set( _Node** nodes, int nb )
573 { _nodes.assign( nodes, nodes + nb ); }
576 // topology of a hexahedron
579 _Link _hexLinks [12];
582 // faces resulted from hexahedron intersection
583 vector< _Face > _polygons;
585 // intresections with EDGEs
586 vector< const E_IntersectPoint* > _eIntPoints;
588 // additional nodes created at intersection points
589 vector< _Node > _intNodes;
591 // nodes inside the hexahedron (at VERTEXes)
592 vector< _Node* > _vIntNodes;
594 // computed volume elements
595 //vector< _volumeDef::Ptr > _volumeDefs;
596 _volumeDef _volumeDefs;
599 double _sizeThreshold, _sideLength[3];
600 int _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
601 int _origNodeInd; // index of _hexNodes[0] node within the _grid
605 Hexahedron(const double sizeThreshold, Grid* grid);
606 int MakeElements(SMESH_MesherHelper& helper,
607 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
608 void ComputeElements();
609 void Init() { init( _i, _j, _k ); }
612 Hexahedron(const Hexahedron& other );
613 void init( size_t i, size_t j, size_t k );
614 void init( size_t i );
615 void addEdges(SMESH_MesherHelper& helper,
616 vector< Hexahedron* >& intersectedHex,
617 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
618 gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
619 double proj, BRepAdaptor_Curve& curve,
620 const gp_XYZ& axis, const gp_XYZ& origin );
621 int getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
622 bool addIntersection( const E_IntersectPoint& ip,
623 vector< Hexahedron* >& hexes,
624 int ijk[], int dIJK[] );
625 bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
626 bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
627 bool findChainOnEdge( const vector< _OrientedLink >& splits,
628 const _OrientedLink& prevSplit,
629 const _OrientedLink& avoidSplit,
632 vector<_Node*>& chn);
633 int addElements(SMESH_MesherHelper& helper);
634 bool isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper ) const;
635 void sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID face);
636 bool isInHole() const;
637 bool checkPolyhedronSize() const;
642 bool debugDumpLink( _Link* link );
643 _Node* findEqualNode( vector< _Node* >& nodes,
644 const E_IntersectPoint* ip,
647 for ( size_t i = 0; i < nodes.size(); ++i )
648 if ( nodes[i]->EdgeIntPnt() == ip ||
649 nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
653 bool isImplementEdges() const { return !_grid->_edgeIntP.empty(); }
654 bool isOutParam(const double uvw[3]) const;
658 // --------------------------------------------------------------------------
660 * \brief Hexahedron computing volumes in one thread
662 struct ParallelHexahedron
664 vector< Hexahedron* >& _hexVec;
665 ParallelHexahedron( vector< Hexahedron* >& hv ): _hexVec(hv) {}
666 void operator() ( const tbb::blocked_range<size_t>& r ) const
668 for ( size_t i = r.begin(); i != r.end(); ++i )
669 if ( Hexahedron* hex = _hexVec[ i ] )
670 hex->ComputeElements();
673 // --------------------------------------------------------------------------
675 * \brief Structure intersecting certain nb of faces with GridLine's in one thread
677 struct ParallelIntersector
679 vector< FaceGridIntersector >& _faceVec;
680 ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
681 void operator() ( const tbb::blocked_range<size_t>& r ) const
683 for ( size_t i = r.begin(); i != r.end(); ++i )
684 _faceVec[i].Intersect();
689 //=============================================================================
690 // Implementation of internal utils
691 //=============================================================================
693 * \brief adjust \a i to have \a val between values[i] and values[i+1]
695 inline void locateValue( int & i, double val, const vector<double>& values,
696 int& di, double tol )
698 //val += values[0]; // input \a val is measured from 0.
699 if ( i > values.size()-2 )
702 while ( i+2 < values.size() && val > values[ i+1 ])
704 while ( i > 0 && val < values[ i ])
707 if ( i > 0 && val - values[ i ] < tol )
709 else if ( i+2 < values.size() && values[ i+1 ] - val < tol )
714 //=============================================================================
716 * Remove coincident intersection points
718 void GridLine::RemoveExcessIntPoints( const double tol )
720 if ( _intPoints.size() < 2 ) return;
722 set< Transition > tranSet;
723 multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
724 while ( ip2 != _intPoints.end() )
728 while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
730 tranSet.insert( ip1->_transition );
731 tranSet.insert( ip2->_transition );
732 ip2->Add( ip1->_faceIDs );
733 _intPoints.erase( ip1 );
736 if ( tranSet.size() > 1 ) // points with different transition coincide
738 bool isIN = tranSet.count( Trans_IN );
739 bool isOUT = tranSet.count( Trans_OUT );
741 (*ip1)._transition = Trans_TANGENT;
743 (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
747 //================================================================================
749 * Return "is OUT" state for nodes before the given intersection point
751 bool GridLine::GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut )
753 if ( ip->_transition == Trans_IN )
755 if ( ip->_transition == Trans_OUT )
757 if ( ip->_transition == Trans_APEX )
759 // singularity point (apex of a cone)
760 if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
762 multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
763 if ( ipAft == _intPoints.end() )
766 if ( ipBef->_transition != ipAft->_transition )
767 return ( ipBef->_transition == Trans_OUT );
768 return ( ipBef->_transition != Trans_OUT );
770 // _transition == Trans_TANGENT
773 //================================================================================
777 void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
778 const SMDS_MeshNode* n) const
780 if ( _faceIDs.empty() )
783 for ( size_t i = 0; i < fIDs.size(); ++i )
785 vector< TGeomID >::iterator it =
786 std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
787 if ( it == _faceIDs.end() )
788 _faceIDs.push_back( fIDs[i] );
793 //================================================================================
795 * Returns index of a common face if any, else zero
797 int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
800 for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
801 if ( avoidFace != other->_faceIDs[i] &&
802 IsOnFace ( other->_faceIDs[i] ))
803 return other->_faceIDs[i];
806 //================================================================================
808 * Returns \c true if \a faceID in in this->_faceIDs
810 bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
812 vector< TGeomID >::const_iterator it =
813 std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
814 return ( it != _faceIDs.end() );
816 //================================================================================
818 * Return an iterator on GridLine's in a given direction
820 LineIndexer Grid::GetLineIndexer(size_t iDir) const
822 const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
823 const string s [] = { "X", "Y", "Z" };
824 LineIndexer li( _coords[0].size(), _coords[1].size(), _coords[2].size(),
825 indices[iDir*3], indices[iDir*3+1], indices[iDir*3+2],
826 s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
829 //=============================================================================
831 * Creates GridLine's of the grid
833 void Grid::SetCoordinates(const vector<double>& xCoords,
834 const vector<double>& yCoords,
835 const vector<double>& zCoords,
836 const double* axesDirs,
837 const Bnd_Box& shapeBox)
839 _coords[0] = xCoords;
840 _coords[1] = yCoords;
841 _coords[2] = zCoords;
843 _axes[0].SetCoord( axesDirs[0],
846 _axes[1].SetCoord( axesDirs[3],
849 _axes[2].SetCoord( axesDirs[6],
852 _axes[0].Normalize();
853 _axes[1].Normalize();
854 _axes[2].Normalize();
856 _invB.SetCols( _axes[0], _axes[1], _axes[2] );
860 _minCellSize = Precision::Infinite();
861 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
863 for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
865 double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
866 if ( cellLen < _minCellSize )
867 _minCellSize = cellLen;
870 if ( _minCellSize < Precision::Confusion() )
871 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
872 SMESH_Comment("Too small cell size: ") << _minCellSize );
873 _tol = _minCellSize / 1000.;
875 // attune grid extremities to shape bounding box
877 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
878 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
879 double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
880 &_coords[0].back(), &_coords[1].back(), &_coords[2].back() };
881 for ( int i = 0; i < 6; ++i )
882 if ( fabs( sP[i] - *cP[i] ) < _tol )
883 *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
885 for ( int iDir = 0; iDir < 3; ++iDir )
887 if ( _coords[iDir][0] - sP[iDir] > _tol )
889 _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
890 _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
892 if ( sP[iDir+3] - _coords[iDir].back() > _tol )
894 _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
895 _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
898 _tol = _minCellSize / 1000.;
900 _origin = ( _coords[0][0] * _axes[0] +
901 _coords[1][0] * _axes[1] +
902 _coords[2][0] * _axes[2] );
905 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
907 LineIndexer li = GetLineIndexer( iDir );
908 _lines[iDir].resize( li.NbLines() );
909 double len = _coords[ iDir ].back() - _coords[iDir].front();
910 for ( ; li.More(); ++li )
912 GridLine& gl = _lines[iDir][ li.LineIndex() ];
913 gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
914 _coords[1][li.J()] * _axes[1] +
915 _coords[2][li.K()] * _axes[2] );
916 gl._line.SetDirection( _axes[ iDir ]);
921 //================================================================================
923 * Computes coordinates of a point in the grid CS
925 void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
927 gp_XYZ p = P * _invB;
928 p.Coord( UVW[0], UVW[1], UVW[2] );
930 //================================================================================
934 void Grid::ComputeNodes(SMESH_MesherHelper& helper)
936 // state of each node of the grid relative to the geometry
937 const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
938 vector< bool > isNodeOut( nbGridNodes, false );
939 _nodes.resize( nbGridNodes, 0 );
940 _gridIntP.resize( nbGridNodes, NULL );
942 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
944 LineIndexer li = GetLineIndexer( iDir );
946 // find out a shift of node index while walking along a GridLine in this direction
947 li.SetIndexOnLine( 0 );
948 size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
949 li.SetIndexOnLine( 1 );
950 const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
952 const vector<double> & coords = _coords[ iDir ];
953 for ( ; li.More(); ++li ) // loop on lines in iDir
955 li.SetIndexOnLine( 0 );
956 nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
958 GridLine& line = _lines[ iDir ][ li.LineIndex() ];
959 const gp_XYZ lineLoc = line._line.Location().XYZ();
960 const gp_XYZ lineDir = line._line.Direction().XYZ();
961 line.RemoveExcessIntPoints( _tol );
962 multiset< F_IntersectPoint >& intPnts = line._intPoints;
963 multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
966 const double* nodeCoord = & coords[0];
967 const double* coord0 = nodeCoord;
968 const double* coordEnd = coord0 + coords.size();
969 double nodeParam = 0;
970 for ( ; ip != intPnts.end(); ++ip )
972 // set OUT state or just skip IN nodes before ip
973 if ( nodeParam < ip->_paramOnLine - _tol )
975 isOut = line.GetIsOutBefore( ip, isOut );
977 while ( nodeParam < ip->_paramOnLine - _tol )
980 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
981 if ( ++nodeCoord < coordEnd )
982 nodeParam = *nodeCoord - *coord0;
986 if ( nodeCoord == coordEnd ) break;
988 // create a mesh node on a GridLine at ip if it does not coincide with a grid node
989 if ( nodeParam > ip->_paramOnLine + _tol )
991 // li.SetIndexOnLine( 0 );
992 // double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
993 // xyz[ li._iConst ] += ip->_paramOnLine;
994 gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
995 ip->_node = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
996 ip->_indexOnLine = nodeCoord-coord0-1;
998 // create a mesh node at ip concident with a grid node
1001 int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1002 if ( !_nodes[ nodeIndex ] )
1004 //li.SetIndexOnLine( nodeCoord-coord0 );
1005 //double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1006 gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1007 _nodes [ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1008 _gridIntP[ nodeIndex ] = & * ip;
1010 if ( _gridIntP[ nodeIndex ] )
1011 _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1013 _gridIntP[ nodeIndex ] = & * ip;
1014 // ip->_node = _nodes[ nodeIndex ]; -- to differ from ip on links
1015 ip->_indexOnLine = nodeCoord-coord0;
1016 if ( ++nodeCoord < coordEnd )
1017 nodeParam = *nodeCoord - *coord0;
1020 // set OUT state to nodes after the last ip
1021 for ( ; nodeCoord < coordEnd; ++nodeCoord )
1022 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
1026 // Create mesh nodes at !OUT nodes of the grid
1028 for ( size_t z = 0; z < _coords[2].size(); ++z )
1029 for ( size_t y = 0; y < _coords[1].size(); ++y )
1030 for ( size_t x = 0; x < _coords[0].size(); ++x )
1032 size_t nodeIndex = NodeIndex( x, y, z );
1033 if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
1035 //_nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
1036 gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1037 _coords[1][y] * _axes[1] +
1038 _coords[2][z] * _axes[2] );
1039 _nodes[ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1044 // check validity of transitions
1045 const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1046 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1048 LineIndexer li = GetLineIndexer( iDir );
1049 for ( ; li.More(); ++li )
1051 multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1052 if ( intPnts.empty() ) continue;
1053 if ( intPnts.size() == 1 )
1055 if ( intPnts.begin()->_transition != Trans_TANGENT &&
1056 intPnts.begin()->_transition != Trans_APEX )
1057 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1058 SMESH_Comment("Wrong SOLE transition of GridLine (")
1059 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1060 << ") along " << li._nameConst
1061 << ": " << trName[ intPnts.begin()->_transition] );
1065 if ( intPnts.begin()->_transition == Trans_OUT )
1066 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1067 SMESH_Comment("Wrong START transition of GridLine (")
1068 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1069 << ") along " << li._nameConst
1070 << ": " << trName[ intPnts.begin()->_transition ]);
1071 if ( intPnts.rbegin()->_transition == Trans_IN )
1072 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1073 SMESH_Comment("Wrong END transition of GridLine (")
1074 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1075 << ") along " << li._nameConst
1076 << ": " << trName[ intPnts.rbegin()->_transition ]);
1083 //=============================================================================
1085 * Intersects TopoDS_Face with all GridLine's
1087 void FaceGridIntersector::Intersect()
1089 FaceLineIntersector intersector;
1090 intersector._surfaceInt = GetCurveFaceIntersector();
1091 intersector._tol = _grid->_tol;
1092 intersector._transOut = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1093 intersector._transIn = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1095 typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1096 PIntFun interFunction;
1098 bool isDirect = true;
1099 BRepAdaptor_Surface surf( _face );
1100 switch ( surf.GetType() ) {
1102 intersector._plane = surf.Plane();
1103 interFunction = &FaceLineIntersector::IntersectWithPlane;
1104 isDirect = intersector._plane.Direct();
1106 case GeomAbs_Cylinder:
1107 intersector._cylinder = surf.Cylinder();
1108 interFunction = &FaceLineIntersector::IntersectWithCylinder;
1109 isDirect = intersector._cylinder.Direct();
1112 intersector._cone = surf.Cone();
1113 interFunction = &FaceLineIntersector::IntersectWithCone;
1114 //isDirect = intersector._cone.Direct();
1116 case GeomAbs_Sphere:
1117 intersector._sphere = surf.Sphere();
1118 interFunction = &FaceLineIntersector::IntersectWithSphere;
1119 isDirect = intersector._sphere.Direct();
1122 intersector._torus = surf.Torus();
1123 interFunction = &FaceLineIntersector::IntersectWithTorus;
1124 //isDirect = intersector._torus.Direct();
1127 interFunction = &FaceLineIntersector::IntersectWithSurface;
1130 std::swap( intersector._transOut, intersector._transIn );
1132 _intersections.clear();
1133 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1135 if ( surf.GetType() == GeomAbs_Plane )
1137 // check if all lines in this direction are parallel to a plane
1138 if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1139 Precision::Angular()))
1141 // find out a transition, that is the same for all lines of a direction
1142 gp_Dir plnNorm = intersector._plane.Axis().Direction();
1143 gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1144 intersector._transition =
1145 ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1147 if ( surf.GetType() == GeomAbs_Cylinder )
1149 // check if all lines in this direction are parallel to a cylinder
1150 if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1151 Precision::Angular()))
1155 // intersect the grid lines with the face
1156 for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1158 GridLine& gridLine = _grid->_lines[iDir][iL];
1159 if ( _bndBox.IsOut( gridLine._line )) continue;
1161 intersector._intPoints.clear();
1162 (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1163 for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1164 _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1168 //================================================================================
1170 * Return true if (_u,_v) is on the face
1172 bool FaceLineIntersector::UVIsOnFace() const
1174 TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1175 return ( state == TopAbs_IN || state == TopAbs_ON );
1177 //================================================================================
1179 * Store an intersection if it is IN or ON the face
1181 void FaceLineIntersector::addIntPoint(const bool toClassify)
1183 if ( !toClassify || UVIsOnFace() )
1186 p._paramOnLine = _w;
1187 p._transition = _transition;
1188 _intPoints.push_back( p );
1191 //================================================================================
1193 * Intersect a line with a plane
1195 void FaceLineIntersector::IntersectWithPlane(const GridLine& gridLine)
1197 IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1198 _w = linPlane.ParamOnConic(1);
1199 if ( isParamOnLineOK( gridLine._length ))
1201 ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1205 //================================================================================
1207 * Intersect a line with a cylinder
1209 void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1211 IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1212 if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1214 _w = linCylinder.ParamOnConic(1);
1215 if ( linCylinder.NbPoints() == 1 )
1216 _transition = Trans_TANGENT;
1218 _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1219 if ( isParamOnLineOK( gridLine._length ))
1221 ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1224 if ( linCylinder.NbPoints() > 1 )
1226 _w = linCylinder.ParamOnConic(2);
1227 if ( isParamOnLineOK( gridLine._length ))
1229 ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1230 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1236 //================================================================================
1238 * Intersect a line with a cone
1240 void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1242 IntAna_IntConicQuad linCone(gridLine._line,_cone);
1243 if ( !linCone.IsDone() ) return;
1245 gp_Vec du, dv, norm;
1246 for ( int i = 1; i <= linCone.NbPoints(); ++i )
1248 _w = linCone.ParamOnConic( i );
1249 if ( !isParamOnLineOK( gridLine._length )) continue;
1250 ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1253 ElSLib::D1( _u, _v, _cone, P, du, dv );
1255 double normSize2 = norm.SquareMagnitude();
1256 if ( normSize2 > Precision::Angular() * Precision::Angular() )
1258 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1259 cos /= sqrt( normSize2 );
1260 if ( cos < -Precision::Angular() )
1261 _transition = _transIn;
1262 else if ( cos > Precision::Angular() )
1263 _transition = _transOut;
1265 _transition = Trans_TANGENT;
1269 _transition = Trans_APEX;
1271 addIntPoint( /*toClassify=*/false);
1275 //================================================================================
1277 * Intersect a line with a sphere
1279 void FaceLineIntersector::IntersectWithSphere (const GridLine& gridLine)
1281 IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1282 if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1284 _w = linSphere.ParamOnConic(1);
1285 if ( linSphere.NbPoints() == 1 )
1286 _transition = Trans_TANGENT;
1288 _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1289 if ( isParamOnLineOK( gridLine._length ))
1291 ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1294 if ( linSphere.NbPoints() > 1 )
1296 _w = linSphere.ParamOnConic(2);
1297 if ( isParamOnLineOK( gridLine._length ))
1299 ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1300 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1306 //================================================================================
1308 * Intersect a line with a torus
1310 void FaceLineIntersector::IntersectWithTorus (const GridLine& gridLine)
1312 IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1313 if ( !linTorus.IsDone()) return;
1315 gp_Vec du, dv, norm;
1316 for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1318 _w = linTorus.ParamOnLine( i );
1319 if ( !isParamOnLineOK( gridLine._length )) continue;
1320 linTorus.ParamOnTorus( i, _u,_v );
1323 ElSLib::D1( _u, _v, _torus, P, du, dv );
1325 double normSize = norm.Magnitude();
1326 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1328 if ( cos < -Precision::Angular() )
1329 _transition = _transIn;
1330 else if ( cos > Precision::Angular() )
1331 _transition = _transOut;
1333 _transition = Trans_TANGENT;
1334 addIntPoint( /*toClassify=*/false);
1338 //================================================================================
1340 * Intersect a line with a non-analytical surface
1342 void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1344 _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1345 if ( !_surfaceInt->IsDone() ) return;
1346 for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1348 _transition = Transition( _surfaceInt->Transition( i ) );
1349 _w = _surfaceInt->WParameter( i );
1350 addIntPoint(/*toClassify=*/false);
1353 //================================================================================
1355 * check if its face can be safely intersected in a thread
1357 bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1362 TopLoc_Location loc;
1363 Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1364 Handle(Geom_RectangularTrimmedSurface) ts =
1365 Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1366 while( !ts.IsNull() ) {
1367 surf = ts->BasisSurface();
1368 ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1370 if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1371 surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1372 if ( !noSafeTShapes.insert((const Standard_Transient*) _face.TShape() ).second )
1376 TopExp_Explorer exp( _face, TopAbs_EDGE );
1377 for ( ; exp.More(); exp.Next() )
1379 bool edgeIsSafe = true;
1380 const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1383 Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1386 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1387 while( !tc.IsNull() ) {
1388 c = tc->BasisCurve();
1389 tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1391 if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1392 c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1399 Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1402 Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1403 while( !tc.IsNull() ) {
1404 c2 = tc->BasisCurve();
1405 tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1407 if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1408 c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1412 if ( !edgeIsSafe && !noSafeTShapes.insert((const Standard_Transient*) e.TShape() ).second )
1417 //================================================================================
1419 * \brief Creates topology of the hexahedron
1421 Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1422 : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbFaceIntNodes(0)
1424 _polygons.reserve(100); // to avoid reallocation;
1426 //set nodes shift within grid->_nodes from the node 000
1427 size_t dx = _grid->NodeIndexDX();
1428 size_t dy = _grid->NodeIndexDY();
1429 size_t dz = _grid->NodeIndexDZ();
1431 size_t i100 = i000 + dx;
1432 size_t i010 = i000 + dy;
1433 size_t i110 = i010 + dx;
1434 size_t i001 = i000 + dz;
1435 size_t i101 = i100 + dz;
1436 size_t i011 = i010 + dz;
1437 size_t i111 = i110 + dz;
1438 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1439 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1440 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1441 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1442 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1443 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1444 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1445 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1447 vector< int > idVec;
1448 // set nodes to links
1449 for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1451 SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1452 _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1453 link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1454 link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1457 // set links to faces
1458 int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1459 for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1461 SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1462 _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1463 bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1464 faceID == SMESH_Block::ID_Fx1z ||
1465 faceID == SMESH_Block::ID_F0yz );
1466 quad._links.resize(4);
1467 vector<_OrientedLink>::iterator frwLinkIt = quad._links.begin();
1468 vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1469 for ( int i = 0; i < 4; ++i )
1471 bool revLink = revFace;
1472 if ( i > 1 ) // reverse links u1 and v0
1474 _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1475 link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1480 //================================================================================
1482 * \brief Copy constructor
1484 Hexahedron::Hexahedron( const Hexahedron& other )
1485 :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbFaceIntNodes(0)
1487 _polygons.reserve(100); // to avoid reallocation;
1489 for ( int i = 0; i < 8; ++i )
1490 _nodeShift[i] = other._nodeShift[i];
1492 for ( int i = 0; i < 12; ++i )
1494 const _Link& srcLink = other._hexLinks[ i ];
1495 _Link& tgtLink = this->_hexLinks[ i ];
1496 tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1497 tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1500 for ( int i = 0; i < 6; ++i )
1502 const _Face& srcQuad = other._hexQuads[ i ];
1503 _Face& tgtQuad = this->_hexQuads[ i ];
1504 tgtQuad._links.resize(4);
1505 for ( int j = 0; j < 4; ++j )
1507 const _OrientedLink& srcLink = srcQuad._links[ j ];
1508 _OrientedLink& tgtLink = tgtQuad._links[ j ];
1509 tgtLink._reverse = srcLink._reverse;
1510 tgtLink._link = _hexLinks + ( srcLink._link - other._hexLinks );
1515 //================================================================================
1517 * \brief Initializes its data by given grid cell
1519 void Hexahedron::init( size_t i, size_t j, size_t k )
1521 _i = i; _j = j; _k = k;
1522 // set nodes of grid to nodes of the hexahedron and
1523 // count nodes at hexahedron corners located IN and ON geometry
1524 _nbCornerNodes = _nbBndNodes = 0;
1525 _origNodeInd = _grid->NodeIndex( i,j,k );
1526 for ( int iN = 0; iN < 8; ++iN )
1528 _hexNodes[iN]._node = _grid->_nodes [ _origNodeInd + _nodeShift[iN] ];
1529 _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
1530 _nbCornerNodes += bool( _hexNodes[iN]._node );
1531 _nbBndNodes += bool( _hexNodes[iN]._intPoint );
1533 _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1534 _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1535 _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1540 if ( _nbFaceIntNodes + _eIntPoints.size() > 0 &&
1541 _nbFaceIntNodes + _nbCornerNodes + _eIntPoints.size() > 3)
1543 _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
1545 // this method can be called in parallel, so use own helper
1546 SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
1548 // create sub-links (_splits) by splitting links with _fIntPoints
1550 for ( int iLink = 0; iLink < 12; ++iLink )
1552 _Link& link = _hexLinks[ iLink ];
1553 link._fIntNodes.resize( link._fIntPoints.size() );
1554 for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
1556 _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
1557 link._fIntNodes[ i ] = & _intNodes.back();
1560 link._splits.clear();
1561 split._nodes[ 0 ] = link._nodes[0];
1562 bool isOut = ( ! link._nodes[0]->Node() );
1563 bool checkTransition;
1564 for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
1566 const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
1567 if ( !isGridNode ) // intersection non-coincident with a grid node
1569 if ( split._nodes[ 0 ]->Node() && !isOut )
1571 split._nodes[ 1 ] = link._fIntNodes[i];
1572 link._splits.push_back( split );
1574 split._nodes[ 0 ] = link._fIntNodes[i];
1575 checkTransition = true;
1577 else // FACE intersection coincident with a grid node (at link ends)
1579 checkTransition = ( i == 0 && link._nodes[0]->Node() );
1581 if ( checkTransition )
1583 if ( link._fIntPoints[i]->_faceIDs.size() > 1 || _eIntPoints.size() > 0 )
1584 isOut = isOutPoint( link, i, helper );
1586 switch ( link._fIntPoints[i]->_transition ) {
1587 case Trans_OUT: isOut = true; break;
1588 case Trans_IN : isOut = false; break;
1590 isOut = isOutPoint( link, i, helper );
1594 if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
1596 split._nodes[ 1 ] = link._nodes[1];
1597 link._splits.push_back( split );
1601 // Create _Node's at intersections with EDGEs.
1603 const double tol2 = _grid->_tol * _grid->_tol;
1604 int facets[3], nbFacets, subEntity;
1606 for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
1608 nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
1609 _Node* equalNode = 0;
1610 switch( nbFacets ) {
1611 case 1: // in a _Face
1613 _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
1614 equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1616 equalNode->Add( _eIntPoints[ iP ] );
1619 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1620 quad._eIntNodes.push_back( & _intNodes.back() );
1624 case 2: // on a _Link
1626 _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
1627 if ( link._splits.size() > 0 )
1629 equalNode = findEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
1631 equalNode->Add( _eIntPoints[ iP ] );
1635 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1636 for ( int iF = 0; iF < 2; ++iF )
1638 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1639 equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1641 equalNode->Add( _eIntPoints[ iP ] );
1644 quad._eIntNodes.push_back( & _intNodes.back() );
1650 case 3: // at a corner
1652 _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
1653 if ( node.Node() > 0 )
1655 if ( node._intPoint )
1656 node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
1660 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1661 for ( int iF = 0; iF < 3; ++iF )
1663 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1664 equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1666 equalNode->Add( _eIntPoints[ iP ] );
1669 quad._eIntNodes.push_back( & _intNodes.back() );
1675 } // switch( nbFacets )
1677 if ( nbFacets == 0 ||
1678 _grid->_shapes( _eIntPoints[ iP ]->_shapeID ).ShapeType() == TopAbs_VERTEX )
1680 equalNode = findEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
1682 equalNode->Add( _eIntPoints[ iP ] );
1684 else if ( nbFacets == 0 ) {
1685 if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
1686 _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1687 _vIntNodes.push_back( & _intNodes.back() );
1690 } // loop on _eIntPoints
1692 else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
1695 // create sub-links (_splits) of whole links
1696 for ( int iLink = 0; iLink < 12; ++iLink )
1698 _Link& link = _hexLinks[ iLink ];
1699 link._splits.clear();
1700 if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
1702 split._nodes[ 0 ] = link._nodes[0];
1703 split._nodes[ 1 ] = link._nodes[1];
1704 link._splits.push_back( split );
1710 //================================================================================
1712 * \brief Initializes its data by given grid cell (countered from zero)
1714 void Hexahedron::init( size_t iCell )
1716 size_t iNbCell = _grid->_coords[0].size() - 1;
1717 size_t jNbCell = _grid->_coords[1].size() - 1;
1718 _i = iCell % iNbCell;
1719 _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1720 _k = iCell / iNbCell / jNbCell;
1724 //================================================================================
1726 * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1728 void Hexahedron::ComputeElements()
1732 int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
1733 if ( _nbCornerNodes + nbIntersections < 4 )
1736 if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
1740 _polygons.reserve( 20 );
1742 // Create polygons from quadrangles
1743 // --------------------------------
1745 vector< _OrientedLink > splits;
1746 vector<_Node*> chainNodes;
1747 _Face* coplanarPolyg;
1749 bool hasEdgeIntersections = !_eIntPoints.empty();
1751 for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1753 _Face& quad = _hexQuads[ iF ] ;
1755 _polygons.resize( _polygons.size() + 1 );
1756 _Face* polygon = &_polygons.back();
1757 polygon->_polyLinks.reserve( 20 );
1760 for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1761 for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
1762 splits.push_back( quad._links[ iE ].ResultLink( iS ));
1764 // add splits of links to a polygon and add _polyLinks to make
1765 // polygon's boundary closed
1767 int nbSplits = splits.size();
1768 if (( nbSplits == 1 ) &&
1769 ( quad._eIntNodes.empty() ||
1770 splits[0].FirstNode()->IsLinked( splits[0].LastNode()->_intPoint )))
1771 //( quad._eIntNodes.empty() || _nbCornerNodes + nbIntersections > 6 ))
1775 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1776 if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
1777 quad._eIntNodes[ iP ]->_usedInFace = 0;
1779 int nbUsedEdgeNodes = 0;
1780 _Face* prevPolyg = 0; // polygon previously created from this quad
1782 while ( nbSplits > 0 )
1785 while ( !splits[ iS ] )
1788 if ( !polygon->_links.empty() )
1790 _polygons.resize( _polygons.size() + 1 );
1791 polygon = &_polygons.back();
1792 polygon->_polyLinks.reserve( 20 );
1794 polygon->_links.push_back( splits[ iS ] );
1795 splits[ iS++ ]._link = 0;
1798 _Node* nFirst = polygon->_links.back().FirstNode();
1799 _Node *n1,*n2 = polygon->_links.back().LastNode();
1800 for ( ; nFirst != n2 && iS < splits.size(); ++iS )
1802 _OrientedLink& split = splits[ iS ];
1803 if ( !split ) continue;
1805 n1 = split.FirstNode();
1808 n1->_intPoint->_faceIDs.size() > 1 )
1810 // n1 is at intersection with EDGE
1811 if ( findChainOnEdge( splits, polygon->_links.back(), split, iS, quad, chainNodes ))
1813 for ( size_t i = 1; i < chainNodes.size(); ++i )
1814 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
1815 prevPolyg = polygon;
1816 n2 = chainNodes.back();
1820 else if ( n1 != n2 )
1822 // try to connect to intersections with EDGEs
1823 if ( quad._eIntNodes.size() > nbUsedEdgeNodes &&
1824 findChain( n2, n1, quad, chainNodes ))
1826 for ( size_t i = 1; i < chainNodes.size(); ++i )
1828 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i] );
1829 nbUsedEdgeNodes += ( chainNodes[i]->IsUsedInFace( polygon ));
1831 if ( chainNodes.back() != n1 )
1833 n2 = chainNodes.back();
1838 // try to connect to a split ending on the same FACE
1841 _OrientedLink foundSplit;
1842 for ( int i = iS; i < splits.size() && !foundSplit; ++i )
1843 if (( foundSplit = splits[ i ]) &&
1844 ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
1850 foundSplit._link = 0;
1854 if ( n2 != foundSplit.FirstNode() )
1856 polygon->AddPolyLink( n2, foundSplit.FirstNode() );
1857 n2 = foundSplit.FirstNode();
1863 if ( n2->IsLinked( nFirst->_intPoint ))
1865 polygon->AddPolyLink( n2, n1, prevPolyg );
1868 } // if ( n1 != n2 )
1870 polygon->_links.push_back( split );
1873 n2 = polygon->_links.back().LastNode();
1877 if ( nFirst != n2 ) // close a polygon
1879 if ( !findChain( n2, nFirst, quad, chainNodes ))
1881 if ( !closePolygon( polygon, chainNodes ))
1882 if ( !isImplementEdges() )
1883 chainNodes.push_back( nFirst );
1885 for ( size_t i = 1; i < chainNodes.size(); ++i )
1887 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
1888 nbUsedEdgeNodes += bool( chainNodes[i]->IsUsedInFace( polygon ));
1892 if ( polygon->_links.size() < 3 && nbSplits > 0 )
1894 polygon->_polyLinks.clear();
1895 polygon->_links.clear();
1897 } // while ( nbSplits > 0 )
1899 if ( polygon->_links.size() < 3 )
1901 _polygons.pop_back();
1903 } // loop on 6 hexahedron sides
1905 // Create polygons closing holes in a polyhedron
1906 // ----------------------------------------------
1908 // clear _usedInFace
1909 for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
1910 _intNodes[ iN ]._usedInFace = 0;
1912 // add polygons to their links and mark used nodes
1913 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1915 _Face& polygon = _polygons[ iP ];
1916 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1918 polygon._links[ iL ].AddFace( &polygon );
1919 polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
1923 vector< _OrientedLink* > freeLinks;
1924 freeLinks.reserve(20);
1925 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1927 _Face& polygon = _polygons[ iP ];
1928 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1929 if ( polygon._links[ iL ].NbFaces() < 2 )
1930 freeLinks.push_back( & polygon._links[ iL ]);
1932 int nbFreeLinks = freeLinks.size();
1933 if ( nbFreeLinks == 1 ) return;
1935 // put not used intersection nodes to _vIntNodes
1936 int nbVertexNodes = 0; // nb not used vertex nodes
1938 for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
1939 nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
1941 const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
1942 for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
1944 if ( _intNodes[ iN ].IsUsedInFace() ) continue;
1945 if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
1947 findEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
1950 _vIntNodes.push_back( &_intNodes[ iN ]);
1956 set<TGeomID> usedFaceIDs;
1957 vector< TGeomID > faces;
1958 TGeomID curFace = 0;
1959 const size_t nbQuadPolygons = _polygons.size();
1960 E_IntersectPoint ipTmp;
1962 // create polygons by making closed chains of free links
1963 size_t iPolygon = _polygons.size();
1964 while ( nbFreeLinks > 0 )
1966 if ( iPolygon == _polygons.size() )
1968 _polygons.resize( _polygons.size() + 1 );
1969 _polygons[ iPolygon ]._polyLinks.reserve( 20 );
1970 _polygons[ iPolygon ]._links.reserve( 20 );
1972 _Face& polygon = _polygons[ iPolygon ];
1974 _OrientedLink* curLink = 0;
1976 if (( !hasEdgeIntersections ) ||
1977 ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
1979 // get a remaining link to start from
1980 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1981 if (( curLink = freeLinks[ iL ] ))
1982 freeLinks[ iL ] = 0;
1983 polygon._links.push_back( *curLink );
1987 // find all links connected to curLink
1988 curNode = curLink->FirstNode();
1990 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1991 if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
1993 curLink = freeLinks[ iL ];
1994 freeLinks[ iL ] = 0;
1996 polygon._links.push_back( *curLink );
1998 } while ( curLink );
2000 else // there are intersections with EDGEs
2002 // get a remaining link to start from, one lying on minimal nb of FACEs
2004 typedef pair< TGeomID, int > TFaceOfLink;
2005 TFaceOfLink faceOfLink( -1, -1 );
2006 TFaceOfLink facesOfLink[3] = { faceOfLink, faceOfLink, faceOfLink };
2007 for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2008 if ( freeLinks[ iL ] )
2010 faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2011 if ( faces.size() == 1 )
2013 faceOfLink = TFaceOfLink( faces[0], iL );
2014 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2016 facesOfLink[0] = faceOfLink;
2018 else if ( facesOfLink[0].first < 0 )
2020 faceOfLink = TFaceOfLink(( faces.empty() ? -1 : faces[0]), iL );
2021 facesOfLink[ 1 + faces.empty() ] = faceOfLink;
2024 for ( int i = 0; faceOfLink.first < 0 && i < 3; ++i )
2025 faceOfLink = facesOfLink[i];
2027 if ( faceOfLink.first < 0 ) // all faces used
2029 for ( size_t iL = 0; iL < freeLinks.size() && faceOfLink.first < 1; ++iL )
2030 if (( curLink = freeLinks[ iL ]))
2033 curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2034 faceOfLink.second = iL;
2036 usedFaceIDs.clear();
2038 curFace = faceOfLink.first;
2039 curLink = freeLinks[ faceOfLink.second ];
2040 freeLinks[ faceOfLink.second ] = 0;
2042 usedFaceIDs.insert( curFace );
2043 polygon._links.push_back( *curLink );
2046 // find all links lying on a curFace
2049 // go forward from curLink
2050 curNode = curLink->LastNode();
2052 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2053 if ( freeLinks[ iL ] &&
2054 freeLinks[ iL ]->FirstNode() == curNode &&
2055 freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2057 curLink = freeLinks[ iL ];
2058 freeLinks[ iL ] = 0;
2059 polygon._links.push_back( *curLink );
2062 } while ( curLink );
2064 std::reverse( polygon._links.begin(), polygon._links.end() );
2066 curLink = & polygon._links.back();
2069 // go backward from curLink
2070 curNode = curLink->FirstNode();
2072 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2073 if ( freeLinks[ iL ] &&
2074 freeLinks[ iL ]->LastNode() == curNode &&
2075 freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2077 curLink = freeLinks[ iL ];
2078 freeLinks[ iL ] = 0;
2079 polygon._links.push_back( *curLink );
2082 } while ( curLink );
2084 curNode = polygon._links.back().FirstNode();
2086 if ( polygon._links[0].LastNode() != curNode )
2088 if ( nbVertexNodes > 0 )
2090 // add links with _vIntNodes if not already used
2092 for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2093 if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2094 _vIntNodes[ iN ]->IsOnFace( curFace ))
2096 _vIntNodes[ iN ]->_usedInFace = &polygon;
2097 chainNodes.push_back( _vIntNodes[ iN ] );
2099 if ( chainNodes.size() > 1 )
2101 sortVertexNodes( chainNodes, curNode, curFace );
2103 for ( int i = 0; i < chainNodes.size(); ++i )
2105 polygon.AddPolyLink( chainNodes[ i ], curNode );
2106 curNode = chainNodes[ i ];
2107 freeLinks.push_back( &polygon._links.back() );
2110 nbVertexNodes -= chainNodes.size();
2112 // if ( polygon._links.size() > 1 )
2114 polygon.AddPolyLink( polygon._links[0].LastNode(), curNode );
2115 freeLinks.push_back( &polygon._links.back() );
2119 } // if there are intersections with EDGEs
2121 if ( polygon._links.size() < 2 ||
2122 polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
2123 return; // closed polygon not found -> invalid polyhedron
2125 if ( polygon._links.size() == 2 )
2127 if ( freeLinks.back() == &polygon._links.back() )
2129 freeLinks.pop_back();
2132 if ( polygon._links.front().NbFaces() > 0 )
2133 polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
2134 if ( polygon._links.back().NbFaces() > 0 )
2135 polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
2137 if ( iPolygon == _polygons.size()-1 )
2138 _polygons.pop_back();
2140 else // polygon._links.size() >= 2
2142 // add polygon to its links
2143 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2145 polygon._links[ iL ].AddFace( &polygon );
2146 polygon._links[ iL ].Reverse();
2148 if ( /*hasEdgeIntersections &&*/ iPolygon == _polygons.size() - 1 )
2150 // check that a polygon does not lie on a hexa side
2152 for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
2154 if ( polygon._links[ iL ].NbFaces() < 2 )
2155 continue; // it's a just added free link
2156 // look for a polygon made on a hexa side and sharing
2157 // two or more haxa links
2159 coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
2160 for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
2161 if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
2162 !coplanarPolyg->IsPolyLink( polygon._links[ iL ]) &&
2163 !coplanarPolyg->IsPolyLink( polygon._links[ iL2 ]) &&
2164 coplanarPolyg < & _polygons[ nbQuadPolygons ])
2166 if ( iL2 == polygon._links.size() )
2169 if ( coplanarPolyg ) // coplanar polygon found
2171 freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
2172 nbFreeLinks -= polygon._polyLinks.size();
2174 // an E_IntersectPoint used to mark nodes of coplanarPolyg
2175 // as lying on curFace while they are not at intersection with geometry
2176 ipTmp._faceIDs.resize(1);
2177 ipTmp._faceIDs[0] = curFace;
2179 // fill freeLinks with links not shared by coplanarPolyg and polygon
2180 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2181 if ( polygon._links[ iL ]._link->_faces[1] &&
2182 polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
2184 _Face* p = polygon._links[ iL ]._link->_faces[0];
2185 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2186 if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
2188 freeLinks.push_back( & p->_links[ iL2 ] );
2190 freeLinks.back()->RemoveFace( &polygon );
2194 for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
2195 if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
2196 coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
2198 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
2199 if ( p == coplanarPolyg )
2200 p = coplanarPolyg->_links[ iL ]._link->_faces[1];
2201 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2202 if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
2204 // set links of coplanarPolyg in place of used freeLinks
2205 // to re-create coplanarPolyg next
2207 for ( ; iL3 < freeLinks.size() && freeLinks[ iL3 ]; ++iL3 );
2208 if ( iL3 < freeLinks.size() )
2209 freeLinks[ iL3 ] = ( & p->_links[ iL2 ] );
2211 freeLinks.push_back( & p->_links[ iL2 ] );
2213 freeLinks[ iL3 ]->RemoveFace( coplanarPolyg );
2214 // mark nodes of coplanarPolyg as lying on curFace
2215 for ( int iN = 0; iN < 2; ++iN )
2217 _Node* n = freeLinks[ iL3 ]->_link->_nodes[ iN ];
2218 if ( n->_intPoint ) n->_intPoint->Add( ipTmp._faceIDs );
2219 else n->_intPoint = &ipTmp;
2224 // set coplanarPolyg to be re-created next
2225 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2226 if ( coplanarPolyg == & _polygons[ iP ] )
2229 _polygons[ iPolygon ]._links.clear();
2230 _polygons[ iPolygon ]._polyLinks.clear();
2233 _polygons.pop_back();
2234 usedFaceIDs.erase( curFace );
2236 } // if ( coplanarPolyg )
2237 } // if ( hasEdgeIntersections ) - search for coplanarPolyg
2239 iPolygon = _polygons.size();
2241 } // end of case ( polygon._links.size() > 2 )
2242 } // while ( nbFreeLinks > 0 )
2244 if ( ! checkPolyhedronSize() )
2249 for ( size_t i = 0; i < 8; ++i )
2250 if ( _hexNodes[ i ]._intPoint == &ipTmp )
2251 _hexNodes[ i ]._intPoint = 0;
2253 // create a classic cell if possible
2256 for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2257 nbPolygons += (_polygons[ iF ]._links.size() > 0 );
2259 //const int nbNodes = _nbCornerNodes + nbIntersections;
2261 for ( size_t i = 0; i < 8; ++i )
2262 nbNodes += _hexNodes[ i ].IsUsedInFace();
2263 for ( size_t i = 0; i < _intNodes.size(); ++i )
2264 nbNodes += _intNodes[ i ].IsUsedInFace();
2266 bool isClassicElem = false;
2267 if ( nbNodes == 8 && nbPolygons == 6 ) isClassicElem = addHexa();
2268 else if ( nbNodes == 4 && nbPolygons == 4 ) isClassicElem = addTetra();
2269 else if ( nbNodes == 6 && nbPolygons == 5 ) isClassicElem = addPenta();
2270 else if ( nbNodes == 5 && nbPolygons == 5 ) isClassicElem = addPyra ();
2271 if ( !isClassicElem )
2273 _volumeDefs._nodes.clear();
2274 _volumeDefs._quantities.clear();
2276 for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2278 const size_t nbLinks = _polygons[ iF ]._links.size();
2279 if ( nbLinks == 0 ) continue;
2280 _volumeDefs._quantities.push_back( nbLinks );
2281 for ( size_t iL = 0; iL < nbLinks; ++iL )
2282 _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
2286 //================================================================================
2288 * \brief Create elements in the mesh
2290 int Hexahedron::MakeElements(SMESH_MesherHelper& helper,
2291 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2293 SMESHDS_Mesh* mesh = helper.GetMeshDS();
2295 size_t nbCells[3] = { _grid->_coords[0].size() - 1,
2296 _grid->_coords[1].size() - 1,
2297 _grid->_coords[2].size() - 1 };
2298 const size_t nbGridCells = nbCells[0] * nbCells[1] * nbCells[2];
2299 vector< Hexahedron* > allHexa( nbGridCells, 0 );
2302 // set intersection nodes from GridLine's to links of allHexa
2303 int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
2304 for ( int iDir = 0; iDir < 3; ++iDir )
2306 int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
2307 dInd[1][ iDirOther[iDir][0] ] = -1;
2308 dInd[2][ iDirOther[iDir][1] ] = -1;
2309 dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
2310 // loop on GridLine's parallel to iDir
2311 LineIndexer lineInd = _grid->GetLineIndexer( iDir );
2312 for ( ; lineInd.More(); ++lineInd )
2314 GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
2315 multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2316 for ( ; ip != line._intPoints.end(); ++ip )
2318 // if ( !ip->_node ) continue; // intersection at a grid node
2319 lineInd.SetIndexOnLine( ip->_indexOnLine );
2320 for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
2322 i = int(lineInd.I()) + dInd[iL][0];
2323 j = int(lineInd.J()) + dInd[iL][1];
2324 k = int(lineInd.K()) + dInd[iL][2];
2325 if ( i < 0 || i >= nbCells[0] ||
2326 j < 0 || j >= nbCells[1] ||
2327 k < 0 || k >= nbCells[2] ) continue;
2329 const size_t hexIndex = _grid->CellIndex( i,j,k );
2330 Hexahedron *& hex = allHexa[ hexIndex ];
2333 hex = new Hexahedron( *this );
2339 const int iLink = iL + iDir * 4;
2340 hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
2341 hex->_nbFaceIntNodes += bool( ip->_node );
2347 // implement geom edges into the mesh
2348 addEdges( helper, allHexa, edge2faceIDsMap );
2350 // add not split hexadrons to the mesh
2352 vector< Hexahedron* > intHexa( nbIntHex, (Hexahedron*) NULL );
2353 for ( size_t i = 0; i < allHexa.size(); ++i )
2355 Hexahedron * & hex = allHexa[ i ];
2358 intHexa.push_back( hex );
2359 if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
2360 continue; // treat intersected hex later
2361 this->init( hex->_i, hex->_j, hex->_k );
2367 if (( _nbCornerNodes == 8 ) &&
2368 ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
2370 // order of _hexNodes is defined by enum SMESH_Block::TShapeID
2371 SMDS_MeshElement* el =
2372 mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
2373 _hexNodes[3].Node(), _hexNodes[1].Node(),
2374 _hexNodes[4].Node(), _hexNodes[6].Node(),
2375 _hexNodes[7].Node(), _hexNodes[5].Node() );
2376 mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
2381 else if ( _nbCornerNodes > 3 && !hex )
2383 // all intersection of hex with geometry are at grid nodes
2384 hex = new Hexahedron( *this );
2388 intHexa.push_back( hex );
2392 // add elements resulted from hexadron intersection
2394 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, intHexa.size() ),
2395 ParallelHexahedron( intHexa ),
2396 tbb::simple_partitioner()); // ComputeElements() is called here
2397 for ( size_t i = 0; i < intHexa.size(); ++i )
2398 if ( Hexahedron * hex = intHexa[ i ] )
2399 nbAdded += hex->addElements( helper );
2401 for ( size_t i = 0; i < intHexa.size(); ++i )
2402 if ( Hexahedron * hex = intHexa[ i ] )
2404 hex->ComputeElements();
2405 nbAdded += hex->addElements( helper );
2409 for ( size_t i = 0; i < allHexa.size(); ++i )
2411 delete allHexa[ i ];
2416 //================================================================================
2418 * \brief Implements geom edges into the mesh
2420 void Hexahedron::addEdges(SMESH_MesherHelper& helper,
2421 vector< Hexahedron* >& hexes,
2422 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2424 if ( edge2faceIDsMap.empty() ) return;
2426 // Prepare planes for intersecting with EDGEs
2429 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
2431 GridPlanes& planes = pln[ iDirZ ];
2432 int iDirX = ( iDirZ + 1 ) % 3;
2433 int iDirY = ( iDirZ + 2 ) % 3;
2434 planes._zNorm = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
2435 planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
2436 planes._zProjs [0] = 0;
2437 const double zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2438 const vector< double > & u = _grid->_coords[ iDirZ ];
2439 for ( int i = 1; i < planes._zProjs.size(); ++i )
2441 planes._zProjs [i] = zFactor * ( u[i] - u[0] );
2445 const double deflection = _grid->_minCellSize / 20.;
2446 const double tol = _grid->_tol;
2447 E_IntersectPoint ip;
2449 // Intersect EDGEs with the planes
2450 map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
2451 for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
2453 const TGeomID edgeID = e2fIt->first;
2454 const TopoDS_Edge & E = TopoDS::Edge( _grid->_shapes( edgeID ));
2455 BRepAdaptor_Curve curve( E );
2456 TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
2457 TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
2459 ip._faceIDs = e2fIt->second;
2460 ip._shapeID = edgeID;
2462 // discretize the EGDE
2463 GCPnts_UniformDeflection discret( curve, deflection, true );
2464 if ( !discret.IsDone() || discret.NbPoints() < 2 )
2467 // perform intersection
2468 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
2470 GridPlanes& planes = pln[ iDirZ ];
2471 int iDirX = ( iDirZ + 1 ) % 3;
2472 int iDirY = ( iDirZ + 2 ) % 3;
2473 double xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
2474 double yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
2475 double zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
2476 int dIJK[3], d000[3] = { 0,0,0 };
2477 double o[3] = { _grid->_coords[0][0],
2478 _grid->_coords[1][0],
2479 _grid->_coords[2][0] };
2481 // locate the 1st point of a segment within the grid
2482 gp_XYZ p1 = discret.Value( 1 ).XYZ();
2483 double u1 = discret.Parameter( 1 );
2484 double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
2486 _grid->ComputeUVW( p1, ip._uvw );
2487 int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
2488 int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
2489 int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
2490 locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
2491 locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
2492 locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
2494 int ijk[3]; // grid index where a segment intersect a plane
2499 // add the 1st vertex point to a hexahedron
2503 ip._shapeID = _grid->_shapes.Add( v1 );
2504 _grid->_edgeIntP.push_back( ip );
2505 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2506 _grid->_edgeIntP.pop_back();
2507 ip._shapeID = edgeID;
2509 for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
2511 // locate the 2nd point of a segment within the grid
2512 gp_XYZ p2 = discret.Value( iP ).XYZ();
2513 double u2 = discret.Parameter( iP );
2514 double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
2516 if ( Abs( zProj2 - zProj1 ) > std::numeric_limits<double>::min() )
2518 locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
2520 // treat intersections with planes between 2 end points of a segment
2521 int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
2522 int iZ = iZ1 + ( iZ1 < iZ2 );
2523 for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
2525 ip._point = findIntPoint( u1, zProj1, u2, zProj2,
2526 planes._zProjs[ iZ ],
2527 curve, planes._zNorm, _grid->_origin );
2528 _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
2529 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2530 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2533 // add ip to hex "above" the plane
2534 _grid->_edgeIntP.push_back( ip );
2536 bool added = addIntersection(_grid->_edgeIntP.back(), hexes, ijk, dIJK);
2538 // add ip to hex "below" the plane
2539 ijk[ iDirZ ] = iZ-1;
2540 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, dIJK ) &&
2542 _grid->_edgeIntP.pop_back();
2550 // add the 2nd vertex point to a hexahedron
2553 ip._shapeID = _grid->_shapes.Add( v2 );
2555 _grid->ComputeUVW( p1, ip._uvw );
2556 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2557 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2559 _grid->_edgeIntP.push_back( ip );
2560 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2561 _grid->_edgeIntP.pop_back();
2562 ip._shapeID = edgeID;
2564 } // loop on 3 grid directions
2569 //================================================================================
2571 * \brief Finds intersection of a curve with a plane
2572 * \param [in] u1 - parameter of one curve point
2573 * \param [in] proj1 - projection of the curve point to the plane normal
2574 * \param [in] u2 - parameter of another curve point
2575 * \param [in] proj2 - projection of the other curve point to the plane normal
2576 * \param [in] proj - projection of a point where the curve intersects the plane
2577 * \param [in] curve - the curve
2578 * \param [in] axis - the plane normal
2579 * \param [in] origin - the plane origin
2580 * \return gp_Pnt - the found intersection point
2582 gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
2583 double u2, double proj2,
2585 BRepAdaptor_Curve& curve,
2587 const gp_XYZ& origin)
2589 double r = (( proj - proj1 ) / ( proj2 - proj1 ));
2590 double u = u1 * ( 1 - r ) + u2 * r;
2591 gp_Pnt p = curve.Value( u );
2592 double newProj = axis * ( p.XYZ() - origin );
2593 if ( Abs( proj - newProj ) > _grid->_tol / 10. )
2596 return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
2598 return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
2603 //================================================================================
2605 * \brief Returns indices of a hexahedron sub-entities holding a point
2606 * \param [in] ip - intersection point
2607 * \param [out] facets - 0-3 facets holding a point
2608 * \param [out] sub - index of a vertex or an edge holding a point
2609 * \return int - number of facets holding a point
2611 int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
2613 enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
2615 int vertex = 0, egdeMask = 0;
2617 if ( Abs( _grid->_coords[0][ _i ] - ip->_uvw[0] ) < _grid->_tol ) {
2618 facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
2621 else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
2622 facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
2626 if ( Abs( _grid->_coords[1][ _j ] - ip->_uvw[1] ) < _grid->_tol ) {
2627 facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
2630 else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
2631 facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
2635 if ( Abs( _grid->_coords[2][ _k ] - ip->_uvw[2] ) < _grid->_tol ) {
2636 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
2639 else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
2640 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
2647 case 0: sub = 0; break;
2648 case 1: sub = facets[0]; break;
2650 const int edge [3][8] = {
2651 { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
2652 SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
2653 { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
2654 SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
2655 { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
2656 SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
2658 switch ( egdeMask ) {
2659 case X | Y: sub = edge[ 0 ][ vertex ]; break;
2660 case X | Z: sub = edge[ 1 ][ vertex ]; break;
2661 default: sub = edge[ 2 ][ vertex ];
2667 sub = vertex + SMESH_Block::ID_FirstV;
2672 //================================================================================
2674 * \brief Adds intersection with an EDGE
2676 bool Hexahedron::addIntersection( const E_IntersectPoint& ip,
2677 vector< Hexahedron* >& hexes,
2678 int ijk[], int dIJK[] )
2682 size_t hexIndex[4] = {
2683 _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
2684 dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
2685 dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
2686 dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
2688 for ( int i = 0; i < 4; ++i )
2690 if ( /*0 <= hexIndex[i] &&*/ hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
2692 Hexahedron* h = hexes[ hexIndex[i] ];
2693 // check if ip is really inside the hex
2695 if ( h->isOutParam( ip._uvw ))
2696 throw SALOME_Exception("ip outside a hex");
2698 h->_eIntPoints.push_back( & ip );
2704 //================================================================================
2706 * \brief Finds nodes at a path from one node to another via intersections with EDGEs
2708 bool Hexahedron::findChain( _Node* n1,
2711 vector<_Node*>& chn )
2714 chn.push_back( n1 );
2715 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2716 if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2717 n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
2718 n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2720 chn.push_back( quad._eIntNodes[ iP ]);
2721 chn.push_back( n2 );
2722 quad._eIntNodes[ iP ]->_usedInFace = &quad;
2729 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2730 if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2731 chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2733 chn.push_back( quad._eIntNodes[ iP ]);
2734 found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
2737 } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
2739 if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
2740 chn.push_back( n2 );
2742 return chn.size() > 1;
2744 //================================================================================
2746 * \brief Try to heal a polygon whose ends are not connected
2748 bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
2750 int i = -1, nbLinks = polygon->_links.size();
2753 vector< _OrientedLink > newLinks;
2754 // find a node lying on the same FACE as the last one
2755 _Node* node = polygon->_links.back().LastNode();
2756 int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
2757 for ( i = nbLinks - 2; i >= 0; --i )
2758 if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
2762 for ( ; i < nbLinks; ++i )
2763 newLinks.push_back( polygon->_links[i] );
2767 // find a node lying on the same FACE as the first one
2768 node = polygon->_links[0].FirstNode();
2769 avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
2770 for ( i = 1; i < nbLinks; ++i )
2771 if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
2774 for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
2775 newLinks.push_back( polygon->_links[i] );
2777 if ( newLinks.size() > 1 )
2779 polygon->_links.swap( newLinks );
2781 chainNodes.push_back( polygon->_links.back().LastNode() );
2782 chainNodes.push_back( polygon->_links[0].FirstNode() );
2787 //================================================================================
2789 * \brief Finds nodes on the same EDGE as the first node of avoidSplit.
2791 * This function is for a case where an EDGE lies on a quad which lies on a FACE
2792 * so that a part of quad in ON and another part in IN
2794 bool Hexahedron::findChainOnEdge( const vector< _OrientedLink >& splits,
2795 const _OrientedLink& prevSplit,
2796 const _OrientedLink& avoidSplit,
2799 vector<_Node*>& chn )
2801 if ( !isImplementEdges() )
2804 _Node* pn1 = prevSplit.FirstNode();
2805 _Node* pn2 = prevSplit.LastNode();
2806 int avoidFace = pn1->IsLinked( pn2->_intPoint ); // FACE under the quad
2807 if ( avoidFace < 1 && pn1->_intPoint )
2810 _Node* n, *stopNode = avoidSplit.LastNode();
2813 if ( !quad._eIntNodes.empty() )
2815 chn.push_back( pn2 );
2820 for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2821 if (( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad )) &&
2822 ( chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint, avoidFace )) &&
2823 ( !avoidFace || quad._eIntNodes[ iP ]->IsOnFace( avoidFace )))
2825 chn.push_back( quad._eIntNodes[ iP ]);
2826 found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
2834 for ( i = splits.size()-1; i >= 0; --i )
2839 n = splits[i].LastNode();
2840 if ( n == stopNode )
2843 ( n->IsLinked( pn2->_intPoint, avoidFace )) &&
2844 ( !avoidFace || n->IsOnFace( avoidFace )))
2847 n = splits[i].FirstNode();
2848 if ( n == stopNode )
2850 if (( n->IsLinked( pn2->_intPoint, avoidFace )) &&
2851 ( !avoidFace || n->IsOnFace( avoidFace )))
2855 if ( n && n != stopNode)
2858 chn.push_back( pn2 );
2865 //================================================================================
2867 * \brief Checks transition at the ginen intersection node of a link
2869 bool Hexahedron::isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper ) const
2873 const bool moreIntPoints = ( iP+1 < link._fIntPoints.size() );
2876 _Node* n1 = link._fIntNodes[ iP ];
2878 n1 = link._nodes[0];
2879 _Node* n2 = moreIntPoints ? link._fIntNodes[ iP+1 ] : 0;
2880 if ( !n2 || !n2->Node() )
2881 n2 = link._nodes[1];
2885 // get all FACEs under n1 and n2
2886 set< TGeomID > faceIDs;
2887 if ( moreIntPoints ) faceIDs.insert( link._fIntPoints[iP+1]->_faceIDs.begin(),
2888 link._fIntPoints[iP+1]->_faceIDs.end() );
2889 if ( n2->_intPoint ) faceIDs.insert( n2->_intPoint->_faceIDs.begin(),
2890 n2->_intPoint->_faceIDs.end() );
2891 if ( faceIDs.empty() )
2892 return false; // n2 is inside
2893 if ( n1->_intPoint ) faceIDs.insert( n1->_intPoint->_faceIDs.begin(),
2894 n1->_intPoint->_faceIDs.end() );
2895 faceIDs.insert( link._fIntPoints[iP]->_faceIDs.begin(),
2896 link._fIntPoints[iP]->_faceIDs.end() );
2898 // get a point between 2 nodes
2899 gp_Pnt p1 = n1->Point();
2900 gp_Pnt p2 = n2->Point();
2901 gp_Pnt pOnLink = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
2903 TopLoc_Location loc;
2905 set< TGeomID >::iterator faceID = faceIDs.begin();
2906 for ( ; faceID != faceIDs.end(); ++faceID )
2908 // project pOnLink on a FACE
2909 if ( *faceID < 1 ) continue;
2910 const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( *faceID ));
2911 GeomAPI_ProjectPointOnSurf& proj =
2912 helper.GetProjector( face, loc, 0.1*_grid->_tol );
2913 gp_Pnt testPnt = pOnLink.Transformed( loc.Transformation().Inverted() );
2914 proj.Perform( testPnt );
2915 if ( proj.IsDone() && proj.NbPoints() > 0 )
2917 Quantity_Parameter u,v;
2918 proj.LowerDistanceParameters( u,v );
2920 if ( proj.LowerDistance() <= 0.1 * _grid->_tol )
2926 // find isOut by normals
2928 if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
2933 if ( face.Orientation() == TopAbs_REVERSED )
2935 gp_Vec v( proj.NearestPoint(), testPnt );
2936 isOut = ( v * normal > 0 );
2941 // classify a projection
2942 if ( !n1->IsOnFace( *faceID ) || !n2->IsOnFace( *faceID ))
2944 BRepTopAdaptor_FClass2d cls( face, Precision::Confusion() );
2945 TopAbs_State state = cls.Perform( gp_Pnt2d( u,v ));
2946 if ( state == TopAbs_OUT )
2958 //================================================================================
2960 * \brief Sort nodes on a FACE
2962 void Hexahedron::sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID faceID)
2964 if ( nodes.size() > 20 ) return;
2966 // get shapes under nodes
2967 TGeomID nShapeIds[20], *nShapeIdsEnd = &nShapeIds[0] + nodes.size();
2968 for ( size_t i = 0; i < nodes.size(); ++i )
2969 if ( !( nShapeIds[i] = nodes[i]->ShapeID() ))
2972 // get shapes of the FACE
2973 const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( faceID ));
2974 list< TopoDS_Edge > edges;
2975 list< int > nbEdges;
2976 int nbW = SMESH_Block::GetOrderedEdges (face, edges, nbEdges);
2978 // select a WIRE - remove EDGEs of irrelevant WIREs from edges
2979 list< TopoDS_Edge >::iterator e = edges.begin(), eEnd = e;
2980 list< int >::iterator nE = nbEdges.begin();
2981 for ( ; nbW > 0; ++nE, --nbW )
2983 std::advance( eEnd, *nE );
2984 for ( ; e != eEnd; ++e )
2985 for ( int i = 0; i < 2; ++i )
2988 _grid->_shapes.FindIndex( *e ) :
2989 _grid->_shapes.FindIndex( SMESH_MesherHelper::IthVertex( 0, *e ));
2991 ( std::find( &nShapeIds[0], nShapeIdsEnd, id ) != nShapeIdsEnd ))
2993 edges.erase( eEnd, edges.end() ); // remove rest wires
2994 e = eEnd = edges.end();
3001 edges.erase( edges.begin(), eEnd ); // remove a current irrelevant wire
3004 // rotate edges to have the first one at least partially out of the hexa
3005 list< TopoDS_Edge >::iterator e = edges.begin(), eMidOut = edges.end();
3006 for ( ; e != edges.end(); ++e )
3008 if ( !_grid->_shapes.FindIndex( *e ))
3013 for ( int i = 0; i < 2 && !isOut; ++i )
3017 TopoDS_Vertex v = SMESH_MesherHelper::IthVertex( 0, *e );
3018 p = BRep_Tool::Pnt( v );
3020 else if ( eMidOut == edges.end() )
3022 TopLoc_Location loc;
3023 Handle(Geom_Curve) c = BRep_Tool::Curve( *e, loc, f, l);
3024 if ( c.IsNull() ) break;
3025 p = c->Value( 0.5 * ( f + l )).Transformed( loc );
3032 _grid->ComputeUVW( p.XYZ(), uvw );
3033 if ( isOutParam( uvw ))
3044 if ( e != edges.end() )
3045 edges.splice( edges.end(), edges, edges.begin(), e );
3046 else if ( eMidOut != edges.end() )
3047 edges.splice( edges.end(), edges, edges.begin(), eMidOut );
3049 // sort nodes accoring to the order of edges
3050 _Node* orderNodes [20];
3051 TGeomID orderShapeIDs[20];
3054 for ( e = edges.begin(); e != edges.end(); ++e )
3056 if (( id = _grid->_shapes.FindIndex( SMESH_MesherHelper::IthVertex( 0, *e ))) &&
3057 (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
3059 orderShapeIDs[ nbN ] = id;
3060 orderNodes [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
3063 if (( id = _grid->_shapes.FindIndex( *e )) &&
3064 (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
3066 orderShapeIDs[ nbN ] = id;
3067 orderNodes [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
3071 if ( nbN != nodes.size() )
3074 bool reverse = ( orderNodes[0 ]->Point().SquareDistance( curNode->Point() ) >
3075 orderNodes[nbN-1]->Point().SquareDistance( curNode->Point() ));
3077 for ( size_t i = 0; i < nodes.size(); ++i )
3078 nodes[ i ] = orderNodes[ reverse ? nbN-1-i : i ];
3081 //================================================================================
3083 * \brief Adds computed elements to the mesh
3085 int Hexahedron::addElements(SMESH_MesherHelper& helper)
3088 // add elements resulted from hexahedron intersection
3089 //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
3091 vector< const SMDS_MeshNode* > nodes( _volumeDefs._nodes.size() );
3092 for ( size_t iN = 0; iN < nodes.size(); ++iN )
3093 if ( !( nodes[iN] = _volumeDefs._nodes[iN]->Node() ))
3095 if ( const E_IntersectPoint* eip = _volumeDefs._nodes[iN]->EdgeIntPnt() )
3096 nodes[iN] = _volumeDefs._nodes[iN]->_intPoint->_node =
3097 helper.AddNode( eip->_point.X(),
3101 throw SALOME_Exception("Bug: no node at intersection point");
3104 if ( !_volumeDefs._quantities.empty() )
3106 helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
3110 switch ( nodes.size() )
3112 case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
3113 nodes[4],nodes[5],nodes[6],nodes[7] );
3115 case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
3117 case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
3120 helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
3124 nbAdded += int ( _volumeDefs._nodes.size() > 0 );
3129 //================================================================================
3131 * \brief Return true if the element is in a hole
3133 bool Hexahedron::isInHole() const
3135 if ( !_vIntNodes.empty() )
3138 const size_t ijk[3] = { _i, _j, _k };
3139 F_IntersectPoint curIntPnt;
3141 // consider a cell to be in a hole if all links in any direction
3142 // comes OUT of geometry
3143 for ( int iDir = 0; iDir < 3; ++iDir )
3145 const vector<double>& coords = _grid->_coords[ iDir ];
3146 LineIndexer li = _grid->GetLineIndexer( iDir );
3147 li.SetIJK( _i,_j,_k );
3148 size_t lineIndex[4] = { li.LineIndex (),
3152 bool allLinksOut = true, hasLinks = false;
3153 for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
3155 const _Link& link = _hexLinks[ iL + 4*iDir ];
3156 // check transition of the first node of a link
3157 const F_IntersectPoint* firstIntPnt = 0;
3158 if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
3160 curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0];
3161 const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
3162 multiset< F_IntersectPoint >::const_iterator ip =
3163 line._intPoints.upper_bound( curIntPnt );
3165 firstIntPnt = &(*ip);
3167 else if ( !link._fIntPoints.empty() )
3169 firstIntPnt = link._fIntPoints[0];
3175 allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
3178 if ( hasLinks && allLinksOut )
3184 //================================================================================
3186 * \brief Return true if a polyhedron passes _sizeThreshold criterion
3188 bool Hexahedron::checkPolyhedronSize() const
3191 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3193 const _Face& polygon = _polygons[iP];
3194 if ( polygon._links.empty() )
3196 gp_XYZ area (0,0,0);
3197 gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
3198 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3200 gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
3204 volume += p1 * area;
3208 double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
3210 return volume > initVolume / _sizeThreshold;
3212 //================================================================================
3214 * \brief Tries to create a hexahedron
3216 bool Hexahedron::addHexa()
3218 int nbQuad = 0, iQuad = -1;
3219 for ( size_t i = 0; i < _polygons.size(); ++i )
3221 if ( _polygons[i]._links.empty() )
3223 if ( _polygons[i]._links.size() != 4 )
3234 for ( int iL = 0; iL < 4; ++iL )
3237 nodes[iL] = _polygons[iQuad]._links[iL].FirstNode();
3240 // find a top node above the base node
3241 _Link* link = _polygons[iQuad]._links[iL]._link;
3242 if ( !link->_faces[0] || !link->_faces[1] )
3243 return debugDumpLink( link );
3244 // a quadrangle sharing <link> with _polygons[iQuad]
3245 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[iQuad] )];
3246 for ( int i = 0; i < 4; ++i )
3247 if ( quad->_links[i]._link == link )
3249 // 1st node of a link opposite to <link> in <quad>
3250 nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
3256 _volumeDefs.set( &nodes[0], 8 );
3260 //================================================================================
3262 * \brief Tries to create a tetrahedron
3264 bool Hexahedron::addTetra()
3267 for ( size_t i = 0; i < _polygons.size() && iTria < 0; ++i )
3268 if ( _polygons[i]._links.size() == 3 )
3274 nodes[0] = _polygons[iTria]._links[0].FirstNode();
3275 nodes[1] = _polygons[iTria]._links[1].FirstNode();
3276 nodes[2] = _polygons[iTria]._links[2].FirstNode();
3278 _Link* link = _polygons[iTria]._links[0]._link;
3279 if ( !link->_faces[0] || !link->_faces[1] )
3280 return debugDumpLink( link );
3282 // a triangle sharing <link> with _polygons[0]
3283 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[iTria] )];
3284 for ( int i = 0; i < 3; ++i )
3285 if ( tria->_links[i]._link == link )
3287 nodes[3] = tria->_links[(i+1)%3].LastNode();
3288 _volumeDefs.set( &nodes[0], 4 );
3294 //================================================================================
3296 * \brief Tries to create a pentahedron
3298 bool Hexahedron::addPenta()
3300 // find a base triangular face
3302 for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
3303 if ( _polygons[ iF ]._links.size() == 3 )
3305 if ( iTri < 0 ) return false;
3310 for ( int iL = 0; iL < 3; ++iL )
3313 nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
3316 // find a top node above the base node
3317 _Link* link = _polygons[ iTri ]._links[iL]._link;
3318 if ( !link->_faces[0] || !link->_faces[1] )
3319 return debugDumpLink( link );
3320 // a quadrangle sharing <link> with a base triangle
3321 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
3322 if ( quad->_links.size() != 4 ) return false;
3323 for ( int i = 0; i < 4; ++i )
3324 if ( quad->_links[i]._link == link )
3326 // 1st node of a link opposite to <link> in <quad>
3327 nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
3333 _volumeDefs.set( &nodes[0], 6 );
3335 return ( nbN == 6 );
3337 //================================================================================
3339 * \brief Tries to create a pyramid
3341 bool Hexahedron::addPyra()
3343 // find a base quadrangle
3345 for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
3346 if ( _polygons[ iF ]._links.size() == 4 )
3348 if ( iQuad < 0 ) return false;
3352 nodes[0] = _polygons[iQuad]._links[0].FirstNode();
3353 nodes[1] = _polygons[iQuad]._links[1].FirstNode();
3354 nodes[2] = _polygons[iQuad]._links[2].FirstNode();
3355 nodes[3] = _polygons[iQuad]._links[3].FirstNode();
3357 _Link* link = _polygons[iQuad]._links[0]._link;
3358 if ( !link->_faces[0] || !link->_faces[1] )
3359 return debugDumpLink( link );
3361 // a triangle sharing <link> with a base quadrangle
3362 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
3363 if ( tria->_links.size() != 3 ) return false;
3364 for ( int i = 0; i < 3; ++i )
3365 if ( tria->_links[i]._link == link )
3367 nodes[4] = tria->_links[(i+1)%3].LastNode();
3368 _volumeDefs.set( &nodes[0], 5 );
3374 //================================================================================
3376 * \brief Dump a link and return \c false
3378 bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
3381 gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
3382 cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
3383 << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
3384 << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
3388 //================================================================================
3390 * \brief Classify a point by grid paremeters
3392 bool Hexahedron::isOutParam(const double uvw[3]) const
3394 return (( _grid->_coords[0][ _i ] - _grid->_tol > uvw[0] ) ||
3395 ( _grid->_coords[0][ _i+1 ] + _grid->_tol < uvw[0] ) ||
3396 ( _grid->_coords[1][ _j ] - _grid->_tol > uvw[1] ) ||
3397 ( _grid->_coords[1][ _j+1 ] + _grid->_tol < uvw[1] ) ||
3398 ( _grid->_coords[2][ _k ] - _grid->_tol > uvw[2] ) ||
3399 ( _grid->_coords[2][ _k+1 ] + _grid->_tol < uvw[2] ));
3402 //================================================================================
3404 * \brief computes exact bounding box with axes parallel to given ones
3406 //================================================================================
3408 void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
3409 const double* axesDirs,
3413 TopoDS_Compound allFacesComp;
3414 b.MakeCompound( allFacesComp );
3415 for ( size_t iF = 0; iF < faceVec.size(); ++iF )
3416 b.Add( allFacesComp, faceVec[ iF ] );
3418 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
3419 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
3421 for ( int i = 0; i < 6; ++i )
3422 farDist = Max( farDist, 10 * sP[i] );
3424 gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
3425 gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
3426 gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
3427 axis[0].Normalize();
3428 axis[1].Normalize();
3429 axis[2].Normalize();
3431 gp_Mat basis( axis[0], axis[1], axis[2] );
3432 gp_Mat bi = basis.Inverted();
3435 for ( int iDir = 0; iDir < 3; ++iDir )
3437 gp_XYZ axis0 = axis[ iDir ];
3438 gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
3439 gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
3440 for ( int isMax = 0; isMax < 2; ++isMax )
3442 double shift = isMax ? farDist : -farDist;
3443 gp_XYZ orig = shift * axis0;
3444 gp_XYZ norm = axis1 ^ axis2;
3445 gp_Pln pln( orig, norm );
3446 norm = pln.Axis().Direction().XYZ();
3447 BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
3449 gp_Pnt& pAxis = isMax ? pMax : pMin;
3450 gp_Pnt pPlane, pFaces;
3451 double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
3456 for ( int i = 0; i < 2; ++i ) {
3457 corner.SetCoord( 1, sP[ i*3 ]);
3458 for ( int j = 0; j < 2; ++j ) {
3459 corner.SetCoord( 2, sP[ i*3 + 1 ]);
3460 for ( int k = 0; k < 2; ++k )
3462 corner.SetCoord( 3, sP[ i*3 + 2 ]);
3468 corner = isMax ? bb.CornerMax() : bb.CornerMin();
3469 pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
3473 gp_XYZ pf = pFaces.XYZ() * bi;
3474 pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
3480 shapeBox.Add( pMin );
3481 shapeBox.Add( pMax );
3488 //=============================================================================
3490 * \brief Generates 3D structured Cartesian mesh in the internal part of
3491 * solid shapes and polyhedral volumes near the shape boundary.
3492 * \param theMesh - mesh to fill in
3493 * \param theShape - a compound of all SOLIDs to mesh
3494 * \retval bool - true in case of success
3496 //=============================================================================
3498 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh & theMesh,
3499 const TopoDS_Shape & theShape)
3501 // The algorithm generates the mesh in following steps:
3503 // 1) Intersection of grid lines with the geometry boundary.
3504 // This step allows to find out if a given node of the initial grid is
3505 // inside or outside the geometry.
3507 // 2) For each cell of the grid, check how many of it's nodes are outside
3508 // of the geometry boundary. Depending on a result of this check
3509 // - skip a cell, if all it's nodes are outside
3510 // - skip a cell, if it is too small according to the size threshold
3511 // - add a hexahedron in the mesh, if all nodes are inside
3512 // - add a polyhedron in the mesh, if some nodes are inside and some outside
3514 _computeCanceled = false;
3516 SMESH_MesherHelper helper( theMesh );
3521 grid._helper = &helper;
3523 vector< TopoDS_Shape > faceVec;
3525 TopTools_MapOfShape faceMap;
3526 TopExp_Explorer fExp;
3527 for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
3528 if ( !faceMap.Add( fExp.Current() ))
3529 faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
3531 for ( fExp.ReInit(); fExp.More(); fExp.Next() )
3532 if ( faceMap.Contains( fExp.Current() ))
3533 faceVec.push_back( fExp.Current() );
3535 vector<FaceGridIntersector> facesItersectors( faceVec.size() );
3536 map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
3537 TopExp_Explorer eExp;
3539 for ( int i = 0; i < faceVec.size(); ++i )
3541 facesItersectors[i]._face = TopoDS::Face ( faceVec[i] );
3542 facesItersectors[i]._faceID = grid._shapes.Add( faceVec[i] );
3543 facesItersectors[i]._grid = &grid;
3544 shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
3546 if ( _hyp->GetToAddEdges() )
3548 helper.SetSubShape( faceVec[i] );
3549 for ( eExp.Init( faceVec[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
3551 const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
3552 if ( !SMESH_Algo::isDegenerated( edge ) &&
3553 !helper.IsRealSeam( edge ))
3554 edge2faceIDsMap[ grid._shapes.Add( edge )].push_back( facesItersectors[i]._faceID );
3559 getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
3561 vector<double> xCoords, yCoords, zCoords;
3562 _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
3564 grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
3566 if ( _computeCanceled ) return false;
3569 { // copy partner faces and curves of not thread-safe types
3570 set< const Standard_Transient* > tshapes;
3571 BRepBuilderAPI_Copy copier;
3572 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3574 if ( !facesItersectors[i].IsThreadSafe(tshapes) )
3576 copier.Perform( facesItersectors[i]._face );
3577 facesItersectors[i]._face = TopoDS::Face( copier );
3581 // Intersection of grid lines with the geometry boundary.
3582 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
3583 ParallelIntersector( facesItersectors ),
3584 tbb::simple_partitioner());
3586 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3587 facesItersectors[i].Intersect();
3590 // put interesection points onto the GridLine's; this is done after intersection
3591 // to avoid contention of facesItersectors for writing into the same GridLine
3592 // in case of parallel work of facesItersectors
3593 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3594 facesItersectors[i].StoreIntersections();
3596 TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
3597 helper.SetSubShape( solidExp.Current() );
3598 helper.SetElementsOnShape( true );
3600 if ( _computeCanceled ) return false;
3602 // create nodes on the geometry
3603 grid.ComputeNodes(helper);
3605 if ( _computeCanceled ) return false;
3607 // create volume elements
3608 Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
3609 int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
3611 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
3614 // make all SOLIDs computed
3615 if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
3617 SMDS_ElemIteratorPtr volIt = sm1->GetElements();
3618 for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
3620 const SMDS_MeshElement* vol = volIt->next();
3621 sm1->RemoveElement( vol, /*isElemDeleted=*/false );
3622 meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
3625 // make other sub-shapes computed
3626 setSubmeshesComputed( theMesh, theShape );
3629 // remove free nodes
3630 if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
3632 TIDSortedNodeSet nodesToRemove;
3633 // get intersection nodes
3634 for ( int iDir = 0; iDir < 3; ++iDir )
3636 vector< GridLine >& lines = grid._lines[ iDir ];
3637 for ( size_t i = 0; i < lines.size(); ++i )
3639 multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
3640 for ( ; ip != lines[i]._intPoints.end(); ++ip )
3641 if ( ip->_node && ip->_node->NbInverseElements() == 0 )
3642 nodesToRemove.insert( nodesToRemove.end(), ip->_node );
3646 for ( size_t i = 0; i < grid._nodes.size(); ++i )
3647 if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
3648 nodesToRemove.insert( nodesToRemove.end(), grid._nodes[i] );
3651 TIDSortedNodeSet::iterator n = nodesToRemove.begin();
3652 for ( ; n != nodesToRemove.end(); ++n )
3653 meshDS->RemoveFreeNode( *n, smDS, /*fromGroups=*/false );
3659 // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
3660 catch ( SMESH_ComputeError& e)
3662 return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
3667 //=============================================================================
3671 //=============================================================================
3673 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh & theMesh,
3674 const TopoDS_Shape & theShape,
3675 MapShapeNbElems& theResMap)
3678 // std::vector<int> aResVec(SMDSEntity_Last);
3679 // for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
3680 // if(IsQuadratic) {
3681 // aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3682 // int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
3683 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
3686 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
3687 // aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3689 // SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
3690 // aResMap.insert(std::make_pair(sm,aResVec));
3695 //=============================================================================
3699 * \brief Event listener setting/unsetting _alwaysComputed flag to
3700 * submeshes of inferior levels to prevent their computing
3702 struct _EventListener : public SMESH_subMeshEventListener
3706 _EventListener(const string& algoName):
3707 SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
3710 // --------------------------------------------------------------------------------
3711 // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
3713 static void setAlwaysComputed( const bool isComputed,
3714 SMESH_subMesh* subMeshOfSolid)
3716 SMESH_subMeshIteratorPtr smIt =
3717 subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3718 while ( smIt->more() )
3720 SMESH_subMesh* sm = smIt->next();
3721 sm->SetIsAlwaysComputed( isComputed );
3723 subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
3726 // --------------------------------------------------------------------------------
3727 // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
3729 virtual void ProcessEvent(const int event,
3730 const int eventType,
3731 SMESH_subMesh* subMeshOfSolid,
3732 SMESH_subMeshEventListenerData* data,
3733 const SMESH_Hypothesis* hyp = 0)
3735 if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
3737 setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
3742 SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
3743 if ( !algo3D || _algoName != algo3D->GetName() )
3744 setAlwaysComputed( false, subMeshOfSolid );
3748 // --------------------------------------------------------------------------------
3749 // set the event listener
3751 static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
3753 subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
3758 }; // struct _EventListener
3762 //================================================================================
3764 * \brief Sets event listener to submeshes if necessary
3765 * \param subMesh - submesh where algo is set
3766 * This method is called when a submesh gets HYP_OK algo_state.
3767 * After being set, event listener is notified on each event of a submesh.
3769 //================================================================================
3771 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
3773 _EventListener::SetOn( subMesh, GetName() );
3776 //================================================================================
3778 * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
3780 //================================================================================
3782 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh& theMesh,
3783 const TopoDS_Shape& theShape)
3785 for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
3786 _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));