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