1 // Copyright (C) 2007-2010 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
23 // File: SMESH_MesherHelper.cxx
24 // Created: 15.02.06 15:22:41
25 // Author: Sergey KUUL
27 #include "SMESH_MesherHelper.hxx"
29 #include "SMDS_FacePosition.hxx"
30 #include "SMDS_EdgePosition.hxx"
31 #include "SMDS_VolumeTool.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_ProxyMesh.hxx"
35 #include <BRepAdaptor_Surface.hxx>
36 #include <BRepTools.hxx>
37 #include <BRepTools_WireExplorer.hxx>
38 #include <BRep_Tool.hxx>
39 #include <Geom2d_Curve.hxx>
40 #include <GeomAPI_ProjectPointOnCurve.hxx>
41 #include <GeomAPI_ProjectPointOnSurf.hxx>
42 #include <Geom_Curve.hxx>
43 #include <Geom_Surface.hxx>
44 #include <ShapeAnalysis.hxx>
46 #include <TopExp_Explorer.hxx>
47 #include <TopTools_ListIteratorOfListOfShape.hxx>
48 #include <TopTools_MapIteratorOfMapOfShape.hxx>
49 #include <TopTools_MapOfShape.hxx>
52 #include <gp_Pnt2d.hxx>
53 #include <gp_Trsf.hxx>
55 #include <Standard_Failure.hxx>
56 #include <Standard_ErrorHandler.hxx>
58 #include <utilities.h>
64 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
68 gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
70 enum { U_periodic = 1, V_periodic = 2 };
73 //================================================================================
77 //================================================================================
79 SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
80 : myParIndex(0), myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false)
82 myPar1[0] = myPar2[0] = myPar1[1] = myPar2[1] = 0;
83 mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
86 //=======================================================================
87 //function : ~SMESH_MesherHelper
89 //=======================================================================
91 SMESH_MesherHelper::~SMESH_MesherHelper()
94 TID2ProjectorOnSurf::iterator i_proj = myFace2Projector.begin();
95 for ( ; i_proj != myFace2Projector.end(); ++i_proj )
96 delete i_proj->second;
99 TID2ProjectorOnCurve::iterator i_proj = myEdge2Projector.begin();
100 for ( ; i_proj != myEdge2Projector.end(); ++i_proj )
101 delete i_proj->second;
105 //=======================================================================
106 //function : IsQuadraticSubMesh
107 //purpose : Check submesh for given shape: if all elements on this shape
108 // are quadratic, quadratic elements will be created.
109 // Also fill myTLinkNodeMap
110 //=======================================================================
112 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
114 SMESHDS_Mesh* meshDS = GetMeshDS();
115 // we can create quadratic elements only if all elements
116 // created on subshapes of given shape are quadratic
117 // also we have to fill myTLinkNodeMap
118 myCreateQuadratic = true;
119 mySeamShapeIds.clear();
120 myDegenShapeIds.clear();
121 TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
122 SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
124 int nbOldLinks = myTLinkNodeMap.size();
126 if ( !myMesh->HasShapeToMesh() )
128 if (( myCreateQuadratic = myMesh->NbFaces( ORDER_QUADRATIC )))
130 SMDS_FaceIteratorPtr fIt = meshDS->facesIterator();
131 while ( fIt->more() )
132 AddTLinks( static_cast< const SMDS_MeshFace* >( fIt->next() ));
137 TopExp_Explorer exp( aSh, subType );
138 TopTools_MapOfShape checkedSubShapes;
139 for (; exp.More() && myCreateQuadratic; exp.Next()) {
140 if ( !checkedSubShapes.Add( exp.Current() ))
141 continue; // needed if aSh is compound of solids
142 if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
143 if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
145 const SMDS_MeshElement* e = it->next();
146 if ( e->GetType() != elemType || !e->IsQuadratic() ) {
147 myCreateQuadratic = false;
152 switch ( e->NbNodes() ) {
154 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
156 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
157 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
158 AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
160 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
161 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
162 AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
163 AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
166 myCreateQuadratic = false;
176 if ( nbOldLinks == myTLinkNodeMap.size() )
177 myCreateQuadratic = false;
179 if(!myCreateQuadratic) {
180 myTLinkNodeMap.clear();
184 return myCreateQuadratic;
187 //=======================================================================
188 //function : SetSubShape
189 //purpose : Set geomerty to make elements on
190 //=======================================================================
192 void SMESH_MesherHelper::SetSubShape(const int aShID)
194 if ( aShID == myShapeID )
197 SetSubShape( GetMeshDS()->IndexToShape( aShID ));
199 SetSubShape( TopoDS_Shape() );
202 //=======================================================================
203 //function : SetSubShape
204 //purpose : Set geomerty to create elements on
205 //=======================================================================
207 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
209 if ( myShape.IsSame( aSh ))
213 mySeamShapeIds.clear();
214 myDegenShapeIds.clear();
216 if ( myShape.IsNull() ) {
220 SMESHDS_Mesh* meshDS = GetMeshDS();
221 myShapeID = meshDS->ShapeToIndex(aSh);
224 // treatment of periodic faces
225 for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
227 const TopoDS_Face& face = TopoDS::Face( eF.Current() );
228 BRepAdaptor_Surface surface( face );
229 if ( surface.IsUPeriodic() || surface.IsVPeriodic() )
231 for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
233 // look for a seam edge
234 const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
235 if ( BRep_Tool::IsClosed( edge, face )) {
236 // initialize myPar1, myPar2 and myParIndex
238 BRep_Tool::UVPoints( edge, face, uv1, uv2 );
239 if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
241 myParIndex |= U_periodic;
242 myPar1[0] = surface.FirstUParameter();
243 myPar2[0] = surface.LastUParameter();
246 myParIndex |= V_periodic;
247 myPar1[1] = surface.FirstVParameter();
248 myPar2[1] = surface.LastVParameter();
250 // store seam shape indices, negative if shape encounters twice
251 int edgeID = meshDS->ShapeToIndex( edge );
252 mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
253 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
254 int vertexID = meshDS->ShapeToIndex( v.Current() );
255 mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
259 // look for a degenerated edge
260 if ( BRep_Tool::Degenerated( edge )) {
261 myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
262 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
263 myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
270 //=======================================================================
271 //function : GetNodeUVneedInFaceNode
272 //purpose : Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
273 // Return true if the face is periodic.
274 // If F is Null, answer about subshape set through IsQuadraticSubMesh() or
276 //=======================================================================
278 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
280 if ( F.IsNull() ) return !mySeamShapeIds.empty();
282 if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
283 return !mySeamShapeIds.empty();
286 Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
287 if ( !aSurface.IsNull() )
288 return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
293 //=======================================================================
294 //function : IsMedium
296 //=======================================================================
298 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode* node,
299 const SMDSAbs_ElementType typeToCheck)
301 return SMESH_MeshEditor::IsMedium( node, typeToCheck );
304 //=======================================================================
305 //function : GetSubShapeByNode
306 //purpose : Return support shape of a node
307 //=======================================================================
309 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
310 const SMESHDS_Mesh* meshDS)
312 int shapeID = node->getshapeId();
313 if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
314 return meshDS->IndexToShape( shapeID );
316 return TopoDS_Shape();
320 //=======================================================================
321 //function : AddTLinkNode
322 //purpose : add a link in my data structure
323 //=======================================================================
325 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
326 const SMDS_MeshNode* n2,
327 const SMDS_MeshNode* n12)
329 // add new record to map
330 SMESH_TLink link( n1, n2 );
331 myTLinkNodeMap.insert( make_pair(link,n12));
334 //================================================================================
336 * \brief Add quadratic links of edge to own data structure
338 //================================================================================
340 void SMESH_MesherHelper::AddTLinks(const SMDS_MeshEdge* edge)
342 if ( edge->IsQuadratic() )
343 AddTLinkNode(edge->GetNode(0), edge->GetNode(1), edge->GetNode(2));
346 //================================================================================
348 * \brief Add quadratic links of face to own data structure
350 //================================================================================
352 void SMESH_MesherHelper::AddTLinks(const SMDS_MeshFace* f)
355 switch ( f->NbNodes() ) {
357 AddTLinkNode(f->GetNode(0),f->GetNode(1),f->GetNode(3));
358 AddTLinkNode(f->GetNode(1),f->GetNode(2),f->GetNode(4));
359 AddTLinkNode(f->GetNode(2),f->GetNode(0),f->GetNode(5)); break;
361 AddTLinkNode(f->GetNode(0),f->GetNode(1),f->GetNode(4));
362 AddTLinkNode(f->GetNode(1),f->GetNode(2),f->GetNode(5));
363 AddTLinkNode(f->GetNode(2),f->GetNode(3),f->GetNode(6));
364 AddTLinkNode(f->GetNode(3),f->GetNode(0),f->GetNode(7));
369 //================================================================================
371 * \brief Add quadratic links of volume to own data structure
373 //================================================================================
375 void SMESH_MesherHelper::AddTLinks(const SMDS_MeshVolume* volume)
377 if ( volume->IsQuadratic() )
379 SMDS_VolumeTool vTool( volume );
380 const SMDS_MeshNode** nodes = vTool.GetNodes();
382 for ( int iF = 1; iF < vTool.NbFaces(); ++iF )
384 const int nbN = vTool.NbFaceNodes( iF );
385 const int* iNodes = vTool.GetFaceNodesIndices( iF );
386 for ( int i = 0; i < nbN; )
388 int iN1 = iNodes[i++];
389 int iN12 = iNodes[i++];
390 int iN2 = iNodes[i++];
391 if ( iN1 > iN2 ) std::swap( iN1, iN2 );
392 int linkID = iN1 * vTool.NbNodes() + iN2;
393 pair< set<int>::iterator, bool > it_isNew = addedLinks.insert( linkID );
394 if ( it_isNew.second )
395 AddTLinkNode( nodes[iN1], nodes[iN2], nodes[iN12] );
397 addedLinks.erase( it_isNew.first ); // each link encounters only twice
403 //================================================================================
405 * \brief Return true if position of nodes on the shape hasn't yet been checked or
406 * the positions proved to be invalid
408 //================================================================================
410 bool SMESH_MesherHelper::toCheckPosOnShape(int shapeID ) const
412 map< int,bool >::const_iterator id_ok = myNodePosShapesValidity.find( shapeID );
413 return ( id_ok == myNodePosShapesValidity.end() || !id_ok->second );
416 //================================================================================
418 * \brief Set validity of positions of nodes on the shape.
419 * Once set, validity is not changed
421 //================================================================================
423 void SMESH_MesherHelper::setPosOnShapeValidity(int shapeID, bool ok ) const
425 ((SMESH_MesherHelper*)this)->myNodePosShapesValidity.insert( make_pair( shapeID, ok));
428 //=======================================================================
429 //function : GetUVOnSeam
430 //purpose : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
431 //=======================================================================
433 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
435 gp_Pnt2d result = uv1;
436 for ( int i = U_periodic; i <= V_periodic ; ++i )
438 if ( myParIndex & i )
440 double p1 = uv1.Coord( i );
441 double dp1 = Abs( p1-myPar1[i-1]), dp2 = Abs( p1-myPar2[i-1]);
442 if ( myParIndex == i ||
443 dp1 < ( myPar2[i-1] - myPar2[i-1] ) / 100. ||
444 dp2 < ( myPar2[i-1] - myPar2[i-1] ) / 100. )
446 double p2 = uv2.Coord( i );
447 double p1Alt = ( dp1 < dp2 ) ? myPar2[i-1] : myPar1[i-1];
448 if ( Abs( p2 - p1 ) > Abs( p2 - p1Alt ))
449 result.SetCoord( i, p1Alt );
456 //=======================================================================
457 //function : GetNodeUV
458 //purpose : Return node UV on face
459 //=======================================================================
461 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
462 const SMDS_MeshNode* n,
463 const SMDS_MeshNode* n2,
466 gp_Pnt2d uv( Precision::Infinite(), Precision::Infinite() );
467 const SMDS_PositionPtr Pos = n->GetPosition();
469 if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
471 // node has position on face
472 const SMDS_FacePosition* fpos =
473 static_cast<const SMDS_FacePosition*>(n->GetPosition());
474 uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
476 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ));
478 else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
480 // node has position on edge => it is needed to find
481 // corresponding edge from face, get pcurve for this
482 // edge and retrieve value from this pcurve
483 const SMDS_EdgePosition* epos =
484 static_cast<const SMDS_EdgePosition*>(n->GetPosition());
485 int edgeID = n->getshapeId();
486 TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
487 double f, l, u = epos->GetUParameter();
488 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
489 bool validU = ( f < u && u < l );
491 uv = C2d->Value( u );
493 uv.SetCoord( Precision::Infinite(),0.);
494 if ( check || !validU )
495 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ),/*force=*/ !validU );
497 // for a node on a seam edge select one of UVs on 2 pcurves
498 if ( n2 && IsSeamShape( edgeID ) )
500 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0, check ));
503 { // adjust uv to period
505 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
506 Standard_Boolean isUPeriodic = S->IsUPeriodic();
507 Standard_Boolean isVPeriodic = S->IsVPeriodic();
508 if ( isUPeriodic || isVPeriodic ) {
509 Standard_Real UF,UL,VF,VL;
510 S->Bounds(UF,UL,VF,VL);
512 uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
514 uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
518 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
520 if ( int vertexID = n->getshapeId() ) {
521 const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
523 uv = BRep_Tool::Parameters( V, F );
526 catch (Standard_Failure& exc) {
529 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() )
530 uvOK = ( V == vert.Current() );
533 MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
534 << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
536 // get UV of a vertex closest to the node
538 gp_Pnt pn = XYZ( n );
539 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) {
540 TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
541 gp_Pnt p = BRep_Tool::Pnt( curV );
542 double curDist = p.SquareDistance( pn );
543 if ( curDist < dist ) {
545 uv = BRep_Tool::Parameters( curV, F );
546 uvOK = ( dist < DBL_MIN );
552 TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
553 for ( ; it.More(); it.Next() ) {
554 if ( it.Value().ShapeType() == TopAbs_EDGE ) {
555 const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
557 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
558 if ( !C2d.IsNull() ) {
559 double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
560 uv = C2d->Value( u );
568 if ( n2 && IsSeamShape( vertexID ) )
569 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
574 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ));
583 //=======================================================================
584 //function : CheckNodeUV
585 //purpose : Check and fix node UV on a face
586 //=======================================================================
588 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
589 const SMDS_MeshNode* n,
593 double distXYZ[4]) const
595 int shapeID = n->getshapeId();
596 bool infinit = ( Precision::IsInfinite( uv.X() ) || Precision::IsInfinite( uv.Y() ));
597 if ( force || toCheckPosOnShape( shapeID ) || infinit )
599 // check that uv is correct
601 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
602 gp_Pnt nodePnt = XYZ( n ), surfPnt(0,0,0);
604 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
606 (dist = nodePnt.Distance( surfPnt = surface->Value( uv.X(), uv.Y() ))) > tol )
608 setPosOnShapeValidity( shapeID, false );
609 if ( !infinit && distXYZ ) {
610 surfPnt.Transform( loc );
612 distXYZ[1] = surfPnt.X(); distXYZ[2] = surfPnt.Y(); distXYZ[3]=surfPnt.Z();
614 // uv incorrect, project the node to surface
615 GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
616 projector.Perform( nodePnt );
617 if ( !projector.IsDone() || projector.NbPoints() < 1 )
619 MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
622 Quantity_Parameter U,V;
623 projector.LowerDistanceParameters(U,V);
625 surfPnt = surface->Value( U, V );
626 dist = nodePnt.Distance( surfPnt );
628 surfPnt.Transform( loc );
630 distXYZ[1] = surfPnt.X(); distXYZ[2] = surfPnt.Y(); distXYZ[3]=surfPnt.Z();
634 MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
637 // store the fixed UV on the face
638 if ( myShape.IsSame(F) && shapeID == myShapeID )
639 const_cast<SMDS_MeshNode*>(n)->SetPosition
640 ( SMDS_PositionPtr( new SMDS_FacePosition( U, V )));
642 else if ( uv.Modulus() > numeric_limits<double>::min() )
644 setPosOnShapeValidity( shapeID, true );
650 //=======================================================================
651 //function : GetProjector
652 //purpose : Return projector intitialized by given face without location, which is returned
653 //=======================================================================
655 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
656 TopLoc_Location& loc,
659 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
660 int faceID = GetMeshDS()->ShapeToIndex( F );
661 TID2ProjectorOnSurf& i2proj = const_cast< TID2ProjectorOnSurf&>( myFace2Projector );
662 TID2ProjectorOnSurf::iterator i_proj = i2proj.find( faceID );
663 if ( i_proj == i2proj.end() )
665 if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
666 double U1, U2, V1, V2;
667 surface->Bounds(U1, U2, V1, V2);
668 GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
669 proj->Init( surface, U1, U2, V1, V2, tol );
670 i_proj = i2proj.insert( make_pair( faceID, proj )).first;
672 return *( i_proj->second );
677 gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
678 gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
679 gp_XY_FunPtr(Subtracted);
682 //=======================================================================
683 //function : applyIn2D
684 //purpose : Perform given operation on two 2d points in parameric space of given surface.
685 // It takes into account period of the surface. Use gp_XY_FunPtr macro
686 // to easily define pointer to function of gp_XY class.
687 //=======================================================================
689 gp_XY SMESH_MesherHelper::applyIn2D(const Handle(Geom_Surface)& surface,
693 const bool resultInPeriod)
695 Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
696 Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
697 if ( !isUPeriodic && !isVPeriodic )
700 // move uv2 not far than half-period from uv1
702 uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
704 uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
707 gp_XY res = fun( uv1, gp_XY(u2,v2) );
709 // move result within period
710 if ( resultInPeriod )
712 Standard_Real UF,UL,VF,VL;
713 surface->Bounds(UF,UL,VF,VL);
715 res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
717 res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
722 //=======================================================================
723 //function : GetMiddleUV
724 //purpose : Return middle UV taking in account surface period
725 //=======================================================================
727 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
731 return applyIn2D( surface, p1, p2, & AverageUV );
734 //=======================================================================
735 //function : GetNodeU
736 //purpose : Return node U on edge
737 //=======================================================================
739 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
740 const SMDS_MeshNode* n,
741 const SMDS_MeshNode* inEdgeNode,
745 const SMDS_PositionPtr pos = n->GetPosition();
746 if ( pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
748 const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( pos );
749 param = epos->GetUParameter();
751 else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
753 if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
756 BRep_Tool::Range( E, f,l );
757 double uInEdge = GetNodeU( E, inEdgeNode );
758 param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
762 SMESHDS_Mesh * meshDS = GetMeshDS();
763 int vertexID = n->getshapeId();
764 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
765 param = BRep_Tool::Parameter( V, E );
770 double tol = BRep_Tool::Tolerance( E );
771 double f,l; BRep_Tool::Range( E, f,l );
772 bool force = ( param < f-tol || param > l+tol );
773 if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
774 force = ( GetMeshDS()->ShapeToIndex( E ) != n->getshapeId() );
776 *check = CheckNodeU( E, n, param, 2*tol, force );
781 //=======================================================================
782 //function : CheckNodeU
783 //purpose : Check and fix node U on an edge
784 // Return false if U is bad and could not be fixed
785 //=======================================================================
787 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
788 const SMDS_MeshNode* n,
792 double distXYZ[4]) const
794 int shapeID = n->getshapeId();
795 if ( force || toCheckPosOnShape( shapeID ))
797 TopLoc_Location loc; double f,l;
798 Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
799 if ( curve.IsNull() ) // degenerated edge
801 if ( u+tol < f || u-tol > l )
803 double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
809 gp_Pnt nodePnt = SMESH_TNodeXYZ( n );
810 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
811 gp_Pnt curvPnt = curve->Value( u );
812 double dist = nodePnt.Distance( curvPnt );
814 curvPnt.Transform( loc );
816 distXYZ[1] = curvPnt.X(); distXYZ[2] = curvPnt.Y(); distXYZ[3]=curvPnt.Z();
820 setPosOnShapeValidity( shapeID, false );
821 // u incorrect, project the node to the curve
822 int edgeID = GetMeshDS()->ShapeToIndex( E );
823 TID2ProjectorOnCurve& i2proj = const_cast< TID2ProjectorOnCurve&>( myEdge2Projector );
824 TID2ProjectorOnCurve::iterator i_proj =
825 i2proj.insert( make_pair( edgeID, (GeomAPI_ProjectPointOnCurve*) 0 )).first;
826 if ( !i_proj->second )
828 i_proj->second = new GeomAPI_ProjectPointOnCurve();
829 i_proj->second->Init( curve, f, l );
831 GeomAPI_ProjectPointOnCurve* projector = i_proj->second;
832 projector->Perform( nodePnt );
833 if ( projector->NbPoints() < 1 )
835 MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
838 Quantity_Parameter U = projector->LowerDistanceParameter();
840 curvPnt = curve->Value( u );
841 dist = nodePnt.Distance( curvPnt );
843 curvPnt.Transform( loc );
845 distXYZ[1] = curvPnt.X(); distXYZ[2] = curvPnt.Y(); distXYZ[3]=curvPnt.Z();
849 MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
850 MESSAGE("distance " << dist << " " << tol );
853 // store the fixed U on the edge
854 if ( myShape.IsSame(E) && shapeID == myShapeID )
855 const_cast<SMDS_MeshNode*>(n)->SetPosition
856 ( SMDS_PositionPtr( new SMDS_EdgePosition( U )));
858 else if ( fabs( u ) > numeric_limits<double>::min() )
860 setPosOnShapeValidity( shapeID, true );
862 if (( u < f-tol || u > l+tol ) && force )
864 // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
867 // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
868 double period = curve->Period();
869 u = ( u < f ) ? u + period : u - period;
871 catch (Standard_Failure& exc)
881 //=======================================================================
882 //function : GetMediumNode
883 //purpose : Return existing or create new medium nodes between given ones
884 //=======================================================================
886 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
887 const SMDS_MeshNode* n2,
890 // Find existing node
892 SMESH_TLink link(n1,n2);
893 ItTLinkNode itLN = myTLinkNodeMap.find( link );
894 if ( itLN != myTLinkNodeMap.end() ) {
895 return (*itLN).second;
898 // Create medium node
901 SMESHDS_Mesh* meshDS = GetMeshDS();
903 if ( IsSeamShape( n1->getshapeId() ))
904 // to get a correct UV of a node on seam, the second node must have checked UV
907 // get type of shape for the new medium node
908 int faceID = -1, edgeID = -1;
909 const SMDS_PositionPtr Pos1 = n1->GetPosition();
910 const SMDS_PositionPtr Pos2 = n2->GetPosition();
912 TopoDS_Edge E; double u [2];
913 TopoDS_Face F; gp_XY uv[2];
914 bool uvOK[2] = { false, false };
916 if( myShape.IsNull() )
918 if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
919 faceID = n1->getshapeId();
921 else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
922 faceID = n2->getshapeId();
925 if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
926 edgeID = n1->getshapeId();
928 if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
929 edgeID = n2->getshapeId();
932 // get positions of the given nodes on shapes
933 TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
934 if ( faceID>0 || shapeType == TopAbs_FACE)
936 if( myShape.IsNull() )
937 F = TopoDS::Face(meshDS->IndexToShape(faceID));
939 F = TopoDS::Face(myShape);
942 uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
943 uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
945 else if (edgeID>0 || shapeType == TopAbs_EDGE)
947 if ( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE &&
948 Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE &&
949 n1->getshapeId() != n2->getshapeId() ) // issue 0021006
950 return getMediumNodeOnComposedWire(n1,n2,force3d);
952 if( myShape.IsNull() )
953 E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
955 E = TopoDS::Edge(myShape);
958 u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
959 u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
963 // we try to create medium node using UV parameters of
964 // nodes, else - medium between corresponding 3d points
967 if ( uvOK[0] && uvOK[1] )
969 if ( IsDegenShape( n1->getshapeId() )) {
970 if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
971 else uv[0].SetCoord( 2, uv[1].Coord( 2 ));
973 else if ( IsDegenShape( n2->getshapeId() )) {
974 if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
975 else uv[1].SetCoord( 2, uv[0].Coord( 2 ));
979 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
980 gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
981 gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
982 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
983 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
984 myTLinkNodeMap.insert(make_pair(link,n12));
988 else if ( !E.IsNull() )
991 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
994 Standard_Boolean isPeriodic = C->IsPeriodic();
997 Standard_Real Period = C->Period();
998 Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
999 Standard_Real pmid = (u[0]+p)/2.;
1000 U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
1005 gp_Pnt P = C->Value( U );
1006 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
1007 meshDS->SetNodeOnEdge(n12, edgeID, U);
1008 myTLinkNodeMap.insert(make_pair(link,n12));
1015 double x = ( n1->X() + n2->X() )/2.;
1016 double y = ( n1->Y() + n2->Y() )/2.;
1017 double z = ( n1->Z() + n2->Z() )/2.;
1018 n12 = meshDS->AddNode(x,y,z);
1022 gp_XY UV = ( uv[0] + uv[1] ) / 2.;
1023 CheckNodeUV( F, n12, UV, 2*BRep_Tool::Tolerance( F ), /*force=*/true);
1024 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
1026 else if ( !E.IsNull() )
1028 double U = ( u[0] + u[1] ) / 2.;
1029 CheckNodeU( E, n12, U, 2*BRep_Tool::Tolerance( E ), /*force=*/true);
1030 meshDS->SetNodeOnEdge(n12, edgeID, U);
1032 else if ( myShapeID > 0 )
1034 meshDS->SetNodeInVolume(n12, myShapeID);
1037 myTLinkNodeMap.insert( make_pair( link, n12 ));
1041 //================================================================================
1043 * \brief Makes a medium node if nodes reside different edges
1045 //================================================================================
1047 const SMDS_MeshNode* SMESH_MesherHelper::getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
1048 const SMDS_MeshNode* n2,
1051 gp_Pnt middle = 0.5 * XYZ(n1) + 0.5 * XYZ(n2);
1052 SMDS_MeshNode* n12 = AddNode( middle.X(), middle.Y(), middle.Z() );
1054 // To find position on edge and 3D position for n12,
1055 // project <middle> to 2 edges and select projection most close to <middle>
1057 double u = 0, distMiddleProj = Precision::Infinite(), distXYZ[4];
1059 TopoDS_Edge edges[2];
1060 for ( int is2nd = 0; is2nd < 2; ++is2nd )
1063 const SMDS_MeshNode* n = is2nd ? n2 : n1;
1064 TopoDS_Shape shape = GetSubShapeByNode( n, GetMeshDS() );
1065 if ( shape.IsNull() || shape.ShapeType() != TopAbs_EDGE )
1068 // project to get U of projection and distance from middle to projection
1069 TopoDS_Edge edge = edges[ is2nd ] = TopoDS::Edge( shape );
1070 double node2MiddleDist = middle.Distance( XYZ(n) );
1071 double foundU = GetNodeU( edge, n );
1072 CheckNodeU( edge, n12, foundU, 2*BRep_Tool::Tolerance(edge), /*force=*/true, distXYZ );
1073 if ( distXYZ[0] < node2MiddleDist )
1075 distMiddleProj = distXYZ[0];
1080 if ( Precision::IsInfinite( distMiddleProj ))
1082 // both projections failed; set n12 on the edge of n1 with U of a common vertex
1083 TopoDS_Vertex vCommon;
1084 if ( TopExp::CommonVertex( edges[0], edges[1], vCommon ))
1085 u = BRep_Tool::Parameter( vCommon, edges[0] );
1088 double f,l, u0 = GetNodeU( edges[0], n1 );
1089 BRep_Tool::Range( edges[0],f,l );
1090 u = ( fabs(u0-f) < fabs(u0-l) ) ? f : l;
1096 // move n12 to position of a successfull projection
1097 double tol = BRep_Tool::Tolerance(edges[ iOkEdge ]);
1098 if ( !force3d && distMiddleProj > 2*tol )
1100 TopLoc_Location loc; double f,l;
1101 Handle(Geom_Curve) curve = BRep_Tool::Curve( edges[iOkEdge],loc,f,l );
1102 gp_Pnt p = curve->Value( u );
1103 GetMeshDS()->MoveNode( n12, p.X(), p.Y(), p.Z() );
1106 GetMeshDS()->SetNodeOnEdge(n12, edges[iOkEdge], u);
1108 myTLinkNodeMap.insert( make_pair( SMESH_TLink(n1,n2), n12 ));
1113 //=======================================================================
1114 //function : AddNode
1115 //purpose : Creates a node
1116 //=======================================================================
1118 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
1120 SMESHDS_Mesh * meshDS = GetMeshDS();
1121 SMDS_MeshNode* node = 0;
1123 node = meshDS->AddNodeWithID( x, y, z, ID );
1125 node = meshDS->AddNode( x, y, z );
1126 if ( mySetElemOnShape && myShapeID > 0 ) {
1127 switch ( myShape.ShapeType() ) {
1128 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
1129 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
1130 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID); break;
1131 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID); break;
1132 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
1139 //=======================================================================
1140 //function : AddEdge
1141 //purpose : Creates quadratic or linear edge
1142 //=======================================================================
1144 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
1145 const SMDS_MeshNode* n2,
1149 SMESHDS_Mesh * meshDS = GetMeshDS();
1151 SMDS_MeshEdge* edge = 0;
1152 if (myCreateQuadratic) {
1153 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1155 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
1157 edge = meshDS->AddEdge(n1, n2, n12);
1161 edge = meshDS->AddEdgeWithID(n1, n2, id);
1163 edge = meshDS->AddEdge(n1, n2);
1166 if ( mySetElemOnShape && myShapeID > 0 )
1167 meshDS->SetMeshElementOnShape( edge, myShapeID );
1172 //=======================================================================
1173 //function : AddFace
1174 //purpose : Creates quadratic or linear triangle
1175 //=======================================================================
1177 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1178 const SMDS_MeshNode* n2,
1179 const SMDS_MeshNode* n3,
1183 SMESHDS_Mesh * meshDS = GetMeshDS();
1184 SMDS_MeshFace* elem = 0;
1186 if( n1==n2 || n2==n3 || n3==n1 )
1189 if(!myCreateQuadratic) {
1191 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
1193 elem = meshDS->AddFace(n1, n2, n3);
1196 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1197 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1198 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1201 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
1203 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
1205 if ( mySetElemOnShape && myShapeID > 0 )
1206 meshDS->SetMeshElementOnShape( elem, myShapeID );
1211 //=======================================================================
1212 //function : AddFace
1213 //purpose : Creates quadratic or linear quadrangle
1214 //=======================================================================
1216 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1217 const SMDS_MeshNode* n2,
1218 const SMDS_MeshNode* n3,
1219 const SMDS_MeshNode* n4,
1223 SMESHDS_Mesh * meshDS = GetMeshDS();
1224 SMDS_MeshFace* elem = 0;
1227 return AddFace(n1,n3,n4,id,force3d);
1230 return AddFace(n1,n2,n4,id,force3d);
1233 return AddFace(n1,n2,n3,id,force3d);
1236 return AddFace(n1,n2,n4,id,force3d);
1239 return AddFace(n1,n2,n3,id,force3d);
1242 return AddFace(n1,n2,n3,id,force3d);
1245 if(!myCreateQuadratic) {
1247 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
1249 elem = meshDS->AddFace(n1, n2, n3, n4);
1252 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1253 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1254 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1255 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1258 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
1260 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
1262 if ( mySetElemOnShape && myShapeID > 0 )
1263 meshDS->SetMeshElementOnShape( elem, myShapeID );
1268 //=======================================================================
1269 //function : AddPolygonalFace
1270 //purpose : Creates polygon, with additional nodes in quadratic mesh
1271 //=======================================================================
1273 SMDS_MeshFace* SMESH_MesherHelper::AddPolygonalFace (const vector<const SMDS_MeshNode*>& nodes,
1277 SMESHDS_Mesh * meshDS = GetMeshDS();
1278 SMDS_MeshFace* elem = 0;
1280 if(!myCreateQuadratic) {
1282 elem = meshDS->AddPolygonalFaceWithID(nodes, id);
1284 elem = meshDS->AddPolygonalFace(nodes);
1287 vector<const SMDS_MeshNode*> newNodes;
1288 for ( int i = 0; i < nodes.size(); ++i )
1290 const SMDS_MeshNode* n1 = nodes[i];
1291 const SMDS_MeshNode* n2 = nodes[(i+1)/nodes.size()];
1292 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1293 newNodes.push_back( n1 );
1294 newNodes.push_back( n12 );
1297 elem = meshDS->AddPolygonalFaceWithID(newNodes, id);
1299 elem = meshDS->AddPolygonalFace(newNodes);
1301 if ( mySetElemOnShape && myShapeID > 0 )
1302 meshDS->SetMeshElementOnShape( elem, myShapeID );
1307 //=======================================================================
1308 //function : AddVolume
1309 //purpose : Creates quadratic or linear prism
1310 //=======================================================================
1312 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1313 const SMDS_MeshNode* n2,
1314 const SMDS_MeshNode* n3,
1315 const SMDS_MeshNode* n4,
1316 const SMDS_MeshNode* n5,
1317 const SMDS_MeshNode* n6,
1321 SMESHDS_Mesh * meshDS = GetMeshDS();
1322 SMDS_MeshVolume* elem = 0;
1323 if(!myCreateQuadratic) {
1325 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
1327 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
1330 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1331 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1332 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1334 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1335 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1336 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1338 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1339 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1340 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1343 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
1344 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1346 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1347 n12, n23, n31, n45, n56, n64, n14, n25, n36);
1349 if ( mySetElemOnShape && myShapeID > 0 )
1350 meshDS->SetMeshElementOnShape( elem, myShapeID );
1355 //=======================================================================
1356 //function : AddVolume
1357 //purpose : Creates quadratic or linear tetrahedron
1358 //=======================================================================
1360 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1361 const SMDS_MeshNode* n2,
1362 const SMDS_MeshNode* n3,
1363 const SMDS_MeshNode* n4,
1367 SMESHDS_Mesh * meshDS = GetMeshDS();
1368 SMDS_MeshVolume* elem = 0;
1369 if(!myCreateQuadratic) {
1371 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1373 elem = meshDS->AddVolume(n1, n2, n3, n4);
1376 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1377 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1378 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1380 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1381 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1382 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1385 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1387 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1389 if ( mySetElemOnShape && myShapeID > 0 )
1390 meshDS->SetMeshElementOnShape( elem, myShapeID );
1395 //=======================================================================
1396 //function : AddVolume
1397 //purpose : Creates quadratic or linear pyramid
1398 //=======================================================================
1400 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1401 const SMDS_MeshNode* n2,
1402 const SMDS_MeshNode* n3,
1403 const SMDS_MeshNode* n4,
1404 const SMDS_MeshNode* n5,
1408 SMDS_MeshVolume* elem = 0;
1409 if(!myCreateQuadratic) {
1411 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1413 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1416 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1417 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1418 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1419 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1421 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1422 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1423 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1424 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1427 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
1432 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
1434 n15, n25, n35, n45);
1436 if ( mySetElemOnShape && myShapeID > 0 )
1437 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1442 //=======================================================================
1443 //function : AddVolume
1444 //purpose : Creates quadratic or linear hexahedron
1445 //=======================================================================
1447 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1448 const SMDS_MeshNode* n2,
1449 const SMDS_MeshNode* n3,
1450 const SMDS_MeshNode* n4,
1451 const SMDS_MeshNode* n5,
1452 const SMDS_MeshNode* n6,
1453 const SMDS_MeshNode* n7,
1454 const SMDS_MeshNode* n8,
1458 SMESHDS_Mesh * meshDS = GetMeshDS();
1459 SMDS_MeshVolume* elem = 0;
1460 if(!myCreateQuadratic) {
1462 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1464 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1467 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1468 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1469 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1470 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1472 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1473 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1474 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1475 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1477 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1478 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1479 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1480 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1483 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1484 n12, n23, n34, n41, n56, n67,
1485 n78, n85, n15, n26, n37, n48, id);
1487 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1488 n12, n23, n34, n41, n56, n67,
1489 n78, n85, n15, n26, n37, n48);
1491 if ( mySetElemOnShape && myShapeID > 0 )
1492 meshDS->SetMeshElementOnShape( elem, myShapeID );
1497 //=======================================================================
1498 //function : AddPolyhedralVolume
1499 //purpose : Creates polyhedron. In quadratic mesh, adds medium nodes
1500 //=======================================================================
1503 SMESH_MesherHelper::AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
1504 const std::vector<int>& quantities,
1508 SMESHDS_Mesh * meshDS = GetMeshDS();
1509 SMDS_MeshVolume* elem = 0;
1510 if(!myCreateQuadratic)
1513 elem = meshDS->AddPolyhedralVolumeWithID(nodes, quantities, id);
1515 elem = meshDS->AddPolyhedralVolume(nodes, quantities);
1519 vector<const SMDS_MeshNode*> newNodes;
1520 vector<int> newQuantities;
1521 for ( int iFace=0, iN=0; iFace < quantities.size(); ++iFace)
1523 int nbNodesInFace = quantities[iFace];
1524 newQuantities.push_back(0);
1525 for ( int i = 0; i < nbNodesInFace; ++i )
1527 const SMDS_MeshNode* n1 = nodes[ iN + i ];
1528 newNodes.push_back( n1 );
1529 newQuantities.back()++;
1531 const SMDS_MeshNode* n2 = nodes[ iN + ( i+1==nbNodesInFace ? 0 : i+1 )];
1532 // if ( n1->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE &&
1533 // n2->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE )
1535 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1536 newNodes.push_back( n12 );
1537 newQuantities.back()++;
1540 iN += nbNodesInFace;
1543 elem = meshDS->AddPolyhedralVolumeWithID( newNodes, newQuantities, id );
1545 elem = meshDS->AddPolyhedralVolume( newNodes, newQuantities );
1547 if ( mySetElemOnShape && myShapeID > 0 )
1548 meshDS->SetMeshElementOnShape( elem, myShapeID );
1553 //=======================================================================
1554 //function : LoadNodeColumns
1555 //purpose : Load nodes bound to face into a map of node columns
1556 //=======================================================================
1558 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1559 const TopoDS_Face& theFace,
1560 const TopoDS_Edge& theBaseEdge,
1561 SMESHDS_Mesh* theMesh,
1562 SMESH_ProxyMesh* theProxyMesh)
1564 const SMESHDS_SubMesh* faceSubMesh = 0;
1567 faceSubMesh = theProxyMesh->GetSubMesh( theFace );
1568 if ( !faceSubMesh ||
1569 faceSubMesh->NbElements() == 0 ||
1570 theProxyMesh->IsTemporary( faceSubMesh->GetElements()->next() ))
1572 // can use a proxy sub-mesh with not temporary elements only
1578 faceSubMesh = theMesh->MeshElements( theFace );
1579 if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
1582 // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
1584 map< double, const SMDS_MeshNode*> sortedBaseNodes;
1585 if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
1586 || sortedBaseNodes.size() < 2 )
1589 int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
1590 map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
1591 double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
1592 for ( ; u_n != sortedBaseNodes.end(); u_n++ )
1594 double par = ( u_n->first - f ) / range;
1595 vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
1596 nCol.resize( nbRows );
1597 nCol[0] = u_n->second;
1599 TParam2ColumnMap::iterator par_nVec_2, par_nVec_1 = theParam2ColumnMap.begin();
1602 for ( ; par_nVec_1 != theParam2ColumnMap.end(); ++par_nVec_1 )
1604 const SMDS_MeshNode* & n = par_nVec_1->second[0];
1605 n = theProxyMesh->GetProxyNode( n );
1609 // fill theParam2ColumnMap column by column by passing from nodes on
1610 // theBaseEdge up via mesh faces on theFace
1612 par_nVec_2 = theParam2ColumnMap.begin();
1613 par_nVec_1 = par_nVec_2++;
1614 TIDSortedElemSet emptySet, avoidSet;
1615 for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
1617 vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
1618 vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
1620 int i1, i2, iRow = 0;
1621 const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
1622 // find face sharing node n1 and n2 and belonging to faceSubMesh
1623 while ( const SMDS_MeshElement* face =
1624 SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
1626 if ( faceSubMesh->Contains( face ))
1628 int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
1631 n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
1632 n2 = face->GetNode( (i1+2) % 4 );
1633 if ( ++iRow >= nbRows )
1639 avoidSet.insert( face );
1641 if ( iRow + 1 < nbRows ) // compact if necessary
1642 nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
1644 return theParam2ColumnMap.size() > 1 && theParam2ColumnMap.begin()->second.size() > 1;
1647 //=======================================================================
1648 //function : NbAncestors
1649 //purpose : Return number of unique ancestors of the shape
1650 //=======================================================================
1652 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1653 const SMESH_Mesh& mesh,
1654 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
1656 TopTools_MapOfShape ancestors;
1657 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1658 for ( ; ansIt.More(); ansIt.Next() ) {
1659 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1660 ancestors.Add( ansIt.Value() );
1662 return ancestors.Extent();
1665 //=======================================================================
1666 //function : GetSubShapeOri
1667 //purpose : Return orientation of sub-shape in the main shape
1668 //=======================================================================
1670 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
1671 const TopoDS_Shape& subShape)
1673 TopAbs_Orientation ori = TopAbs_Orientation(-1);
1674 if ( !shape.IsNull() && !subShape.IsNull() )
1676 TopExp_Explorer e( shape, subShape.ShapeType() );
1677 if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
1678 e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
1679 for ( ; e.More(); e.Next())
1680 if ( subShape.IsSame( e.Current() ))
1683 ori = e.Current().Orientation();
1688 //=======================================================================
1689 //function : IsSubShape
1691 //=======================================================================
1693 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
1694 const TopoDS_Shape& mainShape )
1696 if ( !shape.IsNull() && !mainShape.IsNull() )
1698 for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
1701 if ( shape.IsSame( exp.Current() ))
1704 SCRUTE((shape.IsNull()));
1705 SCRUTE((mainShape.IsNull()));
1709 //=======================================================================
1710 //function : IsSubShape
1712 //=======================================================================
1714 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
1716 if ( shape.IsNull() || !aMesh )
1719 aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
1721 (shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape ));
1724 //================================================================================
1726 * \brief Return maximal tolerance of shape
1728 //================================================================================
1730 double SMESH_MesherHelper::MaxTolerance( const TopoDS_Shape& shape )
1732 double tol = Precision::Confusion();
1733 TopExp_Explorer exp;
1734 for ( exp.Init( shape, TopAbs_FACE ); exp.More(); exp.Next() )
1735 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Face( exp.Current())));
1736 for ( exp.Init( shape, TopAbs_EDGE ); exp.More(); exp.Next() )
1737 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Edge( exp.Current())));
1738 for ( exp.Init( shape, TopAbs_VERTEX ); exp.More(); exp.Next() )
1739 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Vertex( exp.Current())));
1744 //================================================================================
1746 * \brief Check if the first and last vertices of an edge are the same
1747 * \param anEdge - the edge to check
1748 * \retval bool - true if same
1750 //================================================================================
1752 bool SMESH_MesherHelper::IsClosedEdge( const TopoDS_Edge& anEdge )
1754 if ( anEdge.Orientation() >= TopAbs_INTERNAL )
1755 return IsClosedEdge( TopoDS::Edge( anEdge.Oriented( TopAbs_FORWARD )));
1756 return TopExp::FirstVertex( anEdge ).IsSame( TopExp::LastVertex( anEdge ));
1759 //=======================================================================
1760 //function : IsQuadraticMesh
1761 //purpose : Check mesh without geometry for: if all elements on this shape are quadratic,
1762 // quadratic elements will be created.
1763 // Used then generated 3D mesh without geometry.
1764 //=======================================================================
1766 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1768 int NbAllEdgsAndFaces=0;
1769 int NbQuadFacesAndEdgs=0;
1770 int NbFacesAndEdges=0;
1771 //All faces and edges
1772 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1774 //Quadratic faces and edges
1775 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1777 //Linear faces and edges
1778 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1780 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1782 return SMESH_MesherHelper::QUADRATIC;
1784 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1786 return SMESH_MesherHelper::LINEAR;
1789 //Mesh with both type of elements
1790 return SMESH_MesherHelper::COMP;
1793 //=======================================================================
1794 //function : GetOtherParam
1795 //purpose : Return an alternative parameter for a node on seam
1796 //=======================================================================
1798 double SMESH_MesherHelper::GetOtherParam(const double param) const
1800 int i = myParIndex & U_periodic ? 0 : 1;
1801 return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1804 //#include <Perf_Meter.hxx>
1806 //=======================================================================
1807 namespace { // Structures used by FixQuadraticElements()
1808 //=======================================================================
1810 #define __DMP__(txt) \
1812 #define MSG(txt) __DMP__(txt<<endl)
1813 #define MSGBEG(txt) __DMP__(txt)
1815 //const double straightTol2 = 1e-33; // to detect straing links
1816 bool isStraightLink(double linkLen2, double middleNodeMove2)
1818 // straight if <node move> < 1/15 * <link length>
1819 return middleNodeMove2 < 1/15./15. * linkLen2;
1823 // ---------------------------------------
1825 * \brief Quadratic link knowing its faces
1827 struct QLink: public SMESH_TLink
1829 const SMDS_MeshNode* _mediumNode;
1830 mutable vector<const QFace* > _faces;
1831 mutable gp_Vec _nodeMove;
1832 mutable int _nbMoves;
1834 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1835 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1837 //if ( MediumPos() != SMDS_TOP_3DSPACE )
1838 _nodeMove = MediumPnt() - MiddlePnt();
1840 void SetContinuesFaces() const;
1841 const QFace* GetContinuesFace( const QFace* face ) const;
1842 bool OnBoundary() const;
1843 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1844 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1846 SMDS_TypeOfPosition MediumPos() const
1847 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1848 SMDS_TypeOfPosition EndPos(bool isSecond) const
1849 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1850 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1851 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1853 void Move(const gp_Vec& move, bool sum=false) const
1854 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1855 gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1856 bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1857 bool IsStraight() const
1858 { return isStraightLink( (XYZ(node1())-XYZ(node2())).SquareModulus(),
1859 _nodeMove.SquareMagnitude());
1861 bool operator<(const QLink& other) const {
1862 return (node1()->GetID() == other.node1()->GetID() ?
1863 node2()->GetID() < other.node2()->GetID() :
1864 node1()->GetID() < other.node1()->GetID());
1866 // struct PtrComparator {
1867 // bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1870 // ---------------------------------------------------------
1872 * \brief Link in the chain of links; it connects two faces
1876 const QLink* _qlink;
1877 mutable const QFace* _qfaces[2];
1879 TChainLink(const QLink* qlink=0):_qlink(qlink) {
1880 _qfaces[0] = _qfaces[1] = 0;
1882 void SetFace(const QFace* face) const { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1884 bool IsBoundary() const { return !_qfaces[1]; }
1886 void RemoveFace( const QFace* face ) const
1887 { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1889 const QFace* NextFace( const QFace* f ) const
1890 { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1892 const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1893 { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1895 bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1897 operator bool() const { return (_qlink); }
1899 const QLink* operator->() const { return _qlink; }
1901 gp_Vec Normal() const;
1903 // --------------------------------------------------------------------
1904 typedef list< TChainLink > TChain;
1905 typedef set < TChainLink > TLinkSet;
1906 typedef TLinkSet::const_iterator TLinkInSet;
1908 const int theFirstStep = 5;
1910 enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1911 // --------------------------------------------------------------------
1913 * \brief Face shared by two volumes and bound by QLinks
1915 struct QFace: public TIDSortedNodeSet
1917 mutable const SMDS_MeshElement* _volumes[2];
1918 mutable vector< const QLink* > _sides;
1919 mutable bool _sideIsAdded[4]; // added in chain of links
1922 mutable const SMDS_MeshElement* _face;
1925 QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1927 void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1929 int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1931 void AddSelfToLinks() const {
1932 for ( int i = 0; i < _sides.size(); ++i )
1933 _sides[i]->_faces.push_back( this );
1935 int LinkIndex( const QLink* side ) const {
1936 for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1939 bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1941 bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1943 int i = LinkIndex( link._qlink );
1944 if ( i < 0 ) return true;
1945 _sideIsAdded[i] = true;
1946 link.SetFace( this );
1947 // continue from opposite link
1948 return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1950 bool IsBoundary() const { return !_volumes[1]; }
1952 bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1954 bool IsSpoiled(const QLink* bentLink ) const;
1956 TLinkInSet GetBoundaryLink( const TLinkSet& links,
1957 const TChainLink& avoidLink,
1958 TLinkInSet * notBoundaryLink = 0,
1959 const SMDS_MeshNode* nodeToContain = 0,
1960 bool * isAdjacentUsed = 0,
1961 int nbRecursionsLeft = -1) const;
1963 TLinkInSet GetLinkByNode( const TLinkSet& links,
1964 const TChainLink& avoidLink,
1965 const SMDS_MeshNode* nodeToContain) const;
1967 const SMDS_MeshNode* GetNodeInFace() const {
1968 for ( int iL = 0; iL < _sides.size(); ++iL )
1969 if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1973 gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1975 double MoveByBoundary( const TChainLink& theLink,
1976 const gp_Vec& theRefVec,
1977 const TLinkSet& theLinks,
1978 SMESH_MesherHelper* theFaceHelper=0,
1979 const double thePrevLen=0,
1980 const int theStep=theFirstStep,
1981 gp_Vec* theLinkNorm=0,
1982 double theSign=1.0) const;
1985 //================================================================================
1987 * \brief Dump QLink and QFace
1989 ostream& operator << (ostream& out, const QLink& l)
1991 out <<"QLink nodes: "
1992 << l.node1()->GetID() << " - "
1993 << l._mediumNode->GetID() << " - "
1994 << l.node2()->GetID() << endl;
1997 ostream& operator << (ostream& out, const QFace& f)
1999 out <<"QFace nodes: "/*<< &f << " "*/;
2000 for ( TIDSortedNodeSet::const_iterator n = f.begin(); n != f.end(); ++n )
2001 out << (*n)->GetID() << " ";
2002 out << " \tvolumes: "
2003 << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
2004 << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
2005 out << " \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
2009 //================================================================================
2011 * \brief Construct QFace from QLinks
2013 //================================================================================
2015 QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
2017 _volumes[0] = _volumes[1] = 0;
2019 _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
2020 _normal.SetCoord(0,0,0);
2021 for ( int i = 1; i < _sides.size(); ++i ) {
2022 const QLink *l1 = _sides[i-1], *l2 = _sides[i];
2023 insert( l1->node1() ); insert( l1->node2() );
2025 gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
2026 gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
2027 if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
2031 double normSqSize = _normal.SquareMagnitude();
2032 if ( normSqSize > numeric_limits<double>::min() )
2033 _normal /= sqrt( normSqSize );
2035 _normal.SetCoord(1e-33,0,0);
2041 //================================================================================
2043 * \brief Make up a chain of links
2044 * \param iSide - link to add first
2045 * \param chain - chain to fill in
2046 * \param pos - postion of medium nodes the links should have
2047 * \param error - out, specifies what is wrong
2048 * \retval bool - false if valid chain can't be built; "valid" means that links
2049 * of the chain belongs to rectangles bounding hexahedrons
2051 //================================================================================
2053 bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
2055 if ( iSide >= _sides.size() ) // wrong argument iSide
2057 if ( _sideIsAdded[ iSide ]) // already in chain
2060 if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
2063 list< const QFace* > faces( 1, this );
2064 while ( !faces.empty() ) {
2065 const QFace* face = faces.front();
2066 for ( int i = 0; i < face->_sides.size(); ++i ) {
2067 if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
2068 face->_sideIsAdded[i] = true;
2069 // find a face side in the chain
2070 TLinkInSet chLink = links.insert( TChainLink(face->_sides[i])).first;
2071 // TChain::iterator chLink = chain.begin();
2072 // for ( ; chLink != chain.end(); ++chLink )
2073 // if ( chLink->_qlink == face->_sides[i] )
2075 // if ( chLink == chain.end() )
2076 // chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
2077 // add a face to a chained link and put a continues face in the queue
2078 chLink->SetFace( face );
2079 if ( face->_sides[i]->MediumPos() >= pos )
2080 if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
2081 faces.push_back( contFace );
2086 if ( error < ERR_TRI )
2088 chain.insert( chain.end(), links.begin(),links.end() );
2091 _sideIsAdded[iSide] = true; // not to add this link to chain again
2092 const QLink* link = _sides[iSide];
2096 // add link into chain
2097 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
2098 chLink->SetFace( this );
2101 // propagate from quadrangle to neighbour faces
2102 if ( link->MediumPos() >= pos ) {
2103 int nbLinkFaces = link->_faces.size();
2104 if ( nbLinkFaces == 4 || (nbLinkFaces < 4 && link->OnBoundary())) {
2105 // hexahedral mesh or boundary quadrangles - goto a continous face
2106 if ( const QFace* f = link->GetContinuesFace( this ))
2107 return f->GetLinkChain( *chLink, chain, pos, error );
2110 TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
2111 for ( int i = 0; i < nbLinkFaces; ++i )
2112 if ( link->_faces[i] )
2113 link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
2114 if ( error < ERR_PRISM )
2122 //================================================================================
2124 * \brief Return a boundary link of the triangle face
2125 * \param links - set of all links
2126 * \param avoidLink - link not to return
2127 * \param notBoundaryLink - out, neither the returned link nor avoidLink
2128 * \param nodeToContain - node the returned link must contain; if provided, search
2129 * also performed on adjacent faces
2130 * \param isAdjacentUsed - returns true if link is found in adjacent faces
2131 * \param nbRecursionsLeft - to limit recursion
2133 //================================================================================
2135 TLinkInSet QFace::GetBoundaryLink( const TLinkSet& links,
2136 const TChainLink& avoidLink,
2137 TLinkInSet * notBoundaryLink,
2138 const SMDS_MeshNode* nodeToContain,
2139 bool * isAdjacentUsed,
2140 int nbRecursionsLeft) const
2142 TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
2144 typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
2145 TFaceLinkList adjacentFaces;
2147 for ( int iL = 0; iL < _sides.size(); ++iL )
2149 if ( avoidLink._qlink == _sides[iL] )
2151 TLinkInSet link = links.find( _sides[iL] );
2152 if ( link == linksEnd ) continue;
2153 if ( (*link)->MediumPos() > SMDS_TOP_FACE )
2154 continue; // We work on faces here, don't go inside a solid
2157 if ( link->IsBoundary() ) {
2158 if ( !nodeToContain ||
2159 (*link)->node1() == nodeToContain ||
2160 (*link)->node2() == nodeToContain )
2162 boundaryLink = link;
2163 if ( !notBoundaryLink ) break;
2166 else if ( notBoundaryLink ) {
2167 *notBoundaryLink = link;
2168 if ( boundaryLink != linksEnd ) break;
2171 if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
2172 if ( const QFace* adj = link->NextFace( this ))
2173 if ( adj->Contains( nodeToContain ))
2174 adjacentFaces.push_back( make_pair( adj, link ));
2177 if ( isAdjacentUsed ) *isAdjacentUsed = false;
2178 if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
2180 if ( nbRecursionsLeft < 0 )
2181 nbRecursionsLeft = nodeToContain->NbInverseElements();
2182 TFaceLinkList::iterator adj = adjacentFaces.begin();
2183 for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
2184 boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
2185 isAdjacentUsed, nbRecursionsLeft-1);
2186 if ( isAdjacentUsed ) *isAdjacentUsed = true;
2188 return boundaryLink;
2190 //================================================================================
2192 * \brief Return a link ending at the given node but not avoidLink
2194 //================================================================================
2196 TLinkInSet QFace::GetLinkByNode( const TLinkSet& links,
2197 const TChainLink& avoidLink,
2198 const SMDS_MeshNode* nodeToContain) const
2200 for ( int i = 0; i < _sides.size(); ++i )
2201 if ( avoidLink._qlink != _sides[i] &&
2202 (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
2203 return links.find( _sides[ i ]);
2207 //================================================================================
2209 * \brief Return normal to the i-th side pointing outside the face
2211 //================================================================================
2213 gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
2215 gp_Vec norm, vecOut;
2216 // if ( uvHelper ) {
2217 // TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
2218 // const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
2219 // gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
2220 // gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
2221 // norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
2223 // const QLink* otherLink = _sides[(i + 1) % _sides.size()];
2224 // const SMDS_MeshNode* otherNode =
2225 // otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
2226 // gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
2227 // vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
2230 norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
2231 gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
2232 XYZ( _sides[0]->node2() ) +
2233 XYZ( _sides[1]->node1() )) / 3.;
2234 vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
2236 if ( norm * vecOut < 0 )
2238 double mag2 = norm.SquareMagnitude();
2239 if ( mag2 > numeric_limits<double>::min() )
2240 norm /= sqrt( mag2 );
2243 //================================================================================
2245 * \brief Move medium node of theLink according to its distance from boundary
2246 * \param theLink - link to fix
2247 * \param theRefVec - movement of boundary
2248 * \param theLinks - all adjacent links of continous triangles
2249 * \param theFaceHelper - helper is not used so far
2250 * \param thePrevLen - distance from the boundary
2251 * \param theStep - number of steps till movement propagation limit
2252 * \param theLinkNorm - out normal to theLink
2253 * \param theSign - 1 or -1 depending on movement of boundary
2254 * \retval double - distance from boundary to propagation limit or other boundary
2256 //================================================================================
2258 double QFace::MoveByBoundary( const TChainLink& theLink,
2259 const gp_Vec& theRefVec,
2260 const TLinkSet& theLinks,
2261 SMESH_MesherHelper* theFaceHelper,
2262 const double thePrevLen,
2264 gp_Vec* theLinkNorm,
2265 double theSign) const
2268 return thePrevLen; // propagation limit reached
2270 int iL; // index of theLink
2271 for ( iL = 0; iL < _sides.size(); ++iL )
2272 if ( theLink._qlink == _sides[ iL ])
2275 MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
2276 <<" thePrevLen " << thePrevLen);
2277 MSG(string(theStep,'.')<<" "<<*theLink._qlink);
2279 gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
2280 double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
2281 if ( theStep == theFirstStep )
2282 theSign = refProj < 0. ? -1. : 1.;
2283 else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
2284 return thePrevLen; // to propagate movement forward only, not in side dir or backward
2286 int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
2287 TLinkInSet link1 = theLinks.find( _sides[iL1] );
2288 TLinkInSet link2 = theLinks.find( _sides[iL2] );
2289 if ( link1 == theLinks.end() || link2 == theLinks.end() )
2291 const QFace* f1 = link1->NextFace( this ); // adjacent faces
2292 const QFace* f2 = link2->NextFace( this );
2294 // propagate to adjacent faces till limit step or boundary
2295 double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
2296 double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
2297 gp_Vec linkDir1(0,0,0); // initialize to avoid valgrind error ("Conditional jump...")
2298 gp_Vec linkDir2(0,0,0);
2302 len1 = f1->MoveByBoundary
2303 ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
2305 linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
2307 MSG( " --------------- EXCEPTION");
2313 len2 = f2->MoveByBoundary
2314 ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
2316 linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
2318 MSG( " --------------- EXCEPTION");
2323 if ( theStep != theFirstStep )
2325 // choose chain length by direction of propagation most codirected with theRefVec
2326 bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
2327 fullLen = choose1 ? len1 : len2;
2328 double r = thePrevLen / fullLen;
2330 gp_Vec move = linkNorm * refProj * ( 1 - r );
2331 theLink->Move( move, true );
2333 MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
2334 " by " << refProj * ( 1 - r ) << " following " <<
2335 (choose1 ? *link1->_qlink : *link2->_qlink));
2337 if ( theLinkNorm ) *theLinkNorm = linkNorm;
2342 //================================================================================
2344 * \brief Checks if the face is distorted due to bentLink
2346 //================================================================================
2348 bool QFace::IsSpoiled(const QLink* bentLink ) const
2350 // code is valid for convex faces only
2352 for ( TIDSortedNodeSet::const_iterator n = begin(); n!=end(); ++n)
2353 gc += XYZ( *n ) / size();
2354 for (unsigned i = 0; i < _sides.size(); ++i )
2356 if ( _sides[i] == bentLink ) continue;
2357 gp_Vec linkNorm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
2358 gp_Vec vecOut( gc, _sides[i]->MiddlePnt() );
2359 if ( linkNorm * vecOut < 0 )
2361 double mag2 = linkNorm.SquareMagnitude();
2362 if ( mag2 > numeric_limits<double>::min() )
2363 linkNorm /= sqrt( mag2 );
2364 gp_Vec vecBent ( _sides[i]->MiddlePnt(), bentLink->MediumPnt());
2365 gp_Vec vecStraight( _sides[i]->MiddlePnt(), bentLink->MiddlePnt());
2366 if ( vecBent * linkNorm > -0.1*vecStraight.Magnitude() )
2373 //================================================================================
2375 * \brief Find pairs of continues faces
2377 //================================================================================
2379 void QLink::SetContinuesFaces() const
2381 // x0 x - QLink, [-|] - QFace, v - volume
2383 // | Between _faces of link x2 two vertical faces are continues
2384 // x1----x2-----x3 and two horizontal faces are continues. We set vertical faces
2385 // | to _faces[0] and _faces[1] and horizontal faces to
2386 // v2 | v3 _faces[2] and _faces[3] (or vise versa).
2389 if ( _faces.empty() )
2392 for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
2394 // look for a face bounding none of volumes bound by _faces[0]
2395 bool sameVol = false;
2396 int nbVol = _faces[iF]->NbVolumes();
2397 for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
2398 sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
2399 _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
2403 if ( iFaceCont > 0 ) // continues faces found, set one by the other
2405 if ( iFaceCont != 1 )
2406 std::swap( _faces[1], _faces[iFaceCont] );
2408 else if ( _faces.size() > 1 ) // not found, set NULL by the first face
2410 _faces.insert( ++_faces.begin(), 0 );
2413 //================================================================================
2415 * \brief Return a face continues to the given one
2417 //================================================================================
2419 const QFace* QLink::GetContinuesFace( const QFace* face ) const
2421 for ( int i = 0; i < _faces.size(); ++i ) {
2422 if ( _faces[i] == face ) {
2423 int iF = i < 2 ? 1-i : 5-i;
2424 return iF < _faces.size() ? _faces[iF] : 0;
2429 //================================================================================
2431 * \brief True if link is on mesh boundary
2433 //================================================================================
2435 bool QLink::OnBoundary() const
2437 for ( int i = 0; i < _faces.size(); ++i )
2438 if (_faces[i] && _faces[i]->IsBoundary()) return true;
2441 //================================================================================
2443 * \brief Return normal of link of the chain
2445 //================================================================================
2447 gp_Vec TChainLink::Normal() const {
2449 if (_qfaces[0]) norm = _qfaces[0]->_normal;
2450 if (_qfaces[1]) norm += _qfaces[1]->_normal;
2453 //================================================================================
2455 * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
2457 //================================================================================
2459 void fixPrism( TChain& allLinks )
2461 // separate boundary links from internal ones
2462 typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
2463 QLinkSet interLinks, bndLinks1, bndLink2;
2465 bool isCurved = false;
2466 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2467 if ( (*lnk)->OnBoundary() )
2468 bndLinks1.insert( lnk->_qlink );
2470 interLinks.insert( lnk->_qlink );
2471 isCurved = isCurved || !(*lnk)->IsStraight();
2474 return; // no need to move
2476 QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
2478 while ( !interLinks.empty() && !curBndLinks->empty() )
2480 // propagate movement from boundary links to connected internal links
2481 QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
2482 for ( ; bnd != bndEnd; ++bnd )
2484 const QLink* bndLink = *bnd;
2485 for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
2487 const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
2488 if ( !face ) continue;
2489 // find and move internal link opposite to bndLink within the face
2490 int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
2491 const QLink* interLink = face->_sides[ interInd ];
2492 QLinkSet::iterator pInterLink = interLinks.find( interLink );
2493 if ( pInterLink == interLinks.end() ) continue; // not internal link
2494 interLink->Move( bndLink->_nodeMove );
2495 // treated internal links become new boundary ones
2496 interLinks. erase( pInterLink );
2497 newBndLinks->insert( interLink );
2500 curBndLinks->clear();
2501 std::swap( curBndLinks, newBndLinks );
2505 //================================================================================
2507 * \brief Fix links of continues triangles near curved boundary
2509 //================================================================================
2511 void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2513 if ( allLinks.empty() ) return;
2515 TLinkSet linkSet( allLinks.begin(), allLinks.end());
2516 TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2518 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2520 if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2522 // move iff a boundary link is bent towards inside of a face (issue 0021084)
2523 const QFace* face = linkIt->_qfaces[0];
2524 gp_XYZ pIn = ( face->_sides[0]->MiddlePnt() +
2525 face->_sides[1]->MiddlePnt() +
2526 face->_sides[2]->MiddlePnt() ) / 3.;
2527 gp_XYZ insideDir( pIn - (*linkIt)->MiddlePnt());
2528 bool linkBentInside = ((*linkIt)->_nodeMove.Dot( insideDir ) > 0 );
2529 //if ( face->IsSpoiled( linkIt->_qlink ))
2530 if ( linkBentInside )
2531 face->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2536 //================================================================================
2538 * \brief Detect rectangular structure of links and build chains from them
2540 //================================================================================
2542 enum TSplitTriaResult {
2543 _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2544 _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK, _TWISTED_CHAIN };
2546 TSplitTriaResult splitTrianglesIntoChains( TChain & allLinks,
2547 vector< TChain> & resultChains,
2548 SMDS_TypeOfPosition pos )
2550 // put links in the set and evalute number of result chains by number of boundary links
2553 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2554 linkSet.insert( *lnk );
2555 nbBndLinks += lnk->IsBoundary();
2557 resultChains.clear();
2558 resultChains.reserve( nbBndLinks / 2 );
2560 TLinkInSet linkIt, linksEnd = linkSet.end();
2562 // find a boundary link with corner node; corner node has position pos-2
2563 // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2565 SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2566 const SMDS_MeshNode* corner = 0;
2567 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2568 if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2573 TLinkInSet startLink = linkIt;
2574 const SMDS_MeshNode* startCorner = corner;
2575 vector< TChain* > rowChains;
2578 while ( startLink != linksEnd) // loop on columns
2580 // We suppose we have a rectangular structure like shown here. We have found a
2581 // corner of the rectangle (startCorner) and a boundary link sharing
2582 // |/ |/ | the startCorner (startLink). We are going to loop on rows of the
2583 // --o---o---o structure making several chains at once. One chain (columnChain)
2584 // |\ | /| starts at startLink and continues upward (we look at the structure
2585 // \ | \ | / | from such point that startLink is on the bottom of the structure).
2586 // \| \|/ | While going upward we also fill horizontal chains (rowChains) we
2587 // --o---o---o encounter.
2589 // / | \ | \ | startCorner
2594 if ( resultChains.size() == nbBndLinks / 2 )
2596 resultChains.push_back( TChain() );
2597 TChain& columnChain = resultChains.back();
2599 TLinkInSet botLink = startLink; // current horizontal link to go up from
2600 corner = startCorner; // current corner the botLink ends at
2602 while ( botLink != linksEnd ) // loop on rows
2604 // add botLink to the columnChain
2605 columnChain.push_back( *botLink );
2607 const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2609 { // the column ends
2610 if ( botLink == startLink )
2611 return _TWISTED_CHAIN; // issue 0020951
2612 linkSet.erase( botLink );
2613 if ( iRow != rowChains.size() )
2614 return _FEW_ROWS; // different nb of rows in columns
2617 // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2618 // link ending at <corner> (sideLink); there are two cases:
2619 // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2620 // since midQuadLink is not at boundary while sideLink is.
2621 // 2) midQuadLink ends at <corner>
2623 TLinkInSet midQuadLink = linksEnd;
2624 TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2626 if ( isCase2 ) { // find midQuadLink among links of botTria
2627 midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2628 if ( midQuadLink->IsBoundary() )
2629 return _BAD_MIDQUAD;
2631 if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2632 return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2635 columnChain.push_back( *midQuadLink );
2636 if ( iRow >= rowChains.size() ) {
2638 return _MANY_ROWS; // different nb of rows in columns
2639 if ( resultChains.size() == nbBndLinks / 2 )
2641 resultChains.push_back( TChain() );
2642 rowChains.push_back( & resultChains.back() );
2644 rowChains[iRow]->push_back( *sideLink );
2645 rowChains[iRow]->push_back( *midQuadLink );
2647 const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2651 // prepare startCorner and startLink for the next column
2652 startCorner = startLink->NextNode( startCorner );
2654 startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2656 startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2657 // check if no more columns remains
2658 if ( startLink != linksEnd ) {
2659 const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2660 if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2661 startLink = linksEnd; // startLink bounds upTria or botTria
2662 else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2666 // find bottom link and corner for the next row
2667 corner = sideLink->NextNode( corner );
2668 // next bottom link ends at the new corner
2669 linkSet.erase( botLink );
2670 botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2671 if ( botLink == linksEnd || botLink == midQuadLink || botLink == sideLink)
2673 if ( midQuadLink == startLink || sideLink == startLink )
2674 return _TWISTED_CHAIN; // issue 0020951
2675 linkSet.erase( midQuadLink );
2676 linkSet.erase( sideLink );
2678 // make faces neighboring the found ones be boundary
2679 if ( startLink != linksEnd ) {
2680 const QFace* tria = isCase2 ? botTria : upTria;
2681 for ( int iL = 0; iL < 3; ++iL ) {
2682 linkIt = linkSet.find( tria->_sides[iL] );
2683 if ( linkIt != linksEnd )
2684 linkIt->RemoveFace( tria );
2687 if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2688 botLink->RemoveFace( upTria ); // make next botTria first in vector
2695 // In the linkSet, there must remain the last links of rowChains; add them
2696 if ( linkSet.size() != rowChains.size() )
2697 return _BAD_SET_SIZE;
2698 for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2699 // find the link (startLink) ending at startCorner
2701 for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2702 if ( (*startLink)->node1() == startCorner ) {
2703 corner = (*startLink)->node2(); break;
2705 else if ( (*startLink)->node2() == startCorner) {
2706 corner = (*startLink)->node1(); break;
2709 if ( startLink == linksEnd )
2711 rowChains[ iRow ]->push_back( *startLink );
2712 linkSet.erase( startLink );
2713 startCorner = corner;
2720 //=======================================================================
2722 * \brief Move medium nodes of faces and volumes to fix distorted elements
2723 * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2725 * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2727 //=======================================================================
2729 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2731 // 0. Apply algorithm to solids or geom faces
2732 // ----------------------------------------------
2733 if ( myShape.IsNull() ) {
2734 if ( !myMesh->HasShapeToMesh() ) return;
2735 SetSubShape( myMesh->GetShapeToMesh() );
2739 TopTools_IndexedMapOfShape solids;
2740 TopExp::MapShapes(myShape,TopAbs_SOLID,solids);
2741 nbSolids = solids.Extent();
2743 TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2744 for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2745 faces.Add( f.Current() ); // not in solid
2747 for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2748 if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2749 for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2750 faces.Add( f.Current() ); // in not meshed solid
2752 else { // fix nodes in the solid and its faces
2754 MSG("FIX SOLID " << nbSolids-- << " #" << GetMeshDS()->ShapeToIndex(s.Current()));
2756 SMESH_MesherHelper h(*myMesh);
2757 h.SetSubShape( s.Current() );
2758 h.FixQuadraticElements(false);
2761 // fix nodes on geom faces
2763 //int nbfaces = faces.Extent();
2765 for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2766 MSG("FIX FACE " << nbfaces-- << " #" << GetMeshDS()->ShapeToIndex(fIt.Key()));
2767 SMESH_MesherHelper h(*myMesh);
2768 h.SetSubShape( fIt.Key() );
2769 h.FixQuadraticElements(true);
2771 //perf_print_all_meters(1);
2775 // 1. Find out type of elements and get iterator on them
2776 // ---------------------------------------------------
2778 SMDS_ElemIteratorPtr elemIt;
2779 SMDSAbs_ElementType elemType = SMDSAbs_All;
2781 SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2784 if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2785 elemIt = smDS->GetElements();
2786 if ( elemIt->more() ) {
2787 elemType = elemIt->next()->GetType();
2788 elemIt = smDS->GetElements();
2791 if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2794 // 2. Fill in auxiliary data structures
2795 // ----------------------------------
2799 set< QLink >::iterator pLink;
2800 set< QFace >::iterator pFace;
2802 bool isCurved = false;
2803 //bool hasRectFaces = false;
2804 //set<int> nbElemNodeSet;
2806 if ( elemType == SMDSAbs_Volume )
2808 SMDS_VolumeTool volTool;
2809 while ( elemIt->more() ) // loop on volumes
2811 const SMDS_MeshElement* vol = elemIt->next();
2812 if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2814 for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2816 int nbN = volTool.NbFaceNodes( iF );
2817 //nbElemNodeSet.insert( nbN );
2818 const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2819 vector< const QLink* > faceLinks( nbN/2 );
2820 for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2823 QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2824 pLink = links.insert( link ).first;
2825 faceLinks[ iN/2 ] = & *pLink;
2827 isCurved = !link.IsStraight();
2828 if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2829 return; // already fixed
2832 pFace = faces.insert( QFace( faceLinks )).first;
2833 if ( pFace->NbVolumes() == 0 )
2834 pFace->AddSelfToLinks();
2835 pFace->SetVolume( vol );
2836 // hasRectFaces = hasRectFaces ||
2837 // ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2838 // volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2841 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2843 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2844 faceNodes[4],faceNodes[6] );
2848 set< QLink >::iterator pLink = links.begin();
2849 for ( ; pLink != links.end(); ++pLink )
2850 pLink->SetContinuesFaces();
2854 while ( elemIt->more() ) // loop on faces
2856 const SMDS_MeshElement* face = elemIt->next();
2857 if ( !face->IsQuadratic() )
2859 //nbElemNodeSet.insert( face->NbNodes() );
2860 int nbN = face->NbNodes()/2;
2861 vector< const QLink* > faceLinks( nbN );
2862 for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2865 QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2866 pLink = links.insert( link ).first;
2867 faceLinks[ iN ] = & *pLink;
2869 isCurved = !link.IsStraight();
2872 pFace = faces.insert( QFace( faceLinks )).first;
2873 pFace->AddSelfToLinks();
2874 //hasRectFaces = ( hasRectFaces || nbN == 4 );
2878 return; // no curved edges of faces
2880 // 3. Compute displacement of medium nodes
2881 // -------------------------------------
2883 // two loops on faces: the first is to treat boundary links, the second is for internal ones
2884 TopLoc_Location loc;
2885 // not treat boundary of volumic submesh
2886 int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2887 for ( ; isInside < 2; ++isInside ) {
2888 MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2889 SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2890 SMDS_TypeOfPosition bndPos = isInside ? SMDS_TOP_FACE : SMDS_TOP_EDGE;
2892 for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2893 if ( bool(isInside) == pFace->IsBoundary() )
2895 for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2898 // make chain of links connected via continues faces
2901 if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2903 if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2905 vector< TChain > chains;
2906 if ( error == ERR_OK ) { // chain contains continues rectangles
2908 chains[0].splice( chains[0].begin(), rawChain );
2910 else if ( error == ERR_TRI ) { // chain contains continues triangles
2911 TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2912 if ( res != _OK ) { // not quadrangles split into triangles
2913 fixTriaNearBoundary( rawChain, *this );
2917 else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2918 fixPrism( rawChain );
2924 for ( int iC = 0; iC < chains.size(); ++iC )
2926 TChain& chain = chains[iC];
2927 if ( chain.empty() ) continue;
2928 if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2929 MSG("3D straight - ignore");
2932 if ( chain.front()->MediumPos() > bndPos ||
2933 chain.back()->MediumPos() > bndPos ) {
2934 MSG("Internal chain - ignore");
2937 // mesure chain length and compute link position along the chain
2938 double chainLen = 0;
2939 vector< double > linkPos;
2940 MSGBEG( "Link medium nodes: ");
2941 TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2942 for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2943 MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2944 double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2945 while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2946 link1 = chain.erase( link1 );
2947 if ( link1 == chain.end() )
2949 len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2952 linkPos.push_back( chainLen );
2955 if ( linkPos.size() < 2 )
2958 gp_Vec move0 = chain.front()->_nodeMove;
2959 gp_Vec move1 = chain.back ()->_nodeMove;
2962 bool checkUV = true;
2964 // compute node displacement of end links in parametric space of face
2965 const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2966 TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2967 if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2969 face = TopoDS::Face( f );
2970 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2972 for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2974 TChainLink& link = is1 ? chain.back() : chain.front();
2975 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2976 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2977 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2978 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2979 // uvMove = uvm - uv12
2980 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2981 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2982 if ( !is1 ) // correct nodeOnFace for move1 (issue 0020919)
2983 nodeOnFace = (*(++chain.rbegin()))->_mediumNode;
2984 isStraight[is1] = isStraightLink( (uv2-uv1).SquareModulus(),uvMove.SquareModulus());
2986 // if ( move0.SquareMagnitude() < straightTol2 &&
2987 // move1.SquareMagnitude() < straightTol2 ) {
2988 if ( isStraight[0] && isStraight[1] ) {
2989 MSG("2D straight - ignore");
2990 continue; // straight - no need to move nodes of internal links
2995 if ( isInside || face.IsNull() )
2997 // compute node displacement of end links in their local coord systems
2999 TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
3000 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
3001 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
3002 move0.Transform(trsf);
3005 TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
3006 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
3007 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
3008 move1.Transform(trsf);
3011 // compute displacement of medium nodes
3012 link2 = chain.begin();
3015 for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
3017 double r = linkPos[i] / chainLen;
3018 // displacement in local coord system
3019 gp_Vec move = (1. - r) * move0 + r * move1;
3020 if ( isInside || face.IsNull()) {
3021 // transform to global
3022 gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
3023 gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
3024 gp_Vec x = x01.Normalized() + x12.Normalized();
3025 trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
3026 move.Transform(trsf);
3029 // compute 3D displacement by 2D one
3030 Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
3031 gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
3032 gp_XY newUV = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
3033 gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
3034 move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
3036 if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
3037 move.SquareMagnitude())
3039 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
3040 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
3041 MSG( "TOO LONG MOVE \t" <<
3042 "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
3043 "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
3044 "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
3045 "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
3049 (*link1)->Move( move );
3050 MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
3051 << chain.front()->_mediumNode->GetID() <<"-"
3052 << chain.back ()->_mediumNode->GetID() <<
3053 " by " << move.Magnitude());
3055 } // loop on chains of links
3056 } // loop on 2 directions of propagation from quadrangle
3063 for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
3064 if ( pLink->IsMoved() ) {
3065 //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
3066 gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
3067 GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
3072 //=======================================================================
3074 * \brief Iterator on ancestors of the given type
3076 //=======================================================================
3078 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
3080 TopTools_ListIteratorOfListOfShape _ancIter;
3081 TopAbs_ShapeEnum _type;
3082 TopTools_MapOfShape _encountered;
3083 TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
3084 : _ancIter( ancestors ), _type( type )
3086 if ( _ancIter.More() ) {
3087 if ( _ancIter.Value().ShapeType() != _type ) next();
3088 else _encountered.Add( _ancIter.Value() );
3093 return _ancIter.More();
3095 virtual const TopoDS_Shape* next()
3097 const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
3098 if ( _ancIter.More() )
3099 for ( _ancIter.Next(); _ancIter.More(); _ancIter.Next())
3100 if ( _ancIter.Value().ShapeType() == _type && _encountered.Add( _ancIter.Value() ))
3106 //=======================================================================
3108 * \brief Return iterator on ancestors of the given type
3110 //=======================================================================
3112 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
3113 const SMESH_Mesh& mesh,
3114 TopAbs_ShapeEnum ancestorType)
3116 return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));