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