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()
93 TID2ProjectorOnSurf::iterator i_proj = myFace2Projector.begin();
94 for ( ; i_proj != myFace2Projector.end(); ++i_proj )
95 delete i_proj->second;
98 TID2ProjectorOnCurve::iterator i_proj = myEdge2Projector.begin();
99 for ( ; i_proj != myEdge2Projector.end(); ++i_proj )
100 delete i_proj->second;
104 //=======================================================================
105 //function : IsQuadraticSubMesh
106 //purpose : Check submesh for given shape: if all elements on this shape
107 // are quadratic, quadratic elements will be created.
108 // Also fill myTLinkNodeMap
109 //=======================================================================
111 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
113 SMESHDS_Mesh* meshDS = GetMeshDS();
114 // we can create quadratic elements only if all elements
115 // created on subshapes of given shape are quadratic
116 // also we have to fill myTLinkNodeMap
117 myCreateQuadratic = true;
118 mySeamShapeIds.clear();
119 myDegenShapeIds.clear();
120 TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
121 SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
123 int nbOldLinks = myTLinkNodeMap.size();
125 TopExp_Explorer exp( aSh, subType );
126 for (; exp.More() && myCreateQuadratic; exp.Next()) {
127 if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
128 if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
130 const SMDS_MeshElement* e = it->next();
131 if ( e->GetType() != elemType || !e->IsQuadratic() ) {
132 myCreateQuadratic = false;
137 switch ( e->NbNodes() ) {
139 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
141 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
142 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
143 AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
145 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
146 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
147 AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
148 AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
151 myCreateQuadratic = false;
160 if ( nbOldLinks == myTLinkNodeMap.size() )
161 myCreateQuadratic = false;
163 if(!myCreateQuadratic) {
164 myTLinkNodeMap.clear();
168 return myCreateQuadratic;
171 //=======================================================================
172 //function : SetSubShape
173 //purpose : Set geomerty to make elements on
174 //=======================================================================
176 void SMESH_MesherHelper::SetSubShape(const int aShID)
178 if ( aShID == myShapeID )
181 SetSubShape( GetMeshDS()->IndexToShape( aShID ));
183 SetSubShape( TopoDS_Shape() );
186 //=======================================================================
187 //function : SetSubShape
188 //purpose : Set geomerty to create elements on
189 //=======================================================================
191 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
193 if ( myShape.IsSame( aSh ))
197 mySeamShapeIds.clear();
198 myDegenShapeIds.clear();
200 if ( myShape.IsNull() ) {
204 SMESHDS_Mesh* meshDS = GetMeshDS();
205 myShapeID = meshDS->ShapeToIndex(aSh);
208 // treatment of periodic faces
209 for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
211 const TopoDS_Face& face = TopoDS::Face( eF.Current() );
212 BRepAdaptor_Surface surface( face );
213 if ( surface.IsUPeriodic() || surface.IsVPeriodic() )
215 for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
217 // look for a seam edge
218 const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
219 if ( BRep_Tool::IsClosed( edge, face )) {
220 // initialize myPar1, myPar2 and myParIndex
222 BRep_Tool::UVPoints( edge, face, uv1, uv2 );
223 if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
225 myParIndex |= U_periodic;
226 myPar1[0] = surface.FirstUParameter();
227 myPar2[0] = surface.LastUParameter();
230 myParIndex |= V_periodic;
231 myPar1[1] = surface.FirstVParameter();
232 myPar2[1] = surface.LastVParameter();
234 // store seam shape indices, negative if shape encounters twice
235 int edgeID = meshDS->ShapeToIndex( edge );
236 mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
237 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
238 int vertexID = meshDS->ShapeToIndex( v.Current() );
239 mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
243 // look for a degenerated edge
244 if ( BRep_Tool::Degenerated( edge )) {
245 myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
246 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
247 myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
254 //=======================================================================
255 //function : GetNodeUVneedInFaceNode
256 //purpose : Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
257 // Return true if the face is periodic.
258 // If F is Null, answer about subshape set through IsQuadraticSubMesh() or
260 //=======================================================================
262 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
264 if ( F.IsNull() ) return !mySeamShapeIds.empty();
266 if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
267 return !mySeamShapeIds.empty();
270 Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
271 if ( !aSurface.IsNull() )
272 return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
277 //=======================================================================
278 //function : IsMedium
280 //=======================================================================
282 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode* node,
283 const SMDSAbs_ElementType typeToCheck)
285 return SMESH_MeshEditor::IsMedium( node, typeToCheck );
288 //=======================================================================
289 //function : GetSubShapeByNode
290 //purpose : Return support shape of a node
291 //=======================================================================
293 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
294 SMESHDS_Mesh* meshDS)
296 int shapeID = node->getshapeId();
297 if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
298 return meshDS->IndexToShape( shapeID );
300 return TopoDS_Shape();
304 //=======================================================================
305 //function : AddTLinkNode
306 //purpose : add a link in my data structure
307 //=======================================================================
309 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
310 const SMDS_MeshNode* n2,
311 const SMDS_MeshNode* n12)
313 // add new record to map
314 SMESH_TLink link( n1, n2 );
315 myTLinkNodeMap.insert( make_pair(link,n12));
318 //================================================================================
320 * \brief Return true if position of nodes on the shape hasn't yet been checked or
321 * the positions proved to be invalid
323 //================================================================================
325 bool SMESH_MesherHelper::toCheckPosOnShape(int shapeID ) const
327 map< int,bool >::const_iterator id_ok = myNodePosShapesValidity.find( shapeID );
328 return ( id_ok == myNodePosShapesValidity.end() || !id_ok->second );
331 //================================================================================
333 * \brief Set validity of positions of nodes on the shape.
334 * Once set, validity is not changed
336 //================================================================================
338 void SMESH_MesherHelper::setPosOnShapeValidity(int shapeID, bool ok ) const
340 ((SMESH_MesherHelper*)this)->myNodePosShapesValidity.insert( make_pair( shapeID, ok));
343 //=======================================================================
344 //function : GetUVOnSeam
345 //purpose : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
346 //=======================================================================
348 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
350 gp_Pnt2d result = uv1;
351 for ( int i = U_periodic; i <= V_periodic ; ++i )
353 if ( myParIndex & i )
355 double p1 = uv1.Coord( i );
356 double dp1 = Abs( p1-myPar1[i-1]), dp2 = Abs( p1-myPar2[i-1]);
357 if ( myParIndex == i ||
358 dp1 < ( myPar2[i-1] - myPar2[i-1] ) / 100. ||
359 dp2 < ( myPar2[i-1] - myPar2[i-1] ) / 100. )
361 double p2 = uv2.Coord( i );
362 double p1Alt = ( dp1 < dp2 ) ? myPar2[i-1] : myPar1[i-1];
363 if ( Abs( p2 - p1 ) > Abs( p2 - p1Alt ))
364 result.SetCoord( i, p1Alt );
371 //=======================================================================
372 //function : GetNodeUV
373 //purpose : Return node UV on face
374 //=======================================================================
376 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
377 const SMDS_MeshNode* n,
378 const SMDS_MeshNode* n2,
381 gp_Pnt2d uv( Precision::Infinite(), Precision::Infinite() );
382 const SMDS_PositionPtr Pos = n->GetPosition();
384 if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
386 // node has position on face
387 const SMDS_FacePosition* fpos =
388 static_cast<const SMDS_FacePosition*>(n->GetPosition());
389 uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
391 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ));
393 else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
395 // node has position on edge => it is needed to find
396 // corresponding edge from face, get pcurve for this
397 // edge and retrieve value from this pcurve
398 const SMDS_EdgePosition* epos =
399 static_cast<const SMDS_EdgePosition*>(n->GetPosition());
400 int edgeID = n->getshapeId();
401 TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
402 double f, l, u = epos->GetUParameter();
403 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
404 bool validU = ( f < u && u < l );
406 uv = C2d->Value( u );
409 if ( check || !validU )
410 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ),/*force=*/ !validU );
412 // for a node on a seam edge select one of UVs on 2 pcurves
413 if ( n2 && IsSeamShape( edgeID ) )
415 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0, check ));
418 { // adjust uv to period
420 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
421 Standard_Boolean isUPeriodic = S->IsUPeriodic();
422 Standard_Boolean isVPeriodic = S->IsVPeriodic();
423 if ( isUPeriodic || isVPeriodic ) {
424 Standard_Real UF,UL,VF,VL;
425 S->Bounds(UF,UL,VF,VL);
427 uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
429 uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
433 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
435 if ( int vertexID = n->getshapeId() ) {
436 const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
438 uv = BRep_Tool::Parameters( V, F );
441 catch (Standard_Failure& exc) {
444 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() )
445 uvOK = ( V == vert.Current() );
448 MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
449 << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
451 // get UV of a vertex closest to the node
453 gp_Pnt pn = XYZ( n );
454 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) {
455 TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
456 gp_Pnt p = BRep_Tool::Pnt( curV );
457 double curDist = p.SquareDistance( pn );
458 if ( curDist < dist ) {
460 uv = BRep_Tool::Parameters( curV, F );
461 uvOK = ( dist < DBL_MIN );
467 TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
468 for ( ; it.More(); it.Next() ) {
469 if ( it.Value().ShapeType() == TopAbs_EDGE ) {
470 const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
472 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
473 if ( !C2d.IsNull() ) {
474 double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
475 uv = C2d->Value( u );
483 if ( n2 && IsSeamShape( vertexID ) )
484 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
494 //=======================================================================
495 //function : CheckNodeUV
496 //purpose : Check and fix node UV on a face
497 //=======================================================================
499 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
500 const SMDS_MeshNode* n,
503 const bool force) const
505 int shapeID = n->getshapeId();
506 if ( force || toCheckPosOnShape( shapeID ))
509 double tolmin = 1.e-7*myMesh->GetMeshDS()->getMaxDim(); // nodes coordinates are stored in float format
510 if (toldis < tolmin) toldis = tolmin;
511 // check that uv is correct
513 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
514 gp_Pnt nodePnt = XYZ( n );
515 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
516 if ( Precision::IsInfinite( uv.X() ) ||
517 Precision::IsInfinite( uv.Y() ) ||
518 nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > toldis )
520 setPosOnShapeValidity( shapeID, false );
521 // uv incorrect, project the node to surface
522 GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
523 projector.Perform( nodePnt );
524 if ( !projector.IsDone() || projector.NbPoints() < 1 )
526 MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
529 Quantity_Parameter U,V;
530 projector.LowerDistanceParameters(U,V);
532 if ( nodePnt.Distance( surface->Value( U, V )) > toldis )
534 MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
537 // store the fixed UV on the face
538 if ( myShape.IsSame(F) && shapeID == myShapeID )
539 const_cast<SMDS_MeshNode*>(n)->SetPosition
540 ( SMDS_PositionPtr( new SMDS_FacePosition( U, V )));
542 else if ( uv.Modulus() > numeric_limits<double>::min() )
544 setPosOnShapeValidity( shapeID, true );
550 //=======================================================================
551 //function : GetProjector
552 //purpose : Return projector intitialized by given face without location, which is returned
553 //=======================================================================
555 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
556 TopLoc_Location& loc,
559 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
560 int faceID = GetMeshDS()->ShapeToIndex( F );
561 TID2ProjectorOnSurf& i2proj = const_cast< TID2ProjectorOnSurf&>( myFace2Projector );
562 TID2ProjectorOnSurf::iterator i_proj = i2proj.find( faceID );
563 if ( i_proj == i2proj.end() )
565 if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
566 double U1, U2, V1, V2;
567 surface->Bounds(U1, U2, V1, V2);
568 GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
569 proj->Init( surface, U1, U2, V1, V2, tol );
570 i_proj = i2proj.insert( make_pair( faceID, proj )).first;
572 return *( i_proj->second );
577 gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
578 gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
579 gp_XY_FunPtr(Subtracted);
582 //=======================================================================
583 //function : applyIn2D
584 //purpose : Perform given operation on two 2d points in parameric space of given surface.
585 // It takes into account period of the surface. Use gp_XY_FunPtr macro
586 // to easily define pointer to function of gp_XY class.
587 //=======================================================================
589 gp_XY SMESH_MesherHelper::applyIn2D(const Handle(Geom_Surface)& surface,
593 const bool resultInPeriod)
595 Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
596 Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
597 if ( !isUPeriodic && !isVPeriodic )
600 // move uv2 not far than half-period from uv1
602 uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
604 uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
607 gp_XY res = fun( uv1, gp_XY(u2,v2) );
609 // move result within period
610 if ( resultInPeriod )
612 Standard_Real UF,UL,VF,VL;
613 surface->Bounds(UF,UL,VF,VL);
615 res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
617 res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
622 //=======================================================================
623 //function : GetMiddleUV
624 //purpose : Return middle UV taking in account surface period
625 //=======================================================================
627 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
631 return applyIn2D( surface, p1, p2, & AverageUV );
634 //=======================================================================
635 //function : GetNodeU
636 //purpose : Return node U on edge
637 //=======================================================================
639 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
640 const SMDS_MeshNode* n,
641 const SMDS_MeshNode* inEdgeNode,
645 const SMDS_PositionPtr pos = n->GetPosition();
646 if ( pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
648 const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( pos );
649 param = epos->GetUParameter();
651 else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
653 if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
656 BRep_Tool::Range( E, f,l );
657 double uInEdge = GetNodeU( E, inEdgeNode );
658 param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
662 SMESHDS_Mesh * meshDS = GetMeshDS();
663 int vertexID = n->getshapeId();
664 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
665 param = BRep_Tool::Parameter( V, E );
670 double tol = BRep_Tool::Tolerance( E );
671 double f,l; BRep_Tool::Range( E, f,l );
672 bool force = ( param < f-tol || param > l+tol );
673 if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
674 force = ( GetMeshDS()->ShapeToIndex( E ) != n->getshapeId() );
676 *check = CheckNodeU( E, n, param, 2*tol, force );
681 //=======================================================================
682 //function : CheckNodeU
683 //purpose : Check and fix node U on an edge
684 // Return false if U is bad and could not be fixed
685 //=======================================================================
687 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
688 const SMDS_MeshNode* n,
692 double* distance) const
694 int shapeID = n->getshapeId();
695 if ( force || toCheckPosOnShape( shapeID ))
698 double tolmin = 1.e-7*myMesh->GetMeshDS()->getMaxDim(); // nodes coordinates are stored in float format
699 if (toldis < tolmin) toldis = tolmin;
700 // check that u is correct
701 TopLoc_Location loc; double f,l;
702 Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
703 if ( curve.IsNull() ) // degenerated edge
705 if ( u+tol < f || u-tol > l )
707 double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
713 gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
714 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
715 double dist = nodePnt.Distance( curve->Value( u ));
716 if ( distance ) *distance = dist;
719 setPosOnShapeValidity( shapeID, false );
720 // u incorrect, project the node to the curve
721 int edgeID = GetMeshDS()->ShapeToIndex( E );
722 TID2ProjectorOnCurve& i2proj = const_cast< TID2ProjectorOnCurve&>( myEdge2Projector );
723 TID2ProjectorOnCurve::iterator i_proj =
724 i2proj.insert( make_pair( edgeID, (GeomAPI_ProjectPointOnCurve*) 0 )).first;
725 if ( !i_proj->second )
727 i_proj->second = new GeomAPI_ProjectPointOnCurve();
728 i_proj->second->Init( curve, f, l );
730 GeomAPI_ProjectPointOnCurve* projector = i_proj->second;
731 projector->Perform( nodePnt );
732 if ( projector->NbPoints() < 1 )
734 MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
737 Quantity_Parameter U = projector->LowerDistanceParameter();
739 dist = nodePnt.Distance( curve->Value( U ));
740 if ( distance ) *distance = dist;
743 MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
744 MESSAGE("distance " << nodePnt.Distance(curve->Value( U )) << " " << toldis);
747 // store the fixed U on the edge
748 if ( myShape.IsSame(E) && shapeID == myShapeID )
749 const_cast<SMDS_MeshNode*>(n)->SetPosition
750 ( SMDS_PositionPtr( new SMDS_EdgePosition( U )));
752 else if ( fabs( u ) > numeric_limits<double>::min() )
754 setPosOnShapeValidity( shapeID, true );
756 if (( u < f-tol || u > l+tol ) && force )
758 // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
761 // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
762 double period = curve->Period();
763 u = ( u < f ) ? u + period : u - period;
765 catch (Standard_Failure& exc)
775 //=======================================================================
776 //function : GetMediumNode
777 //purpose : Return existing or create new medium nodes between given ones
778 //=======================================================================
780 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
781 const SMDS_MeshNode* n2,
784 // Find existing node
786 SMESH_TLink link(n1,n2);
787 ItTLinkNode itLN = myTLinkNodeMap.find( link );
788 if ( itLN != myTLinkNodeMap.end() ) {
789 return (*itLN).second;
792 // Create medium node
795 SMESHDS_Mesh* meshDS = GetMeshDS();
797 if ( IsSeamShape( n1->getshapeId() ))
798 // to get a correct UV of a node on seam, the second node must have checked UV
801 // get type of shape for the new medium node
802 int faceID = -1, edgeID = -1;
803 const SMDS_PositionPtr Pos1 = n1->GetPosition();
804 const SMDS_PositionPtr Pos2 = n2->GetPosition();
806 TopoDS_Edge E; double u [2];
807 TopoDS_Face F; gp_XY uv[2];
808 bool uvOK[2] = { false, false };
810 if( myShape.IsNull() )
812 if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
813 faceID = n1->getshapeId();
815 else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
816 faceID = n2->getshapeId();
819 if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
820 edgeID = n1->getshapeId();
822 if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
823 edgeID = n2->getshapeId();
826 // get positions of the given nodes on shapes
827 TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
828 if ( faceID>0 || shapeType == TopAbs_FACE)
830 if( myShape.IsNull() )
831 F = TopoDS::Face(meshDS->IndexToShape(faceID));
833 F = TopoDS::Face(myShape);
836 uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
837 uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
839 else if (edgeID>0 || shapeType == TopAbs_EDGE)
841 if ( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE &&
842 Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE &&
843 n1->getshapeId() != n2->getshapeId() ) // issue 0021006
844 return getMediumNodeOnComposedWire(n1,n2,force3d);
846 if( myShape.IsNull() )
847 E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
849 E = TopoDS::Edge(myShape);
852 u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
853 u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
857 // we try to create medium node using UV parameters of
858 // nodes, else - medium between corresponding 3d points
861 if ( uvOK[0] && uvOK[1] )
863 if ( IsDegenShape( n1->getshapeId() )) {
864 if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
865 else uv[0].SetCoord( 2, uv[1].Coord( 2 ));
867 else if ( IsDegenShape( n2->getshapeId() )) {
868 if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
869 else uv[1].SetCoord( 2, uv[0].Coord( 2 ));
873 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
874 gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
875 gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
876 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
877 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
878 myTLinkNodeMap.insert(make_pair(link,n12));
882 else if ( !E.IsNull() )
885 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
888 Standard_Boolean isPeriodic = C->IsPeriodic();
891 Standard_Real Period = C->Period();
892 Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
893 Standard_Real pmid = (u[0]+p)/2.;
894 U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
899 gp_Pnt P = C->Value( U );
900 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
901 meshDS->SetNodeOnEdge(n12, edgeID, U);
902 myTLinkNodeMap.insert(make_pair(link,n12));
909 double x = ( n1->X() + n2->X() )/2.;
910 double y = ( n1->Y() + n2->Y() )/2.;
911 double z = ( n1->Z() + n2->Z() )/2.;
912 n12 = meshDS->AddNode(x,y,z);
916 gp_XY UV = ( uv[0] + uv[1] ) / 2.;
917 CheckNodeUV( F, n12, UV, 2*BRep_Tool::Tolerance( F ), /*force=*/true);
918 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
920 else if ( !E.IsNull() )
922 double U = ( u[0] + u[1] ) / 2.;
923 CheckNodeU( E, n12, U, 2*BRep_Tool::Tolerance( E ), /*force=*/true);
924 meshDS->SetNodeOnEdge(n12, edgeID, U);
926 else if ( myShapeID > 0 )
928 meshDS->SetNodeInVolume(n12, myShapeID);
931 myTLinkNodeMap.insert( make_pair( link, n12 ));
935 //================================================================================
937 * \brief Makes a medium node if nodes reside different edges
939 //================================================================================
941 const SMDS_MeshNode* SMESH_MesherHelper::getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
942 const SMDS_MeshNode* n2,
945 gp_Pnt middle = 0.5 * XYZ(n1) + 0.5 * XYZ(n2);
946 SMDS_MeshNode* n12 = AddNode( middle.X(), middle.Y(), middle.Z() );
948 // To find position on edge and 3D position for n12,
949 // project <middle> to 2 edges and select projection most close to <middle>
951 double u = 0, distMiddleProj = Precision::Infinite();
953 TopoDS_Edge edges[2];
954 for ( int is2nd = 0; is2nd < 2; ++is2nd )
957 const SMDS_MeshNode* n = is2nd ? n2 : n1;
958 TopoDS_Shape shape = GetSubShapeByNode( n, GetMeshDS() );
959 if ( shape.IsNull() || shape.ShapeType() != TopAbs_EDGE )
962 // project to get U of projection and distance from middle to projection
963 TopoDS_Edge edge = edges[ is2nd ] = TopoDS::Edge( shape );
964 double node2MiddleDist = middle.Distance( XYZ(n) );
965 double foundU = GetNodeU( edge, n ), foundDist = node2MiddleDist;
966 CheckNodeU( edge, n12, foundU, 2*BRep_Tool::Tolerance(edge), /*force=*/true, &foundDist );
967 if ( foundDist < node2MiddleDist )
969 distMiddleProj = foundDist;
974 if ( Precision::IsInfinite( distMiddleProj ))
976 // both projections failed; set n12 on the edge of n1 with U of a common vertex
977 TopoDS_Vertex vCommon;
978 if ( TopExp::CommonVertex( edges[0], edges[1], vCommon ))
979 u = BRep_Tool::Parameter( vCommon, edges[0] );
982 double f,l, u0 = GetNodeU( edges[0], n1 );
983 BRep_Tool::Range( edges[0],f,l );
984 u = ( fabs(u0-f) < fabs(u0-l) ) ? f : l;
990 // move n12 to position of a successfull projection
991 double tol = BRep_Tool::Tolerance(edges[ iOkEdge ]);
992 if ( !force3d && distMiddleProj > 2*tol )
994 TopLoc_Location loc; double f,l;
995 Handle(Geom_Curve) curve = BRep_Tool::Curve( edges[iOkEdge],loc,f,l );
996 gp_Pnt p = curve->Value( u );
997 GetMeshDS()->MoveNode( n12, p.X(), p.Y(), p.Z() );
1000 GetMeshDS()->SetNodeOnEdge(n12, edges[iOkEdge], u);
1002 myTLinkNodeMap.insert( make_pair( SMESH_TLink(n1,n2), n12 ));
1007 //=======================================================================
1008 //function : AddNode
1009 //purpose : Creates a node
1010 //=======================================================================
1012 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
1014 SMESHDS_Mesh * meshDS = GetMeshDS();
1015 SMDS_MeshNode* node = 0;
1017 node = meshDS->AddNodeWithID( x, y, z, ID );
1019 node = meshDS->AddNode( x, y, z );
1020 if ( mySetElemOnShape && myShapeID > 0 ) {
1021 switch ( myShape.ShapeType() ) {
1022 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
1023 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
1024 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID); break;
1025 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID); break;
1026 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
1033 //=======================================================================
1034 //function : AddEdge
1035 //purpose : Creates quadratic or linear edge
1036 //=======================================================================
1038 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
1039 const SMDS_MeshNode* n2,
1043 SMESHDS_Mesh * meshDS = GetMeshDS();
1045 SMDS_MeshEdge* edge = 0;
1046 if (myCreateQuadratic) {
1047 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1049 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
1051 edge = meshDS->AddEdge(n1, n2, n12);
1055 edge = meshDS->AddEdgeWithID(n1, n2, id);
1057 edge = meshDS->AddEdge(n1, n2);
1060 if ( mySetElemOnShape && myShapeID > 0 )
1061 meshDS->SetMeshElementOnShape( edge, myShapeID );
1066 //=======================================================================
1067 //function : AddFace
1068 //purpose : Creates quadratic or linear triangle
1069 //=======================================================================
1071 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1072 const SMDS_MeshNode* n2,
1073 const SMDS_MeshNode* n3,
1077 SMESHDS_Mesh * meshDS = GetMeshDS();
1078 SMDS_MeshFace* elem = 0;
1080 if( n1==n2 || n2==n3 || n3==n1 )
1083 if(!myCreateQuadratic) {
1085 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
1087 elem = meshDS->AddFace(n1, n2, n3);
1090 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1091 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1092 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1095 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
1097 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
1099 if ( mySetElemOnShape && myShapeID > 0 )
1100 meshDS->SetMeshElementOnShape( elem, myShapeID );
1105 //=======================================================================
1106 //function : AddFace
1107 //purpose : Creates quadratic or linear quadrangle
1108 //=======================================================================
1110 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1111 const SMDS_MeshNode* n2,
1112 const SMDS_MeshNode* n3,
1113 const SMDS_MeshNode* n4,
1117 SMESHDS_Mesh * meshDS = GetMeshDS();
1118 SMDS_MeshFace* elem = 0;
1121 return AddFace(n1,n3,n4,id,force3d);
1124 return AddFace(n1,n2,n4,id,force3d);
1127 return AddFace(n1,n2,n3,id,force3d);
1130 return AddFace(n1,n2,n4,id,force3d);
1133 return AddFace(n1,n2,n3,id,force3d);
1136 return AddFace(n1,n2,n3,id,force3d);
1139 if(!myCreateQuadratic) {
1141 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
1143 elem = meshDS->AddFace(n1, n2, n3, n4);
1146 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1147 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1148 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1149 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1152 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
1154 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
1156 if ( mySetElemOnShape && myShapeID > 0 )
1157 meshDS->SetMeshElementOnShape( elem, myShapeID );
1162 //=======================================================================
1163 //function : AddPolygonalFace
1164 //purpose : Creates polygon, with additional nodes in quadratic mesh
1165 //=======================================================================
1167 SMDS_MeshFace* SMESH_MesherHelper::AddPolygonalFace (const vector<const SMDS_MeshNode*>& nodes,
1171 SMESHDS_Mesh * meshDS = GetMeshDS();
1172 SMDS_MeshFace* elem = 0;
1174 if(!myCreateQuadratic) {
1176 elem = meshDS->AddPolygonalFaceWithID(nodes, id);
1178 elem = meshDS->AddPolygonalFace(nodes);
1181 vector<const SMDS_MeshNode*> newNodes;
1182 for ( int i = 0; i < nodes.size(); ++i )
1184 const SMDS_MeshNode* n1 = nodes[i];
1185 const SMDS_MeshNode* n2 = nodes[(i+1)/nodes.size()];
1186 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1187 newNodes.push_back( n1 );
1188 newNodes.push_back( n12 );
1191 elem = meshDS->AddPolygonalFaceWithID(newNodes, id);
1193 elem = meshDS->AddPolygonalFace(newNodes);
1195 if ( mySetElemOnShape && myShapeID > 0 )
1196 meshDS->SetMeshElementOnShape( elem, myShapeID );
1201 //=======================================================================
1202 //function : AddVolume
1203 //purpose : Creates quadratic or linear prism
1204 //=======================================================================
1206 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1207 const SMDS_MeshNode* n2,
1208 const SMDS_MeshNode* n3,
1209 const SMDS_MeshNode* n4,
1210 const SMDS_MeshNode* n5,
1211 const SMDS_MeshNode* n6,
1215 SMESHDS_Mesh * meshDS = GetMeshDS();
1216 SMDS_MeshVolume* elem = 0;
1217 if(!myCreateQuadratic) {
1219 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
1221 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
1224 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1225 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1226 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1228 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1229 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1230 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1232 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1233 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1234 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1237 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
1238 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1240 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1241 n12, n23, n31, n45, n56, n64, n14, n25, n36);
1243 if ( mySetElemOnShape && myShapeID > 0 )
1244 meshDS->SetMeshElementOnShape( elem, myShapeID );
1249 //=======================================================================
1250 //function : AddVolume
1251 //purpose : Creates quadratic or linear tetrahedron
1252 //=======================================================================
1254 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1255 const SMDS_MeshNode* n2,
1256 const SMDS_MeshNode* n3,
1257 const SMDS_MeshNode* n4,
1261 SMESHDS_Mesh * meshDS = GetMeshDS();
1262 SMDS_MeshVolume* elem = 0;
1263 if(!myCreateQuadratic) {
1265 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1267 elem = meshDS->AddVolume(n1, n2, n3, n4);
1270 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1271 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1272 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1274 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1275 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1276 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1279 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1281 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1283 if ( mySetElemOnShape && myShapeID > 0 )
1284 meshDS->SetMeshElementOnShape( elem, myShapeID );
1289 //=======================================================================
1290 //function : AddVolume
1291 //purpose : Creates quadratic or linear pyramid
1292 //=======================================================================
1294 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1295 const SMDS_MeshNode* n2,
1296 const SMDS_MeshNode* n3,
1297 const SMDS_MeshNode* n4,
1298 const SMDS_MeshNode* n5,
1302 SMDS_MeshVolume* elem = 0;
1303 if(!myCreateQuadratic) {
1305 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1307 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1310 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1311 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1312 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1313 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1315 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1316 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1317 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1318 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1321 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
1326 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
1328 n15, n25, n35, n45);
1330 if ( mySetElemOnShape && myShapeID > 0 )
1331 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1336 //=======================================================================
1337 //function : AddVolume
1338 //purpose : Creates quadratic or linear hexahedron
1339 //=======================================================================
1341 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1342 const SMDS_MeshNode* n2,
1343 const SMDS_MeshNode* n3,
1344 const SMDS_MeshNode* n4,
1345 const SMDS_MeshNode* n5,
1346 const SMDS_MeshNode* n6,
1347 const SMDS_MeshNode* n7,
1348 const SMDS_MeshNode* n8,
1352 SMESHDS_Mesh * meshDS = GetMeshDS();
1353 SMDS_MeshVolume* elem = 0;
1354 if(!myCreateQuadratic) {
1356 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1358 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1361 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1362 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1363 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1364 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1366 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1367 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1368 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1369 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1371 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1372 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1373 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1374 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1377 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1378 n12, n23, n34, n41, n56, n67,
1379 n78, n85, n15, n26, n37, n48, id);
1381 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1382 n12, n23, n34, n41, n56, n67,
1383 n78, n85, n15, n26, n37, n48);
1385 if ( mySetElemOnShape && myShapeID > 0 )
1386 meshDS->SetMeshElementOnShape( elem, myShapeID );
1391 //=======================================================================
1392 //function : AddPolyhedralVolume
1393 //purpose : Creates polyhedron. In quadratic mesh, adds medium nodes
1394 //=======================================================================
1397 SMESH_MesherHelper::AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
1398 const std::vector<int>& quantities,
1402 SMESHDS_Mesh * meshDS = GetMeshDS();
1403 SMDS_MeshVolume* elem = 0;
1404 if(!myCreateQuadratic)
1407 elem = meshDS->AddPolyhedralVolumeWithID(nodes, quantities, id);
1409 elem = meshDS->AddPolyhedralVolume(nodes, quantities);
1413 vector<const SMDS_MeshNode*> newNodes;
1414 vector<int> newQuantities;
1415 for ( int iFace=0, iN=0; iFace < quantities.size(); ++iFace)
1417 int nbNodesInFace = quantities[iFace];
1418 newQuantities.push_back(0);
1419 for ( int i = 0; i < nbNodesInFace; ++i )
1421 const SMDS_MeshNode* n1 = nodes[ iN + i ];
1422 newNodes.push_back( n1 );
1423 newQuantities.back()++;
1425 const SMDS_MeshNode* n2 = nodes[ iN + ( i+1==nbNodesInFace ? 0 : i+1 )];
1426 // if ( n1->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE &&
1427 // n2->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE )
1429 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1430 newNodes.push_back( n12 );
1431 newQuantities.back()++;
1434 iN += nbNodesInFace;
1437 elem = meshDS->AddPolyhedralVolumeWithID( newNodes, newQuantities, id );
1439 elem = meshDS->AddPolyhedralVolume( newNodes, newQuantities );
1441 if ( mySetElemOnShape && myShapeID > 0 )
1442 meshDS->SetMeshElementOnShape( elem, myShapeID );
1447 //=======================================================================
1448 //function : LoadNodeColumns
1449 //purpose : Load nodes bound to face into a map of node columns
1450 //=======================================================================
1452 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1453 const TopoDS_Face& theFace,
1454 const TopoDS_Edge& theBaseEdge,
1455 SMESHDS_Mesh* theMesh)
1457 SMESHDS_SubMesh* faceSubMesh = theMesh->MeshElements( theFace );
1458 if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
1461 // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
1463 map< double, const SMDS_MeshNode*> sortedBaseNodes;
1464 if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
1465 || sortedBaseNodes.size() < 2 )
1468 int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
1469 map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
1470 double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
1471 for ( ; u_n != sortedBaseNodes.end(); u_n++ )
1473 double par = ( u_n->first - f ) / range;
1474 vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
1475 nCol.resize( nbRows );
1476 nCol[0] = u_n->second;
1479 // fill theParam2ColumnMap column by column by passing from nodes on
1480 // theBaseEdge up via mesh faces on theFace
1482 TParam2ColumnMap::iterator par_nVec_2 = theParam2ColumnMap.begin();
1483 TParam2ColumnMap::iterator par_nVec_1 = par_nVec_2++;
1484 TIDSortedElemSet emptySet, avoidSet;
1485 for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
1487 vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
1488 vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
1490 int i1, i2, iRow = 0;
1491 const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
1492 // find face sharing node n1 and n2 and belonging to faceSubMesh
1493 while ( const SMDS_MeshElement* face =
1494 SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
1496 if ( faceSubMesh->Contains( face ))
1498 int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
1501 n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
1502 n2 = face->GetNode( (i1+2) % 4 );
1503 if ( ++iRow >= nbRows )
1509 avoidSet.insert( face );
1511 if ( iRow + 1 < nbRows ) // compact if necessary
1512 nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
1517 //=======================================================================
1518 //function : NbAncestors
1519 //purpose : Return number of unique ancestors of the shape
1520 //=======================================================================
1522 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1523 const SMESH_Mesh& mesh,
1524 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
1526 TopTools_MapOfShape ancestors;
1527 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1528 for ( ; ansIt.More(); ansIt.Next() ) {
1529 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1530 ancestors.Add( ansIt.Value() );
1532 return ancestors.Extent();
1535 //=======================================================================
1536 //function : GetSubShapeOri
1537 //purpose : Return orientation of sub-shape in the main shape
1538 //=======================================================================
1540 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
1541 const TopoDS_Shape& subShape)
1543 TopAbs_Orientation ori = TopAbs_Orientation(-1);
1544 if ( !shape.IsNull() && !subShape.IsNull() )
1546 TopExp_Explorer e( shape, subShape.ShapeType() );
1547 if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
1548 e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
1549 for ( ; e.More(); e.Next())
1550 if ( subShape.IsSame( e.Current() ))
1553 ori = e.Current().Orientation();
1558 //=======================================================================
1559 //function : IsSubShape
1561 //=======================================================================
1563 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
1564 const TopoDS_Shape& mainShape )
1566 if ( !shape.IsNull() && !mainShape.IsNull() )
1568 for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
1571 if ( shape.IsSame( exp.Current() ))
1574 SCRUTE((shape.IsNull()));
1575 SCRUTE((mainShape.IsNull()));
1579 //=======================================================================
1580 //function : IsSubShape
1582 //=======================================================================
1584 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
1586 if ( shape.IsNull() || !aMesh )
1589 aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
1591 (shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape ));
1594 //================================================================================
1596 * \brief Return maximal tolerance of shape
1598 //================================================================================
1600 double SMESH_MesherHelper::MaxTolerance( const TopoDS_Shape& shape )
1602 double tol = Precision::Confusion();
1603 TopExp_Explorer exp;
1604 for ( exp.Init( shape, TopAbs_FACE ); exp.More(); exp.Next() )
1605 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Face( exp.Current())));
1606 for ( exp.Init( shape, TopAbs_EDGE ); exp.More(); exp.Next() )
1607 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Edge( exp.Current())));
1608 for ( exp.Init( shape, TopAbs_VERTEX ); exp.More(); exp.Next() )
1609 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Vertex( exp.Current())));
1614 //=======================================================================
1615 //function : IsQuadraticMesh
1616 //purpose : Check mesh without geometry for: if all elements on this shape are quadratic,
1617 // quadratic elements will be created.
1618 // Used then generated 3D mesh without geometry.
1619 //=======================================================================
1621 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1623 int NbAllEdgsAndFaces=0;
1624 int NbQuadFacesAndEdgs=0;
1625 int NbFacesAndEdges=0;
1626 //All faces and edges
1627 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1629 //Quadratic faces and edges
1630 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1632 //Linear faces and edges
1633 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1635 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1637 return SMESH_MesherHelper::QUADRATIC;
1639 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1641 return SMESH_MesherHelper::LINEAR;
1644 //Mesh with both type of elements
1645 return SMESH_MesherHelper::COMP;
1648 //=======================================================================
1649 //function : GetOtherParam
1650 //purpose : Return an alternative parameter for a node on seam
1651 //=======================================================================
1653 double SMESH_MesherHelper::GetOtherParam(const double param) const
1655 int i = myParIndex & U_periodic ? 0 : 1;
1656 return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1659 //#include <Perf_Meter.hxx>
1661 //=======================================================================
1662 namespace { // Structures used by FixQuadraticElements()
1663 //=======================================================================
1665 #define __DMP__(txt) \
1667 #define MSG(txt) __DMP__(txt<<endl)
1668 #define MSGBEG(txt) __DMP__(txt)
1670 //const double straightTol2 = 1e-33; // to detect straing links
1671 bool isStraightLink(double linkLen2, double middleNodeMove2)
1673 // straight if <node move> < 1/15 * <link length>
1674 return middleNodeMove2 < 1/15./15. * linkLen2;
1678 // ---------------------------------------
1680 * \brief Quadratic link knowing its faces
1682 struct QLink: public SMESH_TLink
1684 const SMDS_MeshNode* _mediumNode;
1685 mutable vector<const QFace* > _faces;
1686 mutable gp_Vec _nodeMove;
1687 mutable int _nbMoves;
1689 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1690 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1692 //if ( MediumPos() != SMDS_TOP_3DSPACE )
1693 _nodeMove = MediumPnt() - MiddlePnt();
1695 void SetContinuesFaces() const;
1696 const QFace* GetContinuesFace( const QFace* face ) const;
1697 bool OnBoundary() const;
1698 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1699 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1701 SMDS_TypeOfPosition MediumPos() const
1702 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1703 SMDS_TypeOfPosition EndPos(bool isSecond) const
1704 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1705 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1706 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1708 void Move(const gp_Vec& move, bool sum=false) const
1709 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1710 gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1711 bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1712 bool IsStraight() const
1713 { return isStraightLink( (XYZ(node1())-XYZ(node2())).SquareModulus(),
1714 _nodeMove.SquareMagnitude());
1716 bool operator<(const QLink& other) const {
1717 return (node1()->GetID() == other.node1()->GetID() ?
1718 node2()->GetID() < other.node2()->GetID() :
1719 node1()->GetID() < other.node1()->GetID());
1721 // struct PtrComparator {
1722 // bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1725 // ---------------------------------------------------------
1727 * \brief Link in the chain of links; it connects two faces
1731 const QLink* _qlink;
1732 mutable const QFace* _qfaces[2];
1734 TChainLink(const QLink* qlink=0):_qlink(qlink) {
1735 _qfaces[0] = _qfaces[1] = 0;
1737 void SetFace(const QFace* face) const { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1739 bool IsBoundary() const { return !_qfaces[1]; }
1741 void RemoveFace( const QFace* face ) const
1742 { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1744 const QFace* NextFace( const QFace* f ) const
1745 { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1747 const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1748 { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1750 bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1752 operator bool() const { return (_qlink); }
1754 const QLink* operator->() const { return _qlink; }
1756 gp_Vec Normal() const;
1758 // --------------------------------------------------------------------
1759 typedef list< TChainLink > TChain;
1760 typedef set < TChainLink > TLinkSet;
1761 typedef TLinkSet::const_iterator TLinkInSet;
1763 const int theFirstStep = 5;
1765 enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1766 // --------------------------------------------------------------------
1768 * \brief Face shared by two volumes and bound by QLinks
1770 struct QFace: public TIDSortedNodeSet
1772 mutable const SMDS_MeshElement* _volumes[2];
1773 mutable vector< const QLink* > _sides;
1774 mutable bool _sideIsAdded[4]; // added in chain of links
1777 mutable const SMDS_MeshElement* _face;
1780 QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1782 void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1784 int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1786 void AddSelfToLinks() const {
1787 for ( int i = 0; i < _sides.size(); ++i )
1788 _sides[i]->_faces.push_back( this );
1790 int LinkIndex( const QLink* side ) const {
1791 for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1794 bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1796 bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1798 int i = LinkIndex( link._qlink );
1799 if ( i < 0 ) return true;
1800 _sideIsAdded[i] = true;
1801 link.SetFace( this );
1802 // continue from opposite link
1803 return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1805 bool IsBoundary() const { return !_volumes[1]; }
1807 bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1809 bool IsSpoiled(const QLink* bentLink ) const;
1811 TLinkInSet GetBoundaryLink( const TLinkSet& links,
1812 const TChainLink& avoidLink,
1813 TLinkInSet * notBoundaryLink = 0,
1814 const SMDS_MeshNode* nodeToContain = 0,
1815 bool * isAdjacentUsed = 0,
1816 int nbRecursionsLeft = -1) const;
1818 TLinkInSet GetLinkByNode( const TLinkSet& links,
1819 const TChainLink& avoidLink,
1820 const SMDS_MeshNode* nodeToContain) const;
1822 const SMDS_MeshNode* GetNodeInFace() const {
1823 for ( int iL = 0; iL < _sides.size(); ++iL )
1824 if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1828 gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1830 double MoveByBoundary( const TChainLink& theLink,
1831 const gp_Vec& theRefVec,
1832 const TLinkSet& theLinks,
1833 SMESH_MesherHelper* theFaceHelper=0,
1834 const double thePrevLen=0,
1835 const int theStep=theFirstStep,
1836 gp_Vec* theLinkNorm=0,
1837 double theSign=1.0) const;
1840 //================================================================================
1842 * \brief Dump QLink and QFace
1844 ostream& operator << (ostream& out, const QLink& l)
1846 out <<"QLink nodes: "
1847 << l.node1()->GetID() << " - "
1848 << l._mediumNode->GetID() << " - "
1849 << l.node2()->GetID() << endl;
1852 ostream& operator << (ostream& out, const QFace& f)
1854 out <<"QFace nodes: "/*<< &f << " "*/;
1855 for ( TIDSortedNodeSet::const_iterator n = f.begin(); n != f.end(); ++n )
1856 out << (*n)->GetID() << " ";
1857 out << " \tvolumes: "
1858 << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1859 << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1860 out << " \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1864 //================================================================================
1866 * \brief Construct QFace from QLinks
1868 //================================================================================
1870 QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
1872 _volumes[0] = _volumes[1] = 0;
1874 _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1875 _normal.SetCoord(0,0,0);
1876 for ( int i = 1; i < _sides.size(); ++i ) {
1877 const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1878 insert( l1->node1() ); insert( l1->node2() );
1880 gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1881 gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1882 if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1886 double normSqSize = _normal.SquareMagnitude();
1887 if ( normSqSize > numeric_limits<double>::min() )
1888 _normal /= sqrt( normSqSize );
1890 _normal.SetCoord(1e-33,0,0);
1896 //================================================================================
1898 * \brief Make up a chain of links
1899 * \param iSide - link to add first
1900 * \param chain - chain to fill in
1901 * \param pos - postion of medium nodes the links should have
1902 * \param error - out, specifies what is wrong
1903 * \retval bool - false if valid chain can't be built; "valid" means that links
1904 * of the chain belongs to rectangles bounding hexahedrons
1906 //================================================================================
1908 bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1910 if ( iSide >= _sides.size() ) // wrong argument iSide
1912 if ( _sideIsAdded[ iSide ]) // already in chain
1915 if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1918 list< const QFace* > faces( 1, this );
1919 while ( !faces.empty() ) {
1920 const QFace* face = faces.front();
1921 for ( int i = 0; i < face->_sides.size(); ++i ) {
1922 if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
1923 face->_sideIsAdded[i] = true;
1924 // find a face side in the chain
1925 TLinkInSet chLink = links.insert( TChainLink(face->_sides[i])).first;
1926 // TChain::iterator chLink = chain.begin();
1927 // for ( ; chLink != chain.end(); ++chLink )
1928 // if ( chLink->_qlink == face->_sides[i] )
1930 // if ( chLink == chain.end() )
1931 // chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
1932 // add a face to a chained link and put a continues face in the queue
1933 chLink->SetFace( face );
1934 if ( face->_sides[i]->MediumPos() >= pos )
1935 if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
1936 faces.push_back( contFace );
1941 if ( error < ERR_TRI )
1943 chain.insert( chain.end(), links.begin(),links.end() );
1946 _sideIsAdded[iSide] = true; // not to add this link to chain again
1947 const QLink* link = _sides[iSide];
1951 // add link into chain
1952 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1953 chLink->SetFace( this );
1956 // propagate from quadrangle to neighbour faces
1957 if ( link->MediumPos() >= pos ) {
1958 int nbLinkFaces = link->_faces.size();
1959 if ( nbLinkFaces == 4 || (nbLinkFaces < 4 && link->OnBoundary())) {
1960 // hexahedral mesh or boundary quadrangles - goto a continous face
1961 if ( const QFace* f = link->GetContinuesFace( this ))
1962 return f->GetLinkChain( *chLink, chain, pos, error );
1965 TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1966 for ( int i = 0; i < nbLinkFaces; ++i )
1967 if ( link->_faces[i] )
1968 link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1969 if ( error < ERR_PRISM )
1977 //================================================================================
1979 * \brief Return a boundary link of the triangle face
1980 * \param links - set of all links
1981 * \param avoidLink - link not to return
1982 * \param notBoundaryLink - out, neither the returned link nor avoidLink
1983 * \param nodeToContain - node the returned link must contain; if provided, search
1984 * also performed on adjacent faces
1985 * \param isAdjacentUsed - returns true if link is found in adjacent faces
1986 * \param nbRecursionsLeft - to limit recursion
1988 //================================================================================
1990 TLinkInSet QFace::GetBoundaryLink( const TLinkSet& links,
1991 const TChainLink& avoidLink,
1992 TLinkInSet * notBoundaryLink,
1993 const SMDS_MeshNode* nodeToContain,
1994 bool * isAdjacentUsed,
1995 int nbRecursionsLeft) const
1997 TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1999 typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
2000 TFaceLinkList adjacentFaces;
2002 for ( int iL = 0; iL < _sides.size(); ++iL )
2004 if ( avoidLink._qlink == _sides[iL] )
2006 TLinkInSet link = links.find( _sides[iL] );
2007 if ( link == linksEnd ) continue;
2008 if ( (*link)->MediumPos() > SMDS_TOP_FACE )
2009 continue; // We work on faces here, don't go inside a solid
2012 if ( link->IsBoundary() ) {
2013 if ( !nodeToContain ||
2014 (*link)->node1() == nodeToContain ||
2015 (*link)->node2() == nodeToContain )
2017 boundaryLink = link;
2018 if ( !notBoundaryLink ) break;
2021 else if ( notBoundaryLink ) {
2022 *notBoundaryLink = link;
2023 if ( boundaryLink != linksEnd ) break;
2026 if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
2027 if ( const QFace* adj = link->NextFace( this ))
2028 if ( adj->Contains( nodeToContain ))
2029 adjacentFaces.push_back( make_pair( adj, link ));
2032 if ( isAdjacentUsed ) *isAdjacentUsed = false;
2033 if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
2035 if ( nbRecursionsLeft < 0 )
2036 nbRecursionsLeft = nodeToContain->NbInverseElements();
2037 TFaceLinkList::iterator adj = adjacentFaces.begin();
2038 for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
2039 boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
2040 isAdjacentUsed, nbRecursionsLeft-1);
2041 if ( isAdjacentUsed ) *isAdjacentUsed = true;
2043 return boundaryLink;
2045 //================================================================================
2047 * \brief Return a link ending at the given node but not avoidLink
2049 //================================================================================
2051 TLinkInSet QFace::GetLinkByNode( const TLinkSet& links,
2052 const TChainLink& avoidLink,
2053 const SMDS_MeshNode* nodeToContain) const
2055 for ( int i = 0; i < _sides.size(); ++i )
2056 if ( avoidLink._qlink != _sides[i] &&
2057 (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
2058 return links.find( _sides[ i ]);
2062 //================================================================================
2064 * \brief Return normal to the i-th side pointing outside the face
2066 //================================================================================
2068 gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
2070 gp_Vec norm, vecOut;
2071 // if ( uvHelper ) {
2072 // TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
2073 // const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
2074 // gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
2075 // gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
2076 // norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
2078 // const QLink* otherLink = _sides[(i + 1) % _sides.size()];
2079 // const SMDS_MeshNode* otherNode =
2080 // otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
2081 // gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
2082 // vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
2085 norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
2086 gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
2087 XYZ( _sides[0]->node2() ) +
2088 XYZ( _sides[1]->node1() )) / 3.;
2089 vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
2091 if ( norm * vecOut < 0 )
2093 double mag2 = norm.SquareMagnitude();
2094 if ( mag2 > numeric_limits<double>::min() )
2095 norm /= sqrt( mag2 );
2098 //================================================================================
2100 * \brief Move medium node of theLink according to its distance from boundary
2101 * \param theLink - link to fix
2102 * \param theRefVec - movement of boundary
2103 * \param theLinks - all adjacent links of continous triangles
2104 * \param theFaceHelper - helper is not used so far
2105 * \param thePrevLen - distance from the boundary
2106 * \param theStep - number of steps till movement propagation limit
2107 * \param theLinkNorm - out normal to theLink
2108 * \param theSign - 1 or -1 depending on movement of boundary
2109 * \retval double - distance from boundary to propagation limit or other boundary
2111 //================================================================================
2113 double QFace::MoveByBoundary( const TChainLink& theLink,
2114 const gp_Vec& theRefVec,
2115 const TLinkSet& theLinks,
2116 SMESH_MesherHelper* theFaceHelper,
2117 const double thePrevLen,
2119 gp_Vec* theLinkNorm,
2120 double theSign) const
2123 return thePrevLen; // propagation limit reached
2125 int iL; // index of theLink
2126 for ( iL = 0; iL < _sides.size(); ++iL )
2127 if ( theLink._qlink == _sides[ iL ])
2130 MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
2131 <<" thePrevLen " << thePrevLen);
2132 MSG(string(theStep,'.')<<" "<<*theLink._qlink);
2134 gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
2135 double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
2136 if ( theStep == theFirstStep )
2137 theSign = refProj < 0. ? -1. : 1.;
2138 else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
2139 return thePrevLen; // to propagate movement forward only, not in side dir or backward
2141 int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
2142 TLinkInSet link1 = theLinks.find( _sides[iL1] );
2143 TLinkInSet link2 = theLinks.find( _sides[iL2] );
2144 if ( link1 == theLinks.end() || link2 == theLinks.end() )
2146 const QFace* f1 = link1->NextFace( this ); // adjacent faces
2147 const QFace* f2 = link2->NextFace( this );
2149 // propagate to adjacent faces till limit step or boundary
2150 double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
2151 double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
2152 gp_Vec linkDir1(0,0,0); // initialize to avoid valgrind error ("Conditional jump...")
2153 gp_Vec linkDir2(0,0,0);
2157 len1 = f1->MoveByBoundary
2158 ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
2160 linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
2162 MSG( " --------------- EXCEPTION");
2168 len2 = f2->MoveByBoundary
2169 ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
2171 linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
2173 MSG( " --------------- EXCEPTION");
2178 if ( theStep != theFirstStep )
2180 // choose chain length by direction of propagation most codirected with theRefVec
2181 bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
2182 fullLen = choose1 ? len1 : len2;
2183 double r = thePrevLen / fullLen;
2185 gp_Vec move = linkNorm * refProj * ( 1 - r );
2186 theLink->Move( move, true );
2188 MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
2189 " by " << refProj * ( 1 - r ) << " following " <<
2190 (choose1 ? *link1->_qlink : *link2->_qlink));
2192 if ( theLinkNorm ) *theLinkNorm = linkNorm;
2197 //================================================================================
2199 * \brief Checks if the face is distorted due to bentLink
2201 //================================================================================
2203 bool QFace::IsSpoiled(const QLink* bentLink ) const
2205 // code is valid for convex faces only
2207 for ( TIDSortedNodeSet::const_iterator n = begin(); n!=end(); ++n)
2208 gc += XYZ( *n ) / size();
2209 for (unsigned i = 0; i < _sides.size(); ++i )
2211 if ( _sides[i] == bentLink ) continue;
2212 gp_Vec linkNorm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
2213 gp_Vec vecOut( gc, _sides[i]->MiddlePnt() );
2214 if ( linkNorm * vecOut < 0 )
2216 double mag2 = linkNorm.SquareMagnitude();
2217 if ( mag2 > numeric_limits<double>::min() )
2218 linkNorm /= sqrt( mag2 );
2219 gp_Vec vecBent ( _sides[i]->MiddlePnt(), bentLink->MediumPnt());
2220 gp_Vec vecStraight( _sides[i]->MiddlePnt(), bentLink->MiddlePnt());
2221 if ( vecBent * linkNorm > -0.1*vecStraight.Magnitude() )
2228 //================================================================================
2230 * \brief Find pairs of continues faces
2232 //================================================================================
2234 void QLink::SetContinuesFaces() const
2236 // x0 x - QLink, [-|] - QFace, v - volume
2238 // | Between _faces of link x2 two vertical faces are continues
2239 // x1----x2-----x3 and two horizontal faces are continues. We set vertical faces
2240 // | to _faces[0] and _faces[1] and horizontal faces to
2241 // v2 | v3 _faces[2] and _faces[3] (or vise versa).
2244 if ( _faces.empty() )
2247 for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
2249 // look for a face bounding none of volumes bound by _faces[0]
2250 bool sameVol = false;
2251 int nbVol = _faces[iF]->NbVolumes();
2252 for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
2253 sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
2254 _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
2258 if ( iFaceCont > 0 ) // continues faces found, set one by the other
2260 if ( iFaceCont != 1 )
2261 std::swap( _faces[1], _faces[iFaceCont] );
2263 else if ( _faces.size() > 1 ) // not found, set NULL by the first face
2265 _faces.insert( ++_faces.begin(), 0 );
2268 //================================================================================
2270 * \brief Return a face continues to the given one
2272 //================================================================================
2274 const QFace* QLink::GetContinuesFace( const QFace* face ) const
2276 for ( int i = 0; i < _faces.size(); ++i ) {
2277 if ( _faces[i] == face ) {
2278 int iF = i < 2 ? 1-i : 5-i;
2279 return iF < _faces.size() ? _faces[iF] : 0;
2284 //================================================================================
2286 * \brief True if link is on mesh boundary
2288 //================================================================================
2290 bool QLink::OnBoundary() const
2292 for ( int i = 0; i < _faces.size(); ++i )
2293 if (_faces[i] && _faces[i]->IsBoundary()) return true;
2296 //================================================================================
2298 * \brief Return normal of link of the chain
2300 //================================================================================
2302 gp_Vec TChainLink::Normal() const {
2304 if (_qfaces[0]) norm = _qfaces[0]->_normal;
2305 if (_qfaces[1]) norm += _qfaces[1]->_normal;
2308 //================================================================================
2310 * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
2312 //================================================================================
2314 void fixPrism( TChain& allLinks )
2316 // separate boundary links from internal ones
2317 typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
2318 QLinkSet interLinks, bndLinks1, bndLink2;
2320 bool isCurved = false;
2321 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2322 if ( (*lnk)->OnBoundary() )
2323 bndLinks1.insert( lnk->_qlink );
2325 interLinks.insert( lnk->_qlink );
2326 isCurved = isCurved || !(*lnk)->IsStraight();
2329 return; // no need to move
2331 QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
2333 while ( !interLinks.empty() && !curBndLinks->empty() )
2335 // propagate movement from boundary links to connected internal links
2336 QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
2337 for ( ; bnd != bndEnd; ++bnd )
2339 const QLink* bndLink = *bnd;
2340 for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
2342 const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
2343 if ( !face ) continue;
2344 // find and move internal link opposite to bndLink within the face
2345 int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
2346 const QLink* interLink = face->_sides[ interInd ];
2347 QLinkSet::iterator pInterLink = interLinks.find( interLink );
2348 if ( pInterLink == interLinks.end() ) continue; // not internal link
2349 interLink->Move( bndLink->_nodeMove );
2350 // treated internal links become new boundary ones
2351 interLinks. erase( pInterLink );
2352 newBndLinks->insert( interLink );
2355 curBndLinks->clear();
2356 std::swap( curBndLinks, newBndLinks );
2360 //================================================================================
2362 * \brief Fix links of continues triangles near curved boundary
2364 //================================================================================
2366 void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2368 if ( allLinks.empty() ) return;
2370 TLinkSet linkSet( allLinks.begin(), allLinks.end());
2371 TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2373 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2375 if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2377 // move iff a boundary link is bent towards inside of a face (issue 0021084)
2378 const QFace* face = linkIt->_qfaces[0];
2379 gp_XYZ pIn = ( face->_sides[0]->MiddlePnt() +
2380 face->_sides[1]->MiddlePnt() +
2381 face->_sides[2]->MiddlePnt() ) / 3.;
2382 gp_XYZ insideDir( pIn - (*linkIt)->MiddlePnt());
2383 bool linkBentInside = ((*linkIt)->_nodeMove.Dot( insideDir ) > 0 );
2384 //if ( face->IsSpoiled( linkIt->_qlink ))
2385 if ( linkBentInside )
2386 face->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2391 //================================================================================
2393 * \brief Detect rectangular structure of links and build chains from them
2395 //================================================================================
2397 enum TSplitTriaResult {
2398 _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2399 _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK, _TWISTED_CHAIN };
2401 TSplitTriaResult splitTrianglesIntoChains( TChain & allLinks,
2402 vector< TChain> & resultChains,
2403 SMDS_TypeOfPosition pos )
2405 // put links in the set and evalute number of result chains by number of boundary links
2408 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2409 linkSet.insert( *lnk );
2410 nbBndLinks += lnk->IsBoundary();
2412 resultChains.clear();
2413 resultChains.reserve( nbBndLinks / 2 );
2415 TLinkInSet linkIt, linksEnd = linkSet.end();
2417 // find a boundary link with corner node; corner node has position pos-2
2418 // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2420 SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2421 const SMDS_MeshNode* corner = 0;
2422 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2423 if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2428 TLinkInSet startLink = linkIt;
2429 const SMDS_MeshNode* startCorner = corner;
2430 vector< TChain* > rowChains;
2433 while ( startLink != linksEnd) // loop on columns
2435 // We suppose we have a rectangular structure like shown here. We have found a
2436 // corner of the rectangle (startCorner) and a boundary link sharing
2437 // |/ |/ | the startCorner (startLink). We are going to loop on rows of the
2438 // --o---o---o structure making several chains at once. One chain (columnChain)
2439 // |\ | /| starts at startLink and continues upward (we look at the structure
2440 // \ | \ | / | from such point that startLink is on the bottom of the structure).
2441 // \| \|/ | While going upward we also fill horizontal chains (rowChains) we
2442 // --o---o---o encounter.
2444 // / | \ | \ | startCorner
2449 if ( resultChains.size() == nbBndLinks / 2 )
2451 resultChains.push_back( TChain() );
2452 TChain& columnChain = resultChains.back();
2454 TLinkInSet botLink = startLink; // current horizontal link to go up from
2455 corner = startCorner; // current corner the botLink ends at
2457 while ( botLink != linksEnd ) // loop on rows
2459 // add botLink to the columnChain
2460 columnChain.push_back( *botLink );
2462 const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2464 { // the column ends
2465 if ( botLink == startLink )
2466 return _TWISTED_CHAIN; // issue 0020951
2467 linkSet.erase( botLink );
2468 if ( iRow != rowChains.size() )
2469 return _FEW_ROWS; // different nb of rows in columns
2472 // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2473 // link ending at <corner> (sideLink); there are two cases:
2474 // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2475 // since midQuadLink is not at boundary while sideLink is.
2476 // 2) midQuadLink ends at <corner>
2478 TLinkInSet midQuadLink = linksEnd;
2479 TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2481 if ( isCase2 ) { // find midQuadLink among links of botTria
2482 midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2483 if ( midQuadLink->IsBoundary() )
2484 return _BAD_MIDQUAD;
2486 if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2487 return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2490 columnChain.push_back( *midQuadLink );
2491 if ( iRow >= rowChains.size() ) {
2493 return _MANY_ROWS; // different nb of rows in columns
2494 if ( resultChains.size() == nbBndLinks / 2 )
2496 resultChains.push_back( TChain() );
2497 rowChains.push_back( & resultChains.back() );
2499 rowChains[iRow]->push_back( *sideLink );
2500 rowChains[iRow]->push_back( *midQuadLink );
2502 const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2506 // prepare startCorner and startLink for the next column
2507 startCorner = startLink->NextNode( startCorner );
2509 startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2511 startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2512 // check if no more columns remains
2513 if ( startLink != linksEnd ) {
2514 const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2515 if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2516 startLink = linksEnd; // startLink bounds upTria or botTria
2517 else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2521 // find bottom link and corner for the next row
2522 corner = sideLink->NextNode( corner );
2523 // next bottom link ends at the new corner
2524 linkSet.erase( botLink );
2525 botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2526 if ( botLink == linksEnd || botLink == midQuadLink || botLink == sideLink)
2528 if ( midQuadLink == startLink || sideLink == startLink )
2529 return _TWISTED_CHAIN; // issue 0020951
2530 linkSet.erase( midQuadLink );
2531 linkSet.erase( sideLink );
2533 // make faces neighboring the found ones be boundary
2534 if ( startLink != linksEnd ) {
2535 const QFace* tria = isCase2 ? botTria : upTria;
2536 for ( int iL = 0; iL < 3; ++iL ) {
2537 linkIt = linkSet.find( tria->_sides[iL] );
2538 if ( linkIt != linksEnd )
2539 linkIt->RemoveFace( tria );
2542 if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2543 botLink->RemoveFace( upTria ); // make next botTria first in vector
2550 // In the linkSet, there must remain the last links of rowChains; add them
2551 if ( linkSet.size() != rowChains.size() )
2552 return _BAD_SET_SIZE;
2553 for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2554 // find the link (startLink) ending at startCorner
2556 for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2557 if ( (*startLink)->node1() == startCorner ) {
2558 corner = (*startLink)->node2(); break;
2560 else if ( (*startLink)->node2() == startCorner) {
2561 corner = (*startLink)->node1(); break;
2564 if ( startLink == linksEnd )
2566 rowChains[ iRow ]->push_back( *startLink );
2567 linkSet.erase( startLink );
2568 startCorner = corner;
2575 //=======================================================================
2577 * \brief Move medium nodes of faces and volumes to fix distorted elements
2578 * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2580 * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2582 //=======================================================================
2584 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2586 // 0. Apply algorithm to solids or geom faces
2587 // ----------------------------------------------
2588 if ( myShape.IsNull() ) {
2589 if ( !myMesh->HasShapeToMesh() ) return;
2590 SetSubShape( myMesh->GetShapeToMesh() );
2594 TopTools_IndexedMapOfShape solids;
2595 TopExp::MapShapes(myShape,TopAbs_SOLID,solids);
2596 nbSolids = solids.Extent();
2598 TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2599 for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2600 faces.Add( f.Current() ); // not in solid
2602 for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2603 if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2604 for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2605 faces.Add( f.Current() ); // in not meshed solid
2607 else { // fix nodes in the solid and its faces
2608 MSG("FIX SOLID " << nbSolids-- << " #" << GetMeshDS()->ShapeToIndex(s.Current()));
2609 SMESH_MesherHelper h(*myMesh);
2610 h.SetSubShape( s.Current() );
2611 h.FixQuadraticElements(false);
2614 // fix nodes on geom faces
2616 int nbfaces = faces.Extent();
2618 for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2619 MSG("FIX FACE " << nbfaces-- << " #" << GetMeshDS()->ShapeToIndex(fIt.Key()));
2620 SMESH_MesherHelper h(*myMesh);
2621 h.SetSubShape( fIt.Key() );
2622 h.FixQuadraticElements(true);
2624 //perf_print_all_meters(1);
2628 // 1. Find out type of elements and get iterator on them
2629 // ---------------------------------------------------
2631 SMDS_ElemIteratorPtr elemIt;
2632 SMDSAbs_ElementType elemType = SMDSAbs_All;
2634 SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2637 if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2638 elemIt = smDS->GetElements();
2639 if ( elemIt->more() ) {
2640 elemType = elemIt->next()->GetType();
2641 elemIt = smDS->GetElements();
2644 if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2647 // 2. Fill in auxiliary data structures
2648 // ----------------------------------
2652 set< QLink >::iterator pLink;
2653 set< QFace >::iterator pFace;
2655 bool isCurved = false;
2656 //bool hasRectFaces = false;
2657 //set<int> nbElemNodeSet;
2659 if ( elemType == SMDSAbs_Volume )
2661 SMDS_VolumeTool volTool;
2662 while ( elemIt->more() ) // loop on volumes
2664 const SMDS_MeshElement* vol = elemIt->next();
2665 if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2667 for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2669 int nbN = volTool.NbFaceNodes( iF );
2670 //nbElemNodeSet.insert( nbN );
2671 const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2672 vector< const QLink* > faceLinks( nbN/2 );
2673 for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2676 QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2677 pLink = links.insert( link ).first;
2678 faceLinks[ iN/2 ] = & *pLink;
2680 isCurved = !link.IsStraight();
2681 if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2682 return; // already fixed
2685 pFace = faces.insert( QFace( faceLinks )).first;
2686 if ( pFace->NbVolumes() == 0 )
2687 pFace->AddSelfToLinks();
2688 pFace->SetVolume( vol );
2689 // hasRectFaces = hasRectFaces ||
2690 // ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2691 // volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2694 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2696 pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2697 faceNodes[4],faceNodes[6] );
2701 set< QLink >::iterator pLink = links.begin();
2702 for ( ; pLink != links.end(); ++pLink )
2703 pLink->SetContinuesFaces();
2707 while ( elemIt->more() ) // loop on faces
2709 const SMDS_MeshElement* face = elemIt->next();
2710 if ( !face->IsQuadratic() )
2712 //nbElemNodeSet.insert( face->NbNodes() );
2713 int nbN = face->NbNodes()/2;
2714 vector< const QLink* > faceLinks( nbN );
2715 for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2718 QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2719 pLink = links.insert( link ).first;
2720 faceLinks[ iN ] = & *pLink;
2722 isCurved = !link.IsStraight();
2725 pFace = faces.insert( QFace( faceLinks )).first;
2726 pFace->AddSelfToLinks();
2727 //hasRectFaces = ( hasRectFaces || nbN == 4 );
2731 return; // no curved edges of faces
2733 // 3. Compute displacement of medium nodes
2734 // -------------------------------------
2736 // two loops on faces: the first is to treat boundary links, the second is for internal ones
2737 TopLoc_Location loc;
2738 // not treat boundary of volumic submesh
2739 int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2740 for ( ; isInside < 2; ++isInside ) {
2741 MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2742 SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2743 SMDS_TypeOfPosition bndPos = isInside ? SMDS_TOP_FACE : SMDS_TOP_EDGE;
2745 for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2746 if ( bool(isInside) == pFace->IsBoundary() )
2748 for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2751 // make chain of links connected via continues faces
2754 if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2756 if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2758 vector< TChain > chains;
2759 if ( error == ERR_OK ) { // chain contains continues rectangles
2761 chains[0].splice( chains[0].begin(), rawChain );
2763 else if ( error == ERR_TRI ) { // chain contains continues triangles
2764 TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2765 if ( res != _OK ) { // not quadrangles split into triangles
2766 fixTriaNearBoundary( rawChain, *this );
2770 else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2771 fixPrism( rawChain );
2777 for ( int iC = 0; iC < chains.size(); ++iC )
2779 TChain& chain = chains[iC];
2780 if ( chain.empty() ) continue;
2781 if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2782 MSG("3D straight - ignore");
2785 if ( chain.front()->MediumPos() > bndPos ||
2786 chain.back()->MediumPos() > bndPos ) {
2787 MSG("Internal chain - ignore");
2790 // mesure chain length and compute link position along the chain
2791 double chainLen = 0;
2792 vector< double > linkPos;
2793 MSGBEG( "Link medium nodes: ");
2794 TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2795 for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2796 MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2797 double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2798 while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2799 link1 = chain.erase( link1 );
2800 if ( link1 == chain.end() )
2802 len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2805 linkPos.push_back( chainLen );
2808 if ( linkPos.size() < 2 )
2811 gp_Vec move0 = chain.front()->_nodeMove;
2812 gp_Vec move1 = chain.back ()->_nodeMove;
2815 bool checkUV = true;
2817 // compute node displacement of end links in parametric space of face
2818 const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2819 TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2820 if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2822 face = TopoDS::Face( f );
2823 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2825 for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2827 TChainLink& link = is1 ? chain.back() : chain.front();
2828 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2829 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2830 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2831 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2832 // uvMove = uvm - uv12
2833 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2834 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2835 if ( !is1 ) // correct nodeOnFace for move1 (issue 0020919)
2836 nodeOnFace = (*(++chain.rbegin()))->_mediumNode;
2837 isStraight[is1] = isStraightLink( (uv2-uv1).SquareModulus(),uvMove.SquareModulus());
2839 // if ( move0.SquareMagnitude() < straightTol2 &&
2840 // move1.SquareMagnitude() < straightTol2 ) {
2841 if ( isStraight[0] && isStraight[1] ) {
2842 MSG("2D straight - ignore");
2843 continue; // straight - no need to move nodes of internal links
2848 if ( isInside || face.IsNull() )
2850 // compute node displacement of end links in their local coord systems
2852 TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2853 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2854 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2855 move0.Transform(trsf);
2858 TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2859 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2860 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2861 move1.Transform(trsf);
2864 // compute displacement of medium nodes
2865 link2 = chain.begin();
2868 for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2870 double r = linkPos[i] / chainLen;
2871 // displacement in local coord system
2872 gp_Vec move = (1. - r) * move0 + r * move1;
2873 if ( isInside || face.IsNull()) {
2874 // transform to global
2875 gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2876 gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2877 gp_Vec x = x01.Normalized() + x12.Normalized();
2878 trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2879 move.Transform(trsf);
2882 // compute 3D displacement by 2D one
2883 Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
2884 gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2885 gp_XY newUV = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
2886 gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
2887 move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2889 if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2890 move.SquareMagnitude())
2892 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2893 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2894 MSG( "TOO LONG MOVE \t" <<
2895 "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2896 "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2897 "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2898 "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2902 (*link1)->Move( move );
2903 MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2904 << chain.front()->_mediumNode->GetID() <<"-"
2905 << chain.back ()->_mediumNode->GetID() <<
2906 " by " << move.Magnitude());
2908 } // loop on chains of links
2909 } // loop on 2 directions of propagation from quadrangle
2916 for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2917 if ( pLink->IsMoved() ) {
2918 //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2919 gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2920 GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
2925 //=======================================================================
2927 * \brief Iterator on ancestors of the given type
2929 //=======================================================================
2931 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2933 TopTools_ListIteratorOfListOfShape _ancIter;
2934 TopAbs_ShapeEnum _type;
2935 TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2936 : _ancIter( ancestors ), _type( type )
2938 if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
2942 return _ancIter.More();
2944 virtual const TopoDS_Shape* next()
2946 const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2947 if ( _ancIter.More() )
2948 for ( _ancIter.Next(); _ancIter.More(); _ancIter.Next())
2949 if ( _ancIter.Value().ShapeType() == _type )
2955 //=======================================================================
2957 * \brief Return iterator on ancestors of the given type
2959 //=======================================================================
2961 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2962 const SMESH_Mesh& mesh,
2963 TopAbs_ShapeEnum ancestorType)
2965 return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));