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