Salome HOME
0022523: [CEA 1096] Add a colomn "biquadratic" in "Mesh Information"
[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           if ( existingNodeOcTr.NodesAround( SMESH_TNodeXYZ( *node ), dist2foundNodes, groupTol ))
284             (*n2nIt).second = dist2foundNodes.begin()->second;
285         }
286         if ( !n2nIt->second )
287         {
288           // find out if node lies on theShape
289           tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
290           uv.SetCoord( Precision::Infinite(), Precision::Infinite() );
291           if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
292           {
293             SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
294             n2nIt->second = newNode;
295             tgtMesh->SetNodeOnFace( newNode, shapeID, uv.X(), uv.Y() );
296             nbCreatedNodes++;
297           }
298         }
299         if ( !(newNodes[i] = n2nIt->second ))
300           break;
301       }
302       if ( !newNodes.back() )
303         continue; // not all nodes of the face lie on theShape
304
305       // try to find already created face
306       SMDS_MeshElement * newFace = 0;
307       if ( nbCreatedNodes == 0 &&
308            tgtMesh->FindElement(newNodes, SMDSAbs_Face, /*noMedium=*/false))
309         continue; // repeated face in source groups already created 
310
311       // check future face orientation
312       const int nbCorners = face->NbCornerNodes();
313       const bool isQuad   = ( nbCorners != (int) newNodes.size() );
314       if ( toCheckOri )
315       {
316         int iNode = -1;
317         gp_Vec geomNorm;
318         do
319         {
320           uv = helper.GetNodeUV( geomFace, newNodes[++iNode] );
321           surface->D1( uv.X(),uv.Y(), p, du,dv );
322           geomNorm = reverse ? dv^du : du^dv;
323         }
324         while ( geomNorm.SquareMagnitude() < 1e-6 && iNode+1 < nbCorners );
325
326         int iNext = helper.WrapIndex( iNode+1, nbCorners );
327         int iPrev = helper.WrapIndex( iNode-1, nbCorners );
328
329         SMESH_TNodeXYZ prevNode( newNodes[iPrev] );
330         SMESH_TNodeXYZ curNode ( newNodes[iNode] );
331         SMESH_TNodeXYZ nextNode( newNodes[iNext] );
332         gp_Vec n1n0( prevNode - curNode);
333         gp_Vec n1n2( nextNode - curNode );
334         gp_Vec meshNorm = n1n2 ^ n1n0;
335
336         if ( geomNorm * meshNorm < 0 )
337           SMDS_MeshCell::applyInterlace
338             ( SMDS_MeshCell::reverseSmdsOrder( face->GetEntityType() ), newNodes );
339       }
340
341       // make a new face
342       if ( face->IsPoly() )
343         newFace = tgtMesh->AddPolygonalFace( newNodes );
344       else
345         switch ( newNodes.size() )
346         {
347         case 3:
348           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2] );
349           break;
350         case 4:
351           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3] );
352           break;
353         case 6:
354           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2],
355                                       newNodes[3], newNodes[4], newNodes[5]);
356           break;
357         case 8:
358           newFace = tgtMesh->AddFace( newNodes[0], newNodes[1], newNodes[2], newNodes[3],
359                                       newNodes[4], newNodes[5], newNodes[6], newNodes[7]);
360           break;
361         default: continue;
362         }
363       tgtMesh->SetMeshElementOnShape( newFace, shapeID );
364       e2e->insert( make_pair( face, newFace ));
365
366       // collect links
367       const SMDS_MeshNode* medium = 0;
368       for ( int i = 0; i < nbCorners; ++i )
369       {
370         const SMDS_MeshNode* n1 = newNodes[i];
371         const SMDS_MeshNode* n2 = newNodes[ (i+1)%nbCorners ];
372         if ( isQuad ) // quadratic face
373           medium = newNodes[i+nbCorners];
374         link2Nb = linkCount.insert( make_pair( TLink( n1, n2, medium ), 0)).first;
375         ++link2Nb->second;
376         // if ( link2Nb->second == 1 )
377         // {
378         //   // measure link length
379         //   double len2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
380         //   if ( len2 < minGroupTol )
381         //     minGroupTol = len2;
382         // }
383       }
384     }
385     helper.GetMeshDS()->RemoveNode(tmpNode);
386   }
387
388   // ==========================================================
389   // Put nodes on geom edges and create edges on them;
390   // check if the whole geom face is covered by imported faces
391   // ==========================================================
392
393   // use large tolerance for projection of nodes to edges because of
394   // BLSURF mesher specifics (issue 0020918, Study2.hdf)
395   const double projTol = minGroupTol;
396
397   bool isFaceMeshed = false;
398   SMESHDS_SubMesh* tgtFaceSM = tgtMesh->MeshElements( theShape );
399   if ( tgtFaceSM )
400   {
401     // the imported mesh is valid if all external links (encountered once)
402     // lie on geom edges
403     subShapeIDs.erase( shapeID ); // to contain edges and vertices only
404     double u, f, l;
405     for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); ++link2Nb)
406     {
407       const TLink& link = (*link2Nb).first;
408       int nbFaces = link2Nb->second;
409       if ( nbFaces == 1 )
410       {
411         // check if a not shared link lies on face boundary
412         bool nodesOnBoundary = true;
413         list< TopoDS_Shape > bndShapes;
414         for ( int is1stN = 0; is1stN < 2 && nodesOnBoundary; ++is1stN )
415         {
416           const SMDS_MeshNode* n = is1stN ? link.node1() : link.node2();
417           if ( !subShapeIDs.count( n->getshapeId() )) // n is assigned to FACE
418           {
419             for ( size_t iE = 0; iE < edges.size(); ++iE )
420               if ( helper.CheckNodeU( edges[iE], n, u=0, projTol, /*force=*/true ))
421               {
422                 BRep_Tool::Range(edges[iE],f,l);
423                 if ( Abs(u-f) < 2 * faceTol || Abs(u-l) < 2 * faceTol )
424                   // duplicated node on vertex
425                   return error("Source elements overlap one another");
426                 tgtFaceSM->RemoveNode( n, /*isNodeDeleted=*/false );
427                 tgtMesh->SetNodeOnEdge( n, edges[iE], u );
428                 break;
429               }
430             nodesOnBoundary = subShapeIDs.count( n->getshapeId());
431           }
432           if ( nodesOnBoundary )
433           {
434             TopoDS_Shape s = helper.GetSubShapeByNode( n, tgtMesh );
435             if ( s.ShapeType() == TopAbs_VERTEX )
436               bndShapes.push_front( s ); // vertex first
437             else
438               bndShapes.push_back( s ); // edges last
439           }
440         }
441         if ( !nodesOnBoundary )
442         {
443           error("free internal link"); // just for an easier debug
444           break;
445         }
446         if ( bndShapes.front().ShapeType() == TopAbs_EDGE && // all link nodes are on EDGEs
447              bndShapes.front() != bndShapes.back() )
448           // link nodes on different geom edges
449           return error(COMPERR_BAD_INPUT_MESH, "Source nodes mismatch target vertices");
450
451         // find geom edge the link is on
452         if ( bndShapes.back().ShapeType() != TopAbs_EDGE ) // all link nodes are on VERTEXes
453         {
454           // find geom edge by two vertices
455           TopoDS_Shape geomEdge = helper.GetCommonAncestor( bndShapes.back(),
456                                                             bndShapes.front(),
457                                                             theMesh, TopAbs_EDGE );
458           if ( geomEdge.IsNull() )
459           {
460             error("free internal link");
461             break; // vertices belong to different edges
462           }
463           bndShapes.push_back( geomEdge );
464         }
465
466         // create an edge if not yet exists
467         newNodes.resize(2);
468         newNodes[0] = link.node1(), newNodes[1] = link.node2();
469         const SMDS_MeshElement* edge = tgtMesh->FindElement( newNodes, SMDSAbs_Edge );
470         if ( edge ) continue;
471
472         if ( link._reversed ) std::swap( newNodes[0], newNodes[1] );
473         if ( link._medium )
474         {
475           edge = tgtMesh->AddEdge( newNodes[0], newNodes[1], link._medium );
476
477           TopoDS_Edge geomEdge = TopoDS::Edge(bndShapes.back());
478           helper.CheckNodeU( geomEdge, link._medium, u, projTol, /*force=*/true );
479           tgtFaceSM->RemoveNode( link._medium, /*isNodeDeleted=*/false );
480           tgtMesh->SetNodeOnEdge( (SMDS_MeshNode*)link._medium, geomEdge, u );
481         }
482         else
483         {
484           edge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
485         }
486         if ( !edge )
487           return false;
488
489         tgtMesh->SetMeshElementOnShape( edge, bndShapes.back() );
490       }
491       else if ( nbFaces > 2 )
492       {
493         return error( COMPERR_BAD_INPUT_MESH, "Non-manifold source mesh");
494       }
495     }
496     isFaceMeshed = ( link2Nb == linkCount.end() && !linkCount.empty());
497     if ( isFaceMeshed )
498     {
499       // check that source faces do not overlap:
500       // there must be only two edges sharing each vertex and bound to sub-edges of theShape
501       SMESH_MeshEditor editor( &theMesh );
502       set<int>::iterator subID = subShapeIDs.begin();
503       for ( ; subID != subShapeIDs.end(); ++subID )
504       {
505         const TopoDS_Shape& s = tgtMesh->IndexToShape( *subID );
506         if ( s.ShapeType() != TopAbs_VERTEX ) continue;
507         const SMDS_MeshNode* n = SMESH_Algo::VertexNode( TopoDS::Vertex(s), tgtMesh );
508         SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(SMDSAbs_Edge);
509         int nbEdges = 0;
510         while ( eIt->more() )
511         {
512           const SMDS_MeshElement* edge = eIt->next();
513           int sId = editor.FindShape( edge );
514           nbEdges += subShapeIDs.count( sId );
515         }
516         if ( nbEdges < 2 && !helper.IsRealSeam( s ))
517           return false; // weird
518         if ( nbEdges > 2 )
519           return error( COMPERR_BAD_INPUT_MESH, "Source elements overlap one another");
520       }
521     }
522   }
523   if ( !isFaceMeshed )
524     return error( COMPERR_BAD_INPUT_MESH,
525                   "Source elements don't cover totally the geometrical face" );
526
527   if ( helper.HasSeam() )
528   {
529     // links on seam edges are shared by two faces, so no edges were created on them
530     // by the previous detection of 2D mesh boundary
531     for ( size_t iE = 0; iE < edges.size(); ++iE )
532     {
533       if ( !helper.IsRealSeam( edges[iE] )) continue;
534       const TopoDS_Edge& seamEdge = edges[iE];
535       // to find nodes lying on the seamEdge we check nodes of mesh faces sharing a node on one
536       // of its vertices; after finding another node on seamEdge we continue the same way
537       // until finding all nodes.
538       TopoDS_Vertex      seamVertex = helper.IthVertex( 0, seamEdge );
539       const SMDS_MeshNode* vertNode = SMESH_Algo::VertexNode( seamVertex, tgtMesh );
540       set< const SMDS_MeshNode* > checkedNodes; checkedNodes.insert( vertNode );
541       set< const SMDS_MeshElement* > checkedFaces;
542       // as a face can have more than one node on the seamEdge, there is a difficulty in selecting
543       // one of those nodes to treat next; so we simply find all nodes on the seamEdge and
544       // then sort them by U on edge
545       typedef list< pair< double, const SMDS_MeshNode* > > TUNodeList;
546       TUNodeList nodesOnSeam;
547       double u = helper.GetNodeU( seamEdge, vertNode );
548       nodesOnSeam.push_back( make_pair( u, vertNode ));
549       TUNodeList::iterator u2nIt = nodesOnSeam.begin();
550       for ( ; u2nIt != nodesOnSeam.end(); ++u2nIt )
551       {
552         const SMDS_MeshNode* startNode = (*u2nIt).second;
553         SMDS_ElemIteratorPtr faceIt = startNode->GetInverseElementIterator( SMDSAbs_Face );
554         while ( faceIt->more() )
555         {
556           const SMDS_MeshElement* face = faceIt->next();
557           if ( !checkedFaces.insert( face ).second ) continue;
558           for ( int i = 0, nbNodes = face->NbCornerNodes(); i < nbNodes; ++i )
559           {
560             const SMDS_MeshNode* n = face->GetNode( i );
561             if ( n == startNode || !checkedNodes.insert( n ).second ) continue;
562             if ( helper.CheckNodeU( seamEdge, n, u=0, projTol, /*force=*/true ))
563               nodesOnSeam.push_back( make_pair( u, n ));
564           }
565         }
566       }
567       // sort the found nodes by U on the seamEdge; most probably they are in a good order,
568       // so we can use the hint to spead-up map filling
569       map< double, const SMDS_MeshNode* > u2nodeMap;
570       for ( u2nIt = nodesOnSeam.begin(); u2nIt != nodesOnSeam.end(); ++u2nIt )
571         u2nodeMap.insert( u2nodeMap.end(), *u2nIt );
572
573       // create edges
574       {
575         SMESH_MesherHelper seamHelper( theMesh );
576         seamHelper.SetSubShape( edges[ iE ]);
577         seamHelper.SetElementsOnShape( true );
578
579         if ( (*checkedFaces.begin())->IsQuadratic() )
580           for ( set< const SMDS_MeshElement* >::iterator fIt = checkedFaces.begin();
581                 fIt != checkedFaces.end(); ++fIt )
582             seamHelper.AddTLinks( static_cast<const SMDS_MeshFace*>( *fIt ));
583
584         map< double, const SMDS_MeshNode* >::iterator n1, n2, u2nEnd = u2nodeMap.end();
585         for ( n2 = u2nodeMap.begin(), n1 = n2++; n2 != u2nEnd; ++n1, ++n2 )
586         {
587           const SMDS_MeshNode* node1 = n1->second;
588           const SMDS_MeshNode* node2 = n2->second;
589           seamHelper.AddEdge( node1, node2 );
590           if ( node2->getshapeId() == helper.GetSubShapeID() )
591           {
592             tgtFaceSM->RemoveNode( node2, /*isNodeDeleted=*/false );
593             tgtMesh->SetNodeOnEdge( const_cast<SMDS_MeshNode*>( node2 ), seamEdge, n2->first );
594           }
595         }
596       }
597     } // loop on edges to find seam ones
598   } // if ( helper.HasSeam() )
599
600   // notify sub-meshes of edges on computation
601   for ( size_t iE = 0; iE < edges.size(); ++iE )
602   {
603     SMESH_subMesh * sm = theMesh.GetSubMesh( edges[iE] );
604     // if ( SMESH_Algo::isDegenerated( edges[iE] ))
605     //   sm->SetIsAlwaysComputed( true );
606     sm->ComputeStateEngine(SMESH_subMesh::CHECK_COMPUTE_STATE);
607     if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
608       return error(SMESH_Comment("Failed to create segments on the edge ")
609                    << tgtMesh->ShapeToIndex( edges[iE ]));
610   }
611
612   // ============
613   // Copy meshes
614   // ============
615
616   vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
617   for ( size_t i = 0; i < srcMeshes.size(); ++i )
618     StdMeshers_Import_1D::importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
619
620   return true;
621 }
622
623 //=============================================================================
624 /*!
625  * \brief Set needed event listeners and create a submesh for a copied mesh
626  *
627  * This method is called only if a submesh has HYP_OK algo_state.
628  */
629 //=============================================================================
630
631 void StdMeshers_Import_1D2D::SetEventListener(SMESH_subMesh* subMesh)
632 {
633   if ( !_sourceHyp )
634   {
635     const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
636     SMESH_Mesh*         tgtMesh  = subMesh->GetFather();
637     Hypothesis_Status aStatus;
638     CheckHypothesis( *tgtMesh, tgtShape, aStatus );
639   }
640   StdMeshers_Import_1D::setEventListener( subMesh, _sourceHyp );
641 }
642 void StdMeshers_Import_1D2D::SubmeshRestored(SMESH_subMesh* subMesh)
643 {
644   SetEventListener(subMesh);
645 }
646
647 //=============================================================================
648 /*!
649  * Predict nb of mesh entities created by Compute()
650  */
651 //=============================================================================
652
653 bool StdMeshers_Import_1D2D::Evaluate(SMESH_Mesh &         theMesh,
654                                       const TopoDS_Shape & theShape,
655                                       MapShapeNbElems&     aResMap)
656 {
657   if ( !_sourceHyp ) return false;
658
659   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
660   if ( srcGroups.empty() )
661     return error("Invalid source groups");
662
663   vector<int> aVec(SMDSEntity_Last,0);
664
665   bool toCopyMesh, toCopyGroups;
666   _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
667   if ( toCopyMesh ) // the whole mesh is copied
668   {
669     vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
670     for ( unsigned i = 0; i < srcMeshes.size(); ++i )
671     {
672       SMESH_subMesh* sm = StdMeshers_Import_1D::getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
673       if ( !sm || aResMap.count( sm )) continue; // already counted
674       const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
675       for (int i = 0; i < SMDSEntity_Last; i++)
676         aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
677     }
678   }
679   else
680   {
681     // std-like iterator used to get coordinates of nodes of mesh element
682     typedef SMDS_StdIterator< SMESH_TNodeXYZ, SMDS_ElemIteratorPtr > TXyzIterator;
683
684     SMESH_MesherHelper helper(theMesh);
685     helper.SetSubShape(theShape);
686
687     const TopoDS_Face& geomFace = TopoDS::Face( theShape );
688
689     // take into account nodes on vertices
690     TopExp_Explorer exp( theShape, TopAbs_VERTEX );
691     for ( ; exp.More(); exp.Next() )
692       theMesh.GetSubMesh( exp.Current())->Evaluate( aResMap );
693
694     // to count now many times a link between nodes encounters,
695     // negative nb additionally means that a link is quadratic
696     map<SMESH_TLink, int> linkCount;
697     map<SMESH_TLink, int>::iterator link2Nb;
698
699     // count faces and nodes imported from groups
700     set<const SMDS_MeshNode* > allNodes;
701     gp_XY uv;
702     double minGroupTol = 1e100;
703     for ( int iG = 0; iG < srcGroups.size(); ++iG )
704     {
705       const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
706       const double groupTol = 0.5 * sqrt( getMinElemSize2( srcGroup ));
707       minGroupTol = std::min( groupTol, minGroupTol );
708       SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
709       SMDS_MeshNode *tmpNode =helper.AddNode(0,0,0);
710       while ( srcElems->more() ) // loop on group contents
711       {
712         const SMDS_MeshElement* face = srcElems->next();
713         // find out if face is located on geomEdge by projecting
714         // a gravity center of face to geomFace
715         gp_XYZ gc(0,0,0);
716         gc = accumulate( TXyzIterator(face->nodesIterator()), TXyzIterator(), gc)/face->NbNodes();
717         tmpNode->setXYZ( gc.X(), gc.Y(), gc.Z());
718         if ( helper.CheckNodeUV( geomFace, tmpNode, uv, groupTol, /*force=*/true ))
719         {
720           ++aVec[ face->GetEntityType() ];
721
722           // collect links
723           int nbConers = face->NbCornerNodes();
724           for ( int i = 0; i < face->NbNodes(); ++i )
725           {
726             const SMDS_MeshNode* n1 = face->GetNode(i);
727             allNodes.insert( n1 );
728             if ( i < nbConers )
729             {
730               const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbConers );
731               link2Nb = linkCount.insert( make_pair( SMESH_TLink( n1, n2 ), 0)).first;
732               if ( (*link2Nb).second )
733                 link2Nb->second += (link2Nb->second < 0 ) ? -1 : 1;
734               else
735                 link2Nb->second += ( face->IsQuadratic() ) ? -1 : 1;
736             }
737           }
738         }
739       }
740       helper.GetMeshDS()->RemoveNode(tmpNode);
741     }
742
743     int nbNodes = allNodes.size();
744     allNodes.clear();
745
746     // count nodes and edges on geom edges
747
748     double u;
749     for ( exp.Init(theShape, TopAbs_EDGE); exp.More(); exp.Next() )
750     {
751       TopoDS_Edge geomEdge = TopoDS::Edge( exp.Current() );
752       SMESH_subMesh* sm = theMesh.GetSubMesh( geomEdge );
753       vector<int>& edgeVec = aResMap[sm];
754       if ( edgeVec.empty() )
755       {
756         edgeVec.resize(SMDSEntity_Last,0);
757         for ( link2Nb = linkCount.begin(); link2Nb != linkCount.end(); )
758         {
759           const SMESH_TLink& link = (*link2Nb).first;
760           int nbFacesOfLink = Abs( link2Nb->second );
761           bool eraseLink = ( nbFacesOfLink != 1 );
762           if ( nbFacesOfLink == 1 )
763           {
764             if ( helper.CheckNodeU( geomEdge, link.node1(), u, minGroupTol, /*force=*/true )&&
765                  helper.CheckNodeU( geomEdge, link.node2(), u, minGroupTol, /*force=*/true ))
766             {
767               bool isQuadratic = ( link2Nb->second < 0 );
768               ++edgeVec[ isQuadratic ? SMDSEntity_Quad_Edge : SMDSEntity_Edge ];
769               ++edgeVec[ SMDSEntity_Node ];
770               --nbNodes;
771               eraseLink = true;
772             }
773           }
774           if ( eraseLink )
775             linkCount.erase(link2Nb++);
776           else
777             link2Nb++;
778         }
779         if ( edgeVec[ SMDSEntity_Node] > 0 )
780           --edgeVec[ SMDSEntity_Node ]; // for one node on vertex
781       }
782       else if ( !helper.IsSeamShape( geomEdge ) ||
783                 geomEdge.Orientation() == TopAbs_FORWARD )
784       {
785         nbNodes -= 1+edgeVec[ SMDSEntity_Node ];
786       }
787     }
788
789     aVec[SMDSEntity_Node] = nbNodes;
790   }
791
792   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
793   aResMap.insert(make_pair(sm,aVec));
794
795   return true;
796 }