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"
34 #include <BRepAdaptor_Surface.hxx>
35 #include <BRepTools.hxx>
36 #include <BRepTools_WireExplorer.hxx>
37 #include <BRep_Tool.hxx>
38 #include <Geom2d_Curve.hxx>
39 #include <GeomAPI_ProjectPointOnCurve.hxx>
40 #include <GeomAPI_ProjectPointOnSurf.hxx>
41 #include <Geom_Curve.hxx>
42 #include <Geom_Surface.hxx>
43 #include <ShapeAnalysis.hxx>
45 #include <TopExp_Explorer.hxx>
46 #include <TopTools_ListIteratorOfListOfShape.hxx>
47 #include <TopTools_MapIteratorOfMapOfShape.hxx>
48 #include <TopTools_MapOfShape.hxx>
51 #include <gp_Pnt2d.hxx>
52 #include <gp_Trsf.hxx>
54 #include <Standard_Failure.hxx>
55 #include <Standard_ErrorHandler.hxx>
57 #include <utilities.h>
63 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
67 gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
69 enum { U_periodic = 1, V_periodic = 2 };
72 //================================================================================
76 //================================================================================
78 SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
79 : myParIndex(0), myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false)
81 myPar1[0] = myPar2[0] = myPar1[1] = myPar2[1] = 0;
82 mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
85 //=======================================================================
86 //function : ~SMESH_MesherHelper
88 //=======================================================================
90 SMESH_MesherHelper::~SMESH_MesherHelper()
92 TID2Projector::iterator i_proj = myFace2Projector.begin();
93 for ( ; i_proj != myFace2Projector.end(); ++i_proj )
94 delete i_proj->second;
97 //=======================================================================
98 //function : IsQuadraticSubMesh
99 //purpose : Check submesh for given shape: if all elements on this shape
100 // are quadratic, quadratic elements will be created.
101 // Also fill myTLinkNodeMap
102 //=======================================================================
104 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
106 SMESHDS_Mesh* meshDS = GetMeshDS();
107 // we can create quadratic elements only if all elements
108 // created on subshapes of given shape are quadratic
109 // also we have to fill myTLinkNodeMap
110 myCreateQuadratic = true;
111 mySeamShapeIds.clear();
112 myDegenShapeIds.clear();
113 TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
114 SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
116 int nbOldLinks = myTLinkNodeMap.size();
118 TopExp_Explorer exp( aSh, subType );
119 for (; exp.More() && myCreateQuadratic; exp.Next()) {
120 if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
121 if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
123 const SMDS_MeshElement* e = it->next();
124 if ( e->GetType() != elemType || !e->IsQuadratic() ) {
125 myCreateQuadratic = false;
130 switch ( e->NbNodes() ) {
132 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
134 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
135 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
136 AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
138 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
139 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
140 AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
141 AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
144 myCreateQuadratic = false;
153 if ( nbOldLinks == myTLinkNodeMap.size() )
154 myCreateQuadratic = false;
156 if(!myCreateQuadratic) {
157 myTLinkNodeMap.clear();
161 return myCreateQuadratic;
164 //=======================================================================
165 //function : SetSubShape
166 //purpose : Set geomerty to make elements on
167 //=======================================================================
169 void SMESH_MesherHelper::SetSubShape(const int aShID)
171 if ( aShID == myShapeID )
174 SetSubShape( GetMeshDS()->IndexToShape( aShID ));
176 SetSubShape( TopoDS_Shape() );
179 //=======================================================================
180 //function : SetSubShape
181 //purpose : Set geomerty to create elements on
182 //=======================================================================
184 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
186 if ( myShape.IsSame( aSh ))
190 mySeamShapeIds.clear();
191 myDegenShapeIds.clear();
193 if ( myShape.IsNull() ) {
197 SMESHDS_Mesh* meshDS = GetMeshDS();
198 myShapeID = meshDS->ShapeToIndex(aSh);
201 // treatment of periodic faces
202 for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
204 const TopoDS_Face& face = TopoDS::Face( eF.Current() );
205 BRepAdaptor_Surface surface( face );
206 if ( surface.IsUPeriodic() || surface.IsVPeriodic() )
208 for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
210 // look for a seam edge
211 const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
212 if ( BRep_Tool::IsClosed( edge, face )) {
213 // initialize myPar1, myPar2 and myParIndex
215 BRep_Tool::UVPoints( edge, face, uv1, uv2 );
216 if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
218 myParIndex |= U_periodic;
219 myPar1[0] = surface.FirstUParameter();
220 myPar2[0] = surface.LastUParameter();
223 myParIndex |= V_periodic;
224 myPar1[1] = surface.FirstVParameter();
225 myPar2[1] = surface.LastVParameter();
227 // store seam shape indices, negative if shape encounters twice
228 int edgeID = meshDS->ShapeToIndex( edge );
229 mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
230 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
231 int vertexID = meshDS->ShapeToIndex( v.Current() );
232 mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
236 // look for a degenerated edge
237 if ( BRep_Tool::Degenerated( edge )) {
238 myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
239 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
240 myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
247 //=======================================================================
248 //function : GetNodeUVneedInFaceNode
249 //purpose : Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
250 // Return true if the face is periodic.
251 // If F is Null, answer about subshape set through IsQuadraticSubMesh() or
253 //=======================================================================
255 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
257 if ( F.IsNull() ) return !mySeamShapeIds.empty();
259 if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
260 return !mySeamShapeIds.empty();
263 Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
264 if ( !aSurface.IsNull() )
265 return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
270 //=======================================================================
271 //function : IsMedium
273 //=======================================================================
275 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode* node,
276 const SMDSAbs_ElementType typeToCheck)
278 return SMESH_MeshEditor::IsMedium( node, typeToCheck );
281 //=======================================================================
282 //function : GetSubShapeByNode
283 //purpose : Return support shape of a node
284 //=======================================================================
286 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
287 SMESHDS_Mesh* meshDS)
289 int shapeID = node->GetPosition()->GetShapeId();
290 if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
291 return meshDS->IndexToShape( shapeID );
293 return TopoDS_Shape();
297 //=======================================================================
298 //function : AddTLinkNode
299 //purpose : add a link in my data structure
300 //=======================================================================
302 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
303 const SMDS_MeshNode* n2,
304 const SMDS_MeshNode* n12)
306 // add new record to map
307 SMESH_TLink link( n1, n2 );
308 myTLinkNodeMap.insert( make_pair(link,n12));
311 //=======================================================================
312 //function : GetUVOnSeam
313 //purpose : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
314 //=======================================================================
316 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
318 gp_Pnt2d result = uv1;
319 for ( int i = U_periodic; i <= V_periodic ; ++i )
321 if ( myParIndex & i )
323 double p1 = uv1.Coord( i );
324 double dp1 = Abs( p1-myPar1[i-1]), dp2 = Abs( p1-myPar2[i-1]);
325 if ( myParIndex == i ||
326 dp1 < ( myPar2[i-1] - myPar2[i-1] ) / 100. ||
327 dp2 < ( myPar2[i-1] - myPar2[i-1] ) / 100. )
329 double p2 = uv2.Coord( i );
330 double p1Alt = ( dp1 < dp2 ) ? myPar2[i-1] : myPar1[i-1];
331 if ( Abs( p2 - p1 ) > Abs( p2 - p1Alt ))
332 result.SetCoord( i, p1Alt );
339 //=======================================================================
340 //function : GetNodeUV
341 //purpose : Return node UV on face
342 //=======================================================================
344 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
345 const SMDS_MeshNode* n,
346 const SMDS_MeshNode* n2,
349 gp_Pnt2d uv( Precision::Infinite(), Precision::Infinite() );
350 const SMDS_PositionPtr Pos = n->GetPosition();
352 if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
354 // node has position on face
355 const SMDS_FacePosition* fpos =
356 static_cast<const SMDS_FacePosition*>(n->GetPosition().get());
357 uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
359 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( F ));
361 else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
363 // node has position on edge => it is needed to find
364 // corresponding edge from face, get pcurve for this
365 // edge and retrieve value from this pcurve
366 const SMDS_EdgePosition* epos =
367 static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
368 int edgeID = Pos->GetShapeId();
369 TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
370 double f, l, u = epos->GetUParameter();
371 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
372 bool validU = ( f < u && u < l );
374 uv = C2d->Value( u );
377 if ( check || !validU )
378 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( E ),/*force=*/ !validU );
380 // for a node on a seam edge select one of UVs on 2 pcurves
381 if ( n2 && IsSeamShape( edgeID ) )
383 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0, check ));
386 { // adjust uv to period
388 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
389 Standard_Boolean isUPeriodic = S->IsUPeriodic();
390 Standard_Boolean isVPeriodic = S->IsVPeriodic();
391 if ( isUPeriodic || isVPeriodic ) {
392 Standard_Real UF,UL,VF,VL;
393 S->Bounds(UF,UL,VF,VL);
395 uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
397 uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
401 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
403 if ( int vertexID = n->GetPosition()->GetShapeId() ) {
404 const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
406 uv = BRep_Tool::Parameters( V, F );
409 catch (Standard_Failure& exc) {
412 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() )
413 uvOK = ( V == vert.Current() );
416 MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
417 << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
419 // get UV of a vertex closest to the node
421 gp_Pnt pn = XYZ( n );
422 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) {
423 TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
424 gp_Pnt p = BRep_Tool::Pnt( curV );
425 double curDist = p.SquareDistance( pn );
426 if ( curDist < dist ) {
428 uv = BRep_Tool::Parameters( curV, F );
429 uvOK = ( dist < DBL_MIN );
435 TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
436 for ( ; it.More(); it.Next() ) {
437 if ( it.Value().ShapeType() == TopAbs_EDGE ) {
438 const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
440 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
441 if ( !C2d.IsNull() ) {
442 double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
443 uv = C2d->Value( u );
451 if ( n2 && IsSeamShape( vertexID ) )
452 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
462 //=======================================================================
463 //function : CheckNodeUV
464 //purpose : Check and fix node UV on a face
465 //=======================================================================
467 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
468 const SMDS_MeshNode* n,
471 const bool force) const
473 if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
475 // check that uv is correct
477 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
478 gp_Pnt nodePnt = XYZ( n );
479 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
480 if ( Precision::IsInfinite( uv.X() ) ||
481 Precision::IsInfinite( uv.Y() ) ||
482 nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > tol )
484 // uv incorrect, project the node to surface
485 GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
486 projector.Perform( nodePnt );
487 if ( !projector.IsDone() || projector.NbPoints() < 1 )
489 MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
492 Quantity_Parameter U,V;
493 projector.LowerDistanceParameters(U,V);
495 if ( nodePnt.Distance( surface->Value( U, V )) > tol )
497 MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
500 //uv.SetCoord( U,V );
502 else if ( uv.Modulus() > numeric_limits<double>::min() )
504 ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
510 //=======================================================================
511 //function : GetProjector
512 //purpose : Return projector intitialized by given face without location, which is returned
513 //=======================================================================
515 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
516 TopLoc_Location& loc,
519 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
520 int faceID = GetMeshDS()->ShapeToIndex( F );
521 TID2Projector& i2proj = const_cast< TID2Projector&>( myFace2Projector );
522 TID2Projector::iterator i_proj = i2proj.find( faceID );
523 if ( i_proj == i2proj.end() )
525 if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
526 double U1, U2, V1, V2;
527 surface->Bounds(U1, U2, V1, V2);
528 GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
529 proj->Init( surface, U1, U2, V1, V2, tol );
530 i_proj = i2proj.insert( make_pair( faceID, proj )).first;
532 return *( i_proj->second );
537 gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
538 gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
539 gp_XY_FunPtr(Subtracted);
542 //=======================================================================
543 //function : applyIn2D
544 //purpose : Perform given operation on two 2d points in parameric space of given surface.
545 // It takes into account period of the surface. Use gp_XY_FunPtr macro
546 // to easily define pointer to function of gp_XY class.
547 //=======================================================================
549 gp_XY SMESH_MesherHelper::applyIn2D(const Handle(Geom_Surface)& surface,
553 const bool resultInPeriod)
555 Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
556 Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
557 if ( !isUPeriodic && !isVPeriodic )
560 // move uv2 not far than half-period from uv1
562 uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
564 uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
567 gp_XY res = fun( uv1, gp_XY(u2,v2) );
569 // move result within period
570 if ( resultInPeriod )
572 Standard_Real UF,UL,VF,VL;
573 surface->Bounds(UF,UL,VF,VL);
575 res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
577 res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
582 //=======================================================================
583 //function : GetMiddleUV
584 //purpose : Return middle UV taking in account surface period
585 //=======================================================================
587 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
591 return applyIn2D( surface, p1, p2, & AverageUV );
594 //=======================================================================
595 //function : GetNodeU
596 //purpose : Return node U on edge
597 //=======================================================================
599 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
600 const SMDS_MeshNode* n,
601 const SMDS_MeshNode* inEdgeNode,
605 const SMDS_PositionPtr pos = n->GetPosition();
606 if ( pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
608 const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( pos.get() );
609 param = epos->GetUParameter();
611 else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
613 if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
616 BRep_Tool::Range( E, f,l );
617 double uInEdge = GetNodeU( E, inEdgeNode );
618 param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
622 SMESHDS_Mesh * meshDS = GetMeshDS();
623 int vertexID = pos->GetShapeId();
624 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
625 param = BRep_Tool::Parameter( V, E );
630 double tol = BRep_Tool::Tolerance( E );
631 double f,l; BRep_Tool::Range( E, f,l );
632 bool force = ( param < f-tol || param > l+tol );
633 if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
634 force = ( GetMeshDS()->ShapeToIndex( E ) != pos->GetShapeId() );
636 *check = CheckNodeU( E, n, param, tol, force );
641 //=======================================================================
642 //function : CheckNodeU
643 //purpose : Check and fix node U on an edge
644 // Return false if U is bad and could not be fixed
645 //=======================================================================
647 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
648 const SMDS_MeshNode* n,
651 const bool force) const
653 if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
655 // check that u is correct
656 TopLoc_Location loc; double f,l;
657 Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
658 if ( curve.IsNull() ) // degenerated edge
660 if ( u+tol < f || u-tol > l )
662 double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
668 gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
669 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
670 if ( nodePnt.Distance( curve->Value( u )) > tol )
672 // u incorrect, project the node to the curve
673 GeomAPI_ProjectPointOnCurve projector( nodePnt, curve, f, l );
674 if ( projector.NbPoints() < 1 )
676 MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
679 Quantity_Parameter U = projector.LowerDistanceParameter();
681 if ( nodePnt.Distance( curve->Value( U )) > tol )
683 MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
688 else if ( fabs( u ) > numeric_limits<double>::min() )
690 ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
692 if (( u < f-tol || u > l+tol ) && force )
694 // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
697 // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
698 double period = curve->Period();
699 u = ( u < f ) ? u + period : u - period;
701 catch (Standard_Failure& exc)
711 //=======================================================================
712 //function : GetMediumNode
713 //purpose : Return existing or create new medium nodes between given ones
714 //=======================================================================
716 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
717 const SMDS_MeshNode* n2,
720 // Find existing node
722 SMESH_TLink link(n1,n2);
723 ItTLinkNode itLN = myTLinkNodeMap.find( link );
724 if ( itLN != myTLinkNodeMap.end() ) {
725 return (*itLN).second;
728 // Create medium node
731 SMESHDS_Mesh* meshDS = GetMeshDS();
733 // get type of shape for the new medium node
734 int faceID = -1, edgeID = -1;
735 const SMDS_PositionPtr Pos1 = n1->GetPosition();
736 const SMDS_PositionPtr Pos2 = n2->GetPosition();
738 if( myShape.IsNull() )
740 if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
741 faceID = Pos1->GetShapeId();
743 else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
744 faceID = Pos2->GetShapeId();
747 if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
748 edgeID = Pos1->GetShapeId();
750 if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
751 edgeID = Pos2->GetShapeId();
754 // get positions of the given nodes on shapes
755 TopoDS_Edge E; double u [2];
756 TopoDS_Face F; gp_XY uv[2];
757 bool uvOK[2] = { false, false };
758 TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
759 if ( faceID>0 || shapeType == TopAbs_FACE)
761 if( myShape.IsNull() )
762 F = TopoDS::Face(meshDS->IndexToShape(faceID));
764 F = TopoDS::Face(myShape);
767 uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
768 uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
770 else if (edgeID>0 || shapeType == TopAbs_EDGE)
772 if( myShape.IsNull() )
773 E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
775 E = TopoDS::Edge(myShape);
778 u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
779 u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
783 // we try to create medium node using UV parameters of
784 // nodes, else - medium between corresponding 3d points
787 if ( uvOK[0] && uvOK[1] )
789 if ( IsDegenShape( Pos1->GetShapeId() ))
790 if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
791 else uv[0].SetCoord( 2, uv[1].Coord( 2 ));
792 else if ( IsDegenShape( Pos2->GetShapeId() ))
793 if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
794 else uv[1].SetCoord( 2, uv[0].Coord( 2 ));
797 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
798 gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
799 gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
800 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
801 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
802 myTLinkNodeMap.insert(make_pair(link,n12));
806 else if ( !E.IsNull() )
809 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
812 Standard_Boolean isPeriodic = C->IsPeriodic();
815 Standard_Real Period = C->Period();
816 Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
817 Standard_Real pmid = (u[0]+p)/2.;
818 U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
823 gp_Pnt P = C->Value( U );
824 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
825 meshDS->SetNodeOnEdge(n12, edgeID, U);
826 myTLinkNodeMap.insert(make_pair(link,n12));
832 double x = ( n1->X() + n2->X() )/2.;
833 double y = ( n1->Y() + n2->Y() )/2.;
834 double z = ( n1->Z() + n2->Z() )/2.;
835 n12 = meshDS->AddNode(x,y,z);
838 gp_XY UV = ( uv[0] + uv[1] ) / 2.;
839 CheckNodeUV( F, n12, UV, BRep_Tool::Tolerance( F ), /*force=*/true);
840 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
842 else if ( !E.IsNull() )
844 double U = ( u[0] + u[1] ) / 2.;
845 CheckNodeU( E, n12, U, BRep_Tool::Tolerance( E ), /*force=*/true);
846 meshDS->SetNodeOnEdge(n12, edgeID, U);
848 else if ( myShapeID > 0 )
850 meshDS->SetNodeInVolume(n12, myShapeID);
852 myTLinkNodeMap.insert( make_pair( link, n12 ));
856 //=======================================================================
858 //purpose : Creates a node
859 //=======================================================================
861 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
863 SMESHDS_Mesh * meshDS = GetMeshDS();
864 SMDS_MeshNode* node = 0;
866 node = meshDS->AddNodeWithID( x, y, z, ID );
868 node = meshDS->AddNode( x, y, z );
869 if ( mySetElemOnShape && myShapeID > 0 ) {
870 switch ( myShape.ShapeType() ) {
871 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
872 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
873 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID); break;
874 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID); break;
875 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
882 //=======================================================================
884 //purpose : Creates quadratic or linear edge
885 //=======================================================================
887 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
888 const SMDS_MeshNode* n2,
892 SMESHDS_Mesh * meshDS = GetMeshDS();
894 SMDS_MeshEdge* edge = 0;
895 if (myCreateQuadratic) {
896 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
898 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
900 edge = meshDS->AddEdge(n1, n2, n12);
904 edge = meshDS->AddEdgeWithID(n1, n2, id);
906 edge = meshDS->AddEdge(n1, n2);
909 if ( mySetElemOnShape && myShapeID > 0 )
910 meshDS->SetMeshElementOnShape( edge, myShapeID );
915 //=======================================================================
917 //purpose : Creates quadratic or linear triangle
918 //=======================================================================
920 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
921 const SMDS_MeshNode* n2,
922 const SMDS_MeshNode* n3,
926 SMESHDS_Mesh * meshDS = GetMeshDS();
927 SMDS_MeshFace* elem = 0;
929 if( n1==n2 || n2==n3 || n3==n1 )
932 if(!myCreateQuadratic) {
934 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
936 elem = meshDS->AddFace(n1, n2, n3);
939 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
940 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
941 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
944 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
946 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
948 if ( mySetElemOnShape && myShapeID > 0 )
949 meshDS->SetMeshElementOnShape( elem, myShapeID );
954 //=======================================================================
956 //purpose : Creates quadratic or linear quadrangle
957 //=======================================================================
959 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
960 const SMDS_MeshNode* n2,
961 const SMDS_MeshNode* n3,
962 const SMDS_MeshNode* n4,
966 SMESHDS_Mesh * meshDS = GetMeshDS();
967 SMDS_MeshFace* elem = 0;
970 return AddFace(n1,n3,n4,id,force3d);
973 return AddFace(n1,n2,n4,id,force3d);
976 return AddFace(n1,n2,n3,id,force3d);
979 return AddFace(n1,n2,n4,id,force3d);
982 return AddFace(n1,n2,n3,id,force3d);
985 return AddFace(n1,n2,n3,id,force3d);
988 if(!myCreateQuadratic) {
990 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
992 elem = meshDS->AddFace(n1, n2, n3, n4);
995 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
996 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
997 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
998 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1001 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
1003 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
1005 if ( mySetElemOnShape && myShapeID > 0 )
1006 meshDS->SetMeshElementOnShape( elem, myShapeID );
1011 //=======================================================================
1012 //function : AddVolume
1013 //purpose : Creates quadratic or linear prism
1014 //=======================================================================
1016 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1017 const SMDS_MeshNode* n2,
1018 const SMDS_MeshNode* n3,
1019 const SMDS_MeshNode* n4,
1020 const SMDS_MeshNode* n5,
1021 const SMDS_MeshNode* n6,
1025 SMESHDS_Mesh * meshDS = GetMeshDS();
1026 SMDS_MeshVolume* elem = 0;
1027 if(!myCreateQuadratic) {
1029 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
1031 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
1034 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1035 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1036 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1038 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1039 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1040 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1042 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1043 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1044 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1047 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
1048 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1050 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1051 n12, n23, n31, n45, n56, n64, n14, n25, n36);
1053 if ( mySetElemOnShape && myShapeID > 0 )
1054 meshDS->SetMeshElementOnShape( elem, myShapeID );
1059 //=======================================================================
1060 //function : AddVolume
1061 //purpose : Creates quadratic or linear tetrahedron
1062 //=======================================================================
1064 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1065 const SMDS_MeshNode* n2,
1066 const SMDS_MeshNode* n3,
1067 const SMDS_MeshNode* n4,
1071 SMESHDS_Mesh * meshDS = GetMeshDS();
1072 SMDS_MeshVolume* elem = 0;
1073 if(!myCreateQuadratic) {
1075 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1077 elem = meshDS->AddVolume(n1, n2, n3, n4);
1080 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1081 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1082 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1084 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1085 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1086 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1089 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1091 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1093 if ( mySetElemOnShape && myShapeID > 0 )
1094 meshDS->SetMeshElementOnShape( elem, myShapeID );
1099 //=======================================================================
1100 //function : AddVolume
1101 //purpose : Creates quadratic or linear pyramid
1102 //=======================================================================
1104 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1105 const SMDS_MeshNode* n2,
1106 const SMDS_MeshNode* n3,
1107 const SMDS_MeshNode* n4,
1108 const SMDS_MeshNode* n5,
1112 SMDS_MeshVolume* elem = 0;
1113 if(!myCreateQuadratic) {
1115 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1117 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1120 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1121 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1122 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1123 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1125 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1126 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1127 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1128 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1131 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
1136 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
1138 n15, n25, n35, n45);
1140 if ( mySetElemOnShape && myShapeID > 0 )
1141 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1146 //=======================================================================
1147 //function : AddVolume
1148 //purpose : Creates quadratic or linear hexahedron
1149 //=======================================================================
1151 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1152 const SMDS_MeshNode* n2,
1153 const SMDS_MeshNode* n3,
1154 const SMDS_MeshNode* n4,
1155 const SMDS_MeshNode* n5,
1156 const SMDS_MeshNode* n6,
1157 const SMDS_MeshNode* n7,
1158 const SMDS_MeshNode* n8,
1162 SMESHDS_Mesh * meshDS = GetMeshDS();
1163 SMDS_MeshVolume* elem = 0;
1164 if(!myCreateQuadratic) {
1166 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1168 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1171 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1172 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1173 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1174 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1176 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1177 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1178 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1179 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1181 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1182 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1183 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1184 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1187 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1188 n12, n23, n34, n41, n56, n67,
1189 n78, n85, n15, n26, n37, n48, id);
1191 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1192 n12, n23, n34, n41, n56, n67,
1193 n78, n85, n15, n26, n37, n48);
1195 if ( mySetElemOnShape && myShapeID > 0 )
1196 meshDS->SetMeshElementOnShape( elem, myShapeID );
1201 //=======================================================================
1202 //function : LoadNodeColumns
1203 //purpose : Load nodes bound to face into a map of node columns
1204 //=======================================================================
1206 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1207 const TopoDS_Face& theFace,
1208 const TopoDS_Edge& theBaseEdge,
1209 SMESHDS_Mesh* theMesh)
1211 SMESHDS_SubMesh* faceSubMesh = theMesh->MeshElements( theFace );
1212 if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
1215 // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
1217 map< double, const SMDS_MeshNode*> sortedBaseNodes;
1218 if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
1219 || sortedBaseNodes.size() < 2 )
1222 int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
1223 map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
1224 double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
1225 for ( ; u_n != sortedBaseNodes.end(); u_n++ )
1227 double par = ( u_n->first - f ) / range;
1228 vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
1229 nCol.resize( nbRows );
1230 nCol[0] = u_n->second;
1233 // fill theParam2ColumnMap column by column by passing from nodes on
1234 // theBaseEdge up via mesh faces on theFace
1236 TParam2ColumnMap::iterator par_nVec_2 = theParam2ColumnMap.begin();
1237 TParam2ColumnMap::iterator par_nVec_1 = par_nVec_2++;
1238 TIDSortedElemSet emptySet, avoidSet;
1239 for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
1241 vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
1242 vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
1244 int i1, i2, iRow = 0;
1245 const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
1246 // find face sharing node n1 and n2 and belonging to faceSubMesh
1247 while ( const SMDS_MeshElement* face =
1248 SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
1250 if ( faceSubMesh->Contains( face ))
1252 int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
1255 n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
1256 n2 = face->GetNode( (i1+2) % 4 );
1257 if ( ++iRow >= nbRows )
1263 avoidSet.insert( face );
1265 if ( iRow + 1 < nbRows ) // compact if necessary
1266 nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
1271 //=======================================================================
1272 //function : NbAncestors
1273 //purpose : Return number of unique ancestors of the shape
1274 //=======================================================================
1276 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1277 const SMESH_Mesh& mesh,
1278 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
1280 TopTools_MapOfShape ancestors;
1281 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1282 for ( ; ansIt.More(); ansIt.Next() ) {
1283 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1284 ancestors.Add( ansIt.Value() );
1286 return ancestors.Extent();
1289 //=======================================================================
1290 //function : GetSubShapeOri
1291 //purpose : Return orientation of sub-shape in the main shape
1292 //=======================================================================
1294 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
1295 const TopoDS_Shape& subShape)
1297 TopAbs_Orientation ori = TopAbs_Orientation(-1);
1298 if ( !shape.IsNull() && !subShape.IsNull() )
1300 TopExp_Explorer e( shape, subShape.ShapeType() );
1301 if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
1302 e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
1303 for ( ; e.More(); e.Next())
1304 if ( subShape.IsSame( e.Current() ))
1307 ori = e.Current().Orientation();
1312 //=======================================================================
1313 //function : IsSubShape
1315 //=======================================================================
1317 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
1318 const TopoDS_Shape& mainShape )
1320 if ( !shape.IsNull() && !mainShape.IsNull() )
1322 for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
1325 if ( shape.IsSame( exp.Current() ))
1328 SCRUTE((shape.IsNull()));
1329 SCRUTE((mainShape.IsNull()));
1333 //=======================================================================
1334 //function : IsSubShape
1336 //=======================================================================
1338 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
1340 if ( shape.IsNull() || !aMesh )
1343 aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
1345 shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape );
1348 //=======================================================================
1349 //function : IsQuadraticMesh
1350 //purpose : Check mesh without geometry for: if all elements on this shape are quadratic,
1351 // quadratic elements will be created.
1352 // Used then generated 3D mesh without geometry.
1353 //=======================================================================
1355 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1357 int NbAllEdgsAndFaces=0;
1358 int NbQuadFacesAndEdgs=0;
1359 int NbFacesAndEdges=0;
1360 //All faces and edges
1361 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1363 //Quadratic faces and edges
1364 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1366 //Linear faces and edges
1367 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1369 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1371 return SMESH_MesherHelper::QUADRATIC;
1373 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1375 return SMESH_MesherHelper::LINEAR;
1378 //Mesh with both type of elements
1379 return SMESH_MesherHelper::COMP;
1382 //=======================================================================
1383 //function : GetOtherParam
1384 //purpose : Return an alternative parameter for a node on seam
1385 //=======================================================================
1387 double SMESH_MesherHelper::GetOtherParam(const double param) const
1389 int i = myParIndex & U_periodic ? 0 : 1;
1390 return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1393 //=======================================================================
1394 namespace { // Structures used by FixQuadraticElements()
1395 //=======================================================================
1397 #define __DMP__(txt) \
1399 #define MSG(txt) __DMP__(txt<<endl)
1400 #define MSGBEG(txt) __DMP__(txt)
1402 const double straightTol2 = 1e-33; // to detect straing links
1405 // ---------------------------------------
1407 * \brief Quadratic link knowing its faces
1409 struct QLink: public SMESH_TLink
1411 const SMDS_MeshNode* _mediumNode;
1412 mutable vector<const QFace* > _faces;
1413 mutable gp_Vec _nodeMove;
1414 mutable int _nbMoves;
1416 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1417 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1419 //if ( MediumPos() != SMDS_TOP_3DSPACE )
1420 _nodeMove = MediumPnt() - MiddlePnt();
1422 void SetContinuesFaces() const;
1423 const QFace* GetContinuesFace( const QFace* face ) const;
1424 bool OnBoundary() const;
1425 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1426 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1428 SMDS_TypeOfPosition MediumPos() const
1429 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1430 SMDS_TypeOfPosition EndPos(bool isSecond) const
1431 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1432 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1433 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1435 void Move(const gp_Vec& move, bool sum=false) const
1436 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1437 gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1438 bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1439 bool IsStraight() const { return _nodeMove.SquareMagnitude() <= straightTol2; }
1441 bool operator<(const QLink& other) const {
1442 return (node1()->GetID() == other.node1()->GetID() ?
1443 node2()->GetID() < other.node2()->GetID() :
1444 node1()->GetID() < other.node1()->GetID());
1446 struct PtrComparator {
1447 bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1450 // ---------------------------------------------------------
1452 * \brief Link in the chain of links; it connects two faces
1456 const QLink* _qlink;
1457 mutable const QFace* _qfaces[2];
1459 TChainLink(const QLink* qlink=0):_qlink(qlink) {
1460 _qfaces[0] = _qfaces[1] = 0;
1462 void SetFace(const QFace* face) { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1464 bool IsBoundary() const { return !_qfaces[1]; }
1466 void RemoveFace( const QFace* face ) const
1467 { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1469 const QFace* NextFace( const QFace* f ) const
1470 { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1472 const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1473 { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1475 bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1477 operator bool() const { return (_qlink); }
1479 const QLink* operator->() const { return _qlink; }
1481 gp_Vec Normal() const;
1483 // --------------------------------------------------------------------
1484 typedef list< TChainLink > TChain;
1485 typedef set < TChainLink > TLinkSet;
1486 typedef TLinkSet::const_iterator TLinkInSet;
1488 const int theFirstStep = 5;
1490 enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1491 // --------------------------------------------------------------------
1493 * \brief Face shared by two volumes and bound by QLinks
1495 struct QFace: public TIDSortedElemSet
1497 mutable const SMDS_MeshElement* _volumes[2];
1498 mutable vector< const QLink* > _sides;
1499 mutable bool _sideIsAdded[4]; // added in chain of links
1502 mutable const SMDS_MeshElement* _face;
1505 QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1507 void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1509 int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1511 void AddSelfToLinks() const {
1512 for ( int i = 0; i < _sides.size(); ++i )
1513 _sides[i]->_faces.push_back( this );
1515 int LinkIndex( const QLink* side ) const {
1516 for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1519 bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1521 bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1523 int i = LinkIndex( link._qlink );
1524 if ( i < 0 ) return true;
1525 _sideIsAdded[i] = true;
1526 link.SetFace( this );
1527 // continue from opposite link
1528 return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1530 bool IsBoundary() const { return !_volumes[1]; }
1532 bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1534 TLinkInSet GetBoundaryLink( const TLinkSet& links,
1535 const TChainLink& avoidLink,
1536 TLinkInSet * notBoundaryLink = 0,
1537 const SMDS_MeshNode* nodeToContain = 0,
1538 bool * isAdjacentUsed = 0,
1539 int nbRecursionsLeft = -1) const;
1541 TLinkInSet GetLinkByNode( const TLinkSet& links,
1542 const TChainLink& avoidLink,
1543 const SMDS_MeshNode* nodeToContain) const;
1545 const SMDS_MeshNode* GetNodeInFace() const {
1546 for ( int iL = 0; iL < _sides.size(); ++iL )
1547 if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1551 gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1553 double MoveByBoundary( const TChainLink& theLink,
1554 const gp_Vec& theRefVec,
1555 const TLinkSet& theLinks,
1556 SMESH_MesherHelper* theFaceHelper=0,
1557 const double thePrevLen=0,
1558 const int theStep=theFirstStep,
1559 gp_Vec* theLinkNorm=0,
1560 double theSign=1.0) const;
1563 //================================================================================
1565 * \brief Dump QLink and QFace
1567 ostream& operator << (ostream& out, const QLink& l)
1569 out <<"QLink nodes: "
1570 << l.node1()->GetID() << " - "
1571 << l._mediumNode->GetID() << " - "
1572 << l.node2()->GetID() << endl;
1575 ostream& operator << (ostream& out, const QFace& f)
1577 out <<"QFace nodes: "/*<< &f << " "*/;
1578 for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
1579 out << (*n)->GetID() << " ";
1580 out << " \tvolumes: "
1581 << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1582 << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1583 out << " \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1587 //================================================================================
1589 * \brief Construct QFace from QLinks
1591 //================================================================================
1593 QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
1595 _volumes[0] = _volumes[1] = 0;
1597 _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1598 _normal.SetCoord(0,0,0);
1599 for ( int i = 1; i < _sides.size(); ++i ) {
1600 const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1601 insert( l1->node1() ); insert( l1->node2() );
1603 gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1604 gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1605 if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1609 double normSqSize = _normal.SquareMagnitude();
1610 if ( normSqSize > numeric_limits<double>::min() )
1611 _normal /= sqrt( normSqSize );
1613 _normal.SetCoord(1e-33,0,0);
1619 //================================================================================
1621 * \brief Make up a chain of links
1622 * \param iSide - link to add first
1623 * \param chain - chain to fill in
1624 * \param pos - postion of medium nodes the links should have
1625 * \param error - out, specifies what is wrong
1626 * \retval bool - false if valid chain can't be built; "valid" means that links
1627 * of the chain belongs to rectangles bounding hexahedrons
1629 //================================================================================
1631 bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1633 if ( iSide >= _sides.size() ) // wrong argument iSide
1635 if ( _sideIsAdded[ iSide ]) // already in chain
1638 if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1640 list< const QFace* > faces( 1, this );
1641 while ( !faces.empty() ) {
1642 const QFace* face = faces.front();
1643 for ( int i = 0; i < face->_sides.size(); ++i ) {
1644 if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
1645 face->_sideIsAdded[i] = true;
1646 // find a face side in the chain
1647 TChain::iterator chLink = chain.begin();
1648 for ( ; chLink != chain.end(); ++chLink )
1649 if ( chLink->_qlink == face->_sides[i] )
1651 if ( chLink == chain.end() )
1652 chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
1653 // add a face to a chained link and put a continues face in the queue
1654 chLink->SetFace( face );
1655 if ( face->_sides[i]->MediumPos() >= pos )
1656 if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
1657 faces.push_back( contFace );
1662 if ( error < ERR_TRI )
1666 _sideIsAdded[iSide] = true; // not to add this link to chain again
1667 const QLink* link = _sides[iSide];
1671 // add link into chain
1672 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1673 chLink->SetFace( this );
1676 // propagate from quadrangle to neighbour faces
1677 if ( link->MediumPos() >= pos ) {
1678 int nbLinkFaces = link->_faces.size();
1679 if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1680 // hexahedral mesh or boundary quadrangles - goto a continous face
1681 if ( const QFace* f = link->GetContinuesFace( this ))
1682 return f->GetLinkChain( *chLink, chain, pos, error );
1685 TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1686 for ( int i = 0; i < nbLinkFaces; ++i )
1687 if ( link->_faces[i] )
1688 link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1689 if ( error < ERR_PRISM )
1697 //================================================================================
1699 * \brief Return a boundary link of the triangle face
1700 * \param links - set of all links
1701 * \param avoidLink - link not to return
1702 * \param notBoundaryLink - out, neither the returned link nor avoidLink
1703 * \param nodeToContain - node the returned link must contain; if provided, search
1704 * also performed on adjacent faces
1705 * \param isAdjacentUsed - returns true if link is found in adjacent faces
1706 * \param nbRecursionsLeft - to limit recursion
1708 //================================================================================
1710 TLinkInSet QFace::GetBoundaryLink( const TLinkSet& links,
1711 const TChainLink& avoidLink,
1712 TLinkInSet * notBoundaryLink,
1713 const SMDS_MeshNode* nodeToContain,
1714 bool * isAdjacentUsed,
1715 int nbRecursionsLeft) const
1717 TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1719 typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1720 TFaceLinkList adjacentFaces;
1722 for ( int iL = 0; iL < _sides.size(); ++iL )
1724 if ( avoidLink._qlink == _sides[iL] )
1726 TLinkInSet link = links.find( _sides[iL] );
1727 if ( link == linksEnd ) continue;
1728 if ( (*link)->MediumPos() > SMDS_TOP_FACE )
1729 continue; // We work on faces here, don't go inside a solid
1732 if ( link->IsBoundary() ) {
1733 if ( !nodeToContain ||
1734 (*link)->node1() == nodeToContain ||
1735 (*link)->node2() == nodeToContain )
1737 boundaryLink = link;
1738 if ( !notBoundaryLink ) break;
1741 else if ( notBoundaryLink ) {
1742 *notBoundaryLink = link;
1743 if ( boundaryLink != linksEnd ) break;
1746 if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
1747 if ( const QFace* adj = link->NextFace( this ))
1748 if ( adj->Contains( nodeToContain ))
1749 adjacentFaces.push_back( make_pair( adj, link ));
1752 if ( isAdjacentUsed ) *isAdjacentUsed = false;
1753 if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
1755 if ( nbRecursionsLeft < 0 )
1756 nbRecursionsLeft = nodeToContain->NbInverseElements();
1757 TFaceLinkList::iterator adj = adjacentFaces.begin();
1758 for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1759 boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
1760 isAdjacentUsed, nbRecursionsLeft-1);
1761 if ( isAdjacentUsed ) *isAdjacentUsed = true;
1763 return boundaryLink;
1765 //================================================================================
1767 * \brief Return a link ending at the given node but not avoidLink
1769 //================================================================================
1771 TLinkInSet QFace::GetLinkByNode( const TLinkSet& links,
1772 const TChainLink& avoidLink,
1773 const SMDS_MeshNode* nodeToContain) const
1775 for ( int i = 0; i < _sides.size(); ++i )
1776 if ( avoidLink._qlink != _sides[i] &&
1777 (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1778 return links.find( _sides[ i ]);
1782 //================================================================================
1784 * \brief Return normal to the i-th side pointing outside the face
1786 //================================================================================
1788 gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1790 gp_Vec norm, vecOut;
1791 // if ( uvHelper ) {
1792 // TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
1793 // const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
1794 // gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
1795 // gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
1796 // norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
1798 // const QLink* otherLink = _sides[(i + 1) % _sides.size()];
1799 // const SMDS_MeshNode* otherNode =
1800 // otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
1801 // gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
1802 // vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
1805 norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
1806 gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
1807 XYZ( _sides[0]->node2() ) +
1808 XYZ( _sides[1]->node1() )) / 3.;
1809 vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
1811 if ( norm * vecOut < 0 )
1813 double mag2 = norm.SquareMagnitude();
1814 if ( mag2 > numeric_limits<double>::min() )
1815 norm /= sqrt( mag2 );
1818 //================================================================================
1820 * \brief Move medium node of theLink according to its distance from boundary
1821 * \param theLink - link to fix
1822 * \param theRefVec - movement of boundary
1823 * \param theLinks - all adjacent links of continous triangles
1824 * \param theFaceHelper - helper is not used so far
1825 * \param thePrevLen - distance from the boundary
1826 * \param theStep - number of steps till movement propagation limit
1827 * \param theLinkNorm - out normal to theLink
1828 * \param theSign - 1 or -1 depending on movement of boundary
1829 * \retval double - distance from boundary to propagation limit or other boundary
1831 //================================================================================
1833 double QFace::MoveByBoundary( const TChainLink& theLink,
1834 const gp_Vec& theRefVec,
1835 const TLinkSet& theLinks,
1836 SMESH_MesherHelper* theFaceHelper,
1837 const double thePrevLen,
1839 gp_Vec* theLinkNorm,
1840 double theSign) const
1843 return thePrevLen; // propagation limit reached
1845 int iL; // index of theLink
1846 for ( iL = 0; iL < _sides.size(); ++iL )
1847 if ( theLink._qlink == _sides[ iL ])
1850 MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
1851 <<" thePrevLen " << thePrevLen);
1852 MSG(string(theStep,'.')<<" "<<*theLink._qlink);
1854 gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
1855 double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
1856 if ( theStep == theFirstStep )
1857 theSign = refProj < 0. ? -1. : 1.;
1858 else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
1859 return thePrevLen; // to propagate movement forward only, not in side dir or backward
1861 int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
1862 TLinkInSet link1 = theLinks.find( _sides[iL1] );
1863 TLinkInSet link2 = theLinks.find( _sides[iL2] );
1864 if ( link1 == theLinks.end() || link2 == theLinks.end() )
1866 const QFace* f1 = link1->NextFace( this ); // adjacent faces
1867 const QFace* f2 = link2->NextFace( this );
1869 // propagate to adjacent faces till limit step or boundary
1870 double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
1871 double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
1872 gp_Vec linkDir1(0,0,0); // initialize to avoid valgrind error ("Conditional jump...")
1873 gp_Vec linkDir2(0,0,0);
1877 len1 = f1->MoveByBoundary
1878 ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
1880 linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
1882 MSG( " --------------- EXCEPTION");
1888 len2 = f2->MoveByBoundary
1889 ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
1891 linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
1893 MSG( " --------------- EXCEPTION");
1898 if ( theStep != theFirstStep )
1900 // choose chain length by direction of propagation most codirected with theRefVec
1901 bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
1902 fullLen = choose1 ? len1 : len2;
1903 double r = thePrevLen / fullLen;
1905 gp_Vec move = linkNorm * refProj * ( 1 - r );
1906 theLink->Move( move, true );
1908 MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
1909 " by " << refProj * ( 1 - r ) << " following " <<
1910 (choose1 ? *link1->_qlink : *link2->_qlink));
1912 if ( theLinkNorm ) *theLinkNorm = linkNorm;
1917 //================================================================================
1919 * \brief Find pairs of continues faces
1921 //================================================================================
1923 void QLink::SetContinuesFaces() const
1925 // x0 x - QLink, [-|] - QFace, v - volume
1927 // | Between _faces of link x2 two vertical faces are continues
1928 // x1----x2-----x3 and two horizontal faces are continues. We set vertical faces
1929 // | to _faces[0] and _faces[1] and horizontal faces to
1930 // v2 | v3 _faces[2] and _faces[3] (or vise versa).
1933 if ( _faces.empty() )
1936 for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
1938 // look for a face bounding none of volumes bound by _faces[0]
1939 bool sameVol = false;
1940 int nbVol = _faces[iF]->NbVolumes();
1941 for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
1942 sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
1943 _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
1947 if ( iFaceCont > 0 ) // continues faces found, set one by the other
1949 if ( iFaceCont != 1 )
1950 std::swap( _faces[1], _faces[iFaceCont] );
1952 else if ( _faces.size() > 1 ) // not found, set NULL by the first face
1954 _faces.insert( ++_faces.begin(), 0 );
1957 //================================================================================
1959 * \brief Return a face continues to the given one
1961 //================================================================================
1963 const QFace* QLink::GetContinuesFace( const QFace* face ) const
1965 for ( int i = 0; i < _faces.size(); ++i ) {
1966 if ( _faces[i] == face ) {
1967 int iF = i < 2 ? 1-i : 5-i;
1968 return iF < _faces.size() ? _faces[iF] : 0;
1973 //================================================================================
1975 * \brief True if link is on mesh boundary
1977 //================================================================================
1979 bool QLink::OnBoundary() const
1981 for ( int i = 0; i < _faces.size(); ++i )
1982 if (_faces[i] && _faces[i]->IsBoundary()) return true;
1985 //================================================================================
1987 * \brief Return normal of link of the chain
1989 //================================================================================
1991 gp_Vec TChainLink::Normal() const {
1993 if (_qfaces[0]) norm = _qfaces[0]->_normal;
1994 if (_qfaces[1]) norm += _qfaces[1]->_normal;
1997 //================================================================================
1999 * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
2001 //================================================================================
2003 void fixPrism( TChain& allLinks )
2005 // separate boundary links from internal ones
2006 typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
2007 QLinkSet interLinks, bndLinks1, bndLink2;
2009 bool isCurved = false;
2010 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2011 if ( (*lnk)->OnBoundary() )
2012 bndLinks1.insert( lnk->_qlink );
2014 interLinks.insert( lnk->_qlink );
2015 isCurved = isCurved || !(*lnk)->IsStraight();
2018 return; // no need to move
2020 QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
2022 while ( !interLinks.empty() && !curBndLinks->empty() )
2024 // propagate movement from boundary links to connected internal links
2025 QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
2026 for ( ; bnd != bndEnd; ++bnd )
2028 const QLink* bndLink = *bnd;
2029 for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
2031 const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
2032 if ( !face ) continue;
2033 // find and move internal link opposite to bndLink within the face
2034 int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
2035 const QLink* interLink = face->_sides[ interInd ];
2036 QLinkSet::iterator pInterLink = interLinks.find( interLink );
2037 if ( pInterLink == interLinks.end() ) continue; // not internal link
2038 interLink->Move( bndLink->_nodeMove );
2039 // treated internal links become new boundary ones
2040 interLinks. erase( pInterLink );
2041 newBndLinks->insert( interLink );
2044 curBndLinks->clear();
2045 std::swap( curBndLinks, newBndLinks );
2049 //================================================================================
2051 * \brief Fix links of continues triangles near curved boundary
2053 //================================================================================
2055 void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2057 if ( allLinks.empty() ) return;
2059 TLinkSet linkSet( allLinks.begin(), allLinks.end());
2060 TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2062 // move in 2d if we are on geom face
2063 // TopoDS_Face face;
2064 // TopLoc_Location loc;
2065 // SMESH_MesherHelper faceHelper( *helper.GetMesh());
2066 // while ( linkIt->IsBoundary()) ++linkIt;
2067 // if ( linkIt == linksEnd ) return;
2068 // if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
2069 // bool checkPos = true;
2070 // TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
2071 // if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2072 // face = TopoDS::Face( f );
2073 // helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
2077 // faceHelper.SetSubShape( face );
2080 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2082 if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2084 // if ( !face.IsNull() ) {
2085 // const SMDS_MeshNode* inFaceNode =
2086 // faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
2087 // gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
2088 // gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
2089 // gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
2090 // gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2091 // gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
2092 // linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
2095 linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2101 //================================================================================
2103 * \brief Detect rectangular structure of links and build chains from them
2105 //================================================================================
2107 enum TSplitTriaResult {
2108 _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2109 _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK, _TWISTED_CHAIN };
2111 TSplitTriaResult splitTrianglesIntoChains( TChain & allLinks,
2112 vector< TChain> & resultChains,
2113 SMDS_TypeOfPosition pos )
2115 // put links in the set and evalute number of result chains by number of boundary links
2118 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2119 linkSet.insert( *lnk );
2120 nbBndLinks += lnk->IsBoundary();
2122 resultChains.clear();
2123 resultChains.reserve( nbBndLinks / 2 );
2125 TLinkInSet linkIt, linksEnd = linkSet.end();
2127 // find a boundary link with corner node; corner node has position pos-2
2128 // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2130 SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2131 const SMDS_MeshNode* corner = 0;
2132 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2133 if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2138 TLinkInSet startLink = linkIt;
2139 const SMDS_MeshNode* startCorner = corner;
2140 vector< TChain* > rowChains;
2143 while ( startLink != linksEnd) // loop on columns
2145 // We suppose we have a rectangular structure like shown here. We have found a
2146 // corner of the rectangle (startCorner) and a boundary link sharing
2147 // |/ |/ | the startCorner (startLink). We are going to loop on rows of the
2148 // --o---o---o structure making several chains at once. One chain (columnChain)
2149 // |\ | /| starts at startLink and continues upward (we look at the structure
2150 // \ | \ | / | from such point that startLink is on the bottom of the structure).
2151 // \| \|/ | While going upward we also fill horizontal chains (rowChains) we
2152 // --o---o---o encounter.
2154 // / | \ | \ | startCorner
2159 if ( resultChains.size() == nbBndLinks / 2 )
2161 resultChains.push_back( TChain() );
2162 TChain& columnChain = resultChains.back();
2164 TLinkInSet botLink = startLink; // current horizontal link to go up from
2165 corner = startCorner; // current corner the botLink ends at
2167 while ( botLink != linksEnd ) // loop on rows
2169 // add botLink to the columnChain
2170 columnChain.push_back( *botLink );
2172 const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2174 { // the column ends
2175 if ( botLink == startLink )
2176 return _TWISTED_CHAIN; // issue 0020951
2177 linkSet.erase( botLink );
2178 if ( iRow != rowChains.size() )
2179 return _FEW_ROWS; // different nb of rows in columns
2182 // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2183 // link ending at <corner> (sideLink); there are two cases:
2184 // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2185 // since midQuadLink is not at boundary while sideLink is.
2186 // 2) midQuadLink ends at <corner>
2188 TLinkInSet midQuadLink = linksEnd;
2189 TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2191 if ( isCase2 ) { // find midQuadLink among links of botTria
2192 midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2193 if ( midQuadLink->IsBoundary() )
2194 return _BAD_MIDQUAD;
2196 if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2197 return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2200 columnChain.push_back( *midQuadLink );
2201 if ( iRow >= rowChains.size() ) {
2203 return _MANY_ROWS; // different nb of rows in columns
2204 if ( resultChains.size() == nbBndLinks / 2 )
2206 resultChains.push_back( TChain() );
2207 rowChains.push_back( & resultChains.back() );
2209 rowChains[iRow]->push_back( *sideLink );
2210 rowChains[iRow]->push_back( *midQuadLink );
2212 const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2216 // prepare startCorner and startLink for the next column
2217 startCorner = startLink->NextNode( startCorner );
2219 startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2221 startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2222 // check if no more columns remains
2223 if ( startLink != linksEnd ) {
2224 const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2225 if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2226 startLink = linksEnd; // startLink bounds upTria or botTria
2227 else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2231 // find bottom link and corner for the next row
2232 corner = sideLink->NextNode( corner );
2233 // next bottom link ends at the new corner
2234 linkSet.erase( botLink );
2235 botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2236 if ( botLink == linksEnd || botLink == midQuadLink || botLink == sideLink)
2238 if ( midQuadLink == startLink || sideLink == startLink )
2239 return _TWISTED_CHAIN; // issue 0020951
2240 linkSet.erase( midQuadLink );
2241 linkSet.erase( sideLink );
2243 // make faces neighboring the found ones be boundary
2244 if ( startLink != linksEnd ) {
2245 const QFace* tria = isCase2 ? botTria : upTria;
2246 for ( int iL = 0; iL < 3; ++iL ) {
2247 linkIt = linkSet.find( tria->_sides[iL] );
2248 if ( linkIt != linksEnd )
2249 linkIt->RemoveFace( tria );
2252 if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2253 botLink->RemoveFace( upTria ); // make next botTria first in vector
2260 // In the linkSet, there must remain the last links of rowChains; add them
2261 if ( linkSet.size() != rowChains.size() )
2262 return _BAD_SET_SIZE;
2263 for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2264 // find the link (startLink) ending at startCorner
2266 for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2267 if ( (*startLink)->node1() == startCorner ) {
2268 corner = (*startLink)->node2(); break;
2270 else if ( (*startLink)->node2() == startCorner) {
2271 corner = (*startLink)->node1(); break;
2274 if ( startLink == linksEnd )
2276 rowChains[ iRow ]->push_back( *startLink );
2277 linkSet.erase( startLink );
2278 startCorner = corner;
2285 //=======================================================================
2287 * \brief Move medium nodes of faces and volumes to fix distorted elements
2288 * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2290 * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2292 //=======================================================================
2294 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2296 // apply algorithm to solids or geom faces
2297 // ----------------------------------------------
2298 if ( myShape.IsNull() ) {
2299 if ( !myMesh->HasShapeToMesh() ) return;
2300 SetSubShape( myMesh->GetShapeToMesh() );
2302 TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2303 for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2304 faces.Add( f.Current() );
2306 for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2307 if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2308 for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2309 faces.Add( f.Current() );
2311 else { // fix nodes in the solid and its faces
2312 SMESH_MesherHelper h(*myMesh);
2313 h.SetSubShape( s.Current() );
2314 h.FixQuadraticElements(false);
2317 // fix nodes on geom faces
2318 for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2319 SMESH_MesherHelper h(*myMesh);
2320 h.SetSubShape( fIt.Key() );
2321 h.FixQuadraticElements(true);
2326 // Find out type of elements and get iterator on them
2327 // ---------------------------------------------------
2329 SMDS_ElemIteratorPtr elemIt;
2330 SMDSAbs_ElementType elemType = SMDSAbs_All;
2332 SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2335 if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2336 elemIt = smDS->GetElements();
2337 if ( elemIt->more() ) {
2338 elemType = elemIt->next()->GetType();
2339 elemIt = smDS->GetElements();
2342 if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2345 // Fill in auxiliary data structures
2346 // ----------------------------------
2350 set< QLink >::iterator pLink;
2351 set< QFace >::iterator pFace;
2353 bool isCurved = false;
2354 bool hasRectFaces = false;
2355 set<int> nbElemNodeSet;
2357 if ( elemType == SMDSAbs_Volume )
2359 SMDS_VolumeTool volTool;
2360 while ( elemIt->more() ) // loop on volumes
2362 const SMDS_MeshElement* vol = elemIt->next();
2363 if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2365 for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2367 int nbN = volTool.NbFaceNodes( iF );
2368 nbElemNodeSet.insert( nbN );
2369 const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2370 vector< const QLink* > faceLinks( nbN/2 );
2371 for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2374 QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2375 pLink = links.insert( link ).first;
2376 faceLinks[ iN/2 ] = & *pLink;
2378 isCurved = !link.IsStraight();
2379 if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2380 return; // already fixed
2383 pFace = faces.insert( QFace( faceLinks )).first;
2384 if ( pFace->NbVolumes() == 0 )
2385 pFace->AddSelfToLinks();
2386 pFace->SetVolume( vol );
2387 hasRectFaces = hasRectFaces ||
2388 ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2389 volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2392 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2394 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2395 faceNodes[4],faceNodes[6] );
2399 set< QLink >::iterator pLink = links.begin();
2400 for ( ; pLink != links.end(); ++pLink )
2401 pLink->SetContinuesFaces();
2405 while ( elemIt->more() ) // loop on faces
2407 const SMDS_MeshElement* face = elemIt->next();
2408 if ( !face->IsQuadratic() )
2410 nbElemNodeSet.insert( face->NbNodes() );
2411 int nbN = face->NbNodes()/2;
2412 vector< const QLink* > faceLinks( nbN );
2413 for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2416 QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2417 pLink = links.insert( link ).first;
2418 faceLinks[ iN ] = & *pLink;
2420 isCurved = !link.IsStraight();
2423 pFace = faces.insert( QFace( faceLinks )).first;
2424 pFace->AddSelfToLinks();
2425 hasRectFaces = ( hasRectFaces || nbN == 4 );
2429 return; // no curved edges of faces
2431 // Compute displacement of medium nodes
2432 // -------------------------------------
2434 // two loops on faces: the first is to treat boundary links, the second is for internal ones
2435 TopLoc_Location loc;
2436 // not treat boundary of volumic submesh
2437 int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2438 for ( ; isInside < 2; ++isInside ) {
2439 MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2440 SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2441 SMDS_TypeOfPosition bndPos = isInside ? SMDS_TOP_FACE : SMDS_TOP_EDGE;
2443 for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2444 if ( bool(isInside) == pFace->IsBoundary() )
2446 for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2449 // make chain of links connected via continues faces
2452 if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2454 if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2456 vector< TChain > chains;
2457 if ( error == ERR_OK ) { // chain contains continues rectangles
2459 chains[0].splice( chains[0].begin(), rawChain );
2461 else if ( error == ERR_TRI ) { // chain contains continues triangles
2462 TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2463 if ( res != _OK ) { // not quadrangles split into triangles
2464 fixTriaNearBoundary( rawChain, *this );
2468 else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2469 fixPrism( rawChain );
2475 for ( int iC = 0; iC < chains.size(); ++iC )
2477 TChain& chain = chains[iC];
2478 if ( chain.empty() ) continue;
2479 if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2480 MSG("3D straight - ignore");
2483 if ( chain.front()->MediumPos() > bndPos ||
2484 chain.back()->MediumPos() > bndPos ) {
2485 MSG("Internal chain - ignore");
2488 // mesure chain length and compute link position along the chain
2489 double chainLen = 0;
2490 vector< double > linkPos;
2491 MSGBEG( "Link medium nodes: ");
2492 TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2493 for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2494 MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2495 double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2496 while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2497 link1 = chain.erase( link1 );
2498 if ( link1 == chain.end() )
2500 len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2503 linkPos.push_back( chainLen );
2506 if ( linkPos.size() < 2 )
2509 gp_Vec move0 = chain.front()->_nodeMove;
2510 gp_Vec move1 = chain.back ()->_nodeMove;
2513 bool checkUV = true;
2515 // compute node displacement of end links in parametric space of face
2516 const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2517 TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2518 if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2520 face = TopoDS::Face( f );
2521 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2522 for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2524 TChainLink& link = is1 ? chain.back() : chain.front();
2525 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2526 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2527 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2528 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2529 // uvMove = uvm - uv12
2530 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2531 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2532 if ( !is1 ) // correct nodeOnFace for move1 (issue 0020919)
2533 nodeOnFace = (*(++chain.rbegin()))->_mediumNode;
2535 if ( move0.SquareMagnitude() < straightTol2 &&
2536 move1.SquareMagnitude() < straightTol2 ) {
2537 MSG("2D straight - ignore");
2538 continue; // straight - no need to move nodes of internal links
2543 if ( isInside || face.IsNull() )
2545 // compute node displacement of end links in their local coord systems
2547 TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2548 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2549 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2550 move0.Transform(trsf);
2553 TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2554 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2555 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2556 move1.Transform(trsf);
2559 // compute displacement of medium nodes
2560 link2 = chain.begin();
2563 for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2565 double r = linkPos[i] / chainLen;
2566 // displacement in local coord system
2567 gp_Vec move = (1. - r) * move0 + r * move1;
2568 if ( isInside || face.IsNull()) {
2569 // transform to global
2570 gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2571 gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2572 gp_Vec x = x01.Normalized() + x12.Normalized();
2573 trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2574 move.Transform(trsf);
2577 // compute 3D displacement by 2D one
2578 Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
2579 gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2580 gp_XY newUV = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
2581 gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
2582 move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2584 if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2585 move.SquareMagnitude())
2587 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2588 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2589 MSG( "TOO LONG MOVE \t" <<
2590 "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2591 "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2592 "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2593 "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2597 (*link1)->Move( move );
2598 MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2599 << chain.front()->_mediumNode->GetID() <<"-"
2600 << chain.back ()->_mediumNode->GetID() <<
2601 " by " << move.Magnitude());
2603 } // loop on chains of links
2604 } // loop on 2 directions of propagation from quadrangle
2611 for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2612 if ( pLink->IsMoved() ) {
2613 //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2614 gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2615 GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
2620 //=======================================================================
2622 * \brief Iterator on ancestors of the given type
2624 //=======================================================================
2626 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2628 TopTools_ListIteratorOfListOfShape _ancIter;
2629 TopAbs_ShapeEnum _type;
2630 TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2631 : _ancIter( ancestors ), _type( type )
2633 if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
2637 return _ancIter.More();
2639 virtual const TopoDS_Shape* next()
2641 const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2642 if ( _ancIter.More() )
2643 for ( _ancIter.Next(); _ancIter.More(); _ancIter.Next())
2644 if ( _ancIter.Value().ShapeType() == _type )
2650 //=======================================================================
2652 * \brief Return iterator on ancestors of the given type
2654 //=======================================================================
2656 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2657 const SMESH_Mesh& mesh,
2658 TopAbs_ShapeEnum ancestorType)
2660 return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));