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