1 // Copyright (C) 2007-2008 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 // File: SMESH_MesherHelper.cxx
23 // Created: 15.02.06 15:22:41
24 // Author: Sergey KUUL
26 #include "SMESH_MesherHelper.hxx"
28 #include "SMDS_FacePosition.hxx"
29 #include "SMDS_EdgePosition.hxx"
30 #include "SMDS_VolumeTool.hxx"
31 #include "SMESH_subMesh.hxx"
33 #include <BRepAdaptor_Surface.hxx>
34 #include <BRepTools.hxx>
35 #include <BRepTools_WireExplorer.hxx>
36 #include <BRep_Tool.hxx>
37 #include <Geom2d_Curve.hxx>
38 #include <GeomAPI_ProjectPointOnSurf.hxx>
39 #include <Geom_Curve.hxx>
40 #include <Geom_Surface.hxx>
41 #include <ShapeAnalysis.hxx>
43 #include <TopExp_Explorer.hxx>
44 #include <TopTools_ListIteratorOfListOfShape.hxx>
45 #include <TopTools_MapIteratorOfMapOfShape.hxx>
46 #include <TopTools_MapOfShape.hxx>
49 #include <gp_Pnt2d.hxx>
50 #include <gp_Trsf.hxx>
52 #include <Standard_Failure.hxx>
53 #include <Standard_ErrorHandler.hxx>
55 #include <utilities.h>
59 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
63 gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
67 //================================================================================
71 //================================================================================
73 SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
74 : myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false)
76 mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
79 //=======================================================================
80 //function : CheckShape
82 //=======================================================================
84 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
86 SMESHDS_Mesh* meshDS = GetMeshDS();
87 // we can create quadratic elements only if all elements
88 // created on subshapes of given shape are quadratic
89 // also we have to fill myTLinkNodeMap
90 myCreateQuadratic = true;
91 mySeamShapeIds.clear();
92 myDegenShapeIds.clear();
93 TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
94 SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
96 int nbOldLinks = myTLinkNodeMap.size();
98 TopExp_Explorer exp( aSh, subType );
99 for (; exp.More() && myCreateQuadratic; exp.Next()) {
100 if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
101 if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
103 const SMDS_MeshElement* e = it->next();
104 if ( e->GetType() != elemType || !e->IsQuadratic() ) {
105 myCreateQuadratic = false;
110 switch ( e->NbNodes() ) {
112 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
114 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
115 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
116 AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
118 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
119 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
120 AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
121 AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
124 myCreateQuadratic = false;
133 if ( nbOldLinks == myTLinkNodeMap.size() )
134 myCreateQuadratic = false;
136 if(!myCreateQuadratic) {
137 myTLinkNodeMap.clear();
141 return myCreateQuadratic;
144 //================================================================================
146 * \brief Set geomerty to make elements on
147 * \param aSh - geomertic shape
149 //================================================================================
151 void SMESH_MesherHelper::SetSubShape(const int aShID)
153 if ( aShID == myShapeID )
156 SetSubShape( GetMeshDS()->IndexToShape( aShID ));
158 SetSubShape( TopoDS_Shape() );
161 //================================================================================
163 * \brief Set geomerty to make elements on
164 * \param aSh - geomertic shape
166 //================================================================================
168 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
170 if ( myShape.IsSame( aSh ))
174 mySeamShapeIds.clear();
175 myDegenShapeIds.clear();
177 if ( myShape.IsNull() ) {
181 SMESHDS_Mesh* meshDS = GetMeshDS();
182 myShapeID = meshDS->ShapeToIndex(aSh);
184 // treatment of periodic faces
185 for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
187 const TopoDS_Face& face = TopoDS::Face( eF.Current() );
188 BRepAdaptor_Surface surface( face );
189 if ( surface.IsUPeriodic() || surface.IsVPeriodic() )
191 for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
193 // look for a seam edge
194 const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
195 if ( BRep_Tool::IsClosed( edge, face )) {
196 // initialize myPar1, myPar2 and myParIndex
197 if ( mySeamShapeIds.empty() ) {
199 BRep_Tool::UVPoints( edge, face, uv1, uv2 );
200 if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
202 myParIndex = 1; // U periodic
203 myPar1 = surface.FirstUParameter();
204 myPar2 = surface.LastUParameter();
207 myParIndex = 2; // V periodic
208 myPar1 = surface.FirstVParameter();
209 myPar2 = surface.LastVParameter();
212 // store seam shape indices, negative if shape encounters twice
213 int edgeID = meshDS->ShapeToIndex( edge );
214 mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
215 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
216 int vertexID = meshDS->ShapeToIndex( v.Current() );
217 mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
221 // look for a degenerated edge
222 if ( BRep_Tool::Degenerated( edge )) {
223 myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
224 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
225 myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
232 //================================================================================
234 * \brief Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
235 * \param F - the face
236 * \retval bool - return true if the face is periodic
238 //================================================================================
240 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
242 if ( F.IsNull() ) return !mySeamShapeIds.empty();
244 if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
245 return !mySeamShapeIds.empty();
248 Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
249 if ( !aSurface.IsNull() )
250 return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
255 //=======================================================================
256 //function : IsMedium
258 //=======================================================================
260 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode* node,
261 const SMDSAbs_ElementType typeToCheck)
263 return SMESH_MeshEditor::IsMedium( node, typeToCheck );
266 //=======================================================================
268 * \brief Return support shape of a node
269 * \param node - the node
270 * \param meshDS - mesh DS
271 * \retval TopoDS_Shape - found support shape
273 //=======================================================================
275 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
276 SMESHDS_Mesh* meshDS)
278 int shapeID = node->GetPosition()->GetShapeId();
279 if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
280 return meshDS->IndexToShape( shapeID );
282 return TopoDS_Shape();
286 //=======================================================================
287 //function : AddTLinkNode
289 //=======================================================================
291 * Auxilary function for filling myTLinkNodeMap
293 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
294 const SMDS_MeshNode* n2,
295 const SMDS_MeshNode* n12)
297 // add new record to map
298 SMESH_TLink link( n1, n2 );
299 myTLinkNodeMap.insert( make_pair(link,n12));
302 //=======================================================================
304 * \brief Select UV on either of 2 pcurves of a seam edge, closest to the given UV
305 * \param uv1 - UV on the seam
306 * \param uv2 - UV within a face
307 * \retval gp_Pnt2d - selected UV
309 //=======================================================================
311 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
313 double p1 = uv1.Coord( myParIndex );
314 double p2 = uv2.Coord( myParIndex );
315 double p3 = ( Abs( p1 - myPar1 ) < Abs( p1 - myPar2 )) ? myPar2 : myPar1;
316 if ( Abs( p2 - p1 ) > Abs( p2 - p3 ))
318 gp_Pnt2d result = uv1;
319 result.SetCoord( myParIndex, p1 );
323 //=======================================================================
325 * \brief Return node UV on face
326 * \param F - the face
327 * \param n - the node
328 * \param n2 - a node of element being created located inside a face
329 * \param check - optional flag returing false if found UV are invalid
330 * \retval gp_XY - resulting UV
332 //=======================================================================
334 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
335 const SMDS_MeshNode* n,
336 const SMDS_MeshNode* n2,
339 gp_Pnt2d uv( 1e100, 1e100 );
340 const SMDS_PositionPtr Pos = n->GetPosition();
342 if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
344 // node has position on face
345 const SMDS_FacePosition* fpos =
346 static_cast<const SMDS_FacePosition*>(n->GetPosition().get());
347 uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
348 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( F ));
350 else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
352 // node has position on edge => it is needed to find
353 // corresponding edge from face, get pcurve for this
354 // edge and retrieve value from this pcurve
355 const SMDS_EdgePosition* epos =
356 static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
357 int edgeID = Pos->GetShapeId();
358 TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
360 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
361 uv = C2d->Value( epos->GetUParameter() );
362 uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( E ));
364 // for a node on a seam edge select one of UVs on 2 pcurves
365 if ( n2 && IsSeamShape( edgeID ) )
367 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
370 { // adjust uv to period
372 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
373 Standard_Boolean isUPeriodic = S->IsUPeriodic();
374 Standard_Boolean isVPeriodic = S->IsVPeriodic();
375 if ( isUPeriodic || isVPeriodic ) {
376 Standard_Real UF,UL,VF,VL;
377 S->Bounds(UF,UL,VF,VL);
379 uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
381 uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
385 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
387 if ( int vertexID = n->GetPosition()->GetShapeId() ) {
388 const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
390 uv = BRep_Tool::Parameters( V, F );
393 catch (Standard_Failure& exc) {
396 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() )
397 uvOK = ( V == vert.Current() );
400 MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
401 << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
403 // get UV of a vertex closest to the node
405 gp_Pnt pn = XYZ( n );
406 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) {
407 TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
408 gp_Pnt p = BRep_Tool::Pnt( curV );
409 double curDist = p.SquareDistance( pn );
410 if ( curDist < dist ) {
412 uv = BRep_Tool::Parameters( curV, F );
413 uvOK = ( dist < DBL_MIN );
418 TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
419 for ( ; it.More(); it.Next() ) {
420 if ( it.Value().ShapeType() == TopAbs_EDGE ) {
421 const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
423 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
424 if ( !C2d.IsNull() ) {
425 double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
426 uv = C2d->Value( u );
433 if ( n2 && IsSeamShape( vertexID ) )
434 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
444 //=======================================================================
446 * \brief Check and fix node UV on a face
447 * \retval bool - false if UV is bad and could not be fixed
449 //=======================================================================
451 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face& F,
452 const SMDS_MeshNode* n,
454 const double tol) const
456 if ( !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
458 // check that uv is correct
460 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
461 gp_Pnt nodePnt = XYZ( n );
462 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
463 if ( nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > tol )
465 // uv incorrect, project the node to surface
466 GeomAPI_ProjectPointOnSurf projector( nodePnt, surface, tol );
467 if ( !projector.IsDone() || projector.NbPoints() < 1 )
469 MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
472 Quantity_Parameter U,V;
473 projector.LowerDistanceParameters(U,V);
474 if ( nodePnt.Distance( surface->Value( U, V )) > tol )
476 MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
481 else if ( uv.Modulus() > numeric_limits<double>::min() )
483 ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
489 //=======================================================================
491 * \brief Return middle UV taking in account surface period
493 //=======================================================================
495 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
499 if ( surface.IsNull() )
500 return 0.5 * ( p1 + p2 );
501 //checking if surface is periodic
502 Standard_Real UF,UL,VF,VL;
503 surface->Bounds(UF,UL,VF,VL);
506 Standard_Boolean isUPeriodic = surface->IsUPeriodic();
508 Standard_Real UPeriod = surface->UPeriod();
509 Standard_Real p2x = p2.X()+ShapeAnalysis::AdjustByPeriod(p2.X(),p1.X(),UPeriod);
510 Standard_Real pmid = (p1.X()+p2x)/2.;
511 u = pmid+ShapeAnalysis::AdjustToPeriod(pmid,UF,UL);
514 u= (p1.X()+p2.X())/2.;
516 Standard_Boolean isVPeriodic = surface->IsVPeriodic();
518 Standard_Real VPeriod = surface->VPeriod();
519 Standard_Real p2y = p2.Y()+ShapeAnalysis::AdjustByPeriod(p2.Y(),p1.Y(),VPeriod);
520 Standard_Real pmid = (p1.Y()+p2y)/2.;
521 v = pmid+ShapeAnalysis::AdjustToPeriod(pmid,VF,VL);
524 v = (p1.Y()+p2.Y())/2.;
529 //=======================================================================
531 * \brief Return node U on edge
532 * \param E - the Edge
533 * \param n - the node
534 * \retval double - resulting U
536 * Auxilary function called form GetMediumNode()
538 //=======================================================================
540 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
541 const SMDS_MeshNode* n,
545 const SMDS_PositionPtr Pos = n->GetPosition();
546 if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE) {
547 const SMDS_EdgePosition* epos =
548 static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
549 param = epos->GetUParameter();
551 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX) {
552 SMESHDS_Mesh * meshDS = GetMeshDS();
553 int vertexID = n->GetPosition()->GetShapeId();
554 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
555 param = BRep_Tool::Parameter( V, E );
560 //================================================================================
562 * \brief Return existing or create new medium nodes between given ones
563 * \param force3d - if true, new node is the middle of n1 and n2,
564 * else is located on geom face or geom edge
566 //================================================================================
568 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
569 const SMDS_MeshNode* n2,
572 SMESH_TLink link(n1,n2);
573 ItTLinkNode itLN = myTLinkNodeMap.find( link );
574 if ( itLN != myTLinkNodeMap.end() ) {
575 return (*itLN).second;
577 // create medium node
579 SMESHDS_Mesh* meshDS = GetMeshDS();
580 int faceID = -1, edgeID = -1;
581 const SMDS_PositionPtr Pos1 = n1->GetPosition();
582 const SMDS_PositionPtr Pos2 = n2->GetPosition();
584 if( myShape.IsNull() )
586 if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
587 faceID = Pos1->GetShapeId();
589 else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
590 faceID = Pos2->GetShapeId();
593 if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
594 edgeID = Pos1->GetShapeId();
596 if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
597 edgeID = Pos2->GetShapeId();
602 // we try to create medium node using UV parameters of
603 // nodes, else - medium between corresponding 3d points
605 TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
606 if(faceID>0 || shapeType == TopAbs_FACE) {
607 // obtaining a face and 2d points for nodes
609 if( myShape.IsNull() )
610 F = TopoDS::Face(meshDS->IndexToShape(faceID));
612 F = TopoDS::Face(myShape);
616 gp_XY p1 = GetNodeUV(F,n1,n2, &uvOK1);
617 gp_XY p2 = GetNodeUV(F,n2,n1, &uvOK2);
619 if ( uvOK1 && uvOK2 )
621 if ( IsDegenShape( Pos1->GetShapeId() ))
622 p1.SetCoord( myParIndex, p2.Coord( myParIndex ));
623 else if ( IsDegenShape( Pos2->GetShapeId() ))
624 p2.SetCoord( myParIndex, p1.Coord( myParIndex ));
627 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
628 gp_XY uv = GetMiddleUV( S, p1, p2 );
629 gp_Pnt P = S->Value( uv.X(), uv.Y() ).Transformed(loc);
630 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
631 meshDS->SetNodeOnFace(n12, faceID, uv.X(), uv.Y());
632 myTLinkNodeMap.insert(make_pair(link,n12));
636 if (edgeID>0 || shapeType == TopAbs_EDGE) {
639 if( myShape.IsNull() )
640 E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
642 E = TopoDS::Edge(myShape);
646 double p1 = GetNodeU(E,n1);
647 double p2 = GetNodeU(E,n2);
650 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
653 Standard_Boolean isPeriodic = C->IsPeriodic();
656 Standard_Real Period = C->Period();
657 Standard_Real p = p2+ShapeAnalysis::AdjustByPeriod(p2,p1,Period);
658 Standard_Real pmid = (p1+p)/2.;
659 u = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
664 gp_Pnt P = C->Value( u );
665 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
666 meshDS->SetNodeOnEdge(n12, edgeID, u);
667 myTLinkNodeMap.insert(make_pair(link,n12));
673 double x = ( n1->X() + n2->X() )/2.;
674 double y = ( n1->Y() + n2->Y() )/2.;
675 double z = ( n1->Z() + n2->Z() )/2.;
676 n12 = meshDS->AddNode(x,y,z);
678 meshDS->SetNodeOnEdge(n12, edgeID);
680 meshDS->SetNodeOnFace(n12, faceID);
682 meshDS->SetNodeInVolume(n12, myShapeID);
683 myTLinkNodeMap.insert( make_pair( link, n12 ));
687 //=======================================================================
691 //=======================================================================
693 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
695 SMESHDS_Mesh * meshDS = GetMeshDS();
696 SMDS_MeshNode* node = 0;
698 node = meshDS->AddNodeWithID( x, y, z, ID );
700 node = meshDS->AddNode( x, y, z );
701 if ( mySetElemOnShape && myShapeID > 0 ) {
702 switch ( myShape.ShapeType() ) {
703 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
704 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
705 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID); break;
706 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID); break;
707 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
714 //=======================================================================
716 * Creates quadratic or linear edge
718 //=======================================================================
720 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
721 const SMDS_MeshNode* n2,
725 SMESHDS_Mesh * meshDS = GetMeshDS();
727 SMDS_MeshEdge* edge = 0;
728 if (myCreateQuadratic) {
729 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
731 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
733 edge = meshDS->AddEdge(n1, n2, n12);
737 edge = meshDS->AddEdgeWithID(n1, n2, id);
739 edge = meshDS->AddEdge(n1, n2);
742 if ( mySetElemOnShape && myShapeID > 0 )
743 meshDS->SetMeshElementOnShape( edge, myShapeID );
748 //=======================================================================
750 * Creates quadratic or linear triangle
752 //=======================================================================
754 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
755 const SMDS_MeshNode* n2,
756 const SMDS_MeshNode* n3,
760 SMESHDS_Mesh * meshDS = GetMeshDS();
761 SMDS_MeshFace* elem = 0;
763 if( n1==n2 || n2==n3 || n3==n1 )
766 if(!myCreateQuadratic) {
768 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
770 elem = meshDS->AddFace(n1, n2, n3);
773 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
774 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
775 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
778 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
780 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
782 if ( mySetElemOnShape && myShapeID > 0 )
783 meshDS->SetMeshElementOnShape( elem, myShapeID );
788 //=======================================================================
790 * Creates quadratic or linear quadrangle
792 //=======================================================================
794 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
795 const SMDS_MeshNode* n2,
796 const SMDS_MeshNode* n3,
797 const SMDS_MeshNode* n4,
801 SMESHDS_Mesh * meshDS = GetMeshDS();
802 SMDS_MeshFace* elem = 0;
805 return AddFace(n1,n3,n4,id,force3d);
808 return AddFace(n1,n2,n4,id,force3d);
811 return AddFace(n1,n2,n3,id,force3d);
814 return AddFace(n1,n2,n4,id,force3d);
817 return AddFace(n1,n2,n3,id,force3d);
820 return AddFace(n1,n2,n3,id,force3d);
823 if(!myCreateQuadratic) {
825 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
827 elem = meshDS->AddFace(n1, n2, n3, n4);
830 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
831 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
832 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
833 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
836 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
838 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
840 if ( mySetElemOnShape && myShapeID > 0 )
841 meshDS->SetMeshElementOnShape( elem, myShapeID );
846 //=======================================================================
848 * Creates quadratic or linear volume
850 //=======================================================================
852 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
853 const SMDS_MeshNode* n2,
854 const SMDS_MeshNode* n3,
855 const SMDS_MeshNode* n4,
856 const SMDS_MeshNode* n5,
857 const SMDS_MeshNode* n6,
861 SMESHDS_Mesh * meshDS = GetMeshDS();
862 SMDS_MeshVolume* elem = 0;
863 if(!myCreateQuadratic) {
865 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
867 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
870 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
871 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
872 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
874 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
875 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
876 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
878 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
879 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
880 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
883 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
884 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
886 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
887 n12, n23, n31, n45, n56, n64, n14, n25, n36);
889 if ( mySetElemOnShape && myShapeID > 0 )
890 meshDS->SetMeshElementOnShape( elem, myShapeID );
895 //=======================================================================
897 * Creates quadratic or linear volume
899 //=======================================================================
901 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
902 const SMDS_MeshNode* n2,
903 const SMDS_MeshNode* n3,
904 const SMDS_MeshNode* n4,
908 SMESHDS_Mesh * meshDS = GetMeshDS();
909 SMDS_MeshVolume* elem = 0;
910 if(!myCreateQuadratic) {
912 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
914 elem = meshDS->AddVolume(n1, n2, n3, n4);
917 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
918 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
919 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
921 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
922 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
923 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
926 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
928 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
930 if ( mySetElemOnShape && myShapeID > 0 )
931 meshDS->SetMeshElementOnShape( elem, myShapeID );
936 //=======================================================================
938 * Creates quadratic or linear pyramid
940 //=======================================================================
942 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
943 const SMDS_MeshNode* n2,
944 const SMDS_MeshNode* n3,
945 const SMDS_MeshNode* n4,
946 const SMDS_MeshNode* n5,
950 SMDS_MeshVolume* elem = 0;
951 if(!myCreateQuadratic) {
953 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
955 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
958 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
959 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
960 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
961 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
963 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
964 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
965 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
966 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
969 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
974 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
978 if ( mySetElemOnShape && myShapeID > 0 )
979 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
984 //=======================================================================
986 * Creates quadratic or linear hexahedron
988 //=======================================================================
990 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
991 const SMDS_MeshNode* n2,
992 const SMDS_MeshNode* n3,
993 const SMDS_MeshNode* n4,
994 const SMDS_MeshNode* n5,
995 const SMDS_MeshNode* n6,
996 const SMDS_MeshNode* n7,
997 const SMDS_MeshNode* n8,
1001 SMESHDS_Mesh * meshDS = GetMeshDS();
1002 SMDS_MeshVolume* elem = 0;
1003 if(!myCreateQuadratic) {
1005 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1007 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1010 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1011 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1012 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1013 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1015 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1016 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1017 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1018 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1020 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1021 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1022 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1023 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1026 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1027 n12, n23, n34, n41, n56, n67,
1028 n78, n85, n15, n26, n37, n48, id);
1030 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1031 n12, n23, n34, n41, n56, n67,
1032 n78, n85, n15, n26, n37, n48);
1034 if ( mySetElemOnShape && myShapeID > 0 )
1035 meshDS->SetMeshElementOnShape( elem, myShapeID );
1040 //=======================================================================
1042 * \brief Load nodes bound to face into a map of node columns
1043 * \param theParam2ColumnMap - map of node columns to fill
1044 * \param theFace - the face on which nodes are searched for
1045 * \param theBaseEdge - the edge nodes of which are columns' bases
1046 * \param theMesh - the mesh containing nodes
1047 * \retval bool - false if something is wrong
1049 * The key of the map is a normalized parameter of each
1050 * base node on theBaseEdge.
1051 * This method works in supposition that nodes on the face
1052 * forms a rectangular grid and elements can be quardrangles or triangles
1054 //=======================================================================
1056 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1057 const TopoDS_Face& theFace,
1058 const TopoDS_Edge& theBaseEdge,
1059 SMESHDS_Mesh* theMesh)
1061 // get vertices of theBaseEdge
1062 TopoDS_Vertex vfb, vlb, vft; // first and last, bottom and top vertices
1063 TopoDS_Edge eFrw = TopoDS::Edge( theBaseEdge.Oriented( TopAbs_FORWARD ));
1064 TopExp::Vertices( eFrw, vfb, vlb );
1066 // find the other edges of theFace and orientation of e1
1067 TopoDS_Edge e1, e2, eTop;
1068 bool rev1, CumOri = false;
1069 TopExp_Explorer exp( theFace, TopAbs_EDGE );
1071 for ( ; exp.More(); exp.Next() ) {
1072 if ( ++nbEdges > 4 ) {
1073 return false; // more than 4 edges in theFace
1075 TopoDS_Edge e = TopoDS::Edge( exp.Current() );
1076 if ( theBaseEdge.IsSame( e ))
1078 TopoDS_Vertex vCommon;
1079 if ( !TopExp::CommonVertex( theBaseEdge, e, vCommon ))
1081 else if ( vCommon.IsSame( vfb )) {
1083 vft = TopExp::LastVertex( e1, CumOri );
1084 rev1 = vfb.IsSame( vft );
1086 vft = TopExp::FirstVertex( e1, CumOri );
1091 if ( nbEdges < 4 ) {
1092 return false; // less than 4 edges in theFace
1094 if ( e2.IsNull() && vfb.IsSame( vlb ))
1097 // submeshes corresponding to shapes
1098 SMESHDS_SubMesh* smFace = theMesh->MeshElements( theFace );
1099 SMESHDS_SubMesh* smb = theMesh->MeshElements( theBaseEdge );
1100 SMESHDS_SubMesh* smt = theMesh->MeshElements( eTop );
1101 SMESHDS_SubMesh* sm1 = theMesh->MeshElements( e1 );
1102 SMESHDS_SubMesh* sm2 = theMesh->MeshElements( e2 );
1103 SMESHDS_SubMesh* smVfb = theMesh->MeshElements( vfb );
1104 SMESHDS_SubMesh* smVlb = theMesh->MeshElements( vlb );
1105 SMESHDS_SubMesh* smVft = theMesh->MeshElements( vft );
1106 if (!smFace || !smb || !smt || !sm1 || !sm2 || !smVfb || !smVlb || !smVft ) {
1107 RETURN_BAD_RESULT( "NULL submesh " <<smFace<<" "<<smb<<" "<<smt<<" "<<
1108 sm1<<" "<<sm2<<" "<<smVfb<<" "<<smVlb<<" "<<smVft);
1110 if ( smb->NbNodes() != smt->NbNodes() || sm1->NbNodes() != sm2->NbNodes() ) {
1111 RETURN_BAD_RESULT(" Diff nb of nodes on opposite edges" );
1113 if (smVfb->NbNodes() != 1 || smVlb->NbNodes() != 1 || smVft->NbNodes() != 1) {
1114 RETURN_BAD_RESULT("Empty submesh of vertex");
1116 // define whether mesh is quadratic
1117 bool isQuadraticMesh = false;
1118 SMDS_ElemIteratorPtr eIt = smFace->GetElements();
1119 if ( !eIt->more() ) {
1120 RETURN_BAD_RESULT("No elements on the face");
1122 const SMDS_MeshElement* e = eIt->next();
1123 isQuadraticMesh = e->IsQuadratic();
1125 if ( sm1->NbNodes() * smb->NbNodes() != smFace->NbNodes() ) {
1126 // check quadratic case
1127 if ( isQuadraticMesh ) {
1128 // what if there are quadrangles and triangles mixed?
1129 // int n1 = sm1->NbNodes()/2;
1130 // int n2 = smb->NbNodes()/2;
1131 // int n3 = sm1->NbNodes() - n1;
1132 // int n4 = smb->NbNodes() - n2;
1133 // int nf = sm1->NbNodes()*smb->NbNodes() - n3*n4;
1134 // if( nf != smFace->NbNodes() ) {
1135 // MESSAGE( "Wrong nb face nodes: " <<
1136 // sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
1141 RETURN_BAD_RESULT( "Wrong nb face nodes: " <<
1142 sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
1146 int vsize = sm1->NbNodes() + 2;
1147 int hsize = smb->NbNodes() + 2;
1148 if(isQuadraticMesh) {
1149 vsize = vsize - sm1->NbNodes()/2 -1;
1150 hsize = hsize - smb->NbNodes()/2 -1;
1153 // load nodes from theBaseEdge
1155 std::set<const SMDS_MeshNode*> loadedNodes;
1156 const SMDS_MeshNode* nullNode = 0;
1158 std::vector<const SMDS_MeshNode*> & nVecf = theParam2ColumnMap[ 0.];
1159 nVecf.resize( vsize, nullNode );
1160 loadedNodes.insert( nVecf[ 0 ] = smVfb->GetNodes()->next() );
1162 std::vector<const SMDS_MeshNode*> & nVecl = theParam2ColumnMap[ 1.];
1163 nVecl.resize( vsize, nullNode );
1164 loadedNodes.insert( nVecl[ 0 ] = smVlb->GetNodes()->next() );
1167 BRep_Tool::Range( eFrw, f, l );
1168 double range = l - f;
1169 SMDS_NodeIteratorPtr nIt = smb->GetNodes();
1170 const SMDS_MeshNode* node;
1171 while ( nIt->more() ) {
1173 if(IsMedium(node, SMDSAbs_Edge))
1175 const SMDS_EdgePosition* pos =
1176 dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
1180 double u = ( pos->GetUParameter() - f ) / range;
1181 std::vector<const SMDS_MeshNode*> & nVec = theParam2ColumnMap[ u ];
1182 nVec.resize( vsize, nullNode );
1183 loadedNodes.insert( nVec[ 0 ] = node );
1185 if ( theParam2ColumnMap.size() != hsize ) {
1186 RETURN_BAD_RESULT( "Wrong node positions on theBaseEdge" );
1189 // load nodes from e1
1191 std::map< double, const SMDS_MeshNode*> sortedNodes; // sort by param on edge
1192 nIt = sm1->GetNodes();
1193 while ( nIt->more() ) {
1197 const SMDS_EdgePosition* pos =
1198 dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
1202 sortedNodes.insert( std::make_pair( pos->GetUParameter(), node ));
1204 loadedNodes.insert( nVecf[ vsize - 1 ] = smVft->GetNodes()->next() );
1205 std::map< double, const SMDS_MeshNode*>::iterator u_n = sortedNodes.begin();
1206 int row = rev1 ? vsize - 1 : 0;
1207 int dRow = rev1 ? -1 : +1;
1208 for ( ; u_n != sortedNodes.end(); u_n++ ) {
1210 loadedNodes.insert( nVecf[ row ] = u_n->second );
1213 // try to load the rest nodes
1215 // get all faces from theFace
1216 TIDSortedElemSet allFaces, foundFaces;
1217 eIt = smFace->GetElements();
1218 while ( eIt->more() ) {
1219 const SMDS_MeshElement* e = eIt->next();
1220 if ( e->GetType() == SMDSAbs_Face )
1221 allFaces.insert( e );
1223 // Starting from 2 neighbour nodes on theBaseEdge, look for a face
1224 // the nodes belong to, and between the nodes of the found face,
1225 // look for a not loaded node considering this node to be the next
1226 // in a column of the starting second node. Repeat, starting
1227 // from nodes next to the previous starting nodes in their columns,
1228 // and so on while a face can be found. Then go the the next pair
1229 // of nodes on theBaseEdge.
1230 TParam2ColumnMap::iterator par_nVec_1 = theParam2ColumnMap.begin();
1231 TParam2ColumnMap::iterator par_nVec_2 = par_nVec_1;
1234 for ( par_nVec_2++; par_nVec_2 != theParam2ColumnMap.end(); par_nVec_1++, par_nVec_2++ ) {
1237 const SMDS_MeshNode* n1 = par_nVec_1->second[ row ];
1238 const SMDS_MeshNode* n2 = par_nVec_2->second[ row ];
1239 const SMDS_MeshElement* face = 0;
1240 bool lastColOnClosedFace = ( nVecf[ row ] == n2 );
1242 // look for a face by 2 nodes
1243 face = SMESH_MeshEditor::FindFaceInSet( n1, n2, allFaces, foundFaces );
1245 int nbFaceNodes = face->NbNodes();
1246 if ( face->IsQuadratic() )
1248 if ( nbFaceNodes>4 ) {
1249 RETURN_BAD_RESULT(" Too many nodes in a face: " << nbFaceNodes );
1251 // look for a not loaded node of the <face>
1253 const SMDS_MeshNode* n3 = 0; // a node defferent from n1 and n2
1254 for ( int i = 0; i < nbFaceNodes && !found; ++i ) {
1255 node = face->GetNode( i );
1256 found = loadedNodes.insert( node ).second;
1257 if ( !found && node != n1 && node != n2 )
1260 if ( lastColOnClosedFace && row + 1 < vsize ) {
1261 node = nVecf[ row + 1 ];
1262 found = ( face->GetNodeIndex( node ) >= 0 );
1265 if ( ++row > vsize - 1 ) {
1266 RETURN_BAD_RESULT( "Too many nodes in column "<< col <<": "<< row+1);
1268 par_nVec_2->second[ row ] = node;
1269 foundFaces.insert( face );
1271 if ( nbFaceNodes==4 ) {
1272 n1 = par_nVec_1->second[ row ];
1275 else if ( nbFaceNodes==3 && n3 == par_nVec_1->second[ row + 1 ] ) {
1279 RETURN_BAD_RESULT( "Not quad mesh, column "<< col );
1283 while ( face && n1 && n2 );
1285 if ( row < vsize - 1 ) {
1286 MESSAGE( "Too few nodes in column "<< col <<": "<< row+1);
1287 MESSAGE( "Base node 1: "<< par_nVec_1->second[0]);
1288 MESSAGE( "Base node 2: "<< par_nVec_2->second[0]);
1289 if ( n1 ) { MESSAGE( "Current node 1: "<< n1); }
1290 else { MESSAGE( "Current node 1: NULL"); }
1291 if ( n2 ) { MESSAGE( "Current node 2: "<< n2); }
1292 else { MESSAGE( "Current node 2: NULL"); }
1293 MESSAGE( "first base node: "<< theParam2ColumnMap.begin()->second[0]);
1294 MESSAGE( "last base node: "<< theParam2ColumnMap.rbegin()->second[0]);
1297 } // loop on columns
1302 //=======================================================================
1304 * \brief Return number of unique ancestors of the shape
1306 //=======================================================================
1308 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1309 const SMESH_Mesh& mesh,
1310 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
1312 TopTools_MapOfShape ancestors;
1313 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1314 for ( ; ansIt.More(); ansIt.Next() ) {
1315 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1316 ancestors.Add( ansIt.Value() );
1318 return ancestors.Extent();
1321 //=======================================================================
1323 * Check mesh without geometry for: if all elements on this shape are quadratic,
1324 * quadratic elements will be created.
1325 * Used then generated 3D mesh without geometry.
1327 //=======================================================================
1329 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1331 int NbAllEdgsAndFaces=0;
1332 int NbQuadFacesAndEdgs=0;
1333 int NbFacesAndEdges=0;
1334 //All faces and edges
1335 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1337 //Quadratic faces and edges
1338 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1340 //Linear faces and edges
1341 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1343 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1345 return SMESH_MesherHelper::QUADRATIC;
1347 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1349 return SMESH_MesherHelper::LINEAR;
1352 //Mesh with both type of elements
1353 return SMESH_MesherHelper::COMP;
1356 //=======================================================================
1358 * \brief Return an alternative parameter for a node on seam
1360 //=======================================================================
1362 double SMESH_MesherHelper::GetOtherParam(const double param) const
1364 return fabs(param-myPar1) < fabs(param-myPar2) ? myPar2 : myPar1;
1367 //=======================================================================
1368 namespace { // Structures used by FixQuadraticElements()
1369 //=======================================================================
1371 #define __DMP__(txt) \
1373 #define MSG(txt) __DMP__(txt<<endl)
1374 #define MSGBEG(txt) __DMP__(txt)
1376 const double straightTol2 = 1e-33; // to detect straing links
1379 // ---------------------------------------
1381 * \brief Quadratic link knowing its faces
1383 struct QLink: public SMESH_TLink
1385 const SMDS_MeshNode* _mediumNode;
1386 mutable vector<const QFace* > _faces;
1387 mutable gp_Vec _nodeMove;
1388 mutable int _nbMoves;
1390 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1391 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1393 //if ( MediumPos() != SMDS_TOP_3DSPACE )
1394 _nodeMove = MediumPnt() - MiddlePnt();
1396 void SetContinuesFaces() const;
1397 const QFace* GetContinuesFace( const QFace* face ) const;
1398 bool OnBoundary() const;
1399 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1400 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1402 SMDS_TypeOfPosition MediumPos() const
1403 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1404 SMDS_TypeOfPosition EndPos(bool isSecond) const
1405 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1406 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1407 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1409 void Move(const gp_Vec& move, bool sum=false) const
1410 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1411 gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1412 bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1413 bool IsStraight() const { return _nodeMove.SquareMagnitude() <= straightTol2; }
1415 bool operator<(const QLink& other) const {
1416 return (node1()->GetID() == other.node1()->GetID() ?
1417 node2()->GetID() < other.node2()->GetID() :
1418 node1()->GetID() < other.node1()->GetID());
1420 struct PtrComparator {
1421 bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1424 // ---------------------------------------------------------
1426 * \brief Link in the chain of links; it connects two faces
1430 const QLink* _qlink;
1431 mutable const QFace* _qfaces[2];
1433 TChainLink(const QLink* qlink=0):_qlink(qlink) {
1434 _qfaces[0] = _qfaces[1] = 0;
1436 void SetFace(const QFace* face) { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1438 bool IsBoundary() const { return !_qfaces[1]; }
1440 void RemoveFace( const QFace* face ) const
1441 { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1443 const QFace* NextFace( const QFace* f ) const
1444 { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1446 const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1447 { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1449 bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1451 operator bool() const { return (_qlink); }
1453 const QLink* operator->() const { return _qlink; }
1455 gp_Vec Normal() const;
1457 // --------------------------------------------------------------------
1458 typedef list< TChainLink > TChain;
1459 typedef set < TChainLink > TLinkSet;
1460 typedef TLinkSet::const_iterator TLinkInSet;
1462 const int theFirstStep = 5;
1464 enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1465 // --------------------------------------------------------------------
1467 * \brief Face shared by two volumes and bound by QLinks
1469 struct QFace: public TIDSortedElemSet
1471 mutable const SMDS_MeshElement* _volumes[2];
1472 mutable vector< const QLink* > _sides;
1473 mutable bool _sideIsAdded[4]; // added in chain of links
1476 QFace( const vector< const QLink*>& links );
1478 void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1480 int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1482 void AddSelfToLinks() const {
1483 for ( int i = 0; i < _sides.size(); ++i )
1484 _sides[i]->_faces.push_back( this );
1486 int LinkIndex( const QLink* side ) const {
1487 for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1490 bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const;
1492 bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1494 int i = LinkIndex( link._qlink );
1495 if ( i < 0 ) return true;
1496 _sideIsAdded[i] = true;
1497 link.SetFace( this );
1498 // continue from opposite link
1499 return GetLinkChain( (i+2)%_sides.size(), chain, pos, error );
1501 bool IsBoundary() const { return !_volumes[1]; }
1503 bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1505 TLinkInSet GetBoundaryLink( const TLinkSet& links,
1506 const TChainLink& avoidLink,
1507 TLinkInSet * notBoundaryLink = 0,
1508 const SMDS_MeshNode* nodeToContain = 0,
1509 bool * isAdjacentUsed = 0) const;
1511 TLinkInSet GetLinkByNode( const TLinkSet& links,
1512 const TChainLink& avoidLink,
1513 const SMDS_MeshNode* nodeToContain) const;
1515 const SMDS_MeshNode* GetNodeInFace() const {
1516 for ( int iL = 0; iL < _sides.size(); ++iL )
1517 if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1521 gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1523 double MoveByBoundary( const TChainLink& theLink,
1524 const gp_Vec& theRefVec,
1525 const TLinkSet& theLinks,
1526 SMESH_MesherHelper* theFaceHelper=0,
1527 const double thePrevLen=0,
1528 const int theStep=theFirstStep,
1529 gp_Vec* theLinkNorm=0,
1530 double theSign=1.0) const;
1533 //================================================================================
1535 * \brief Dump QLink and QFace
1537 ostream& operator << (ostream& out, const QLink& l)
1539 out <<"QLink nodes: "
1540 << l.node1()->GetID() << " - "
1541 << l._mediumNode->GetID() << " - "
1542 << l.node2()->GetID() << endl;
1545 ostream& operator << (ostream& out, const QFace& f)
1547 out <<"QFace nodes: "/*<< &f << " "*/;
1548 for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
1549 out << (*n)->GetID() << " ";
1550 out << " \tvolumes: "
1551 << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1552 << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1553 out << " \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1557 //================================================================================
1559 * \brief Construct QFace from QLinks
1561 //================================================================================
1563 QFace::QFace( const vector< const QLink*>& links )
1565 _volumes[0] = _volumes[1] = 0;
1567 _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1568 _normal.SetCoord(0,0,0);
1569 for ( int i = 1; i < _sides.size(); ++i ) {
1570 const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1571 insert( l1->node1() ); insert( l1->node2() );
1573 gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1574 gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1575 if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1579 double normSqSize = _normal.SquareMagnitude();
1580 if ( normSqSize > numeric_limits<double>::min() )
1581 _normal /= sqrt( normSqSize );
1583 _normal.SetCoord(1e-33,0,0);
1585 //================================================================================
1587 * \brief Make up chain of links
1588 * \param iSide - link to add first
1589 * \param chain - chain to fill in
1590 * \param pos - postion of medium nodes the links should have
1591 * \param error - out, specifies what is wrong
1592 * \retval bool - false if valid chain can't be built; "valid" means that links
1593 * of the chain belongs to rectangles bounding hexahedrons
1595 //================================================================================
1597 bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1599 if ( iSide >= _sides.size() ) // wrong argument iSide
1601 if ( _sideIsAdded[ iSide ]) // already in chain
1604 if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1606 for ( int i = 0; i < _sides.size(); ++i ) {
1607 if ( !_sideIsAdded[i] && _sides[i] ) {
1608 _sideIsAdded[i]=true;
1609 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(_sides[i]));
1610 chLink->SetFace( this );
1611 if ( _sides[i]->MediumPos() >= pos )
1612 if ( const QFace* f = _sides[i]->GetContinuesFace( this ))
1613 f->GetLinkChain( *chLink, chain, pos, error );
1616 if ( error < ERR_TRI )
1620 _sideIsAdded[iSide] = true; // not to add this link to chain again
1621 const QLink* link = _sides[iSide];
1625 // add link into chain
1626 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1627 chLink->SetFace( this );
1630 // propagate from rectangle to neighbour faces
1631 if ( link->MediumPos() >= pos ) {
1632 int nbLinkFaces = link->_faces.size();
1633 if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1634 // hexahedral mesh or boundary quadrangles - goto a continous face
1635 if ( const QFace* f = link->GetContinuesFace( this ))
1636 return f->GetLinkChain( *chLink, chain, pos, error );
1639 TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1640 for ( int i = 0; i < nbLinkFaces; ++i )
1641 if ( link->_faces[i] )
1642 link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1643 if ( error < ERR_PRISM )
1651 //================================================================================
1653 * \brief Return a boundary link of the triangle face
1654 * \param links - set of all links
1655 * \param avoidLink - link not to return
1656 * \param notBoundaryLink - out, neither the returned link nor avoidLink
1657 * \param nodeToContain - node the returned link must contain; if provided, search
1658 * also performed on adjacent faces
1659 * \param isAdjacentUsed - returns true if link is found in adjacent faces
1661 //================================================================================
1663 TLinkInSet QFace::GetBoundaryLink( const TLinkSet& links,
1664 const TChainLink& avoidLink,
1665 TLinkInSet * notBoundaryLink,
1666 const SMDS_MeshNode* nodeToContain,
1667 bool * isAdjacentUsed) const
1669 TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1671 typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1672 TFaceLinkList adjacentFaces;
1674 for ( int iL = 0; iL < _sides.size(); ++iL )
1676 if ( avoidLink._qlink == _sides[iL] )
1678 TLinkInSet link = links.find( _sides[iL] );
1679 if ( link == linksEnd ) continue;
1682 if ( link->IsBoundary() ) {
1683 if ( !nodeToContain ||
1684 (*link)->node1() == nodeToContain ||
1685 (*link)->node2() == nodeToContain )
1687 boundaryLink = link;
1688 if ( !notBoundaryLink ) break;
1691 else if ( notBoundaryLink ) {
1692 *notBoundaryLink = link;
1693 if ( boundaryLink != linksEnd ) break;
1696 if ( boundaryLink == linksEnd && nodeToContain ) // cellect adjacent faces
1697 if ( const QFace* adj = link->NextFace( this ))
1698 if ( adj->Contains( nodeToContain ))
1699 adjacentFaces.push_back( make_pair( adj, link ));
1702 if ( isAdjacentUsed ) *isAdjacentUsed = false;
1703 if ( boundaryLink == linksEnd && nodeToContain ) // check adjacent faces
1705 TFaceLinkList::iterator adj = adjacentFaces.begin();
1706 for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1707 boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second),
1708 0, nodeToContain, isAdjacentUsed);
1709 if ( isAdjacentUsed ) *isAdjacentUsed = true;
1711 return boundaryLink;
1713 //================================================================================
1715 * \brief Return a link ending at the given node but not avoidLink
1717 //================================================================================
1719 TLinkInSet QFace::GetLinkByNode( const TLinkSet& links,
1720 const TChainLink& avoidLink,
1721 const SMDS_MeshNode* nodeToContain) const
1723 for ( int i = 0; i < _sides.size(); ++i )
1724 if ( avoidLink._qlink != _sides[i] &&
1725 (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1726 return links.find( _sides[ i ]);
1730 //================================================================================
1732 * \brief Return normal to the i-th side pointing outside the face
1734 //================================================================================
1736 gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1738 gp_Vec norm, vecOut;
1739 // if ( uvHelper ) {
1740 // TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
1741 // const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
1742 // gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
1743 // gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
1744 // norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
1746 // const QLink* otherLink = _sides[(i + 1) % _sides.size()];
1747 // const SMDS_MeshNode* otherNode =
1748 // otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
1749 // gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
1750 // vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
1753 norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
1754 gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
1755 XYZ( _sides[0]->node2() ) +
1756 XYZ( _sides[1]->node1() )) / 3.;
1757 vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
1759 if ( norm * vecOut < 0 )
1761 double mag2 = norm.SquareMagnitude();
1762 if ( mag2 > numeric_limits<double>::min() )
1763 norm /= sqrt( mag2 );
1766 //================================================================================
1768 * \brief Move medium node of theLink according to its distance from boundary
1769 * \param theLink - link to fix
1770 * \param theRefVec - movement of boundary
1771 * \param theLinks - all adjacent links of continous triangles
1772 * \param theFaceHelper - helper is not used so far
1773 * \param thePrevLen - distance from the boundary
1774 * \param theStep - number of steps till movement propagation limit
1775 * \param theLinkNorm - out normal to theLink
1776 * \param theSign - 1 or -1 depending on movement of boundary
1777 * \retval double - distance from boundary to propagation limit or other boundary
1779 //================================================================================
1781 double QFace::MoveByBoundary( const TChainLink& theLink,
1782 const gp_Vec& theRefVec,
1783 const TLinkSet& theLinks,
1784 SMESH_MesherHelper* theFaceHelper,
1785 const double thePrevLen,
1787 gp_Vec* theLinkNorm,
1788 double theSign) const
1791 return thePrevLen; // propagation limit reached
1793 int iL; // index of theLink
1794 for ( iL = 0; iL < _sides.size(); ++iL )
1795 if ( theLink._qlink == _sides[ iL ])
1798 MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
1799 <<" thePrevLen " << thePrevLen);
1800 MSG(string(theStep,'.')<<" "<<*theLink._qlink);
1802 gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
1803 double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
1804 if ( theStep == theFirstStep )
1805 theSign = refProj < 0. ? -1. : 1.;
1806 else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
1807 return thePrevLen; // to propagate movement forward only, not in side dir or backward
1809 int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
1810 TLinkInSet link1 = theLinks.find( _sides[iL1] );
1811 TLinkInSet link2 = theLinks.find( _sides[iL2] );
1812 const QFace* f1 = link1->NextFace( this ); // adjacent faces
1813 const QFace* f2 = link2->NextFace( this );
1815 // propagate to adjacent faces till limit step or boundary
1816 double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
1817 double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
1818 gp_Vec linkDir1, linkDir2;
1822 len1 = f1->MoveByBoundary
1823 ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
1825 linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
1827 MSG( " --------------- EXCEPTION");
1833 len2 = f2->MoveByBoundary
1834 ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
1836 linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
1838 MSG( " --------------- EXCEPTION");
1843 if ( theStep != theFirstStep )
1845 // choose chain length by direction of propagation most codirected with theRefVec
1846 bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
1847 fullLen = choose1 ? len1 : len2;
1848 double r = thePrevLen / fullLen;
1850 gp_Vec move = linkNorm * refProj * ( 1 - r );
1851 theLink->Move( move, true );
1853 MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
1854 " by " << refProj * ( 1 - r ) << " following " <<
1855 (choose1 ? *link1->_qlink : *link2->_qlink));
1857 if ( theLinkNorm ) *theLinkNorm = linkNorm;
1862 //================================================================================
1864 * \brief Find pairs of continues faces
1866 //================================================================================
1868 void QLink::SetContinuesFaces() const
1870 // x0 x - QLink, [-|] - QFace, v - volume
1872 // | Between _faces of link x2 two vertical faces are continues
1873 // x1----x2-----x3 and two horizontal faces are continues. We set vertical faces
1874 // | to _faces[0] and _faces[1] and horizontal faces to
1875 // v2 | v3 _faces[2] and _faces[3] (or vise versa).
1878 if ( _faces.empty() )
1881 for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
1883 // look for a face bounding none of volumes bound by _faces[0]
1884 bool sameVol = false;
1885 int nbVol = _faces[iF]->NbVolumes();
1886 for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
1887 sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
1888 _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
1892 if ( iFaceCont > 0 ) // continues faces found, set one by the other
1894 if ( iFaceCont != 1 )
1895 std::swap( _faces[1], _faces[iFaceCont] );
1897 else if ( _faces.size() > 1 ) // not found, set NULL by the first face
1899 _faces.insert( ++_faces.begin(), 0 );
1902 //================================================================================
1904 * \brief Return a face continues to the given one
1906 //================================================================================
1908 const QFace* QLink::GetContinuesFace( const QFace* face ) const
1910 for ( int i = 0; i < _faces.size(); ++i ) {
1911 if ( _faces[i] == face ) {
1912 int iF = i < 2 ? 1-i : 5-i;
1913 return iF < _faces.size() ? _faces[iF] : 0;
1918 //================================================================================
1920 * \brief True if link is on mesh boundary
1922 //================================================================================
1924 bool QLink::OnBoundary() const
1926 for ( int i = 0; i < _faces.size(); ++i )
1927 if (_faces[i] && _faces[i]->IsBoundary()) return true;
1930 //================================================================================
1932 * \brief Return normal of link of the chain
1934 //================================================================================
1936 gp_Vec TChainLink::Normal() const {
1938 if (_qfaces[0]) norm = _qfaces[0]->_normal;
1939 if (_qfaces[1]) norm += _qfaces[1]->_normal;
1942 //================================================================================
1944 * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
1946 //================================================================================
1948 void fixPrism( TChain& allLinks )
1950 // separate boundary links from internal ones
1951 typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
1952 QLinkSet interLinks, bndLinks1, bndLink2;
1954 bool isCurved = false;
1955 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
1956 if ( (*lnk)->OnBoundary() )
1957 bndLinks1.insert( lnk->_qlink );
1959 interLinks.insert( lnk->_qlink );
1960 isCurved = isCurved || !(*lnk)->IsStraight();
1963 return; // no need to move
1965 QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
1967 while ( !interLinks.empty() && !curBndLinks->empty() )
1969 // propagate movement from boundary links to connected internal links
1970 QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
1971 for ( ; bnd != bndEnd; ++bnd )
1973 const QLink* bndLink = *bnd;
1974 for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
1976 const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
1977 if ( !face ) continue;
1978 // find and move internal link opposite to bndLink within the face
1979 int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
1980 const QLink* interLink = face->_sides[ interInd ];
1981 QLinkSet::iterator pInterLink = interLinks.find( interLink );
1982 if ( pInterLink == interLinks.end() ) continue; // not internal link
1983 interLink->Move( bndLink->_nodeMove );
1984 // treated internal links become new boundary ones
1985 interLinks. erase( pInterLink );
1986 newBndLinks->insert( interLink );
1989 curBndLinks->clear();
1990 std::swap( curBndLinks, newBndLinks );
1994 //================================================================================
1996 * \brief Fix links of continues triangles near curved boundary
1998 //================================================================================
2000 void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2002 if ( allLinks.empty() ) return;
2004 TLinkSet linkSet( allLinks.begin(), allLinks.end());
2005 TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2007 // move in 2d if we are on geom face
2008 // TopoDS_Face face;
2009 // TopLoc_Location loc;
2010 // SMESH_MesherHelper faceHelper( *helper.GetMesh());
2011 // while ( linkIt->IsBoundary()) ++linkIt;
2012 // if ( linkIt == linksEnd ) return;
2013 // if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
2014 // bool checkPos = true;
2015 // TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
2016 // if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2017 // face = TopoDS::Face( f );
2018 // helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
2022 // faceHelper.SetSubShape( face );
2025 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2027 if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2029 // if ( !face.IsNull() ) {
2030 // const SMDS_MeshNode* inFaceNode =
2031 // faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
2032 // gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
2033 // gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
2034 // gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
2035 // gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2036 // gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
2037 // linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
2040 linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2046 //================================================================================
2048 * \brief Detect rectangular structure of links and build chains from them
2050 //================================================================================
2052 enum TSplitTriaResult {
2053 _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2054 _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK };
2056 TSplitTriaResult splitTrianglesIntoChains( TChain & allLinks,
2057 vector< TChain> & resultChains,
2058 SMDS_TypeOfPosition pos )
2060 // put links in the set and evalute number of result chains by number of boundary links
2063 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2064 linkSet.insert( *lnk );
2065 nbBndLinks += lnk->IsBoundary();
2067 resultChains.clear();
2068 resultChains.reserve( nbBndLinks / 2 );
2070 TLinkInSet linkIt, linksEnd = linkSet.end();
2072 // find a boundary link with corner node; corner node has position pos-2
2073 // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2075 SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2076 const SMDS_MeshNode* corner = 0;
2077 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2078 if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2083 TLinkInSet startLink = linkIt;
2084 const SMDS_MeshNode* startCorner = corner;
2085 vector< TChain* > rowChains;
2088 while ( startLink != linksEnd) // loop on columns
2090 // We suppose we have a rectangular structure like shown here. We have found a
2091 // corner of the rectangle (startCorner) and a boundary link sharing
2092 // |/ |/ | the startCorner (startLink). We are going to loop on rows of the
2093 // --o---o---o structure making several chains at once. One chain (columnChain)
2094 // |\ | /| starts at startLink and continues upward (we look at the structure
2095 // \ | \ | / | from such point that startLink is on the bottom of the structure).
2096 // \| \|/ | While going upward we also fill horizontal chains (rowChains) we
2097 // --o---o---o encounter.
2099 // / | \ | \ | startCorner
2104 if ( resultChains.size() == nbBndLinks / 2 )
2106 resultChains.push_back( TChain() );
2107 TChain& columnChain = resultChains.back();
2109 TLinkInSet botLink = startLink; // current horizontal link to go up from
2110 corner = startCorner; // current corner the botLink ends at
2112 while ( botLink != linksEnd ) // loop on rows
2114 // add botLink to the columnChain
2115 columnChain.push_back( *botLink );
2117 const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2119 { // the column ends
2120 linkSet.erase( botLink );
2121 if ( iRow != rowChains.size() )
2122 return _FEW_ROWS; // different nb of rows in columns
2125 // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2126 // link ending at <corner> (sideLink); there are two cases:
2127 // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2128 // since midQuadLink is not at boundary while sideLink is.
2129 // 2) midQuadLink ends at <corner>
2131 TLinkInSet midQuadLink = linksEnd;
2132 TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2134 if ( isCase2 ) { // find midQuadLink among links of botTria
2135 midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2136 if ( midQuadLink->IsBoundary() )
2137 return _BAD_MIDQUAD;
2139 if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2140 return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2143 columnChain.push_back( *midQuadLink );
2144 if ( iRow >= rowChains.size() ) {
2146 return _MANY_ROWS; // different nb of rows in columns
2147 if ( resultChains.size() == nbBndLinks / 2 )
2149 resultChains.push_back( TChain() );
2150 rowChains.push_back( & resultChains.back() );
2152 rowChains[iRow]->push_back( *sideLink );
2153 rowChains[iRow]->push_back( *midQuadLink );
2155 const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2159 // prepare startCorner and startLink for the next column
2160 startCorner = startLink->NextNode( startCorner );
2162 startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2164 startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2165 // check if no more columns remains
2166 if ( startLink != linksEnd ) {
2167 const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2168 if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2169 startLink = linksEnd; // startLink bounds upTria or botTria
2170 else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2174 // find bottom link and corner for the next row
2175 corner = sideLink->NextNode( corner );
2176 // next bottom link ends at the new corner
2177 linkSet.erase( botLink );
2178 botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2179 if ( botLink == linksEnd || botLink == (isCase2 ? midQuadLink : sideLink))
2181 linkSet.erase( midQuadLink );
2182 linkSet.erase( sideLink );
2184 // make faces neighboring the found ones be boundary
2185 if ( startLink != linksEnd ) {
2186 const QFace* tria = isCase2 ? botTria : upTria;
2187 for ( int iL = 0; iL < 3; ++iL ) {
2188 linkIt = linkSet.find( tria->_sides[iL] );
2189 if ( linkIt != linksEnd )
2190 linkIt->RemoveFace( tria );
2193 if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2194 botLink->RemoveFace( upTria ); // make next botTria first in vector
2201 // In the linkSet, there must remain the last links of rowChains; add them
2202 if ( linkSet.size() != rowChains.size() )
2203 return _BAD_SET_SIZE;
2204 for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2205 // find the link (startLink) ending at startCorner
2207 for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2208 if ( (*startLink)->node1() == startCorner ) {
2209 corner = (*startLink)->node2(); break;
2211 else if ( (*startLink)->node2() == startCorner) {
2212 corner = (*startLink)->node1(); break;
2215 if ( startLink == linksEnd )
2217 rowChains[ iRow ]->push_back( *startLink );
2218 linkSet.erase( startLink );
2219 startCorner = corner;
2226 //=======================================================================
2228 * \brief Move medium nodes of faces and volumes to fix distorted elements
2229 * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2231 * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2233 //=======================================================================
2235 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2237 // apply algorithm to solids or geom faces
2238 // ----------------------------------------------
2239 if ( myShape.IsNull() ) {
2240 if ( !myMesh->HasShapeToMesh() ) return;
2241 SetSubShape( myMesh->GetShapeToMesh() );
2243 TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2244 for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2245 faces.Add( f.Current() );
2247 for ( TopExp_Explorer v(myShape,TopAbs_SOLID); v.More(); v.Next() ) {
2248 if ( myMesh->GetSubMesh( v.Current() )->IsEmpty() ) { // get faces of solid
2249 for ( TopExp_Explorer f( v.Current(), TopAbs_FACE); f.More(); f.Next() )
2250 faces.Add( f.Current() );
2252 else { // fix nodes in the solid and its faces
2253 SMESH_MesherHelper h(*myMesh);
2254 h.SetSubShape( v.Current() );
2255 h.FixQuadraticElements(false);
2258 // fix nodes on geom faces
2259 for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2260 SMESH_MesherHelper h(*myMesh);
2261 h.SetSubShape( fIt.Key() );
2262 h.FixQuadraticElements();
2267 // Find out type of elements and get iterator on them
2268 // ---------------------------------------------------
2270 SMDS_ElemIteratorPtr elemIt;
2271 SMDSAbs_ElementType elemType = SMDSAbs_All;
2273 SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2276 if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2277 elemIt = smDS->GetElements();
2278 if ( elemIt->more() ) {
2279 elemType = elemIt->next()->GetType();
2280 elemIt = smDS->GetElements();
2283 if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2286 // Fill in auxiliary data structures
2287 // ----------------------------------
2291 set< QLink >::iterator pLink;
2292 set< QFace >::iterator pFace;
2294 bool isCurved = false;
2295 bool hasRectFaces = false;
2296 set<int> nbElemNodeSet;
2298 if ( elemType == SMDSAbs_Volume )
2300 SMDS_VolumeTool volTool;
2301 while ( elemIt->more() ) // loop on volumes
2303 const SMDS_MeshElement* vol = elemIt->next();
2304 if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2306 for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2308 int nbN = volTool.NbFaceNodes( iF );
2309 nbElemNodeSet.insert( nbN );
2310 const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2311 vector< const QLink* > faceLinks( nbN/2 );
2312 for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2315 QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2316 pLink = links.insert( link ).first;
2317 faceLinks[ iN/2 ] = & *pLink;
2319 isCurved = !link.IsStraight();
2320 if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2321 return; // already fixed
2324 pFace = faces.insert( QFace( faceLinks )).first;
2325 if ( pFace->NbVolumes() == 0 )
2326 pFace->AddSelfToLinks();
2327 pFace->SetVolume( vol );
2328 hasRectFaces = hasRectFaces ||
2329 ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2330 volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2333 set< QLink >::iterator pLink = links.begin();
2334 for ( ; pLink != links.end(); ++pLink )
2335 pLink->SetContinuesFaces();
2339 while ( elemIt->more() ) // loop on faces
2341 const SMDS_MeshElement* face = elemIt->next();
2342 if ( !face->IsQuadratic() )
2344 nbElemNodeSet.insert( face->NbNodes() );
2345 int nbN = face->NbNodes()/2;
2346 vector< const QLink* > faceLinks( nbN );
2347 for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2350 QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2351 pLink = links.insert( link ).first;
2352 faceLinks[ iN ] = & *pLink;
2354 isCurved = !link.IsStraight();
2357 pFace = faces.insert( QFace( faceLinks )).first;
2358 pFace->AddSelfToLinks();
2359 hasRectFaces = ( hasRectFaces || nbN == 4 );
2363 return; // no curved edges of faces
2365 // Compute displacement of medium nodes
2366 // -------------------------------------
2368 // two loops on faces: the first is to treat boundary links, the second is for internal ones
2369 TopLoc_Location loc;
2370 // not treat boundary of volumic submesh
2371 int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2372 for ( ; isInside < 2; ++isInside ) {
2373 MSG( "--------------- LOOP " << isInside << " ------------------");
2374 SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2376 for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2377 if ( bool(isInside) == pFace->IsBoundary() )
2379 for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2382 // make chain of links connected via continues faces
2385 if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2387 if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2389 vector< TChain > chains;
2390 if ( error == ERR_OK ) { // chains contains continues rectangles
2392 chains[0].splice( chains[0].begin(), rawChain );
2394 else if ( error == ERR_TRI ) { // chains contains continues triangles
2395 TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2396 if ( res != _OK ) { // not rectangles split into triangles
2397 fixTriaNearBoundary( rawChain, *this );
2401 else if ( error == ERR_PRISM ) { // side faces of prisms
2402 fixPrism( rawChain );
2408 for ( int iC = 0; iC < chains.size(); ++iC )
2410 TChain& chain = chains[iC];
2411 if ( chain.empty() ) continue;
2412 if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2416 // mesure chain length and compute link position along the chain
2417 double chainLen = 0;
2418 vector< double > linkPos;
2419 MSGBEG( "Link medium nodes: ");
2420 TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2421 for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2422 MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2423 double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2424 while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2425 link1 = chain.erase( link1 );
2426 if ( link1 == chain.end() )
2428 len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2431 linkPos.push_back( chainLen );
2434 if ( linkPos.size() < 2 )
2437 gp_Vec move0 = chain.front()->_nodeMove;
2438 gp_Vec move1 = chain.back ()->_nodeMove;
2441 bool checkUV = true;
2443 // compute node displacement of end links in parametric space of face
2444 const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2445 TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2446 if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2447 face = TopoDS::Face( f );
2448 for ( int is1 = 0; is1 < 2; ++is1 ) { // move0 or move1
2449 TChainLink& link = is1 ? chain.back() : chain.front();
2450 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2451 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2452 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2453 gp_XY uvMove = uvm - GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2454 if ( is1 ) move1.SetCoord( uvMove.X(), uvMove.Y(), 0 );
2455 else move0.SetCoord( uvMove.X(), uvMove.Y(), 0 );
2457 if ( move0.SquareMagnitude() < straightTol2 &&
2458 move1.SquareMagnitude() < straightTol2 ) {
2460 continue; // straight - no need to move nodes of internal links
2465 if ( isInside || face.IsNull() )
2467 // compute node displacement of end links in their local coord systems
2469 TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2470 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2471 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2472 move0.Transform(trsf);
2475 TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2476 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2477 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2478 move1.Transform(trsf);
2481 // compute displacement of medium nodes
2482 link2 = chain.begin();
2485 for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2487 double r = linkPos[i] / chainLen;
2488 // displacement in local coord system
2489 gp_Vec move = (1. - r) * move0 + r * move1;
2490 if ( isInside || face.IsNull()) {
2491 // transform to global
2492 gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2493 gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2494 gp_Vec x = x01.Normalized() + x12.Normalized();
2495 trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2496 move.Transform(trsf);
2499 // compute 3D displacement by 2D one
2500 gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2501 gp_XY newUV = oldUV + gp_XY( move.X(), move.Y() );
2502 gp_Pnt newPnt = BRep_Tool::Surface(face,loc)->Value( newUV.X(), newUV.Y());
2503 move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2505 if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2506 move.SquareMagnitude())
2508 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2509 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2510 MSG( "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2511 "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2512 "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2513 "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2517 (*link1)->Move( move );
2518 MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2519 << chain.front()->_mediumNode->GetID() <<"-"
2520 << chain.back ()->_mediumNode->GetID() <<
2521 " by " << move.Magnitude());
2523 } // loop on chains of links
2524 } // loop on 2 directions of propagation from quadrangle
2531 for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2532 if ( pLink->IsMoved() ) {
2533 //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2534 gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2535 GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());