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