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>
62 double getMinElemSize2( const SMESHDS_GroupBase* srcGroup )
64 double minSize2 = 1e100;
65 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
66 while ( srcElems->more() ) // loop on group contents
68 const SMDS_MeshElement* face = srcElems->next();
69 int nbN = face->NbCornerNodes();
71 SMESH_TNodeXYZ prevN( face->GetNode( nbN-1 ));
72 for ( int i = 0; i < nbN; ++i )
74 SMESH_TNodeXYZ n( face->GetNode( i ) );
75 double size2 = ( n - prevN ).SquareModulus();
76 minSize2 = std::min( minSize2, size2 );
84 //=============================================================================
86 * Creates StdMeshers_Import_1D2D
88 //=============================================================================
90 StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, int studyId, SMESH_Gen * gen)
91 :SMESH_2D_Algo(hypId, studyId, gen), _sourceHyp(0)
93 MESSAGE("StdMeshers_Import_1D2D::StdMeshers_Import_1D2D");
94 _name = "Import_1D2D";
95 _shapeType = (1 << TopAbs_FACE);
97 _compatibleHypothesis.push_back("ImportSource2D");
98 _requireDescretBoundary = false;
101 //=============================================================================
103 * Check presence of a hypothesis
105 //=============================================================================
107 bool StdMeshers_Import_1D2D::CheckHypothesis
109 const TopoDS_Shape& aShape,
110 SMESH_Hypothesis::Hypothesis_Status& aStatus)
114 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
115 if ( hyps.size() == 0 )
117 aStatus = SMESH_Hypothesis::HYP_MISSING;
118 return false; // can't work with no hypothesis
121 if ( hyps.size() > 1 )
123 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
127 const SMESHDS_Hypothesis *theHyp = hyps.front();
129 string hypName = theHyp->GetName();
131 if (hypName == _compatibleHypothesis.front())
133 _sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
134 aStatus = SMESH_Hypothesis::HYP_OK;
138 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
145 * \brief OrientedLink additionally storing a medium node
147 struct TLink : public SMESH_OrientedLink
149 const SMDS_MeshNode* _medium;
150 TLink( const SMDS_MeshNode* n1,
151 const SMDS_MeshNode* n2,
152 const SMDS_MeshNode* medium=0)
153 : SMESH_OrientedLink( n1,n2 ), _medium( medium ) {}
157 //=============================================================================
159 * Import elements from the other mesh
161 //=============================================================================
163 bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
165 if ( !_sourceHyp ) return false;
167 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
168 if ( srcGroups.empty() )
169 return error("Invalid source groups");
171 bool allGroupsEmpty = true;
172 for ( size_t iG = 0; iG < srcGroups.size() && allGroupsEmpty; ++iG )
173 allGroupsEmpty = srcGroups[iG]->GetGroupDS()->IsEmpty();
174 if ( allGroupsEmpty )
175 return error("No faces in source groups");
177 SMESH_MesherHelper helper(theMesh);
178 helper.SetSubShape(theShape);
179 SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
181 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
182 const double faceTol = helper.MaxTolerance( geomFace );
183 const int shapeID = tgtMesh->ShapeToIndex( geomFace );
184 const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
186 Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
188 ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
189 gp_Pnt p; gp_Vec du, dv;
191 set<int> subShapeIDs;
192 subShapeIDs.insert( shapeID );
194 // get nodes on vertices
195 list < SMESH_TNodeXYZ > vertexNodes;
196 list < SMESH_TNodeXYZ >::iterator vNIt;
197 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
198 for ( ; exp.More(); exp.Next() )
200 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
201 if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
203 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
206 _gen->Compute(theMesh,v,/*anUpward=*/true);
207 n = SMESH_Algo::VertexNode( v, tgtMesh );
208 if ( !n ) return false; // very strange
210 vertexNodes.push_back( SMESH_TNodeXYZ( n ));
213 // to count now many times a link between nodes encounters
214 map<TLink, int> linkCount;
215 map<TLink, int>::iterator link2Nb;
216 double minGroupTol = Precision::Infinite();
218 // =========================
219 // Import faces from groups
220 // =========================
222 StdMeshers_Import_1D::TNodeNodeMap* n2n;
223 StdMeshers_Import_1D::TElemElemMap* e2e;
224 vector<const SMDS_MeshNode*> newNodes;
225 for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
227 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
229 const int meshID = srcGroup->GetMesh()->GetPersistentId();
230 const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
231 if ( !srcMesh ) continue;
232 StdMeshers_Import_1D::getMaps( srcMesh, &theMesh, n2n, e2e );
234 const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
235 minGroupTol = std::min( groupTol, minGroupTol );
237 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
238 SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
239 gp_XY uv( Precision::Infinite(), Precision::Infinite() );
240 while ( srcElems->more() ) // loop on group contents
242 const SMDS_MeshElement* face = srcElems->next();
243 // find or create nodes of a new face
244 newNodes.resize( face->NbNodes() );
246 int nbCreatedNodes = 0;
247 SMDS_MeshElement::iterator node = face->begin_nodes();
248 for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
250 StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
253 if ( !subShapeIDs.count( n2nIt->second->getshapeId() ))
258 // find an existing vertex node
259 for ( vNIt = vertexNodes.begin(); vNIt != vertexNodes.end(); ++vNIt)
260 if ( vNIt->SquareDistance( *node ) < groupTol * groupTol)
262 (*n2nIt).second = vNIt->_node;
263 vertexNodes.erase( vNIt );
267 if ( !n2nIt->second )
269 // find out if node lies on theShape
270 tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
271 uv.SetCoord( Precision::Infinite(), Precision::Infinite() );
272 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
274 SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
275 n2nIt->second = newNode;
276 tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
280 if ( !(newNodes[i] = n2nIt->second ))
283 if ( !newNodes.back() )
284 continue; // not all nodes of the face lie on theShape
286 // try to find already created face
287 SMDS_MeshElement * newFace = 0;
288 if ( nbCreatedNodes == 0 &&
289 tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
290 continue; // repeated face in source groups already created
292 // check future face orientation
299 uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
300 surface->D1( uv.X(),uv.Y(), p, du,dv );
301 geomNorm = reverse ? dv^du : du^dv;
303 while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < face->NbCornerNodes());
305 int iNext = helper.WrapIndex( iNode+1, face->NbCornerNodes() );
306 int iPrev = helper.WrapIndex( iNode-1, face->NbCornerNodes() );
308 SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
309 SMESH_TNodeXYZ curNode ( newNodes[iNode] );
310 SMESH_TNodeXYZ nextNode( newNodes[iNext] );
311 gp_Vec n1n0( prevNode - curNode);
312 gp_Vec n1n2( nextNode - curNode );
313 gp_Vec meshNorm = n1n2 ^ n1n0;
315 if ( geomNorm * meshNorm < 0 )
316 std::reverse( newNodes.begin(), newNodes.end() );
320 switch ( newNodes.size() )
323 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
326 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
329 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
330 newNodes[3], newNodes[4], newNodes[5]);
333 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
334 newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
338 tgtMesh->SetMeshElementOnShape( newFace, shapeID );
339 e2e->insert( make_pair( face, newFace ));
342 int nbNodes = face->NbCornerNodes();
343 const SMDS_MeshNode* medium = 0;
344 for ( int i = 0; i < nbNodes; ++i )
346 const SMDS_MeshNode* n1 = newNodes[i];
347 const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbNodes ];
348 if ( newFace->IsQuadratic() )
349 medium = newNodes[i+nbNodes];
350 link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
352 // if ( link2Nb->second == 1 )
354 // // measure link length
355 // double len2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
356 // if ( len2 < minGroupTol )
357 // minGroupTol = len2;
361 helper.GetMeshDS()->RemoveNode(tmpNode);
364 // ==========================================================
365 // Put nodes on geom edges and create edges on them;
366 // check if the whole geom face is covered by imported faces
367 // ==========================================================
369 vector< TopoDS_Edge > edges;
370 for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
371 if ( subShapeIDs.insert( tgtMesh->ShapeToIndex( exp.Current() )).second )
372 edges.push_back( TopoDS::Edge( exp.Current() ));
374 // use large tolerance for projection of nodes to edges because of
375 // BLSURF mesher specifics (issue 0020918, Study2.hdf)
376 const double projTol = minGroupTol;
378 bool isFaceMeshed = false;
379 SMESHDS_SubMesh* tgtFaceSM = tgtMesh->MeshElements( theShape );
382 // the imported mesh is valid if all external links (encountered once)
384 subShapeIDs.erase( shapeID ); // to contain edges and vertices only
386 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
388 const TLink& link = (*link2Nb).first;
389 int nbFaces = link2Nb->second;
392 // check if a not shared link lies on face boundary
393 bool nodesOnBoundary = true;
394 list< TopoDS_Shape > bndShapes;
395 for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
397 const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
398 if ( !subShapeIDs.count( n->getshapeId() ))
400 for ( size_t iE = 0; iE < edges.size(); ++iE )
401 if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
403 BRep_Tool::Range(edges[iE],f,l);
404 if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
405 // duplicated node on vertex
406 return error("Source elements overlap one another");
407 tgtFaceSM->RemoveNode( n, /*isNodeDeleted=*/false );
408 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)n, edges[iE], u );
411 nodesOnBoundary = subShapeIDs.count( n->getshapeId());
413 if ( nodesOnBoundary )
415 TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
416 if ( s.ShapeType() == TopAbs_VERTEX )
417 bndShapes.push_front( s ); // vertex first
419 bndShapes.push_back( s ); // edges last
422 if ( !nodesOnBoundary )
423 break; // error: free internal link
424 if ( bndShapes.front().ShapeType() == TopAbs_EDGE &&
425 bndShapes.front() != bndShapes.back() )
426 break; // error: link nodes on different geom edges
428 // find geom edge the link is on
429 if ( bndShapes.back().ShapeType() != TopAbs_EDGE )
431 // find geom edge by two vertices
432 TopoDS_Shape geomEdge = helper.GetCommonAncestor( bndShapes.back(),
434 theMesh, TopAbs_EDGE );
435 if ( geomEdge.IsNull() )
436 break; // vertices belong to different edges -> error: free internal link
437 bndShapes.push_back( geomEdge );
440 // create an edge if not yet exists
442 newNodes[0] = link.node1(), newNodes[1] = link.node2();
443 const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
444 if ( edge ) continue;
446 if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
449 newNodes.push_back( link._medium );
450 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], newNodes[2] );
452 TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
453 helper.CheckNodeU( geomEdge, link._medium, u, projTol, /*force=*/true );
454 tgtFaceSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
455 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
459 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
464 tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
466 else if ( nbFaces > 2 )
468 return error( "Non-manifold source mesh");
471 isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
474 // check that source faces do not overlap:
475 // there must be only two edges sharing each vertex and bound to sub-edges of theShape
476 SMESH_MeshEditor editor( &theMesh );
477 set<int>::iterator subID = subShapeIDs.begin();
478 for ( ; subID != subShapeIDs.end(); ++subID )
480 const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
481 if ( s.ShapeType() != TopAbs_VERTEX ) continue;
482 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
483 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
485 while ( eIt->more() )
487 const SMDS_MeshElement* edge = eIt->next();
488 int sId = editor.FindShape( edge );
489 nbEdges += subShapeIDs.count( sId );
492 return false; // weird
494 return error( "Source elements overlap one another");
499 return error( "Source elements don't cover totally the geometrical face" );
501 if ( helper.HasSeam() )
503 // links on seam edges are shared by two faces, so no edges were created on them
504 // by the previous detection of 2D mesh boundary
505 for ( size_t iE = 0; iE < edges.size(); ++iE )
507 if ( !helper.IsRealSeam( edges[iE] )) continue;
508 const TopoDS_Edge& seamEdge = edges[iE];
509 // to find nodes lying on the seamEdge we check nodes of mesh faces sharing a node on one
510 // of its vertices; after finding another node on seamEdge we continue the same way
511 // until finding all nodes.
512 TopoDS_Vertex seamVertex = helper.IthVertex( 0, seamEdge );
513 const SMDS_MeshNode* vertNode = SMESH_Algo::VertexNode( seamVertex, tgtMesh );
514 set< const SMDS_MeshNode* > checkedNodes; checkedNodes.insert( vertNode );
515 set< const SMDS_MeshElement* > checkedFaces;
516 // as a face can have more than one node on the seamEdge, there is a difficulty in selecting
517 // one of those nodes to treat next; so we simply find all nodes on the seamEdge and
518 // then sort them by U on edge
519 typedef list< pair< double, const SMDS_MeshNode* > > TUNodeList;
520 TUNodeList nodesOnSeam;
521 double u = helper.GetNodeU( seamEdge, vertNode );
522 nodesOnSeam.push_back( make_pair( u, vertNode ));
523 TUNodeList::iterator u2nIt = nodesOnSeam.begin();
524 for ( ; u2nIt != nodesOnSeam.end(); ++u2nIt )
526 const SMDS_MeshNode* startNode = (*u2nIt).second;
527 SMDS_ElemIteratorPtr faceIt = startNode->GetInverseElementIterator( SMDSAbs_Face );
528 while ( faceIt->more() )
530 const SMDS_MeshElement* face = faceIt->next();
531 if ( !checkedFaces.insert( face ).second ) continue;
532 for ( int i = 0, nbNodes = face->NbCornerNodes(); i < nbNodes; ++i )
534 const SMDS_MeshNode* n = face->GetNode( i );
535 if ( n == startNode || !checkedNodes.insert( n ).second ) continue;
536 if ( helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true ))
537 nodesOnSeam.push_back( make_pair( u, n ));
541 // sort the found nodes by U on the seamEdge; most probably they are in a good order,
542 // so we can use the hint to spead-up map filling
543 map< double, const SMDS_MeshNode* > u2nodeMap;
544 for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
545 u2nodeMap.insert( u2nodeMap.end(), *u2nIt );
549 SMESH_MesherHelper seamHelper( theMesh );
550 seamHelper.SetSubShape( edges[ iE ]);
551 seamHelper.SetElementsOnShape( true );
553 if ( (*checkedFaces.begin())->IsQuadratic() )
554 for ( set< const SMDS_MeshElement* >::iterator fIt = checkedFaces.begin();
555 fIt != checkedFaces.end(); ++fIt )
556 seamHelper.AddTLinks( static_cast<const SMDS_MeshFace*>( *fIt ));
558 map< double, const SMDS_MeshNode* >::iterator n1, n2, u2nEnd = u2nodeMap.end();
559 for ( n2 = u2nodeMap.begin(), n1 = n2++; n2 != u2nEnd; ++n1, ++n2 )
561 const SMDS_MeshNode* node1 = n1->second;
562 const SMDS_MeshNode* node2 = n2->second;
563 seamHelper.AddEdge( node1, node2 );
564 if ( node2->getshapeId() == helper.GetSubShapeID() )
566 tgtFaceSM->RemoveNode( node2, /*isNodeDeleted=*/false );
567 tgtMesh->SetNodeOnEdge( const_cast<SMDS_MeshNode*>( node2 ), seamEdge, n2->first );
571 } // loop on edges to find seam ones
572 } // if ( helper.HasSeam() )
574 // notify sub-meshes of edges on computation
575 for ( size_t iE = 0; iE < edges.size(); ++iE )
577 SMESH_subMesh * sm = theMesh.GetSubMesh( edges[iE] );
578 if ( BRep_Tool::Degenerated( edges[iE] ))
579 sm->SetIsAlwaysComputed( true );
580 sm->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
581 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
582 return error(SMESH_Comment("Failed to create segments on the edge ")
583 << tgtMesh->ShapeToIndex( edges[iE ]));
590 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
591 for ( size_t i = 0; i < srcMeshes.size(); ++i )
592 StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
597 //=============================================================================
599 * \brief Set needed event listeners and create a submesh for a copied mesh
601 * This method is called only if a submesh has HYP_OK algo_state.
603 //=============================================================================
605 void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
609 const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
610 SMESH_Mesh* tgtMesh = subMesh->GetFather();
611 Hypothesis_Status aStatus;
612 CheckHypothesis( *tgtMesh, tgtShape, aStatus );
614 StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
616 void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
618 SetEventListener(subMesh);
621 //=============================================================================
623 * Predict nb of mesh entities created by Compute()
625 //=============================================================================
627 bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh & theMesh,
628 const TopoDS_Shape & theShape,
629 MapShapeNbElems& aResMap)
631 if ( !_sourceHyp ) return false;
633 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
634 if ( srcGroups.empty() )
635 return error("Invalid source groups");
637 vector<int> aVec(SMDSEntity_Last,0);
639 bool toCopyMesh, toCopyGroups;
640 _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
641 if ( toCopyMesh ) // the whole mesh is copied
643 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
644 for ( unsigned i = 0; i < srcMeshes.size(); ++i )
646 SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
647 if ( !sm || aResMap.count( sm )) continue; // already counted
648 const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
649 for (int i = 0; i < SMDSEntity_Last; i++)
650 aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
655 // std-like iterator used to get coordinates of nodes of mesh element
656 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
658 SMESH_MesherHelper helper(theMesh);
659 helper.SetSubShape(theShape);
661 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
663 // take into account nodes on vertices
664 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
665 for ( ; exp.More(); exp.Next() )
666 theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
668 // to count now many times a link between nodes encounters,
669 // negative nb additionally means that a link is quadratic
670 map<SMESH_TLink, int> linkCount;
671 map<SMESH_TLink, int>::iterator link2Nb;
673 // count faces and nodes imported from groups
674 set<const SMDS_MeshNode* > allNodes;
676 double minGroupTol = 1e100;
677 for ( int iG = 0; iG < srcGroups.size(); ++iG )
679 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
680 const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
681 minGroupTol = std::min( groupTol, minGroupTol );
682 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
683 SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
684 while ( srcElems->more() ) // loop on group contents
686 const SMDS_MeshElement* face = srcElems->next();
687 // find out if face is located on geomEdge by projecting
688 // a gravity center of face to geomFace
690 gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
691 tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
692 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
694 ++aVec[ face->GetEntityType() ];
697 int nbConers = face->NbCornerNodes();
698 for ( int i = 0; i < face->NbNodes(); ++i )
700 const SMDS_MeshNode* n1 = face->GetNode(i);
701 allNodes.insert( n1 );
704 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
705 link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
706 if ( (*link2Nb).second )
707 link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
709 link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
714 helper.GetMeshDS()->RemoveNode(tmpNode);
717 int nbNodes = allNodes.size();
720 // count nodes and edges on geom edges
723 for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
725 TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
726 SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
727 vector<int>& edgeVec = aResMap[sm];
728 if ( edgeVec.empty() )
730 edgeVec.resize(SMDSEntity_Last,0);
731 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
733 const SMESH_TLink& link = (*link2Nb).first;
734 int nbFacesOfLink = Abs( link2Nb->second );
735 bool eraseLink = ( nbFacesOfLink != 1 );
736 if ( nbFacesOfLink == 1 )
738 if ( helper.CheckNodeU( geomEdge, link.node1(), u, minGroupTol, /*force=*/true )&&
739 helper.CheckNodeU( geomEdge, link.node2(), u, minGroupTol, /*force=*/true ))
741 bool isQuadratic = ( link2Nb->second < 0 );
742 ++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
743 ++edgeVec[ SMDSEntity_Node ];
749 linkCount.erase(link2Nb++);
753 if ( edgeVec[ SMDSEntity_Node] > 0 )
754 --edgeVec[ SMDSEntity_Node ]; // for one node on vertex
756 else if ( !helper.IsSeamShape( geomEdge ) ||
757 geomEdge.Orientation() == TopAbs_FORWARD )
759 nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
763 aVec[SMDSEntity_Node] = nbNodes;
766 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
767 aResMap.insert(make_pair(sm,aVec));