1 // Copyright (C) 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // SMESH SMESH : implementaion of SMESH idl descriptions
24 // File : StdMeshers_Import_1D2D.cxx
27 #include "StdMeshers_Import_1D2D.hxx"
29 #include "StdMeshers_Import_1D.hxx"
30 #include "StdMeshers_ImportSource.hxx"
32 #include "SMDS_MeshElement.hxx"
33 #include "SMDS_MeshNode.hxx"
34 #include "SMESHDS_Group.hxx"
35 #include "SMESHDS_Mesh.hxx"
36 #include "SMESH_Comment.hxx"
37 #include "SMESH_Gen.hxx"
38 #include "SMESH_Group.hxx"
39 #include "SMESH_Mesh.hxx"
40 #include "SMESH_MesherHelper.hxx"
41 #include "SMESH_subMesh.hxx"
43 #include "Utils_SALOME_Exception.hxx"
44 #include "utilities.h"
46 #include <BRep_Builder.hxx>
47 #include <BRep_Tool.hxx>
49 #include <TopExp_Explorer.hxx>
51 #include <TopoDS_Compound.hxx>
52 #include <TopoDS_Edge.hxx>
53 #include <TopoDS_Vertex.hxx>
54 #include <Precision.hxx>
60 //=============================================================================
62 * Creates StdMeshers_Import_1D2D
64 //=============================================================================
66 StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, int studyId, SMESH_Gen * gen)
67 :SMESH_2D_Algo(hypId, studyId, gen), _sourceHyp(0)
69 MESSAGE("StdMeshers_Import_1D2D::StdMeshers_Import_1D2D");
70 _name = "Import_1D2D";
71 _shapeType = (1 << TopAbs_FACE);
73 _compatibleHypothesis.push_back("ImportSource2D");
74 _requireDescretBoundary = false;
77 //=============================================================================
79 * Check presence of a hypothesis
81 //=============================================================================
83 bool StdMeshers_Import_1D2D::CheckHypothesis
85 const TopoDS_Shape& aShape,
86 SMESH_Hypothesis::Hypothesis_Status& aStatus)
90 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
91 if ( hyps.size() == 0 )
93 aStatus = SMESH_Hypothesis::HYP_MISSING;
94 return false; // can't work with no hypothesis
97 if ( hyps.size() > 1 )
99 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
103 const SMESHDS_Hypothesis *theHyp = hyps.front();
105 string hypName = theHyp->GetName();
107 if (hypName == _compatibleHypothesis.front())
109 _sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
110 aStatus = SMESH_Hypothesis::HYP_OK;
114 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
121 * \brief OrientedLink additionally storing a medium node
123 struct TLink : public SMESH_OrientedLink
125 const SMDS_MeshNode* _medium;
126 TLink( const SMDS_MeshNode* n1,
127 const SMDS_MeshNode* n2,
128 const SMDS_MeshNode* medium=0)
129 : SMESH_OrientedLink( n1,n2 ), _medium( medium ) {}
133 //=============================================================================
135 * Import elements from the other mesh
137 //=============================================================================
139 bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
141 if ( !_sourceHyp ) return false;
143 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
144 if ( srcGroups.empty() )
145 return error("Invalid source groups");
147 SMESH_MesherHelper helper(theMesh);
148 helper.SetSubShape(theShape);
149 SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
151 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
152 const double faceTol = helper.MaxTolerance( geomFace );
153 const int shapeID = tgtMesh->ShapeToIndex( geomFace );
154 const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
156 Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
158 ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
159 gp_Pnt p; gp_Vec du, dv;
161 set<int> subShapeIDs;
162 subShapeIDs.insert( shapeID );
164 // get nodes on vertices
165 list < SMESH_TNodeXYZ > vertexNodes;
166 list < SMESH_TNodeXYZ >::iterator vNIt;
167 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
168 for ( ; exp.More(); exp.Next() )
170 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
171 if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
173 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
176 _gen->Compute(theMesh,v,/*anUpward=*/true);
177 n = SMESH_Algo::VertexNode( v, tgtMesh );
178 if ( !n ) return false; // very strange
180 vertexNodes.push_back( SMESH_TNodeXYZ( n ));
183 // to count now many times a link between nodes encounters
184 map<TLink, int> linkCount;
185 map<TLink, int>::iterator link2Nb;
186 double minLinkLen2 = Precision::Infinite();
188 // =========================
189 // Import faces from groups
190 // =========================
192 StdMeshers_Import_1D::TNodeNodeMap* n2n;
193 StdMeshers_Import_1D::TElemElemMap* e2e;
194 vector<const SMDS_MeshNode*> newNodes;
195 for ( int iG = 0; iG < srcGroups.size(); ++iG )
197 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
199 const int meshID = srcGroup->GetMesh()->GetPersistentId();
200 const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
201 if ( !srcMesh ) continue;
202 StdMeshers_Import_1D::getMaps( srcMesh, &theMesh, n2n, e2e );
204 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
205 SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
206 gp_XY uv( Precision::Infinite(), Precision::Infinite() );
207 while ( srcElems->more() ) // loop on group contents
209 const SMDS_MeshElement* face = srcElems->next();
210 // find or create nodes of a new face
211 newNodes.resize( face->NbNodes() );
213 int nbCreatedNodes = 0;
214 SMDS_MeshElement::iterator node = face->begin_nodes();
215 for ( unsigned i = 0; i < newNodes.size(); ++i, ++node )
217 TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
220 if ( !subShapeIDs.count( n2nIt->second->getshapeId() ))
225 // find an existing vertex node
226 for ( vNIt = vertexNodes.begin(); vNIt != vertexNodes.end(); ++vNIt)
227 if ( vNIt->SquareDistance( *node ) < 10 * faceTol * faceTol)
229 (*n2nIt).second = vNIt->_node;
230 vertexNodes.erase( vNIt );
234 if ( !n2nIt->second )
236 // find out if node lies on theShape
237 tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
238 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, 10 * faceTol, /*force=*/true ))
240 SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
241 n2nIt->second = newNode;
242 tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
246 if ( !(newNodes[i] = n2nIt->second ))
249 if ( !newNodes.back() )
250 continue; // not all nodes of the face lie on theShape
252 // try to find already created face
253 SMDS_MeshElement * newFace = 0;
254 if ( nbCreatedNodes == 0 &&
255 tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
256 continue; // repeated face in source groups already created
258 // check future face orientation
265 uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
266 surface->D1( uv.X(),uv.Y(), p, du,dv );
267 geomNorm = reverse ? dv^du : du^dv;
269 while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < face->NbCornerNodes());
271 int iNext = helper.WrapIndex( iNode+1, face->NbCornerNodes() );
272 int iPrev = helper.WrapIndex( iNode-1, face->NbCornerNodes() );
274 SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
275 SMESH_TNodeXYZ curNode ( newNodes[iNode] );
276 SMESH_TNodeXYZ nextNode( newNodes[iNext] );
277 gp_Vec n1n0( prevNode - curNode);
278 gp_Vec n1n2( nextNode - curNode );
279 gp_Vec meshNorm = n1n2 ^ n1n0;
281 if ( geomNorm * meshNorm < 0 )
282 std::reverse( newNodes.begin(), newNodes.end() );
286 switch ( newNodes.size() )
289 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
292 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
295 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
296 newNodes[3], newNodes[4], newNodes[5]);
299 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
300 newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
304 tgtMesh->SetMeshElementOnShape( newFace, shapeID );
305 e2e->insert( make_pair( face, newFace ));
308 int nbNodes = face->NbCornerNodes();
309 const SMDS_MeshNode* medium = 0;
310 for ( int i = 0; i < nbNodes; ++i )
312 const SMDS_MeshNode* n1 = newNodes[i];
313 const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbNodes ];
314 if ( newFace->IsQuadratic() )
315 medium = newNodes[i+nbNodes];
316 link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
318 if ( link2Nb->second == 1 )
320 // measure link length
321 double len2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
322 if ( len2 < minLinkLen2 )
327 helper.GetMeshDS()->RemoveNode(tmpNode);
330 // ==========================================================
331 // Put nodes on geom edges and create edges on them;
332 // check if the whole geom face is covered by imported faces
333 // ==========================================================
335 vector< TopoDS_Edge > edges;
336 for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
337 if ( subShapeIDs.insert( tgtMesh->ShapeToIndex( exp.Current() )).second )
338 edges.push_back( TopoDS::Edge( exp.Current() ));
340 // use large tolerance for projection of nodes to edges because of
341 // BLSURF mesher specifics (issue 0020918, Study2.hdf)
342 const double projTol = 1e-3 * sqrt( minLinkLen2 );
344 bool isFaceMeshed = false;
345 if ( SMESHDS_SubMesh* tgtSM = tgtMesh->MeshElements( theShape ))
347 // the imported mesh is valid if all external links (encountered once)
349 subShapeIDs.erase( shapeID ); // to contain edges and vertices only
351 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
353 const TLink& link = (*link2Nb).first;
354 int nbFaces = link2Nb->second;
357 // check if a not shared link lie on face boundary
358 bool nodesOnBoundary = true;
359 list< TopoDS_Shape > bndShapes;
360 for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
362 const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
363 if ( !subShapeIDs.count( n->getshapeId() ))
365 for ( unsigned iE = 0; iE < edges.size(); ++iE )
366 if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
368 BRep_Tool::Range(edges[iE],f,l);
369 if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
370 // duplicated node on vertex
371 return error("Source elements overlap one another");
372 tgtSM->RemoveNode( n, /*isNodeDeleted=*/false );
373 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)n, edges[iE], u );
376 nodesOnBoundary = subShapeIDs.count( n->getshapeId());
378 if ( nodesOnBoundary )
380 TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
381 if ( s.ShapeType() == TopAbs_VERTEX )
382 bndShapes.push_front( s ); // vertex first
384 bndShapes.push_back( s ); // edges last
387 if ( !nodesOnBoundary )
388 break; // error: free internal link
389 if ( bndShapes.front().ShapeType() == TopAbs_EDGE &&
390 bndShapes.front() != bndShapes.back() )
391 break; // error: link nodes on different geom edges
393 // find geom edge the link is on
394 if ( bndShapes.back().ShapeType() != TopAbs_EDGE )
396 // find geom edge by two vertices
397 TopoDS_Shape geomEdge;
398 PShapeIteratorPtr edgeIt = helper.GetAncestors( bndShapes.back(), theMesh, TopAbs_EDGE );
399 while ( edgeIt->more() )
401 geomEdge = *(edgeIt->next());
402 if ( !helper.IsSubShape( bndShapes.front(), geomEdge ))
405 if ( geomEdge.IsNull() )
406 break; // vertices belong to different edges -> error: free internal link
407 bndShapes.push_back( geomEdge );
410 // create an edge if not yet exists
412 newNodes[0] = link.node1(), newNodes[1] = link.node2();
413 const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
414 if ( edge ) continue;
416 if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
419 newNodes.push_back( link._medium );
420 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], newNodes[2] );
422 TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
423 helper.CheckNodeU( geomEdge, link._medium, u, 10*faceTol, /*force=*/true );
424 tgtSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
425 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
429 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
434 tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
436 else if ( nbFaces > 2 )
438 return error( "Non-manifold source mesh");
441 isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
444 // check that source faces do not overlap:
445 // there must be only two edges sharing each vertex and bound to sub-edges of theShape
446 SMESH_MeshEditor editor( &theMesh );
447 set<int>::iterator subID = subShapeIDs.begin();
448 for ( ; subID != subShapeIDs.end(); ++subID )
450 const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
451 if ( s.ShapeType() != TopAbs_VERTEX ) continue;
452 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
453 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
455 while ( eIt->more() )
457 const SMDS_MeshElement* edge = eIt->next();
458 int sId = editor.FindShape( edge );
459 nbEdges += subShapeIDs.count( sId );
462 return false; // weird
464 return error( "Source elements overlap one another");
469 return error( "Source elements don't cover totally the geometrical face" );
471 // notify sub-meshes of edges on computation
472 for ( unsigned iE = 0; iE < edges.size(); ++iE )
473 theMesh.GetSubMesh( edges[iE] )->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
479 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
480 for ( unsigned i = 0; i < srcMeshes.size(); ++i )
481 StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
486 //=============================================================================
488 * \brief Set needed event listeners and create a submesh for a copied mesh
490 * This method is called only if a submesh has HYP_OK algo_state.
492 //=============================================================================
494 void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
498 const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
499 SMESH_Mesh* tgtMesh = subMesh->GetFather();
500 Hypothesis_Status aStatus;
501 CheckHypothesis( *tgtMesh, tgtShape, aStatus );
503 StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
505 void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
507 SetEventListener(subMesh);
510 //=============================================================================
512 * Predict nb of mesh entities created by Compute()
514 //=============================================================================
516 bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh & theMesh,
517 const TopoDS_Shape & theShape,
518 MapShapeNbElems& aResMap)
520 if ( !_sourceHyp ) return false;
522 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
523 if ( srcGroups.empty() )
524 return error("Invalid source groups");
526 vector<int> aVec(SMDSEntity_Last,0);
528 bool toCopyMesh, toCopyGroups;
529 _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
530 if ( toCopyMesh ) // the whole mesh is copied
532 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
533 for ( unsigned i = 0; i < srcMeshes.size(); ++i )
535 SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
536 if ( !sm || aResMap.count( sm )) continue; // already counted
537 const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
538 for (int i = 0; i < SMDSEntity_Last; i++)
539 aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
544 // std-like iterator used to get coordinates of nodes of mesh element
545 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
547 SMESH_MesherHelper helper(theMesh);
548 helper.SetSubShape(theShape);
550 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
551 const double faceTol = helper.MaxTolerance( geomFace );
553 // take into account nodes on vertices
554 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
555 for ( ; exp.More(); exp.Next() )
556 theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
558 // to count now many times a link between nodes encounters,
559 // negative nb additionally means that a link is quadratic
560 map<SMESH_TLink, int> linkCount;
561 map<SMESH_TLink, int>::iterator link2Nb;
563 // count faces and nodes imported from groups
564 set<const SMDS_MeshNode* > allNodes;
566 for ( int iG = 0; iG < srcGroups.size(); ++iG )
568 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
569 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
570 SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
571 while ( srcElems->more() ) // loop on group contents
573 const SMDS_MeshElement* face = srcElems->next();
574 // find out if face is located on geomEdge by projecting
575 // a gravity center of face to geomFace
577 gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
578 tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
579 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, 10 * faceTol, /*force=*/true ))
581 ++aVec[ face->GetEntityType() ];
584 int nbConers = face->NbCornerNodes();
585 for ( int i = 0; i < face->NbNodes(); ++i )
587 const SMDS_MeshNode* n1 = face->GetNode(i);
588 allNodes.insert( n1 );
591 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
592 link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
593 if ( (*link2Nb).second )
594 link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
596 link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
601 helper.GetMeshDS()->RemoveNode(tmpNode);
604 int nbNodes = allNodes.size();
607 // count nodes and edges on geom edges
610 for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
612 TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
613 SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
614 vector<int>& edgeVec = aResMap[sm];
615 if ( edgeVec.empty() )
617 edgeVec.resize(SMDSEntity_Last,0);
618 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
620 const SMESH_TLink& link = (*link2Nb).first;
621 int nbFacesOfLink = Abs( link2Nb->second );
622 bool eraseLink = ( nbFacesOfLink != 1 );
623 if ( nbFacesOfLink == 1 )
625 if ( helper.CheckNodeU( geomEdge, link.node1(), u, 10*faceTol, /*force=*/true )&&
626 helper.CheckNodeU( geomEdge, link.node2(), u, 10*faceTol, /*force=*/true ))
628 bool isQuadratic = ( link2Nb->second < 0 );
629 ++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
630 ++edgeVec[ SMDSEntity_Node ];
636 linkCount.erase(link2Nb++);
640 if ( edgeVec[ SMDSEntity_Node] > 0 )
641 --edgeVec[ SMDSEntity_Node ]; // for one node on vertex
643 else if ( !helper.IsSeamShape( geomEdge ) ||
644 geomEdge.Orientation() == TopAbs_FORWARD )
646 nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
650 aVec[SMDSEntity_Node] = nbNodes;
653 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
654 aResMap.insert(make_pair(sm,aVec));