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