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