1 // Copyright (C) 2007-2012 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"
39 #include <BRepAdaptor_Surface.hxx>
40 #include <BRepBndLib.hxx>
41 #include <BRepBuilderAPI_Copy.hxx>
42 #include <BRepTools.hxx>
43 #include <BRep_Tool.hxx>
44 #include <Bnd_Box.hxx>
46 #include <Geom2d_BSplineCurve.hxx>
47 #include <Geom2d_BezierCurve.hxx>
48 #include <Geom2d_TrimmedCurve.hxx>
49 #include <Geom_BSplineCurve.hxx>
50 #include <Geom_BSplineSurface.hxx>
51 #include <Geom_BezierCurve.hxx>
52 #include <Geom_BezierSurface.hxx>
53 #include <Geom_RectangularTrimmedSurface.hxx>
54 #include <Geom_TrimmedCurve.hxx>
55 #include <IntAna_IntConicQuad.hxx>
56 #include <IntAna_IntLinTorus.hxx>
57 #include <IntAna_Quadric.hxx>
58 #include <IntCurveSurface_TransitionOnCurve.hxx>
59 #include <IntCurvesFace_Intersector.hxx>
60 #include <Poly_Triangulation.hxx>
61 #include <Precision.hxx>
63 #include <TopExp_Explorer.hxx>
64 #include <TopLoc_Location.hxx>
65 #include <TopTools_MapIteratorOfMapOfShape.hxx>
66 #include <TopTools_MapOfShape.hxx>
68 #include <TopoDS_Face.hxx>
69 #include <TopoDS_TShape.hxx>
70 #include <gp_Cone.hxx>
71 #include <gp_Cylinder.hxx>
74 #include <gp_Pnt2d.hxx>
75 #include <gp_Sphere.hxx>
76 #include <gp_Torus.hxx>
80 #include <tbb/parallel_for.h>
81 //#include <tbb/enumerable_thread_specific.h>
88 #define ELLIPSOLID_WORKAROUND // remove it as soon as http://tracker.dev.opencascade.org/view.php?id=22809 is solved
90 #ifdef ELLIPSOLID_WORKAROUND
91 #include <BRepIntCurveSurface_Inter.hxx>
92 #include <BRepTopAdaptor_TopolTool.hxx>
93 #include <BRepAdaptor_HSurface.hxx>
96 //=============================================================================
100 //=============================================================================
102 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
103 :SMESH_3D_Algo(hypId, studyId, gen)
105 _name = "Cartesian_3D";
106 _shapeType = (1 << TopAbs_SOLID); // 1 bit /shape type
107 _compatibleHypothesis.push_back("CartesianParameters3D");
109 _onlyUnaryInput = false; // to mesh all SOLIDs at once
110 _requireDiscreteBoundary = false; // 2D mesh not needed
111 _supportSubmeshes = false; // do not use any existing mesh
114 //=============================================================================
116 * Check presence of a hypothesis
118 //=============================================================================
120 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh& aMesh,
121 const TopoDS_Shape& aShape,
122 Hypothesis_Status& aStatus)
124 aStatus = SMESH_Hypothesis::HYP_MISSING;
126 const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
127 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
128 if ( h == hyps.end())
133 for ( ; h != hyps.end(); ++h )
135 if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
137 aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
142 return aStatus == HYP_OK;
147 //=============================================================================
148 // Definitions of internal utils
149 // --------------------------------------------------------------------------
151 Trans_TANGENT = IntCurveSurface_Tangent,
152 Trans_IN = IntCurveSurface_In,
153 Trans_OUT = IntCurveSurface_Out,
156 // --------------------------------------------------------------------------
158 * \brief Data of intersection between a GridLine and a TopoDS_Face
160 struct IntersectionPoint
163 mutable Transition _transition;
164 mutable const SMDS_MeshNode* _node;
165 mutable size_t _indexOnLine;
167 IntersectionPoint(): _node(0) {}
168 bool operator< ( const IntersectionPoint& o ) const { return _paramOnLine < o._paramOnLine; }
170 // --------------------------------------------------------------------------
172 * \brief A line of the grid and its intersections with 2D geometry
177 double _length; // line length
178 multiset< IntersectionPoint > _intPoints;
180 void RemoveExcessIntPoints( const double tol );
181 bool GetIsOutBefore( multiset< IntersectionPoint >::iterator ip, bool prevIsOut );
183 // --------------------------------------------------------------------------
185 * \brief Iterator on the parallel grid lines of one direction
191 size_t _iVar1, _iVar2, _iConst;
192 string _name1, _name2, _nameConst;
194 LineIndexer( size_t sz1, size_t sz2, size_t sz3,
195 size_t iv1, size_t iv2, size_t iConst,
196 const string& nv1, const string& nv2, const string& nConst )
198 _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
199 _curInd[0] = _curInd[1] = _curInd[2] = 0;
200 _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
201 _name1 = nv1; _name2 = nv2; _nameConst = nConst;
204 size_t I() const { return _curInd[0]; }
205 size_t J() const { return _curInd[1]; }
206 size_t K() const { return _curInd[2]; }
207 void SetIJK( size_t i, size_t j, size_t k )
209 _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
213 if ( ++_curInd[_iVar1] == _size[_iVar1] )
214 _curInd[_iVar1] = 0, ++_curInd[_iVar2];
216 bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
217 size_t LineIndex () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
218 size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
219 size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
220 size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
221 void SetIndexOnLine (size_t i) { _curInd[ _iConst ] = i; }
222 size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
224 // --------------------------------------------------------------------------
226 * \brief Container of GridLine's
230 vector< double > _coords[3]; // coordinates of grid nodes
231 vector< GridLine > _lines [3]; // in 3 directions
232 double _tol, _minCellSize;
234 vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
235 vector< bool > _isBndNode; // is mesh node at intersection with geometry
237 size_t CellIndex( size_t i, size_t j, size_t k ) const
239 return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
241 size_t NodeIndex( size_t i, size_t j, size_t k ) const
243 return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
245 size_t NodeIndexDX() const { return 1; }
246 size_t NodeIndexDY() const { return _coords[0].size(); }
247 size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
249 LineIndexer GetLineIndexer(size_t iDir) const;
251 void SetCoordinates(const vector<double>& xCoords,
252 const vector<double>& yCoords,
253 const vector<double>& zCoords,
254 const TopoDS_Shape& shape );
255 void ComputeNodes(SMESH_MesherHelper& helper);
257 #ifdef ELLIPSOLID_WORKAROUND
258 // --------------------------------------------------------------------------
260 * \brief struct temporary replacing IntCurvesFace_Intersector until
261 * OCCT bug 0022809 is fixed
262 * http://tracker.dev.opencascade.org/view.php?id=22809
264 struct TMP_IntCurvesFace_Intersector
266 BRepAdaptor_Surface _surf;
268 BRepIntCurveSurface_Inter _intcs;
269 vector<IntCurveSurface_IntersectionPoint> _points;
270 BRepTopAdaptor_TopolTool _clsf;
272 TMP_IntCurvesFace_Intersector(const TopoDS_Face& face, const double tol)
273 :_surf( face ), _tol( tol ), _clsf( new BRepAdaptor_HSurface(_surf) ) {}
274 Bnd_Box Bounding() const { Bnd_Box b; BRepBndLib::Add (_surf.Face(), b); return b; }
275 void Perform( const gp_Lin& line, const double w0, const double w1 )
278 for ( _intcs.Init( _surf.Face(), line, _tol ); _intcs.More(); _intcs.Next() )
279 if ( w0 <= _intcs.W() && _intcs.W() <= w1 )
280 _points.push_back( _intcs.Point() );
282 bool IsDone() const { return true; }
283 int NbPnt() const { return _points.size(); }
284 IntCurveSurface_TransitionOnCurve Transition( const int i ) const { return _points[ i-1 ].Transition(); }
285 double WParameter( const int i ) const { return _points[ i-1 ].W(); }
286 TopAbs_State ClassifyUVPoint(const gp_Pnt2d& p) { return _clsf.Classify( p, _tol ); }
288 #define __IntCurvesFace_Intersector TMP_IntCurvesFace_Intersector
290 #define __IntCurvesFace_Intersector IntCurvesFace_Intersector
292 // --------------------------------------------------------------------------
294 * \brief Intersector of TopoDS_Face with all GridLine's
296 struct FaceGridIntersector
301 __IntCurvesFace_Intersector* _surfaceInt;
302 vector< std::pair< GridLine*, IntersectionPoint > > _intersections;
304 FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
306 bool IsInGrid(const Bnd_Box& gridBox);
308 void StoreIntersections()
310 for ( size_t i = 0; i < _intersections.size(); ++i )
311 _intersections[i].first->_intPoints.insert( _intersections[i].second );
313 const Bnd_Box& GetFaceBndBox()
315 GetCurveFaceIntersector();
318 __IntCurvesFace_Intersector* GetCurveFaceIntersector()
322 _surfaceInt = new __IntCurvesFace_Intersector( _face, Precision::PConfusion() );
323 _bndBox = _surfaceInt->Bounding();
324 if ( _bndBox.IsVoid() )
325 BRepBndLib::Add (_face, _bndBox);
329 bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
331 // --------------------------------------------------------------------------
333 * \brief Intersector of a surface with a GridLine
335 struct FaceLineIntersector
338 double _u, _v, _w; // params on the face and the line
339 Transition _transition; // transition of at intersection (see IntCurveSurface.cdl)
340 Transition _transIn, _transOut; // IN and OUT transitions depending of face orientation
343 gp_Cylinder _cylinder;
347 __IntCurvesFace_Intersector* _surfaceInt;
349 vector< IntersectionPoint > _intPoints;
351 void IntersectWithPlane (const GridLine& gridLine);
352 void IntersectWithCylinder(const GridLine& gridLine);
353 void IntersectWithCone (const GridLine& gridLine);
354 void IntersectWithSphere (const GridLine& gridLine);
355 void IntersectWithTorus (const GridLine& gridLine);
356 void IntersectWithSurface (const GridLine& gridLine);
358 bool UVIsOnFace() const;
359 void addIntPoint(const bool toClassify=true);
360 bool isParamOnLineOK( const double linLength )
362 return -_tol < _w && _w < linLength + _tol;
364 FaceLineIntersector():_surfaceInt(0) {}
365 ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
367 // --------------------------------------------------------------------------
369 * \brief Class representing topology of the hexahedron and creating a mesh
370 * volume basing on analysis of hexahedron intersection with geometry
374 // --------------------------------------------------------------------------------
377 // --------------------------------------------------------------------------------
378 struct _Node //!< node either at a hexahedron corner or at GridLine intersection
380 const SMDS_MeshNode* _node; // mesh node at hexahedron corner
381 const IntersectionPoint* _intPoint;
383 _Node(const SMDS_MeshNode* n=0, const IntersectionPoint* ip=0):_node(n), _intPoint(ip) {}
384 const SMDS_MeshNode* Node() const { return _intPoint ? _intPoint->_node : _node; }
385 //bool IsCorner() const { return _node; }
387 // --------------------------------------------------------------------------------
388 struct _Link // link connecting two _Node's
391 vector< _Node> _intNodes; // _Node's at GridLine intersections
392 vector< _Link > _splits;
393 vector< _Face*> _faces;
395 // --------------------------------------------------------------------------------
400 _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
401 void Reverse() { _reverse = !_reverse; }
402 int NbResultLinks() const { return _link->_splits.size(); }
403 _OrientedLink ResultLink(int i) const
405 return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
407 _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
408 _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
410 // --------------------------------------------------------------------------------
413 vector< _OrientedLink > _links;
414 vector< _Link > _polyLinks; // links added to close a polygonal face
416 // --------------------------------------------------------------------------------
417 struct _volumeDef // holder of nodes of a volume mesh element
419 vector< const SMDS_MeshNode* > _nodes;
420 vector< int > _quantities;
421 typedef boost::shared_ptr<_volumeDef> Ptr;
422 void set( const vector< const SMDS_MeshNode* >& nodes,
423 const vector< int > quant = vector< int >() )
424 { _nodes = nodes; _quantities = quant; }
425 // static Ptr New( const vector< const SMDS_MeshNode* >& nodes,
426 // const vector< int > quant = vector< int >() )
428 // _volumeDef* def = new _volumeDef;
429 // def->_nodes = nodes;
430 // def->_quantities = quant;
431 // return Ptr( def );
435 // topology of a hexahedron
441 // faces resulted from hexahedron intersection
442 vector< _Face > _polygons;
444 // computed volume elements
445 //vector< _volumeDef::Ptr > _volumeDefs;
446 _volumeDef _volumeDefs;
449 double _sizeThreshold, _sideLength[3];
450 int _nbCornerNodes, _nbIntNodes, _nbBndNodes;
451 int _origNodeInd; // index of _hexNodes[0] node within the _grid
455 Hexahedron(const double sizeThreshold, Grid* grid);
456 int MakeElements(SMESH_MesherHelper& helper);
457 void ComputeElements();
458 void Init() { init( _i, _j, _k ); }
461 Hexahedron(const Hexahedron& other );
462 void init( size_t i, size_t j, size_t k );
463 void init( size_t i );
464 int addElements(SMESH_MesherHelper& helper);
465 bool isInHole() const;
466 bool checkPolyhedronSize() const;
474 // --------------------------------------------------------------------------
476 * \brief Hexahedron computing volumes in one thread
478 struct ParallelHexahedron
480 vector< Hexahedron* >& _hexVec;
482 ParallelHexahedron( vector< Hexahedron* >& hv, vector<int>& ind): _hexVec(hv), _index(ind) {}
483 void operator() ( const tbb::blocked_range<size_t>& r ) const
485 for ( size_t i = r.begin(); i != r.end(); ++i )
486 if ( Hexahedron* hex = _hexVec[ _index[i]] )
487 hex->ComputeElements();
490 // --------------------------------------------------------------------------
492 * \brief Structure intersecting certain nb of faces with GridLine's in one thread
494 struct ParallelIntersector
496 vector< FaceGridIntersector >& _faceVec;
497 ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
498 void operator() ( const tbb::blocked_range<size_t>& r ) const
500 for ( size_t i = r.begin(); i != r.end(); ++i )
501 _faceVec[i].Intersect();
506 //=============================================================================
507 // Implementation of internal utils
508 //=============================================================================
510 * Remove coincident intersection points
512 void GridLine::RemoveExcessIntPoints( const double tol )
514 if ( _intPoints.size() < 2 ) return;
516 set< Transition > tranSet;
517 multiset< IntersectionPoint >::iterator ip1, ip2 = _intPoints.begin();
518 while ( ip2 != _intPoints.end() )
522 while ( ip2->_paramOnLine - ip1->_paramOnLine <= tol && ip2 != _intPoints.end())
524 tranSet.insert( ip1->_transition );
525 tranSet.insert( ip2->_transition );
526 _intPoints.erase( ip1 );
529 if ( tranSet.size() > 1 ) // points with different transition coincide
531 bool isIN = tranSet.count( Trans_IN );
532 bool isOUT = tranSet.count( Trans_OUT );
534 (*ip1)._transition = Trans_TANGENT;
536 (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
540 //================================================================================
542 * Return "is OUT" state for nodes before the given intersection point
544 bool GridLine::GetIsOutBefore( multiset< IntersectionPoint >::iterator ip, bool prevIsOut )
546 if ( ip->_transition == Trans_IN )
548 if ( ip->_transition == Trans_OUT )
550 if ( ip->_transition == Trans_APEX )
552 // singularity point (apex of a cone)
553 if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
555 multiset< IntersectionPoint >::iterator ipBef = ip, ipAft = ++ip;
556 if ( ipAft == _intPoints.end() )
559 if ( ipBef->_transition != ipAft->_transition )
560 return ( ipBef->_transition == Trans_OUT );
561 return ( ipBef->_transition != Trans_OUT );
563 return prevIsOut; // _transition == Trans_TANGENT
565 //================================================================================
567 * Return an iterator on GridLine's in a given direction
569 LineIndexer Grid::GetLineIndexer(size_t iDir) const
571 const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
572 const string s[] = { "X", "Y", "Z" };
573 LineIndexer li( _coords[0].size(), _coords[1].size(), _coords[2].size(),
574 indices[iDir*3], indices[iDir*3+1], indices[iDir*3+2],
575 s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
578 //=============================================================================
580 * Creates GridLine's of the grid
582 void Grid::SetCoordinates(const vector<double>& xCoords,
583 const vector<double>& yCoords,
584 const vector<double>& zCoords,
585 const TopoDS_Shape& shape)
587 _coords[0] = xCoords;
588 _coords[1] = yCoords;
589 _coords[2] = zCoords;
592 _minCellSize = Precision::Infinite();
593 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
595 for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
597 double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
598 if ( cellLen < _minCellSize )
599 _minCellSize = cellLen;
602 if ( _minCellSize < Precision::Confusion() )
603 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
604 SMESH_Comment("Too small cell size: ") << _tol );
605 _tol = _minCellSize / 1000.;
607 // attune grid extremities to shape bounding box computed by vertices
609 for ( TopExp_Explorer vExp( shape, TopAbs_VERTEX ); vExp.More(); vExp.Next() )
610 shapeBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vExp.Current() )));
612 double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
613 shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
614 double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
615 &_coords[0].back(), &_coords[1].back(), &_coords[2].back() };
616 for ( int i = 0; i < 6; ++i )
617 if ( fabs( sP[i] - *cP[i] ) < _tol )
618 *cP[i] = sP[i] + _tol/1000. * ( i < 3 ? +1 : -1 );
621 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
623 LineIndexer li = GetLineIndexer( iDir );
624 _lines[iDir].resize( li.NbLines() );
625 double len = _coords[ iDir ].back() - _coords[iDir].front();
626 gp_Vec dir( iDir==0, iDir==1, iDir==2 );
627 for ( ; li.More(); ++li )
629 GridLine& gl = _lines[iDir][ li.LineIndex() ];
630 gl._line.SetLocation(gp_Pnt(_coords[0][li.I()], _coords[1][li.J()], _coords[2][li.K()]));
631 gl._line.SetDirection( dir );
636 //================================================================================
640 void Grid::ComputeNodes(SMESH_MesherHelper& helper)
642 // state of each node of the grid relative to the geomerty
643 const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
644 vector< bool > isNodeOut( nbGridNodes, false );
645 _nodes.resize( nbGridNodes, 0 );
646 _isBndNode.resize( nbGridNodes, false );
648 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
650 LineIndexer li = GetLineIndexer( iDir );
652 // find out a shift of node index while walking along a GridLine in this direction
653 li.SetIndexOnLine( 0 );
654 size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
655 li.SetIndexOnLine( 1 );
656 const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
658 const vector<double> & coords = _coords[ iDir ];
659 for ( ; li.More(); ++li ) // loop on lines in iDir
661 li.SetIndexOnLine( 0 );
662 nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
664 GridLine& line = _lines[ iDir ][ li.LineIndex() ];
665 line.RemoveExcessIntPoints( _tol );
666 multiset< IntersectionPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
667 multiset< IntersectionPoint >::iterator ip = intPnts.begin();
670 const double* nodeCoord = & coords[0], *coord0 = nodeCoord, *coordEnd = coord0 + coords.size();
671 double nodeParam = 0;
672 for ( ; ip != intPnts.end(); ++ip )
674 // set OUT state or just skip IN nodes before ip
675 if ( nodeParam < ip->_paramOnLine - _tol )
677 isOut = line.GetIsOutBefore( ip, isOut );
679 while ( nodeParam < ip->_paramOnLine - _tol )
682 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
683 if ( ++nodeCoord < coordEnd )
684 nodeParam = *nodeCoord - *coord0;
688 if ( nodeCoord == coordEnd ) break;
690 // create a mesh node on a GridLine at ip if it does not coincide with a grid node
691 if ( nodeParam > ip->_paramOnLine + _tol )
693 li.SetIndexOnLine( 0 );
694 double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
695 xyz[ li._iConst ] += ip->_paramOnLine;
696 ip->_node = helper.AddNode( xyz[0], xyz[1], xyz[2] );
697 ip->_indexOnLine = nodeCoord-coord0-1;
699 // create a mesh node at ip concident with a grid node
702 int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
703 if ( ! _nodes[ nodeIndex ] )
705 li.SetIndexOnLine( nodeCoord-coord0 );
706 double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
707 _nodes[ nodeIndex ] = helper.AddNode( xyz[0], xyz[1], xyz[2] );
708 _isBndNode[ nodeIndex ] = true;
710 //ip->_node = _nodes[ nodeIndex ];
711 ip->_indexOnLine = nodeCoord-coord0;
712 if ( ++nodeCoord < coordEnd )
713 nodeParam = *nodeCoord - *coord0;
716 // set OUT state to nodes after the last ip
717 for ( ; nodeCoord < coordEnd; ++nodeCoord )
718 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
722 // Create mesh nodes at !OUT nodes of the grid
724 for ( size_t z = 0; z < _coords[2].size(); ++z )
725 for ( size_t y = 0; y < _coords[1].size(); ++y )
726 for ( size_t x = 0; x < _coords[0].size(); ++x )
728 size_t nodeIndex = NodeIndex( x, y, z );
729 if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
730 _nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
734 // check validity of transitions
735 const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
736 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
738 LineIndexer li = GetLineIndexer( iDir );
739 for ( ; li.More(); ++li )
741 multiset< IntersectionPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
742 if ( intPnts.empty() ) continue;
743 if ( intPnts.size() == 1 )
745 if ( intPnts.begin()->_transition != Trans_TANGENT &&
746 intPnts.begin()->_transition != Trans_APEX )
747 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
748 SMESH_Comment("Wrong SOLE transition of GridLine (")
749 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
750 << ") along " << li._nameConst
751 << ": " << trName[ intPnts.begin()->_transition] );
755 if ( intPnts.begin()->_transition == Trans_OUT )
756 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
757 SMESH_Comment("Wrong START transition of GridLine (")
758 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
759 << ") along " << li._nameConst
760 << ": " << trName[ intPnts.begin()->_transition ]);
761 if ( intPnts.rbegin()->_transition == Trans_IN )
762 throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
763 SMESH_Comment("Wrong END transition of GridLine (")
764 << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
765 << ") along " << li._nameConst
766 << ": " << trName[ intPnts.rbegin()->_transition ]);
773 //=============================================================================
775 * Checks if the face is encosed by the grid
777 bool FaceGridIntersector::IsInGrid(const Bnd_Box& gridBox)
779 double x0,y0,z0, x1,y1,z1;
780 const Bnd_Box& faceBox = GetFaceBndBox();
781 faceBox.Get(x0,y0,z0, x1,y1,z1);
783 if ( !gridBox.IsOut( gp_Pnt( x0,y0,z0 )) &&
784 !gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
787 double X0,Y0,Z0, X1,Y1,Z1;
788 gridBox.Get(X0,Y0,Z0, X1,Y1,Z1);
789 double faceP[6] = { x0,y0,z0, x1,y1,z1 };
790 double gridP[6] = { X0,Y0,Z0, X1,Y1,Z1 };
791 gp_Dir axes[3] = { gp::DX(), gp::DY(), gp::DZ() };
792 for ( int iDir = 0; iDir < 6; ++iDir )
794 if ( iDir < 3 && gridP[ iDir ] <= faceP[ iDir ] ) continue;
795 if ( iDir >= 3 && gridP[ iDir ] >= faceP[ iDir ] ) continue;
797 // check if the face intersects a side of a gridBox
799 gp_Pnt p = iDir < 3 ? gp_Pnt( X0,Y0,Z0 ) : gp_Pnt( X1,Y1,Z1 );
800 gp_Ax1 norm( p, axes[ iDir % 3 ] );
801 if ( iDir < 3 ) norm.Reverse();
803 gp_XYZ O = norm.Location().XYZ(), N = norm.Direction().XYZ();
805 TopLoc_Location loc = _face.Location();
806 Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(_face,loc);
807 if ( !aPoly.IsNull() )
809 if ( !loc.IsIdentity() )
811 norm.Transform( loc.Transformation().Inverted() );
812 O = norm.Location().XYZ(), N = norm.Direction().XYZ();
814 const double deflection = aPoly->Deflection();
816 const TColgp_Array1OfPnt& nodes = aPoly->Nodes();
817 for ( int i = nodes.Lower(); i <= nodes.Upper(); ++i )
818 if (( nodes( i ).XYZ() - O ) * N > _grid->_tol + deflection )
823 BRepAdaptor_Surface surf( _face );
824 double u0, u1, v0, v1, du, dv, u, v;
825 BRepTools::UVBounds( _face, u0, u1, v0, v1);
826 if ( surf.GetType() == GeomAbs_Plane ) {
827 du = u1 - u0, dv = v1 - v0;
830 du = surf.UResolution( _grid->_minCellSize / 10. );
831 dv = surf.VResolution( _grid->_minCellSize / 10. );
833 for ( u = u0, v = v0; u <= u1 && v <= v1; u += du, v += dv )
835 gp_Pnt p = surf.Value( u, v );
836 if (( p.XYZ() - O ) * N > _grid->_tol )
838 TopAbs_State state = GetCurveFaceIntersector()->ClassifyUVPoint(gp_Pnt2d( u, v ));
839 if ( state == TopAbs_IN || state == TopAbs_ON )
847 //=============================================================================
849 * Intersects TopoDS_Face with all GridLine's
851 void FaceGridIntersector::Intersect()
853 FaceLineIntersector intersector;
854 intersector._surfaceInt = GetCurveFaceIntersector();
855 intersector._tol = _grid->_tol;
856 intersector._transOut = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
857 intersector._transIn = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
859 typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
860 PIntFun interFunction;
862 BRepAdaptor_Surface surf( _face );
863 switch ( surf.GetType() ) {
865 intersector._plane = surf.Plane();
866 interFunction = &FaceLineIntersector::IntersectWithPlane;
868 case GeomAbs_Cylinder:
869 intersector._cylinder = surf.Cylinder();
870 interFunction = &FaceLineIntersector::IntersectWithCylinder;
873 intersector._cone = surf.Cone();
874 interFunction = &FaceLineIntersector::IntersectWithCone;
877 intersector._sphere = surf.Sphere();
878 interFunction = &FaceLineIntersector::IntersectWithSphere;
881 intersector._torus = surf.Torus();
882 interFunction = &FaceLineIntersector::IntersectWithTorus;
885 interFunction = &FaceLineIntersector::IntersectWithSurface;
888 _intersections.clear();
889 for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
891 if ( surf.GetType() == GeomAbs_Plane )
893 // check if all lines in this direction are parallel to a plane
894 if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
895 Precision::Angular()))
897 // find out a transition, that is the same for all lines of a direction
898 gp_Dir plnNorm = intersector._plane.Axis().Direction();
899 gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
900 intersector._transition =
901 ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
903 if ( surf.GetType() == GeomAbs_Cylinder )
905 // check if all lines in this direction are parallel to a cylinder
906 if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
907 Precision::Angular()))
911 // intersect the grid lines with the face
912 for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
914 GridLine& gridLine = _grid->_lines[iDir][iL];
915 if ( _bndBox.IsOut( gridLine._line )) continue;
917 intersector._intPoints.clear();
918 (intersector.*interFunction)( gridLine );
919 for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
920 _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
924 //================================================================================
926 * Return true if (_u,_v) is on the face
928 bool FaceLineIntersector::UVIsOnFace() const
930 TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
931 return ( state == TopAbs_IN || state == TopAbs_ON );
933 //================================================================================
935 * Store an intersection if it is IN or ON the face
937 void FaceLineIntersector::addIntPoint(const bool toClassify)
939 if ( !toClassify || UVIsOnFace() )
943 p._transition = _transition;
944 _intPoints.push_back( p );
947 //================================================================================
949 * Intersect a line with a plane
951 void FaceLineIntersector::IntersectWithPlane (const GridLine& gridLine)
953 IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
954 _w = linPlane.ParamOnConic(1);
955 if ( isParamOnLineOK( gridLine._length ))
957 ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
961 //================================================================================
963 * Intersect a line with a cylinder
965 void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
967 IntAna_IntConicQuad linCylinder( gridLine._line,_cylinder);
968 if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
970 _w = linCylinder.ParamOnConic(1);
971 if ( linCylinder.NbPoints() == 1 )
972 _transition = Trans_TANGENT;
974 _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
975 if ( isParamOnLineOK( gridLine._length ))
977 ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
980 if ( linCylinder.NbPoints() > 1 )
982 _w = linCylinder.ParamOnConic(2);
983 if ( isParamOnLineOK( gridLine._length ))
985 ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
986 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
992 //================================================================================
994 * Intersect a line with a cone
996 void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
998 IntAna_IntConicQuad linCone(gridLine._line,_cone);
999 if ( !linCone.IsDone() ) return;
1001 gp_Vec du, dv, norm;
1002 for ( int i = 1; i <= linCone.NbPoints(); ++i )
1004 _w = linCone.ParamOnConic( i );
1005 if ( !isParamOnLineOK( gridLine._length )) continue;
1006 ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1009 ElSLib::D1( _u, _v, _cone, P, du, dv );
1011 double normSize2 = norm.SquareMagnitude();
1012 if ( normSize2 > Precision::Angular() * Precision::Angular() )
1014 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1015 cos /= sqrt( normSize2 );
1016 if ( cos < -Precision::Angular() )
1017 _transition = _transIn;
1018 else if ( cos > Precision::Angular() )
1019 _transition = _transOut;
1021 _transition = Trans_TANGENT;
1025 _transition = Trans_APEX;
1027 addIntPoint( /*toClassify=*/false);
1031 //================================================================================
1033 * Intersect a line with a sphere
1035 void FaceLineIntersector::IntersectWithSphere (const GridLine& gridLine)
1037 IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1038 if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1040 _w = linSphere.ParamOnConic(1);
1041 if ( linSphere.NbPoints() == 1 )
1042 _transition = Trans_TANGENT;
1044 _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1045 if ( isParamOnLineOK( gridLine._length ))
1047 ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1050 if ( linSphere.NbPoints() > 1 )
1052 _w = linSphere.ParamOnConic(2);
1053 if ( isParamOnLineOK( gridLine._length ))
1055 ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1056 _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1062 //================================================================================
1064 * Intersect a line with a torus
1066 void FaceLineIntersector::IntersectWithTorus (const GridLine& gridLine)
1068 IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1069 if ( !linTorus.IsDone()) return;
1071 gp_Vec du, dv, norm;
1072 for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1074 _w = linTorus.ParamOnLine( i );
1075 if ( !isParamOnLineOK( gridLine._length )) continue;
1076 linTorus.ParamOnTorus( i, _u,_v );
1079 ElSLib::D1( _u, _v, _torus, P, du, dv );
1081 double normSize = norm.Magnitude();
1082 double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1084 if ( cos < -Precision::Angular() )
1085 _transition = _transIn;
1086 else if ( cos > Precision::Angular() )
1087 _transition = _transOut;
1089 _transition = Trans_TANGENT;
1090 addIntPoint( /*toClassify=*/false);
1094 //================================================================================
1096 * Intersect a line with a non-analytical surface
1098 void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1100 _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1101 if ( !_surfaceInt->IsDone() ) return;
1102 for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1104 _transition = Transition( _surfaceInt->Transition( i ) );
1105 _w = _surfaceInt->WParameter( i );
1106 addIntPoint(/*toClassify=*/false);
1109 //================================================================================
1111 * check if its face can be safely intersected in a thread
1113 bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1118 TopLoc_Location loc;
1119 Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1120 Handle(Geom_RectangularTrimmedSurface) ts =
1121 Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1122 while( !ts.IsNull() ) {
1123 surf = ts->BasisSurface();
1124 ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1126 if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1127 surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1128 if ( !noSafeTShapes.insert((const Standard_Transient*) _face.TShape() ).second )
1132 TopExp_Explorer exp( _face, TopAbs_EDGE );
1133 for ( ; exp.More(); exp.Next() )
1135 bool edgeIsSafe = true;
1136 const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1139 Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1142 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1143 while( !tc.IsNull() ) {
1144 c = tc->BasisCurve();
1145 tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1147 if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1148 c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1155 Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1158 Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1159 while( !tc.IsNull() ) {
1160 c2 = tc->BasisCurve();
1161 tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1163 if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1164 c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1168 if ( !edgeIsSafe && !noSafeTShapes.insert((const Standard_Transient*) e.TShape() ).second )
1173 //================================================================================
1175 * \brief Creates topology of the hexahedron
1177 Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1178 : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbIntNodes(0)
1180 _polygons.reserve(100); // to avoid reallocation;
1182 //set nodes shift within grid->_nodes from the node 000
1183 size_t dx = _grid->NodeIndexDX();
1184 size_t dy = _grid->NodeIndexDY();
1185 size_t dz = _grid->NodeIndexDZ();
1187 size_t i100 = i000 + dx;
1188 size_t i010 = i000 + dy;
1189 size_t i110 = i010 + dx;
1190 size_t i001 = i000 + dz;
1191 size_t i101 = i100 + dz;
1192 size_t i011 = i010 + dz;
1193 size_t i111 = i110 + dz;
1194 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1195 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1196 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1197 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1198 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1199 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1200 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1201 _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1203 vector< int > idVec;
1204 // set nodes to links
1205 for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1207 SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1208 _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1209 link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1210 link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1211 link._intNodes.reserve( 10 ); // to avoid reallocation
1212 link._splits.reserve( 10 );
1215 // set links to faces
1216 int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1217 for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1219 SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1220 _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1221 bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1222 faceID == SMESH_Block::ID_Fx1z ||
1223 faceID == SMESH_Block::ID_F0yz );
1224 quad._links.resize(4);
1225 vector<_OrientedLink>::iterator frwLinkIt = quad._links.begin();
1226 vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1227 for ( int i = 0; i < 4; ++i )
1229 bool revLink = revFace;
1230 if ( i > 1 ) // reverse links u1 and v0
1232 _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1233 link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1238 //================================================================================
1240 * \brief Copy constructor
1242 Hexahedron::Hexahedron( const Hexahedron& other )
1243 :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbIntNodes(0)
1245 _polygons.reserve(100); // to avoid reallocation;
1247 for ( int i = 0; i < 8; ++i )
1248 _nodeShift[i] = other._nodeShift[i];
1250 for ( int i = 0; i < 12; ++i )
1252 const _Link& srcLink = other._hexLinks[ i ];
1253 _Link& tgtLink = this->_hexLinks[ i ];
1254 tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1255 tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1256 tgtLink._intNodes.reserve( 10 ); // to avoid reallocation
1257 tgtLink._splits.reserve( 10 );
1260 for ( int i = 0; i < 6; ++i )
1262 const _Face& srcQuad = other._hexQuads[ i ];
1263 _Face& tgtQuad = this->_hexQuads[ i ];
1264 tgtQuad._links.resize(4);
1265 for ( int j = 0; j < 4; ++j )
1267 const _OrientedLink& srcLink = srcQuad._links[ j ];
1268 _OrientedLink& tgtLink = tgtQuad._links[ j ];
1269 tgtLink._reverse = srcLink._reverse;
1270 tgtLink._link = _hexLinks + ( srcLink._link - other._hexLinks );
1275 //================================================================================
1277 * \brief Initializes its data by given grid cell
1279 void Hexahedron::init( size_t i, size_t j, size_t k )
1281 _i = i; _j = j; _k = k;
1282 // set nodes of grid to nodes of the hexahedron and
1283 // count nodes at hexahedron corners located IN and ON geometry
1284 _nbCornerNodes = _nbBndNodes = 0;
1285 _origNodeInd = _grid->NodeIndex( i,j,k );
1286 for ( int iN = 0; iN < 8; ++iN )
1288 _hexNodes[iN]._node = _grid->_nodes[ _origNodeInd + _nodeShift[iN] ];
1289 _nbCornerNodes += bool( _hexNodes[iN]._node );
1290 _nbBndNodes += _grid->_isBndNode[ _origNodeInd + _nodeShift[iN] ];
1293 _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1294 _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1295 _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1297 if ( _nbCornerNodes < 8 && _nbIntNodes + _nbCornerNodes > 3)
1300 // create sub-links (_splits) by splitting links with _intNodes
1301 for ( int iLink = 0; iLink < 12; ++iLink )
1303 _Link& link = _hexLinks[ iLink ];
1304 link._splits.clear();
1305 split._nodes[ 0 ] = link._nodes[0];
1306 for ( size_t i = 0; i < link._intNodes.size(); ++ i )
1308 if ( split._nodes[ 0 ]->Node() )
1310 split._nodes[ 1 ] = &link._intNodes[i];
1311 link._splits.push_back( split );
1313 split._nodes[ 0 ] = &link._intNodes[i];
1315 if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() )
1317 split._nodes[ 1 ] = link._nodes[1];
1318 link._splits.push_back( split );
1323 //================================================================================
1325 * \brief Initializes its data by given grid cell (countered from zero)
1327 void Hexahedron::init( size_t iCell )
1329 size_t iNbCell = _grid->_coords[0].size() - 1;
1330 size_t jNbCell = _grid->_coords[1].size() - 1;
1331 _i = iCell % iNbCell;
1332 _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1333 _k = iCell / iNbCell / jNbCell;
1337 //================================================================================
1339 * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1341 void Hexahedron::ComputeElements()
1345 if ( _nbCornerNodes + _nbIntNodes < 4 )
1348 if ( _nbBndNodes == _nbCornerNodes && isInHole() )
1353 vector<const SMDS_MeshNode* > polyhedraNodes;
1354 vector<int> quantities;
1356 // create polygons from quadrangles and get their nodes
1358 vector<_Node*> nodes;
1359 nodes.reserve( _nbCornerNodes + _nbIntNodes );
1362 polyLink._faces.reserve( 1 );
1364 for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1366 const _Face& quad = _hexQuads[ iF ] ;
1368 _polygons.resize( _polygons.size() + 1 );
1369 _Face& polygon = _polygons.back();
1370 polygon._links.clear();
1371 polygon._polyLinks.clear(); polygon._polyLinks.reserve( 10 );
1373 // add splits of a link to a polygon and collect info on nodes
1374 //int nbIn = 0, nbOut = 0, nbCorners = 0;
1376 for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1378 int nbSpits = quad._links[ iE ].NbResultLinks();
1379 for ( int iS = 0; iS < nbSpits; ++iS )
1381 _OrientedLink split = quad._links[ iE ].ResultLink( iS );
1382 _Node* n = split.FirstNode();
1383 if ( !polygon._links.empty() )
1385 _Node* nPrev = polygon._links.back().LastNode();
1388 polyLink._nodes[0] = nPrev;
1389 polyLink._nodes[1] = n;
1390 polygon._polyLinks.push_back( polyLink );
1391 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
1392 nodes.push_back( nPrev );
1395 polygon._links.push_back( split );
1396 nodes.push_back( n );
1399 if ( polygon._links.size() > 1 )
1401 _Node* n1 = polygon._links.back().LastNode();
1402 _Node* n2 = polygon._links.front().FirstNode();
1405 polyLink._nodes[0] = n1;
1406 polyLink._nodes[1] = n2;
1407 polygon._polyLinks.push_back( polyLink );
1408 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
1409 nodes.push_back( n1 );
1411 // add polygon to its links
1412 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1413 polygon._links[ iL ]._link->_faces.push_back( &polygon );
1414 // store polygon nodes
1415 quantities.push_back( nodes.size() );
1416 for ( size_t i = 0; i < nodes.size(); ++i )
1417 polyhedraNodes.push_back( nodes[i]->Node() );
1421 _polygons.resize( _polygons.size() - 1 );
1425 // create polygons closing holes in a polyhedron
1428 vector< _OrientedLink* > freeLinks;
1429 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1431 _Face& polygon = _polygons[ iP ];
1432 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1433 if ( polygon._links[ iL ]._link->_faces.size() < 2 )
1434 freeLinks.push_back( & polygon._links[ iL ]);
1436 // make closed chains of free links
1437 int nbFreeLinks = freeLinks.size();
1438 if ( 0 < nbFreeLinks && nbFreeLinks < 3 ) return;
1439 while ( nbFreeLinks > 0 )
1442 _polygons.resize( _polygons.size() + 1 );
1443 _Face& polygon = _polygons.back();
1444 polygon._links.clear();
1446 // get a remaining link to start from
1447 _OrientedLink* curLink = 0;
1448 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1449 if (( curLink = freeLinks[ iL ] ))
1450 freeLinks[ iL ] = 0;
1451 nodes.push_back( curLink->LastNode() );
1452 polygon._links.push_back( *curLink );
1454 // find all links connected to curLink
1458 curNode = curLink->FirstNode();
1460 for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1461 if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
1463 curLink = freeLinks[ iL ];
1464 freeLinks[ iL ] = 0;
1465 nodes.push_back( curNode );
1466 polygon._links.push_back( *curLink );
1468 } while ( curLink );
1470 nbFreeLinks -= polygon._links.size();
1472 if ( curNode != nodes.front() || polygon._links.size() < 3 )
1473 return; // closed polygon not found -> invalid polyhedron
1475 quantities.push_back( nodes.size() );
1476 for ( size_t i = 0; i < nodes.size(); ++i )
1477 polyhedraNodes.push_back( nodes[i]->Node() );
1479 // add polygon to its links and reverse links
1480 for ( size_t i = 0; i < polygon._links.size(); ++i )
1482 polygon._links[i].Reverse();
1483 polygon._links[i]._link->_faces.push_back( &polygon );
1486 //const size_t firstPoly = _polygons.size();
1489 if ( ! checkPolyhedronSize() )
1494 // create a classic cell if possible
1495 const int nbNodes = _nbCornerNodes + _nbIntNodes;
1496 bool isClassicElem = false;
1497 if ( nbNodes == 8 && _polygons.size() == 6 ) isClassicElem = addHexa();
1498 else if ( nbNodes == 4 && _polygons.size() == 4 ) isClassicElem = addTetra();
1499 else if ( nbNodes == 6 && _polygons.size() == 5 ) isClassicElem = addPenta();
1500 else if ( nbNodes == 5 && _polygons.size() == 5 ) isClassicElem = addPyra ();
1501 if ( !isClassicElem )
1502 _volumeDefs.set( polyhedraNodes, quantities );
1504 //================================================================================
1506 * \brief Create elements in the mesh
1508 int Hexahedron::MakeElements(SMESH_MesherHelper& helper)
1510 SMESHDS_Mesh* mesh = helper.GetMeshDS();
1512 size_t nbCells[3] = { _grid->_coords[0].size() - 1,
1513 _grid->_coords[1].size() - 1,
1514 _grid->_coords[2].size() - 1 };
1515 const size_t nbGridCells = nbCells[0] *nbCells [1] * nbCells[2];
1516 vector< Hexahedron* > intersectedHex( nbGridCells, 0 );
1519 // set intersection nodes from GridLine's to links of intersectedHex
1520 int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
1521 for ( int iDir = 0; iDir < 3; ++iDir )
1523 int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
1524 dInd[1][ iDirOther[iDir][0] ] = -1;
1525 dInd[2][ iDirOther[iDir][1] ] = -1;
1526 dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
1527 // loop on GridLine's parallel to iDir
1528 LineIndexer lineInd = _grid->GetLineIndexer( iDir );
1529 for ( ; lineInd.More(); ++lineInd )
1531 GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
1532 multiset< IntersectionPoint >::const_iterator ip = line._intPoints.begin();
1533 for ( ; ip != line._intPoints.end(); ++ip )
1535 if ( !ip->_node ) continue;
1536 lineInd.SetIndexOnLine( ip->_indexOnLine );
1537 for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
1539 i = int(lineInd.I()) + dInd[iL][0];
1540 j = int(lineInd.J()) + dInd[iL][1];
1541 k = int(lineInd.K()) + dInd[iL][2];
1542 if ( i < 0 || i >= nbCells[0] ||
1543 j < 0 || j >= nbCells[1] ||
1544 k < 0 || k >= nbCells[2] ) continue;
1546 const size_t hexIndex = _grid->CellIndex( i,j,k );
1547 Hexahedron *& hex = intersectedHex[ hexIndex ];
1550 hex = new Hexahedron( *this );
1556 const int iLink = iL + iDir * 4;
1557 hex->_hexLinks[iLink]._intNodes.push_back( _Node( 0, &(*ip) ));
1564 // add not split hexadrons to the mesh
1566 vector<int> intHexInd( nbIntHex );
1568 for ( size_t i = 0; i < intersectedHex.size(); ++i )
1570 Hexahedron * & hex = intersectedHex[ i ];
1573 intHexInd[ nbIntHex++ ] = i;
1574 if ( hex->_nbIntNodes > 0 ) continue;
1575 init( hex->_i, hex->_j, hex->_k );
1581 if ( _nbCornerNodes == 8 && ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
1583 // order of _hexNodes is defined by enum SMESH_Block::TShapeID
1584 SMDS_MeshElement* el =
1585 mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
1586 _hexNodes[3].Node(), _hexNodes[1].Node(),
1587 _hexNodes[4].Node(), _hexNodes[6].Node(),
1588 _hexNodes[7].Node(), _hexNodes[5].Node() );
1589 mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
1594 intersectedHex[ i ] = 0;
1598 else if ( _nbCornerNodes > 3 && !hex )
1600 // all intersection of hex with geometry are at grid nodes
1601 hex = new Hexahedron( *this );
1603 intHexInd.push_back(0);
1604 intHexInd[ nbIntHex++ ] = i;
1608 // add elements resulted from hexadron intersection
1610 intHexInd.resize( nbIntHex );
1611 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, nbIntHex ),
1612 ParallelHexahedron( intersectedHex, intHexInd ),
1613 tbb::simple_partitioner()); // ComputeElements() is called here
1614 for ( size_t i = 0; i < intHexInd.size(); ++i )
1615 if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
1616 nbAdded += hex->addElements( helper );
1618 for ( size_t i = 0; i < intHexInd.size(); ++i )
1619 if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
1621 hex->ComputeElements();
1622 nbAdded += hex->addElements( helper );
1626 for ( size_t i = 0; i < intersectedHex.size(); ++i )
1627 if ( intersectedHex[ i ] )
1628 delete intersectedHex[ i ];
1633 //================================================================================
1635 * \brief Adds computed elements to the mesh
1637 int Hexahedron::addElements(SMESH_MesherHelper& helper)
1640 // add elements resulted from hexahedron intersection
1641 //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
1643 vector< const SMDS_MeshNode* >& nodes = _volumeDefs._nodes;
1645 if ( !_volumeDefs._quantities.empty() )
1647 helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
1651 switch ( nodes.size() )
1653 case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
1654 nodes[4],nodes[5],nodes[6],nodes[7] );
1656 case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
1658 case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
1661 helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
1665 nbAdded += int ( _volumeDefs._nodes.size() > 0 );
1670 //================================================================================
1672 * \brief Return true if the element is in a hole
1674 bool Hexahedron::isInHole() const
1676 const int ijk[3] = { _i, _j, _k };
1677 IntersectionPoint curIntPnt;
1679 // consider a cell to be in a hole if all links in any direction
1680 // comes OUT of geometry
1681 for ( int iDir = 0; iDir < 3; ++iDir )
1683 const vector<double>& coords = _grid->_coords[ iDir ];
1684 LineIndexer li = _grid->GetLineIndexer( iDir );
1685 li.SetIJK( _i,_j,_k );
1686 size_t lineIndex[4] = { li.LineIndex (),
1690 bool allLinksOut = true, hasLinks = false;
1691 for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
1693 const _Link& link = _hexLinks[ iL + 4*iDir ];
1694 // check transition of the first node of a link
1695 const IntersectionPoint* firstIntPnt = 0;
1696 if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
1698 curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0];
1699 const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
1700 multiset< IntersectionPoint >::const_iterator ip =
1701 line._intPoints.upper_bound( curIntPnt );
1703 firstIntPnt = &(*ip);
1705 else if ( !link._intNodes.empty() )
1707 firstIntPnt = link._intNodes[0]._intPoint;
1713 allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
1716 if ( hasLinks && allLinksOut )
1722 //================================================================================
1724 * \brief Return true if a polyhedron passes _sizeThreshold criterion
1726 bool Hexahedron::checkPolyhedronSize() const
1729 for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1731 const _Face& polygon = _polygons[iP];
1732 gp_XYZ area (0,0,0);
1733 SMESH_TNodeXYZ p1 ( polygon._links[ 0 ].FirstNode()->Node() );
1734 for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1736 SMESH_TNodeXYZ p2 ( polygon._links[ iL ].LastNode()->Node() );
1740 volume += p1 * area;
1744 double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
1746 return volume > initVolume / _sizeThreshold;
1748 //================================================================================
1750 * \brief Tries to create a hexahedron
1752 bool Hexahedron::addHexa()
1754 if ( _polygons[0]._links.size() != 4 ||
1755 _polygons[1]._links.size() != 4 ||
1756 _polygons[2]._links.size() != 4 ||
1757 _polygons[3]._links.size() != 4 ||
1758 _polygons[4]._links.size() != 4 ||
1759 _polygons[5]._links.size() != 4 )
1761 const SMDS_MeshNode* nodes[8];
1763 for ( int iL = 0; iL < 4; ++iL )
1766 nodes[iL] = _polygons[0]._links[iL].FirstNode()->Node();
1769 // find a top node above the base node
1770 _Link* link = _polygons[0]._links[iL]._link;
1771 ASSERT( link->_faces.size() > 1 );
1772 // a quadrangle sharing <link> with _polygons[0]
1773 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
1774 for ( int i = 0; i < 4; ++i )
1775 if ( quad->_links[i]._link == link )
1777 // 1st node of a link opposite to <link> in <quad>
1778 nodes[iL+4] = quad->_links[(i+2)%4].FirstNode()->Node();
1784 _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+8 ));
1788 //================================================================================
1790 * \brief Tries to create a tetrahedron
1792 bool Hexahedron::addTetra()
1794 const SMDS_MeshNode* nodes[4];
1795 nodes[0] = _polygons[0]._links[0].FirstNode()->Node();
1796 nodes[1] = _polygons[0]._links[1].FirstNode()->Node();
1797 nodes[2] = _polygons[0]._links[2].FirstNode()->Node();
1799 _Link* link = _polygons[0]._links[0]._link;
1800 ASSERT( link->_faces.size() > 1 );
1802 // a triangle sharing <link> with _polygons[0]
1803 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
1804 for ( int i = 0; i < 3; ++i )
1805 if ( tria->_links[i]._link == link )
1807 nodes[3] = tria->_links[(i+1)%3].LastNode()->Node();
1808 _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+4 ));
1814 //================================================================================
1816 * \brief Tries to create a pentahedron
1818 bool Hexahedron::addPenta()
1820 // find a base triangular face
1822 for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
1823 if ( _polygons[ iF ]._links.size() == 3 )
1825 if ( iTri < 0 ) return false;
1828 const SMDS_MeshNode* nodes[6];
1830 for ( int iL = 0; iL < 3; ++iL )
1833 nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode()->Node();
1836 // find a top node above the base node
1837 _Link* link = _polygons[ iTri ]._links[iL]._link;
1838 ASSERT( link->_faces.size() > 1 );
1839 // a quadrangle sharing <link> with a base triangle
1840 _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
1841 if ( quad->_links.size() != 4 ) return false;
1842 for ( int i = 0; i < 4; ++i )
1843 if ( quad->_links[i]._link == link )
1845 // 1st node of a link opposite to <link> in <quad>
1846 nodes[iL+3] = quad->_links[(i+2)%4].FirstNode()->Node();
1852 _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+6 ));
1854 return ( nbN == 6 );
1856 //================================================================================
1858 * \brief Tries to create a pyramid
1860 bool Hexahedron::addPyra()
1862 // find a base quadrangle
1864 for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
1865 if ( _polygons[ iF ]._links.size() == 4 )
1867 if ( iQuad < 0 ) return false;
1870 const SMDS_MeshNode* nodes[5];
1871 nodes[0] = _polygons[iQuad]._links[0].FirstNode()->Node();
1872 nodes[1] = _polygons[iQuad]._links[1].FirstNode()->Node();
1873 nodes[2] = _polygons[iQuad]._links[2].FirstNode()->Node();
1874 nodes[3] = _polygons[iQuad]._links[3].FirstNode()->Node();
1876 _Link* link = _polygons[iQuad]._links[0]._link;
1877 ASSERT( link->_faces.size() > 1 );
1879 // a triangle sharing <link> with a base quadrangle
1880 _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
1881 if ( tria->_links.size() != 3 ) return false;
1882 for ( int i = 0; i < 3; ++i )
1883 if ( tria->_links[i]._link == link )
1885 nodes[4] = tria->_links[(i+1)%3].LastNode()->Node();
1886 _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+5 ));
1895 //=============================================================================
1897 * \brief Generates 3D structured Cartesian mesh in the internal part of
1898 * solid shapes and polyhedral volumes near the shape boundary.
1899 * \param theMesh - mesh to fill in
1900 * \param theShape - a compound of all SOLIDs to mesh
1901 * \retval bool - true in case of success
1903 //=============================================================================
1905 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh & theMesh,
1906 const TopoDS_Shape & theShape)
1908 // The algorithm generates the mesh in following steps:
1910 // 1) Intersection of grid lines with the geometry boundary.
1911 // This step allows to find out if a given node of the initial grid is
1912 // inside or outside the geometry.
1914 // 2) For each cell of the grid, check how many of it's nodes are outside
1915 // of the geometry boundary. Depending on a result of this check
1916 // - skip a cell, if all it's nodes are outside
1917 // - skip a cell, if it is too small according to the size threshold
1918 // - add a hexahedron in the mesh, if all nodes are inside
1919 // - add a polyhedron in the mesh, if some nodes are inside and some outside
1921 _computeCanceled = false;
1927 TopTools_MapOfShape faceMap;
1928 for ( TopExp_Explorer fExp( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
1929 if ( !faceMap.Add( fExp.Current() ))
1930 faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
1933 vector<FaceGridIntersector> facesItersectors( faceMap.Extent() );
1934 TopTools_MapIteratorOfMapOfShape faceMppIt( faceMap );
1935 for ( int i = 0; faceMppIt.More(); faceMppIt.Next(), ++i )
1937 facesItersectors[i]._face = TopoDS::Face( faceMppIt.Key() );
1938 facesItersectors[i]._grid = &grid;
1939 shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
1942 vector<double> xCoords, yCoords, zCoords;
1943 _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
1945 grid.SetCoordinates( xCoords, yCoords, zCoords, theShape );
1947 // check if the grid encloses the shape
1948 if ( !_hyp->IsGridBySpacing(0) ||
1949 !_hyp->IsGridBySpacing(1) ||
1950 !_hyp->IsGridBySpacing(2) )
1953 gridBox.Add( gp_Pnt( xCoords[0], yCoords[0], zCoords[0] ));
1954 gridBox.Add( gp_Pnt( xCoords.back(), yCoords.back(), zCoords.back() ));
1955 double x0,y0,z0, x1,y1,z1;
1956 shapeBox.Get(x0,y0,z0, x1,y1,z1);
1957 if ( gridBox.IsOut( gp_Pnt( x0,y0,z0 )) ||
1958 gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
1959 for ( size_t i = 0; i < facesItersectors.size(); ++i )
1961 if ( !facesItersectors[i].IsInGrid( gridBox ))
1962 return error("The grid doesn't enclose the geometry");
1963 #ifdef ELLIPSOLID_WORKAROUND
1964 delete facesItersectors[i]._surfaceInt, facesItersectors[i]._surfaceInt = 0;
1968 if ( _computeCanceled ) return false;
1971 { // copy partner faces and curves of not thread-safe types
1972 set< const Standard_Transient* > tshapes;
1973 BRepBuilderAPI_Copy copier;
1974 for ( size_t i = 0; i < facesItersectors.size(); ++i )
1976 if ( !facesItersectors[i].IsThreadSafe(tshapes) )
1978 copier.Perform( facesItersectors[i]._face );
1979 facesItersectors[i]._face = TopoDS::Face( copier );
1983 // Intersection of grid lines with the geometry boundary.
1984 tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
1985 ParallelIntersector( facesItersectors ),
1986 tbb::simple_partitioner());
1988 for ( size_t i = 0; i < facesItersectors.size(); ++i )
1989 facesItersectors[i].Intersect();
1992 // put interesection points onto the GridLine's; this is done after intersection
1993 // to avoid contention of facesItersectors for writing into the same GridLine
1994 // in case of parallel work of facesItersectors
1995 for ( size_t i = 0; i < facesItersectors.size(); ++i )
1996 facesItersectors[i].StoreIntersections();
1998 SMESH_MesherHelper helper( theMesh );
1999 TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
2000 helper.SetSubShape( solidExp.Current() );
2001 helper.SetElementsOnShape( true );
2003 if ( _computeCanceled ) return false;
2005 // create nodes on the geometry
2006 grid.ComputeNodes(helper);
2008 if ( _computeCanceled ) return false;
2010 // create volume elements
2011 Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
2012 int nbAdded = hex.MakeElements( helper );
2014 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
2017 // make all SOLIDS computed
2018 if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
2020 SMDS_ElemIteratorPtr volIt = sm1->GetElements();
2021 for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
2023 const SMDS_MeshElement* vol = volIt->next();
2024 sm1->RemoveElement( vol, /*isElemDeleted=*/false );
2025 meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
2028 // make other sub-shapes computed
2029 setSubmeshesComputed( theMesh, theShape );
2032 // remove free nodes
2033 if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
2035 // intersection nodes
2036 for ( int iDir = 0; iDir < 3; ++iDir )
2038 vector< GridLine >& lines = grid._lines[ iDir ];
2039 for ( size_t i = 0; i < lines.size(); ++i )
2041 multiset< IntersectionPoint >::iterator ip = lines[i]._intPoints.begin();
2042 for ( ; ip != lines[i]._intPoints.end(); ++ip )
2043 if ( ip->_node && ip->_node->NbInverseElements() == 0 )
2044 meshDS->RemoveFreeNode( ip->_node, smDS, /*fromGroups=*/false );
2048 for ( size_t i = 0; i < grid._nodes.size(); ++i )
2049 if ( !grid._isBndNode[i] ) // nodes on boundary are already removed
2050 if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
2051 meshDS->RemoveFreeNode( grid._nodes[i], smDS, /*fromGroups=*/false );
2057 // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
2058 catch ( SMESH_ComputeError& e)
2060 return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
2065 //=============================================================================
2069 //=============================================================================
2071 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh & theMesh,
2072 const TopoDS_Shape & theShape,
2073 MapShapeNbElems& theResMap)
2076 // std::vector<int> aResVec(SMDSEntity_Last);
2077 // for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
2078 // if(IsQuadratic) {
2079 // aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
2080 // int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
2081 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
2084 // aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
2085 // aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
2087 // SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
2088 // aResMap.insert(std::make_pair(sm,aResVec));
2093 //=============================================================================
2097 * \brief Event listener setting/unsetting _alwaysComputed flag to
2098 * submeshes of inferior levels to prevent their computing
2100 struct _EventListener : public SMESH_subMeshEventListener
2104 _EventListener(const string& algoName):
2105 SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
2108 // --------------------------------------------------------------------------------
2109 // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
2111 static void setAlwaysComputed( const bool isComputed,
2112 SMESH_subMesh* subMeshOfSolid)
2114 SMESH_subMeshIteratorPtr smIt =
2115 subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
2116 while ( smIt->more() )
2118 SMESH_subMesh* sm = smIt->next();
2119 sm->SetIsAlwaysComputed( isComputed );
2123 // --------------------------------------------------------------------------------
2124 // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
2126 virtual void ProcessEvent(const int event,
2127 const int eventType,
2128 SMESH_subMesh* subMeshOfSolid,
2129 SMESH_subMeshEventListenerData* data,
2130 const SMESH_Hypothesis* hyp = 0)
2132 if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
2134 setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
2139 SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
2140 if ( !algo3D || _algoName != algo3D->GetName() )
2141 setAlwaysComputed( false, subMeshOfSolid );
2145 // --------------------------------------------------------------------------------
2146 // set the event listener
2148 static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
2150 subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
2155 }; // struct _EventListener
2159 //================================================================================
2161 * \brief Sets event listener to submeshes if necessary
2162 * \param subMesh - submesh where algo is set
2163 * This method is called when a submesh gets HYP_OK algo_state.
2164 * After being set, event listener is notified on each event of a submesh.
2166 //================================================================================
2168 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
2170 _EventListener::SetOn( subMesh, GetName() );
2173 //================================================================================
2175 * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
2177 //================================================================================
2179 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh& theMesh,
2180 const TopoDS_Shape& theShape)
2182 for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
2183 _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));