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 bool allGroupsEmpty = true;
148 for ( size_t iG = 0; iG < srcGroups.size() && allGroupsEmpty; ++iG )
149 allGroupsEmpty = srcGroups[iG]->GetGroupDS()->IsEmpty();
150 if ( allGroupsEmpty )
151 return error("No faces in source groups");
153 SMESH_MesherHelper helper(theMesh);
154 helper.SetSubShape(theShape);
155 SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
157 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
158 const double faceTol = helper.MaxTolerance( geomFace );
159 const int shapeID = tgtMesh->ShapeToIndex( geomFace );
160 const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
162 Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
164 ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
165 gp_Pnt p; gp_Vec du, dv;
167 set<int> subShapeIDs;
168 subShapeIDs.insert( shapeID );
170 // get nodes on vertices
171 list < SMESH_TNodeXYZ > vertexNodes;
172 list < SMESH_TNodeXYZ >::iterator vNIt;
173 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
174 for ( ; exp.More(); exp.Next() )
176 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
177 if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
179 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
182 _gen->Compute(theMesh,v,/*anUpward=*/true);
183 n = SMESH_Algo::VertexNode( v, tgtMesh );
184 if ( !n ) return false; // very strange
186 vertexNodes.push_back( SMESH_TNodeXYZ( n ));
189 // to count now many times a link between nodes encounters
190 map<TLink, int> linkCount;
191 map<TLink, int>::iterator link2Nb;
192 double minLinkLen2 = Precision::Infinite();
194 // =========================
195 // Import faces from groups
196 // =========================
198 StdMeshers_Import_1D::TNodeNodeMap* n2n;
199 StdMeshers_Import_1D::TElemElemMap* e2e;
200 vector<const SMDS_MeshNode*> newNodes;
201 for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
203 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
205 const int meshID = srcGroup->GetMesh()->GetPersistentId();
206 const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
207 if ( !srcMesh ) continue;
208 StdMeshers_Import_1D::getMaps( srcMesh, &theMesh, n2n, e2e );
210 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
211 SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
212 gp_XY uv( Precision::Infinite(), Precision::Infinite() );
213 while ( srcElems->more() ) // loop on group contents
215 const SMDS_MeshElement* face = srcElems->next();
216 // find or create nodes of a new face
217 newNodes.resize( face->NbNodes() );
219 int nbCreatedNodes = 0;
220 SMDS_MeshElement::iterator node = face->begin_nodes();
221 for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
223 StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
226 if ( !subShapeIDs.count( n2nIt->second->getshapeId() ))
231 // find an existing vertex node
232 for ( vNIt = vertexNodes.begin(); vNIt != vertexNodes.end(); ++vNIt)
233 if ( vNIt->SquareDistance( *node ) < 10 * faceTol * faceTol)
235 (*n2nIt).second = vNIt->_node;
236 vertexNodes.erase( vNIt );
240 if ( !n2nIt->second )
242 // find out if node lies on theShape
243 tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
244 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, 10 * faceTol, /*force=*/true ))
246 SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
247 n2nIt->second = newNode;
248 tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
252 if ( !(newNodes[i] = n2nIt->second ))
255 if ( !newNodes.back() )
256 continue; // not all nodes of the face lie on theShape
258 // try to find already created face
259 SMDS_MeshElement * newFace = 0;
260 if ( nbCreatedNodes == 0 &&
261 tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
262 continue; // repeated face in source groups already created
264 // check future face orientation
271 uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
272 surface->D1( uv.X(),uv.Y(), p, du,dv );
273 geomNorm = reverse ? dv^du : du^dv;
275 while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < face->NbCornerNodes());
277 int iNext = helper.WrapIndex( iNode+1, face->NbCornerNodes() );
278 int iPrev = helper.WrapIndex( iNode-1, face->NbCornerNodes() );
280 SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
281 SMESH_TNodeXYZ curNode ( newNodes[iNode] );
282 SMESH_TNodeXYZ nextNode( newNodes[iNext] );
283 gp_Vec n1n0( prevNode - curNode);
284 gp_Vec n1n2( nextNode - curNode );
285 gp_Vec meshNorm = n1n2 ^ n1n0;
287 if ( geomNorm * meshNorm < 0 )
288 std::reverse( newNodes.begin(), newNodes.end() );
292 switch ( newNodes.size() )
295 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
298 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
301 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
302 newNodes[3], newNodes[4], newNodes[5]);
305 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
306 newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
310 tgtMesh->SetMeshElementOnShape( newFace, shapeID );
311 e2e->insert( make_pair( face, newFace ));
314 int nbNodes = face->NbCornerNodes();
315 const SMDS_MeshNode* medium = 0;
316 for ( int i = 0; i < nbNodes; ++i )
318 const SMDS_MeshNode* n1 = newNodes[i];
319 const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbNodes ];
320 if ( newFace->IsQuadratic() )
321 medium = newNodes[i+nbNodes];
322 link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
324 if ( link2Nb->second == 1 )
326 // measure link length
327 double len2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
328 if ( len2 < minLinkLen2 )
333 helper.GetMeshDS()->RemoveNode(tmpNode);
336 // ==========================================================
337 // Put nodes on geom edges and create edges on them;
338 // check if the whole geom face is covered by imported faces
339 // ==========================================================
341 vector< TopoDS_Edge > edges;
342 for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
343 if ( subShapeIDs.insert( tgtMesh->ShapeToIndex( exp.Current() )).second )
344 edges.push_back( TopoDS::Edge( exp.Current() ));
346 // use large tolerance for projection of nodes to edges because of
347 // BLSURF mesher specifics (issue 0020918, Study2.hdf)
348 const double projTol = 1e-3 * sqrt( minLinkLen2 );
350 bool isFaceMeshed = false;
351 SMESHDS_SubMesh* tgtFaceSM = tgtMesh->MeshElements( theShape );
354 // the imported mesh is valid if all external links (encountered once)
356 subShapeIDs.erase( shapeID ); // to contain edges and vertices only
358 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
360 const TLink& link = (*link2Nb).first;
361 int nbFaces = link2Nb->second;
364 // check if a not shared link lies on face boundary
365 bool nodesOnBoundary = true;
366 list< TopoDS_Shape > bndShapes;
367 for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
369 const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
370 if ( !subShapeIDs.count( n->getshapeId() ))
372 for ( size_t iE = 0; iE < edges.size(); ++iE )
373 if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
375 BRep_Tool::Range(edges[iE],f,l);
376 if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
377 // duplicated node on vertex
378 return error("Source elements overlap one another");
379 tgtFaceSM->RemoveNode( n, /*isNodeDeleted=*/false );
380 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)n, edges[iE], u );
383 nodesOnBoundary = subShapeIDs.count( n->getshapeId());
385 if ( nodesOnBoundary )
387 TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
388 if ( s.ShapeType() == TopAbs_VERTEX )
389 bndShapes.push_front( s ); // vertex first
391 bndShapes.push_back( s ); // edges last
394 if ( !nodesOnBoundary )
395 break; // error: free internal link
396 if ( bndShapes.front().ShapeType() == TopAbs_EDGE &&
397 bndShapes.front() != bndShapes.back() )
398 break; // error: link nodes on different geom edges
400 // find geom edge the link is on
401 if ( bndShapes.back().ShapeType() != TopAbs_EDGE )
403 // find geom edge by two vertices
404 TopoDS_Shape geomEdge;
405 PShapeIteratorPtr edgeIt = helper.GetAncestors( bndShapes.back(), theMesh, TopAbs_EDGE );
406 while ( edgeIt->more() )
408 geomEdge = *(edgeIt->next());
409 if ( !helper.IsSubShape( bndShapes.front(), geomEdge ))
412 if ( geomEdge.IsNull() )
413 break; // vertices belong to different edges -> error: free internal link
414 bndShapes.push_back( geomEdge );
417 // create an edge if not yet exists
419 newNodes[0] = link.node1(), newNodes[1] = link.node2();
420 const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
421 if ( edge ) continue;
423 if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
426 newNodes.push_back( link._medium );
427 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], newNodes[2] );
429 TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
430 helper.CheckNodeU( geomEdge, link._medium, u, 10*faceTol, /*force=*/true );
431 tgtFaceSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
432 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
436 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
441 tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
443 else if ( nbFaces > 2 )
445 return error( "Non-manifold source mesh");
448 isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
451 // check that source faces do not overlap:
452 // there must be only two edges sharing each vertex and bound to sub-edges of theShape
453 SMESH_MeshEditor editor( &theMesh );
454 set<int>::iterator subID = subShapeIDs.begin();
455 for ( ; subID != subShapeIDs.end(); ++subID )
457 const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
458 if ( s.ShapeType() != TopAbs_VERTEX ) continue;
459 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
460 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
462 while ( eIt->more() )
464 const SMDS_MeshElement* edge = eIt->next();
465 int sId = editor.FindShape( edge );
466 nbEdges += subShapeIDs.count( sId );
469 return false; // weird
471 return error( "Source elements overlap one another");
476 return error( "Source elements don't cover totally the geometrical face" );
478 if ( helper.HasSeam() )
480 // links on seam edges are shared by two faces, so no edges were created on them
481 // by the previous detection of 2D mesh boundary
482 for ( size_t iE = 0; iE < edges.size(); ++iE )
484 if ( !helper.IsRealSeam( edges[iE] )) continue;
485 const TopoDS_Edge& seamEdge = edges[iE];
486 // to find nodes lying on the seamEdge we check nodes of mesh faces sharing a node on one
487 // of its vertices; after finding another node on seamEdge we continue the same way
488 // until finding all nodes.
489 TopoDS_Vertex seamVertex = helper.IthVertex( 0, seamEdge );
490 const SMDS_MeshNode* vertNode = SMESH_Algo::VertexNode( seamVertex, tgtMesh );
491 set< const SMDS_MeshNode* > checkedNodes; checkedNodes.insert( vertNode );
492 set< const SMDS_MeshElement* > checkedFaces;
493 // as a face can have more than one node on the seamEdge, there is a difficulty in selecting
494 // one of those nodes to treat next; so we simply find all nodes on the seamEdge and
495 // then sort them by U on edge
496 typedef list< pair< double, const SMDS_MeshNode* > > TUNodeList;
497 TUNodeList nodesOnSeam;
498 double u = helper.GetNodeU( seamEdge, vertNode );
499 nodesOnSeam.push_back( make_pair( u, vertNode ));
500 TUNodeList::iterator u2nIt = nodesOnSeam.begin();
501 for ( ; u2nIt != nodesOnSeam.end(); ++u2nIt )
503 const SMDS_MeshNode* startNode = (*u2nIt).second;
504 SMDS_ElemIteratorPtr faceIt = startNode->GetInverseElementIterator( SMDSAbs_Face );
505 while ( faceIt->more() )
507 const SMDS_MeshElement* face = faceIt->next();
508 if ( !checkedFaces.insert( face ).second ) continue;
509 for ( int i = 0, nbNodes = face->NbCornerNodes(); i < nbNodes; ++i )
511 const SMDS_MeshNode* n = face->GetNode( i );
512 if ( n == startNode || !checkedNodes.insert( n ).second ) continue;
513 if ( helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true ))
514 nodesOnSeam.push_back( make_pair( u, n ));
518 // sort the found nodes by U on the seamEdge; most probably they are in a good order,
519 // so we can use the hint to spead-up map filling
520 map< double, const SMDS_MeshNode* > u2nodeMap;
521 for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
522 u2nodeMap.insert( u2nodeMap.end(), *u2nIt );
526 SMESH_MesherHelper seamHelper( theMesh );
527 seamHelper.SetSubShape( edges[ iE ]);
528 seamHelper.SetElementsOnShape( true );
530 if ( (*checkedFaces.begin())->IsQuadratic() )
531 for ( set< const SMDS_MeshElement* >::iterator fIt = checkedFaces.begin();
532 fIt != checkedFaces.end(); ++fIt )
533 seamHelper.AddTLinks( static_cast<const SMDS_MeshFace*>( *fIt ));
535 map< double, const SMDS_MeshNode* >::iterator n1, n2, u2nEnd = u2nodeMap.end();
536 for ( n2 = u2nodeMap.begin(), n1 = n2++; n2 != u2nEnd; ++n1, ++n2 )
538 const SMDS_MeshNode* node1 = n1->second;
539 const SMDS_MeshNode* node2 = n2->second;
540 seamHelper.AddEdge( node1, node2 );
541 if ( node2->getshapeId() == helper.GetSubShapeID() )
543 tgtFaceSM->RemoveNode( node2, /*isNodeDeleted=*/false );
544 tgtMesh->SetNodeOnEdge( const_cast<SMDS_MeshNode*>( node2 ), seamEdge, n2->first );
548 } // loop on edges to find seam ones
549 } // if ( helper.HasSeam() )
551 // notify sub-meshes of edges on computation
552 for ( size_t iE = 0; iE < edges.size(); ++iE )
554 SMESH_subMesh * sm = theMesh.GetSubMesh( edges[iE] );
555 if ( BRep_Tool::Degenerated( edges[iE] ))
556 sm->SetIsAlwaysComputed( true );
557 sm->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
558 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
559 return error(SMESH_Comment("Failed to create segments on the edge ")
560 << tgtMesh->ShapeToIndex( edges[iE ]));
567 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
568 for ( size_t i = 0; i < srcMeshes.size(); ++i )
569 StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
574 //=============================================================================
576 * \brief Set needed event listeners and create a submesh for a copied mesh
578 * This method is called only if a submesh has HYP_OK algo_state.
580 //=============================================================================
582 void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
586 const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
587 SMESH_Mesh* tgtMesh = subMesh->GetFather();
588 Hypothesis_Status aStatus;
589 CheckHypothesis( *tgtMesh, tgtShape, aStatus );
591 StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
593 void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
595 SetEventListener(subMesh);
598 //=============================================================================
600 * Predict nb of mesh entities created by Compute()
602 //=============================================================================
604 bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh & theMesh,
605 const TopoDS_Shape & theShape,
606 MapShapeNbElems& aResMap)
608 if ( !_sourceHyp ) return false;
610 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
611 if ( srcGroups.empty() )
612 return error("Invalid source groups");
614 vector<int> aVec(SMDSEntity_Last,0);
616 bool toCopyMesh, toCopyGroups;
617 _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
618 if ( toCopyMesh ) // the whole mesh is copied
620 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
621 for ( unsigned i = 0; i < srcMeshes.size(); ++i )
623 SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
624 if ( !sm || aResMap.count( sm )) continue; // already counted
625 const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
626 for (int i = 0; i < SMDSEntity_Last; i++)
627 aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
632 // std-like iterator used to get coordinates of nodes of mesh element
633 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
635 SMESH_MesherHelper helper(theMesh);
636 helper.SetSubShape(theShape);
638 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
639 const double faceTol = helper.MaxTolerance( geomFace );
641 // take into account nodes on vertices
642 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
643 for ( ; exp.More(); exp.Next() )
644 theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
646 // to count now many times a link between nodes encounters,
647 // negative nb additionally means that a link is quadratic
648 map<SMESH_TLink, int> linkCount;
649 map<SMESH_TLink, int>::iterator link2Nb;
651 // count faces and nodes imported from groups
652 set<const SMDS_MeshNode* > allNodes;
654 for ( int iG = 0; iG < srcGroups.size(); ++iG )
656 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
657 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
658 SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
659 while ( srcElems->more() ) // loop on group contents
661 const SMDS_MeshElement* face = srcElems->next();
662 // find out if face is located on geomEdge by projecting
663 // a gravity center of face to geomFace
665 gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
666 tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
667 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, 10 * faceTol, /*force=*/true ))
669 ++aVec[ face->GetEntityType() ];
672 int nbConers = face->NbCornerNodes();
673 for ( int i = 0; i < face->NbNodes(); ++i )
675 const SMDS_MeshNode* n1 = face->GetNode(i);
676 allNodes.insert( n1 );
679 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
680 link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
681 if ( (*link2Nb).second )
682 link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
684 link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
689 helper.GetMeshDS()->RemoveNode(tmpNode);
692 int nbNodes = allNodes.size();
695 // count nodes and edges on geom edges
698 for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
700 TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
701 SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
702 vector<int>& edgeVec = aResMap[sm];
703 if ( edgeVec.empty() )
705 edgeVec.resize(SMDSEntity_Last,0);
706 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
708 const SMESH_TLink& link = (*link2Nb).first;
709 int nbFacesOfLink = Abs( link2Nb->second );
710 bool eraseLink = ( nbFacesOfLink != 1 );
711 if ( nbFacesOfLink == 1 )
713 if ( helper.CheckNodeU( geomEdge, link.node1(), u, 10*faceTol, /*force=*/true )&&
714 helper.CheckNodeU( geomEdge, link.node2(), u, 10*faceTol, /*force=*/true ))
716 bool isQuadratic = ( link2Nb->second < 0 );
717 ++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
718 ++edgeVec[ SMDSEntity_Node ];
724 linkCount.erase(link2Nb++);
728 if ( edgeVec[ SMDSEntity_Node] > 0 )
729 --edgeVec[ SMDSEntity_Node ]; // for one node on vertex
731 else if ( !helper.IsSeamShape( geomEdge ) ||
732 geomEdge.Orientation() == TopAbs_FORWARD )
734 nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
738 aVec[SMDSEntity_Node] = nbNodes;
741 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
742 aResMap.insert(make_pair(sm,aVec));