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 ));
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);
494 if ( nodePnt.Distance( surface->Value( U, V )) > tol )
496 MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
501 else if ( uv.Modulus() > numeric_limits<double>::min() )
503 ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
509 //=======================================================================
510 //function : GetProjector
511 //purpose : Return projector intitialized by given face without location, which is returned
512 //=======================================================================
514 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
515 TopLoc_Location& loc,
518 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
519 int faceID = GetMeshDS()->ShapeToIndex( F );
520 TID2Projector& i2proj = const_cast< TID2Projector&>( myFace2Projector );
521 TID2Projector::iterator i_proj = i2proj.find( faceID );
522 if ( i_proj == i2proj.end() )
524 if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
525 double U1, U2, V1, V2;
526 surface->Bounds(U1, U2, V1, V2);
527 GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
528 proj->Init( surface, U1, U2, V1, V2, tol );
529 i_proj = i2proj.insert( make_pair( faceID, proj )).first;
531 return *( i_proj->second );
536 gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
537 gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
538 gp_XY_FunPtr(Subtracted);
541 //=======================================================================
542 //function : applyIn2D
543 //purpose : Perform given operation on two 2d points in parameric space of given surface.
544 // It takes into account period of the surface. Use gp_XY_FunPtr macro
545 // to easily define pointer to function of gp_XY class.
546 //=======================================================================
548 gp_XY SMESH_MesherHelper::applyIn2D(const Handle(Geom_Surface)& surface,
552 const bool resultInPeriod)
554 Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
555 Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
556 if ( !isUPeriodic && !isVPeriodic )
559 // move uv2 not far than half-period from uv1
561 uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
563 uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
566 gp_XY res = fun( uv1, gp_XY(u2,v2) );
568 // move result within period
569 if ( resultInPeriod )
571 Standard_Real UF,UL,VF,VL;
572 surface->Bounds(UF,UL,VF,VL);
574 res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
576 res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
581 //=======================================================================
582 //function : GetMiddleUV
583 //purpose : Return middle UV taking in account surface period
584 //=======================================================================
586 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
590 return applyIn2D( surface, p1, p2, & AverageUV );
593 //=======================================================================
594 //function : GetNodeU
595 //purpose : Return node U on edge
596 //=======================================================================
598 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
599 const SMDS_MeshNode* n,
603 const SMDS_PositionPtr Pos = n->GetPosition();
604 if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE) {
605 const SMDS_EdgePosition* epos =
606 static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
607 param = epos->GetUParameter();
609 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX) {
610 SMESHDS_Mesh * meshDS = GetMeshDS();
611 int vertexID = n->GetPosition()->GetShapeId();
612 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
613 param = BRep_Tool::Parameter( V, E );
616 *check = CheckNodeU( E, n, param, BRep_Tool::Tolerance( E ));
620 //=======================================================================
621 //function : CheckNodeU
622 //purpose : Check and fix node U on an edge
623 // Return false if U is bad and could not be fixed
624 //=======================================================================
626 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
627 const SMDS_MeshNode* n,
630 const bool force) const
632 if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
634 // check that u is correct
635 TopLoc_Location loc; double f,l;
636 Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
637 if ( curve.IsNull() ) // degenerated edge
639 if ( u+tol < f || u-tol > l )
641 double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
647 gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
648 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
649 if ( nodePnt.Distance( curve->Value( u )) > tol )
651 // u incorrect, project the node to the curve
652 GeomAPI_ProjectPointOnCurve projector( nodePnt, curve, f, l );
653 if ( projector.NbPoints() < 1 )
655 MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
658 Quantity_Parameter U = projector.LowerDistanceParameter();
659 if ( nodePnt.Distance( curve->Value( U )) > tol )
661 MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
666 else if ( fabs( u ) > numeric_limits<double>::min() )
668 ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
675 //=======================================================================
676 //function : GetMediumNode
677 //purpose : Return existing or create new medium nodes between given ones
678 //=======================================================================
680 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
681 const SMDS_MeshNode* n2,
684 // Find existing node
686 SMESH_TLink link(n1,n2);
687 ItTLinkNode itLN = myTLinkNodeMap.find( link );
688 if ( itLN != myTLinkNodeMap.end() ) {
689 return (*itLN).second;
692 // Create medium node
695 SMESHDS_Mesh* meshDS = GetMeshDS();
697 // get type of shape for the new medium node
698 int faceID = -1, edgeID = -1;
699 const SMDS_PositionPtr Pos1 = n1->GetPosition();
700 const SMDS_PositionPtr Pos2 = n2->GetPosition();
702 if( myShape.IsNull() )
704 if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
705 faceID = Pos1->GetShapeId();
707 else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
708 faceID = Pos2->GetShapeId();
711 if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
712 edgeID = Pos1->GetShapeId();
714 if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
715 edgeID = Pos2->GetShapeId();
718 // get positions of the given nodes on shapes
719 TopoDS_Edge E; double u [2];
720 TopoDS_Face F; gp_XY uv[2];
721 bool uvOK[2] = { false, false };
722 TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
723 if ( faceID>0 || shapeType == TopAbs_FACE)
725 if( myShape.IsNull() )
726 F = TopoDS::Face(meshDS->IndexToShape(faceID));
728 F = TopoDS::Face(myShape);
731 uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
732 uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
734 else if (edgeID>0 || shapeType == TopAbs_EDGE)
736 if( myShape.IsNull() )
737 E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
739 E = TopoDS::Edge(myShape);
742 u[0] = GetNodeU(E,n1, force3d ? 0 : &uvOK[0]);
743 u[1] = GetNodeU(E,n2, force3d ? 0 : &uvOK[1]);
747 // we try to create medium node using UV parameters of
748 // nodes, else - medium between corresponding 3d points
751 if ( uvOK[0] && uvOK[1] )
753 if ( IsDegenShape( Pos1->GetShapeId() ))
754 if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
755 else uv[0].SetCoord( 2, uv[1].Coord( 2 ));
756 else if ( IsDegenShape( Pos2->GetShapeId() ))
757 if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
758 else uv[1].SetCoord( 2, uv[0].Coord( 2 ));
761 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
762 gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
763 gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
764 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
765 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
766 myTLinkNodeMap.insert(make_pair(link,n12));
770 else if ( !E.IsNull() )
773 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
776 Standard_Boolean isPeriodic = C->IsPeriodic();
779 Standard_Real Period = C->Period();
780 Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
781 Standard_Real pmid = (u[0]+p)/2.;
782 U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
787 gp_Pnt P = C->Value( U );
788 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
789 meshDS->SetNodeOnEdge(n12, edgeID, U);
790 myTLinkNodeMap.insert(make_pair(link,n12));
796 double x = ( n1->X() + n2->X() )/2.;
797 double y = ( n1->Y() + n2->Y() )/2.;
798 double z = ( n1->Z() + n2->Z() )/2.;
799 n12 = meshDS->AddNode(x,y,z);
802 gp_XY UV = ( uv[0] + uv[1] ) / 2.;
803 CheckNodeUV( F, n12, UV, BRep_Tool::Tolerance( F ), /*force=*/true);
804 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
806 else if ( !E.IsNull() )
808 double U = ( u[0] + u[1] ) / 2.;
809 CheckNodeU( E, n12, U, BRep_Tool::Tolerance( E ), /*force=*/true);
810 meshDS->SetNodeOnEdge(n12, edgeID, U);
814 meshDS->SetNodeInVolume(n12, myShapeID);
816 myTLinkNodeMap.insert( make_pair( link, n12 ));
820 //=======================================================================
822 //purpose : Creates a node
823 //=======================================================================
825 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
827 SMESHDS_Mesh * meshDS = GetMeshDS();
828 SMDS_MeshNode* node = 0;
830 node = meshDS->AddNodeWithID( x, y, z, ID );
832 node = meshDS->AddNode( x, y, z );
833 if ( mySetElemOnShape && myShapeID > 0 ) {
834 switch ( myShape.ShapeType() ) {
835 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
836 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
837 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID); break;
838 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID); break;
839 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
846 //=======================================================================
848 //purpose : Creates quadratic or linear edge
849 //=======================================================================
851 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
852 const SMDS_MeshNode* n2,
856 SMESHDS_Mesh * meshDS = GetMeshDS();
858 SMDS_MeshEdge* edge = 0;
859 if (myCreateQuadratic) {
860 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
862 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
864 edge = meshDS->AddEdge(n1, n2, n12);
868 edge = meshDS->AddEdgeWithID(n1, n2, id);
870 edge = meshDS->AddEdge(n1, n2);
873 if ( mySetElemOnShape && myShapeID > 0 )
874 meshDS->SetMeshElementOnShape( edge, myShapeID );
879 //=======================================================================
881 //purpose : Creates quadratic or linear triangle
882 //=======================================================================
884 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
885 const SMDS_MeshNode* n2,
886 const SMDS_MeshNode* n3,
890 SMESHDS_Mesh * meshDS = GetMeshDS();
891 SMDS_MeshFace* elem = 0;
893 if( n1==n2 || n2==n3 || n3==n1 )
896 if(!myCreateQuadratic) {
898 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
900 elem = meshDS->AddFace(n1, n2, n3);
903 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
904 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
905 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
908 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
910 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
912 if ( mySetElemOnShape && myShapeID > 0 )
913 meshDS->SetMeshElementOnShape( elem, myShapeID );
918 //=======================================================================
920 //purpose : Creates quadratic or linear quadrangle
921 //=======================================================================
923 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
924 const SMDS_MeshNode* n2,
925 const SMDS_MeshNode* n3,
926 const SMDS_MeshNode* n4,
930 SMESHDS_Mesh * meshDS = GetMeshDS();
931 SMDS_MeshFace* elem = 0;
934 return AddFace(n1,n3,n4,id,force3d);
937 return AddFace(n1,n2,n4,id,force3d);
940 return AddFace(n1,n2,n3,id,force3d);
943 return AddFace(n1,n2,n4,id,force3d);
946 return AddFace(n1,n2,n3,id,force3d);
949 return AddFace(n1,n2,n3,id,force3d);
952 if(!myCreateQuadratic) {
954 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
956 elem = meshDS->AddFace(n1, n2, n3, n4);
959 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
960 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
961 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
962 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
965 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
967 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
969 if ( mySetElemOnShape && myShapeID > 0 )
970 meshDS->SetMeshElementOnShape( elem, myShapeID );
975 //=======================================================================
976 //function : AddVolume
977 //purpose : Creates quadratic or linear prism
978 //=======================================================================
980 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
981 const SMDS_MeshNode* n2,
982 const SMDS_MeshNode* n3,
983 const SMDS_MeshNode* n4,
984 const SMDS_MeshNode* n5,
985 const SMDS_MeshNode* n6,
989 SMESHDS_Mesh * meshDS = GetMeshDS();
990 SMDS_MeshVolume* elem = 0;
991 if(!myCreateQuadratic) {
993 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
995 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
998 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
999 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1000 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1002 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1003 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1004 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1006 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1007 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1008 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1011 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
1012 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1014 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1015 n12, n23, n31, n45, n56, n64, n14, n25, n36);
1017 if ( mySetElemOnShape && myShapeID > 0 )
1018 meshDS->SetMeshElementOnShape( elem, myShapeID );
1023 //=======================================================================
1024 //function : AddVolume
1025 //purpose : Creates quadratic or linear tetrahedron
1026 //=======================================================================
1028 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1029 const SMDS_MeshNode* n2,
1030 const SMDS_MeshNode* n3,
1031 const SMDS_MeshNode* n4,
1035 SMESHDS_Mesh * meshDS = GetMeshDS();
1036 SMDS_MeshVolume* elem = 0;
1037 if(!myCreateQuadratic) {
1039 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1041 elem = meshDS->AddVolume(n1, n2, n3, n4);
1044 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1045 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1046 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1048 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1049 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1050 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1053 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1055 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1057 if ( mySetElemOnShape && myShapeID > 0 )
1058 meshDS->SetMeshElementOnShape( elem, myShapeID );
1063 //=======================================================================
1064 //function : AddVolume
1065 //purpose : Creates quadratic or linear pyramid
1066 //=======================================================================
1068 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1069 const SMDS_MeshNode* n2,
1070 const SMDS_MeshNode* n3,
1071 const SMDS_MeshNode* n4,
1072 const SMDS_MeshNode* n5,
1076 SMDS_MeshVolume* elem = 0;
1077 if(!myCreateQuadratic) {
1079 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1081 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1084 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1085 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1086 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1087 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1089 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1090 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1091 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1092 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1095 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
1100 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
1102 n15, n25, n35, n45);
1104 if ( mySetElemOnShape && myShapeID > 0 )
1105 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1110 //=======================================================================
1111 //function : AddVolume
1112 //purpose : Creates quadratic or linear hexahedron
1113 //=======================================================================
1115 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1116 const SMDS_MeshNode* n2,
1117 const SMDS_MeshNode* n3,
1118 const SMDS_MeshNode* n4,
1119 const SMDS_MeshNode* n5,
1120 const SMDS_MeshNode* n6,
1121 const SMDS_MeshNode* n7,
1122 const SMDS_MeshNode* n8,
1126 SMESHDS_Mesh * meshDS = GetMeshDS();
1127 SMDS_MeshVolume* elem = 0;
1128 if(!myCreateQuadratic) {
1130 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1132 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1135 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1136 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1137 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1138 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1140 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1141 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1142 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1143 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1145 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1146 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1147 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1148 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1151 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1152 n12, n23, n34, n41, n56, n67,
1153 n78, n85, n15, n26, n37, n48, id);
1155 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1156 n12, n23, n34, n41, n56, n67,
1157 n78, n85, n15, n26, n37, n48);
1159 if ( mySetElemOnShape && myShapeID > 0 )
1160 meshDS->SetMeshElementOnShape( elem, myShapeID );
1165 //=======================================================================
1166 //function : LoadNodeColumns
1167 //purpose : Load nodes bound to face into a map of node columns
1168 //=======================================================================
1170 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1171 const TopoDS_Face& theFace,
1172 const TopoDS_Edge& theBaseEdge,
1173 SMESHDS_Mesh* theMesh)
1175 SMESHDS_SubMesh* faceSubMesh = theMesh->MeshElements( theFace );
1176 if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
1179 // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
1181 map< double, const SMDS_MeshNode*> sortedBaseNodes;
1182 if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
1183 || sortedBaseNodes.size() < 2 )
1186 int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
1187 map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
1188 double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
1189 for ( ; u_n != sortedBaseNodes.end(); u_n++ )
1191 double par = ( u_n->first - f ) / range;
1192 vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
1193 nCol.resize( nbRows );
1194 nCol[0] = u_n->second;
1197 // fill theParam2ColumnMap column by column by passing from nodes on
1198 // theBaseEdge up via mesh faces on theFace
1200 TParam2ColumnMap::iterator par_nVec_2 = theParam2ColumnMap.begin();
1201 TParam2ColumnMap::iterator par_nVec_1 = par_nVec_2++;
1202 TIDSortedElemSet emptySet, avoidSet;
1203 for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
1205 vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
1206 vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
1208 int i1, i2, iRow = 0;
1209 const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
1210 // find face sharing node n1 and n2 and belonging to faceSubMesh
1211 while ( const SMDS_MeshElement* face =
1212 SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
1214 if ( faceSubMesh->Contains( face ))
1216 int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
1219 n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
1220 n2 = face->GetNode( (i1+2) % 4 );
1221 if ( ++iRow >= nbRows )
1227 avoidSet.insert( face );
1229 if ( iRow + 1 < nbRows ) // compact if necessary
1230 nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
1235 //=======================================================================
1236 //function : NbAncestors
1237 //purpose : Return number of unique ancestors of the shape
1238 //=======================================================================
1240 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1241 const SMESH_Mesh& mesh,
1242 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
1244 TopTools_MapOfShape ancestors;
1245 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1246 for ( ; ansIt.More(); ansIt.Next() ) {
1247 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1248 ancestors.Add( ansIt.Value() );
1250 return ancestors.Extent();
1253 //=======================================================================
1254 //function : GetSubShapeOri
1255 //purpose : Return orientation of sub-shape in the main shape
1256 //=======================================================================
1258 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
1259 const TopoDS_Shape& subShape)
1261 TopAbs_Orientation ori = TopAbs_Orientation(-1);
1262 if ( !shape.IsNull() && !subShape.IsNull() )
1264 TopExp_Explorer e( shape, subShape.ShapeType() );
1265 if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
1266 e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
1267 for ( ; e.More(); e.Next())
1268 if ( subShape.IsSame( e.Current() ))
1271 ori = e.Current().Orientation();
1276 //=======================================================================
1277 //function : IsSubShape
1279 //=======================================================================
1281 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
1282 const TopoDS_Shape& mainShape )
1284 if ( !shape.IsNull() && !mainShape.IsNull() )
1286 for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
1289 if ( shape.IsSame( exp.Current() ))
1292 SCRUTE((shape.IsNull()));
1293 SCRUTE((mainShape.IsNull()));
1297 //=======================================================================
1298 //function : IsSubShape
1300 //=======================================================================
1302 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
1304 if ( shape.IsNull() || !aMesh )
1307 aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
1309 shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape );
1312 //=======================================================================
1313 //function : IsQuadraticMesh
1314 //purpose : Check mesh without geometry for: if all elements on this shape are quadratic,
1315 // quadratic elements will be created.
1316 // Used then generated 3D mesh without geometry.
1317 //=======================================================================
1319 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1321 int NbAllEdgsAndFaces=0;
1322 int NbQuadFacesAndEdgs=0;
1323 int NbFacesAndEdges=0;
1324 //All faces and edges
1325 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1327 //Quadratic faces and edges
1328 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1330 //Linear faces and edges
1331 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1333 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1335 return SMESH_MesherHelper::QUADRATIC;
1337 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1339 return SMESH_MesherHelper::LINEAR;
1342 //Mesh with both type of elements
1343 return SMESH_MesherHelper::COMP;
1346 //=======================================================================
1347 //function : GetOtherParam
1348 //purpose : Return an alternative parameter for a node on seam
1349 //=======================================================================
1351 double SMESH_MesherHelper::GetOtherParam(const double param) const
1353 int i = myParIndex & U_periodic ? 0 : 1;
1354 return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1357 //=======================================================================
1358 namespace { // Structures used by FixQuadraticElements()
1359 //=======================================================================
1361 #define __DMP__(txt) \
1363 #define MSG(txt) __DMP__(txt<<endl)
1364 #define MSGBEG(txt) __DMP__(txt)
1366 const double straightTol2 = 1e-33; // to detect straing links
1369 // ---------------------------------------
1371 * \brief Quadratic link knowing its faces
1373 struct QLink: public SMESH_TLink
1375 const SMDS_MeshNode* _mediumNode;
1376 mutable vector<const QFace* > _faces;
1377 mutable gp_Vec _nodeMove;
1378 mutable int _nbMoves;
1380 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1381 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1383 //if ( MediumPos() != SMDS_TOP_3DSPACE )
1384 _nodeMove = MediumPnt() - MiddlePnt();
1386 void SetContinuesFaces() const;
1387 const QFace* GetContinuesFace( const QFace* face ) const;
1388 bool OnBoundary() const;
1389 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1390 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1392 SMDS_TypeOfPosition MediumPos() const
1393 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1394 SMDS_TypeOfPosition EndPos(bool isSecond) const
1395 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1396 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1397 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1399 void Move(const gp_Vec& move, bool sum=false) const
1400 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1401 gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1402 bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1403 bool IsStraight() const { return _nodeMove.SquareMagnitude() <= straightTol2; }
1405 bool operator<(const QLink& other) const {
1406 return (node1()->GetID() == other.node1()->GetID() ?
1407 node2()->GetID() < other.node2()->GetID() :
1408 node1()->GetID() < other.node1()->GetID());
1410 struct PtrComparator {
1411 bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1414 // ---------------------------------------------------------
1416 * \brief Link in the chain of links; it connects two faces
1420 const QLink* _qlink;
1421 mutable const QFace* _qfaces[2];
1423 TChainLink(const QLink* qlink=0):_qlink(qlink) {
1424 _qfaces[0] = _qfaces[1] = 0;
1426 void SetFace(const QFace* face) { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1428 bool IsBoundary() const { return !_qfaces[1]; }
1430 void RemoveFace( const QFace* face ) const
1431 { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1433 const QFace* NextFace( const QFace* f ) const
1434 { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1436 const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1437 { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1439 bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1441 operator bool() const { return (_qlink); }
1443 const QLink* operator->() const { return _qlink; }
1445 gp_Vec Normal() const;
1447 // --------------------------------------------------------------------
1448 typedef list< TChainLink > TChain;
1449 typedef set < TChainLink > TLinkSet;
1450 typedef TLinkSet::const_iterator TLinkInSet;
1452 const int theFirstStep = 5;
1454 enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1455 // --------------------------------------------------------------------
1457 * \brief Face shared by two volumes and bound by QLinks
1459 struct QFace: public TIDSortedElemSet
1461 mutable const SMDS_MeshElement* _volumes[2];
1462 mutable vector< const QLink* > _sides;
1463 mutable bool _sideIsAdded[4]; // added in chain of links
1466 mutable const SMDS_MeshElement* _face;
1469 QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1471 void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1473 int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1475 void AddSelfToLinks() const {
1476 for ( int i = 0; i < _sides.size(); ++i )
1477 _sides[i]->_faces.push_back( this );
1479 int LinkIndex( const QLink* side ) const {
1480 for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1483 bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1485 bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1487 int i = LinkIndex( link._qlink );
1488 if ( i < 0 ) return true;
1489 _sideIsAdded[i] = true;
1490 link.SetFace( this );
1491 // continue from opposite link
1492 return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1494 bool IsBoundary() const { return !_volumes[1]; }
1496 bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1498 TLinkInSet GetBoundaryLink( const TLinkSet& links,
1499 const TChainLink& avoidLink,
1500 TLinkInSet * notBoundaryLink = 0,
1501 const SMDS_MeshNode* nodeToContain = 0,
1502 bool * isAdjacentUsed = 0,
1503 int nbRecursionsLeft = -1) const;
1505 TLinkInSet GetLinkByNode( const TLinkSet& links,
1506 const TChainLink& avoidLink,
1507 const SMDS_MeshNode* nodeToContain) const;
1509 const SMDS_MeshNode* GetNodeInFace() const {
1510 for ( int iL = 0; iL < _sides.size(); ++iL )
1511 if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1515 gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1517 double MoveByBoundary( const TChainLink& theLink,
1518 const gp_Vec& theRefVec,
1519 const TLinkSet& theLinks,
1520 SMESH_MesherHelper* theFaceHelper=0,
1521 const double thePrevLen=0,
1522 const int theStep=theFirstStep,
1523 gp_Vec* theLinkNorm=0,
1524 double theSign=1.0) const;
1527 //================================================================================
1529 * \brief Dump QLink and QFace
1531 ostream& operator << (ostream& out, const QLink& l)
1533 out <<"QLink nodes: "
1534 << l.node1()->GetID() << " - "
1535 << l._mediumNode->GetID() << " - "
1536 << l.node2()->GetID() << endl;
1539 ostream& operator << (ostream& out, const QFace& f)
1541 out <<"QFace nodes: "/*<< &f << " "*/;
1542 for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
1543 out << (*n)->GetID() << " ";
1544 out << " \tvolumes: "
1545 << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1546 << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1547 out << " \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1551 //================================================================================
1553 * \brief Construct QFace from QLinks
1555 //================================================================================
1557 QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
1559 _volumes[0] = _volumes[1] = 0;
1561 _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1562 _normal.SetCoord(0,0,0);
1563 for ( int i = 1; i < _sides.size(); ++i ) {
1564 const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1565 insert( l1->node1() ); insert( l1->node2() );
1567 gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1568 gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1569 if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1573 double normSqSize = _normal.SquareMagnitude();
1574 if ( normSqSize > numeric_limits<double>::min() )
1575 _normal /= sqrt( normSqSize );
1577 _normal.SetCoord(1e-33,0,0);
1583 //================================================================================
1585 * \brief Make up chain of links
1586 * \param iSide - link to add first
1587 * \param chain - chain to fill in
1588 * \param pos - postion of medium nodes the links should have
1589 * \param error - out, specifies what is wrong
1590 * \retval bool - false if valid chain can't be built; "valid" means that links
1591 * of the chain belongs to rectangles bounding hexahedrons
1593 //================================================================================
1595 bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1597 if ( iSide >= _sides.size() ) // wrong argument iSide
1599 if ( _sideIsAdded[ iSide ]) // already in chain
1602 if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1604 list< const QFace* > faces( 1, this );
1605 for (list< const QFace* >::iterator fIt = faces.begin(); fIt != faces.end(); ++fIt ) {
1606 const QFace* face = *fIt;
1607 for ( int i = 0; i < face->_sides.size(); ++i ) {
1608 if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
1609 face->_sideIsAdded[i] = true;
1610 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
1611 chLink->SetFace( face );
1612 if ( face->_sides[i]->MediumPos() >= pos )
1613 if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
1614 faces.push_back( contFace );
1618 if ( error < ERR_TRI )
1622 _sideIsAdded[iSide] = true; // not to add this link to chain again
1623 const QLink* link = _sides[iSide];
1627 // add link into chain
1628 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1629 chLink->SetFace( this );
1632 // propagate from quadrangle to neighbour faces
1633 if ( link->MediumPos() >= pos ) {
1634 int nbLinkFaces = link->_faces.size();
1635 if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1636 // hexahedral mesh or boundary quadrangles - goto a continous face
1637 if ( const QFace* f = link->GetContinuesFace( this ))
1638 return f->GetLinkChain( *chLink, chain, pos, error );
1641 TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1642 for ( int i = 0; i < nbLinkFaces; ++i )
1643 if ( link->_faces[i] )
1644 link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1645 if ( error < ERR_PRISM )
1653 //================================================================================
1655 * \brief Return a boundary link of the triangle face
1656 * \param links - set of all links
1657 * \param avoidLink - link not to return
1658 * \param notBoundaryLink - out, neither the returned link nor avoidLink
1659 * \param nodeToContain - node the returned link must contain; if provided, search
1660 * also performed on adjacent faces
1661 * \param isAdjacentUsed - returns true if link is found in adjacent faces
1662 * \param nbRecursionsLeft - to limit recursion
1664 //================================================================================
1666 TLinkInSet QFace::GetBoundaryLink( const TLinkSet& links,
1667 const TChainLink& avoidLink,
1668 TLinkInSet * notBoundaryLink,
1669 const SMDS_MeshNode* nodeToContain,
1670 bool * isAdjacentUsed,
1671 int nbRecursionsLeft) const
1673 TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1675 typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1676 TFaceLinkList adjacentFaces;
1678 for ( int iL = 0; iL < _sides.size(); ++iL )
1680 if ( avoidLink._qlink == _sides[iL] )
1682 TLinkInSet link = links.find( _sides[iL] );
1683 if ( link == linksEnd ) continue;
1684 if ( (*link)->MediumPos() > SMDS_TOP_FACE )
1685 continue; // We work on faces here, don't go into a volume
1688 if ( link->IsBoundary() ) {
1689 if ( !nodeToContain ||
1690 (*link)->node1() == nodeToContain ||
1691 (*link)->node2() == nodeToContain )
1693 boundaryLink = link;
1694 if ( !notBoundaryLink ) break;
1697 else if ( notBoundaryLink ) {
1698 *notBoundaryLink = link;
1699 if ( boundaryLink != linksEnd ) break;
1702 if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
1703 if ( const QFace* adj = link->NextFace( this ))
1704 if ( adj->Contains( nodeToContain ))
1705 adjacentFaces.push_back( make_pair( adj, link ));
1708 if ( isAdjacentUsed ) *isAdjacentUsed = false;
1709 if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
1711 if ( nbRecursionsLeft < 0 )
1712 nbRecursionsLeft = nodeToContain->NbInverseElements();
1713 TFaceLinkList::iterator adj = adjacentFaces.begin();
1714 for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1715 boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
1716 isAdjacentUsed, nbRecursionsLeft-1);
1717 if ( isAdjacentUsed ) *isAdjacentUsed = true;
1719 return boundaryLink;
1721 //================================================================================
1723 * \brief Return a link ending at the given node but not avoidLink
1725 //================================================================================
1727 TLinkInSet QFace::GetLinkByNode( const TLinkSet& links,
1728 const TChainLink& avoidLink,
1729 const SMDS_MeshNode* nodeToContain) const
1731 for ( int i = 0; i < _sides.size(); ++i )
1732 if ( avoidLink._qlink != _sides[i] &&
1733 (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1734 return links.find( _sides[ i ]);
1738 //================================================================================
1740 * \brief Return normal to the i-th side pointing outside the face
1742 //================================================================================
1744 gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1746 gp_Vec norm, vecOut;
1747 // if ( uvHelper ) {
1748 // TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
1749 // const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
1750 // gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
1751 // gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
1752 // norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
1754 // const QLink* otherLink = _sides[(i + 1) % _sides.size()];
1755 // const SMDS_MeshNode* otherNode =
1756 // otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
1757 // gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
1758 // vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
1761 norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
1762 gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
1763 XYZ( _sides[0]->node2() ) +
1764 XYZ( _sides[1]->node1() )) / 3.;
1765 vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
1767 if ( norm * vecOut < 0 )
1769 double mag2 = norm.SquareMagnitude();
1770 if ( mag2 > numeric_limits<double>::min() )
1771 norm /= sqrt( mag2 );
1774 //================================================================================
1776 * \brief Move medium node of theLink according to its distance from boundary
1777 * \param theLink - link to fix
1778 * \param theRefVec - movement of boundary
1779 * \param theLinks - all adjacent links of continous triangles
1780 * \param theFaceHelper - helper is not used so far
1781 * \param thePrevLen - distance from the boundary
1782 * \param theStep - number of steps till movement propagation limit
1783 * \param theLinkNorm - out normal to theLink
1784 * \param theSign - 1 or -1 depending on movement of boundary
1785 * \retval double - distance from boundary to propagation limit or other boundary
1787 //================================================================================
1789 double QFace::MoveByBoundary( const TChainLink& theLink,
1790 const gp_Vec& theRefVec,
1791 const TLinkSet& theLinks,
1792 SMESH_MesherHelper* theFaceHelper,
1793 const double thePrevLen,
1795 gp_Vec* theLinkNorm,
1796 double theSign) const
1799 return thePrevLen; // propagation limit reached
1801 int iL; // index of theLink
1802 for ( iL = 0; iL < _sides.size(); ++iL )
1803 if ( theLink._qlink == _sides[ iL ])
1806 MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
1807 <<" thePrevLen " << thePrevLen);
1808 MSG(string(theStep,'.')<<" "<<*theLink._qlink);
1810 gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
1811 double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
1812 if ( theStep == theFirstStep )
1813 theSign = refProj < 0. ? -1. : 1.;
1814 else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
1815 return thePrevLen; // to propagate movement forward only, not in side dir or backward
1817 int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
1818 TLinkInSet link1 = theLinks.find( _sides[iL1] );
1819 TLinkInSet link2 = theLinks.find( _sides[iL2] );
1820 if ( link1 == theLinks.end() || link2 == theLinks.end() )
1822 const QFace* f1 = link1->NextFace( this ); // adjacent faces
1823 const QFace* f2 = link2->NextFace( this );
1825 // propagate to adjacent faces till limit step or boundary
1826 double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
1827 double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
1828 gp_Vec linkDir1, linkDir2;
1832 len1 = f1->MoveByBoundary
1833 ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
1835 linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
1837 MSG( " --------------- EXCEPTION");
1843 len2 = f2->MoveByBoundary
1844 ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
1846 linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
1848 MSG( " --------------- EXCEPTION");
1853 if ( theStep != theFirstStep )
1855 // choose chain length by direction of propagation most codirected with theRefVec
1856 bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
1857 fullLen = choose1 ? len1 : len2;
1858 double r = thePrevLen / fullLen;
1860 gp_Vec move = linkNorm * refProj * ( 1 - r );
1861 theLink->Move( move, true );
1863 MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
1864 " by " << refProj * ( 1 - r ) << " following " <<
1865 (choose1 ? *link1->_qlink : *link2->_qlink));
1867 if ( theLinkNorm ) *theLinkNorm = linkNorm;
1872 //================================================================================
1874 * \brief Find pairs of continues faces
1876 //================================================================================
1878 void QLink::SetContinuesFaces() const
1880 // x0 x - QLink, [-|] - QFace, v - volume
1882 // | Between _faces of link x2 two vertical faces are continues
1883 // x1----x2-----x3 and two horizontal faces are continues. We set vertical faces
1884 // | to _faces[0] and _faces[1] and horizontal faces to
1885 // v2 | v3 _faces[2] and _faces[3] (or vise versa).
1888 if ( _faces.empty() )
1891 for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
1893 // look for a face bounding none of volumes bound by _faces[0]
1894 bool sameVol = false;
1895 int nbVol = _faces[iF]->NbVolumes();
1896 for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
1897 sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
1898 _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
1902 if ( iFaceCont > 0 ) // continues faces found, set one by the other
1904 if ( iFaceCont != 1 )
1905 std::swap( _faces[1], _faces[iFaceCont] );
1907 else if ( _faces.size() > 1 ) // not found, set NULL by the first face
1909 _faces.insert( ++_faces.begin(), 0 );
1912 //================================================================================
1914 * \brief Return a face continues to the given one
1916 //================================================================================
1918 const QFace* QLink::GetContinuesFace( const QFace* face ) const
1920 for ( int i = 0; i < _faces.size(); ++i ) {
1921 if ( _faces[i] == face ) {
1922 int iF = i < 2 ? 1-i : 5-i;
1923 return iF < _faces.size() ? _faces[iF] : 0;
1928 //================================================================================
1930 * \brief True if link is on mesh boundary
1932 //================================================================================
1934 bool QLink::OnBoundary() const
1936 for ( int i = 0; i < _faces.size(); ++i )
1937 if (_faces[i] && _faces[i]->IsBoundary()) return true;
1940 //================================================================================
1942 * \brief Return normal of link of the chain
1944 //================================================================================
1946 gp_Vec TChainLink::Normal() const {
1948 if (_qfaces[0]) norm = _qfaces[0]->_normal;
1949 if (_qfaces[1]) norm += _qfaces[1]->_normal;
1952 //================================================================================
1954 * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
1956 //================================================================================
1958 void fixPrism( TChain& allLinks )
1960 // separate boundary links from internal ones
1961 typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
1962 QLinkSet interLinks, bndLinks1, bndLink2;
1964 bool isCurved = false;
1965 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
1966 if ( (*lnk)->OnBoundary() )
1967 bndLinks1.insert( lnk->_qlink );
1969 interLinks.insert( lnk->_qlink );
1970 isCurved = isCurved || !(*lnk)->IsStraight();
1973 return; // no need to move
1975 QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
1977 while ( !interLinks.empty() && !curBndLinks->empty() )
1979 // propagate movement from boundary links to connected internal links
1980 QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
1981 for ( ; bnd != bndEnd; ++bnd )
1983 const QLink* bndLink = *bnd;
1984 for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
1986 const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
1987 if ( !face ) continue;
1988 // find and move internal link opposite to bndLink within the face
1989 int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
1990 const QLink* interLink = face->_sides[ interInd ];
1991 QLinkSet::iterator pInterLink = interLinks.find( interLink );
1992 if ( pInterLink == interLinks.end() ) continue; // not internal link
1993 interLink->Move( bndLink->_nodeMove );
1994 // treated internal links become new boundary ones
1995 interLinks. erase( pInterLink );
1996 newBndLinks->insert( interLink );
1999 curBndLinks->clear();
2000 std::swap( curBndLinks, newBndLinks );
2004 //================================================================================
2006 * \brief Fix links of continues triangles near curved boundary
2008 //================================================================================
2010 void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2012 if ( allLinks.empty() ) return;
2014 TLinkSet linkSet( allLinks.begin(), allLinks.end());
2015 TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2017 // move in 2d if we are on geom face
2018 // TopoDS_Face face;
2019 // TopLoc_Location loc;
2020 // SMESH_MesherHelper faceHelper( *helper.GetMesh());
2021 // while ( linkIt->IsBoundary()) ++linkIt;
2022 // if ( linkIt == linksEnd ) return;
2023 // if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
2024 // bool checkPos = true;
2025 // TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
2026 // if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2027 // face = TopoDS::Face( f );
2028 // helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
2032 // faceHelper.SetSubShape( face );
2035 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2037 if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2039 // if ( !face.IsNull() ) {
2040 // const SMDS_MeshNode* inFaceNode =
2041 // faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
2042 // gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
2043 // gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
2044 // gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
2045 // gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2046 // gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
2047 // linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
2050 linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2056 //================================================================================
2058 * \brief Detect rectangular structure of links and build chains from them
2060 //================================================================================
2062 enum TSplitTriaResult {
2063 _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2064 _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK };
2066 TSplitTriaResult splitTrianglesIntoChains( TChain & allLinks,
2067 vector< TChain> & resultChains,
2068 SMDS_TypeOfPosition pos )
2070 // put links in the set and evalute number of result chains by number of boundary links
2073 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2074 linkSet.insert( *lnk );
2075 nbBndLinks += lnk->IsBoundary();
2077 resultChains.clear();
2078 resultChains.reserve( nbBndLinks / 2 );
2080 TLinkInSet linkIt, linksEnd = linkSet.end();
2082 // find a boundary link with corner node; corner node has position pos-2
2083 // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2085 SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2086 const SMDS_MeshNode* corner = 0;
2087 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2088 if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2093 TLinkInSet startLink = linkIt;
2094 const SMDS_MeshNode* startCorner = corner;
2095 vector< TChain* > rowChains;
2098 while ( startLink != linksEnd) // loop on columns
2100 // We suppose we have a rectangular structure like shown here. We have found a
2101 // corner of the rectangle (startCorner) and a boundary link sharing
2102 // |/ |/ | the startCorner (startLink). We are going to loop on rows of the
2103 // --o---o---o structure making several chains at once. One chain (columnChain)
2104 // |\ | /| starts at startLink and continues upward (we look at the structure
2105 // \ | \ | / | from such point that startLink is on the bottom of the structure).
2106 // \| \|/ | While going upward we also fill horizontal chains (rowChains) we
2107 // --o---o---o encounter.
2109 // / | \ | \ | startCorner
2114 if ( resultChains.size() == nbBndLinks / 2 )
2116 resultChains.push_back( TChain() );
2117 TChain& columnChain = resultChains.back();
2119 TLinkInSet botLink = startLink; // current horizontal link to go up from
2120 corner = startCorner; // current corner the botLink ends at
2122 while ( botLink != linksEnd ) // loop on rows
2124 // add botLink to the columnChain
2125 columnChain.push_back( *botLink );
2127 const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2129 { // the column ends
2130 linkSet.erase( botLink );
2131 if ( iRow != rowChains.size() )
2132 return _FEW_ROWS; // different nb of rows in columns
2135 // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2136 // link ending at <corner> (sideLink); there are two cases:
2137 // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2138 // since midQuadLink is not at boundary while sideLink is.
2139 // 2) midQuadLink ends at <corner>
2141 TLinkInSet midQuadLink = linksEnd;
2142 TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2144 if ( isCase2 ) { // find midQuadLink among links of botTria
2145 midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2146 if ( midQuadLink->IsBoundary() )
2147 return _BAD_MIDQUAD;
2149 if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2150 return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2153 columnChain.push_back( *midQuadLink );
2154 if ( iRow >= rowChains.size() ) {
2156 return _MANY_ROWS; // different nb of rows in columns
2157 if ( resultChains.size() == nbBndLinks / 2 )
2159 resultChains.push_back( TChain() );
2160 rowChains.push_back( & resultChains.back() );
2162 rowChains[iRow]->push_back( *sideLink );
2163 rowChains[iRow]->push_back( *midQuadLink );
2165 const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2169 // prepare startCorner and startLink for the next column
2170 startCorner = startLink->NextNode( startCorner );
2172 startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2174 startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2175 // check if no more columns remains
2176 if ( startLink != linksEnd ) {
2177 const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2178 if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2179 startLink = linksEnd; // startLink bounds upTria or botTria
2180 else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2184 // find bottom link and corner for the next row
2185 corner = sideLink->NextNode( corner );
2186 // next bottom link ends at the new corner
2187 linkSet.erase( botLink );
2188 botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2189 if ( botLink == linksEnd || botLink == (isCase2 ? midQuadLink : sideLink))
2191 linkSet.erase( midQuadLink );
2192 linkSet.erase( sideLink );
2194 // make faces neighboring the found ones be boundary
2195 if ( startLink != linksEnd ) {
2196 const QFace* tria = isCase2 ? botTria : upTria;
2197 for ( int iL = 0; iL < 3; ++iL ) {
2198 linkIt = linkSet.find( tria->_sides[iL] );
2199 if ( linkIt != linksEnd )
2200 linkIt->RemoveFace( tria );
2203 if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2204 botLink->RemoveFace( upTria ); // make next botTria first in vector
2211 // In the linkSet, there must remain the last links of rowChains; add them
2212 if ( linkSet.size() != rowChains.size() )
2213 return _BAD_SET_SIZE;
2214 for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2215 // find the link (startLink) ending at startCorner
2217 for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2218 if ( (*startLink)->node1() == startCorner ) {
2219 corner = (*startLink)->node2(); break;
2221 else if ( (*startLink)->node2() == startCorner) {
2222 corner = (*startLink)->node1(); break;
2225 if ( startLink == linksEnd )
2227 rowChains[ iRow ]->push_back( *startLink );
2228 linkSet.erase( startLink );
2229 startCorner = corner;
2236 //=======================================================================
2238 * \brief Move medium nodes of faces and volumes to fix distorted elements
2239 * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2241 * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2243 //=======================================================================
2245 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2247 // apply algorithm to solids or geom faces
2248 // ----------------------------------------------
2249 if ( myShape.IsNull() ) {
2250 if ( !myMesh->HasShapeToMesh() ) return;
2251 SetSubShape( myMesh->GetShapeToMesh() );
2253 TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2254 for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2255 faces.Add( f.Current() );
2257 for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2258 if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2259 for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2260 faces.Add( f.Current() );
2262 else { // fix nodes in the solid and its faces
2263 SMESH_MesherHelper h(*myMesh);
2264 h.SetSubShape( s.Current() );
2265 h.FixQuadraticElements(false);
2268 // fix nodes on geom faces
2269 for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2270 SMESH_MesherHelper h(*myMesh);
2271 h.SetSubShape( fIt.Key() );
2272 h.FixQuadraticElements(true);
2277 // Find out type of elements and get iterator on them
2278 // ---------------------------------------------------
2280 SMDS_ElemIteratorPtr elemIt;
2281 SMDSAbs_ElementType elemType = SMDSAbs_All;
2283 SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2286 if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2287 elemIt = smDS->GetElements();
2288 if ( elemIt->more() ) {
2289 elemType = elemIt->next()->GetType();
2290 elemIt = smDS->GetElements();
2293 if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2296 // Fill in auxiliary data structures
2297 // ----------------------------------
2301 set< QLink >::iterator pLink;
2302 set< QFace >::iterator pFace;
2304 bool isCurved = false;
2305 bool hasRectFaces = false;
2306 set<int> nbElemNodeSet;
2308 if ( elemType == SMDSAbs_Volume )
2310 SMDS_VolumeTool volTool;
2311 while ( elemIt->more() ) // loop on volumes
2313 const SMDS_MeshElement* vol = elemIt->next();
2314 if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2316 for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2318 int nbN = volTool.NbFaceNodes( iF );
2319 nbElemNodeSet.insert( nbN );
2320 const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2321 vector< const QLink* > faceLinks( nbN/2 );
2322 for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2325 QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2326 pLink = links.insert( link ).first;
2327 faceLinks[ iN/2 ] = & *pLink;
2329 isCurved = !link.IsStraight();
2330 if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2331 return; // already fixed
2334 pFace = faces.insert( QFace( faceLinks )).first;
2335 if ( pFace->NbVolumes() == 0 )
2336 pFace->AddSelfToLinks();
2337 pFace->SetVolume( vol );
2338 hasRectFaces = hasRectFaces ||
2339 ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2340 volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2343 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2345 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2346 faceNodes[4],faceNodes[6] );
2350 set< QLink >::iterator pLink = links.begin();
2351 for ( ; pLink != links.end(); ++pLink )
2352 pLink->SetContinuesFaces();
2356 while ( elemIt->more() ) // loop on faces
2358 const SMDS_MeshElement* face = elemIt->next();
2359 if ( !face->IsQuadratic() )
2361 nbElemNodeSet.insert( face->NbNodes() );
2362 int nbN = face->NbNodes()/2;
2363 vector< const QLink* > faceLinks( nbN );
2364 for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2367 QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2368 pLink = links.insert( link ).first;
2369 faceLinks[ iN ] = & *pLink;
2371 isCurved = !link.IsStraight();
2374 pFace = faces.insert( QFace( faceLinks )).first;
2375 pFace->AddSelfToLinks();
2376 hasRectFaces = ( hasRectFaces || nbN == 4 );
2380 return; // no curved edges of faces
2382 // Compute displacement of medium nodes
2383 // -------------------------------------
2385 // two loops on faces: the first is to treat boundary links, the second is for internal ones
2386 TopLoc_Location loc;
2387 // not treat boundary of volumic submesh
2388 int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2389 for ( ; isInside < 2; ++isInside ) {
2390 MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2391 SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2393 for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2394 if ( bool(isInside) == pFace->IsBoundary() )
2396 for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2399 // make chain of links connected via continues faces
2402 if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2404 if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2406 vector< TChain > chains;
2407 if ( error == ERR_OK ) { // chain contains continues rectangles
2409 chains[0].splice( chains[0].begin(), rawChain );
2411 else if ( error == ERR_TRI ) { // chain contains continues triangles
2412 TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2413 if ( res != _OK ) { // not quadrangles split into triangles
2414 fixTriaNearBoundary( rawChain, *this );
2418 else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2419 fixPrism( rawChain );
2425 for ( int iC = 0; iC < chains.size(); ++iC )
2427 TChain& chain = chains[iC];
2428 if ( chain.empty() ) continue;
2429 if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2433 // mesure chain length and compute link position along the chain
2434 double chainLen = 0;
2435 vector< double > linkPos;
2436 MSGBEG( "Link medium nodes: ");
2437 TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2438 for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2439 MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2440 double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2441 while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2442 link1 = chain.erase( link1 );
2443 if ( link1 == chain.end() )
2445 len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2448 linkPos.push_back( chainLen );
2451 if ( linkPos.size() < 2 )
2454 gp_Vec move0 = chain.front()->_nodeMove;
2455 gp_Vec move1 = chain.back ()->_nodeMove;
2458 bool checkUV = true;
2460 // compute node displacement of end links in parametric space of face
2461 const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2462 TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2463 if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2465 face = TopoDS::Face( f );
2466 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2467 for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2469 TChainLink& link = is1 ? chain.back() : chain.front();
2470 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2471 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2472 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2473 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2474 // uvMove = uvm - uv12
2475 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2476 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2478 if ( move0.SquareMagnitude() < straightTol2 &&
2479 move1.SquareMagnitude() < straightTol2 ) {
2481 continue; // straight - no need to move nodes of internal links
2486 if ( isInside || face.IsNull() )
2488 // compute node displacement of end links in their local coord systems
2490 TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2491 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2492 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2493 move0.Transform(trsf);
2496 TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2497 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2498 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2499 move1.Transform(trsf);
2502 // compute displacement of medium nodes
2503 link2 = chain.begin();
2506 for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2508 double r = linkPos[i] / chainLen;
2509 // displacement in local coord system
2510 gp_Vec move = (1. - r) * move0 + r * move1;
2511 if ( isInside || face.IsNull()) {
2512 // transform to global
2513 gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2514 gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2515 gp_Vec x = x01.Normalized() + x12.Normalized();
2516 trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2517 move.Transform(trsf);
2520 // compute 3D displacement by 2D one
2521 Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
2522 gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2523 gp_XY newUV = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
2524 gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
2525 move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2527 if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2528 move.SquareMagnitude())
2530 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2531 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2532 MSG( "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2533 "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2534 "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2535 "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2539 (*link1)->Move( move );
2540 MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2541 << chain.front()->_mediumNode->GetID() <<"-"
2542 << chain.back ()->_mediumNode->GetID() <<
2543 " by " << move.Magnitude());
2545 } // loop on chains of links
2546 } // loop on 2 directions of propagation from quadrangle
2553 for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2554 if ( pLink->IsMoved() ) {
2555 //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2556 gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2557 GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
2562 //=======================================================================
2564 * \brief Iterator on ancestors of the given type
2566 //=======================================================================
2568 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2570 TopTools_ListIteratorOfListOfShape _ancIter;
2571 TopAbs_ShapeEnum _type;
2572 TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2573 : _ancIter( ancestors ), _type( type )
2575 if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
2579 return _ancIter.More();
2581 virtual const TopoDS_Shape* next()
2583 const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2584 if ( _ancIter.More() )
2585 for ( _ancIter.Next(); _ancIter.More(); _ancIter.Next())
2586 if ( _ancIter.Value().ShapeType() == _type )
2592 //=======================================================================
2594 * \brief Return iterator on ancestors of the given type
2596 //=======================================================================
2598 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2599 const SMESH_Mesh& mesh,
2600 TopAbs_ShapeEnum ancestorType)
2602 return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));