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>
57 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
61 gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
65 //================================================================================
69 //================================================================================
71 SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
72 : myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false), myCheckNodePos(false)
74 mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
77 //=======================================================================
78 //function : CheckShape
80 //=======================================================================
82 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
84 SMESHDS_Mesh* meshDS = GetMeshDS();
85 // we can create quadratic elements only if all elements
86 // created on subshapes of given shape are quadratic
87 // also we have to fill myTLinkNodeMap
88 myCreateQuadratic = true;
89 mySeamShapeIds.clear();
90 myDegenShapeIds.clear();
91 TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
92 SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
94 int nbOldLinks = myTLinkNodeMap.size();
96 TopExp_Explorer exp( aSh, subType );
97 for (; exp.More() && myCreateQuadratic; exp.Next()) {
98 if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
99 if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
101 const SMDS_MeshElement* e = it->next();
102 if ( e->GetType() != elemType || !e->IsQuadratic() ) {
103 myCreateQuadratic = false;
108 switch ( e->NbNodes() ) {
110 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
112 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
113 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
114 AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
116 AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
117 AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
118 AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
119 AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
122 myCreateQuadratic = false;
131 if ( nbOldLinks == myTLinkNodeMap.size() )
132 myCreateQuadratic = false;
134 if(!myCreateQuadratic) {
135 myTLinkNodeMap.clear();
139 return myCreateQuadratic;
142 //================================================================================
144 * \brief Set geomerty to make elements on
145 * \param aSh - geomertic shape
147 //================================================================================
149 void SMESH_MesherHelper::SetSubShape(const int aShID)
151 if ( aShID == myShapeID )
154 SetSubShape( GetMeshDS()->IndexToShape( aShID ));
156 SetSubShape( TopoDS_Shape() );
159 //================================================================================
161 * \brief Set geomerty to make elements on
162 * \param aSh - geomertic shape
164 //================================================================================
166 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
168 if ( myShape.IsSame( aSh ))
172 mySeamShapeIds.clear();
173 myDegenShapeIds.clear();
175 if ( myShape.IsNull() ) {
179 SMESHDS_Mesh* meshDS = GetMeshDS();
180 myShapeID = meshDS->ShapeToIndex(aSh);
182 // treatment of periodic faces
183 for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
185 const TopoDS_Face& face = TopoDS::Face( eF.Current() );
186 BRepAdaptor_Surface surface( face );
187 if ( surface.IsUPeriodic() || surface.IsVPeriodic() )
189 for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
191 // look for a seam edge
192 const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
193 if ( BRep_Tool::IsClosed( edge, face )) {
194 // initialize myPar1, myPar2 and myParIndex
195 if ( mySeamShapeIds.empty() ) {
197 BRep_Tool::UVPoints( edge, face, uv1, uv2 );
198 if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
200 myParIndex = 1; // U periodic
201 myPar1 = surface.FirstUParameter();
202 myPar2 = surface.LastUParameter();
205 myParIndex = 2; // V periodic
206 myPar1 = surface.FirstVParameter();
207 myPar2 = surface.LastVParameter();
210 // store seam shape indices, negative if shape encounters twice
211 int edgeID = meshDS->ShapeToIndex( edge );
212 mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
213 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
214 int vertexID = meshDS->ShapeToIndex( v.Current() );
215 mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
219 // look for a degenerated edge
220 if ( BRep_Tool::Degenerated( edge )) {
221 myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
222 for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
223 myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
230 //================================================================================
232 * \brief Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
233 * \param F - the face
234 * \retval bool - return true if the face is periodic
236 //================================================================================
238 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
240 if ( F.IsNull() ) return !mySeamShapeIds.empty();
242 if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
243 return !mySeamShapeIds.empty();
246 Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
247 if ( !aSurface.IsNull() )
248 return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
253 //=======================================================================
254 //function : IsMedium
256 //=======================================================================
258 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode* node,
259 const SMDSAbs_ElementType typeToCheck)
261 return SMESH_MeshEditor::IsMedium( node, typeToCheck );
264 //=======================================================================
266 * \brief Return support shape of a node
267 * \param node - the node
268 * \param meshDS - mesh DS
269 * \retval TopoDS_Shape - found support shape
271 //=======================================================================
273 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
274 SMESHDS_Mesh* meshDS)
276 int shapeID = node->GetPosition()->GetShapeId();
277 if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
278 return meshDS->IndexToShape( shapeID );
280 return TopoDS_Shape();
284 //=======================================================================
285 //function : AddTLinkNode
287 //=======================================================================
289 * Auxilary function for filling myTLinkNodeMap
291 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
292 const SMDS_MeshNode* n2,
293 const SMDS_MeshNode* n12)
295 // add new record to map
296 SMESH_TLink link( n1, n2 );
297 myTLinkNodeMap.insert( make_pair(link,n12));
300 //=======================================================================
302 * \brief Select UV on either of 2 pcurves of a seam edge, closest to the given UV
303 * \param uv1 - UV on the seam
304 * \param uv2 - UV within a face
305 * \retval gp_Pnt2d - selected UV
307 //=======================================================================
309 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
311 double p1 = uv1.Coord( myParIndex );
312 double p2 = uv2.Coord( myParIndex );
313 double p3 = ( Abs( p1 - myPar1 ) < Abs( p1 - myPar2 )) ? myPar2 : myPar1;
314 if ( Abs( p2 - p1 ) > Abs( p2 - p3 ))
316 gp_Pnt2d result = uv1;
317 result.SetCoord( myParIndex, p1 );
321 //=======================================================================
323 * \brief Return node UV on face
324 * \param F - the face
325 * \param n - the node
326 * \param n2 - a node of element being created located inside a face
327 * \retval gp_XY - resulting UV
329 * Auxilary function called form GetMediumNode()
331 //=======================================================================
333 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face& F,
334 const SMDS_MeshNode* n,
335 const SMDS_MeshNode* n2,
338 gp_Pnt2d uv( 1e100, 1e100 );
339 const SMDS_PositionPtr Pos = n->GetPosition();
340 if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
342 // node has position on face
343 const SMDS_FacePosition* fpos =
344 static_cast<const SMDS_FacePosition*>(n->GetPosition().get());
345 uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
346 if ( check && *check )
348 // check that uv is correct
350 Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
351 double tol = 2 * BRep_Tool::Tolerance( F );
352 gp_Pnt nodePnt = XYZ( n );
353 if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
354 if ( nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > tol ) {
355 // uv incorrect, project the node to surface
356 GeomAPI_ProjectPointOnSurf projector( nodePnt, surface, tol );
357 if ( !projector.IsDone() || projector.NbPoints() < 1 ) {
358 MESSAGE( "SMESH_MesherHelper::GetNodeUV() failed to project" )
361 Quantity_Parameter U,V;
362 projector.LowerDistanceParameters(U,V);
363 if ( nodePnt.Distance( surface->Value( U, V )) > tol )
364 MESSAGE( "SMESH_MesherHelper::GetNodeUV(), invalid projection" );
367 else if ( uv.XY().Modulus() > numeric_limits<double>::min() ) {
368 *check = false; // parameters are OK, do not check further more
372 else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
374 // node has position on edge => it is needed to find
375 // corresponding edge from face, get pcurve for this
376 // edge and retrieve value from this pcurve
377 const SMDS_EdgePosition* epos =
378 static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
379 int edgeID = Pos->GetShapeId();
380 TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
382 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
383 uv = C2d->Value( epos->GetUParameter() );
384 // for a node on a seam edge select one of UVs on 2 pcurves
385 if ( n2 && IsSeamShape( edgeID ) )
386 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
388 // adjust uv to period
390 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
391 Standard_Boolean isUPeriodic = S->IsUPeriodic();
392 Standard_Boolean isVPeriodic = S->IsVPeriodic();
393 if ( isUPeriodic || isVPeriodic ) {
394 Standard_Real UF,UL,VF,VL;
395 S->Bounds(UF,UL,VF,VL);
397 uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
399 uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
402 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
404 if ( int vertexID = n->GetPosition()->GetShapeId() ) {
406 const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
408 uv = BRep_Tool::Parameters( V, F );
410 catch (Standard_Failure& exc) {
414 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !ok && vert.More(); vert.Next() )
415 ok = ( V == vert.Current() );
418 MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
419 << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
421 // get UV of a vertex closest to the node
423 gp_Pnt pn = XYZ( n );
424 for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !ok && vert.More(); vert.Next() ) {
425 TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
426 gp_Pnt p = BRep_Tool::Pnt( curV );
427 double curDist = p.SquareDistance( pn );
428 if ( curDist < dist ) {
430 uv = BRep_Tool::Parameters( curV, F );
431 if ( dist < DBL_MIN ) break;
436 TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
437 for ( ; it.More(); it.Next() ) {
438 if ( it.Value().ShapeType() == TopAbs_EDGE ) {
439 const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
441 Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
442 if ( !C2d.IsNull() ) {
443 double u = ( V == TopExp::FirstVertex( edge ) ) ? f : l;
444 uv = C2d->Value( u );
451 if ( n2 && IsSeamShape( vertexID ) )
452 uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
458 //=======================================================================
460 * \brief Return middle UV taking in account surface period
462 //=======================================================================
464 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
468 if ( surface.IsNull() )
469 return 0.5 * ( p1 + p2 );
470 //checking if surface is periodic
471 Standard_Real UF,UL,VF,VL;
472 surface->Bounds(UF,UL,VF,VL);
475 Standard_Boolean isUPeriodic = surface->IsUPeriodic();
477 Standard_Real UPeriod = surface->UPeriod();
478 Standard_Real p2x = p2.X()+ShapeAnalysis::AdjustByPeriod(p2.X(),p1.X(),UPeriod);
479 Standard_Real pmid = (p1.X()+p2x)/2.;
480 u = pmid+ShapeAnalysis::AdjustToPeriod(pmid,UF,UL);
483 u= (p1.X()+p2.X())/2.;
485 Standard_Boolean isVPeriodic = surface->IsVPeriodic();
487 Standard_Real VPeriod = surface->VPeriod();
488 Standard_Real p2y = p2.Y()+ShapeAnalysis::AdjustByPeriod(p2.Y(),p1.Y(),VPeriod);
489 Standard_Real pmid = (p1.Y()+p2y)/2.;
490 v = pmid+ShapeAnalysis::AdjustToPeriod(pmid,VF,VL);
493 v = (p1.Y()+p2.Y())/2.;
498 //=======================================================================
500 * \brief Return node U on edge
501 * \param E - the Edge
502 * \param n - the node
503 * \retval double - resulting U
505 * Auxilary function called form GetMediumNode()
507 //=======================================================================
509 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge& E,
510 const SMDS_MeshNode* n,
514 const SMDS_PositionPtr Pos = n->GetPosition();
515 if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE) {
516 const SMDS_EdgePosition* epos =
517 static_cast<const SMDS_EdgePosition*>(n->GetPosition().get());
518 param = epos->GetUParameter();
520 else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX) {
521 SMESHDS_Mesh * meshDS = GetMeshDS();
522 int vertexID = n->GetPosition()->GetShapeId();
523 const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
524 param = BRep_Tool::Parameter( V, E );
529 //================================================================================
531 * \brief Return existing or create new medium nodes between given ones
532 * \param force3d - if true, new node is the middle of n1 and n2,
533 * else is located on geom face or geom edge
535 //================================================================================
537 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
538 const SMDS_MeshNode* n2,
541 SMESH_TLink link(n1,n2);
542 ItTLinkNode itLN = myTLinkNodeMap.find( link );
543 if ( itLN != myTLinkNodeMap.end() ) {
544 return (*itLN).second;
546 // create medium node
548 SMESHDS_Mesh* meshDS = GetMeshDS();
549 int faceID = -1, edgeID = -1;
550 const SMDS_PositionPtr Pos1 = n1->GetPosition();
551 const SMDS_PositionPtr Pos2 = n2->GetPosition();
553 if( myShape.IsNull() )
555 if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
556 faceID = Pos1->GetShapeId();
558 else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
559 faceID = Pos2->GetShapeId();
562 if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
563 edgeID = Pos1->GetShapeId();
565 if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
566 edgeID = Pos2->GetShapeId();
571 // we try to create medium node using UV parameters of
572 // nodes, else - medium between corresponding 3d points
574 TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
575 if(faceID>0 || shapeType == TopAbs_FACE) {
576 // obtaining a face and 2d points for nodes
578 if( myShape.IsNull() )
579 F = TopoDS::Face(meshDS->IndexToShape(faceID));
581 F = TopoDS::Face(myShape);
585 gp_XY p1 = GetNodeUV(F,n1,n2, &myCheckNodePos);
586 gp_XY p2 = GetNodeUV(F,n2,n1, &myCheckNodePos);
588 if ( IsDegenShape( Pos1->GetShapeId() ))
589 p1.SetCoord( myParIndex, p2.Coord( myParIndex ));
590 else if ( IsDegenShape( Pos2->GetShapeId() ))
591 p2.SetCoord( myParIndex, p1.Coord( myParIndex ));
594 Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
595 gp_XY uv = GetMiddleUV( S, p1, p2 );
596 gp_Pnt P = S->Value( uv.X(), uv.Y() ).Transformed(loc);
597 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
598 meshDS->SetNodeOnFace(n12, faceID, uv.X(), uv.Y());
599 myTLinkNodeMap.insert(make_pair(link,n12));
602 if (edgeID>0 || shapeType == TopAbs_EDGE) {
605 if( myShape.IsNull() )
606 E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
608 E = TopoDS::Edge(myShape);
612 double p1 = GetNodeU(E,n1, &myCheckNodePos);
613 double p2 = GetNodeU(E,n2, &myCheckNodePos);
616 Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
619 Standard_Boolean isPeriodic = C->IsPeriodic();
622 Standard_Real Period = C->Period();
623 Standard_Real p = p2+ShapeAnalysis::AdjustByPeriod(p2,p1,Period);
624 Standard_Real pmid = (p1+p)/2.;
625 u = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
630 gp_Pnt P = C->Value( u );
631 n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
632 meshDS->SetNodeOnEdge(n12, edgeID, u);
633 myTLinkNodeMap.insert(make_pair(link,n12));
639 double x = ( n1->X() + n2->X() )/2.;
640 double y = ( n1->Y() + n2->Y() )/2.;
641 double z = ( n1->Z() + n2->Z() )/2.;
642 n12 = meshDS->AddNode(x,y,z);
644 meshDS->SetNodeOnEdge(n12, edgeID);
646 meshDS->SetNodeOnFace(n12, faceID);
648 meshDS->SetNodeInVolume(n12, myShapeID);
649 myTLinkNodeMap.insert( make_pair( link, n12 ));
653 //=======================================================================
657 //=======================================================================
659 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
661 SMESHDS_Mesh * meshDS = GetMeshDS();
662 SMDS_MeshNode* node = 0;
664 node = meshDS->AddNodeWithID( x, y, z, ID );
666 node = meshDS->AddNode( x, y, z );
667 if ( mySetElemOnShape && myShapeID > 0 ) {
668 switch ( myShape.ShapeType() ) {
669 case TopAbs_SOLID: meshDS->SetNodeInVolume( node, myShapeID); break;
670 case TopAbs_SHELL: meshDS->SetNodeInVolume( node, myShapeID); break;
671 case TopAbs_FACE: meshDS->SetNodeOnFace( node, myShapeID); break;
672 case TopAbs_EDGE: meshDS->SetNodeOnEdge( node, myShapeID); break;
673 case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
680 //=======================================================================
682 * Creates quadratic or linear edge
684 //=======================================================================
686 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
687 const SMDS_MeshNode* n2,
691 SMESHDS_Mesh * meshDS = GetMeshDS();
693 SMDS_MeshEdge* edge = 0;
694 if (myCreateQuadratic) {
695 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
697 edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
699 edge = meshDS->AddEdge(n1, n2, n12);
703 edge = meshDS->AddEdgeWithID(n1, n2, id);
705 edge = meshDS->AddEdge(n1, n2);
708 if ( mySetElemOnShape && myShapeID > 0 )
709 meshDS->SetMeshElementOnShape( edge, myShapeID );
714 //=======================================================================
716 * Creates quadratic or linear triangle
718 //=======================================================================
720 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
721 const SMDS_MeshNode* n2,
722 const SMDS_MeshNode* n3,
726 SMESHDS_Mesh * meshDS = GetMeshDS();
727 SMDS_MeshFace* elem = 0;
728 if(!myCreateQuadratic) {
730 elem = meshDS->AddFaceWithID(n1, n2, n3, id);
732 elem = meshDS->AddFace(n1, n2, n3);
735 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
736 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
737 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
740 elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
742 elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
744 if ( mySetElemOnShape && myShapeID > 0 )
745 meshDS->SetMeshElementOnShape( elem, myShapeID );
750 //=======================================================================
752 * Creates quadratic or linear quadrangle
754 //=======================================================================
756 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
757 const SMDS_MeshNode* n2,
758 const SMDS_MeshNode* n3,
759 const SMDS_MeshNode* n4,
763 SMESHDS_Mesh * meshDS = GetMeshDS();
764 SMDS_MeshFace* elem = 0;
765 if(!myCreateQuadratic) {
767 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
769 elem = meshDS->AddFace(n1, n2, n3, n4);
772 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
773 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
774 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
775 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
778 elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
780 elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
782 if ( mySetElemOnShape && myShapeID > 0 )
783 meshDS->SetMeshElementOnShape( elem, myShapeID );
788 //=======================================================================
790 * Creates quadratic or linear volume
792 //=======================================================================
794 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
795 const SMDS_MeshNode* n2,
796 const SMDS_MeshNode* n3,
797 const SMDS_MeshNode* n4,
798 const SMDS_MeshNode* n5,
799 const SMDS_MeshNode* n6,
803 SMESHDS_Mesh * meshDS = GetMeshDS();
804 SMDS_MeshVolume* elem = 0;
805 if(!myCreateQuadratic) {
807 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
809 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
812 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
813 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
814 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
816 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
817 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
818 const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
820 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
821 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
822 const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
825 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
826 n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
828 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
829 n12, n23, n31, n45, n56, n64, n14, n25, n36);
831 if ( mySetElemOnShape && myShapeID > 0 )
832 meshDS->SetMeshElementOnShape( elem, myShapeID );
837 //=======================================================================
839 * Creates quadratic or linear volume
841 //=======================================================================
843 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
844 const SMDS_MeshNode* n2,
845 const SMDS_MeshNode* n3,
846 const SMDS_MeshNode* n4,
850 SMESHDS_Mesh * meshDS = GetMeshDS();
851 SMDS_MeshVolume* elem = 0;
852 if(!myCreateQuadratic) {
854 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
856 elem = meshDS->AddVolume(n1, n2, n3, n4);
859 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
860 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
861 const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
863 const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
864 const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
865 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
868 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
870 elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
872 if ( mySetElemOnShape && myShapeID > 0 )
873 meshDS->SetMeshElementOnShape( elem, myShapeID );
878 //=======================================================================
880 * Creates quadratic or linear pyramid
882 //=======================================================================
884 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
885 const SMDS_MeshNode* n2,
886 const SMDS_MeshNode* n3,
887 const SMDS_MeshNode* n4,
888 const SMDS_MeshNode* n5,
892 SMDS_MeshVolume* elem = 0;
893 if(!myCreateQuadratic) {
895 elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
897 elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
900 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
901 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
902 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
903 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
905 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
906 const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
907 const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
908 const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
911 elem = GetMeshDS()->AddVolumeWithID ( n1, n2, n3, n4, n5,
916 elem = GetMeshDS()->AddVolume( n1, n2, n3, n4, n5,
920 if ( mySetElemOnShape && myShapeID > 0 )
921 GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
926 //=======================================================================
928 * Creates quadratic or linear hexahedron
930 //=======================================================================
932 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
933 const SMDS_MeshNode* n2,
934 const SMDS_MeshNode* n3,
935 const SMDS_MeshNode* n4,
936 const SMDS_MeshNode* n5,
937 const SMDS_MeshNode* n6,
938 const SMDS_MeshNode* n7,
939 const SMDS_MeshNode* n8,
943 SMESHDS_Mesh * meshDS = GetMeshDS();
944 SMDS_MeshVolume* elem = 0;
945 if(!myCreateQuadratic) {
947 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
949 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
952 const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
953 const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
954 const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
955 const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
957 const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
958 const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
959 const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
960 const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
962 const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
963 const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
964 const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
965 const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
968 elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
969 n12, n23, n34, n41, n56, n67,
970 n78, n85, n15, n26, n37, n48, id);
972 elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
973 n12, n23, n34, n41, n56, n67,
974 n78, n85, n15, n26, n37, n48);
976 if ( mySetElemOnShape && myShapeID > 0 )
977 meshDS->SetMeshElementOnShape( elem, myShapeID );
982 //=======================================================================
984 * \brief Load nodes bound to face into a map of node columns
985 * \param theParam2ColumnMap - map of node columns to fill
986 * \param theFace - the face on which nodes are searched for
987 * \param theBaseEdge - the edge nodes of which are columns' bases
988 * \param theMesh - the mesh containing nodes
989 * \retval bool - false if something is wrong
991 * The key of the map is a normalized parameter of each
992 * base node on theBaseEdge.
993 * This method works in supposition that nodes on the face
994 * forms a rectangular grid and elements can be quardrangles or triangles
996 //=======================================================================
998 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
999 const TopoDS_Face& theFace,
1000 const TopoDS_Edge& theBaseEdge,
1001 SMESHDS_Mesh* theMesh)
1003 // get vertices of theBaseEdge
1004 TopoDS_Vertex vfb, vlb, vft; // first and last, bottom and top vertices
1005 TopoDS_Edge eFrw = TopoDS::Edge( theBaseEdge.Oriented( TopAbs_FORWARD ));
1006 TopExp::Vertices( eFrw, vfb, vlb );
1008 // find the other edges of theFace and orientation of e1
1009 TopoDS_Edge e1, e2, eTop;
1010 bool rev1, CumOri = false;
1011 TopExp_Explorer exp( theFace, TopAbs_EDGE );
1013 for ( ; exp.More(); exp.Next() ) {
1014 if ( ++nbEdges > 4 ) {
1015 return false; // more than 4 edges in theFace
1017 TopoDS_Edge e = TopoDS::Edge( exp.Current() );
1018 if ( theBaseEdge.IsSame( e ))
1020 TopoDS_Vertex vCommon;
1021 if ( !TopExp::CommonVertex( theBaseEdge, e, vCommon ))
1023 else if ( vCommon.IsSame( vfb )) {
1025 vft = TopExp::LastVertex( e1, CumOri );
1026 rev1 = vfb.IsSame( vft );
1028 vft = TopExp::FirstVertex( e1, CumOri );
1033 if ( nbEdges < 4 ) {
1034 return false; // less than 4 edges in theFace
1036 if ( e2.IsNull() && vfb.IsSame( vlb ))
1039 // submeshes corresponding to shapes
1040 SMESHDS_SubMesh* smFace = theMesh->MeshElements( theFace );
1041 SMESHDS_SubMesh* smb = theMesh->MeshElements( theBaseEdge );
1042 SMESHDS_SubMesh* smt = theMesh->MeshElements( eTop );
1043 SMESHDS_SubMesh* sm1 = theMesh->MeshElements( e1 );
1044 SMESHDS_SubMesh* sm2 = theMesh->MeshElements( e2 );
1045 SMESHDS_SubMesh* smVfb = theMesh->MeshElements( vfb );
1046 SMESHDS_SubMesh* smVlb = theMesh->MeshElements( vlb );
1047 SMESHDS_SubMesh* smVft = theMesh->MeshElements( vft );
1048 if (!smFace || !smb || !smt || !sm1 || !sm2 || !smVfb || !smVlb || !smVft ) {
1049 RETURN_BAD_RESULT( "NULL submesh " <<smFace<<" "<<smb<<" "<<smt<<" "<<
1050 sm1<<" "<<sm2<<" "<<smVfb<<" "<<smVlb<<" "<<smVft);
1052 if ( smb->NbNodes() != smt->NbNodes() || sm1->NbNodes() != sm2->NbNodes() ) {
1053 RETURN_BAD_RESULT(" Diff nb of nodes on opposite edges" );
1055 if (smVfb->NbNodes() != 1 || smVlb->NbNodes() != 1 || smVft->NbNodes() != 1) {
1056 RETURN_BAD_RESULT("Empty submesh of vertex");
1058 // define whether mesh is quadratic
1059 bool isQuadraticMesh = false;
1060 SMDS_ElemIteratorPtr eIt = smFace->GetElements();
1061 if ( !eIt->more() ) {
1062 RETURN_BAD_RESULT("No elements on the face");
1064 const SMDS_MeshElement* e = eIt->next();
1065 isQuadraticMesh = e->IsQuadratic();
1067 if ( sm1->NbNodes() * smb->NbNodes() != smFace->NbNodes() ) {
1068 // check quadratic case
1069 if ( isQuadraticMesh ) {
1070 // what if there are quadrangles and triangles mixed?
1071 // int n1 = sm1->NbNodes()/2;
1072 // int n2 = smb->NbNodes()/2;
1073 // int n3 = sm1->NbNodes() - n1;
1074 // int n4 = smb->NbNodes() - n2;
1075 // int nf = sm1->NbNodes()*smb->NbNodes() - n3*n4;
1076 // if( nf != smFace->NbNodes() ) {
1077 // MESSAGE( "Wrong nb face nodes: " <<
1078 // sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
1083 RETURN_BAD_RESULT( "Wrong nb face nodes: " <<
1084 sm1->NbNodes()<<" "<<smb->NbNodes()<<" "<<smFace->NbNodes());
1088 int vsize = sm1->NbNodes() + 2;
1089 int hsize = smb->NbNodes() + 2;
1090 if(isQuadraticMesh) {
1091 vsize = vsize - sm1->NbNodes()/2 -1;
1092 hsize = hsize - smb->NbNodes()/2 -1;
1095 // load nodes from theBaseEdge
1097 std::set<const SMDS_MeshNode*> loadedNodes;
1098 const SMDS_MeshNode* nullNode = 0;
1100 std::vector<const SMDS_MeshNode*> & nVecf = theParam2ColumnMap[ 0.];
1101 nVecf.resize( vsize, nullNode );
1102 loadedNodes.insert( nVecf[ 0 ] = smVfb->GetNodes()->next() );
1104 std::vector<const SMDS_MeshNode*> & nVecl = theParam2ColumnMap[ 1.];
1105 nVecl.resize( vsize, nullNode );
1106 loadedNodes.insert( nVecl[ 0 ] = smVlb->GetNodes()->next() );
1109 BRep_Tool::Range( eFrw, f, l );
1110 double range = l - f;
1111 SMDS_NodeIteratorPtr nIt = smb->GetNodes();
1112 const SMDS_MeshNode* node;
1113 while ( nIt->more() ) {
1115 if(IsMedium(node, SMDSAbs_Edge))
1117 const SMDS_EdgePosition* pos =
1118 dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
1122 double u = ( pos->GetUParameter() - f ) / range;
1123 std::vector<const SMDS_MeshNode*> & nVec = theParam2ColumnMap[ u ];
1124 nVec.resize( vsize, nullNode );
1125 loadedNodes.insert( nVec[ 0 ] = node );
1127 if ( theParam2ColumnMap.size() != hsize ) {
1128 RETURN_BAD_RESULT( "Wrong node positions on theBaseEdge" );
1131 // load nodes from e1
1133 std::map< double, const SMDS_MeshNode*> sortedNodes; // sort by param on edge
1134 nIt = sm1->GetNodes();
1135 while ( nIt->more() ) {
1139 const SMDS_EdgePosition* pos =
1140 dynamic_cast<const SMDS_EdgePosition*>( node->GetPosition().get() );
1144 sortedNodes.insert( std::make_pair( pos->GetUParameter(), node ));
1146 loadedNodes.insert( nVecf[ vsize - 1 ] = smVft->GetNodes()->next() );
1147 std::map< double, const SMDS_MeshNode*>::iterator u_n = sortedNodes.begin();
1148 int row = rev1 ? vsize - 1 : 0;
1149 int dRow = rev1 ? -1 : +1;
1150 for ( ; u_n != sortedNodes.end(); u_n++ ) {
1152 loadedNodes.insert( nVecf[ row ] = u_n->second );
1155 // try to load the rest nodes
1157 // get all faces from theFace
1158 TIDSortedElemSet allFaces, foundFaces;
1159 eIt = smFace->GetElements();
1160 while ( eIt->more() ) {
1161 const SMDS_MeshElement* e = eIt->next();
1162 if ( e->GetType() == SMDSAbs_Face )
1163 allFaces.insert( e );
1165 // Starting from 2 neighbour nodes on theBaseEdge, look for a face
1166 // the nodes belong to, and between the nodes of the found face,
1167 // look for a not loaded node considering this node to be the next
1168 // in a column of the starting second node. Repeat, starting
1169 // from nodes next to the previous starting nodes in their columns,
1170 // and so on while a face can be found. Then go the the next pair
1171 // of nodes on theBaseEdge.
1172 TParam2ColumnMap::iterator par_nVec_1 = theParam2ColumnMap.begin();
1173 TParam2ColumnMap::iterator par_nVec_2 = par_nVec_1;
1176 for ( par_nVec_2++; par_nVec_2 != theParam2ColumnMap.end(); par_nVec_1++, par_nVec_2++ ) {
1179 const SMDS_MeshNode* n1 = par_nVec_1->second[ row ];
1180 const SMDS_MeshNode* n2 = par_nVec_2->second[ row ];
1181 const SMDS_MeshElement* face = 0;
1182 bool lastColOnClosedFace = ( nVecf[ row ] == n2 );
1184 // look for a face by 2 nodes
1185 face = SMESH_MeshEditor::FindFaceInSet( n1, n2, allFaces, foundFaces );
1187 int nbFaceNodes = face->NbNodes();
1188 if ( face->IsQuadratic() )
1190 if ( nbFaceNodes>4 ) {
1191 RETURN_BAD_RESULT(" Too many nodes in a face: " << nbFaceNodes );
1193 // look for a not loaded node of the <face>
1195 const SMDS_MeshNode* n3 = 0; // a node defferent from n1 and n2
1196 for ( int i = 0; i < nbFaceNodes && !found; ++i ) {
1197 node = face->GetNode( i );
1198 found = loadedNodes.insert( node ).second;
1199 if ( !found && node != n1 && node != n2 )
1202 if ( lastColOnClosedFace && row + 1 < vsize ) {
1203 node = nVecf[ row + 1 ];
1204 found = ( face->GetNodeIndex( node ) >= 0 );
1207 if ( ++row > vsize - 1 ) {
1208 RETURN_BAD_RESULT( "Too many nodes in column "<< col <<": "<< row+1);
1210 par_nVec_2->second[ row ] = node;
1211 foundFaces.insert( face );
1213 if ( nbFaceNodes==4 ) {
1214 n1 = par_nVec_1->second[ row ];
1217 else if ( nbFaceNodes==3 && n3 == par_nVec_1->second[ row + 1 ] ) {
1221 RETURN_BAD_RESULT( "Not quad mesh, column "<< col );
1225 while ( face && n1 && n2 );
1227 if ( row < vsize - 1 ) {
1228 MESSAGE( "Too few nodes in column "<< col <<": "<< row+1);
1229 MESSAGE( "Base node 1: "<< par_nVec_1->second[0]);
1230 MESSAGE( "Base node 2: "<< par_nVec_2->second[0]);
1231 if ( n1 ) { MESSAGE( "Current node 1: "<< n1); }
1232 else { MESSAGE( "Current node 1: NULL"); }
1233 if ( n2 ) { MESSAGE( "Current node 2: "<< n2); }
1234 else { MESSAGE( "Current node 2: NULL"); }
1235 MESSAGE( "first base node: "<< theParam2ColumnMap.begin()->second[0]);
1236 MESSAGE( "last base node: "<< theParam2ColumnMap.rbegin()->second[0]);
1239 } // loop on columns
1244 //=======================================================================
1246 * \brief Return number of unique ancestors of the shape
1248 //=======================================================================
1250 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1251 const SMESH_Mesh& mesh,
1252 TopAbs_ShapeEnum ancestorType/*=TopAbs_SHAPE*/)
1254 TopTools_MapOfShape ancestors;
1255 TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1256 for ( ; ansIt.More(); ansIt.Next() ) {
1257 if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1258 ancestors.Add( ansIt.Value() );
1260 return ancestors.Extent();
1263 //=======================================================================
1265 * Check mesh without geometry for: if all elements on this shape are quadratic,
1266 * quadratic elements will be created.
1267 * Used then generated 3D mesh without geometry.
1269 //=======================================================================
1271 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1273 int NbAllEdgsAndFaces=0;
1274 int NbQuadFacesAndEdgs=0;
1275 int NbFacesAndEdges=0;
1276 //All faces and edges
1277 NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1279 //Quadratic faces and edges
1280 NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1282 //Linear faces and edges
1283 NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1285 if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1287 return SMESH_MesherHelper::QUADRATIC;
1289 else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1291 return SMESH_MesherHelper::LINEAR;
1294 //Mesh with both type of elements
1295 return SMESH_MesherHelper::COMP;
1298 //=======================================================================
1300 * \brief Return an alternative parameter for a node on seam
1302 //=======================================================================
1304 double SMESH_MesherHelper::GetOtherParam(const double param) const
1306 return fabs(param-myPar1) < fabs(param-myPar2) ? myPar2 : myPar1;
1309 //=======================================================================
1310 namespace { // Structures used by FixQuadraticElements()
1311 //=======================================================================
1313 #define __DMP__(txt) \
1315 #define MSG(txt) __DMP__(txt<<endl)
1316 #define MSGBEG(txt) __DMP__(txt)
1318 const double straightTol2 = 1e-33; // to detect straing links
1321 // ---------------------------------------
1323 * \brief Quadratic link knowing its faces
1325 struct QLink: public SMESH_TLink
1327 const SMDS_MeshNode* _mediumNode;
1328 mutable vector<const QFace* > _faces;
1329 mutable gp_Vec _nodeMove;
1330 mutable int _nbMoves;
1332 QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1333 SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1335 //if ( MediumPos() != SMDS_TOP_3DSPACE )
1336 _nodeMove = MediumPnt() - MiddlePnt();
1338 void SetContinuesFaces() const;
1339 const QFace* GetContinuesFace( const QFace* face ) const;
1340 bool OnBoundary() const;
1341 gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1342 gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1344 SMDS_TypeOfPosition MediumPos() const
1345 { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1346 SMDS_TypeOfPosition EndPos(bool isSecond) const
1347 { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1348 const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1349 { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1351 void Move(const gp_Vec& move, bool sum=false) const
1352 { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1353 gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1354 bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1355 bool IsStraight() const { return _nodeMove.SquareMagnitude() <= straightTol2; }
1357 bool operator<(const QLink& other) const {
1358 return (node1()->GetID() == other.node1()->GetID() ?
1359 node2()->GetID() < other.node2()->GetID() :
1360 node1()->GetID() < other.node1()->GetID());
1362 struct PtrComparator {
1363 bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1366 // ---------------------------------------------------------
1368 * \brief Link in the chain of links; it connects two faces
1372 const QLink* _qlink;
1373 mutable const QFace* _qfaces[2];
1375 TChainLink(const QLink* qlink=0):_qlink(qlink) {
1376 _qfaces[0] = _qfaces[1] = 0;
1378 void SetFace(const QFace* face) { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1380 bool IsBoundary() const { return !_qfaces[1]; }
1382 void RemoveFace( const QFace* face ) const
1383 { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) swap(_qfaces[0],_qfaces[1]); }
1385 const QFace* NextFace( const QFace* f ) const
1386 { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1388 const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1389 { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1391 bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1393 operator bool() const { return (_qlink); }
1395 const QLink* operator->() const { return _qlink; }
1397 gp_Vec Normal() const;
1399 // --------------------------------------------------------------------
1400 typedef list< TChainLink > TChain;
1401 typedef set < TChainLink > TLinkSet;
1402 typedef TLinkSet::iterator TLinkInSet;
1404 const int theFirstStep = 5;
1406 enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1407 // --------------------------------------------------------------------
1409 * \brief Face shared by two volumes and bound by QLinks
1411 struct QFace: public TIDSortedElemSet
1413 mutable const SMDS_MeshElement* _volumes[2];
1414 mutable vector< const QLink* > _sides;
1415 mutable bool _sideIsAdded[4]; // added in chain of links
1418 QFace( const vector< const QLink*>& links );
1420 void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1422 int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1424 void AddSelfToLinks() const {
1425 for ( int i = 0; i < _sides.size(); ++i )
1426 _sides[i]->_faces.push_back( this );
1428 int LinkIndex( const QLink* side ) const {
1429 for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1432 bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const;
1434 bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1436 int i = LinkIndex( link._qlink );
1437 if ( i < 0 ) return true;
1438 _sideIsAdded[i] = true;
1439 link.SetFace( this );
1440 // continue from opposite link
1441 return GetLinkChain( (i+2)%_sides.size(), chain, pos, error );
1443 bool IsBoundary() const { return !_volumes[1]; }
1445 bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1447 TLinkInSet GetBoundaryLink( const TLinkSet& links,
1448 const TChainLink& avoidLink,
1449 TLinkInSet * notBoundaryLink = 0,
1450 const SMDS_MeshNode* nodeToContain = 0,
1451 bool * isAdjacentUsed = 0) const;
1453 TLinkInSet GetLinkByNode( const TLinkSet& links,
1454 const TChainLink& avoidLink,
1455 const SMDS_MeshNode* nodeToContain) const;
1457 const SMDS_MeshNode* GetNodeInFace() const {
1458 for ( int iL = 0; iL < _sides.size(); ++iL )
1459 if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1463 gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1465 double MoveByBoundary( const TChainLink& theLink,
1466 const gp_Vec& theRefVec,
1467 const TLinkSet& theLinks,
1468 SMESH_MesherHelper* theFaceHelper=0,
1469 const double thePrevLen=0,
1470 const int theStep=theFirstStep,
1471 gp_Vec* theLinkNorm=0,
1472 double theSign=1.0) const;
1475 //================================================================================
1477 * \brief Dump QLink and QFace
1479 ostream& operator << (ostream& out, const QLink& l)
1481 out <<"QLink nodes: "
1482 << l.node1()->GetID() << " - "
1483 << l._mediumNode->GetID() << " - "
1484 << l.node2()->GetID() << endl;
1487 ostream& operator << (ostream& out, const QFace& f)
1489 out <<"QFace nodes: "/*<< &f << " "*/;
1490 for ( TIDSortedElemSet::iterator n = f.begin(); n != f.end(); ++n )
1491 out << (*n)->GetID() << " ";
1492 out << " \tvolumes: "
1493 << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1494 << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1495 out << " \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1499 //================================================================================
1501 * \brief Construct QFace from QLinks
1503 //================================================================================
1505 QFace::QFace( const vector< const QLink*>& links )
1507 _volumes[0] = _volumes[1] = 0;
1509 _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1510 _normal.SetCoord(0,0,0);
1511 for ( int i = 1; i < _sides.size(); ++i ) {
1512 const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1513 insert( l1->node1() ); insert( l1->node2() );
1515 gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1516 gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1517 if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1521 double normSqSize = _normal.SquareMagnitude();
1522 if ( normSqSize > numeric_limits<double>::min() )
1523 _normal /= sqrt( normSqSize );
1525 _normal.SetCoord(1e-33,0,0);
1527 //================================================================================
1529 * \brief Make up chain of links
1530 * \param iSide - link to add first
1531 * \param chain - chain to fill in
1532 * \param pos - postion of medium nodes the links should have
1533 * \param error - out, specifies what is wrong
1534 * \retval bool - false if valid chain can't be built; "valid" means that links
1535 * of the chain belongs to rectangles bounding hexahedrons
1537 //================================================================================
1539 bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1541 if ( iSide >= _sides.size() ) // wrong argument iSide
1543 if ( _sideIsAdded[ iSide ]) // already in chain
1546 if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1548 for ( int i = 0; i < _sides.size(); ++i ) {
1549 if ( !_sideIsAdded[i] && _sides[i] ) {
1550 _sideIsAdded[i]=true;
1551 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(_sides[i]));
1552 chLink->SetFace( this );
1553 if ( _sides[i]->MediumPos() >= pos )
1554 if ( const QFace* f = _sides[i]->GetContinuesFace( this ))
1555 f->GetLinkChain( *chLink, chain, pos, error );
1558 if ( error < ERR_TRI )
1562 _sideIsAdded[iSide] = true; // not to add this link to chain again
1563 const QLink* link = _sides[iSide];
1567 // add link into chain
1568 TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1569 chLink->SetFace( this );
1572 // propagate from rectangle to neighbour faces
1573 if ( link->MediumPos() >= pos ) {
1574 int nbLinkFaces = link->_faces.size();
1575 if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1576 // hexahedral mesh or boundary quadrangles - goto a continous face
1577 if ( const QFace* f = link->GetContinuesFace( this ))
1578 return f->GetLinkChain( *chLink, chain, pos, error );
1581 TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1582 for ( int i = 0; i < nbLinkFaces; ++i )
1583 if ( link->_faces[i] )
1584 link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1585 if ( error < ERR_PRISM )
1593 //================================================================================
1595 * \brief Return a boundary link of the triangle face
1596 * \param links - set of all links
1597 * \param avoidLink - link not to return
1598 * \param notBoundaryLink - out, neither the returned link nor avoidLink
1599 * \param nodeToContain - node the returned link must contain; if provided, search
1600 * also performed on adjacent faces
1601 * \param isAdjacentUsed - returns true if link is found in adjacent faces
1603 //================================================================================
1605 TLinkInSet QFace::GetBoundaryLink( const TLinkSet& links,
1606 const TChainLink& avoidLink,
1607 TLinkInSet * notBoundaryLink,
1608 const SMDS_MeshNode* nodeToContain,
1609 bool * isAdjacentUsed) const
1611 TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1613 typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1614 TFaceLinkList adjacentFaces;
1616 for ( int iL = 0; iL < _sides.size(); ++iL )
1618 if ( avoidLink._qlink == _sides[iL] )
1620 TLinkInSet link = links.find( _sides[iL] );
1621 if ( link == linksEnd ) continue;
1624 if ( link->IsBoundary() ) {
1625 if ( !nodeToContain ||
1626 (*link)->node1() == nodeToContain ||
1627 (*link)->node2() == nodeToContain )
1629 boundaryLink = link;
1630 if ( !notBoundaryLink ) break;
1633 else if ( notBoundaryLink ) {
1634 *notBoundaryLink = link;
1635 if ( boundaryLink != linksEnd ) break;
1638 if ( boundaryLink == linksEnd && nodeToContain ) // cellect adjacent faces
1639 if ( const QFace* adj = link->NextFace( this ))
1640 if ( adj->Contains( nodeToContain ))
1641 adjacentFaces.push_back( make_pair( adj, link ));
1644 if ( isAdjacentUsed ) *isAdjacentUsed = false;
1645 if ( boundaryLink == linksEnd && nodeToContain ) // check adjacent faces
1647 TFaceLinkList::iterator adj = adjacentFaces.begin();
1648 for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1649 boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second),
1650 0, nodeToContain, isAdjacentUsed);
1651 if ( isAdjacentUsed ) *isAdjacentUsed = true;
1653 return boundaryLink;
1655 //================================================================================
1657 * \brief Return a link ending at the given node but not avoidLink
1659 //================================================================================
1661 TLinkInSet QFace::GetLinkByNode( const TLinkSet& links,
1662 const TChainLink& avoidLink,
1663 const SMDS_MeshNode* nodeToContain) const
1665 for ( int i = 0; i < _sides.size(); ++i )
1666 if ( avoidLink._qlink != _sides[i] &&
1667 (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1668 return links.find( _sides[ i ]);
1672 //================================================================================
1674 * \brief Return normal to the i-th side pointing outside the face
1676 //================================================================================
1678 gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1680 gp_Vec norm, vecOut;
1681 // if ( uvHelper ) {
1682 // TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
1683 // const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
1684 // gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
1685 // gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
1686 // norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
1688 // const QLink* otherLink = _sides[(i + 1) % _sides.size()];
1689 // const SMDS_MeshNode* otherNode =
1690 // otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
1691 // gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
1692 // vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
1695 norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
1696 gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
1697 XYZ( _sides[0]->node2() ) +
1698 XYZ( _sides[1]->node1() )) / 3.;
1699 vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
1701 if ( norm * vecOut < 0 )
1703 double mag2 = norm.SquareMagnitude();
1704 if ( mag2 > numeric_limits<double>::min() )
1705 norm /= sqrt( mag2 );
1708 //================================================================================
1710 * \brief Move medium node of theLink according to its distance from boundary
1711 * \param theLink - link to fix
1712 * \param theRefVec - movement of boundary
1713 * \param theLinks - all adjacent links of continous triangles
1714 * \param theFaceHelper - helper is not used so far
1715 * \param thePrevLen - distance from the boundary
1716 * \param theStep - number of steps till movement propagation limit
1717 * \param theLinkNorm - out normal to theLink
1718 * \param theSign - 1 or -1 depending on movement of boundary
1719 * \retval double - distance from boundary to propagation limit or other boundary
1721 //================================================================================
1723 double QFace::MoveByBoundary( const TChainLink& theLink,
1724 const gp_Vec& theRefVec,
1725 const TLinkSet& theLinks,
1726 SMESH_MesherHelper* theFaceHelper,
1727 const double thePrevLen,
1729 gp_Vec* theLinkNorm,
1730 double theSign) const
1733 return thePrevLen; // propagation limit reached
1735 int iL; // index of theLink
1736 for ( iL = 0; iL < _sides.size(); ++iL )
1737 if ( theLink._qlink == _sides[ iL ])
1740 MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
1741 <<" thePrevLen " << thePrevLen);
1742 MSG(string(theStep,'.')<<" "<<*theLink._qlink);
1744 gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
1745 double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
1746 if ( theStep == theFirstStep )
1747 theSign = refProj < 0. ? -1. : 1.;
1748 else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
1749 return thePrevLen; // to propagate movement forward only, not in side dir or backward
1751 int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
1752 TLinkInSet link1 = theLinks.find( _sides[iL1] );
1753 TLinkInSet link2 = theLinks.find( _sides[iL2] );
1754 const QFace* f1 = link1->NextFace( this ); // adjacent faces
1755 const QFace* f2 = link2->NextFace( this );
1757 // propagate to adjacent faces till limit step or boundary
1758 double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
1759 double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
1760 gp_Vec linkDir1, linkDir2;
1764 len1 = f1->MoveByBoundary
1765 ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
1767 linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
1769 MSG( " --------------- EXCEPTION");
1775 len2 = f2->MoveByBoundary
1776 ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
1778 linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
1780 MSG( " --------------- EXCEPTION");
1785 if ( theStep != theFirstStep )
1787 // choose chain length by direction of propagation most codirected with theRefVec
1788 bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
1789 fullLen = choose1 ? len1 : len2;
1790 double r = thePrevLen / fullLen;
1792 gp_Vec move = linkNorm * refProj * ( 1 - r );
1793 theLink->Move( move, true );
1795 MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
1796 " by " << refProj * ( 1 - r ) << " following " <<
1797 (choose1 ? *link1->_qlink : *link2->_qlink));
1799 if ( theLinkNorm ) *theLinkNorm = linkNorm;
1804 //================================================================================
1806 * \brief Find pairs of continues faces
1808 //================================================================================
1810 void QLink::SetContinuesFaces() const
1812 // x0 x - QLink, [-|] - QFace, v - volume
1814 // | Between _faces of link x2 two vertical faces are continues
1815 // x1----x2-----x3 and two horizontal faces are continues. We set vertical faces
1816 // | to _faces[0] and _faces[1] and horizontal faces to
1817 // v2 | v3 _faces[2] and _faces[3] (or vise versa).
1820 if ( _faces.empty() )
1823 for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
1825 // look for a face bounding none of volumes bound by _faces[0]
1826 bool sameVol = false;
1827 int nbVol = _faces[iF]->NbVolumes();
1828 for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
1829 sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
1830 _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
1834 if ( iFaceCont > 0 ) // continues faces found, set one by the other
1836 if ( iFaceCont != 1 )
1837 swap( _faces[1], _faces[iFaceCont] );
1839 else if ( _faces.size() > 1 ) // not found, set NULL by the first face
1841 _faces.insert( ++_faces.begin(), 0 );
1844 //================================================================================
1846 * \brief Return a face continues to the given one
1848 //================================================================================
1850 const QFace* QLink::GetContinuesFace( const QFace* face ) const
1852 for ( int i = 0; i < _faces.size(); ++i ) {
1853 if ( _faces[i] == face ) {
1854 int iF = i < 2 ? 1-i : 5-i;
1855 return iF < _faces.size() ? _faces[iF] : 0;
1860 //================================================================================
1862 * \brief True if link is on mesh boundary
1864 //================================================================================
1866 bool QLink::OnBoundary() const
1868 for ( int i = 0; i < _faces.size(); ++i )
1869 if (_faces[i] && _faces[i]->IsBoundary()) return true;
1872 //================================================================================
1874 * \brief Return normal of link of the chain
1876 //================================================================================
1878 gp_Vec TChainLink::Normal() const {
1880 if (_qfaces[0]) norm = _qfaces[0]->_normal;
1881 if (_qfaces[1]) norm += _qfaces[1]->_normal;
1884 //================================================================================
1886 * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
1888 //================================================================================
1890 void fixPrism( TChain& allLinks )
1892 // separate boundary links from internal ones
1893 typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
1894 QLinkSet interLinks, bndLinks1, bndLink2;
1896 bool isCurved = false;
1897 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
1898 if ( (*lnk)->OnBoundary() )
1899 bndLinks1.insert( lnk->_qlink );
1901 interLinks.insert( lnk->_qlink );
1902 isCurved = isCurved || !(*lnk)->IsStraight();
1905 return; // no need to move
1907 QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
1909 while ( !interLinks.empty() && !curBndLinks->empty() )
1911 // propagate movement from boundary links to connected internal links
1912 QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
1913 for ( ; bnd != bndEnd; ++bnd )
1915 const QLink* bndLink = *bnd;
1916 for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
1918 const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
1919 if ( !face ) continue;
1920 // find and move internal link opposite to bndLink within the face
1921 int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
1922 const QLink* interLink = face->_sides[ interInd ];
1923 QLinkSet::iterator pInterLink = interLinks.find( interLink );
1924 if ( pInterLink == interLinks.end() ) continue; // not internal link
1925 interLink->Move( bndLink->_nodeMove );
1926 // treated internal links become new boundary ones
1927 interLinks. erase( pInterLink );
1928 newBndLinks->insert( interLink );
1931 curBndLinks->clear();
1932 swap( curBndLinks, newBndLinks );
1936 //================================================================================
1938 * \brief Fix links of continues triangles near curved boundary
1940 //================================================================================
1942 void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
1944 if ( allLinks.empty() ) return;
1946 TLinkSet linkSet( allLinks.begin(), allLinks.end());
1947 TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
1949 // move in 2d if we are on geom face
1950 // TopoDS_Face face;
1951 // TopLoc_Location loc;
1952 // SMESH_MesherHelper faceHelper( *helper.GetMesh());
1953 // while ( linkIt->IsBoundary()) ++linkIt;
1954 // if ( linkIt == linksEnd ) return;
1955 // if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
1956 // bool checkPos = true;
1957 // TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
1958 // if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
1959 // face = TopoDS::Face( f );
1960 // helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
1964 // faceHelper.SetSubShape( face );
1967 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
1969 if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
1971 // if ( !face.IsNull() ) {
1972 // const SMDS_MeshNode* inFaceNode =
1973 // faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
1974 // gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
1975 // gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
1976 // gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
1977 // gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
1978 // gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
1979 // linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
1982 linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
1988 //================================================================================
1990 * \brief Detect rectangular structure of links and build chains from them
1992 //================================================================================
1994 enum TSplitTriaResult {
1995 _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
1996 _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK };
1998 TSplitTriaResult splitTrianglesIntoChains( TChain & allLinks,
1999 vector< TChain> & resultChains,
2000 SMDS_TypeOfPosition pos )
2002 // put links in the set and evalute number of result chains by number of boundary links
2005 for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2006 linkSet.insert( *lnk );
2007 nbBndLinks += lnk->IsBoundary();
2009 resultChains.clear();
2010 resultChains.reserve( nbBndLinks / 2 );
2012 TLinkInSet linkIt, linksEnd = linkSet.end();
2014 // find a boundary link with corner node; corner node has position pos-2
2015 // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2017 SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2018 const SMDS_MeshNode* corner = 0;
2019 for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2020 if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2025 TLinkInSet startLink = linkIt;
2026 const SMDS_MeshNode* startCorner = corner;
2027 vector< TChain* > rowChains;
2030 while ( startLink != linksEnd) // loop on columns
2032 // We suppose we have a rectangular structure like shown here. We have found a
2033 // corner of the rectangle (startCorner) and a boundary link sharing
2034 // |/ |/ | the startCorner (startLink). We are going to loop on rows of the
2035 // --o---o---o structure making several chains at once. One chain (columnChain)
2036 // |\ | /| starts at startLink and continues upward (we look at the structure
2037 // \ | \ | / | from such point that startLink is on the bottom of the structure).
2038 // \| \|/ | While going upward we also fill horizontal chains (rowChains) we
2039 // --o---o---o encounter.
2041 // / | \ | \ | startCorner
2046 if ( resultChains.size() == nbBndLinks / 2 )
2048 resultChains.push_back( TChain() );
2049 TChain& columnChain = resultChains.back();
2051 TLinkInSet botLink = startLink; // current horizontal link to go up from
2052 corner = startCorner; // current corner the botLink ends at
2054 while ( botLink != linksEnd ) // loop on rows
2056 // add botLink to the columnChain
2057 columnChain.push_back( *botLink );
2059 const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2061 { // the column ends
2062 linkSet.erase( botLink );
2063 if ( iRow != rowChains.size() )
2064 return _FEW_ROWS; // different nb of rows in columns
2067 // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2068 // link ending at <corner> (sideLink); there are two cases:
2069 // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2070 // since midQuadLink is not at boundary while sideLink is.
2071 // 2) midQuadLink ends at <corner>
2073 TLinkInSet midQuadLink = linksEnd;
2074 TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2076 if ( isCase2 ) { // find midQuadLink among links of botTria
2077 midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2078 if ( midQuadLink->IsBoundary() )
2079 return _BAD_MIDQUAD;
2081 if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2082 return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2085 columnChain.push_back( *midQuadLink );
2086 if ( iRow >= rowChains.size() ) {
2088 return _MANY_ROWS; // different nb of rows in columns
2089 if ( resultChains.size() == nbBndLinks / 2 )
2091 resultChains.push_back( TChain() );
2092 rowChains.push_back( & resultChains.back() );
2094 rowChains[iRow]->push_back( *sideLink );
2095 rowChains[iRow]->push_back( *midQuadLink );
2097 const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2101 // prepare startCorner and startLink for the next column
2102 startCorner = startLink->NextNode( startCorner );
2104 startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2106 startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2107 // check if no more columns remains
2108 if ( startLink != linksEnd ) {
2109 const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2110 if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2111 startLink = linksEnd; // startLink bounds upTria or botTria
2112 else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2116 // find bottom link and corner for the next row
2117 corner = sideLink->NextNode( corner );
2118 // next bottom link ends at the new corner
2119 linkSet.erase( botLink );
2120 botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2121 if ( botLink == linksEnd || botLink == (isCase2 ? midQuadLink : sideLink))
2123 linkSet.erase( midQuadLink );
2124 linkSet.erase( sideLink );
2126 // make faces neighboring the found ones be boundary
2127 if ( startLink != linksEnd ) {
2128 const QFace* tria = isCase2 ? botTria : upTria;
2129 for ( int iL = 0; iL < 3; ++iL ) {
2130 linkIt = linkSet.find( tria->_sides[iL] );
2131 if ( linkIt != linksEnd )
2132 linkIt->RemoveFace( tria );
2135 if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2136 botLink->RemoveFace( upTria ); // make next botTria first in vector
2143 // In the linkSet, there must remain the last links of rowChains; add them
2144 if ( linkSet.size() != rowChains.size() )
2145 return _BAD_SET_SIZE;
2146 for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2147 // find the link (startLink) ending at startCorner
2149 for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2150 if ( (*startLink)->node1() == startCorner ) {
2151 corner = (*startLink)->node2(); break;
2153 else if ( (*startLink)->node2() == startCorner) {
2154 corner = (*startLink)->node1(); break;
2157 if ( startLink == linksEnd )
2159 rowChains[ iRow ]->push_back( *startLink );
2160 linkSet.erase( startLink );
2161 startCorner = corner;
2168 //=======================================================================
2170 * \brief Move medium nodes of faces and volumes to fix distorted elements
2171 * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2173 * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2175 //=======================================================================
2177 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2179 // apply algorithm to solids or geom faces
2180 // ----------------------------------------------
2181 if ( myShape.IsNull() ) {
2182 if ( !myMesh->HasShapeToMesh() ) return;
2183 SetSubShape( myMesh->GetShapeToMesh() );
2185 TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2186 for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2187 faces.Add( f.Current() );
2189 for ( TopExp_Explorer v(myShape,TopAbs_SOLID); v.More(); v.Next() ) {
2190 if ( myMesh->GetSubMesh( v.Current() )->IsEmpty() ) { // get faces of solid
2191 for ( TopExp_Explorer f( v.Current(), TopAbs_FACE); f.More(); f.Next() )
2192 faces.Add( f.Current() );
2194 else { // fix nodes in the solid and its faces
2195 SMESH_MesherHelper h(*myMesh);
2196 h.SetSubShape( v.Current() );
2197 h.FixQuadraticElements(false);
2200 // fix nodes on geom faces
2201 for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2202 SMESH_MesherHelper h(*myMesh);
2203 h.SetSubShape( fIt.Key() );
2204 h.FixQuadraticElements();
2209 // Find out type of elements and get iterator on them
2210 // ---------------------------------------------------
2212 SMDS_ElemIteratorPtr elemIt;
2213 SMDSAbs_ElementType elemType = SMDSAbs_All;
2215 SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2218 if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2219 elemIt = smDS->GetElements();
2220 if ( elemIt->more() ) {
2221 elemType = elemIt->next()->GetType();
2222 elemIt = smDS->GetElements();
2225 if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2228 // Fill in auxiliary data structures
2229 // ----------------------------------
2233 set< QLink >::iterator pLink;
2234 set< QFace >::iterator pFace;
2236 bool isCurved = false;
2237 bool hasRectFaces = false;
2238 set<int> nbElemNodeSet;
2240 if ( elemType == SMDSAbs_Volume )
2242 SMDS_VolumeTool volTool;
2243 while ( elemIt->more() ) // loop on volumes
2245 const SMDS_MeshElement* vol = elemIt->next();
2246 if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2248 for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2250 int nbN = volTool.NbFaceNodes( iF );
2251 nbElemNodeSet.insert( nbN );
2252 const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2253 vector< const QLink* > faceLinks( nbN/2 );
2254 for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2257 QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2258 pLink = links.insert( link ).first;
2259 faceLinks[ iN/2 ] = & *pLink;
2261 isCurved = !link.IsStraight();
2262 if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2263 return; // already fixed
2266 pFace = faces.insert( QFace( faceLinks )).first;
2267 if ( pFace->NbVolumes() == 0 )
2268 pFace->AddSelfToLinks();
2269 pFace->SetVolume( vol );
2270 hasRectFaces = hasRectFaces ||
2271 ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2272 volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2275 set< QLink >::iterator pLink = links.begin();
2276 for ( ; pLink != links.end(); ++pLink )
2277 pLink->SetContinuesFaces();
2281 while ( elemIt->more() ) // loop on faces
2283 const SMDS_MeshElement* face = elemIt->next();
2284 if ( !face->IsQuadratic() )
2286 nbElemNodeSet.insert( face->NbNodes() );
2287 int nbN = face->NbNodes()/2;
2288 vector< const QLink* > faceLinks( nbN );
2289 for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2292 QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2293 pLink = links.insert( link ).first;
2294 faceLinks[ iN ] = & *pLink;
2296 isCurved = !link.IsStraight();
2299 pFace = faces.insert( QFace( faceLinks )).first;
2300 pFace->AddSelfToLinks();
2301 hasRectFaces = ( hasRectFaces || nbN == 4 );
2305 return; // no curved edges of faces
2307 // Compute displacement of medium nodes
2308 // -------------------------------------
2310 // two loops on faces: the first is to treat boundary links, the second is for internal ones
2311 TopLoc_Location loc;
2312 // not treat boundary of volumic submesh
2313 int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2314 for ( ; isInside < 2; ++isInside ) {
2315 MSG( "--------------- LOOP " << isInside << " ------------------");
2316 SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2318 for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2319 if ( bool(isInside) == pFace->IsBoundary() )
2321 for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2324 // make chain of links connected via continues faces
2327 if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2329 if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2331 vector< TChain > chains;
2332 if ( error == ERR_OK ) { // chains contains continues rectangles
2334 chains[0].splice( chains[0].begin(), rawChain );
2336 else if ( error == ERR_TRI ) { // chains contains continues triangles
2337 TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2338 if ( res != _OK ) { // not rectangles split into triangles
2339 fixTriaNearBoundary( rawChain, *this );
2343 else if ( error == ERR_PRISM ) { // side faces of prisms
2344 fixPrism( rawChain );
2350 for ( int iC = 0; iC < chains.size(); ++iC )
2352 TChain& chain = chains[iC];
2353 if ( chain.empty() ) continue;
2354 if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2358 // mesure chain length and compute link position along the chain
2359 double chainLen = 0;
2360 vector< double > linkPos;
2361 MSGBEG( "Link medium nodes: ");
2362 TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2363 for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2364 MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2365 double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2366 while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2367 link1 = chain.erase( link1 );
2368 if ( link1 == chain.end() )
2370 len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2373 linkPos.push_back( chainLen );
2376 if ( linkPos.size() < 2 )
2379 gp_Vec move0 = chain.front()->_nodeMove;
2380 gp_Vec move1 = chain.back ()->_nodeMove;
2383 bool checkUV = true;
2385 // compute node displacement of end links in parametric space of face
2386 const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2387 TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2388 if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2389 face = TopoDS::Face( f );
2390 for ( int is1 = 0; is1 < 2; ++is1 ) { // move0 or move1
2391 TChainLink& link = is1 ? chain.back() : chain.front();
2392 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2393 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2394 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2395 gp_XY uvMove = uvm - GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2396 if ( is1 ) move1.SetCoord( uvMove.X(), uvMove.Y(), 0 );
2397 else move0.SetCoord( uvMove.X(), uvMove.Y(), 0 );
2399 if ( move0.SquareMagnitude() < straightTol2 &&
2400 move1.SquareMagnitude() < straightTol2 ) {
2402 continue; // straight - no need to move nodes of internal links
2407 if ( isInside || face.IsNull() )
2409 // compute node displacement of end links in their local coord systems
2411 TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2412 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2413 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2414 move0.Transform(trsf);
2417 TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2418 trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2419 gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2420 move1.Transform(trsf);
2423 // compute displacement of medium nodes
2424 link2 = chain.begin();
2427 for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2429 double r = linkPos[i] / chainLen;
2430 // displacement in local coord system
2431 gp_Vec move = (1. - r) * move0 + r * move1;
2432 if ( isInside || face.IsNull()) {
2433 // transform to global
2434 gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2435 gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2436 gp_Vec x = x01.Normalized() + x12.Normalized();
2437 trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2438 move.Transform(trsf);
2441 // compute 3D displacement by 2D one
2442 gp_XY oldUV = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2443 gp_XY newUV = oldUV + gp_XY( move.X(), move.Y() );
2444 gp_Pnt newPnt = BRep_Tool::Surface(face,loc)->Value( newUV.X(), newUV.Y());
2445 move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2447 if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2448 move.SquareMagnitude())
2450 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2451 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2452 MSG( "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2453 "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2454 "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2455 "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2459 (*link1)->Move( move );
2460 MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2461 << chain.front()->_mediumNode->GetID() <<"-"
2462 << chain.back ()->_mediumNode->GetID() <<
2463 " by " << move.Magnitude());
2465 } // loop on chains of links
2466 } // loop on 2 directions of propagation from quadrangle
2473 for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2474 if ( pLink->IsMoved() ) {
2475 //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2476 gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2477 GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());