1 // Copyright (C) 2007-2013 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.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 // File : StdMeshers_Cartesian_3D.cxx
25 #include "StdMeshers_Cartesian_3D.hxx"
27 #include "SMDS_MeshNode.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Comment.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_CartesianParameters3D.hxx"
36 #include <utilities.h>
37 #include <Utils_ExceptHandlers.hxx>
38 #include <Basics_OCCTVersion.hxx>
40 #include <GEOMUtils.hxx>
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepBuilderAPI_Copy.hxx>
46 #include <BRepBuilderAPI_MakeFace.hxx>
47 #include <BRepTools.hxx>
48 #include <BRep_Builder.hxx>
49 #include <BRep_Tool.hxx>
50 #include <Bnd_B3d.hxx>
51 #include <Bnd_Box.hxx>
53 #include <GCPnts_UniformDeflection.hxx>
54 #include <Geom2d_BSplineCurve.hxx>
55 #include <Geom2d_BezierCurve.hxx>
56 #include <Geom2d_TrimmedCurve.hxx>
57 #include <Geom_BSplineCurve.hxx>
58 #include <Geom_BSplineSurface.hxx>
59 #include <Geom_BezierCurve.hxx>
60 #include <Geom_BezierSurface.hxx>
61 #include <Geom_RectangularTrimmedSurface.hxx>
62 #include <Geom_TrimmedCurve.hxx>
63 #include <IntAna_IntConicQuad.hxx>
64 #include <IntAna_IntLinTorus.hxx>
65 #include <IntAna_Quadric.hxx>
66 #include <IntCurveSurface_TransitionOnCurve.hxx>
67 #include <IntCurvesFace_Intersector.hxx>
68 #include <Poly_Triangulation.hxx>
69 #include <Precision.hxx>
71 #include <TopExp_Explorer.hxx>
72 #include <TopLoc_Location.hxx>
73 #include <TopTools_MapOfShape.hxx>
75 #include <TopoDS_Compound.hxx>
76 #include <TopoDS_Face.hxx>
77 #include <TopoDS_TShape.hxx>
78 #include <gp_Cone.hxx>
79 #include <gp_Cylinder.hxx>
82 #include <gp_Pnt2d.hxx>
83 #include <gp_Sphere.hxx>
84 #include <gp_Torus.hxx>
88 #include <tbb/parallel_for.h>
89 //#include <tbb/enumerable_thread_specific.h>
98 #if OCC_VERSION_LARGE <= 0x06050300
99 // workaround is required only for OCCT6.5.3 and older (see OCC22809)
100 #define ELLIPSOLID_WORKAROUND
103 #ifdef ELLIPSOLID_WORKAROUND
104 #include <BRepIntCurveSurface_Inter.hxx>
105 #include <BRepTopAdaptor_TopolTool.hxx>
106 #include <BRepAdaptor_HSurface.hxx>
109 //=============================================================================
113 //=============================================================================
115 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
116 :SMESH_3D_Algo(hypId, studyId, gen)
118 _name = "Cartesian_3D";
119 _shapeType = (1 << TopAbs_SOLID); // 1 bit /shape type
120 _compatibleHypothesis.push_back("CartesianParameters3D");
122 _onlyUnaryInput = false; // to mesh all SOLIDs at once
123 _requireDiscreteBoundary = false; // 2D mesh not needed
124 _supportSubmeshes = false; // do not use any existing mesh
127 //=============================================================================
129 * Check presence of a hypothesis
131 //=============================================================================
133 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh& aMesh,
134 const TopoDS_Shape& aShape,
135 Hypothesis_Status& aStatus)
137 aStatus = SMESH_Hypothesis::HYP_MISSING;
139 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
140 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
141 if ( h == hyps.end())
146 for ( ; h != hyps.end(); ++h )
148 if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
150 aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
155 return aStatus == HYP_OK;
162 //=============================================================================
163 // Definitions of internal utils
164 // --------------------------------------------------------------------------
166 Trans_TANGENT = IntCurveSurface_Tangent,
167 Trans_IN = IntCurveSurface_In,
168 Trans_OUT = IntCurveSurface_Out,
171 // --------------------------------------------------------------------------
173 * \brief Common data of any intersection between a Grid and a shape
175 struct B_IntersectPoint
177 mutable const SMDS_MeshNode* _node;
178 mutable vector< TGeomID > _faceIDs;
180 B_IntersectPoint(): _node(NULL) {}
181 void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
182 int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
183 bool IsOnFace( int faceID ) const;
184 virtual ~B_IntersectPoint() {}
186 // --------------------------------------------------------------------------
188 * \brief Data of intersection between a GridLine and a TopoDS_Face
190 struct F_IntersectPoint : public B_IntersectPoint
193 mutable Transition _transition;
194 mutable size_t _indexOnLine;
196 bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
198 // --------------------------------------------------------------------------
200 * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
202 struct E_IntersectPoint : public B_IntersectPoint
208 // --------------------------------------------------------------------------
210 * \brief A line of the grid and its intersections with 2D geometry
215 double _length; // line length
216 multiset< F_IntersectPoint > _intPoints;
218 void RemoveExcessIntPoints( const double tol );
219 bool GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut );
221 // --------------------------------------------------------------------------
223 * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
227 gp_XYZ _uNorm, _vNorm, _zNorm;
228 vector< gp_XYZ > _origins; // origin points of all planes in one direction
229 vector< double > _zProjs; // projections of origins to _zNorm
231 // --------------------------------------------------------------------------
233 * \brief Iterator on the parallel grid lines of one direction
239 size_t _iVar1, _iVar2, _iConst;
240 string _name1, _name2, _nameConst;
242 LineIndexer( size_t sz1, size_t sz2, size_t sz3,
243 size_t iv1, size_t iv2, size_t iConst,
244 const string& nv1, const string& nv2, const string& nConst )
246 _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
247 _curInd[0] = _curInd[1] = _curInd[2] = 0;
248 _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
249 _name1 = nv1; _name2 = nv2; _nameConst = nConst;
252 size_t I() const { return _curInd[0]; }
253 size_t J() const { return _curInd[1]; }
254 size_t K() const { return _curInd[2]; }
255 void SetIJK( size_t i, size_t j, size_t k )
257 _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
261 if ( ++_curInd[_iVar1] == _size[_iVar1] )
262 _curInd[_iVar1] = 0, ++_curInd[_iVar2];
264 bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
265 size_t LineIndex () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
266 size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
267 size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
268 size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
269 void SetIndexOnLine (size_t i) { _curInd[ _iConst ] = i; }
270 size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
272 // --------------------------------------------------------------------------
274 * \brief Container of GridLine's
278 vector< double > _coords[3]; // coordinates of grid nodes
279 gp_XYZ _axes [3]; // axis directions
280 vector< GridLine > _lines [3]; // in 3 directions
281 double _tol, _minCellSize;
283 gp_Mat _invB; // inverted basis of _axes
284 //bool _isOrthogonalAxes;
286 vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
287 vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
289 list< E_IntersectPoint > _edgeIntP; // intersections with EDGEs
290 TopTools_IndexedMapOfShape _shapes;
292 size_t CellIndex( size_t i, size_t j, size_t k ) const
294 return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
296 size_t NodeIndex( size_t i, size_t j, size_t k ) const
298 return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
300 size_t NodeIndexDX() const { return 1; }
301 size_t NodeIndexDY() const { return _coords[0].size(); }
302 size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
304 LineIndexer GetLineIndexer(size_t iDir) const;
306 void SetCoordinates(const vector<double>& xCoords,
307 const vector<double>& yCoords,
308 const vector<double>& zCoords,
309 const double* axesDirs,
310 const Bnd_Box& bndBox );
311 void ComputeUVW(const gp_XYZ& p, double uvw[3]);
312 void ComputeNodes(SMESH_MesherHelper& helper);
314 #ifdef ELLIPSOLID_WORKAROUND
315 // --------------------------------------------------------------------------
317 * \brief struct temporary replacing IntCurvesFace_Intersector until
318 * OCCT bug 0022809 is fixed
319 * http://tracker.dev.opencascade.org/view.php?id=22809
321 struct TMP_IntCurvesFace_Intersector
323 BRepAdaptor_Surface _surf;
325 BRepIntCurveSurface_Inter _intcs;
326 vector<IntCurveSurface_IntersectionPoint> _points;
327 BRepTopAdaptor_TopolTool _clsf;
329 TMP_IntCurvesFace_Intersector(const TopoDS_Face& face, const double tol)
330 :_surf( face ), _tol( tol ), _clsf( new BRepAdaptor_HSurface(_surf) ) {}
331 Bnd_Box Bounding() const { Bnd_Box b; BRepBndLib::Add (_surf.Face(), b); return b; }
332 void Perform( const gp_Lin& line, const double w0, const double w1 )
335 for ( _intcs.Init( _surf.Face(), line, _tol ); _intcs.More(); _intcs.Next() )
336 if ( w0 <= _intcs.W() && _intcs.W() <= w1 )
337 _points.push_back( _intcs.Point() );
339 bool IsDone() const { return true; }
340 int NbPnt() const { return _points.size(); }
341 IntCurveSurface_TransitionOnCurve Transition( const int i ) const { return _points[ i-1 ].Transition(); }
342 double WParameter( const int i ) const { return _points[ i-1 ].W(); }
343 TopAbs_State ClassifyUVPoint(const gp_Pnt2d& p) { return _clsf.Classify( p, _tol ); }
345 #define __IntCurvesFace_Intersector TMP_IntCurvesFace_Intersector
347 #define __IntCurvesFace_Intersector IntCurvesFace_Intersector
349 // --------------------------------------------------------------------------
351 * \brief Intersector of TopoDS_Face with all GridLine's
353 struct FaceGridIntersector
359 __IntCurvesFace_Intersector* _surfaceInt;
360 vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
362 FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
364 bool IsInGrid(const Bnd_Box& gridBox);
366 void StoreIntersections()
368 for ( size_t i = 0; i < _intersections.size(); ++i )
370 multiset< F_IntersectPoint >::iterator ip =
371 _intersections[i].first->_intPoints.insert( _intersections[i].second );
372 ip->_faceIDs.reserve( 1 );
373 ip->_faceIDs.push_back( _faceID );
376 const Bnd_Box& GetFaceBndBox()
378 GetCurveFaceIntersector();
381 __IntCurvesFace_Intersector* GetCurveFaceIntersector()
385 _surfaceInt = new __IntCurvesFace_Intersector( _face, Precision::PConfusion() );
386 _bndBox = _surfaceInt->Bounding();
387 if ( _bndBox.IsVoid() )
388 BRepBndLib::Add (_face, _bndBox);
392 bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
394 // --------------------------------------------------------------------------
396 * \brief Intersector of a surface with a GridLine
398 struct FaceLineIntersector
401 double _u, _v, _w; // params on the face and the line
402 Transition _transition; // transition of at intersection (see IntCurveSurface.cdl)
403 Transition _transIn, _transOut; // IN and OUT transitions depending of face orientation
406 gp_Cylinder _cylinder;
410 __IntCurvesFace_Intersector* _surfaceInt;
412 vector< F_IntersectPoint > _intPoints;
414 void IntersectWithPlane (const GridLine& gridLine);
415 void IntersectWithCylinder(const GridLine& gridLine);
416 void IntersectWithCone (const GridLine& gridLine);
417 void IntersectWithSphere (const GridLine& gridLine);
418 void IntersectWithTorus (const GridLine& gridLine);
419 void IntersectWithSurface (const GridLine& gridLine);
421 bool UVIsOnFace() const;
422 void addIntPoint(const bool toClassify=true);
423 bool isParamOnLineOK( const double linLength )
425 return -_tol < _w && _w < linLength + _tol;
427 FaceLineIntersector():_surfaceInt(0) {}
428 ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
430 // --------------------------------------------------------------------------
432 * \brief Class representing topology of the hexahedron and creating a mesh
433 * volume basing on analysis of hexahedron intersection with geometry
437 // --------------------------------------------------------------------------------
440 // --------------------------------------------------------------------------------
441 struct _Node //!< node either at a hexahedron corner or at intersection
443 const SMDS_MeshNode* _node; // mesh node at hexahedron corner
444 const B_IntersectPoint* _intPoint;
447 _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
448 :_node(n), _intPoint(ip), _isUsedInFace(0) {}
449 const SMDS_MeshNode* Node() const
450 { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
451 const F_IntersectPoint* FaceIntPnt() const
452 { return static_cast< const F_IntersectPoint* >( _intPoint ); }
453 const E_IntersectPoint* EdgeIntPnt() const
454 { return static_cast< const E_IntersectPoint* >( _intPoint ); }
455 void Add( const E_IntersectPoint* ip )
460 else if ( !_intPoint->_node ) {
461 ip->Add( _intPoint->_faceIDs );
465 _intPoint->Add( ip->_faceIDs );
468 int IsLinked( const B_IntersectPoint* other,
469 int avoidFace=-1 ) const // returns id of a common face
471 return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
473 bool IsOnFace( int faceID ) const // returns true if faceID is found
475 return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
479 if ( const SMDS_MeshNode* n = Node() )
480 return SMESH_TNodeXYZ( n );
481 if ( const E_IntersectPoint* eip =
482 dynamic_cast< const E_IntersectPoint* >( _intPoint ))
484 return gp_Pnt( 1e100, 0, 0 );
487 // --------------------------------------------------------------------------------
488 struct _Link // link connecting two _Node's
491 vector< _Node > _intNodes; // _Node's at GridLine intersections
492 vector< _Link > _splits;
493 vector< _Face*> _faces;
495 // --------------------------------------------------------------------------------
500 _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
501 void Reverse() { _reverse = !_reverse; }
502 int NbResultLinks() const { return _link->_splits.size(); }
503 _OrientedLink ResultLink(int i) const
505 return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
507 _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
508 _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
509 operator bool() const { return _link; }
510 vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
512 vector< TGeomID > faces;
513 const B_IntersectPoint *ip0, *ip1;
514 if (( ip0 = _link->_nodes[0]->_intPoint ) &&
515 ( ip1 = _link->_nodes[1]->_intPoint ))
517 for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
518 if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
519 !usedIDs.count( ip0->_faceIDs[i] ) )
520 faces.push_back( ip0->_faceIDs[i] );
524 bool HasEdgeNodes() const
526 return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
527 dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
530 // --------------------------------------------------------------------------------
533 vector< _OrientedLink > _links; // links on GridLine's
534 vector< _Link > _polyLinks; // links added to close a polygonal face
535 vector< _Node > _edgeNodes; // nodes at intersection with EDGEs
537 // --------------------------------------------------------------------------------
538 struct _volumeDef // holder of nodes of a volume mesh element
540 //vector< const SMDS_MeshNode* > _nodes;
541 vector< _Node* > _nodes;
542 vector< int > _quantities;
543 typedef boost::shared_ptr<_volumeDef> Ptr;
544 void set( const vector< _Node* >& nodes,
545 const vector< int >& quant = vector< int >() )
546 { _nodes = nodes; _quantities = quant; }
547 // static Ptr New( const vector< const SMDS_MeshNode* >& nodes,
548 // const vector< int > quant = vector< int >() )
550 // _volumeDef* def = new _volumeDef;
551 // def->_nodes = nodes;
552 // def->_quantities = quant;
553 // return Ptr( def );
557 // topology of a hexahedron
560 _Link _hexLinks [12];
563 // faces resulted from hexahedron intersection
564 vector< _Face > _polygons;
566 // intresections with EDGEs
567 vector< const E_IntersectPoint* > _edgeIntPnts;
569 // nodes inside the hexahedron (at VERTEXes)
570 vector< _Node > _vertexNodes;
572 // computed volume elements
573 //vector< _volumeDef::Ptr > _volumeDefs;
574 _volumeDef _volumeDefs;
577 double _sizeThreshold, _sideLength[3];
578 int _nbCornerNodes, _nbIntNodes, _nbBndNodes;
579 int _origNodeInd; // index of _hexNodes[0] node within the _grid
583 Hexahedron(const double sizeThreshold, Grid* grid);
584 int MakeElements(SMESH_MesherHelper& helper,
585 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
586 void ComputeElements();
587 void Init() { init( _i, _j, _k ); }
590 Hexahedron(const Hexahedron& other );
591 void init( size_t i, size_t j, size_t k );
592 void init( size_t i );
593 void addEdges(SMESH_MesherHelper& helper,
594 vector< Hexahedron* >& intersectedHex,
595 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
596 gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
597 double proj, BRepAdaptor_Curve& curve,
598 const gp_XYZ& axis, const gp_XYZ& origin );
599 int getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
600 bool addIntersection( const E_IntersectPoint& ip,
601 vector< Hexahedron* >& hexes,
602 int ijk[], int dIJK[] );
603 bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
604 bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
605 int addElements(SMESH_MesherHelper& helper);
606 bool is1stNodeOut( int iLink ) const;
607 bool isInHole() const;
608 bool checkPolyhedronSize() const;
613 bool debugDumpLink( _Link* link );
614 _Node* FindEqualNode( vector< _Node >& nodes,
615 const E_IntersectPoint* ip,
618 for ( size_t i = 0; i < nodes.size(); ++i )
619 if ( nodes[i].Point().SquareDistance( ip->_point ) <= tol2 )
626 // --------------------------------------------------------------------------
628 * \brief Hexahedron computing volumes in one thread
630 struct ParallelHexahedron
632 vector< Hexahedron* >& _hexVec;
634 ParallelHexahedron( vector< Hexahedron* >& hv, vector<int>& ind): _hexVec(hv), _index(ind) {}
635 void operator() ( const tbb::blocked_range<size_t>& r ) const
637 for ( size_t i = r.begin(); i != r.end(); ++i )
638 if ( Hexahedron* hex = _hexVec[ _index[i]] )
639 hex->ComputeElements();
642 // --------------------------------------------------------------------------
644 * \brief Structure intersecting certain nb of faces with GridLine's in one thread
646 struct ParallelIntersector
648 vector< FaceGridIntersector >& _faceVec;
649 ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
650 void operator() ( const tbb::blocked_range<size_t>& r ) const
652 for ( size_t i = r.begin(); i != r.end(); ++i )
653 _faceVec[i].Intersect();
658 //=============================================================================
659 // Implementation of internal utils
660 //=============================================================================
662 * \brief adjust \a i to have \a val between values[i] and values[i+1]
664 inline void locateValue( int & i, double val, const vector<double>& values,
665 int& di, double tol )
667 //val += values[0]; // input \a val is measured from 0.
668 if ( i > values.size()-2 )
671 while ( i+2 < values.size() && val > values[ i+1 ])
673 while ( i > 0 && val < values[ i ])
676 if ( i > 0 && val - values[ i ] < tol )
678 else if ( i+2 < values.size() && values[ i+1 ] - val < tol )
683 //=============================================================================
685 * Remove coincident intersection points
687 void GridLine::RemoveExcessIntPoints( const double tol )
689 if ( _intPoints.size() < 2 ) return;
691 set< Transition > tranSet;
692 multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
693 while ( ip2 != _intPoints.end() )
697 while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
699 tranSet.insert( ip1->_transition );
700 tranSet.insert( ip2->_transition );
701 ip2->Add( ip1->_faceIDs );
702 _intPoints.erase( ip1 );
705 if ( tranSet.size() > 1 ) // points with different transition coincide
707 bool isIN = tranSet.count( Trans_IN );
708 bool isOUT = tranSet.count( Trans_OUT );
710 (*ip1)._transition = Trans_TANGENT;
712 (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
716 //================================================================================
718 * Return "is OUT" state for nodes before the given intersection point
720 bool GridLine::GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut )
722 if ( ip->_transition == Trans_IN )
724 if ( ip->_transition == Trans_OUT )
726 if ( ip->_transition == Trans_APEX )
728 // singularity point (apex of a cone)
729 if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
731 multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
732 if ( ipAft == _intPoints.end() )
735 if ( ipBef->_transition != ipAft->_transition )
736 return ( ipBef->_transition == Trans_OUT );
737 return ( ipBef->_transition != Trans_OUT );
739 // _transition == Trans_TANGENT
742 //================================================================================
746 void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
747 const SMDS_MeshNode* n) const
749 if ( _faceIDs.empty() )
752 for ( size_t i = 0; i < fIDs.size(); ++i )
754 vector< TGeomID >::iterator it =
755 std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
756 if ( it == _faceIDs.end() )
757 _faceIDs.push_back( fIDs[i] );
762 //================================================================================
764 * Returns index of a common face if any, else zero
766 int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
769 for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
770 if ( avoidFace != other->_faceIDs[i] &&
771 IsOnFace ( other->_faceIDs[i] ))
772 return other->_faceIDs[i];
775 //================================================================================
777 * Returns \c true if \a faceID in in this->_faceIDs
779 bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
781 vector< TGeomID >::const_iterator it =
782 std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
783 return ( it != _faceIDs.end() );
785 //================================================================================
787 * Return an iterator on GridLine's in a given direction
789 LineIndexer Grid::GetLineIndexer(size_t iDir) const
791 const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
792 const string s [] = { "X", "Y", "Z" };
793 LineIndexer li( _coords[0].size(), _coords[1].size(), _coords[2].size(),
794 indices[iDir*3], indices[iDir*3+1], indices[iDir*3+2],
795 s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
798 //=============================================================================
800 * Creates GridLine's of the grid
802 void Grid::SetCoordinates(const vector<double>& xCoords,
803 const vector<double>& yCoords,
804 const vector<double>& zCoords,
805 const double* axesDirs,
806 const Bnd_Box& shapeBox)
808 _coords[0] = xCoords;
809 _coords[1] = yCoords;
810 _coords[2] = zCoords;
812 _axes[0].SetCoord( axesDirs[0],
815 _axes[1].SetCoord( axesDirs[3],
818 _axes[2].SetCoord( axesDirs[6],
821 _axes[0].Normalize();
822 _axes[1].Normalize();
823 _axes[2].Normalize();
825 _invB.SetCols( _axes[0], _axes[1], _axes[2] );
828 // _isOrthogonalAxes = ( Abs( _axes[0] * _axes[1] ) < 1e-20 &&
829 // Abs( _axes[1] * _axes[2] ) < 1e-20 &&
830 // Abs( _axes[2] * _axes[0] ) < 1e-20 );
833 _minCellSize = Precision::Infinite();
834 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
836 for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
838 double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
839 if ( cellLen < _minCellSize )
840 _minCellSize = cellLen;
843 if ( _minCellSize < Precision::Confusion() )
844 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
845 SMESH_Comment("Too small cell size: ") << _minCellSize );
846 _tol = _minCellSize / 1000.;
848 // attune grid extremities to shape bounding box
850 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
851 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
852 double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
853 &_coords[0].back(), &_coords[1].back(), &_coords[2].back() };
854 for ( int i = 0; i < 6; ++i )
855 if ( fabs( sP[i] - *cP[i] ) < _tol )
856 *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
858 for ( int iDir = 0; iDir < 3; ++iDir )
860 if ( _coords[iDir][0] - sP[iDir] > _tol )
862 _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
863 _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
865 if ( sP[iDir+3] - _coords[iDir].back() > _tol )
867 _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
868 _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
871 _tol = _minCellSize / 1000.;
873 _origin = ( _coords[0][0] * _axes[0] +
874 _coords[1][0] * _axes[1] +
875 _coords[2][0] * _axes[2] );
878 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
880 LineIndexer li = GetLineIndexer( iDir );
881 _lines[iDir].resize( li.NbLines() );
882 double len = _coords[ iDir ].back() - _coords[iDir].front();
883 for ( ; li.More(); ++li )
885 GridLine& gl = _lines[iDir][ li.LineIndex() ];
886 gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
887 _coords[1][li.J()] * _axes[1] +
888 _coords[2][li.K()] * _axes[2] );
889 gl._line.SetDirection( _axes[ iDir ]);
894 //================================================================================
896 * Computes coordinates of a point in the grid CS
898 void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
900 // gp_XYZ p = P - _origin;
901 // UVW[ 0 ] = p.X() * _invB( 1, 1 ) + p.Y() * _invB( 1, 2 ) + p.Z() * _invB( 1, 3 );
902 // UVW[ 1 ] = p.X() * _invB( 2, 1 ) + p.Y() * _invB( 2, 2 ) + p.Z() * _invB( 2, 3 );
903 // UVW[ 2 ] = p.X() * _invB( 3, 1 ) + p.Y() * _invB( 3, 2 ) + p.Z() * _invB( 3, 3 );
904 // UVW[ 0 ] += _coords[0][0];
905 // UVW[ 1 ] += _coords[1][0];
906 // UVW[ 2 ] += _coords[2][0];
907 gp_XYZ p = P * _invB;
908 p.Coord( UVW[0], UVW[1], UVW[2] );
910 //================================================================================
914 void Grid::ComputeNodes(SMESH_MesherHelper& helper)
916 // state of each node of the grid relative to the geometry
917 const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
918 vector< bool > isNodeOut( nbGridNodes, false );
919 _nodes.resize( nbGridNodes, 0 );
920 _gridIntP.resize( nbGridNodes, NULL );
922 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
924 LineIndexer li = GetLineIndexer( iDir );
926 // find out a shift of node index while walking along a GridLine in this direction
927 li.SetIndexOnLine( 0 );
928 size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
929 li.SetIndexOnLine( 1 );
930 const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
932 const vector<double> & coords = _coords[ iDir ];
933 for ( ; li.More(); ++li ) // loop on lines in iDir
935 li.SetIndexOnLine( 0 );
936 nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
938 GridLine& line = _lines[ iDir ][ li.LineIndex() ];
939 const gp_XYZ lineLoc = line._line.Location().XYZ();
940 const gp_XYZ lineDir = line._line.Direction().XYZ();
941 line.RemoveExcessIntPoints( _tol );
942 multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
943 multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
946 const double* nodeCoord = & coords[0];
947 const double* coord0 = nodeCoord;
948 const double* coordEnd = coord0 + coords.size();
949 double nodeParam = 0;
950 for ( ; ip != intPnts.end(); ++ip )
952 // set OUT state or just skip IN nodes before ip
953 if ( nodeParam < ip->_paramOnLine - _tol )
955 isOut = line.GetIsOutBefore( ip, isOut );
957 while ( nodeParam < ip->_paramOnLine - _tol )
960 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
961 if ( ++nodeCoord < coordEnd )
962 nodeParam = *nodeCoord - *coord0;
966 if ( nodeCoord == coordEnd ) break;
968 // create a mesh node on a GridLine at ip if it does not coincide with a grid node
969 if ( nodeParam > ip->_paramOnLine + _tol )
971 // li.SetIndexOnLine( 0 );
972 // double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
973 // xyz[ li._iConst ] += ip->_paramOnLine;
974 gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
975 ip->_node = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
976 ip->_indexOnLine = nodeCoord-coord0-1;
978 // create a mesh node at ip concident with a grid node
981 int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
982 if ( !_nodes[ nodeIndex ] )
984 //li.SetIndexOnLine( nodeCoord-coord0 );
985 //double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
986 gp_XYZ xyz = lineLoc + nodeParam * lineDir;
987 _nodes [ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
988 _gridIntP[ nodeIndex ] = & * ip;
990 if ( _gridIntP[ nodeIndex ] )
991 _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
993 _gridIntP[ nodeIndex ] = & * ip;
994 // ip->_node = _nodes[ nodeIndex ]; -- to differ from ip on links
995 ip->_indexOnLine = nodeCoord-coord0;
996 if ( ++nodeCoord < coordEnd )
997 nodeParam = *nodeCoord - *coord0;
1000 // set OUT state to nodes after the last ip
1001 for ( ; nodeCoord < coordEnd; ++nodeCoord )
1002 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
1006 // Create mesh nodes at !OUT nodes of the grid
1008 for ( size_t z = 0; z < _coords[2].size(); ++z )
1009 for ( size_t y = 0; y < _coords[1].size(); ++y )
1010 for ( size_t x = 0; x < _coords[0].size(); ++x )
1012 size_t nodeIndex = NodeIndex( x, y, z );
1013 if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
1015 //_nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
1016 gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1017 _coords[1][y] * _axes[1] +
1018 _coords[2][z] * _axes[2] );
1019 _nodes[ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1024 // check validity of transitions
1025 const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1026 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1028 LineIndexer li = GetLineIndexer( iDir );
1029 for ( ; li.More(); ++li )
1031 multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1032 if ( intPnts.empty() ) continue;
1033 if ( intPnts.size() == 1 )
1035 if ( intPnts.begin()->_transition != Trans_TANGENT &&
1036 intPnts.begin()->_transition != Trans_APEX )
1037 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1038 SMESH_Comment("Wrong SOLE transition of GridLine (")
1039 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1040 << ") along " << li._nameConst
1041 << ": " << trName[ intPnts.begin()->_transition] );
1045 if ( intPnts.begin()->_transition == Trans_OUT )
1046 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1047 SMESH_Comment("Wrong START transition of GridLine (")
1048 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1049 << ") along " << li._nameConst
1050 << ": " << trName[ intPnts.begin()->_transition ]);
1051 if ( intPnts.rbegin()->_transition == Trans_IN )
1052 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1053 SMESH_Comment("Wrong END transition of GridLine (")
1054 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1055 << ") along " << li._nameConst
1056 << ": " << trName[ intPnts.rbegin()->_transition ]);
1063 //=============================================================================
1065 * Checks if the face is encosed by the grid
1067 bool FaceGridIntersector::IsInGrid(const Bnd_Box& gridBox)
1069 // double x0,y0,z0, x1,y1,z1;
1070 // const Bnd_Box& faceBox = GetFaceBndBox();
1071 // faceBox.Get(x0,y0,z0, x1,y1,z1);
1073 // if ( !gridBox.IsOut( gp_Pnt( x0,y0,z0 )) &&
1074 // !gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
1077 // double X0,Y0,Z0, X1,Y1,Z1;
1078 // gridBox.Get(X0,Y0,Z0, X1,Y1,Z1);
1079 // double faceP[6] = { x0,y0,z0, x1,y1,z1 };
1080 // double gridP[6] = { X0,Y0,Z0, X1,Y1,Z1 };
1081 // gp_Dir axes[3] = { gp::DX(), gp::DY(), gp::DZ() };
1082 // for ( int iDir = 0; iDir < 6; ++iDir )
1084 // if ( iDir < 3 && gridP[ iDir ] <= faceP[ iDir ] ) continue;
1085 // if ( iDir >= 3 && gridP[ iDir ] >= faceP[ iDir ] ) continue;
1087 // // check if the face intersects a side of a gridBox
1089 // gp_Pnt p = iDir < 3 ? gp_Pnt( X0,Y0,Z0 ) : gp_Pnt( X1,Y1,Z1 );
1090 // gp_Ax1 norm( p, axes[ iDir % 3 ] );
1091 // if ( iDir < 3 ) norm.Reverse();
1093 // gp_XYZ O = norm.Location().XYZ(), N = norm.Direction().XYZ();
1095 // TopLoc_Location loc = _face.Location();
1096 // Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(_face,loc);
1097 // if ( !aPoly.IsNull() )
1099 // if ( !loc.IsIdentity() )
1101 // norm.Transform( loc.Transformation().Inverted() );
1102 // O = norm.Location().XYZ(), N = norm.Direction().XYZ();
1104 // const double deflection = aPoly->Deflection();
1106 // const TColgp_Array1OfPnt& nodes = aPoly->Nodes();
1107 // for ( int i = nodes.Lower(); i <= nodes.Upper(); ++i )
1108 // if (( nodes( i ).XYZ() - O ) * N > _grid->_tol + deflection )
1113 // BRepAdaptor_Surface surf( _face );
1114 // double u0, u1, v0, v1, du, dv, u, v;
1115 // BRepTools::UVBounds( _face, u0, u1, v0, v1);
1116 // if ( surf.GetType() == GeomAbs_Plane ) {
1117 // du = u1 - u0, dv = v1 - v0;
1120 // du = surf.UResolution( _grid->_minCellSize / 10. );
1121 // dv = surf.VResolution( _grid->_minCellSize / 10. );
1123 // for ( u = u0, v = v0; u <= u1 && v <= v1; u += du, v += dv )
1125 // gp_Pnt p = surf.Value( u, v );
1126 // if (( p.XYZ() - O ) * N > _grid->_tol )
1128 // TopAbs_State state = GetCurveFaceIntersector()->ClassifyUVPoint(gp_Pnt2d( u, v ));
1129 // if ( state == TopAbs_IN || state == TopAbs_ON )
1137 //=============================================================================
1139 * Intersects TopoDS_Face with all GridLine's
1141 void FaceGridIntersector::Intersect()
1143 FaceLineIntersector intersector;
1144 intersector._surfaceInt = GetCurveFaceIntersector();
1145 intersector._tol = _grid->_tol;
1146 intersector._transOut = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1147 intersector._transIn = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1149 typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1150 PIntFun interFunction;
1152 BRepAdaptor_Surface surf( _face );
1153 switch ( surf.GetType() ) {
1155 intersector._plane = surf.Plane();
1156 interFunction = &FaceLineIntersector::IntersectWithPlane;
1158 case GeomAbs_Cylinder:
1159 intersector._cylinder = surf.Cylinder();
1160 interFunction = &FaceLineIntersector::IntersectWithCylinder;
1163 intersector._cone = surf.Cone();
1164 interFunction = &FaceLineIntersector::IntersectWithCone;
1166 case GeomAbs_Sphere:
1167 intersector._sphere = surf.Sphere();
1168 interFunction = &FaceLineIntersector::IntersectWithSphere;
1171 intersector._torus = surf.Torus();
1172 interFunction = &FaceLineIntersector::IntersectWithTorus;
1175 interFunction = &FaceLineIntersector::IntersectWithSurface;
1178 _intersections.clear();
1179 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1181 if ( surf.GetType() == GeomAbs_Plane )
1183 // check if all lines in this direction are parallel to a plane
1184 if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1185 Precision::Angular()))
1187 // find out a transition, that is the same for all lines of a direction
1188 gp_Dir plnNorm = intersector._plane.Axis().Direction();
1189 gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1190 intersector._transition =
1191 ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1193 if ( surf.GetType() == GeomAbs_Cylinder )
1195 // check if all lines in this direction are parallel to a cylinder
1196 if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1197 Precision::Angular()))
1201 // intersect the grid lines with the face
1202 for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1204 GridLine& gridLine = _grid->_lines[iDir][iL];
1205 if ( _bndBox.IsOut( gridLine._line )) continue;
1207 intersector._intPoints.clear();
1208 (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1209 for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1210 _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1214 //================================================================================
1216 * Return true if (_u,_v) is on the face
1218 bool FaceLineIntersector::UVIsOnFace() const
1220 TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1221 return ( state == TopAbs_IN || state == TopAbs_ON );
1223 //================================================================================
1225 * Store an intersection if it is IN or ON the face
1227 void FaceLineIntersector::addIntPoint(const bool toClassify)
1229 if ( !toClassify || UVIsOnFace() )
1232 p._paramOnLine = _w;
1233 p._transition = _transition;
1234 _intPoints.push_back( p );
1237 //================================================================================
1239 * Intersect a line with a plane
1241 void FaceLineIntersector::IntersectWithPlane (const GridLine& gridLine)
1243 IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1244 _w = linPlane.ParamOnConic(1);
1245 if ( isParamOnLineOK( gridLine._length ))
1247 ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1251 //================================================================================
1253 * Intersect a line with a cylinder
1255 void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1257 IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1258 if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1260 _w = linCylinder.ParamOnConic(1);
1261 if ( linCylinder.NbPoints() == 1 )
1262 _transition = Trans_TANGENT;
1264 _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1265 if ( isParamOnLineOK( gridLine._length ))
1267 ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1270 if ( linCylinder.NbPoints() > 1 )
1272 _w = linCylinder.ParamOnConic(2);
1273 if ( isParamOnLineOK( gridLine._length ))
1275 ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1276 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1282 //================================================================================
1284 * Intersect a line with a cone
1286 void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1288 IntAna_IntConicQuad linCone(gridLine._line,_cone);
1289 if ( !linCone.IsDone() ) return;
1291 gp_Vec du, dv, norm;
1292 for ( int i = 1; i <= linCone.NbPoints(); ++i )
1294 _w = linCone.ParamOnConic( i );
1295 if ( !isParamOnLineOK( gridLine._length )) continue;
1296 ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1299 ElSLib::D1( _u, _v, _cone, P, du, dv );
1301 double normSize2 = norm.SquareMagnitude();
1302 if ( normSize2 > Precision::Angular() * Precision::Angular() )
1304 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1305 cos /= sqrt( normSize2 );
1306 if ( cos < -Precision::Angular() )
1307 _transition = _transIn;
1308 else if ( cos > Precision::Angular() )
1309 _transition = _transOut;
1311 _transition = Trans_TANGENT;
1315 _transition = Trans_APEX;
1317 addIntPoint( /*toClassify=*/false);
1321 //================================================================================
1323 * Intersect a line with a sphere
1325 void FaceLineIntersector::IntersectWithSphere (const GridLine& gridLine)
1327 IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1328 if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1330 _w = linSphere.ParamOnConic(1);
1331 if ( linSphere.NbPoints() == 1 )
1332 _transition = Trans_TANGENT;
1334 _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1335 if ( isParamOnLineOK( gridLine._length ))
1337 ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1340 if ( linSphere.NbPoints() > 1 )
1342 _w = linSphere.ParamOnConic(2);
1343 if ( isParamOnLineOK( gridLine._length ))
1345 ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1346 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1352 //================================================================================
1354 * Intersect a line with a torus
1356 void FaceLineIntersector::IntersectWithTorus (const GridLine& gridLine)
1358 IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1359 if ( !linTorus.IsDone()) return;
1361 gp_Vec du, dv, norm;
1362 for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1364 _w = linTorus.ParamOnLine( i );
1365 if ( !isParamOnLineOK( gridLine._length )) continue;
1366 linTorus.ParamOnTorus( i, _u,_v );
1369 ElSLib::D1( _u, _v, _torus, P, du, dv );
1371 double normSize = norm.Magnitude();
1372 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1374 if ( cos < -Precision::Angular() )
1375 _transition = _transIn;
1376 else if ( cos > Precision::Angular() )
1377 _transition = _transOut;
1379 _transition = Trans_TANGENT;
1380 addIntPoint( /*toClassify=*/false);
1384 //================================================================================
1386 * Intersect a line with a non-analytical surface
1388 void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1390 _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1391 if ( !_surfaceInt->IsDone() ) return;
1392 for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1394 _transition = Transition( _surfaceInt->Transition( i ) );
1395 _w = _surfaceInt->WParameter( i );
1396 addIntPoint(/*toClassify=*/false);
1399 //================================================================================
1401 * check if its face can be safely intersected in a thread
1403 bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1408 TopLoc_Location loc;
1409 Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1410 Handle(Geom_RectangularTrimmedSurface) ts =
1411 Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1412 while( !ts.IsNull() ) {
1413 surf = ts->BasisSurface();
1414 ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1416 if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1417 surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1418 if ( !noSafeTShapes.insert((const Standard_Transient*) _face.TShape() ).second )
1422 TopExp_Explorer exp( _face, TopAbs_EDGE );
1423 for ( ; exp.More(); exp.Next() )
1425 bool edgeIsSafe = true;
1426 const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1429 Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1432 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1433 while( !tc.IsNull() ) {
1434 c = tc->BasisCurve();
1435 tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1437 if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1438 c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1445 Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1448 Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1449 while( !tc.IsNull() ) {
1450 c2 = tc->BasisCurve();
1451 tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1453 if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1454 c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1458 if ( !edgeIsSafe && !noSafeTShapes.insert((const Standard_Transient*) e.TShape() ).second )
1463 //================================================================================
1465 * \brief Creates topology of the hexahedron
1467 Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1468 : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbIntNodes(0)
1470 _polygons.reserve(100); // to avoid reallocation;
1472 //set nodes shift within grid->_nodes from the node 000
1473 size_t dx = _grid->NodeIndexDX();
1474 size_t dy = _grid->NodeIndexDY();
1475 size_t dz = _grid->NodeIndexDZ();
1477 size_t i100 = i000 + dx;
1478 size_t i010 = i000 + dy;
1479 size_t i110 = i010 + dx;
1480 size_t i001 = i000 + dz;
1481 size_t i101 = i100 + dz;
1482 size_t i011 = i010 + dz;
1483 size_t i111 = i110 + dz;
1484 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1485 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1486 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1487 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1488 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1489 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1490 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1491 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1493 vector< int > idVec;
1494 // set nodes to links
1495 for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1497 SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1498 _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1499 link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1500 link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1501 link._intNodes.reserve( 10 ); // to avoid reallocation
1502 link._splits.reserve( 10 );
1505 // set links to faces
1506 int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1507 for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1509 SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1510 _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1511 bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1512 faceID == SMESH_Block::ID_Fx1z ||
1513 faceID == SMESH_Block::ID_F0yz );
1514 quad._links.resize(4);
1515 vector<_OrientedLink>::iterator frwLinkIt = quad._links.begin();
1516 vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1517 for ( int i = 0; i < 4; ++i )
1519 bool revLink = revFace;
1520 if ( i > 1 ) // reverse links u1 and v0
1522 _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1523 link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1528 //================================================================================
1530 * \brief Copy constructor
1532 Hexahedron::Hexahedron( const Hexahedron& other )
1533 :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbIntNodes(0)
1535 _polygons.reserve(100); // to avoid reallocation;
1537 for ( int i = 0; i < 8; ++i )
1538 _nodeShift[i] = other._nodeShift[i];
1540 for ( int i = 0; i < 12; ++i )
1542 const _Link& srcLink = other._hexLinks[ i ];
1543 _Link& tgtLink = this->_hexLinks[ i ];
1544 tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1545 tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1546 tgtLink._intNodes.reserve( 10 ); // to avoid reallocation
1547 tgtLink._splits.reserve( 10 );
1550 for ( int i = 0; i < 6; ++i )
1552 const _Face& srcQuad = other._hexQuads[ i ];
1553 _Face& tgtQuad = this->_hexQuads[ i ];
1554 tgtQuad._links.resize(4);
1555 for ( int j = 0; j < 4; ++j )
1557 const _OrientedLink& srcLink = srcQuad._links[ j ];
1558 _OrientedLink& tgtLink = tgtQuad._links[ j ];
1559 tgtLink._reverse = srcLink._reverse;
1560 tgtLink._link = _hexLinks + ( srcLink._link - other._hexLinks );
1565 //================================================================================
1567 * \brief Initializes its data by given grid cell
1569 void Hexahedron::init( size_t i, size_t j, size_t k )
1571 _i = i; _j = j; _k = k;
1572 // set nodes of grid to nodes of the hexahedron and
1573 // count nodes at hexahedron corners located IN and ON geometry
1574 _nbCornerNodes = _nbBndNodes = 0;
1575 _origNodeInd = _grid->NodeIndex( i,j,k );
1576 for ( int iN = 0; iN < 8; ++iN )
1578 _hexNodes[iN]._node = _grid->_nodes [ _origNodeInd + _nodeShift[iN] ];
1579 _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
1580 _nbCornerNodes += bool( _hexNodes[iN]._node );
1581 _nbBndNodes += bool( _hexNodes[iN]._intPoint );
1584 _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1585 _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1586 _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1588 if ( _nbIntNodes + _edgeIntPnts.size() > 0 &&
1589 _nbIntNodes + _nbCornerNodes + _edgeIntPnts.size() > 3)
1592 // create sub-links (_splits) by splitting links with _intNodes
1593 for ( int iLink = 0; iLink < 12; ++iLink )
1595 _Link& link = _hexLinks[ iLink ];
1596 link._splits.clear();
1597 split._nodes[ 0 ] = link._nodes[0];
1598 bool isOut = ( ! link._nodes[0]->Node() ); // is1stNodeOut( iLink );
1599 bool checkTransition;
1600 for ( size_t i = 0; i < link._intNodes.size(); ++i )
1602 if ( link._intNodes[i].Node() ) // intersection non-coinsident with a grid node
1604 if ( split._nodes[ 0 ]->Node() && !isOut )
1606 split._nodes[ 1 ] = &link._intNodes[i];
1607 link._splits.push_back( split );
1609 split._nodes[ 0 ] = &link._intNodes[i];
1610 checkTransition = true;
1612 else // FACE intersection coinsident with a grid node
1614 checkTransition = ( link._nodes[0]->Node() );
1616 if ( checkTransition )
1618 switch ( link._intNodes[i].FaceIntPnt()->_transition ) {
1619 case Trans_OUT: isOut = true; break;
1620 case Trans_IN : isOut = false; break;
1621 default:; // isOut remains the same
1625 if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
1627 split._nodes[ 1 ] = link._nodes[1];
1628 link._splits.push_back( split );
1632 // Create _Node's at intersections with EDGEs.
1634 const double tol2 = _grid->_tol * _grid->_tol;
1635 int facets[3], nbFacets, subEntity;
1637 for ( size_t iP = 0; iP < _edgeIntPnts.size(); ++iP )
1639 nbFacets = getEntity( _edgeIntPnts[iP], facets, subEntity );
1640 _Node* equalNode = 0;
1641 switch( nbFacets ) {
1642 case 1: // in a _Face
1644 _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
1645 equalNode = FindEqualNode( quad._edgeNodes, _edgeIntPnts[ iP ], tol2 );
1647 equalNode->Add( _edgeIntPnts[ iP ] );
1650 quad._edgeNodes.push_back( _Node( 0, _edgeIntPnts[ iP ]));
1655 case 2: // on a _Link
1657 _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
1658 if ( link._splits.size() > 0 )
1660 equalNode = FindEqualNode( link._intNodes, _edgeIntPnts[ iP ], tol2 );
1662 equalNode->Add( _edgeIntPnts[ iP ] );
1666 for ( int iF = 0; iF < 2; ++iF )
1668 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1669 equalNode = FindEqualNode( quad._edgeNodes, _edgeIntPnts[ iP ], tol2 );
1671 equalNode->Add( _edgeIntPnts[ iP ] );
1674 quad._edgeNodes.push_back( _Node( 0, _edgeIntPnts[ iP ]));
1681 case 3: // at a corner
1683 _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
1684 if ( node.Node() > 0 )
1686 if ( node._intPoint )
1687 node._intPoint->Add( _edgeIntPnts[ iP ]->_faceIDs, _edgeIntPnts[ iP ]->_node );
1691 for ( int iF = 0; iF < 3; ++iF )
1693 _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1694 equalNode = FindEqualNode( quad._edgeNodes, _edgeIntPnts[ iP ], tol2 );
1696 equalNode->Add( _edgeIntPnts[ iP ] );
1699 quad._edgeNodes.push_back( _Node( 0, _edgeIntPnts[ iP ]));
1706 default: // inside a hex
1708 equalNode = FindEqualNode( _vertexNodes, _edgeIntPnts[ iP ], tol2 );
1710 equalNode->Add( _edgeIntPnts[ iP ] );
1713 _vertexNodes.push_back( _Node( 0, _edgeIntPnts[iP] ));
1717 } // switch( nbFacets )
1719 } // loop on _edgeIntPnts
1721 else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbIntNodes == 0
1724 // create sub-links (_splits) of whole links
1725 for ( int iLink = 0; iLink < 12; ++iLink )
1727 _Link& link = _hexLinks[ iLink ];
1728 link._splits.clear();
1729 if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
1731 split._nodes[ 0 ] = link._nodes[0];
1732 split._nodes[ 1 ] = link._nodes[1];
1733 link._splits.push_back( split );
1739 //================================================================================
1741 * \brief Initializes its data by given grid cell (countered from zero)
1743 void Hexahedron::init( size_t iCell )
1745 size_t iNbCell = _grid->_coords[0].size() - 1;
1746 size_t jNbCell = _grid->_coords[1].size() - 1;
1747 _i = iCell % iNbCell;
1748 _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1749 _k = iCell / iNbCell / jNbCell;
1753 //================================================================================
1755 * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1757 void Hexahedron::ComputeElements()
1761 if ( _nbCornerNodes + _nbIntNodes < 4 )
1764 if ( _nbBndNodes == _nbCornerNodes && _nbIntNodes == 0 && isInHole() )
1768 _polygons.reserve( 20 );
1770 // Create polygons from quadrangles
1771 // --------------------------------
1774 vector< _OrientedLink > splits;
1775 vector<_Node*> chainNodes;
1777 bool hasEdgeIntersections = !_edgeIntPnts.empty();
1779 for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1781 _Face& quad = _hexQuads[ iF ] ;
1783 _polygons.resize( _polygons.size() + 1 );
1784 _Face* polygon = &_polygons.back();
1785 polygon->_polyLinks.reserve( 20 );
1788 for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1789 for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
1790 splits.push_back( quad._links[ iE ].ResultLink( iS ));
1792 // add splits of links to a polygon and add _polyLinks to make
1793 // polygon's boundary closed
1795 int nbSplits = splits.size();
1796 if ( nbSplits < 2 && quad._edgeNodes.empty() )
1799 if ( nbSplits == 0 && !quad._edgeNodes.empty() )
1801 // make _vertexNodes from _edgeNodes of an empty quad
1802 const double tol2 = _grid->_tol * _grid->_tol;
1803 for ( size_t iP = 0; iP < quad._edgeNodes.size(); ++iP )
1806 FindEqualNode( _vertexNodes, quad._edgeNodes[ iP ].EdgeIntPnt(), tol2 );
1808 equalNode->Add( quad._edgeNodes[ iP ].EdgeIntPnt() );
1810 _vertexNodes.push_back( quad._edgeNodes[ iP ]);
1814 for ( size_t iP = 0; iP < quad._edgeNodes.size(); ++iP )
1815 quad._edgeNodes[ iP ]._isUsedInFace = false;
1817 int nbUsedEdgeNodes = 0;
1819 while ( nbSplits > 0 )
1822 while ( !splits[ iS ] )
1825 if ( !polygon->_links.empty() )
1827 _polygons.resize( _polygons.size() + 1 );
1828 polygon = &_polygons.back();
1829 polygon->_polyLinks.reserve( 20 );
1831 polygon->_links.push_back( splits[ iS ] );
1832 splits[ iS++ ]._link = 0;
1835 _Node* nFirst = polygon->_links.back().FirstNode();
1836 _Node *n1,*n2 = polygon->_links.back().LastNode();
1837 for ( ; nFirst != n2 && iS < splits.size(); ++iS )
1839 _OrientedLink& split = splits[ iS ];
1840 if ( !split ) continue;
1842 n1 = split.FirstNode();
1845 // try to connect to intersections with EDGEs
1846 if ( quad._edgeNodes.size() > nbUsedEdgeNodes &&
1847 findChain( n2, n1, quad, chainNodes ))
1849 for ( size_t i = 1; i < chainNodes.size(); ++i )
1851 polyLink._nodes[0] = chainNodes[i-1];
1852 polyLink._nodes[1] = chainNodes[i];
1853 polygon->_polyLinks.push_back( polyLink );
1854 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1855 nbUsedEdgeNodes += polyLink._nodes[1]->_isUsedInFace;
1857 if ( chainNodes.back() != n1 )
1859 n2 = chainNodes.back();
1864 // try to connect to a split ending on the same FACE
1867 _OrientedLink foundSplit;
1868 for ( int i = iS; i < splits.size() && !foundSplit; ++i )
1869 if (( foundSplit = splits[ i ]) &&
1870 ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
1872 polyLink._nodes[0] = n2;
1873 polyLink._nodes[1] = foundSplit.FirstNode();
1874 polygon->_polyLinks.push_back( polyLink );
1875 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1880 foundSplit._link = 0;
1884 n2 = foundSplit.FirstNode();
1889 if ( n2->IsLinked( nFirst->_intPoint ))
1891 polyLink._nodes[0] = n2;
1892 polyLink._nodes[1] = n1;
1893 polygon->_polyLinks.push_back( polyLink );
1894 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1898 polygon->_links.push_back( split );
1901 n2 = polygon->_links.back().LastNode();
1905 if ( nFirst != n2 ) // close a polygon
1907 if ( !findChain( n2, nFirst, quad, chainNodes ))
1909 if ( !closePolygon( polygon, chainNodes ))
1910 ;//chainNodes.push_back( nFirst );
1912 for ( size_t i = 1; i < chainNodes.size(); ++i )
1914 polyLink._nodes[0] = chainNodes[i-1];
1915 polyLink._nodes[1] = chainNodes[i];
1916 polygon->_polyLinks.push_back( polyLink );
1917 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1921 if ( polygon->_links.size() < 3 && nbSplits > 0 )
1923 polygon->_polyLinks.clear();
1924 polygon->_links.clear();
1926 } // while ( nbSplits > 0 )
1928 if ( polygon->_links.size() < 3 )
1929 _polygons.pop_back();
1931 } // loop on 6 sides of a hexahedron
1933 // Create polygons closing holes in a polyhedron
1934 // ----------------------------------------------
1936 // add polygons to their links
1937 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1939 _Face& polygon = _polygons[ iP ];
1940 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1942 polygon._links[ iL ]._link->_faces.reserve( 2 );
1943 polygon._links[ iL ]._link->_faces.push_back( &polygon );
1947 vector< _OrientedLink* > freeLinks;
1948 freeLinks.reserve(20);
1949 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1951 _Face& polygon = _polygons[ iP ];
1952 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1953 if ( polygon._links[ iL ]._link->_faces.size() < 2 )
1954 freeLinks.push_back( & polygon._links[ iL ]);
1956 int nbFreeLinks = freeLinks.size();
1957 if ( nbFreeLinks < 3 ) return;
1959 set<TGeomID> usedFaceIDs;
1961 // make closed chains of free links
1962 while ( nbFreeLinks > 0 )
1964 _polygons.resize( _polygons.size() + 1 );
1965 _Face& polygon = _polygons.back();
1966 polygon._polyLinks.reserve( 20 );
1967 polygon._links.reserve( 20 );
1969 _OrientedLink* curLink = 0;
1971 if ( !hasEdgeIntersections )
1973 // get a remaining link to start from
1974 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1975 if (( curLink = freeLinks[ iL ] ))
1976 freeLinks[ iL ] = 0;
1977 polygon._links.push_back( *curLink );
1981 // find all links connected to curLink
1982 curNode = curLink->FirstNode();
1984 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1985 if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
1987 curLink = freeLinks[ iL ];
1988 freeLinks[ iL ] = 0;
1990 polygon._links.push_back( *curLink );
1992 } while ( curLink );
1994 else // there are intersections with EDGEs
1997 // get a remaining link to start from, one lying on minimal
2000 vector< pair< TGeomID, int > > facesOfLink[3];
2001 pair< TGeomID, int > faceOfLink( -1, -1 );
2002 vector< TGeomID > faces;
2003 for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2004 if ( freeLinks[ iL ] )
2006 faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2007 if ( faces.size() == 1 )
2009 faceOfLink = make_pair( faces[0], iL );
2010 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2012 facesOfLink[0].push_back( faceOfLink );
2014 else if ( facesOfLink[0].empty() )
2016 faceOfLink = make_pair(( faces.empty() ? -1 : faces[0]), iL );
2017 facesOfLink[ 1 + faces.empty() ].push_back( faceOfLink );
2020 for ( int i = 0; faceOfLink.second < 0 && i < 3; ++i )
2021 if ( !facesOfLink[i].empty() )
2022 faceOfLink = facesOfLink[i][0];
2024 if ( faceOfLink.first < 0 ) // all faces used
2026 for ( size_t i = 0; i < facesOfLink[2].size() && faceOfLink.first < 1; ++i )
2028 curLink = freeLinks[ facesOfLink[2][i].second ];
2029 faceOfLink.first = curLink->FirstNode()->IsLinked( curLink->FirstNode()->_intPoint );
2031 usedFaceIDs.clear();
2033 curFace = faceOfLink.first;
2034 curLink = freeLinks[ faceOfLink.second ];
2035 freeLinks[ faceOfLink.second ] = 0;
2037 usedFaceIDs.insert( curFace );
2038 polygon._links.push_back( *curLink );
2041 // find all links bounding a FACE of curLink
2044 // go forward from curLink
2045 curNode = curLink->LastNode();
2047 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2048 if ( freeLinks[ iL ] &&
2049 freeLinks[ iL ]->FirstNode() == curNode &&
2050 freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2052 curLink = freeLinks[ iL ];
2053 freeLinks[ iL ] = 0;
2054 polygon._links.push_back( *curLink );
2057 } while ( curLink );
2059 std::reverse( polygon._links.begin(), polygon._links.end() );
2061 curLink = & polygon._links.back();
2064 // go backward from curLink
2065 curNode = curLink->FirstNode();
2067 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2068 if ( freeLinks[ iL ] &&
2069 freeLinks[ iL ]->LastNode() == curNode &&
2070 freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2072 curLink = freeLinks[ iL ];
2073 freeLinks[ iL ] = 0;
2074 polygon._links.push_back( *curLink );
2077 } while ( curLink );
2079 curNode = polygon._links.back().FirstNode();
2081 if ( polygon._links[0].LastNode() != curNode )
2083 if ( !_vertexNodes.empty() )
2085 // add links with _vertexNodes if not already used
2086 for ( size_t iN = 0; iN < _vertexNodes.size(); ++iN )
2087 if ( _vertexNodes[ iN ].IsOnFace( curFace ))
2089 bool used = ( curNode == &_vertexNodes[ iN ] );
2090 for ( size_t iL = 0; iL < polygon._links.size() && !used; ++iL )
2091 used = ( &_vertexNodes[ iN ] == polygon._links[ iL ].LastNode() );
2094 polyLink._nodes[0] = &_vertexNodes[ iN ];
2095 polyLink._nodes[1] = curNode;
2096 polygon._polyLinks.push_back( polyLink );
2097 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
2098 freeLinks.push_back( &polygon._links.back() );
2100 curNode = &_vertexNodes[ iN ];
2102 // TODO: to reorder _vertexNodes within polygon, if there are several ones
2105 polyLink._nodes[0] = polygon._links[0].LastNode();
2106 polyLink._nodes[1] = curNode;
2107 polygon._polyLinks.push_back( polyLink );
2108 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
2109 freeLinks.push_back( &polygon._links.back() );
2113 } // if there are intersections with EDGEs
2115 if ( polygon._links.size() < 2 ||
2116 polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
2117 return; // closed polygon not found -> invalid polyhedron
2119 if ( polygon._links.size() == 2 )
2121 if ( freeLinks.back() == &polygon._links.back() )
2123 freeLinks.back() = 0;
2126 vector< _Face*>& polygs1 = polygon._links.front()._link->_faces;
2127 vector< _Face*>& polygs2 = polygon._links.back()._link->_faces;
2128 _Face* polyg1 = ( polygs1.empty() ? 0 : polygs1[0] );
2129 _Face* polyg2 = ( polygs2.empty() ? 0 : polygs2[0] );
2130 if ( polyg1 ) polygs2.push_back( polyg1 );
2131 if ( polyg2 ) polygs1.push_back( polyg2 );
2132 _polygons.pop_back();
2136 // add polygon to its links
2137 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2139 polygon._links[ iL ]._link->_faces.reserve( 2 );
2140 polygon._links[ iL ]._link->_faces.push_back( &polygon );
2141 polygon._links[ iL ].Reverse();
2144 } // while ( nbFreeLinks > 0 )
2146 if ( ! checkPolyhedronSize() )
2151 // create a classic cell if possible
2152 const int nbNodes = _nbCornerNodes + _nbIntNodes;
2153 bool isClassicElem = false;
2154 if ( nbNodes == 8 && _polygons.size() == 6 ) isClassicElem = addHexa();
2155 else if ( nbNodes == 4 && _polygons.size() == 4 ) isClassicElem = addTetra();
2156 else if ( nbNodes == 6 && _polygons.size() == 5 ) isClassicElem = addPenta();
2157 else if ( nbNodes == 5 && _polygons.size() == 5 ) isClassicElem = addPyra ();
2158 if ( !isClassicElem )
2160 _volumeDefs._nodes.clear();
2161 _volumeDefs._quantities.clear();
2163 for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2165 const size_t nbLinks = _polygons[ iF ]._links.size();
2166 _volumeDefs._quantities.push_back( nbLinks );
2167 for ( size_t iL = 0; iL < nbLinks; ++iL )
2168 _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
2172 //================================================================================
2174 * \brief Create elements in the mesh
2176 int Hexahedron::MakeElements(SMESH_MesherHelper& helper,
2177 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2179 SMESHDS_Mesh* mesh = helper.GetMeshDS();
2181 size_t nbCells[3] = { _grid->_coords[0].size() - 1,
2182 _grid->_coords[1].size() - 1,
2183 _grid->_coords[2].size() - 1 };
2184 const size_t nbGridCells = nbCells[0] * nbCells[1] * nbCells[2];
2185 vector< Hexahedron* > intersectedHex( nbGridCells, 0 );
2188 // set intersection nodes from GridLine's to links of intersectedHex
2189 int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
2190 for ( int iDir = 0; iDir < 3; ++iDir )
2192 int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
2193 dInd[1][ iDirOther[iDir][0] ] = -1;
2194 dInd[2][ iDirOther[iDir][1] ] = -1;
2195 dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
2196 // loop on GridLine's parallel to iDir
2197 LineIndexer lineInd = _grid->GetLineIndexer( iDir );
2198 for ( ; lineInd.More(); ++lineInd )
2200 GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
2201 multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2202 for ( ; ip != line._intPoints.end(); ++ip )
2204 // if ( !ip->_node ) continue; // intersection at a grid node
2205 lineInd.SetIndexOnLine( ip->_indexOnLine );
2206 for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
2208 i = int(lineInd.I()) + dInd[iL][0];
2209 j = int(lineInd.J()) + dInd[iL][1];
2210 k = int(lineInd.K()) + dInd[iL][2];
2211 if ( i < 0 || i >= nbCells[0] ||
2212 j < 0 || j >= nbCells[1] ||
2213 k < 0 || k >= nbCells[2] ) continue;
2215 const size_t hexIndex = _grid->CellIndex( i,j,k );
2216 Hexahedron *& hex = intersectedHex[ hexIndex ];
2219 hex = new Hexahedron( *this );
2225 const int iLink = iL + iDir * 4;
2226 hex->_hexLinks[iLink]._intNodes.push_back( _Node( 0, &(*ip) ));
2227 hex->_nbIntNodes += bool( ip->_node );
2233 // implement geom edges into the mesh
2234 addEdges( helper, intersectedHex, edge2faceIDsMap );
2236 // add not split hexadrons to the mesh
2238 vector<int> intHexInd( nbIntHex );
2240 for ( size_t i = 0; i < intersectedHex.size(); ++i )
2242 Hexahedron * & hex = intersectedHex[ i ];
2245 intHexInd[ nbIntHex++ ] = i;
2246 if ( hex->_nbIntNodes > 0 || ! hex->_edgeIntPnts.empty())
2247 continue; // treat intersected hex later
2248 this->init( hex->_i, hex->_j, hex->_k );
2254 if (( _nbCornerNodes == 8 ) &&
2255 ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
2257 // order of _hexNodes is defined by enum SMESH_Block::TShapeID
2258 SMDS_MeshElement* el =
2259 mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
2260 _hexNodes[3].Node(), _hexNodes[1].Node(),
2261 _hexNodes[4].Node(), _hexNodes[6].Node(),
2262 _hexNodes[7].Node(), _hexNodes[5].Node() );
2263 mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
2268 intersectedHex[ i ] = 0;
2272 else if ( _nbCornerNodes > 3 && !hex )
2274 // all intersection of hex with geometry are at grid nodes
2275 hex = new Hexahedron( *this );
2280 intHexInd.push_back(0);
2281 intHexInd[ nbIntHex++ ] = i;
2285 // add elements resulted from hexadron intersection
2287 intHexInd.resize( nbIntHex );
2288 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, nbIntHex ),
2289 ParallelHexahedron( intersectedHex, intHexInd ),
2290 tbb::simple_partitioner()); // ComputeElements() is called here
2291 for ( size_t i = 0; i < intHexInd.size(); ++i )
2292 if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
2293 nbAdded += hex->addElements( helper );
2295 for ( size_t i = 0; i < intHexInd.size(); ++i )
2296 if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
2298 hex->ComputeElements();
2299 nbAdded += hex->addElements( helper );
2303 for ( size_t i = 0; i < intersectedHex.size(); ++i )
2304 if ( intersectedHex[ i ] )
2305 delete intersectedHex[ i ];
2310 //================================================================================
2312 * \brief Implements geom edges into the mesh
2314 void Hexahedron::addEdges(SMESH_MesherHelper& helper,
2315 vector< Hexahedron* >& hexes,
2316 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2318 if ( edge2faceIDsMap.empty() ) return;
2320 // Prepare planes for intersecting with EDGEs
2323 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
2325 GridPlanes& planes = pln[ iDirZ ];
2326 int iDirX = ( iDirZ + 1 ) % 3;
2327 int iDirY = ( iDirZ + 2 ) % 3;
2328 // planes._uNorm = ( _grid->_axes[ iDirY ] ^ _grid->_axes[ iDirZ ] ).Normalized();
2329 // planes._vNorm = ( _grid->_axes[ iDirZ ] ^ _grid->_axes[ iDirX ] ).Normalized();
2330 planes._zNorm = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
2331 planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
2332 planes._zProjs [0] = 0;
2333 const double zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2334 const vector< double > & u = _grid->_coords[ iDirZ ];
2335 for ( int i = 1; i < planes._zProjs.size(); ++i )
2337 planes._zProjs [i] = zFactor * ( u[i] - u[0] );
2341 const double deflection = _grid->_minCellSize / 20.;
2342 const double tol = _grid->_tol;
2343 E_IntersectPoint ip;
2345 // Intersect EDGEs with the planes
2346 map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
2347 for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
2349 const TGeomID edgeID = e2fIt->first;
2350 const TopoDS_Edge & E = TopoDS::Edge( _grid->_shapes( edgeID ));
2351 BRepAdaptor_Curve curve( E );
2352 TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
2353 TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
2355 ip._faceIDs = e2fIt->second;
2356 ip._shapeID = edgeID;
2358 // discretize the EGDE
2359 GCPnts_UniformDeflection discret( curve, deflection, true );
2360 if ( !discret.IsDone() || discret.NbPoints() < 2 )
2363 // perform intersection
2364 for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
2366 GridPlanes& planes = pln[ iDirZ ];
2367 int iDirX = ( iDirZ + 1 ) % 3;
2368 int iDirY = ( iDirZ + 2 ) % 3;
2369 double xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
2370 double yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
2371 double zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
2372 //double zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2373 int dIJK[3], d000[3] = { 0,0,0 };
2374 double o[3] = { _grid->_coords[0][0],
2375 _grid->_coords[1][0],
2376 _grid->_coords[2][0] };
2378 // locate the 1st point of a segment within the grid
2379 gp_XYZ p1 = discret.Value( 1 ).XYZ();
2380 double u1 = discret.Parameter( 1 );
2381 double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
2383 _grid->ComputeUVW( p1, ip._uvw );
2384 int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
2385 int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
2386 int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
2387 locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
2388 locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
2389 locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
2391 int ijk[3]; // grid index where a segment intersect a plane
2396 // add the 1st vertex point to a hexahedron
2400 _grid->_edgeIntP.push_back( ip );
2401 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2402 _grid->_edgeIntP.pop_back();
2404 for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
2406 // locate the 2nd point of a segment within the grid
2407 gp_XYZ p2 = discret.Value( iP ).XYZ();
2408 double u2 = discret.Parameter( iP );
2409 double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
2411 locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
2413 // treat intersections with planes between 2 end points of a segment
2414 int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
2415 int iZ = iZ1 + ( iZ1 < iZ2 );
2416 for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
2418 ip._point = findIntPoint( u1, zProj1, u2, zProj2,
2419 planes._zProjs[ iZ ],
2420 curve, planes._zNorm, _grid->_origin );
2421 _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
2422 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2423 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2426 // add ip to hex "above" the plane
2427 _grid->_edgeIntP.push_back( ip );
2429 bool added = addIntersection(_grid->_edgeIntP.back(), hexes, ijk, dIJK);
2431 // add ip to hex "below" the plane
2432 ijk[ iDirZ ] = iZ-1;
2433 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, dIJK ) &&
2435 _grid->_edgeIntP.pop_back();
2442 // add the 2nd vertex point to a hexahedron
2446 _grid->ComputeUVW( p1, ip._uvw );
2447 locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2448 locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2450 _grid->_edgeIntP.push_back( ip );
2451 if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2452 _grid->_edgeIntP.pop_back();
2454 } // loop on 3 grid directions
2457 // Create nodes at found intersections
2458 // const E_IntersectPoint* eip;
2459 // for ( size_t i = 0; i < hexes.size(); ++i )
2461 // Hexahedron* h = hexes[i];
2462 // if ( !h ) continue;
2463 // for ( int iF = 0; iF < 6; ++iF )
2465 // _Face& quad = h->_hexQuads[ iF ];
2466 // for ( size_t iP = 0; iP < quad._edgeNodes.size(); ++iP )
2467 // if ( !quad._edgeNodes[ iP ]._node )
2468 // if (( eip = quad._edgeNodes[ iP ].EdgeIntPnt() ))
2469 // quad._edgeNodes[ iP ]._intPoint->_node = helper.AddNode( eip->_point.X(),
2471 // eip->_point.Z() );
2473 // for ( size_t iP = 0; iP < hexes[i]->_vertexNodes.size(); ++iP )
2474 // if (( eip = h->_vertexNodes[ iP ].EdgeIntPnt() ))
2475 // h->_vertexNodes[ iP ]._intPoint->_node = helper.AddNode( eip->_point.X(),
2477 // eip->_point.Z() );
2481 //================================================================================
2483 * \brief Finds intersection of a curve with a plane
2484 * \param [in] u1 - parameter of one curve point
2485 * \param [in] proj1 - projection of the curve point to the plane normal
2486 * \param [in] u2 - parameter of another curve point
2487 * \param [in] proj2 - projection of the other curve point to the plane normal
2488 * \param [in] proj - projection of a point where the curve intersects the plane
2489 * \param [in] curve - the curve
2490 * \param [in] axis - the plane normal
2491 * \param [in] origin - the plane origin
2492 * \return gp_Pnt - the found intersection point
2494 gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
2495 double u2, double proj2,
2497 BRepAdaptor_Curve& curve,
2499 const gp_XYZ& origin)
2501 double r = (( proj - proj1 ) / ( proj2 - proj1 ));
2502 double u = u1 * ( 1 - r ) + u2 * r;
2503 gp_Pnt p = curve.Value( u );
2504 double newProj = axis * ( p.XYZ() - origin );
2505 if ( Abs( proj - newProj ) > _grid->_tol / 10. )
2508 return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
2510 return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
2515 //================================================================================
2517 * \brief Returns indices of a hexahedron sub-entities holding a point
2518 * \param [in] ip - intersection point
2519 * \param [out] facets - 0-3 facets holding a point
2520 * \param [out] sub - index of a vertex or an edge holding a point
2521 * \return int - number of facets holding a point
2523 int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
2525 enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
2527 int vertex = 0, egdeMask = 0;
2529 if ( Abs( _grid->_coords[0][ _i ] - ip->_uvw[0] ) < _grid->_tol ) {
2530 facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
2533 else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
2534 facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
2538 if ( Abs( _grid->_coords[1][ _j ] - ip->_uvw[1] ) < _grid->_tol ) {
2539 facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
2542 else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
2543 facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
2547 if ( Abs( _grid->_coords[2][ _k ] - ip->_uvw[2] ) < _grid->_tol ) {
2548 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
2551 else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
2552 facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
2559 case 0: sub = 0; break;
2560 case 1: sub = facets[0]; break;
2562 const int edge [3][8] = {
2563 { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
2564 SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
2565 { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
2566 SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
2567 { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
2568 SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
2570 switch ( egdeMask ) {
2571 case X | Y: sub = edge[ 0 ][ vertex ]; break;
2572 case X | Z: sub = edge[ 1 ][ vertex ]; break;
2573 default: sub = edge[ 2 ][ vertex ];
2579 sub = vertex + SMESH_Block::ID_FirstV;
2584 //================================================================================
2586 * \brief Adds intersection with an EDGE
2588 bool Hexahedron::addIntersection( const E_IntersectPoint& ip,
2589 vector< Hexahedron* >& hexes,
2590 int ijk[], int dIJK[] )
2594 size_t hexIndex[4] = {
2595 _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
2596 dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
2597 dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
2598 dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
2600 for ( int i = 0; i < 4; ++i )
2602 if ( /*0 <= hexIndex[i] &&*/ hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
2604 Hexahedron* h = hexes[ hexIndex[i] ];
2605 // check if ip is really inside the hex
2607 if (( _grid->_coords[0][ h->_i ] - _grid->_tol > ip._uvw[0] ) ||
2608 ( _grid->_coords[0][ h->_i+1 ] + _grid->_tol < ip._uvw[0] ) ||
2609 ( _grid->_coords[1][ h->_j ] - _grid->_tol > ip._uvw[1] ) ||
2610 ( _grid->_coords[1][ h->_j+1 ] + _grid->_tol < ip._uvw[1] ) ||
2611 ( _grid->_coords[2][ h->_k ] - _grid->_tol > ip._uvw[2] ) ||
2612 ( _grid->_coords[2][ h->_k+1 ] + _grid->_tol < ip._uvw[2] ))
2613 throw SALOME_Exception("ip outside a hex");
2615 h->_edgeIntPnts.push_back( & ip );
2621 //================================================================================
2623 * \brief Finds nodes at a path from one node to another via intersections with EDGEs
2625 bool Hexahedron::findChain( _Node* n1,
2628 vector<_Node*>& chn )
2631 chn.push_back( n1 );
2632 for ( size_t iP = 0; iP < quad._edgeNodes.size(); ++iP )
2633 if ( !quad._edgeNodes[ iP ]._isUsedInFace &&
2634 n1->IsLinked( quad._edgeNodes[ iP ]._intPoint ) &&
2635 n2->IsLinked( quad._edgeNodes[ iP ]._intPoint ))
2637 chn.push_back( & quad._edgeNodes[ iP ]);
2638 chn.push_back( n2 );
2639 quad._edgeNodes[ iP ]._isUsedInFace = true;
2646 for ( size_t iP = 0; iP < quad._edgeNodes.size(); ++iP )
2647 if ( !quad._edgeNodes[ iP ]._isUsedInFace &&
2648 chn.back()->IsLinked( quad._edgeNodes[ iP ]._intPoint ))
2650 chn.push_back( & quad._edgeNodes[ iP ]);
2651 found = quad._edgeNodes[ iP ]._isUsedInFace = true;
2654 } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
2656 if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
2657 chn.push_back( n2 );
2659 return chn.size() > 1;
2661 //================================================================================
2663 * \brief Try to heal a polygon whose ends are not connected
2665 bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
2667 int i = -1, nbLinks = polygon->_links.size();
2670 vector< _OrientedLink > newLinks;
2671 // find a node lying on the same FACE as the last one
2672 _Node* node = polygon->_links.back().LastNode();
2673 int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
2674 for ( i = nbLinks - 2; i >= 0; --i )
2675 if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
2679 for ( ; i < nbLinks; ++i )
2680 newLinks.push_back( polygon->_links[i] );
2684 // find a node lying on the same FACE as the first one
2685 node = polygon->_links[0].FirstNode();
2686 avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
2687 for ( i = 1; i < nbLinks; ++i )
2688 if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
2691 for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
2692 newLinks.push_back( polygon->_links[i] );
2694 if ( newLinks.size() > 1 )
2696 polygon->_links.swap( newLinks );
2698 chainNodes.push_back( polygon->_links.back().LastNode() );
2699 chainNodes.push_back( polygon->_links[0].FirstNode() );
2704 //================================================================================
2706 * \brief Checks transition at the 1st node of a link
2708 bool Hexahedron::is1stNodeOut( int iLink ) const
2710 if ( !_hexLinks[ iLink ]._nodes[0]->Node() ) // no node
2712 if ( !_hexLinks[ iLink ]._nodes[0]->_intPoint ) // no intersection with geometry
2714 switch ( _hexLinks[ iLink ]._nodes[0]->FaceIntPnt()->_transition ) {
2715 case Trans_OUT: return true;
2716 case Trans_IN : return false;
2717 default: ; // tangent transition
2720 // ijk of a GridLine corresponding to the link
2721 int iDir = iLink / 4;
2722 int indSub = iLink % 4;
2723 LineIndexer li = _grid->GetLineIndexer( iDir );
2724 li.SetIJK( _i,_j,_k );
2725 size_t lineIndex[4] = { li.LineIndex (),
2729 GridLine& line = _grid->_lines[ iDir ][ lineIndex[ indSub ]];
2731 // analyze transition of previous ip
2733 multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2734 for ( ; ip != line._intPoints.end(); ++ip )
2736 if ( &(*ip) == _hexLinks[ iLink ]._nodes[0]->_intPoint )
2738 switch ( ip->_transition ) {
2739 case Trans_OUT: isOut = true;
2740 case Trans_IN : isOut = false;
2745 if ( ip == line._intPoints.end() )
2746 cout << "BUG: Wrong GridLine. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl;
2750 //================================================================================
2752 * \brief Adds computed elements to the mesh
2754 int Hexahedron::addElements(SMESH_MesherHelper& helper)
2757 // add elements resulted from hexahedron intersection
2758 //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
2760 vector< const SMDS_MeshNode* > nodes( _volumeDefs._nodes.size() );
2761 for ( size_t iN = 0; iN < nodes.size(); ++iN )
2762 if ( !( nodes[iN] = _volumeDefs._nodes[iN]->Node() ))
2764 if ( const E_IntersectPoint* eip = _volumeDefs._nodes[iN]->EdgeIntPnt() )
2765 nodes[iN] = _volumeDefs._nodes[iN]->_intPoint->_node =
2766 helper.AddNode( eip->_point.X(),
2770 throw SALOME_Exception("Bug: no node at intersection point");
2773 if ( !_volumeDefs._quantities.empty() )
2775 helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
2779 switch ( nodes.size() )
2781 case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
2782 nodes[4],nodes[5],nodes[6],nodes[7] );
2784 case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
2786 case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
2789 helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
2793 nbAdded += int ( _volumeDefs._nodes.size() > 0 );
2798 //================================================================================
2800 * \brief Return true if the element is in a hole
2802 bool Hexahedron::isInHole() const
2804 if ( !_vertexNodes.empty() )
2807 const int ijk[3] = { _i, _j, _k };
2808 F_IntersectPoint curIntPnt;
2810 // consider a cell to be in a hole if all links in any direction
2811 // comes OUT of geometry
2812 for ( int iDir = 0; iDir < 3; ++iDir )
2814 const vector<double>& coords = _grid->_coords[ iDir ];
2815 LineIndexer li = _grid->GetLineIndexer( iDir );
2816 li.SetIJK( _i,_j,_k );
2817 size_t lineIndex[4] = { li.LineIndex (),
2821 bool allLinksOut = true, hasLinks = false;
2822 for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
2824 const _Link& link = _hexLinks[ iL + 4*iDir ];
2825 // check transition of the first node of a link
2826 const F_IntersectPoint* firstIntPnt = 0;
2827 if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
2829 curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0];
2830 const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
2831 multiset< F_IntersectPoint >::const_iterator ip =
2832 line._intPoints.upper_bound( curIntPnt );
2834 firstIntPnt = &(*ip);
2836 else if ( !link._intNodes.empty() )
2838 firstIntPnt = link._intNodes[0].FaceIntPnt();
2844 allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
2847 if ( hasLinks && allLinksOut )
2853 //================================================================================
2855 * \brief Return true if a polyhedron passes _sizeThreshold criterion
2857 bool Hexahedron::checkPolyhedronSize() const
2860 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2862 const _Face& polygon = _polygons[iP];
2863 gp_XYZ area (0,0,0);
2864 gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
2865 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2867 gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
2871 volume += p1 * area;
2875 double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
2877 return volume > initVolume / _sizeThreshold;
2879 //================================================================================
2881 * \brief Tries to create a hexahedron
2883 bool Hexahedron::addHexa()
2885 if ( _polygons[0]._links.size() != 4 ||
2886 _polygons[1]._links.size() != 4 ||
2887 _polygons[2]._links.size() != 4 ||
2888 _polygons[3]._links.size() != 4 ||
2889 _polygons[4]._links.size() != 4 ||
2890 _polygons[5]._links.size() != 4 )
2894 for ( int iL = 0; iL < 4; ++iL )
2897 nodes[iL] = _polygons[0]._links[iL].FirstNode();
2900 // find a top node above the base node
2901 _Link* link = _polygons[0]._links[iL]._link;
2902 //ASSERT( link->_faces.size() > 1 );
2903 if ( link->_faces.size() < 2 )
2904 return debugDumpLink( link );
2905 // a quadrangle sharing <link> with _polygons[0]
2906 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
2907 for ( int i = 0; i < 4; ++i )
2908 if ( quad->_links[i]._link == link )
2910 // 1st node of a link opposite to <link> in <quad>
2911 nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
2917 _volumeDefs.set( vector< _Node* >( nodes, nodes+8 ));
2921 //================================================================================
2923 * \brief Tries to create a tetrahedron
2925 bool Hexahedron::addTetra()
2928 nodes[0] = _polygons[0]._links[0].FirstNode();
2929 nodes[1] = _polygons[0]._links[1].FirstNode();
2930 nodes[2] = _polygons[0]._links[2].FirstNode();
2932 _Link* link = _polygons[0]._links[0]._link;
2933 //ASSERT( link->_faces.size() > 1 );
2934 if ( link->_faces.size() < 2 )
2935 return debugDumpLink( link );
2937 // a triangle sharing <link> with _polygons[0]
2938 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
2939 for ( int i = 0; i < 3; ++i )
2940 if ( tria->_links[i]._link == link )
2942 nodes[3] = tria->_links[(i+1)%3].LastNode();
2943 _volumeDefs.set( vector< _Node* >( nodes, nodes+4 ));
2949 //================================================================================
2951 * \brief Tries to create a pentahedron
2953 bool Hexahedron::addPenta()
2955 // find a base triangular face
2957 for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
2958 if ( _polygons[ iF ]._links.size() == 3 )
2960 if ( iTri < 0 ) return false;
2965 for ( int iL = 0; iL < 3; ++iL )
2968 nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
2971 // find a top node above the base node
2972 _Link* link = _polygons[ iTri ]._links[iL]._link;
2973 //ASSERT( link->_faces.size() > 1 );
2974 if ( link->_faces.size() < 2 )
2975 return debugDumpLink( link );
2976 // a quadrangle sharing <link> with a base triangle
2977 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
2978 if ( quad->_links.size() != 4 ) return false;
2979 for ( int i = 0; i < 4; ++i )
2980 if ( quad->_links[i]._link == link )
2982 // 1st node of a link opposite to <link> in <quad>
2983 nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
2989 _volumeDefs.set( vector< _Node* >( nodes, nodes+6 ));
2991 return ( nbN == 6 );
2993 //================================================================================
2995 * \brief Tries to create a pyramid
2997 bool Hexahedron::addPyra()
2999 // find a base quadrangle
3001 for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
3002 if ( _polygons[ iF ]._links.size() == 4 )
3004 if ( iQuad < 0 ) return false;
3008 nodes[0] = _polygons[iQuad]._links[0].FirstNode();
3009 nodes[1] = _polygons[iQuad]._links[1].FirstNode();
3010 nodes[2] = _polygons[iQuad]._links[2].FirstNode();
3011 nodes[3] = _polygons[iQuad]._links[3].FirstNode();
3013 _Link* link = _polygons[iQuad]._links[0]._link;
3014 ASSERT( link->_faces.size() > 1 );
3015 if ( link->_faces.size() < 2 )
3016 return debugDumpLink( link );
3018 // a triangle sharing <link> with a base quadrangle
3019 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
3020 if ( tria->_links.size() != 3 ) return false;
3021 for ( int i = 0; i < 3; ++i )
3022 if ( tria->_links[i]._link == link )
3024 nodes[4] = tria->_links[(i+1)%3].LastNode();
3025 _volumeDefs.set( vector< _Node* >( nodes, nodes+5 ));
3031 //================================================================================
3033 * \brief Dump a link and return \c false
3035 bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
3038 gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
3039 cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
3040 << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
3041 << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
3046 //================================================================================
3048 * \brief computes exact bounding box with axes parallel to given ones
3050 //================================================================================
3052 void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
3053 const double* axesDirs,
3057 TopoDS_Compound allFacesComp;
3058 b.MakeCompound( allFacesComp );
3059 for ( size_t iF = 0; iF < faceVec.size(); ++iF )
3060 b.Add( allFacesComp, faceVec[ iF ] );
3062 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
3063 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
3065 for ( int i = 0; i < 6; ++i )
3066 farDist = Max( farDist, 10 * sP[i] );
3068 gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
3069 gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
3070 gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
3071 axis[0].Normalize();
3072 axis[1].Normalize();
3073 axis[2].Normalize();
3075 gp_Mat basis( axis[0], axis[1], axis[2] );
3076 gp_Mat bi = basis.Inverted();
3079 for ( int iDir = 0; iDir < 3; ++iDir )
3081 gp_XYZ axis0 = axis[ iDir ];
3082 gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
3083 gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
3084 for ( int isMax = 0; isMax < 2; ++isMax )
3086 double shift = isMax ? farDist : -farDist;
3087 gp_XYZ orig = shift * axis0;
3088 gp_XYZ norm = axis1 ^ axis2;
3089 gp_Pln pln( orig, norm );
3090 norm = pln.Axis().Direction().XYZ();
3091 BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
3093 gp_Pnt& pAxis = isMax ? pMax : pMin;
3094 gp_Pnt pPlane, pFaces;
3095 double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
3100 for ( int i = 0; i < 2; ++i ) {
3101 corner.SetCoord( 1, sP[ i*3 ]);
3102 for ( int j = 0; j < 2; ++j ) {
3103 corner.SetCoord( 2, sP[ i*3 + 1 ]);
3104 for ( int k = 0; k < 2; ++k )
3106 corner.SetCoord( 3, sP[ i*3 + 2 ]);
3112 corner = isMax ? bb.CornerMax() : bb.CornerMin();
3113 pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
3117 gp_XYZ pf = pFaces.XYZ() * bi;
3118 pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
3124 shapeBox.Add( pMin );
3125 shapeBox.Add( pMax );
3132 //=============================================================================
3134 * \brief Generates 3D structured Cartesian mesh in the internal part of
3135 * solid shapes and polyhedral volumes near the shape boundary.
3136 * \param theMesh - mesh to fill in
3137 * \param theShape - a compound of all SOLIDs to mesh
3138 * \retval bool - true in case of success
3140 //=============================================================================
3142 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh & theMesh,
3143 const TopoDS_Shape & theShape)
3145 // The algorithm generates the mesh in following steps:
3147 // 1) Intersection of grid lines with the geometry boundary.
3148 // This step allows to find out if a given node of the initial grid is
3149 // inside or outside the geometry.
3151 // 2) For each cell of the grid, check how many of it's nodes are outside
3152 // of the geometry boundary. Depending on a result of this check
3153 // - skip a cell, if all it's nodes are outside
3154 // - skip a cell, if it is too small according to the size threshold
3155 // - add a hexahedron in the mesh, if all nodes are inside
3156 // - add a polyhedron in the mesh, if some nodes are inside and some outside
3158 _computeCanceled = false;
3160 SMESH_MesherHelper helper( theMesh );
3166 vector< TopoDS_Shape > faceVec;
3168 TopTools_MapOfShape faceMap;
3169 TopExp_Explorer fExp;
3170 for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
3171 if ( !faceMap.Add( fExp.Current() ))
3172 faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
3174 for ( fExp.ReInit(); fExp.More(); fExp.Next() )
3175 if ( faceMap.Contains( fExp.Current() ))
3176 faceVec.push_back( fExp.Current() );
3178 vector<FaceGridIntersector> facesItersectors( faceVec.size() );
3179 map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
3180 TopExp_Explorer eExp;
3182 for ( int i = 0; i < faceVec.size(); ++i )
3184 facesItersectors[i]._face = TopoDS::Face ( faceVec[i] );
3185 facesItersectors[i]._faceID = grid._shapes.Add( faceVec[i] );
3186 facesItersectors[i]._grid = &grid;
3187 shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
3189 if ( _hyp->GetToAddEdges() )
3191 helper.SetSubShape( faceVec[i] );
3192 for ( eExp.Init( faceVec[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
3194 const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
3195 if ( !SMESH_Algo::isDegenerated( edge ) &&
3196 !helper.IsRealSeam( edge ))
3197 edge2faceIDsMap[ grid._shapes.Add( edge )].push_back( facesItersectors[i]._faceID );
3202 getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
3204 vector<double> xCoords, yCoords, zCoords;
3205 _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
3207 grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
3209 if ( _computeCanceled ) return false;
3212 { // copy partner faces and curves of not thread-safe types
3213 set< const Standard_Transient* > tshapes;
3214 BRepBuilderAPI_Copy copier;
3215 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3217 if ( !facesItersectors[i].IsThreadSafe(tshapes) )
3219 copier.Perform( facesItersectors[i]._face );
3220 facesItersectors[i]._face = TopoDS::Face( copier );
3224 // Intersection of grid lines with the geometry boundary.
3225 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
3226 ParallelIntersector( facesItersectors ),
3227 tbb::simple_partitioner());
3229 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3230 facesItersectors[i].Intersect();
3233 // put interesection points onto the GridLine's; this is done after intersection
3234 // to avoid contention of facesItersectors for writing into the same GridLine
3235 // in case of parallel work of facesItersectors
3236 for ( size_t i = 0; i < facesItersectors.size(); ++i )
3237 facesItersectors[i].StoreIntersections();
3239 TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
3240 helper.SetSubShape( solidExp.Current() );
3241 helper.SetElementsOnShape( true );
3243 if ( _computeCanceled ) return false;
3245 // create nodes on the geometry
3246 grid.ComputeNodes(helper);
3248 if ( _computeCanceled ) return false;
3250 // create volume elements
3251 Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
3252 int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
3254 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
3257 // make all SOLIDs computed
3258 if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
3260 SMDS_ElemIteratorPtr volIt = sm1->GetElements();
3261 for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
3263 const SMDS_MeshElement* vol = volIt->next();
3264 sm1->RemoveElement( vol, /*isElemDeleted=*/false );
3265 meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
3268 // make other sub-shapes computed
3269 setSubmeshesComputed( theMesh, theShape );
3272 // remove free nodes
3273 if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
3275 TIDSortedNodeSet nodesToRemove;
3276 // get intersection nodes
3277 for ( int iDir = 0; iDir < 3; ++iDir )
3279 vector< GridLine >& lines = grid._lines[ iDir ];
3280 for ( size_t i = 0; i < lines.size(); ++i )
3282 multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
3283 for ( ; ip != lines[i]._intPoints.end(); ++ip )
3284 if ( ip->_node && ip->_node->NbInverseElements() == 0 )
3285 nodesToRemove.insert( nodesToRemove.end(), ip->_node );
3289 for ( size_t i = 0; i < grid._nodes.size(); ++i )
3290 if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
3291 nodesToRemove.insert( nodesToRemove.end(), grid._nodes[i] );
3294 TIDSortedNodeSet::iterator n = nodesToRemove.begin();
3295 for ( ; n != nodesToRemove.end(); ++n )
3296 meshDS->RemoveFreeNode( *n, smDS, /*fromGroups=*/false );
3302 // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
3303 catch ( SMESH_ComputeError& e)
3305 return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
3310 //=============================================================================
3314 //=============================================================================
3316 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh & theMesh,
3317 const TopoDS_Shape & theShape,
3318 MapShapeNbElems& theResMap)
3321 // std::vector<int> aResVec(SMDSEntity_Last);
3322 // for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
3323 // if(IsQuadratic) {
3324 // aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3325 // int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
3326 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
3329 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
3330 // aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3332 // SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
3333 // aResMap.insert(std::make_pair(sm,aResVec));
3338 //=============================================================================
3342 * \brief Event listener setting/unsetting _alwaysComputed flag to
3343 * submeshes of inferior levels to prevent their computing
3345 struct _EventListener : public SMESH_subMeshEventListener
3349 _EventListener(const string& algoName):
3350 SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
3353 // --------------------------------------------------------------------------------
3354 // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
3356 static void setAlwaysComputed( const bool isComputed,
3357 SMESH_subMesh* subMeshOfSolid)
3359 SMESH_subMeshIteratorPtr smIt =
3360 subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3361 while ( smIt->more() )
3363 SMESH_subMesh* sm = smIt->next();
3364 sm->SetIsAlwaysComputed( isComputed );
3366 subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
3369 // --------------------------------------------------------------------------------
3370 // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
3372 virtual void ProcessEvent(const int event,
3373 const int eventType,
3374 SMESH_subMesh* subMeshOfSolid,
3375 SMESH_subMeshEventListenerData* data,
3376 const SMESH_Hypothesis* hyp = 0)
3378 if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
3380 setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
3385 SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
3386 if ( !algo3D || _algoName != algo3D->GetName() )
3387 setAlwaysComputed( false, subMeshOfSolid );
3391 // --------------------------------------------------------------------------------
3392 // set the event listener
3394 static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
3396 subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
3401 }; // struct _EventListener
3405 //================================================================================
3407 * \brief Sets event listener to submeshes if necessary
3408 * \param subMesh - submesh where algo is set
3409 * This method is called when a submesh gets HYP_OK algo_state.
3410 * After being set, event listener is notified on each event of a submesh.
3412 //================================================================================
3414 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
3416 _EventListener::SetOn( subMesh, GetName() );
3419 //================================================================================
3421 * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
3423 //================================================================================
3425 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh& theMesh,
3426 const TopoDS_Shape& theShape)
3428 for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
3429 _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));