1 // Copyright (C) 2007-2019 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, or (at your option) any later version.
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 : implementation 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_ControlsDef.hxx"
38 #include "SMESH_Gen.hxx"
39 #include "SMESH_Group.hxx"
40 #include "SMESH_Mesh.hxx"
41 #include "SMESH_MesherHelper.hxx"
42 #include "SMESH_OctreeNode.hxx"
43 #include "SMESH_MeshEditor.hxx"
44 #include "SMESH_subMesh.hxx"
46 #include "Utils_SALOME_Exception.hxx"
47 #include "utilities.h"
49 #include <BRepBndLib.hxx>
50 #include <BRepClass_FaceClassifier.hxx>
51 #include <BRepTools.hxx>
52 #include <BRepTopAdaptor_FClass2d.hxx>
53 #include <BRep_Builder.hxx>
54 #include <BRep_Tool.hxx>
55 #include <Bnd_B2d.hxx>
56 #include <Bnd_Box.hxx>
57 #include <GeomAPI_ProjectPointOnSurf.hxx>
58 #include <GeomAdaptor_Surface.hxx>
59 #include <Precision.hxx>
61 #include <TopExp_Explorer.hxx>
63 #include <TopoDS_Compound.hxx>
64 #include <TopoDS_Edge.hxx>
65 #include <TopoDS_Vertex.hxx>
73 double getMinElemSize2( const SMESHDS_GroupBase* srcGroup )
75 double minSize2 = 1e100;
76 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
77 while ( srcElems->more() ) // loop on group contents
79 const SMDS_MeshElement* face = srcElems->next();
80 int nbN = face->NbCornerNodes();
82 SMESH_TNodeXYZ prevN( face->GetNode( nbN-1 ));
83 for ( int i = 0; i < nbN; ++i )
85 SMESH_TNodeXYZ n( face->GetNode( i ) );
86 double size2 = ( n - prevN ).SquareModulus();
87 minSize2 = std::min( minSize2, size2 );
95 //=============================================================================
97 * Creates StdMeshers_Import_1D2D
99 //=============================================================================
101 StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, SMESH_Gen * gen)
102 :SMESH_2D_Algo(hypId, gen), _sourceHyp(0)
104 _name = "Import_1D2D";
105 _shapeType = (1 << TopAbs_FACE);
107 _compatibleHypothesis.push_back("ImportSource2D");
108 _requireDiscreteBoundary = false;
109 _supportSubmeshes = true;
112 //=============================================================================
114 * Check presence of a hypothesis
116 //=============================================================================
118 bool StdMeshers_Import_1D2D::CheckHypothesis
120 const TopoDS_Shape& aShape,
121 SMESH_Hypothesis::Hypothesis_Status& aStatus)
125 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
126 if ( hyps.size() == 0 )
128 aStatus = SMESH_Hypothesis::HYP_MISSING;
129 return false; // can't work with no hypothesis
132 if ( hyps.size() > 1 )
134 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
138 const SMESHDS_Hypothesis *theHyp = hyps.front();
140 string hypName = theHyp->GetName();
142 if (hypName == _compatibleHypothesis.front())
144 _sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
145 aStatus = SMESH_Hypothesis::HYP_OK;
149 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
156 * \brief OrientedLink additionally storing a medium node
158 struct TLink : public SMESH_OrientedLink
160 const SMDS_MeshNode* _medium;
161 TLink( const SMDS_MeshNode* n1,
162 const SMDS_MeshNode* n2,
163 const SMDS_MeshNode* medium=0)
164 : SMESH_OrientedLink( n1,n2 ), _medium( medium ) {}
168 //=============================================================================
170 * Import elements from the other mesh
172 //=============================================================================
174 bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
176 if ( !_sourceHyp ) return false;
178 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups(/*loaded=*/true);
179 if ( srcGroups.empty() )
180 return error("Invalid source groups");
182 bool allGroupsEmpty = true;
183 for ( size_t iG = 0; iG < srcGroups.size() && allGroupsEmpty; ++iG )
184 allGroupsEmpty = srcGroups[iG]->GetGroupDS()->IsEmpty();
185 if ( allGroupsEmpty )
186 return error("No faces in source groups");
188 SMESH_MesherHelper helper(theMesh);
189 helper.SetSubShape(theShape);
190 SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
192 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
193 const double faceTol = helper.MaxTolerance( geomFace );
194 const int shapeID = tgtMesh->ShapeToIndex( geomFace );
195 const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
198 Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
200 ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace ) == TopAbs_REVERSED );
201 gp_Pnt p; gp_Vec du, dv;
203 // BRepClass_FaceClassifier is most time consuming, so minimize its usage
204 const double clsfTol = 10 * BRep_Tool::Tolerance( geomFace );
205 BRepTopAdaptor_FClass2d classifier( geomFace, clsfTol ); //Brimless_FaceClassifier classifier;
209 Standard_Real umin,umax,vmin,vmax;
210 BRepTools::UVBounds(geomFace,umin,umax,vmin,vmax);
211 gp_XY pmin( umin,vmin ), pmax( umax,vmax );
212 bndBox2d.Add( pmin );
213 bndBox2d.Add( pmax );
214 if ( helper.HasSeam() )
216 const int i = helper.GetPeriodicIndex();
217 pmin.SetCoord( i, helper.GetOtherParam( pmin.Coord( i )));
218 pmax.SetCoord( i, helper.GetOtherParam( pmax.Coord( i )));
219 bndBox2d.Add( pmin );
220 bndBox2d.Add( pmax );
222 bndBox2d.Enlarge( 1e-2 * Sqrt( bndBox2d.SquareExtent() ));
224 BRepBndLib::Add( geomFace, bndBox3d );
225 bndBox3d.Enlarge( 1e-2 * sqrt( bndBox3d.SquareExtent() ));
228 set<int> subShapeIDs;
229 subShapeIDs.insert( shapeID );
231 // nodes already existing on sub-shapes of the FACE
232 TIDSortedNodeSet existingNodes;
234 // get/make nodes on vertices and add them to existingNodes
235 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
236 for ( ; exp.More(); exp.Next() )
238 const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
239 if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
241 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
244 _gen->Compute(theMesh,v,/*anUpward=*/true);
245 n = SMESH_Algo::VertexNode( v, tgtMesh );
246 if ( !n ) return false; // very strange
248 existingNodes.insert( n );
251 // get EDGEs and their ids and get existing nodes on EDGEs
252 vector< TopoDS_Edge > edges;
253 for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
255 const TopoDS_Edge & edge = TopoDS::Edge( exp.Current() );
256 if ( !SMESH_Algo::isDegenerated( edge ))
257 if ( subShapeIDs.insert( tgtMesh->ShapeToIndex( edge )).second )
259 edges.push_back( edge );
260 if ( SMESHDS_SubMesh* eSM = tgtMesh->MeshElements( edge ))
262 typedef SMDS_StdIterator< const SMDS_MeshNode*, SMDS_NodeIteratorPtr > iterator;
263 existingNodes.insert( iterator( eSM->GetNodes() ), iterator() );
267 // octree to find existing nodes
268 SMESH_OctreeNode existingNodeOcTr( existingNodes );
269 std::map<double, const SMDS_MeshNode*> dist2foundNodes;
271 // to count now many times a link between nodes encounters
272 map<TLink, int> linkCount;
273 map<TLink, int>::iterator link2Nb;
274 double minGroupTol = Precision::Infinite();
276 SMESH::Controls::ElementsOnShape onEdgeClassifier;
277 if ( helper.HasSeam() )
279 TopoDS_Compound edgesCompound;
280 BRep_Builder builder;
281 builder.MakeCompound( edgesCompound );
282 for ( size_t iE = 0; iE < edges.size(); ++iE )
283 builder.Add( edgesCompound, edges[ iE ]);
284 onEdgeClassifier.SetShape( edgesCompound, SMDSAbs_Node );
287 // =========================
288 // Import faces from groups
289 // =========================
291 StdMeshers_Import_1D::TNodeNodeMap* n2n;
292 StdMeshers_Import_1D::TElemElemMap* e2e;
293 StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt;
294 pair< StdMeshers_Import_1D::TNodeNodeMap::iterator, bool > it_isnew;
295 vector<TopAbs_State> nodeState;
296 vector<const SMDS_MeshNode*> newNodes; // of a face
297 set <const SMDS_MeshNode*> bndNodes; // nodes classified ON
298 vector<bool> isNodeIn; // nodes classified IN, by node ID
300 for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
302 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
304 const int meshID = srcGroup->GetMesh()->GetPersistentId();
305 const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
306 if ( !srcMesh ) continue;
307 StdMeshers_Import_1D::getMaps( srcMesh, &theMesh, n2n, e2e );
309 const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
310 minGroupTol = std::min( groupTol, minGroupTol );
312 // clsfTol is 2D tolerance of a probe line
313 //GeomAdaptor_Surface S( surface );
314 // const double clsfTol = Min( S.UResolution( 0.1 * groupTol ), -- issue 0023092
315 // S.VResolution( 0.1 * groupTol ));
316 // another idea: try to use max tol of all edges
317 //const double clsfTol = 10 * BRep_Tool::Tolerance( geomFace ); // 0.1 * groupTol;
319 if ( helper.HasSeam() )
320 onEdgeClassifier.SetMesh( srcMesh->GetMeshDS() );
322 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
323 while ( srcElems->more() ) // loop on group contents
325 const SMDS_MeshElement* face = srcElems->next();
327 SMDS_MeshElement::iterator node = face->begin_nodes();
328 if ( bndBox3d.IsOut( SMESH_TNodeXYZ( *node )))
331 // find or create nodes of a new face
332 nodeState.resize( face->NbNodes() );
333 newNodes.resize( nodeState.size() );
335 int nbCreatedNodes = 0;
336 bool isOut = false, isIn = false; // if at least one node isIn - do not classify other nodes
337 for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
339 SMESH_TNodeXYZ nXYZ = *node;
340 nodeState[ i ] = TopAbs_UNKNOWN;
343 it_isnew = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 ));
344 n2nIt = it_isnew.first;
346 const SMDS_MeshNode* & newNode = n2nIt->second;
347 if ( !it_isnew.second && !newNode )
348 break; // a node is mapped to NULL - it is OUT of the FACE
352 if ( !subShapeIDs.count( newNode->getshapeId() ))
353 break; // node is Imported onto other FACE
354 if ( newNode->GetID() < (int) isNodeIn.size() &&
355 isNodeIn[ newNode->GetID() ])
357 if ( !isIn && bndNodes.count( *node ))
358 nodeState[ i ] = TopAbs_ON;
362 // find a node pre-existing on EDGE or VERTEX
363 dist2foundNodes.clear();
364 existingNodeOcTr.NodesAround( nXYZ, dist2foundNodes, groupTol );
365 if ( !dist2foundNodes.empty() )
367 newNode = dist2foundNodes.begin()->second;
368 nodeState[ i ] = TopAbs_ON;
374 // find out if node lies on the surface of theShape
375 gp_XY uv( Precision::Infinite(), 0 );
376 isOut = ( !helper.CheckNodeUV( geomFace, *node, uv, groupTol, /*force=*/true ) ||
377 bndBox2d.IsOut( uv ));
379 if ( !isOut && !isIn ) // classify
381 nodeState[i] = classifier.Perform( uv ); //classifier.Perform( geomFace, uv, clsfTol );
382 //nodeState[i] = classifier.State();
383 isOut = ( nodeState[i] == TopAbs_OUT );
384 if ( isOut && helper.IsOnSeam( uv ) && onEdgeClassifier.IsSatisfy( (*node)->GetID() ))
386 // uv.SetCoord( iCoo, helper.GetOtherParam( uv.Coord( iCoo )));
387 // classifier.Perform( geomFace, uv, clsfTol );
388 // nodeState[i] = classifier.State();
389 // isOut = ( nodeState[i] == TopAbs_OUT );
390 nodeState[i] = TopAbs_ON;
394 if ( !isOut ) // create a new node
396 newNode = tgtMesh->AddNode( nXYZ.X(), nXYZ.Y(), nXYZ.Z());
397 tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
399 if ( newNode->GetID() >= (int) isNodeIn.size() )
401 isNodeIn.push_back( false ); // allow allocate more than newNode->GetID()
402 isNodeIn.resize( newNode->GetID() + 1, false );
404 if ( nodeState[i] == TopAbs_ON )
405 bndNodes.insert( *node );
406 else if ( nodeState[i] != TopAbs_UNKNOWN )
407 isNodeIn[ newNode->GetID() ] = isIn = true;
410 if ( !(newNodes[i] = newNode ) || isOut )
413 } // loop on face nodes
415 if ( !newNodes.back() )
416 continue; // not all nodes of the face lie on theShape
418 if ( !isIn ) // if all nodes are on FACE boundary, a mesh face can be OUT
420 // check state of nodes created for other faces
421 for ( size_t i = 0; i < nodeState.size() && !isIn; ++i )
423 if ( nodeState[i] != TopAbs_UNKNOWN ) continue;
424 gp_XY uv = helper.GetNodeUV( geomFace, newNodes[i] );
425 nodeState[i] = classifier.Perform( uv ); //geomFace, uv, clsfTol );
426 //nodeState[i] = classifier.State();
427 isIn = ( nodeState[i] == TopAbs_IN );
429 if ( !isIn ) // classify face center
431 gp_XYZ gc( 0., 0., 0 );
432 for ( size_t i = 0; i < newNodes.size(); ++i )
433 gc += SMESH_TNodeXYZ( newNodes[i] );
434 gc /= newNodes.size();
437 GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( geomFace,
439 helper.MaxTolerance( geomFace ));
440 if ( !loc.IsIdentity() ) loc.Transformation().Inverted().Transforms( gc );
442 if ( !proj.IsDone() || proj.NbPoints() < 1 )
445 proj.LowerDistanceParameters(U,V);
447 //classifier.Perform( geomFace, uv, clsfTol );
448 TopAbs_State state = classifier.Perform( uv );
449 if ( state != TopAbs_IN )
454 // try to find already created face
455 SMDS_MeshElement * newFace = 0;
456 if ( nbCreatedNodes == 0 &&
457 tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
458 continue; // repeated face in source groups already created
460 // check future face orientation
461 const int nbCorners = face->NbCornerNodes();
462 const bool isQuad = ( nbCorners != (int) newNodes.size() );
469 gp_XY uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
470 surface->D1( uv.X(),uv.Y(), p, du,dv );
471 geomNorm = reverse ? dv^du : du^dv;
473 while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < nbCorners );
475 int iNext = helper.WrapIndex( iNode+1, nbCorners );
476 int iPrev = helper.WrapIndex( iNode-1, nbCorners );
478 SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
479 SMESH_TNodeXYZ curNode ( newNodes[iNode] );
480 SMESH_TNodeXYZ nextNode( newNodes[iNext] );
481 gp_Vec n1n0( prevNode - curNode);
482 gp_Vec n1n2( nextNode - curNode );
483 gp_Vec meshNorm = n1n2 ^ n1n0;
485 if ( geomNorm * meshNorm < 0 )
486 SMDS_MeshCell::applyInterlace
487 ( SMDS_MeshCell::reverseSmdsOrder( face->GetEntityType(), newNodes.size() ), newNodes );
491 if ( face->IsPoly() )
492 newFace = tgtMesh->AddPolygonalFace( newNodes );
494 switch ( newNodes.size() )
497 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
500 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
503 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
504 newNodes[3], newNodes[4], newNodes[5]);
507 newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
508 newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
512 tgtMesh->SetMeshElementOnShape( newFace, shapeID );
513 e2e->insert( make_pair( face, newFace ));
516 const SMDS_MeshNode* medium = 0;
517 for ( int i = 0; i < nbCorners; ++i )
519 const SMDS_MeshNode* n1 = newNodes[i];
520 const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbCorners ];
521 if ( isQuad ) // quadratic face
522 medium = newNodes[i+nbCorners];
523 link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
526 } // loop on group contents
528 // Remove OUT nodes from n2n map
529 for ( n2nIt = n2n->begin(); n2nIt != n2n->end(); )
530 if ( !n2nIt->second )
531 n2n->erase( n2nIt++ );
535 } // loop on src groups
537 // remove free nodes created on EDGEs
539 set<const SMDS_MeshNode*>::iterator node = bndNodes.begin();
540 for ( ; node != bndNodes.end(); ++node )
542 n2nIt = n2n->find( *node );
543 const SMDS_MeshNode* newNode = n2nIt->second;
544 if ( newNode && newNode->NbInverseElements() == 0 )
546 tgtMesh->RemoveFreeNode( newNode, 0, false );
552 // ==========================================================
553 // Put nodes on geom edges and create edges on them;
554 // check if the whole geom face is covered by imported faces
555 // ==========================================================
557 // use large tolerance for projection of nodes to edges because of
558 // BLSURF mesher specifics (issue 0020918, Study2.hdf)
559 const double projTol = minGroupTol;
561 bool isFaceMeshed = false;
562 SMESHDS_SubMesh* tgtFaceSM = tgtMesh->MeshElements( theShape );
565 // the imported mesh is valid if all external links (encountered once)
567 subShapeIDs.erase( shapeID ); // to contain edges and vertices only
569 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
571 const TLink& link = (*link2Nb).first;
572 int nbFaces = link2Nb->second;
575 // check if a not shared link lies on face boundary
576 bool nodesOnBoundary = true;
577 list< TopoDS_Shape > bndShapes;
578 for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
580 const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
581 if ( !subShapeIDs.count( n->getshapeId() )) // n is assigned to FACE
583 for ( size_t iE = 0; iE < edges.size(); ++iE )
584 if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
586 BRep_Tool::Range(edges[iE],f,l);
587 if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
588 // duplicated node on vertex
589 return error("Source elements overlap one another");
590 tgtFaceSM->RemoveNode( n );
591 tgtMesh->SetNodeOnEdge( n, edges[iE], u );
594 nodesOnBoundary = subShapeIDs.count( n->getshapeId());
596 if ( nodesOnBoundary )
598 TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
599 if ( s.ShapeType() == TopAbs_VERTEX )
600 bndShapes.push_front( s ); // vertex first
602 bndShapes.push_back( s ); // edges last
605 if ( !nodesOnBoundary )
607 error("free internal link"); // just for an easier debug
610 if ( bndShapes.front().ShapeType() == TopAbs_EDGE && // all link nodes are on EDGEs
611 bndShapes.front() != bndShapes.back() )
612 // link nodes on different geom edges
613 return error(COMPERR_BAD_INPUT_MESH, "Source nodes mismatch target vertices");
615 // find geom edge the link is on
616 if ( bndShapes.back().ShapeType() != TopAbs_EDGE ) // all link nodes are on VERTEXes
618 // find geom edge by two vertices
619 TopoDS_Shape geomEdge = helper.GetCommonAncestor( bndShapes.back(),
621 theMesh, TopAbs_EDGE );
622 if ( geomEdge.IsNull() )
624 error("free internal link");
625 break; // vertices belong to different edges
627 bndShapes.push_back( geomEdge );
630 // create an edge if not yet exists
632 newNodes[0] = link.node1(), newNodes[1] = link.node2();
633 const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
634 if ( edge ) continue;
636 if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
639 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], link._medium );
641 TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
642 helper.CheckNodeU( geomEdge, link._medium, u, projTol, /*force=*/true );
643 tgtFaceSM->RemoveNode( link._medium );
644 tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
648 edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
653 tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
655 else if ( nbFaces > 2 )
657 return error( COMPERR_BAD_INPUT_MESH, "Non-manifold source mesh");
660 isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
663 // check that source faces do not overlap:
664 // there must be only two edges sharing each vertex and bound to sub-edges of theShape
665 SMESH_MeshEditor editor( &theMesh );
666 set<int>::iterator subID = subShapeIDs.begin();
667 for ( ; subID != subShapeIDs.end(); ++subID )
669 const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
670 if ( s.ShapeType() != TopAbs_VERTEX ) continue;
671 const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
672 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
674 while ( eIt->more() )
676 const SMDS_MeshElement* edge = eIt->next();
677 int sId = editor.FindShape( edge );
678 nbEdges += subShapeIDs.count( sId );
680 if ( nbEdges < 2 && !helper.IsRealSeam( s ))
681 return false; // weird
683 return error( COMPERR_BAD_INPUT_MESH, "Source elements overlap one another");
688 return error( COMPERR_BAD_INPUT_MESH,
689 "Source elements don't cover totally the geometrical face" );
691 if ( helper.HasSeam() )
693 // links on seam edges are shared by two faces, so no edges were created on them
694 // by the previous detection of 2D mesh boundary
695 for ( size_t iE = 0; iE < edges.size(); ++iE )
697 if ( !helper.IsRealSeam( edges[iE] )) continue;
698 const TopoDS_Edge& seamEdge = edges[iE];
699 // to find nodes lying on the seamEdge we check nodes of mesh faces sharing a node on one
700 // of its vertices; after finding another node on seamEdge we continue the same way
701 // until finding all nodes.
702 TopoDS_Vertex seamVertex = helper.IthVertex( 0, seamEdge );
703 const SMDS_MeshNode* vertNode = SMESH_Algo::VertexNode( seamVertex, tgtMesh );
704 set< const SMDS_MeshNode* > checkedNodes; checkedNodes.insert( vertNode );
705 set< const SMDS_MeshElement* > checkedFaces;
706 // as a face can have more than one node on the seamEdge, there is a difficulty in selecting
707 // one of those nodes to treat next; so we simply find all nodes on the seamEdge and
708 // then sort them by U on edge
709 typedef list< pair< double, const SMDS_MeshNode* > > TUNodeList;
710 TUNodeList nodesOnSeam;
711 double u = helper.GetNodeU( seamEdge, vertNode );
712 nodesOnSeam.push_back( make_pair( u, vertNode ));
713 TUNodeList::iterator u2nIt = nodesOnSeam.begin();
714 for ( ; u2nIt != nodesOnSeam.end(); ++u2nIt )
716 const SMDS_MeshNode* startNode = (*u2nIt).second;
717 SMDS_ElemIteratorPtr faceIt = startNode->GetInverseElementIterator( SMDSAbs_Face );
718 while ( faceIt->more() )
720 const SMDS_MeshElement* face = faceIt->next();
721 if ( !checkedFaces.insert( face ).second ) continue;
722 for ( int i = 0, nbNodes = face->NbCornerNodes(); i < nbNodes; ++i )
724 const SMDS_MeshNode* n = face->GetNode( i );
725 if ( n == startNode || !checkedNodes.insert( n ).second ) continue;
726 if ( helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true ))
727 nodesOnSeam.push_back( make_pair( u, n ));
731 // sort the found nodes by U on the seamEdge; most probably they are in a good order,
732 // so we can use the hint to spead-up map filling
733 map< double, const SMDS_MeshNode* > u2nodeMap;
734 for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
735 u2nodeMap.insert( u2nodeMap.end(), *u2nIt );
739 SMESH_MesherHelper seamHelper( theMesh );
740 seamHelper.SetSubShape( edges[ iE ]);
741 seamHelper.SetElementsOnShape( true );
743 if ( !checkedFaces.empty() && (*checkedFaces.begin())->IsQuadratic() )
744 for ( set< const SMDS_MeshElement* >::iterator fIt = checkedFaces.begin();
745 fIt != checkedFaces.end(); ++fIt )
746 seamHelper.AddTLinks( static_cast<const SMDS_MeshFace*>( *fIt ));
748 map< double, const SMDS_MeshNode* >::iterator n1, n2, u2nEnd = u2nodeMap.end();
749 for ( n2 = u2nodeMap.begin(), n1 = n2++; n2 != u2nEnd; ++n1, ++n2 )
751 const SMDS_MeshNode* node1 = n1->second;
752 const SMDS_MeshNode* node2 = n2->second;
753 seamHelper.AddEdge( node1, node2 );
754 if ( node2->getshapeId() == helper.GetSubShapeID() )
756 tgtFaceSM->RemoveNode( node2 );
757 tgtMesh->SetNodeOnEdge( const_cast<SMDS_MeshNode*>( node2 ), seamEdge, n2->first );
761 } // loop on edges to find seam ones
762 } // if ( helper.HasSeam() )
764 // notify sub-meshes of edges on computation
765 for ( size_t iE = 0; iE < edges.size(); ++iE )
767 SMESH_subMesh * sm = theMesh.GetSubMesh( edges[iE] );
768 // if ( SMESH_Algo::isDegenerated( edges[iE] ))
769 // sm->SetIsAlwaysComputed( true );
770 sm->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
771 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
772 return error(SMESH_Comment("Failed to create segments on the edge #") << sm->GetId());
779 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
780 for ( size_t i = 0; i < srcMeshes.size(); ++i )
781 StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
786 //=============================================================================
788 * \brief Set needed event listeners and create a submesh for a copied mesh
790 * This method is called only if a submesh has HYP_OK algo_state.
792 //=============================================================================
794 void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
798 const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
799 SMESH_Mesh* tgtMesh = subMesh->GetFather();
800 Hypothesis_Status aStatus;
801 CheckHypothesis( *tgtMesh, tgtShape, aStatus );
803 StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
805 void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
807 SetEventListener(subMesh);
810 //=============================================================================
812 * Predict nb of mesh entities created by Compute()
814 //=============================================================================
816 bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh & theMesh,
817 const TopoDS_Shape & theShape,
818 MapShapeNbElems& aResMap)
820 if ( !_sourceHyp ) return false;
822 const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
823 if ( srcGroups.empty() )
824 return error("Invalid source groups");
826 vector<int> aVec(SMDSEntity_Last,0);
828 bool toCopyMesh, toCopyGroups;
829 _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
830 if ( toCopyMesh ) // the whole mesh is copied
832 vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
833 for ( unsigned i = 0; i < srcMeshes.size(); ++i )
835 SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
836 if ( !sm || aResMap.count( sm )) continue; // already counted
837 const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
838 for (int i = 0; i < SMDSEntity_Last; i++)
839 aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
844 // std-like iterator used to get coordinates of nodes of mesh element
845 typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
847 SMESH_MesherHelper helper(theMesh);
848 helper.SetSubShape(theShape);
850 const TopoDS_Face& geomFace = TopoDS::Face( theShape );
852 // take into account nodes on vertices
853 TopExp_Explorer exp( theShape, TopAbs_VERTEX );
854 for ( ; exp.More(); exp.Next() )
855 theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
857 // to count now many times a link between nodes encounters,
858 // negative nb additionally means that a link is quadratic
859 map<SMESH_TLink, int> linkCount;
860 map<SMESH_TLink, int>::iterator link2Nb;
862 // count faces and nodes imported from groups
863 set<const SMDS_MeshNode* > allNodes;
865 double minGroupTol = 1e100;
866 for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
868 const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
869 const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
870 minGroupTol = std::min( groupTol, minGroupTol );
871 SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
872 SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
873 while ( srcElems->more() ) // loop on group contents
875 const SMDS_MeshElement* face = srcElems->next();
876 // find out if face is located on geomEdge by projecting
877 // a gravity center of face to geomFace
879 gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
880 tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
881 if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
883 ++aVec[ face->GetEntityType() ];
886 int nbConers = face->NbCornerNodes();
887 for ( int i = 0; i < face->NbNodes(); ++i )
889 const SMDS_MeshNode* n1 = face->GetNode(i);
890 allNodes.insert( n1 );
893 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
894 link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
895 if ( (*link2Nb).second )
896 link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
898 link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
903 helper.GetMeshDS()->RemoveNode(tmpNode);
906 int nbNodes = allNodes.size();
909 // count nodes and edges on geom edges
912 for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
914 TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
915 SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
916 vector<int>& edgeVec = aResMap[sm];
917 if ( edgeVec.empty() )
919 edgeVec.resize(SMDSEntity_Last,0);
920 for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
922 const SMESH_TLink& link = (*link2Nb).first;
923 int nbFacesOfLink = Abs( link2Nb->second );
924 bool eraseLink = ( nbFacesOfLink != 1 );
925 if ( nbFacesOfLink == 1 )
927 if ( helper.CheckNodeU( geomEdge, link.node1(), u, minGroupTol, /*force=*/true )&&
928 helper.CheckNodeU( geomEdge, link.node2(), u, minGroupTol, /*force=*/true ))
930 bool isQuadratic = ( link2Nb->second < 0 );
931 ++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
932 ++edgeVec[ SMDSEntity_Node ];
938 linkCount.erase(link2Nb++);
942 if ( edgeVec[ SMDSEntity_Node] > 0 )
943 --edgeVec[ SMDSEntity_Node ]; // for one node on vertex
945 else if ( !helper.IsSeamShape( geomEdge ) ||
946 geomEdge.Orientation() == TopAbs_FORWARD )
948 nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
952 aVec[SMDSEntity_Node] = nbNodes;
955 SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
956 aResMap.insert(make_pair(sm,aVec));