Salome HOME
22564: EDF SMESH: Regression, FindNodeClosestTo does not return the right node on...
[modules/smesh.git] / src / StdMeshers / StdMeshers_Import_1D2D.cxx
1 // Copyright (C) 2007-2014  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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.
10 //
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.
15 //
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
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  SMESH SMESH : implementaion of SMESH idl descriptions
24 //  File   : StdMeshers_Import_1D2D.cxx
25 //  Module : SMESH
26 //
27 #include "StdMeshers_Import_1D2D.hxx"
28
29 #include "StdMeshers_Import_1D.hxx"
30 #include "StdMeshers_ImportSource.hxx"
31
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_OctreeNode.hxx"
42 #include "SMESH_subMesh.hxx"
43
44 #include "Utils_SALOME_Exception.hxx"
45 #include "utilities.h"
46
47 #include <BRep_Builder.hxx>
48 #include <BRep_Tool.hxx>
49 #include <TopExp.hxx>
50 #include <TopExp_Explorer.hxx>
51 #include <TopoDS.hxx>
52 #include <TopoDS_Compound.hxx>
53 #include <TopoDS_Edge.hxx>
54 #include <TopoDS_Vertex.hxx>
55 #include <Precision.hxx>
56
57 #include <numeric>
58
59 using namespace std;
60
61 namespace
62 {
63   double getMinElemSize2( const SMESHDS_GroupBase* srcGroup )
64   {
65     double minSize2 = 1e100;
66     SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
67     while ( srcElems->more() ) // loop on group contents
68     {
69       const SMDS_MeshElement* face = srcElems->next();
70       int nbN = face->NbCornerNodes();
71
72       SMESH_TNodeXYZ prevN( face->GetNode( nbN-1 ));
73       for ( int i = 0; i < nbN; ++i )
74       {
75         SMESH_TNodeXYZ n( face->GetNode( i ) );
76         double size2 = ( n - prevN ).SquareModulus();
77         minSize2 = std::min( minSize2, size2 );
78         prevN = n;
79       }
80     }
81     return minSize2;
82   }
83 }
84
85 //=============================================================================
86 /*!
87  * Creates StdMeshers_Import_1D2D
88  */
89 //=============================================================================
90
91 StdMeshers_Import_1D2D::StdMeshers_Import_1D2D(int hypId, int studyId, SMESH_Gen * gen)
92   :SMESH_2D_Algo(hypId, studyId, gen), _sourceHyp(0)
93 {
94   MESSAGE("StdMeshers_Import_1D2D::StdMeshers_Import_1D2D");
95   _name = "Import_1D2D";
96   _shapeType = (1 << TopAbs_FACE);
97
98   _compatibleHypothesis.push_back("ImportSource2D");
99   _requireDiscreteBoundary = false;
100 }
101
102 //=============================================================================
103 /*!
104  * Check presence of a hypothesis
105  */
106 //=============================================================================
107
108 bool StdMeshers_Import_1D2D::CheckHypothesis
109                          (SMESH_Mesh&                          aMesh,
110                           const TopoDS_Shape&                  aShape,
111                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
112 {
113   _sourceHyp = 0;
114
115   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
116   if ( hyps.size() == 0 )
117   {
118     aStatus = SMESH_Hypothesis::HYP_MISSING;
119     return false;  // can't work with no hypothesis
120   }
121
122   if ( hyps.size() > 1 )
123   {
124     aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
125     return false;
126   }
127
128   const SMESHDS_Hypothesis *theHyp = hyps.front();
129
130   string hypName = theHyp->GetName();
131
132   if (hypName == _compatibleHypothesis.front())
133   {
134     _sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
135     aStatus = SMESH_Hypothesis::HYP_OK;
136     return true;
137   }
138
139   aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
140   return true;
141 }
142
143 namespace
144 {
145   /*!
146    * \brief OrientedLink additionally storing a medium node
147    */
148   struct TLink : public SMESH_OrientedLink
149   {
150     const SMDS_MeshNode* _medium;
151     TLink( const SMDS_MeshNode* n1,
152            const SMDS_MeshNode* n2,
153            const SMDS_MeshNode* medium=0)
154       : SMESH_OrientedLink( n1,n2 ), _medium( medium ) {}
155   };
156 }
157
158 //=============================================================================
159 /*!
160  * Import elements from the other mesh 
161  */
162 //=============================================================================
163
164 bool StdMeshers_Import_1D2D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
165 {
166   if ( !_sourceHyp ) return false;
167
168   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups(/*loaded=*/true);
169   if ( srcGroups.empty() )
170     return error("Invalid source groups");
171
172   bool allGroupsEmpty = true;
173   for ( size_t iG = 0; iG < srcGroups.size() && allGroupsEmpty; ++iG )
174     allGroupsEmpty = srcGroups[iG]->GetGroupDS()->IsEmpty();
175   if ( allGroupsEmpty )
176     return error("No faces in source groups");
177
178   SMESH_MesherHelper helper(theMesh);
179   helper.SetSubShape(theShape);
180   SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
181
182   const TopoDS_Face& geomFace = TopoDS::Face( theShape );
183   const double faceTol = helper.MaxTolerance( geomFace );
184   const int shapeID = tgtMesh->ShapeToIndex( geomFace );
185   const bool toCheckOri = (helper.NbAncestors( geomFace, theMesh, TopAbs_SOLID ) == 1 );
186
187   Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
188   const bool reverse = 
189     ( helper.GetSubShapeOri( tgtMesh->ShapeToMesh(), geomFace) == TopAbs_REVERSED );
190   gp_Pnt p; gp_Vec du, dv;
191
192   set<int> subShapeIDs;
193   subShapeIDs.insert( shapeID );
194
195   // nodes already existing on sub-shapes of the FACE
196   TIDSortedNodeSet existingNodes;
197
198   // get/make nodes on vertices and add them to existingNodes
199   TopExp_Explorer exp( theShape, TopAbs_VERTEX );
200   for ( ; exp.More(); exp.Next() )
201   {
202     const TopoDS_Vertex& v = TopoDS::Vertex( exp.Current() );
203     if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
204       continue;
205     const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
206     if ( !n )
207     {
208       _gen->Compute(theMesh,v,/*anUpward=*/true);
209       n = SMESH_Algo::VertexNode( v, tgtMesh );
210       if ( !n ) return false; // very strange
211     }
212     existingNodes.insert( n );
213   }
214
215   // get EDGESs and their ids and get existing nodes on EDGEs
216   vector< TopoDS_Edge > edges;
217   for ( exp.Init( theShape, TopAbs_EDGE ); exp.More(); exp.Next() )
218   {
219     const TopoDS_Edge & edge = TopoDS::Edge( exp.Current() );
220     if ( !SMESH_Algo::isDegenerated( edge ))
221       if ( subShapeIDs.insert( tgtMesh->ShapeToIndex( edge )).second )
222       {
223         edges.push_back( edge );
224         if ( SMESHDS_SubMesh* eSM = tgtMesh->MeshElements( edge ))
225         {
226           typedef SMDS_StdIterator< const SMDS_MeshNode*, SMDS_NodeIteratorPtr > iterator;
227           existingNodes.insert( iterator( eSM->GetNodes() ), iterator() );
228         }
229       }
230   }
231   // octree to find existing nodes
232   SMESH_OctreeNode existingNodeOcTr( existingNodes );
233   std::map<double, const SMDS_MeshNode*> dist2foundNodes;
234
235   // to count now many times a link between nodes encounters
236   map<TLink, int> linkCount;
237   map<TLink, int>::iterator link2Nb;
238   double minGroupTol = Precision::Infinite();
239
240   // =========================
241   // Import faces from groups
242   // =========================
243
244   StdMeshers_Import_1D::TNodeNodeMap* n2n;
245   StdMeshers_Import_1D::TElemElemMap* e2e;
246   vector<const SMDS_MeshNode*> newNodes;
247   for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
248   {
249     const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
250
251     const int meshID = srcGroup->GetMesh()->GetPersistentId();
252     const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
253     if ( !srcMesh ) continue;
254     StdMeshers_Import_1D::getMaps( srcMesh, &theMesh, n2n, e2e );
255
256     const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
257     minGroupTol = std::min( groupTol, minGroupTol );
258
259     SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
260     SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
261     gp_XY uv( Precision::Infinite(), Precision::Infinite() );
262     while ( srcElems->more() ) // loop on group contents
263     {
264       const SMDS_MeshElement* face = srcElems->next();
265       // find or create nodes of a new face
266       newNodes.resize( face->NbNodes() );
267       newNodes.back() = 0;
268       int nbCreatedNodes = 0;
269       SMDS_MeshElement::iterator node = face->begin_nodes();
270       for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
271       {
272         StdMeshers_Import_1D::TNodeNodeMap::iterator n2nIt =
273           n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
274         if ( n2nIt->second )
275         {
276           if ( !subShapeIDs.count( n2nIt->second->getshapeId() )) // node already on an EDGE
277             break;
278         }
279         else
280         {
281           // find a pre-existing node
282           dist2foundNodes.clear();
283           existingNodeOcTr.NodesAround( SMESH_TNodeXYZ( *node ), dist2foundNodes, groupTol );
284           if ( !dist2foundNodes.empty() )
285             (*n2nIt).second = dist2foundNodes.begin()->second;
286         }
287         if ( !n2nIt->second )
288         {
289           // find out if node lies on theShape
290           tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
291           uv.SetCoord( Precision::Infinite(), Precision::Infinite() );
292           if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
293           {
294             SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
295             n2nIt->second = newNode;
296             tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
297             nbCreatedNodes++;
298           }
299         }
300         if ( !(newNodes[i] = n2nIt->second ))
301           break;
302       }
303       if ( !newNodes.back() )
304         continue; // not all nodes of the face lie on theShape
305
306       // try to find already created face
307       SMDS_MeshElement * newFace = 0;
308       if ( nbCreatedNodes == 0 &&
309            tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
310         continue; // repeated face in source groups already created 
311
312       // check future face orientation
313       const int nbCorners = face->NbCornerNodes();
314       const bool isQuad   = ( nbCorners != (int) newNodes.size() );
315       if ( toCheckOri )
316       {
317         int iNode = -1;
318         gp_Vec geomNorm;
319         do
320         {
321           uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
322           surface->D1( uv.X(),uv.Y(), p, du,dv );
323           geomNorm = reverse ? dv^du : du^dv;
324         }
325         while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < nbCorners );
326
327         int iNext = helper.WrapIndex( iNode+1, nbCorners );
328         int iPrev = helper.WrapIndex( iNode-1, nbCorners );
329
330         SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
331         SMESH_TNodeXYZ curNode ( newNodes[iNode] );
332         SMESH_TNodeXYZ nextNode( newNodes[iNext] );
333         gp_Vec n1n0( prevNode - curNode);
334         gp_Vec n1n2( nextNode - curNode );
335         gp_Vec meshNorm = n1n2 ^ n1n0;
336
337         if ( geomNorm * meshNorm < 0 )
338           SMDS_MeshCell::applyInterlace
339             ( SMDS_MeshCell::reverseSmdsOrder( face->GetEntityType() ), newNodes );
340       }
341
342       // make a new face
343       if ( face->IsPoly() )
344         newFace = tgtMesh->AddPolygonalFace( newNodes );
345       else
346         switch ( newNodes.size() )
347         {
348         case 3:
349           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
350           break;
351         case 4:
352           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
353           break;
354         case 6:
355           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
356                                       newNodes[3], newNodes[4], newNodes[5]);
357           break;
358         case 8:
359           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
360                                       newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
361           break;
362         default: continue;
363         }
364       tgtMesh->SetMeshElementOnShape( newFace, shapeID );
365       e2e->insert( make_pair( face, newFace ));
366
367       // collect links
368       const SMDS_MeshNode* medium = 0;
369       for ( int i = 0; i < nbCorners; ++i )
370       {
371         const SMDS_MeshNode* n1 = newNodes[i];
372         const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbCorners ];
373         if ( isQuad ) // quadratic face
374           medium = newNodes[i+nbCorners];
375         link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
376         ++link2Nb->second;
377         // if ( link2Nb->second == 1 )
378         // {
379         //   // measure link length
380         //   double len2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
381         //   if ( len2 < minGroupTol )
382         //     minGroupTol = len2;
383         // }
384       }
385     }
386     helper.GetMeshDS()->RemoveNode(tmpNode);
387   }
388
389   // ==========================================================
390   // Put nodes on geom edges and create edges on them;
391   // check if the whole geom face is covered by imported faces
392   // ==========================================================
393
394   // use large tolerance for projection of nodes to edges because of
395   // BLSURF mesher specifics (issue 0020918, Study2.hdf)
396   const double projTol = minGroupTol;
397
398   bool isFaceMeshed = false;
399   SMESHDS_SubMesh* tgtFaceSM = tgtMesh->MeshElements( theShape );
400   if ( tgtFaceSM )
401   {
402     // the imported mesh is valid if all external links (encountered once)
403     // lie on geom edges
404     subShapeIDs.erase( shapeID ); // to contain edges and vertices only
405     double u, f, l;
406     for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
407     {
408       const TLink& link = (*link2Nb).first;
409       int nbFaces = link2Nb->second;
410       if ( nbFaces == 1 )
411       {
412         // check if a not shared link lies on face boundary
413         bool nodesOnBoundary = true;
414         list< TopoDS_Shape > bndShapes;
415         for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
416         {
417           const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
418           if ( !subShapeIDs.count( n->getshapeId() )) // n is assigned to FACE
419           {
420             for ( size_t iE = 0; iE < edges.size(); ++iE )
421               if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
422               {
423                 BRep_Tool::Range(edges[iE],f,l);
424                 if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
425                   // duplicated node on vertex
426                   return error("Source elements overlap one another");
427                 tgtFaceSM->RemoveNode( n, /*isNodeDeleted=*/false );
428                 tgtMesh->SetNodeOnEdge( n, edges[iE], u );
429                 break;
430               }
431             nodesOnBoundary = subShapeIDs.count( n->getshapeId());
432           }
433           if ( nodesOnBoundary )
434           {
435             TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
436             if ( s.ShapeType() == TopAbs_VERTEX )
437               bndShapes.push_front( s ); // vertex first
438             else
439               bndShapes.push_back( s ); // edges last
440           }
441         }
442         if ( !nodesOnBoundary )
443         {
444           error("free internal link"); // just for an easier debug
445           break;
446         }
447         if ( bndShapes.front().ShapeType() == TopAbs_EDGE && // all link nodes are on EDGEs
448              bndShapes.front() != bndShapes.back() )
449           // link nodes on different geom edges
450           return error(COMPERR_BAD_INPUT_MESH, "Source nodes mismatch target vertices");
451
452         // find geom edge the link is on
453         if ( bndShapes.back().ShapeType() != TopAbs_EDGE ) // all link nodes are on VERTEXes
454         {
455           // find geom edge by two vertices
456           TopoDS_Shape geomEdge = helper.GetCommonAncestor( bndShapes.back(),
457                                                             bndShapes.front(),
458                                                             theMesh, TopAbs_EDGE );
459           if ( geomEdge.IsNull() )
460           {
461             error("free internal link");
462             break; // vertices belong to different edges
463           }
464           bndShapes.push_back( geomEdge );
465         }
466
467         // create an edge if not yet exists
468         newNodes.resize(2);
469         newNodes[0] = link.node1(), newNodes[1] = link.node2();
470         const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
471         if ( edge ) continue;
472
473         if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
474         if ( link._medium )
475         {
476           edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], link._medium );
477
478           TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
479           helper.CheckNodeU( geomEdge, link._medium, u, projTol, /*force=*/true );
480           tgtFaceSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
481           tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
482         }
483         else
484         {
485           edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
486         }
487         if ( !edge )
488           return false;
489
490         tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
491       }
492       else if ( nbFaces > 2 )
493       {
494         return error( COMPERR_BAD_INPUT_MESH, "Non-manifold source mesh");
495       }
496     }
497     isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
498     if ( isFaceMeshed )
499     {
500       // check that source faces do not overlap:
501       // there must be only two edges sharing each vertex and bound to sub-edges of theShape
502       SMESH_MeshEditor editor( &theMesh );
503       set<int>::iterator subID = subShapeIDs.begin();
504       for ( ; subID != subShapeIDs.end(); ++subID )
505       {
506         const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
507         if ( s.ShapeType() != TopAbs_VERTEX ) continue;
508         const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
509         SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
510         int nbEdges = 0;
511         while ( eIt->more() )
512         {
513           const SMDS_MeshElement* edge = eIt->next();
514           int sId = editor.FindShape( edge );
515           nbEdges += subShapeIDs.count( sId );
516         }
517         if ( nbEdges < 2 && !helper.IsRealSeam( s ))
518           return false; // weird
519         if ( nbEdges > 2 )
520           return error( COMPERR_BAD_INPUT_MESH, "Source elements overlap one another");
521       }
522     }
523   }
524   if ( !isFaceMeshed )
525     return error( COMPERR_BAD_INPUT_MESH,
526                   "Source elements don't cover totally the geometrical face" );
527
528   if ( helper.HasSeam() )
529   {
530     // links on seam edges are shared by two faces, so no edges were created on them
531     // by the previous detection of 2D mesh boundary
532     for ( size_t iE = 0; iE < edges.size(); ++iE )
533     {
534       if ( !helper.IsRealSeam( edges[iE] )) continue;
535       const TopoDS_Edge& seamEdge = edges[iE];
536       // to find nodes lying on the seamEdge we check nodes of mesh faces sharing a node on one
537       // of its vertices; after finding another node on seamEdge we continue the same way
538       // until finding all nodes.
539       TopoDS_Vertex      seamVertex = helper.IthVertex( 0, seamEdge );
540       const SMDS_MeshNode* vertNode = SMESH_Algo::VertexNode( seamVertex, tgtMesh );
541       set< const SMDS_MeshNode* > checkedNodes; checkedNodes.insert( vertNode );
542       set< const SMDS_MeshElement* > checkedFaces;
543       // as a face can have more than one node on the seamEdge, there is a difficulty in selecting
544       // one of those nodes to treat next; so we simply find all nodes on the seamEdge and
545       // then sort them by U on edge
546       typedef list< pair< double, const SMDS_MeshNode* > > TUNodeList;
547       TUNodeList nodesOnSeam;
548       double u = helper.GetNodeU( seamEdge, vertNode );
549       nodesOnSeam.push_back( make_pair( u, vertNode ));
550       TUNodeList::iterator u2nIt = nodesOnSeam.begin();
551       for ( ; u2nIt != nodesOnSeam.end(); ++u2nIt )
552       {
553         const SMDS_MeshNode* startNode = (*u2nIt).second;
554         SMDS_ElemIteratorPtr faceIt = startNode->GetInverseElementIterator( SMDSAbs_Face );
555         while ( faceIt->more() )
556         {
557           const SMDS_MeshElement* face = faceIt->next();
558           if ( !checkedFaces.insert( face ).second ) continue;
559           for ( int i = 0, nbNodes = face->NbCornerNodes(); i < nbNodes; ++i )
560           {
561             const SMDS_MeshNode* n = face->GetNode( i );
562             if ( n == startNode || !checkedNodes.insert( n ).second ) continue;
563             if ( helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true ))
564               nodesOnSeam.push_back( make_pair( u, n ));
565           }
566         }
567       }
568       // sort the found nodes by U on the seamEdge; most probably they are in a good order,
569       // so we can use the hint to spead-up map filling
570       map< double, const SMDS_MeshNode* > u2nodeMap;
571       for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
572         u2nodeMap.insert( u2nodeMap.end(), *u2nIt );
573
574       // create edges
575       {
576         SMESH_MesherHelper seamHelper( theMesh );
577         seamHelper.SetSubShape( edges[ iE ]);
578         seamHelper.SetElementsOnShape( true );
579
580         if ( (*checkedFaces.begin())->IsQuadratic() )
581           for ( set< const SMDS_MeshElement* >::iterator fIt = checkedFaces.begin();
582                 fIt != checkedFaces.end(); ++fIt )
583             seamHelper.AddTLinks( static_cast<const SMDS_MeshFace*>( *fIt ));
584
585         map< double, const SMDS_MeshNode* >::iterator n1, n2, u2nEnd = u2nodeMap.end();
586         for ( n2 = u2nodeMap.begin(), n1 = n2++; n2 != u2nEnd; ++n1, ++n2 )
587         {
588           const SMDS_MeshNode* node1 = n1->second;
589           const SMDS_MeshNode* node2 = n2->second;
590           seamHelper.AddEdge( node1, node2 );
591           if ( node2->getshapeId() == helper.GetSubShapeID() )
592           {
593             tgtFaceSM->RemoveNode( node2, /*isNodeDeleted=*/false );
594             tgtMesh->SetNodeOnEdge( const_cast<SMDS_MeshNode*>( node2 ), seamEdge, n2->first );
595           }
596         }
597       }
598     } // loop on edges to find seam ones
599   } // if ( helper.HasSeam() )
600
601   // notify sub-meshes of edges on computation
602   for ( size_t iE = 0; iE < edges.size(); ++iE )
603   {
604     SMESH_subMesh * sm = theMesh.GetSubMesh( edges[iE] );
605     // if ( SMESH_Algo::isDegenerated( edges[iE] ))
606     //   sm->SetIsAlwaysComputed( true );
607     sm->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
608     if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
609       return error(SMESH_Comment("Failed to create segments on the edge ")
610                    << tgtMesh->ShapeToIndex( edges[iE ]));
611   }
612
613   // ============
614   // Copy meshes
615   // ============
616
617   vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
618   for ( size_t i = 0; i < srcMeshes.size(); ++i )
619     StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
620
621   return true;
622 }
623
624 //=============================================================================
625 /*!
626  * \brief Set needed event listeners and create a submesh for a copied mesh
627  *
628  * This method is called only if a submesh has HYP_OK algo_state.
629  */
630 //=============================================================================
631
632 void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
633 {
634   if ( !_sourceHyp )
635   {
636     const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
637     SMESH_Mesh*         tgtMesh  = subMesh->GetFather();
638     Hypothesis_Status aStatus;
639     CheckHypothesis( *tgtMesh, tgtShape, aStatus );
640   }
641   StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
642 }
643 void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
644 {
645   SetEventListener(subMesh);
646 }
647
648 //=============================================================================
649 /*!
650  * Predict nb of mesh entities created by Compute()
651  */
652 //=============================================================================
653
654 bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh &         theMesh,
655                                       const TopoDS_Shape & theShape,
656                                       MapShapeNbElems&     aResMap)
657 {
658   if ( !_sourceHyp ) return false;
659
660   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
661   if ( srcGroups.empty() )
662     return error("Invalid source groups");
663
664   vector<int> aVec(SMDSEntity_Last,0);
665
666   bool toCopyMesh, toCopyGroups;
667   _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
668   if ( toCopyMesh ) // the whole mesh is copied
669   {
670     vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
671     for ( unsigned i = 0; i < srcMeshes.size(); ++i )
672     {
673       SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
674       if ( !sm || aResMap.count( sm )) continue; // already counted
675       const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
676       for (int i = 0; i < SMDSEntity_Last; i++)
677         aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
678     }
679   }
680   else
681   {
682     // std-like iterator used to get coordinates of nodes of mesh element
683     typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
684
685     SMESH_MesherHelper helper(theMesh);
686     helper.SetSubShape(theShape);
687
688     const TopoDS_Face& geomFace = TopoDS::Face( theShape );
689
690     // take into account nodes on vertices
691     TopExp_Explorer exp( theShape, TopAbs_VERTEX );
692     for ( ; exp.More(); exp.Next() )
693       theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
694
695     // to count now many times a link between nodes encounters,
696     // negative nb additionally means that a link is quadratic
697     map<SMESH_TLink, int> linkCount;
698     map<SMESH_TLink, int>::iterator link2Nb;
699
700     // count faces and nodes imported from groups
701     set<const SMDS_MeshNode* > allNodes;
702     gp_XY uv;
703     double minGroupTol = 1e100;
704     for ( int iG = 0; iG < srcGroups.size(); ++iG )
705     {
706       const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
707       const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
708       minGroupTol = std::min( groupTol, minGroupTol );
709       SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
710       SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
711       while ( srcElems->more() ) // loop on group contents
712       {
713         const SMDS_MeshElement* face = srcElems->next();
714         // find out if face is located on geomEdge by projecting
715         // a gravity center of face to geomFace
716         gp_XYZ gc(0,0,0);
717         gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
718         tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
719         if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
720         {
721           ++aVec[ face->GetEntityType() ];
722
723           // collect links
724           int nbConers = face->NbCornerNodes();
725           for ( int i = 0; i < face->NbNodes(); ++i )
726           {
727             const SMDS_MeshNode* n1 = face->GetNode(i);
728             allNodes.insert( n1 );
729             if ( i < nbConers )
730             {
731               const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
732               link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
733               if ( (*link2Nb).second )
734                 link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
735               else
736                 link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
737             }
738           }
739         }
740       }
741       helper.GetMeshDS()->RemoveNode(tmpNode);
742     }
743
744     int nbNodes = allNodes.size();
745     allNodes.clear();
746
747     // count nodes and edges on geom edges
748
749     double u;
750     for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
751     {
752       TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
753       SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
754       vector<int>& edgeVec = aResMap[sm];
755       if ( edgeVec.empty() )
756       {
757         edgeVec.resize(SMDSEntity_Last,0);
758         for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
759         {
760           const SMESH_TLink& link = (*link2Nb).first;
761           int nbFacesOfLink = Abs( link2Nb->second );
762           bool eraseLink = ( nbFacesOfLink != 1 );
763           if ( nbFacesOfLink == 1 )
764           {
765             if ( helper.CheckNodeU( geomEdge, link.node1(), u, minGroupTol, /*force=*/true )&&
766                  helper.CheckNodeU( geomEdge, link.node2(), u, minGroupTol, /*force=*/true ))
767             {
768               bool isQuadratic = ( link2Nb->second < 0 );
769               ++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
770               ++edgeVec[ SMDSEntity_Node ];
771               --nbNodes;
772               eraseLink = true;
773             }
774           }
775           if ( eraseLink )
776             linkCount.erase(link2Nb++);
777           else
778             link2Nb++;
779         }
780         if ( edgeVec[ SMDSEntity_Node] > 0 )
781           --edgeVec[ SMDSEntity_Node ]; // for one node on vertex
782       }
783       else if ( !helper.IsSeamShape( geomEdge ) ||
784                 geomEdge.Orientation() == TopAbs_FORWARD )
785       {
786         nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
787       }
788     }
789
790     aVec[SMDSEntity_Node] = nbNodes;
791   }
792
793   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
794   aResMap.insert(make_pair(sm,aVec));
795
796   return true;
797 }