1 // Copyright (C) 2007-2014 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, or (at your option) any later version.
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_EdgePosition.hxx"
30 #include "SMDS_FaceOfNodes.hxx"
31 #include "SMDS_FacePosition.hxx"
32 #include "SMDS_IteratorOnIterators.hxx"
33 #include "SMDS_VolumeTool.hxx"
34 #include "SMESH_Block.hxx"
35 #include "SMESH_MeshAlgos.hxx"
36 #include "SMESH_ProxyMesh.hxx"
37 #include "SMESH_subMesh.hxx"
39 #include <BRepAdaptor_Curve.hxx>
40 #include <BRepAdaptor_Surface.hxx>
41 #include <BRepTools.hxx>
42 #include <BRep_Tool.hxx>
43 #include <Geom2d_Curve.hxx>
44 #include <GeomAPI_ProjectPointOnCurve.hxx>
45 #include <GeomAPI_ProjectPointOnSurf.hxx>
46 #include <Geom_Curve.hxx>
47 #include <Geom_RectangularTrimmedSurface.hxx>
48 #include <Geom_Surface.hxx>
49 #include <ShapeAnalysis.hxx>
51 #include <TopExp_Explorer.hxx>
52 #include <TopTools_ListIteratorOfListOfShape.hxx>
53 #include <TopTools_MapIteratorOfMapOfShape.hxx>
54 #include <TopTools_MapOfShape.hxx>
57 #include <gp_Pnt2d.hxx>
58 #include <gp_Trsf.hxx>
60 #include <Standard_Failure.hxx>
61 #include <Standard_ErrorHandler.hxx>
63 #include <utilities.h>
69 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
73 gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
75 enum { U_periodic = 1, V_periodic = 2 };
78 //================================================================================
82 //================================================================================
84 SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
88 myCreateQuadratic(false),
89 myCreateBiQuadratic(false),
90 myFixNodeParameters(false)
92 myPar1[0] = myPar2[0] = myPar1[1] = myPar2[1] = 0;
93 mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
96 //=======================================================================
97 //function : ~SMESH_MesherHelper
99 //=======================================================================
101 SMESH_MesherHelper::~SMESH_MesherHelper()
104 TID2ProjectorOnSurf::iterator i_proj = myFace2Projector.begin();
105 for ( ; i_proj != myFace2Projector.end(); ++i_proj )
106 delete i_proj->second;
109 TID2ProjectorOnCurve::iterator i_proj = myEdge2Projector.begin();
110 for ( ; i_proj != myEdge2Projector.end(); ++i_proj )
111 delete i_proj->second;
115 //=======================================================================
116 //function : IsQuadraticSubMesh
117 //purpose : Check submesh for given shape: if all elements on this shape
118 // are quadratic, quadratic elements will be created.
119 // Also fill myTLinkNodeMap
120 //=======================================================================
122 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
124 SMESHDS_Mesh* meshDS = GetMeshDS();
125 // we can create quadratic elements only if all elements
126 // created on sub-shapes of given shape are quadratic
127 // also we have to fill myTLinkNodeMap
128 myCreateQuadratic = true;
129 mySeamShapeIds.clear();
130 myDegenShapeIds.clear();
131 TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
132 if ( aSh.ShapeType()==TopAbs_COMPOUND )
134 TopoDS_Iterator subIt( aSh );
136 subType = ( subIt.Value().ShapeType()==TopAbs_FACE ) ? TopAbs_EDGE : TopAbs_FACE;
138 SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
141 int nbOldLinks = myTLinkNodeMap.size();
143 if ( !myMesh->HasShapeToMesh() )
145 if (( myCreateQuadratic = myMesh->NbFaces( ORDER_QUADRATIC )))
147 SMDS_FaceIteratorPtr fIt = meshDS->facesIterator();
148 while ( fIt->more() )
149 AddTLinks( static_cast< const SMDS_MeshFace* >( fIt->next() ));
154 TopExp_Explorer exp( aSh, subType );
155 TopTools_MapOfShape checkedSubShapes;
156 for (; exp.More() && myCreateQuadratic; exp.Next()) {
157 if ( !checkedSubShapes.Add( exp.Current() ))
158 continue; // needed if aSh is compound of solids
159 if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
160 if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
162 const SMDS_MeshElement* e = it->next();
163 if ( e->GetType() != elemType || !e->IsQuadratic() ) {
164 myCreateQuadratic = false;
169 switch ( e->NbCornerNodes() ) {
171 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
173 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
174 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
175 AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
177 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
178 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
179 AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
180 AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
183 myCreateQuadratic = false;
193 if ( nbOldLinks == myTLinkNodeMap.size() )
194 myCreateQuadratic = false;
196 if(!myCreateQuadratic) {
197 myTLinkNodeMap.clear();
201 return myCreateQuadratic;
204 //=======================================================================
205 //function : SetSubShape
206 //purpose : Set geometry to make elements on
207 //=======================================================================
209 void SMESH_MesherHelper::SetSubShape(const int aShID)
211 if ( aShID == myShapeID )
214 SetSubShape( GetMeshDS()->IndexToShape( aShID ));
216 SetSubShape( TopoDS_Shape() );
219 //=======================================================================
220 //function : SetSubShape
221 //purpose : Set geometry to create elements on
222 //=======================================================================
224 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
226 if ( myShape.IsSame( aSh ))
230 mySeamShapeIds.clear();
231 myDegenShapeIds.clear();
233 if ( myShape.IsNull() ) {
237 SMESHDS_Mesh* meshDS = GetMeshDS();
238 myShapeID = meshDS->ShapeToIndex(aSh);
241 // treatment of periodic faces
242 for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
244 const TopoDS_Face& face = TopoDS::Face( eF.Current() );
246 Handle(Geom_Surface) surface = BRep_Tool::Surface( face, loc );
248 // if ( surface->IsUPeriodic() || surface->IsVPeriodic() ||
249 // surface->IsUClosed() || surface->IsVClosed() )
251 //while ( surface->IsKind(STANDARD_TYPE(Geom_RectangularTrimmedSurface )))
252 //surface = Handle(Geom_RectangularTrimmedSurface)::DownCast( surface )->BasisSurface();
253 GeomAdaptor_Surface surf( surface );
255 for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
257 // look for a seam edge
258 const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
259 if ( BRep_Tool::IsClosed( edge, face )) {
260 // initialize myPar1, myPar2 and myParIndex
262 BRep_Tool::UVPoints( edge, face, uv1, uv2 );
263 if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
265 myParIndex |= U_periodic;
266 myPar1[0] = surf.FirstUParameter();
267 myPar2[0] = surf.LastUParameter();
270 myParIndex |= V_periodic;
271 myPar1[1] = surf.FirstVParameter();
272 myPar2[1] = surf.LastVParameter();
274 // store seam shape indices, negative if shape encounters twice
275 int edgeID = meshDS->ShapeToIndex( edge );
276 mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
277 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
278 int vertexID = meshDS->ShapeToIndex( v.Current() );
279 mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
283 // look for a degenerated edge
284 if ( SMESH_Algo::isDegenerated( edge )) {
285 myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
286 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
287 myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
290 if ( !myDegenShapeIds.empty() && !myParIndex ) {
291 if ( surface->IsUPeriodic() || surface->IsUClosed() ) {
292 myParIndex |= U_periodic;
293 myPar1[0] = surf.FirstUParameter();
294 myPar2[0] = surf.LastUParameter();
296 else if ( surface->IsVPeriodic() || surface->IsVClosed() ) {
297 myParIndex |= V_periodic;
298 myPar1[1] = surf.FirstVParameter();
299 myPar2[1] = surf.LastVParameter();
306 //=======================================================================
307 //function : GetNodeUVneedInFaceNode
308 //purpose : Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
309 // Return true if the face is periodic.
310 // If F is Null, answer about sub-shape set through IsQuadraticSubMesh() or
312 //=======================================================================
314 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
316 if ( F.IsNull() ) return !mySeamShapeIds.empty();
318 if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
319 return !mySeamShapeIds.empty();
322 Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
323 if ( !aSurface.IsNull() )
324 return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
329 //=======================================================================
330 //function : IsMedium
332 //=======================================================================
334 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode* node,
335 const SMDSAbs_ElementType typeToCheck)
337 return SMESH_MeshEditor::IsMedium( node, typeToCheck );
340 //=======================================================================
341 //function : GetSubShapeByNode
342 //purpose : Return support shape of a node
343 //=======================================================================
345 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
346 const SMESHDS_Mesh* meshDS)
348 int shapeID = node ? node->getshapeId() : 0;
349 if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
350 return meshDS->IndexToShape( shapeID );
352 return TopoDS_Shape();
356 //=======================================================================
357 //function : AddTLinkNode
358 //purpose : add a link in my data structure
359 //=======================================================================
361 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
362 const SMDS_MeshNode* n2,
363 const SMDS_MeshNode* n12)
365 // add new record to map
366 SMESH_TLink link( n1, n2 );
367 myTLinkNodeMap.insert( make_pair(link,n12));
370 //================================================================================
372 * \brief Add quadratic links of edge to own data structure
374 //================================================================================
376 void SMESH_MesherHelper::AddTLinks(const SMDS_MeshEdge* edge)
378 if ( edge->IsQuadratic() )
379 AddTLinkNode(edge->GetNode(0), edge->GetNode(1), edge->GetNode(2));
382 //================================================================================
384 * \brief Add quadratic links of face to own data structure
386 //================================================================================
388 void SMESH_MesherHelper::AddTLinks(const SMDS_MeshFace* f)
391 switch ( f->NbNodes() ) {
393 // myMapWithCentralNode.insert
394 // ( make_pair( TBiQuad( f->GetNode(0),f->GetNode(1),f->GetNode(2) ),
396 // break; -- add medium nodes as well
398 AddTLinkNode(f->GetNode(0),f->GetNode(1),f->GetNode(3));
399 AddTLinkNode(f->GetNode(1),f->GetNode(2),f->GetNode(4));
400 AddTLinkNode(f->GetNode(2),f->GetNode(0),f->GetNode(5)); break;
403 // myMapWithCentralNode.insert
404 // ( make_pair( TBiQuad( f->GetNode(0),f->GetNode(1),f->GetNode(2),f->GetNode(3) ),
406 // break; -- add medium nodes as well
408 AddTLinkNode(f->GetNode(0),f->GetNode(1),f->GetNode(4));
409 AddTLinkNode(f->GetNode(1),f->GetNode(2),f->GetNode(5));
410 AddTLinkNode(f->GetNode(2),f->GetNode(3),f->GetNode(6));
411 AddTLinkNode(f->GetNode(3),f->GetNode(0),f->GetNode(7)); break;
416 //================================================================================
418 * \brief Add quadratic links of volume to own data structure
420 //================================================================================
422 void SMESH_MesherHelper::AddTLinks(const SMDS_MeshVolume* volume)
424 if ( volume->IsQuadratic() )
426 SMDS_VolumeTool vTool( volume );
427 const SMDS_MeshNode** nodes = vTool.GetNodes();
429 for ( int iF = 1; iF < vTool.NbFaces(); ++iF )
431 const int nbN = vTool.NbFaceNodes( iF );
432 const int* iNodes = vTool.GetFaceNodesIndices( iF );
433 for ( int i = 0; i < nbN; )
435 int iN1 = iNodes[i++];
436 int iN12 = iNodes[i++];
438 if ( iN1 > iN2 ) std::swap( iN1, iN2 );
439 int linkID = iN1 * vTool.NbNodes() + iN2;
440 pair< set<int>::iterator, bool > it_isNew = addedLinks.insert( linkID );
441 if ( it_isNew.second )
442 AddTLinkNode( nodes[iN1], nodes[iN2], nodes[iN12] );
444 addedLinks.erase( it_isNew.first ); // each link encounters only twice
446 if ( vTool.NbNodes() == 27 )
448 const SMDS_MeshNode* nFCenter = nodes[ vTool.GetCenterNodeIndex( iF )];
449 if ( nFCenter->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE )
450 myMapWithCentralNode.insert
451 ( make_pair( TBiQuad( nodes[ iNodes[0]], nodes[ iNodes[1]],
452 nodes[ iNodes[2]], nodes[ iNodes[3]] ),
459 //================================================================================
461 * \brief Return true if position of nodes on the shape hasn't yet been checked or
462 * the positions proved to be invalid
464 //================================================================================
466 bool SMESH_MesherHelper::toCheckPosOnShape(int shapeID ) const
468 map< int,bool >::const_iterator id_ok = myNodePosShapesValidity.find( shapeID );
469 return ( id_ok == myNodePosShapesValidity.end() || !id_ok->second );
472 //================================================================================
474 * \brief Set validity of positions of nodes on the shape.
475 * Once set, validity is not changed
477 //================================================================================
479 void SMESH_MesherHelper::setPosOnShapeValidity(int shapeID, bool ok ) const
481 std::map< int,bool >::iterator sh_ok =
482 ((SMESH_MesherHelper*)this)->myNodePosShapesValidity.insert( make_pair( shapeID, ok)).first;
487 //=======================================================================
488 //function : ToFixNodeParameters
489 //purpose : Enables fixing node parameters on EDGEs and FACEs in
490 // GetNodeU(...,check=true), GetNodeUV(...,check=true), CheckNodeUV() and
491 // CheckNodeU() in case if a node lies on a shape set via SetSubShape().
493 //=======================================================================
495 void SMESH_MesherHelper::ToFixNodeParameters(bool toFix)
497 myFixNodeParameters = toFix;
501 //=======================================================================
502 //function : GetUVOnSeam
503 //purpose : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
504 //=======================================================================
506 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
508 gp_Pnt2d result = uv1;
509 for ( int i = U_periodic; i <= V_periodic ; ++i )
511 if ( myParIndex & i )
513 double p1 = uv1.Coord( i );
514 double dp1 = Abs( p1-myPar1[i-1]), dp2 = Abs( p1-myPar2[i-1]);
515 if ( myParIndex == i ||
516 dp1 < ( myPar2[i-1] - myPar2[i-1] ) / 100. ||
517 dp2 < ( myPar2[i-1] - myPar2[i-1] ) / 100. )
519 double p2 = uv2.Coord( i );
520 double p1Alt = ( dp1 < dp2 ) ? myPar2[i-1] : myPar1[i-1];
521 if ( Abs( p2 - p1 ) > Abs( p2 - p1Alt ))
522 result.SetCoord( i, p1Alt );
529 //=======================================================================
530 //function : GetNodeUV
531 //purpose : Return node UV on face
532 //=======================================================================
534 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
535 const SMDS_MeshNode* n,
536 const SMDS_MeshNode* n2,
539 gp_Pnt2d uv( Precision::Infinite(), Precision::Infinite() );
541 const SMDS_PositionPtr Pos = n->GetPosition();
543 if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
545 // node has position on face
546 const SMDS_FacePosition* fpos =
547 static_cast<const SMDS_FacePosition*>(n->GetPosition());
548 uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
550 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ));
552 else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
554 // node has position on edge => it is needed to find
555 // corresponding edge from face, get pcurve for this
556 // edge and retrieve value from this pcurve
557 const SMDS_EdgePosition* epos =
558 static_cast<const SMDS_EdgePosition*>(n->GetPosition());
559 int edgeID = n->getshapeId();
560 TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
561 double f, l, u = epos->GetUParameter();
562 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
563 bool validU = ( f < u && u < l );
565 uv = C2d->Value( u );
567 uv.SetCoord( Precision::Infinite(),0.);
568 if ( check || !validU )
569 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ),/*force=*/ !validU );
571 // for a node on a seam edge select one of UVs on 2 pcurves
572 if ( n2 && IsSeamShape( edgeID ) )
574 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0, check ));
577 { // adjust uv to period
579 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
580 Standard_Boolean isUPeriodic = S->IsUPeriodic();
581 Standard_Boolean isVPeriodic = S->IsVPeriodic();
582 if ( isUPeriodic || isVPeriodic ) {
583 Standard_Real UF,UL,VF,VL;
584 S->Bounds(UF,UL,VF,VL);
586 uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
588 uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
592 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
594 if ( int vertexID = n->getshapeId() ) {
595 const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
597 uv = BRep_Tool::Parameters( V, F );
600 catch (Standard_Failure& exc) {
603 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() )
604 uvOK = ( V == vert.Current() );
606 MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
607 << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
608 // get UV of a vertex closest to the node
610 gp_Pnt pn = XYZ( n );
611 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) {
612 TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
613 gp_Pnt p = BRep_Tool::Pnt( curV );
614 double curDist = p.SquareDistance( pn );
615 if ( curDist < dist ) {
617 uv = BRep_Tool::Parameters( curV, F );
618 uvOK = ( dist < DBL_MIN );
624 TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
625 for ( ; it.More(); it.Next() ) {
626 if ( it.Value().ShapeType() == TopAbs_EDGE ) {
627 const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
629 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
630 if ( !C2d.IsNull() ) {
631 double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
632 uv = C2d->Value( u );
640 if ( n2 && IsSeamShape( vertexID ) )
641 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
646 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 10*MaxTolerance( F ));
655 //=======================================================================
656 //function : CheckNodeUV
657 //purpose : Check and fix node UV on a face
658 //=======================================================================
660 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
661 const SMDS_MeshNode* n,
665 double distXYZ[4]) const
667 int shapeID = n->getshapeId();
668 bool infinit = ( Precision::IsInfinite( uv.X() ) || Precision::IsInfinite( uv.Y() ));
669 bool zero = ( uv.X() == 0. && uv.Y() == 0. );
670 if ( force || toCheckPosOnShape( shapeID ) || infinit || zero )
672 // check that uv is correct
674 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
675 gp_Pnt nodePnt = XYZ( n ), surfPnt(0,0,0);
677 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
679 (dist = nodePnt.Distance( surfPnt = surface->Value( uv.X(), uv.Y() ))) > tol )
681 setPosOnShapeValidity( shapeID, false );
682 if ( !infinit && distXYZ ) {
683 surfPnt.Transform( loc );
685 distXYZ[1] = surfPnt.X(); distXYZ[2] = surfPnt.Y(); distXYZ[3]=surfPnt.Z();
687 // uv incorrect, project the node to surface
688 GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
689 projector.Perform( nodePnt );
690 if ( !projector.IsDone() || projector.NbPoints() < 1 )
692 MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
695 Quantity_Parameter U,V;
696 projector.LowerDistanceParameters(U,V);
698 surfPnt = surface->Value( U, V );
699 dist = nodePnt.Distance( surfPnt );
701 surfPnt.Transform( loc );
703 distXYZ[1] = surfPnt.X(); distXYZ[2] = surfPnt.Y(); distXYZ[3]=surfPnt.Z();
707 MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
710 // store the fixed UV on the face
711 if ( myShape.IsSame(F) && shapeID == myShapeID && myFixNodeParameters )
712 const_cast<SMDS_MeshNode*>(n)->SetPosition
713 ( SMDS_PositionPtr( new SMDS_FacePosition( U, V )));
715 else if ( uv.Modulus() > numeric_limits<double>::min() )
717 setPosOnShapeValidity( shapeID, true );
723 //=======================================================================
724 //function : GetProjector
725 //purpose : Return projector intitialized by given face without location, which is returned
726 //=======================================================================
728 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
729 TopLoc_Location& loc,
732 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
733 int faceID = GetMeshDS()->ShapeToIndex( F );
734 TID2ProjectorOnSurf& i2proj = const_cast< TID2ProjectorOnSurf&>( myFace2Projector );
735 TID2ProjectorOnSurf::iterator i_proj = i2proj.find( faceID );
736 if ( i_proj == i2proj.end() )
738 if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
739 double U1, U2, V1, V2;
740 surface->Bounds(U1, U2, V1, V2);
741 GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
742 proj->Init( surface, U1, U2, V1, V2, tol );
743 i_proj = i2proj.insert( make_pair( faceID, proj )).first;
745 return *( i_proj->second );
750 gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
751 gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
752 gp_XY_FunPtr(Subtracted);
755 //=======================================================================
756 //function : applyIn2D
757 //purpose : Perform given operation on two 2d points in parameric space of given surface.
758 // It takes into account period of the surface. Use gp_XY_FunPtr macro
759 // to easily define pointer to function of gp_XY class.
760 //=======================================================================
762 gp_XY SMESH_MesherHelper::applyIn2D(const Handle(Geom_Surface)& surface,
766 const bool resultInPeriod)
768 Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
769 Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
770 if ( !isUPeriodic && !isVPeriodic )
773 // move uv2 not far than half-period from uv1
775 uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
777 uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
780 gp_XY res = fun( uv1, gp_XY(u2,v2) );
782 // move result within period
783 if ( resultInPeriod )
785 Standard_Real UF,UL,VF,VL;
786 surface->Bounds(UF,UL,VF,VL);
788 res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
790 res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
795 //=======================================================================
796 //function : GetMiddleUV
797 //purpose : Return middle UV taking in account surface period
798 //=======================================================================
800 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
805 // the proper place of getting basic surface seems to be in applyIn2D()
806 // but we put it here to decrease a risk of regressions just before releasing a version
807 Handle(Geom_Surface) surf = surface;
808 while ( !surf.IsNull() && surf->IsKind(STANDARD_TYPE(Geom_RectangularTrimmedSurface )))
809 surf = Handle(Geom_RectangularTrimmedSurface)::DownCast( surf )->BasisSurface();
811 return applyIn2D( surf, p1, p2, & AverageUV );
814 //=======================================================================
815 //function : GetCenterUV
816 //purpose : Return UV for the central node of a biquadratic triangle
817 //=======================================================================
819 gp_XY SMESH_MesherHelper::GetCenterUV(const gp_XY& uv1,
825 bool * isBadTria/*=0*/)
828 gp_XY uvAvg = ( uv12 + uv23 + uv31 ) / 3.;
830 if (( badTria = (( uvAvg - uv1 ) * ( uvAvg - uv23 ) > 0 )))
831 uvAvg = ( uv1 + uv23 ) / 2.;
832 else if (( badTria = (( uvAvg - uv2 ) * ( uvAvg - uv31 ) > 0 )))
833 uvAvg = ( uv2 + uv31 ) / 2.;
834 else if (( badTria = (( uvAvg - uv3 ) * ( uvAvg - uv12 ) > 0 )))
835 uvAvg = ( uv3 + uv12 ) / 2.;
838 *isBadTria = badTria;
842 //=======================================================================
843 //function : GetNodeU
844 //purpose : Return node U on edge
845 //=======================================================================
847 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
848 const SMDS_MeshNode* n,
849 const SMDS_MeshNode* inEdgeNode,
852 double param = Precision::Infinite();
854 const SMDS_PositionPtr pos = n->GetPosition();
855 if ( pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
857 const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( pos );
858 param = epos->GetUParameter();
860 else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
862 if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
865 BRep_Tool::Range( E, f,l );
866 double uInEdge = GetNodeU( E, inEdgeNode );
867 param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
871 SMESHDS_Mesh * meshDS = GetMeshDS();
872 int vertexID = n->getshapeId();
873 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
874 param = BRep_Tool::Parameter( V, E );
879 double tol = BRep_Tool::Tolerance( E );
880 double f,l; BRep_Tool::Range( E, f,l );
881 bool force = ( param < f-tol || param > l+tol );
882 if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
883 force = ( GetMeshDS()->ShapeToIndex( E ) != n->getshapeId() );
885 *check = CheckNodeU( E, n, param, 2*tol, force );
890 //=======================================================================
891 //function : CheckNodeU
892 //purpose : Check and fix node U on an edge
893 // Return false if U is bad and could not be fixed
894 //=======================================================================
896 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge& E,
897 const SMDS_MeshNode* n,
901 double distXYZ[4]) const
903 int shapeID = n->getshapeId();
904 bool infinit = Precision::IsInfinite( u );
905 bool zero = ( u == 0. );
906 if ( force || toCheckPosOnShape( shapeID ) || infinit || zero )
908 TopLoc_Location loc; double f,l;
909 Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
910 if ( curve.IsNull() ) // degenerated edge
912 if ( u+tol < f || u-tol > l )
914 double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
920 gp_Pnt nodePnt = SMESH_TNodeXYZ( n );
921 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
926 curvPnt = curve->Value( u );
927 dist = nodePnt.Distance( curvPnt );
929 curvPnt.Transform( loc );
931 distXYZ[1] = curvPnt.X(); distXYZ[2] = curvPnt.Y(); distXYZ[3]=curvPnt.Z();
936 setPosOnShapeValidity( shapeID, false );
937 // u incorrect, project the node to the curve
938 int edgeID = GetMeshDS()->ShapeToIndex( E );
939 TID2ProjectorOnCurve& i2proj = const_cast< TID2ProjectorOnCurve&>( myEdge2Projector );
940 TID2ProjectorOnCurve::iterator i_proj =
941 i2proj.insert( make_pair( edgeID, (GeomAPI_ProjectPointOnCurve*) 0 )).first;
942 if ( !i_proj->second )
944 i_proj->second = new GeomAPI_ProjectPointOnCurve();
945 i_proj->second->Init( curve, f, l );
947 GeomAPI_ProjectPointOnCurve* projector = i_proj->second;
948 projector->Perform( nodePnt );
949 if ( projector->NbPoints() < 1 )
951 MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
954 Quantity_Parameter U = projector->LowerDistanceParameter();
956 MESSAGE(" f " << f << " l " << l << " u " << u);
957 curvPnt = curve->Value( u );
958 dist = nodePnt.Distance( curvPnt );
960 curvPnt.Transform( loc );
962 distXYZ[1] = curvPnt.X(); distXYZ[2] = curvPnt.Y(); distXYZ[3]=curvPnt.Z();
966 MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
967 MESSAGE("distance " << dist << " " << tol );
970 // store the fixed U on the edge
971 if ( myShape.IsSame(E) && shapeID == myShapeID && myFixNodeParameters )
972 const_cast<SMDS_MeshNode*>(n)->SetPosition
973 ( SMDS_PositionPtr( new SMDS_EdgePosition( U )));
975 else if ( fabs( u ) > numeric_limits<double>::min() )
977 setPosOnShapeValidity( shapeID, true );
979 if (( u < f-tol || u > l+tol ) && force )
981 MESSAGE("u < f-tol || u > l+tol ; u " << u << " f " << f << " l " << l);
982 // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
985 // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
986 double period = curve->Period();
987 u = ( u < f ) ? u + period : u - period;
989 catch (Standard_Failure& exc)
999 //=======================================================================
1000 //function : GetMediumPos
1001 //purpose : Return index and type of the shape (EDGE or FACE only) to
1002 // set a medium node on
1003 //param : useCurSubShape - if true, returns the shape set via SetSubShape()
1005 //=======================================================================
1007 std::pair<int, TopAbs_ShapeEnum>
1008 SMESH_MesherHelper::GetMediumPos(const SMDS_MeshNode* n1,
1009 const SMDS_MeshNode* n2,
1010 const bool useCurSubShape)
1012 if ( useCurSubShape && !myShape.IsNull() )
1013 return std::make_pair( myShapeID, myShape.ShapeType() );
1015 TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
1019 if (( myShapeID == n1->getshapeId() || myShapeID == n2->getshapeId() ) && myShapeID > 0 )
1021 shapeType = myShape.ShapeType();
1022 shapeID = myShapeID;
1024 else if ( n1->getshapeId() == n2->getshapeId() )
1026 shapeID = n2->getshapeId();
1027 shape = GetSubShapeByNode( n1, GetMeshDS() );
1031 const SMDS_TypeOfPosition Pos1 = n1->GetPosition()->GetTypeOfPosition();
1032 const SMDS_TypeOfPosition Pos2 = n2->GetPosition()->GetTypeOfPosition();
1034 if ( Pos1 == SMDS_TOP_3DSPACE || Pos2 == SMDS_TOP_3DSPACE )
1037 else if ( Pos1 == SMDS_TOP_FACE || Pos2 == SMDS_TOP_FACE )
1039 if ( Pos1 != SMDS_TOP_FACE || Pos2 != SMDS_TOP_FACE )
1041 if ( Pos1 != SMDS_TOP_FACE ) std::swap( n1,n2 );
1042 TopoDS_Shape F = GetSubShapeByNode( n1, GetMeshDS() );
1043 TopoDS_Shape S = GetSubShapeByNode( n2, GetMeshDS() );
1044 if ( IsSubShape( S, F ))
1046 shapeType = TopAbs_FACE;
1047 shapeID = n1->getshapeId();
1051 else if ( Pos1 == SMDS_TOP_EDGE && Pos2 == SMDS_TOP_EDGE )
1053 TopoDS_Shape E1 = GetSubShapeByNode( n1, GetMeshDS() );
1054 TopoDS_Shape E2 = GetSubShapeByNode( n2, GetMeshDS() );
1055 shape = GetCommonAncestor( E1, E2, *myMesh, TopAbs_FACE );
1057 else if ( Pos1 == SMDS_TOP_VERTEX && Pos2 == SMDS_TOP_VERTEX )
1059 TopoDS_Shape V1 = GetSubShapeByNode( n1, GetMeshDS() );
1060 TopoDS_Shape V2 = GetSubShapeByNode( n2, GetMeshDS() );
1061 shape = GetCommonAncestor( V1, V2, *myMesh, TopAbs_EDGE );
1062 if ( shape.IsNull() ) shape = GetCommonAncestor( V1, V2, *myMesh, TopAbs_FACE );
1064 else // VERTEX and EDGE
1066 if ( Pos1 != SMDS_TOP_VERTEX ) std::swap( n1,n2 );
1067 TopoDS_Shape V = GetSubShapeByNode( n1, GetMeshDS() );
1068 TopoDS_Shape E = GetSubShapeByNode( n2, GetMeshDS() );
1069 if ( IsSubShape( V, E ))
1072 shape = GetCommonAncestor( V, E, *myMesh, TopAbs_FACE );
1076 if ( !shape.IsNull() )
1079 shapeID = GetMeshDS()->ShapeToIndex( shape );
1080 shapeType = shape.ShapeType();
1082 return make_pair( shapeID, shapeType );
1085 //=======================================================================
1086 //function : GetCentralNode
1087 //purpose : Return existing or create a new central node for a quardilateral
1088 // quadratic face given its 8 nodes.
1089 //@param : force3d - true means node creation in between the given nodes,
1090 // else node position is found on a geometrical face if any.
1091 //=======================================================================
1093 const SMDS_MeshNode* SMESH_MesherHelper::GetCentralNode(const SMDS_MeshNode* n1,
1094 const SMDS_MeshNode* n2,
1095 const SMDS_MeshNode* n3,
1096 const SMDS_MeshNode* n4,
1097 const SMDS_MeshNode* n12,
1098 const SMDS_MeshNode* n23,
1099 const SMDS_MeshNode* n34,
1100 const SMDS_MeshNode* n41,
1103 SMDS_MeshNode *centralNode = 0; // central node to return
1105 // Find an existing central node
1107 TBiQuad keyOfMap(n1,n2,n3,n4);
1108 std::map<TBiQuad, const SMDS_MeshNode* >::iterator itMapCentralNode;
1109 itMapCentralNode = myMapWithCentralNode.find( keyOfMap );
1110 if ( itMapCentralNode != myMapWithCentralNode.end() )
1112 return (*itMapCentralNode).second;
1115 // Get type of shape for the new central node
1117 TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
1121 TopTools_ListIteratorOfListOfShape it;
1123 std::map< int, int > faceId2nbNodes;
1124 std::map< int, int > ::iterator itMapWithIdFace;
1126 SMESHDS_Mesh* meshDS = GetMeshDS();
1128 // check if a face lies on a FACE, i.e. its all corner nodes lie either on the FACE or
1129 // on sub-shapes of the FACE
1130 if ( GetMesh()->HasShapeToMesh() )
1132 const SMDS_MeshNode* nodes[] = { n1, n2, n3, n4 };
1133 for(int i = 0; i < 4; i++)
1135 shape = GetSubShapeByNode( nodes[i], meshDS );
1136 if ( shape.IsNull() ) break;
1137 if ( shape.ShapeType() == TopAbs_SOLID )
1139 solidID = nodes[i]->getshapeId();
1140 shapeType = TopAbs_SOLID;
1143 if ( shape.ShapeType() == TopAbs_FACE )
1145 faceID = nodes[i]->getshapeId();
1146 itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1147 itMapWithIdFace->second++;
1151 PShapeIteratorPtr it = GetAncestors(shape, *GetMesh(), TopAbs_FACE );
1152 while ( const TopoDS_Shape* face = it->next() )
1154 faceID = meshDS->ShapeToIndex( *face );
1155 itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1156 itMapWithIdFace->second++;
1161 if ( solidID < 1 && !faceId2nbNodes.empty() ) // SOLID not found
1163 // find ID of the FACE the four corner nodes belong to
1164 itMapWithIdFace = faceId2nbNodes.begin();
1165 for ( ; itMapWithIdFace != faceId2nbNodes.end(); ++itMapWithIdFace)
1167 if ( itMapWithIdFace->second == 4 )
1169 shapeType = TopAbs_FACE;
1170 faceID = (*itMapWithIdFace).first;
1177 if ( shapeType == TopAbs_FACE )
1179 F = TopoDS::Face( meshDS->IndexToShape( faceID ));
1186 bool toCheck = true;
1187 if ( !F.IsNull() && !force3d )
1189 uvAvg = calcTFI (0.5, 0.5,
1190 GetNodeUV(F,n1,n3,&toCheck), GetNodeUV(F,n2,n4,&toCheck),
1191 GetNodeUV(F,n3,n1,&toCheck), GetNodeUV(F,n4,n2,&toCheck),
1192 GetNodeUV(F,n12,n3), GetNodeUV(F,n23,n4),
1193 GetNodeUV(F,n34,n2), GetNodeUV(F,n41,n2));
1194 TopLoc_Location loc;
1195 Handle( Geom_Surface ) S = BRep_Tool::Surface( F, loc );
1196 P = S->Value( uvAvg.X(), uvAvg.Y() ).Transformed( loc );
1197 centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1198 // if ( mySetElemOnShape ) node is not elem!
1199 meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1201 else // ( force3d || F.IsNull() )
1203 P = calcTFI (0.5, 0.5,
1204 SMESH_TNodeXYZ(n1), SMESH_TNodeXYZ(n2),
1205 SMESH_TNodeXYZ(n3), SMESH_TNodeXYZ(n4),
1206 SMESH_TNodeXYZ(n12), SMESH_TNodeXYZ(n23),
1207 SMESH_TNodeXYZ(n34), SMESH_TNodeXYZ(n41));
1208 centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1210 if ( !F.IsNull() ) // force3d
1212 uvAvg = (GetNodeUV(F,n1,n3,&toCheck) +
1213 GetNodeUV(F,n2,n4,&toCheck) +
1214 GetNodeUV(F,n3,n1,&toCheck) +
1215 GetNodeUV(F,n4,n2,&toCheck)) / 4;
1216 //CheckNodeUV( F, centralNode, uvAvg, 2*BRep_Tool::Tolerance( F ), /*force=*/true);
1217 meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1219 else if ( solidID > 0 )
1221 meshDS->SetNodeInVolume( centralNode, solidID );
1223 else if ( myShapeID > 0 && mySetElemOnShape )
1225 meshDS->SetMeshElementOnShape( centralNode, myShapeID );
1228 myMapWithCentralNode.insert( std::make_pair( keyOfMap, centralNode ) );
1232 //=======================================================================
1233 //function : GetCentralNode
1234 //purpose : Return existing or create a new central node for a
1235 // quadratic triangle given its 6 nodes.
1236 //@param : force3d - true means node creation in between the given nodes,
1237 // else node position is found on a geometrical face if any.
1238 //=======================================================================
1240 const SMDS_MeshNode* SMESH_MesherHelper::GetCentralNode(const SMDS_MeshNode* n1,
1241 const SMDS_MeshNode* n2,
1242 const SMDS_MeshNode* n3,
1243 const SMDS_MeshNode* n12,
1244 const SMDS_MeshNode* n23,
1245 const SMDS_MeshNode* n31,
1248 SMDS_MeshNode *centralNode = 0; // central node to return
1250 // Find an existing central node
1252 TBiQuad keyOfMap(n1,n2,n3);
1253 std::map<TBiQuad, const SMDS_MeshNode* >::iterator itMapCentralNode;
1254 itMapCentralNode = myMapWithCentralNode.find( keyOfMap );
1255 if ( itMapCentralNode != myMapWithCentralNode.end() )
1257 return (*itMapCentralNode).second;
1260 // Get type of shape for the new central node
1262 TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
1266 TopTools_ListIteratorOfListOfShape it;
1268 std::map< int, int > faceId2nbNodes;
1269 std::map< int, int > ::iterator itMapWithIdFace;
1271 SMESHDS_Mesh* meshDS = GetMeshDS();
1273 // check if a face lies on a FACE, i.e. its all corner nodes lie either on the FACE or
1274 // on sub-shapes of the FACE
1275 if ( GetMesh()->HasShapeToMesh() )
1277 const SMDS_MeshNode* nodes[] = { n1, n2, n3 };
1278 for(int i = 0; i < 3; i++)
1280 shape = GetSubShapeByNode( nodes[i], meshDS );
1281 if ( shape.IsNull() ) break;
1282 if ( shape.ShapeType() == TopAbs_SOLID )
1284 solidID = nodes[i]->getshapeId();
1285 shapeType = TopAbs_SOLID;
1288 if ( shape.ShapeType() == TopAbs_FACE )
1290 faceID = nodes[i]->getshapeId();
1291 itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1292 itMapWithIdFace->second++;
1296 PShapeIteratorPtr it = GetAncestors(shape, *GetMesh(), TopAbs_FACE );
1297 while ( const TopoDS_Shape* face = it->next() )
1299 faceID = meshDS->ShapeToIndex( *face );
1300 itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1301 itMapWithIdFace->second++;
1306 if ( solidID < 1 && !faceId2nbNodes.empty() ) // SOLID not found
1308 // find ID of the FACE the four corner nodes belong to
1309 itMapWithIdFace = faceId2nbNodes.begin();
1310 for ( ; itMapWithIdFace != faceId2nbNodes.end(); ++itMapWithIdFace)
1312 if ( itMapWithIdFace->second == 3 )
1314 shapeType = TopAbs_FACE;
1315 faceID = (*itMapWithIdFace).first;
1325 if ( shapeType == TopAbs_FACE )
1327 F = TopoDS::Face( meshDS->IndexToShape( faceID ));
1329 gp_XY uv1 = GetNodeUV( F, n1, n23, &check );
1330 gp_XY uv2 = GetNodeUV( F, n2, n31, &check );
1331 gp_XY uv3 = GetNodeUV( F, n3, n12, &check );
1332 gp_XY uv12 = GetNodeUV( F, n12, n3, &check );
1333 gp_XY uv23 = GetNodeUV( F, n23, n1, &check );
1334 gp_XY uv31 = GetNodeUV( F, n31, n2, &check );
1335 uvAvg = GetCenterUV( uv1,uv2,uv3, uv12,uv23,uv31, &badTria );
1340 // Create a central node
1343 if ( !F.IsNull() && !force3d )
1345 TopLoc_Location loc;
1346 Handle( Geom_Surface ) S = BRep_Tool::Surface( F, loc );
1347 P = S->Value( uvAvg.X(), uvAvg.Y() ).Transformed( loc );
1348 centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1349 // if ( mySetElemOnShape ) node is not elem!
1350 meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1352 else // ( force3d || F.IsNull() )
1354 P = ( SMESH_TNodeXYZ( n12 ) +
1355 SMESH_TNodeXYZ( n23 ) +
1356 SMESH_TNodeXYZ( n31 ) ) / 3;
1357 centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1359 if ( !F.IsNull() ) // force3d
1361 meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1363 else if ( solidID > 0 )
1365 meshDS->SetNodeInVolume( centralNode, solidID );
1367 else if ( myShapeID > 0 && mySetElemOnShape )
1369 meshDS->SetMeshElementOnShape( centralNode, myShapeID );
1372 myMapWithCentralNode.insert( std::make_pair( keyOfMap, centralNode ) );
1376 //=======================================================================
1377 //function : GetMediumNode
1378 //purpose : Return existing or create a new medium node between given ones
1379 //=======================================================================
1381 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
1382 const SMDS_MeshNode* n2,
1385 // Find existing node
1387 SMESH_TLink link(n1,n2);
1388 ItTLinkNode itLN = myTLinkNodeMap.find( link );
1389 if ( itLN != myTLinkNodeMap.end() ) {
1390 return (*itLN).second;
1393 // Create medium node
1396 SMESHDS_Mesh* meshDS = GetMeshDS();
1398 if ( IsSeamShape( n1->getshapeId() ))
1399 // to get a correct UV of a node on seam, the second node must have checked UV
1400 std::swap( n1, n2 );
1402 // get type of shape for the new medium node
1403 int faceID = -1, edgeID = -1;
1404 TopoDS_Edge E; double u [2];
1405 TopoDS_Face F; gp_XY uv[2];
1406 bool uvOK[2] = { false, false };
1408 pair<int, TopAbs_ShapeEnum> pos = GetMediumPos( n1, n2, mySetElemOnShape );
1409 // calling GetMediumPos() with useCurSubShape=mySetElemOnShape is OK only for the
1410 // case where the lower dim mesh is already constructed, else, nodes on EDGEs are
1411 // assigned to FACE, for example.
1413 // get positions of the given nodes on shapes
1414 if ( pos.second == TopAbs_FACE )
1416 F = TopoDS::Face(meshDS->IndexToShape( faceID = pos.first ));
1417 uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
1418 uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
1420 else if ( pos.second == TopAbs_EDGE )
1422 const SMDS_PositionPtr Pos1 = n1->GetPosition();
1423 const SMDS_PositionPtr Pos2 = n2->GetPosition();
1424 if ( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE &&
1425 Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE &&
1426 n1->getshapeId() != n2->getshapeId() )
1429 return getMediumNodeOnComposedWire(n1,n2,force3d);
1431 E = TopoDS::Edge(meshDS->IndexToShape( edgeID = pos.first ));
1433 u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
1434 u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
1436 catch ( Standard_Failure& f )
1438 // issue 22502 / a node is on VERTEX not belonging to E
1439 // issue 22568 / both nodes are on non-connected VERTEXes
1440 return getMediumNodeOnComposedWire(n1,n2,force3d);
1444 if ( !force3d & uvOK[0] && uvOK[1] )
1446 // we try to create medium node using UV parameters of
1447 // nodes, else - medium between corresponding 3d points
1450 //if ( uvOK[0] && uvOK[1] )
1452 if ( IsDegenShape( n1->getshapeId() )) {
1453 if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
1454 else uv[0].SetCoord( 2, uv[1].Coord( 2 ));
1456 else if ( IsDegenShape( n2->getshapeId() )) {
1457 if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
1458 else uv[1].SetCoord( 2, uv[0].Coord( 2 ));
1461 TopLoc_Location loc;
1462 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
1463 gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
1464 gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
1465 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
1466 // if ( mySetElemOnShape ) node is not elem!
1467 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
1468 myTLinkNodeMap.insert(make_pair(link,n12));
1472 else if ( !E.IsNull() )
1475 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
1478 Standard_Boolean isPeriodic = C->IsPeriodic();
1481 Standard_Real Period = C->Period();
1482 Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
1483 Standard_Real pmid = (u[0]+p)/2.;
1484 U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
1489 gp_Pnt P = C->Value( U );
1490 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
1491 //if ( mySetElemOnShape ) node is not elem!
1492 meshDS->SetNodeOnEdge(n12, edgeID, U);
1493 myTLinkNodeMap.insert(make_pair(link,n12));
1500 double x = ( n1->X() + n2->X() )/2.;
1501 double y = ( n1->Y() + n2->Y() )/2.;
1502 double z = ( n1->Z() + n2->Z() )/2.;
1503 n12 = meshDS->AddNode(x,y,z);
1505 //if ( mySetElemOnShape ) node is not elem!
1509 gp_XY UV = ( uv[0] + uv[1] ) / 2.;
1510 CheckNodeUV( F, n12, UV, 2 * BRep_Tool::Tolerance( F ), /*force=*/true);
1511 meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
1513 else if ( !E.IsNull() )
1515 double U = ( u[0] + u[1] ) / 2.;
1516 CheckNodeU( E, n12, U, 2 * BRep_Tool::Tolerance( E ), /*force=*/true);
1517 meshDS->SetNodeOnEdge(n12, edgeID, U);
1519 else if ( myShapeID > 0 && mySetElemOnShape )
1521 meshDS->SetMeshElementOnShape(n12, myShapeID);
1525 myTLinkNodeMap.insert( make_pair( link, n12 ));
1529 //================================================================================
1531 * \brief Makes a medium node if nodes reside different edges
1533 //================================================================================
1535 const SMDS_MeshNode* SMESH_MesherHelper::getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
1536 const SMDS_MeshNode* n2,
1539 SMESH_TNodeXYZ p1( n1 ), p2( n2 );
1540 gp_Pnt middle = 0.5 * p1 + 0.5 * p2;
1541 SMDS_MeshNode* n12 = AddNode( middle.X(), middle.Y(), middle.Z() );
1543 // To find position on edge and 3D position for n12,
1544 // project <middle> to 2 edges and select projection most close to <middle>
1546 TopoDS_Edge bestEdge;
1547 double u = 0, distMiddleProj = Precision::Infinite(), distXYZ[4], f,l;
1549 // get shapes under the nodes
1550 TopoDS_Shape shape[2];
1552 for ( int is2nd = 0; is2nd < 2; ++is2nd )
1554 const SMDS_MeshNode* n = is2nd ? n2 : n1;
1555 TopoDS_Shape S = GetSubShapeByNode( n, GetMeshDS() );
1557 shape[ nbShapes++ ] = S;
1560 vector< TopoDS_Shape > edges;
1561 for ( int iS = 0; iS < nbShapes; ++iS )
1563 switch ( shape[iS].ShapeType() ) {
1566 edges.push_back( shape[iS] );
1572 if ( nbShapes == 2 && iS==0 && shape[1-iS].ShapeType() == TopAbs_VERTEX )
1573 edge = GetCommonAncestor( shape[iS], shape[1-iS], *myMesh, TopAbs_EDGE );
1575 if ( edge.IsNull() )
1577 PShapeIteratorPtr eIt = GetAncestors( shape[iS], *myMesh, TopAbs_EDGE );
1578 while( const TopoDS_Shape* e = eIt->next() )
1579 edges.push_back( *e );
1585 if ( nbShapes == 1 || shape[1-iS].ShapeType() < TopAbs_EDGE )
1586 for ( TopExp_Explorer e( shape[iS], TopAbs_EDGE ); e.More(); e.Next() )
1587 edges.push_back( e.Current() );
1594 // project to get U of projection and distance from middle to projection
1595 for ( size_t iE = 0; iE < edges.size(); ++iE )
1597 const TopoDS_Edge& edge = TopoDS::Edge( edges[ iE ]);
1598 distXYZ[0] = distMiddleProj;
1600 CheckNodeU( edge, n12, testU, 2 * BRep_Tool::Tolerance(edge), /*force=*/true, distXYZ );
1601 if ( distXYZ[0] < distMiddleProj )
1603 distMiddleProj = distXYZ[0];
1609 // // both projections failed; set n12 on the edge of n1 with U of a common vertex
1610 // TopoDS_Vertex vCommon;
1611 // if ( TopExp::CommonVertex( edges[0], edges[1], vCommon ))
1612 // u = BRep_Tool::Parameter( vCommon, edges[0] );
1615 // double f,l, u0 = GetNodeU( edges[0], n1 );
1616 // BRep_Tool::Range( edges[0],f,l );
1617 // u = ( fabs(u0-f) < fabs(u0-l) ) ? f : l;
1620 // distMiddleProj = 0;
1623 if ( !bestEdge.IsNull() )
1625 // move n12 to position of a successfull projection
1626 //double tol = BRep_Tool::Tolerance(edges[ iOkEdge ]);
1627 if ( !force3d /*&& distMiddleProj > 2*tol*/ )
1629 TopLoc_Location loc;
1630 Handle(Geom_Curve) curve = BRep_Tool::Curve( bestEdge,loc,f,l );
1631 gp_Pnt p = curve->Value( u ).Transformed( loc );
1632 GetMeshDS()->MoveNode( n12, p.X(), p.Y(), p.Z() );
1634 //if ( mySetElemOnShape ) node is not elem!
1636 int edgeID = GetMeshDS()->ShapeToIndex( bestEdge );
1637 if ( edgeID != n12->getshapeId() )
1638 GetMeshDS()->UnSetNodeOnShape( n12 );
1639 GetMeshDS()->SetNodeOnEdge(n12, edgeID, u);
1642 myTLinkNodeMap.insert( make_pair( SMESH_TLink(n1,n2), n12 ));
1647 //=======================================================================
1648 //function : AddNode
1649 //purpose : Creates a node
1650 //=======================================================================
1652 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID,
1655 SMESHDS_Mesh * meshDS = GetMeshDS();
1656 SMDS_MeshNode* node = 0;
1658 node = meshDS->AddNodeWithID( x, y, z, ID );
1660 node = meshDS->AddNode( x, y, z );
1661 if ( mySetElemOnShape && myShapeID > 0 ) { // node is not elem ?
1662 switch ( myShape.ShapeType() ) {
1663 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
1664 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
1665 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID, u, v); break;
1666 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID, u); break;
1667 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
1674 //=======================================================================
1675 //function : AddEdge
1676 //purpose : Creates quadratic or linear edge
1677 //=======================================================================
1679 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
1680 const SMDS_MeshNode* n2,
1684 SMESHDS_Mesh * meshDS = GetMeshDS();
1686 SMDS_MeshEdge* edge = 0;
1687 if (myCreateQuadratic) {
1688 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1690 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
1692 edge = meshDS->AddEdge(n1, n2, n12);
1696 edge = meshDS->AddEdgeWithID(n1, n2, id);
1698 edge = meshDS->AddEdge(n1, n2);
1701 if ( mySetElemOnShape && myShapeID > 0 )
1702 meshDS->SetMeshElementOnShape( edge, myShapeID );
1707 //=======================================================================
1708 //function : AddFace
1709 //purpose : Creates quadratic or linear triangle
1710 //=======================================================================
1712 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1713 const SMDS_MeshNode* n2,
1714 const SMDS_MeshNode* n3,
1718 SMESHDS_Mesh * meshDS = GetMeshDS();
1719 SMDS_MeshFace* elem = 0;
1721 if( n1==n2 || n2==n3 || n3==n1 )
1724 if(!myCreateQuadratic) {
1726 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
1728 elem = meshDS->AddFace(n1, n2, n3);
1731 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1732 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1733 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1734 if(myCreateBiQuadratic)
1736 const SMDS_MeshNode* nCenter = GetCentralNode(n1, n2, n3, n12, n23, n31, force3d);
1738 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, nCenter, id);
1740 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31, nCenter);
1745 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
1747 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
1750 if ( mySetElemOnShape && myShapeID > 0 )
1751 meshDS->SetMeshElementOnShape( elem, myShapeID );
1756 //=======================================================================
1757 //function : AddFace
1758 //purpose : Creates bi-quadratic, quadratic or linear quadrangle
1759 //=======================================================================
1761 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1762 const SMDS_MeshNode* n2,
1763 const SMDS_MeshNode* n3,
1764 const SMDS_MeshNode* n4,
1768 SMESHDS_Mesh * meshDS = GetMeshDS();
1769 SMDS_MeshFace* elem = 0;
1772 return AddFace(n1,n3,n4,id,force3d);
1775 return AddFace(n1,n2,n4,id,force3d);
1778 return AddFace(n1,n2,n3,id,force3d);
1781 return AddFace(n1,n2,n4,id,force3d);
1784 return AddFace(n1,n2,n3,id,force3d);
1787 return AddFace(n1,n2,n3,id,force3d);
1790 if(!myCreateQuadratic) {
1792 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
1794 elem = meshDS->AddFace(n1, n2, n3, n4);
1797 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1798 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1799 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1800 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1801 if(myCreateBiQuadratic)
1803 const SMDS_MeshNode* nCenter = GetCentralNode(n1, n2, n3, n4, n12, n23, n34, n41, force3d);
1805 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, nCenter, id);
1807 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41, nCenter);
1812 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
1814 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
1817 if ( mySetElemOnShape && myShapeID > 0 )
1818 meshDS->SetMeshElementOnShape( elem, myShapeID );
1823 //=======================================================================
1824 //function : AddPolygonalFace
1825 //purpose : Creates polygon, with additional nodes in quadratic mesh
1826 //=======================================================================
1828 SMDS_MeshFace* SMESH_MesherHelper::AddPolygonalFace (const vector<const SMDS_MeshNode*>& nodes,
1832 SMESHDS_Mesh * meshDS = GetMeshDS();
1833 SMDS_MeshFace* elem = 0;
1835 if(!myCreateQuadratic) {
1837 elem = meshDS->AddPolygonalFaceWithID(nodes, id);
1839 elem = meshDS->AddPolygonalFace(nodes);
1842 vector<const SMDS_MeshNode*> newNodes;
1843 for ( int i = 0; i < nodes.size(); ++i )
1845 const SMDS_MeshNode* n1 = nodes[i];
1846 const SMDS_MeshNode* n2 = nodes[(i+1)%nodes.size()];
1847 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1848 newNodes.push_back( n1 );
1849 newNodes.push_back( n12 );
1852 elem = meshDS->AddPolygonalFaceWithID(newNodes, id);
1854 elem = meshDS->AddPolygonalFace(newNodes);
1856 if ( mySetElemOnShape && myShapeID > 0 )
1857 meshDS->SetMeshElementOnShape( elem, myShapeID );
1862 //=======================================================================
1863 //function : AddVolume
1864 //purpose : Creates quadratic or linear prism
1865 //=======================================================================
1867 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1868 const SMDS_MeshNode* n2,
1869 const SMDS_MeshNode* n3,
1870 const SMDS_MeshNode* n4,
1871 const SMDS_MeshNode* n5,
1872 const SMDS_MeshNode* n6,
1876 SMESHDS_Mesh * meshDS = GetMeshDS();
1877 SMDS_MeshVolume* elem = 0;
1878 if(!myCreateQuadratic) {
1880 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
1882 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
1885 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1886 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1887 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1889 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1890 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1891 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1893 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1894 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1895 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1898 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
1899 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1901 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1902 n12, n23, n31, n45, n56, n64, n14, n25, n36);
1904 if ( mySetElemOnShape && myShapeID > 0 )
1905 meshDS->SetMeshElementOnShape( elem, myShapeID );
1910 //=======================================================================
1911 //function : AddVolume
1912 //purpose : Creates quadratic or linear tetrahedron
1913 //=======================================================================
1915 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1916 const SMDS_MeshNode* n2,
1917 const SMDS_MeshNode* n3,
1918 const SMDS_MeshNode* n4,
1922 SMESHDS_Mesh * meshDS = GetMeshDS();
1923 SMDS_MeshVolume* elem = 0;
1924 if(!myCreateQuadratic) {
1926 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1928 elem = meshDS->AddVolume(n1, n2, n3, n4);
1931 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1932 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1933 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1935 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1936 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1937 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1940 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1942 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1944 if ( mySetElemOnShape && myShapeID > 0 )
1945 meshDS->SetMeshElementOnShape( elem, myShapeID );
1950 //=======================================================================
1951 //function : AddVolume
1952 //purpose : Creates quadratic or linear pyramid
1953 //=======================================================================
1955 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1956 const SMDS_MeshNode* n2,
1957 const SMDS_MeshNode* n3,
1958 const SMDS_MeshNode* n4,
1959 const SMDS_MeshNode* n5,
1963 SMDS_MeshVolume* elem = 0;
1964 if(!myCreateQuadratic) {
1966 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1968 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1971 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1972 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1973 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1974 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1976 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1977 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1978 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1979 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1982 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
1987 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
1989 n15, n25, n35, n45);
1991 if ( mySetElemOnShape && myShapeID > 0 )
1992 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1997 //=======================================================================
1998 //function : AddVolume
1999 //purpose : Creates bi-quadratic, quadratic or linear hexahedron
2000 //=======================================================================
2002 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2003 const SMDS_MeshNode* n2,
2004 const SMDS_MeshNode* n3,
2005 const SMDS_MeshNode* n4,
2006 const SMDS_MeshNode* n5,
2007 const SMDS_MeshNode* n6,
2008 const SMDS_MeshNode* n7,
2009 const SMDS_MeshNode* n8,
2013 SMESHDS_Mesh * meshDS = GetMeshDS();
2014 SMDS_MeshVolume* elem = 0;
2015 if(!myCreateQuadratic) {
2017 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
2019 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
2022 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
2023 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
2024 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
2025 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
2027 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
2028 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
2029 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
2030 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
2032 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
2033 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
2034 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
2035 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
2036 if(myCreateBiQuadratic)
2038 const SMDS_MeshNode* n1234 = GetCentralNode(n1,n2,n3,n4,n12,n23,n34,n41,force3d);
2039 const SMDS_MeshNode* n1256 = GetCentralNode(n1,n2,n5,n6,n12,n26,n56,n15,force3d);
2040 const SMDS_MeshNode* n2367 = GetCentralNode(n2,n3,n6,n7,n23,n37,n67,n26,force3d);
2041 const SMDS_MeshNode* n3478 = GetCentralNode(n3,n4,n7,n8,n34,n48,n78,n37,force3d);
2042 const SMDS_MeshNode* n1458 = GetCentralNode(n1,n4,n5,n8,n41,n48,n15,n85,force3d);
2043 const SMDS_MeshNode* n5678 = GetCentralNode(n5,n6,n7,n8,n56,n67,n78,n85,force3d);
2045 vector<gp_XYZ> pointsOnShapes( SMESH_Block::ID_Shell );
2047 pointsOnShapes[ SMESH_Block::ID_V000 ] = SMESH_TNodeXYZ( n4 );
2048 pointsOnShapes[ SMESH_Block::ID_V100 ] = SMESH_TNodeXYZ( n8 );
2049 pointsOnShapes[ SMESH_Block::ID_V010 ] = SMESH_TNodeXYZ( n3 );
2050 pointsOnShapes[ SMESH_Block::ID_V110 ] = SMESH_TNodeXYZ( n7 );
2051 pointsOnShapes[ SMESH_Block::ID_V001 ] = SMESH_TNodeXYZ( n1 );
2052 pointsOnShapes[ SMESH_Block::ID_V101 ] = SMESH_TNodeXYZ( n5 );
2053 pointsOnShapes[ SMESH_Block::ID_V011 ] = SMESH_TNodeXYZ( n2 );
2054 pointsOnShapes[ SMESH_Block::ID_V111 ] = SMESH_TNodeXYZ( n6 );
2056 pointsOnShapes[ SMESH_Block::ID_Ex00 ] = SMESH_TNodeXYZ( n48 );
2057 pointsOnShapes[ SMESH_Block::ID_Ex10 ] = SMESH_TNodeXYZ( n37 );
2058 pointsOnShapes[ SMESH_Block::ID_E0y0 ] = SMESH_TNodeXYZ( n15 );
2059 pointsOnShapes[ SMESH_Block::ID_E1y0 ] = SMESH_TNodeXYZ( n26 );
2060 pointsOnShapes[ SMESH_Block::ID_Ex01 ] = SMESH_TNodeXYZ( n34 );
2061 pointsOnShapes[ SMESH_Block::ID_Ex11 ] = SMESH_TNodeXYZ( n78 );
2062 pointsOnShapes[ SMESH_Block::ID_E0y1 ] = SMESH_TNodeXYZ( n12 );
2063 pointsOnShapes[ SMESH_Block::ID_E1y1 ] = SMESH_TNodeXYZ( n56 );
2064 pointsOnShapes[ SMESH_Block::ID_E00z ] = SMESH_TNodeXYZ( n41 );
2065 pointsOnShapes[ SMESH_Block::ID_E10z ] = SMESH_TNodeXYZ( n85 );
2066 pointsOnShapes[ SMESH_Block::ID_E01z ] = SMESH_TNodeXYZ( n23 );
2067 pointsOnShapes[ SMESH_Block::ID_E11z ] = SMESH_TNodeXYZ( n67 );
2069 pointsOnShapes[ SMESH_Block::ID_Fxy0 ] = SMESH_TNodeXYZ( n3478 );
2070 pointsOnShapes[ SMESH_Block::ID_Fxy1 ] = SMESH_TNodeXYZ( n1256 );
2071 pointsOnShapes[ SMESH_Block::ID_Fx0z ] = SMESH_TNodeXYZ( n1458 );
2072 pointsOnShapes[ SMESH_Block::ID_Fx1z ] = SMESH_TNodeXYZ( n2367 );
2073 pointsOnShapes[ SMESH_Block::ID_F0yz ] = SMESH_TNodeXYZ( n1234 );
2074 pointsOnShapes[ SMESH_Block::ID_F1yz ] = SMESH_TNodeXYZ( n5678 );
2076 gp_XYZ centerCube(0.5, 0.5, 0.5);
2078 SMESH_Block::ShellPoint( centerCube, pointsOnShapes, nCenterElem );
2079 const SMDS_MeshNode* nCenter =
2080 meshDS->AddNode( nCenterElem.X(), nCenterElem.Y(), nCenterElem.Z() );
2081 meshDS->SetNodeInVolume( nCenter, myShapeID );
2084 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
2085 n12, n23, n34, n41, n56, n67,
2086 n78, n85, n15, n26, n37, n48,
2087 n1234, n1256, n2367, n3478, n1458, n5678, nCenter, id);
2089 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
2090 n12, n23, n34, n41, n56, n67,
2091 n78, n85, n15, n26, n37, n48,
2092 n1234, n1256, n2367, n3478, n1458, n5678, nCenter);
2097 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
2098 n12, n23, n34, n41, n56, n67,
2099 n78, n85, n15, n26, n37, n48, id);
2101 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
2102 n12, n23, n34, n41, n56, n67,
2103 n78, n85, n15, n26, n37, n48);
2106 if ( mySetElemOnShape && myShapeID > 0 )
2107 meshDS->SetMeshElementOnShape( elem, myShapeID );
2112 //=======================================================================
2113 //function : AddVolume
2114 //purpose : Creates LINEAR!!!!!!!!! octahedron
2115 //=======================================================================
2117 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2118 const SMDS_MeshNode* n2,
2119 const SMDS_MeshNode* n3,
2120 const SMDS_MeshNode* n4,
2121 const SMDS_MeshNode* n5,
2122 const SMDS_MeshNode* n6,
2123 const SMDS_MeshNode* n7,
2124 const SMDS_MeshNode* n8,
2125 const SMDS_MeshNode* n9,
2126 const SMDS_MeshNode* n10,
2127 const SMDS_MeshNode* n11,
2128 const SMDS_MeshNode* n12,
2132 SMESHDS_Mesh * meshDS = GetMeshDS();
2133 SMDS_MeshVolume* elem = 0;
2135 elem = meshDS->AddVolumeWithID(n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,id);
2137 elem = meshDS->AddVolume(n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12);
2138 if ( mySetElemOnShape && myShapeID > 0 )
2139 meshDS->SetMeshElementOnShape( elem, myShapeID );
2143 //=======================================================================
2144 //function : AddPolyhedralVolume
2145 //purpose : Creates polyhedron. In quadratic mesh, adds medium nodes
2146 //=======================================================================
2149 SMESH_MesherHelper::AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
2150 const std::vector<int>& quantities,
2154 SMESHDS_Mesh * meshDS = GetMeshDS();
2155 SMDS_MeshVolume* elem = 0;
2156 if(!myCreateQuadratic)
2159 elem = meshDS->AddPolyhedralVolumeWithID(nodes, quantities, id);
2161 elem = meshDS->AddPolyhedralVolume(nodes, quantities);
2165 vector<const SMDS_MeshNode*> newNodes;
2166 vector<int> newQuantities;
2167 for ( int iFace=0, iN=0; iFace < quantities.size(); ++iFace)
2169 int nbNodesInFace = quantities[iFace];
2170 newQuantities.push_back(0);
2171 for ( int i = 0; i < nbNodesInFace; ++i )
2173 const SMDS_MeshNode* n1 = nodes[ iN + i ];
2174 newNodes.push_back( n1 );
2175 newQuantities.back()++;
2177 const SMDS_MeshNode* n2 = nodes[ iN + ( i+1==nbNodesInFace ? 0 : i+1 )];
2178 // if ( n1->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE &&
2179 // n2->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE )
2181 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
2182 newNodes.push_back( n12 );
2183 newQuantities.back()++;
2186 iN += nbNodesInFace;
2189 elem = meshDS->AddPolyhedralVolumeWithID( newNodes, newQuantities, id );
2191 elem = meshDS->AddPolyhedralVolume( newNodes, newQuantities );
2193 if ( mySetElemOnShape && myShapeID > 0 )
2194 meshDS->SetMeshElementOnShape( elem, myShapeID );
2201 //================================================================================
2203 * \brief Check if a node belongs to any face of sub-mesh
2205 //================================================================================
2207 bool isNodeInSubMesh( const SMDS_MeshNode* n, const SMESHDS_SubMesh* sm )
2209 SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator( SMDSAbs_Face );
2210 while ( fIt->more() )
2211 if ( sm->Contains( fIt->next() ))
2217 //=======================================================================
2218 //function : IsSameElemGeometry
2219 //purpose : Returns true if all elements of a sub-mesh are of same shape
2220 //=======================================================================
2222 bool SMESH_MesherHelper::IsSameElemGeometry(const SMESHDS_SubMesh* smDS,
2223 SMDSAbs_GeometryType shape,
2224 const bool nullSubMeshRes)
2226 if ( !smDS ) return nullSubMeshRes;
2228 SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
2229 while ( elemIt->more() ) {
2230 const SMDS_MeshElement* e = elemIt->next();
2231 if ( e->GetGeomType() != shape )
2237 //=======================================================================
2238 //function : LoadNodeColumns
2239 //purpose : Load nodes bound to face into a map of node columns
2240 //=======================================================================
2242 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
2243 const TopoDS_Face& theFace,
2244 const TopoDS_Edge& theBaseEdge,
2245 SMESHDS_Mesh* theMesh,
2246 SMESH_ProxyMesh* theProxyMesh)
2248 return LoadNodeColumns(theParam2ColumnMap,
2250 std::list<TopoDS_Edge>(1,theBaseEdge),
2255 //=======================================================================
2256 //function : LoadNodeColumns
2257 //purpose : Load nodes bound to face into a map of node columns
2258 //=======================================================================
2260 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
2261 const TopoDS_Face& theFace,
2262 const std::list<TopoDS_Edge>& theBaseSide,
2263 SMESHDS_Mesh* theMesh,
2264 SMESH_ProxyMesh* theProxyMesh)
2266 // get a right sub-mesh of theFace
2268 const SMESHDS_SubMesh* faceSubMesh = 0;
2271 faceSubMesh = theProxyMesh->GetSubMesh( theFace );
2272 if ( !faceSubMesh ||
2273 faceSubMesh->NbElements() == 0 ||
2274 theProxyMesh->IsTemporary( faceSubMesh->GetElements()->next() ))
2276 // can use a proxy sub-mesh with not temporary elements only
2282 faceSubMesh = theMesh->MeshElements( theFace );
2283 if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
2286 if ( theParam2ColumnMap.empty() )
2288 // get data of edges for normalization of params
2289 vector< double > length;
2291 list<TopoDS_Edge>::const_iterator edge;
2293 for ( edge = theBaseSide.begin(); edge != theBaseSide.end(); ++edge )
2295 double len = std::max( 1e-10, SMESH_Algo::EdgeLength( *edge ));
2297 length.push_back( len );
2301 // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
2302 edge = theBaseSide.begin();
2303 for ( int iE = 0; edge != theBaseSide.end(); ++edge, ++iE )
2305 map< double, const SMDS_MeshNode*> sortedBaseNN;
2306 SMESH_Algo::GetSortedNodesOnEdge( theMesh, *edge,/*noMedium=*/true, sortedBaseNN);
2307 if ( sortedBaseNN.empty() ) continue;
2309 map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNN.begin();
2310 if ( theProxyMesh ) // from sortedBaseNN remove nodes not shared by faces of faceSubMesh
2312 const SMDS_MeshNode* n1 = (++sortedBaseNN.begin())->second;
2313 const SMDS_MeshNode* n2 = (++sortedBaseNN.rbegin())->second;
2314 bool allNodesAreProxy = ( n1 != theProxyMesh->GetProxyNode( n1 ) &&
2315 n2 != theProxyMesh->GetProxyNode( n2 ));
2316 if ( allNodesAreProxy )
2317 for ( u_n = sortedBaseNN.begin(); u_n != sortedBaseNN.end(); u_n++ )
2318 u_n->second = theProxyMesh->GetProxyNode( u_n->second );
2320 if ( u_n = sortedBaseNN.begin(), !isNodeInSubMesh( u_n->second, faceSubMesh ))
2322 while ( ++u_n != sortedBaseNN.end() && !isNodeInSubMesh( u_n->second, faceSubMesh ));
2323 sortedBaseNN.erase( sortedBaseNN.begin(), u_n );
2325 if ( u_n = --sortedBaseNN.end(), !isNodeInSubMesh( u_n->second, faceSubMesh ))
2327 while ( u_n != sortedBaseNN.begin() && !isNodeInSubMesh( (--u_n)->second, faceSubMesh ));
2328 sortedBaseNN.erase( ++u_n, sortedBaseNN.end() );
2330 if ( sortedBaseNN.empty() ) continue;
2334 BRep_Tool::Range( *edge, f, l );
2335 if ( edge->Orientation() == TopAbs_REVERSED ) std::swap( f, l );
2336 const double coeff = 1. / ( l - f ) * length[iE] / fullLen;
2337 const double prevPar = theParam2ColumnMap.empty() ? 0 : theParam2ColumnMap.rbegin()->first;
2338 for ( u_n = sortedBaseNN.begin(); u_n != sortedBaseNN.end(); u_n++ )
2340 double par = prevPar + coeff * ( u_n->first - f );
2341 TParam2ColumnMap::iterator u2nn =
2342 theParam2ColumnMap.insert( theParam2ColumnMap.end(), make_pair( par, TNodeColumn()));
2343 u2nn->second.push_back( u_n->second );
2346 if ( theParam2ColumnMap.empty() )
2351 int prevNbRows = theParam2ColumnMap.begin()->second.size(); // current, at least 1 here
2352 int expectedNbRows = faceSubMesh->NbElements() / ( theParam2ColumnMap.size()-1 ); // to be added
2354 // fill theParam2ColumnMap column by column by passing from nodes on
2355 // theBaseEdge up via mesh faces on theFace
2357 TParam2ColumnMap::iterator par_nVec_1, par_nVec_2;
2358 par_nVec_2 = theParam2ColumnMap.begin();
2359 par_nVec_1 = par_nVec_2++;
2360 TIDSortedElemSet emptySet, avoidSet;
2361 for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
2363 vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
2364 vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
2365 nCol1.resize( prevNbRows + expectedNbRows );
2366 nCol2.resize( prevNbRows + expectedNbRows );
2368 int i1, i2, foundNbRows = 0;
2369 const SMDS_MeshNode *n1 = nCol1[ prevNbRows-1 ];
2370 const SMDS_MeshNode *n2 = nCol2[ prevNbRows-1 ];
2371 // find face sharing node n1 and n2 and belonging to faceSubMesh
2372 while ( const SMDS_MeshElement* face =
2373 SMESH_MeshAlgos::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
2375 if ( faceSubMesh->Contains( face ))
2377 int nbNodes = face->NbCornerNodes();
2380 if ( foundNbRows + 1 > expectedNbRows )
2382 n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
2383 n2 = face->GetNode( (i1+2) % 4 );
2384 nCol1[ prevNbRows + foundNbRows] = n1;
2385 nCol2[ prevNbRows + foundNbRows] = n2;
2388 avoidSet.insert( face );
2390 if ( foundNbRows != expectedNbRows )
2394 return ( theParam2ColumnMap.size() > 1 &&
2395 theParam2ColumnMap.begin()->second.size() == prevNbRows + expectedNbRows );
2400 //================================================================================
2402 * \brief Return true if a node is at a corner of a 2D structured mesh of FACE
2404 //================================================================================
2406 bool isCornerOfStructure( const SMDS_MeshNode* n,
2407 const SMESHDS_SubMesh* faceSM,
2408 SMESH_MesherHelper& faceAnalyser )
2410 int nbFacesInSM = 0;
2412 SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator( SMDSAbs_Face );
2413 while ( fIt->more() )
2414 nbFacesInSM += faceSM->Contains( fIt->next() );
2416 if ( nbFacesInSM == 1 )
2419 if ( nbFacesInSM == 2 && n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
2421 return faceAnalyser.IsRealSeam( n->getshapeId() );
2427 //=======================================================================
2428 //function : IsStructured
2429 //purpose : Return true if 2D mesh on FACE is structured
2430 //=======================================================================
2432 bool SMESH_MesherHelper::IsStructured( SMESH_subMesh* faceSM )
2434 SMESHDS_SubMesh* fSM = faceSM->GetSubMeshDS();
2435 if ( !fSM || fSM->NbElements() == 0 )
2438 list< TopoDS_Edge > edges;
2439 list< int > nbEdgesInWires;
2440 int nbWires = SMESH_Block::GetOrderedEdges( TopoDS::Face( faceSM->GetSubShape() ),
2441 edges, nbEdgesInWires );
2442 if ( nbWires != 1 /*|| nbEdgesInWires.front() != 4*/ ) // allow composite sides
2445 // algo: find corners of a structure and then analyze nb of faces and
2446 // length of structure sides
2448 SMESHDS_Mesh* meshDS = faceSM->GetFather()->GetMeshDS();
2449 SMESH_MesherHelper faceAnalyser( *faceSM->GetFather() );
2450 faceAnalyser.SetSubShape( faceSM->GetSubShape() );
2452 // rotate edges to get the first node being at corner
2453 // (in principle it's not necessary but so far none SALOME algo can make
2454 // such a structured mesh that all corner nodes are not on VERTEXes)
2455 bool isCorner = false;
2456 int nbRemainEdges = nbEdgesInWires.front();
2458 TopoDS_Vertex V = IthVertex( 0, edges.front() );
2459 isCorner = isCornerOfStructure( SMESH_Algo::VertexNode( V, meshDS ),
2462 edges.splice( edges.end(), edges, edges.begin() );
2466 while ( !isCorner && nbRemainEdges > 0 );
2471 // get all nodes from EDGEs
2472 list< const SMDS_MeshNode* > nodes;
2473 list< TopoDS_Edge >::iterator edge = edges.begin();
2474 for ( ; edge != edges.end(); ++edge )
2476 map< double, const SMDS_MeshNode* > u2Nodes;
2477 if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, *edge,
2478 /*skipMedium=*/true, u2Nodes ))
2481 list< const SMDS_MeshNode* > edgeNodes;
2482 map< double, const SMDS_MeshNode* >::iterator u2n = u2Nodes.begin();
2483 for ( ; u2n != u2Nodes.end(); ++u2n )
2484 edgeNodes.push_back( u2n->second );
2485 if ( edge->Orientation() == TopAbs_REVERSED )
2486 edgeNodes.reverse();
2488 if ( !nodes.empty() && nodes.back() == edgeNodes.front() )
2489 edgeNodes.pop_front();
2490 nodes.splice( nodes.end(), edgeNodes, edgeNodes.begin(), edgeNodes.end() );
2493 // get length of structured sides
2494 vector<int> nbEdgesInSide;
2496 list< const SMDS_MeshNode* >::iterator n = ++nodes.begin();
2497 for ( ; n != nodes.end(); ++n )
2500 if ( isCornerOfStructure( *n, fSM, faceAnalyser )) {
2501 nbEdgesInSide.push_back( nbEdges );
2507 if ( nbEdgesInSide.size() != 4 )
2509 if ( nbEdgesInSide[0] != nbEdgesInSide[2] )
2511 if ( nbEdgesInSide[1] != nbEdgesInSide[3] )
2513 if ( nbEdgesInSide[0] * nbEdgesInSide[1] != fSM->NbElements() )
2519 //================================================================================
2521 * \brief Find out elements orientation on a geometrical face
2522 * \param theFace - The face correctly oriented in the shape being meshed
2523 * \retval bool - true if the face normal and the normal of first element
2524 * in the correspoding submesh point in different directions
2526 //================================================================================
2528 bool SMESH_MesherHelper::IsReversedSubMesh (const TopoDS_Face& theFace)
2530 if ( theFace.IsNull() )
2533 // find out orientation of a meshed face
2534 int faceID = GetMeshDS()->ShapeToIndex( theFace );
2535 TopoDS_Shape aMeshedFace = GetMeshDS()->IndexToShape( faceID );
2536 bool isReversed = ( theFace.Orientation() != aMeshedFace.Orientation() );
2538 const SMESHDS_SubMesh * aSubMeshDSFace = GetMeshDS()->MeshElements( faceID );
2539 if ( !aSubMeshDSFace )
2542 // find an element with a good normal
2544 bool normalOK = false;
2546 SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
2547 while ( !normalOK && iteratorElem->more() ) // loop on elements on theFace
2549 const SMDS_MeshElement* elem = iteratorElem->next();
2550 if ( elem && elem->NbCornerNodes() > 2 )
2552 SMESH_TNodeXYZ nPnt[3];
2553 SMDS_ElemIteratorPtr nodesIt = elem->nodesIterator();
2554 int iNodeOnFace = 0, iPosDim = SMDS_TOP_VERTEX;
2555 for ( int iN = 0; nodesIt->more() && iN < 3; ++iN) // loop on nodes
2557 nPnt[ iN ] = nodesIt->next();
2558 if ( nPnt[ iN ]._node->GetPosition()->GetTypeOfPosition() > iPosDim )
2561 iPosDim = nPnt[ iN ]._node->GetPosition()->GetTypeOfPosition();
2565 gp_Vec v01( nPnt[0], nPnt[1] ), v02( nPnt[0], nPnt[2] );
2566 if ( v01.SquareMagnitude() > RealSmall() &&
2567 v02.SquareMagnitude() > RealSmall() )
2570 if (( normalOK = ( Ne.SquareMagnitude() > RealSmall() )))
2571 uv = GetNodeUV( theFace, nPnt[iNodeOnFace]._node, 0, &normalOK );
2578 // face normal at node position
2579 TopLoc_Location loc;
2580 Handle(Geom_Surface) surf = BRep_Tool::Surface( theFace, loc );
2581 // if ( surf.IsNull() || surf->Continuity() < GeomAbs_C1 )
2582 // some surfaces not detected as GeomAbs_C1 are nevertheless correct for meshing
2583 if ( surf.IsNull() || surf->Continuity() < GeomAbs_C0 )
2586 MESSAGE("surf->Continuity() < GeomAbs_C1 " << (surf->Continuity() < GeomAbs_C1));
2589 gp_Vec d1u, d1v; gp_Pnt p;
2590 surf->D1( uv.X(), uv.Y(), p, d1u, d1v );
2591 gp_Vec Nf = (d1u ^ d1v).Transformed( loc );
2593 if ( theFace.Orientation() == TopAbs_REVERSED )
2596 return Ne * Nf < 0.;
2599 //=======================================================================
2601 //purpose : Count nb of sub-shapes
2602 //=======================================================================
2604 int SMESH_MesherHelper::Count(const TopoDS_Shape& shape,
2605 const TopAbs_ShapeEnum type,
2606 const bool ignoreSame)
2609 TopTools_IndexedMapOfShape map;
2610 TopExp::MapShapes( shape, type, map );
2611 return map.Extent();
2615 for ( TopExp_Explorer exp( shape, type ); exp.More(); exp.Next() )
2621 //=======================================================================
2622 //function : NbAncestors
2623 //purpose : Return number of unique ancestors of the shape
2624 //=======================================================================
2626 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
2627 const SMESH_Mesh& mesh,
2628 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
2630 TopTools_MapOfShape ancestors;
2631 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
2632 for ( ; ansIt.More(); ansIt.Next() ) {
2633 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
2634 ancestors.Add( ansIt.Value() );
2636 return ancestors.Extent();
2639 //=======================================================================
2640 //function : GetSubShapeOri
2641 //purpose : Return orientation of sub-shape in the main shape
2642 //=======================================================================
2644 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
2645 const TopoDS_Shape& subShape)
2647 TopAbs_Orientation ori = TopAbs_Orientation(-1);
2648 if ( !shape.IsNull() && !subShape.IsNull() )
2650 TopExp_Explorer e( shape, subShape.ShapeType() );
2651 if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
2652 e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
2653 for ( ; e.More(); e.Next())
2654 if ( subShape.IsSame( e.Current() ))
2657 ori = e.Current().Orientation();
2662 //=======================================================================
2663 //function : IsSubShape
2665 //=======================================================================
2667 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
2668 const TopoDS_Shape& mainShape )
2670 if ( !shape.IsNull() && !mainShape.IsNull() )
2672 for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
2675 if ( shape.IsSame( exp.Current() ))
2678 SCRUTE((shape.IsNull()));
2679 SCRUTE((mainShape.IsNull()));
2683 //=======================================================================
2684 //function : IsSubShape
2686 //=======================================================================
2688 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
2690 if ( shape.IsNull() || !aMesh )
2693 aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
2695 (shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape ));
2698 //================================================================================
2700 * \brief Return maximal tolerance of shape
2702 //================================================================================
2704 double SMESH_MesherHelper::MaxTolerance( const TopoDS_Shape& shape )
2706 double tol = Precision::Confusion();
2707 TopExp_Explorer exp;
2708 for ( exp.Init( shape, TopAbs_FACE ); exp.More(); exp.Next() )
2709 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Face( exp.Current())));
2710 for ( exp.Init( shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2711 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Edge( exp.Current())));
2712 for ( exp.Init( shape, TopAbs_VERTEX ); exp.More(); exp.Next() )
2713 tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Vertex( exp.Current())));
2718 //================================================================================
2720 * \brief Return an angle between two EDGEs sharing a common VERTEX with reference
2721 * of the FACE normal
2722 * \return double - the angle (between -Pi and Pi), negative if the angle is concave,
2723 * 1e100 in case of failure
2724 * \waring Care about order of the EDGEs and their orientation to be as they are
2725 * within the FACE! Don't pass degenerated EDGEs neither!
2727 //================================================================================
2729 double SMESH_MesherHelper::GetAngle( const TopoDS_Edge & theE1,
2730 const TopoDS_Edge & theE2,
2731 const TopoDS_Face & theFace)
2733 double angle = 1e100;
2736 TopoDS_Vertex vCommon;
2737 if ( !TopExp::CommonVertex( theE1, theE2, vCommon ))
2740 Handle(Geom_Curve) c1 = BRep_Tool::Curve( theE1, f,l );
2741 Handle(Geom_Curve) c2 = BRep_Tool::Curve( theE2, f,l );
2742 Handle(Geom2d_Curve) c2d1 = BRep_Tool::CurveOnSurface( theE1, theFace, f,l );
2743 Handle(Geom_Surface) surf = BRep_Tool::Surface( theFace );
2744 double p1 = BRep_Tool::Parameter( vCommon, theE1 );
2745 double p2 = BRep_Tool::Parameter( vCommon, theE2 );
2746 if ( c1.IsNull() || c2.IsNull() )
2748 gp_XY uv = c2d1->Value( p1 ).XY();
2749 gp_Vec du, dv; gp_Pnt p;
2750 surf->D1( uv.X(), uv.Y(), p, du, dv );
2751 gp_Vec vec1, vec2, vecRef = du ^ dv;
2754 while ( vecRef.SquareMagnitude() < std::numeric_limits<double>::min() )
2756 double dp = ( l - f ) / 1000.;
2757 p1tmp += dp * (( Abs( p1 - f ) > Abs( p1 - l )) ? +1. : -1.);
2758 uv = c2d1->Value( p1tmp ).XY();
2759 surf->D1( uv.X(), uv.Y(), p, du, dv );
2761 if ( ++nbLoops > 10 )
2764 cout << "SMESH_MesherHelper::GetAngle(): Captured in a sigularity" << endl;
2769 if ( theFace.Orientation() == TopAbs_REVERSED )
2771 c1->D1( p1, p, vec1 );
2772 c2->D1( p2, p, vec2 );
2773 TopoDS_Face F = theFace;
2774 if ( F.Orientation() == TopAbs_INTERNAL )
2775 F.Orientation( TopAbs_FORWARD );
2776 if ( theE1.Orientation() /*GetSubShapeOri( F, theE1 )*/ == TopAbs_REVERSED )
2778 if ( theE2.Orientation() /*GetSubShapeOri( F, theE2 )*/ == TopAbs_REVERSED )
2780 angle = vec1.AngleWithRef( vec2, vecRef );
2788 //================================================================================
2790 * \brief Check if the first and last vertices of an edge are the same
2791 * \param anEdge - the edge to check
2792 * \retval bool - true if same
2794 //================================================================================
2796 bool SMESH_MesherHelper::IsClosedEdge( const TopoDS_Edge& anEdge )
2798 if ( anEdge.Orientation() >= TopAbs_INTERNAL )
2799 return IsClosedEdge( TopoDS::Edge( anEdge.Oriented( TopAbs_FORWARD )));
2800 return TopExp::FirstVertex( anEdge ).IsSame( TopExp::LastVertex( anEdge ));
2803 //================================================================================
2805 * \brief Wrapper over TopExp::FirstVertex() and TopExp::LastVertex() fixing them
2806 * in the case of INTERNAL edge
2808 //================================================================================
2810 TopoDS_Vertex SMESH_MesherHelper::IthVertex( const bool is2nd,
2814 if ( anEdge.Orientation() >= TopAbs_INTERNAL )
2815 anEdge.Orientation( TopAbs_FORWARD );
2817 const TopAbs_Orientation tgtOri = is2nd ? TopAbs_REVERSED : TopAbs_FORWARD;
2818 TopoDS_Iterator vIt( anEdge, CumOri );
2819 while ( vIt.More() && vIt.Value().Orientation() != tgtOri )
2822 return ( vIt.More() ? TopoDS::Vertex(vIt.Value()) : TopoDS_Vertex() );
2825 //================================================================================
2827 * \brief Return type of shape contained in a group
2828 * \param group - a shape of type TopAbs_COMPOUND
2829 * \param avoidCompound - not to return TopAbs_COMPOUND
2831 //================================================================================
2833 TopAbs_ShapeEnum SMESH_MesherHelper::GetGroupType(const TopoDS_Shape& group,
2834 const bool avoidCompound)
2836 if ( !group.IsNull() )
2838 if ( group.ShapeType() != TopAbs_COMPOUND )
2839 return group.ShapeType();
2841 // iterate on a compound
2842 TopoDS_Iterator it( group );
2844 return avoidCompound ? GetGroupType( it.Value() ) : it.Value().ShapeType();
2846 return TopAbs_SHAPE;
2849 //=======================================================================
2850 //function : IsQuadraticMesh
2851 //purpose : Check mesh without geometry for: if all elements on this shape are quadratic,
2852 // quadratic elements will be created.
2853 // Used then generated 3D mesh without geometry.
2854 //=======================================================================
2856 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
2858 int NbAllEdgsAndFaces=0;
2859 int NbQuadFacesAndEdgs=0;
2860 int NbFacesAndEdges=0;
2861 //All faces and edges
2862 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
2863 if ( NbAllEdgsAndFaces == 0 )
2864 return SMESH_MesherHelper::LINEAR;
2866 //Quadratic faces and edges
2867 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
2869 //Linear faces and edges
2870 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
2872 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
2874 return SMESH_MesherHelper::QUADRATIC;
2876 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
2878 return SMESH_MesherHelper::LINEAR;
2881 //Mesh with both type of elements
2882 return SMESH_MesherHelper::COMP;
2885 //=======================================================================
2886 //function : GetOtherParam
2887 //purpose : Return an alternative parameter for a node on seam
2888 //=======================================================================
2890 double SMESH_MesherHelper::GetOtherParam(const double param) const
2892 int i = myParIndex & U_periodic ? 0 : 1;
2893 return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
2898 //=======================================================================
2900 * \brief Iterator on ancestors of the given type
2902 //=======================================================================
2904 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2906 TopTools_ListIteratorOfListOfShape _ancIter;
2907 TopAbs_ShapeEnum _type;
2908 TopTools_MapOfShape _encountered;
2909 TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2910 : _ancIter( ancestors ), _type( type )
2912 if ( _ancIter.More() ) {
2913 if ( _ancIter.Value().ShapeType() != _type ) next();
2914 else _encountered.Add( _ancIter.Value() );
2919 return _ancIter.More();
2921 virtual const TopoDS_Shape* next()
2923 const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2924 if ( _ancIter.More() )
2925 for ( _ancIter.Next(); _ancIter.More(); _ancIter.Next())
2926 if ( _ancIter.Value().ShapeType() == _type && _encountered.Add( _ancIter.Value() ))
2934 //=======================================================================
2936 * \brief Return iterator on ancestors of the given type
2938 //=======================================================================
2940 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2941 const SMESH_Mesh& mesh,
2942 TopAbs_ShapeEnum ancestorType)
2944 return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));
2947 //=======================================================================
2948 //function : GetCommonAncestor
2949 //purpose : Find a common ancestors of two shapes of the given type
2950 //=======================================================================
2952 TopoDS_Shape SMESH_MesherHelper::GetCommonAncestor(const TopoDS_Shape& shape1,
2953 const TopoDS_Shape& shape2,
2954 const SMESH_Mesh& mesh,
2955 TopAbs_ShapeEnum ancestorType)
2957 TopoDS_Shape commonAnc;
2958 if ( !shape1.IsNull() && !shape2.IsNull() )
2960 PShapeIteratorPtr ancIt = GetAncestors( shape1, mesh, ancestorType );
2961 while ( const TopoDS_Shape* anc = ancIt->next() )
2962 if ( IsSubShape( shape2, *anc ))
2971 //#include <Perf_Meter.hxx>
2973 //=======================================================================
2974 namespace { // Structures used by FixQuadraticElements()
2975 //=======================================================================
2977 #define __DMP__(txt) \
2979 #define MSG(txt) __DMP__(txt<<endl)
2980 #define MSGBEG(txt) __DMP__(txt)
2982 //const double straightTol2 = 1e-33; // to detect straing links
2983 bool isStraightLink(double linkLen2, double middleNodeMove2)
2985 // straight if <node move> < 1/15 * <link length>
2986 return middleNodeMove2 < 1/15./15. * linkLen2;
2990 // ---------------------------------------
2992 * \brief Quadratic link knowing its faces
2994 struct QLink: public SMESH_TLink
2996 const SMDS_MeshNode* _mediumNode;
2997 mutable vector<const QFace* > _faces;
2998 mutable gp_Vec _nodeMove;
2999 mutable int _nbMoves;
3001 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
3002 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
3004 //if ( MediumPos() != SMDS_TOP_3DSPACE )
3005 _nodeMove = MediumPnt() - MiddlePnt();
3007 void SetContinuesFaces() const;
3008 const QFace* GetContinuesFace( const QFace* face ) const;
3009 bool OnBoundary() const;
3010 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
3011 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
3013 SMDS_TypeOfPosition MediumPos() const
3014 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
3015 SMDS_TypeOfPosition EndPos(bool isSecond) const
3016 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
3017 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
3018 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
3020 void Move(const gp_Vec& move, bool sum=false) const
3021 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }