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