Salome HOME
bos #20087: [CEA 20080] spacing value constraint not taken into account in Body fitti...
[modules/smesh.git] / src / StdMeshers / StdMeshers_Import_1D.cxx
1 // Copyright (C) 2007-2020  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_1D.cxx
25 //  Module : SMESH
26 //
27 #include "StdMeshers_Import_1D.hxx"
28 #include "StdMeshers_ImportSource.hxx"
29
30 #include "SMDS_MeshElement.hxx"
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESHDS_Group.hxx"
33 #include "SMESHDS_Mesh.hxx"
34 #include "SMESH_Comment.hxx"
35 #include "SMESH_Gen.hxx"
36 #include "SMESH_Group.hxx"
37 #include "SMESH_HypoFilter.hxx"
38 #include "SMESH_Mesh.hxx"
39 #include "SMESH_MeshEditor.hxx"
40 #include "SMESH_MesherHelper.hxx"
41 #include "SMESH_Octree.hxx"
42 #include "SMESH_subMesh.hxx"
43 #include "SMESH_subMeshEventListener.hxx"
44
45 #include "Utils_SALOME_Exception.hxx"
46 #include "utilities.h"
47
48 #include <BRepAdaptor_Curve.hxx>
49 #include <BRep_Builder.hxx>
50 #include <BRep_Tool.hxx>
51 #include <BndLib_Add3dCurve.hxx>
52 #include <GCPnts_TangentialDeflection.hxx>
53 #include <ShapeAnalysis_Curve.hxx>
54 #include <TopExp.hxx>
55 #include <TopExp_Explorer.hxx>
56 #include <TopoDS.hxx>
57 #include <TopoDS_Compound.hxx>
58 #include <TopoDS_Edge.hxx>
59 #include <TopoDS_Vertex.hxx>
60
61 using namespace std;
62
63 //================================================================================
64 namespace // INTERNAL STUFF
65 //================================================================================
66 {
67   /*!
68    * \brief Compute point position on a curve. Use octree to fast reject far points
69    */
70   class CurveProjector : public SMESH_Octree
71   {
72   public:
73     CurveProjector( const TopoDS_Edge& edge, double enlarge );
74
75     bool IsOnCurve( const gp_XYZ& point, double & distance2, double & u );
76
77     bool IsOut( const gp_XYZ& point ) const { return getBox()->IsOut( point ); }
78
79   protected:
80     CurveProjector() {}
81     SMESH_Octree* newChild() const { return new CurveProjector; }
82     void          buildChildrenData();
83     Bnd_B3d*      buildRootBox();
84
85   private:
86     struct CurveSegment : public Bnd_B3d
87     {
88       double _chord, _chord2, _length2;
89       gp_Pnt _pFirst, _pLast;
90       gp_Lin _line;
91       Handle(Geom_Curve) _curve;
92
93       CurveSegment() {}
94       void Init( const gp_Pnt& pf, const gp_Pnt& pl,
95                  double uf, double ul, double tol, Handle(Geom_Curve)& curve );
96       bool IsOn( const gp_XYZ& point, double & distance2, double & u );
97       bool IsInContact( const Bnd_B3d& bb );
98     };
99     std::vector< CurveSegment > _segments;
100   };
101
102   //===============================================================================
103   /*!
104    * \brief Create an octree of curve segments
105    */
106   //================================================================================
107
108   CurveProjector::CurveProjector( const TopoDS_Edge& edge, double enlarge )
109     :SMESH_Octree( 0 )
110   {
111     double f,l;
112     Handle(Geom_Curve) curve = BRep_Tool::Curve( edge, f, l );
113     double curDeflect = 0.3; // Curvature deflection
114     double angDeflect = 1e+100; // Angular deflection - don't control chordal error
115     GCPnts_TangentialDeflection div( BRepAdaptor_Curve( edge ), angDeflect, curDeflect );
116     _segments.resize( div.NbPoints() - 1 );
117     for ( int i = 1; i < div.NbPoints(); ++i )
118       try {
119         _segments[ i - 1 ].Init( div.Value( i ),     div.Value( i+1 ),
120                                  div.Parameter( i ), div.Parameter( i+1 ),
121                                  enlarge, curve );
122       }
123       catch ( Standard_Failure ) {
124         _segments.resize( _segments.size() - 1 );
125         --i;
126       }
127     if ( _segments.size() < 3 )
128       myIsLeaf = true;
129
130     compute();
131
132     if ( _segments.size() == 1 )
133       myBox->Enlarge( enlarge );
134   }
135
136   //================================================================================
137   /*!
138    * \brief Return the maximal box
139    */
140   //================================================================================
141
142   Bnd_B3d* CurveProjector::buildRootBox()
143   {
144     Bnd_B3d* box = new Bnd_B3d;
145     for ( size_t i = 0; i < _segments.size(); ++i )
146       box->Add( _segments[i] );
147     return box;
148   }
149
150   //================================================================================
151   /*!
152    * \brief Redistribute segments among children
153    */
154   //================================================================================
155
156   void CurveProjector::buildChildrenData()
157   {
158     bool allIn = true;
159     for ( size_t i = 0; i < _segments.size(); ++i )
160     {
161       for (int j = 0; j < 8; j++)
162       {
163         if ( _segments[i].IsInContact( *myChildren[j]->getBox() ))
164           ((CurveProjector*)myChildren[j])->_segments.push_back( _segments[i]);
165         else
166           allIn = false;
167       }
168     }
169     if ( allIn && _segments.size() < 3 )
170     {
171       myIsLeaf = true;
172       for (int j = 0; j < 8; j++)
173         static_cast<CurveProjector*>( myChildren[j])->myIsLeaf = true;
174     }
175     else
176     {
177       SMESHUtils::FreeVector( _segments ); // = _segments.clear() + free memory
178
179       for (int j = 0; j < 8; j++)
180       {
181         CurveProjector* child = static_cast<CurveProjector*>( myChildren[j]);
182         if ( child->_segments.size() < 3 )
183           child->myIsLeaf = true;
184       }
185     }
186   }
187
188   //================================================================================
189   /*!
190    * \brief Return true if a point is close to the curve
191    *  \param [in] point - the point
192    *  \param [out] distance2 - distance to the curve
193    *  \param [out] u - parameter on the curve
194    *  \return bool - is the point is close to the curve
195    */
196   //================================================================================
197
198   bool CurveProjector::IsOnCurve( const gp_XYZ& point, double & distance2, double & u )
199   {
200     if ( getBox()->IsOut( point ))
201       return false;
202
203     if ( isLeaf() )
204     {
205       for ( size_t i = 0; i < _segments.size(); ++i )
206         if ( !_segments[i].IsOut( point ) &&
207              _segments[i].IsOn( point, distance2, u ))
208           return true;
209     }
210     else
211     {
212       for (int i = 0; i < 8; i++)
213         if (((CurveProjector*) myChildren[i])->IsOnCurve( point, distance2, u ))
214           return true;
215     }
216     return false;
217   }
218   
219   //================================================================================
220   /*!
221    * \brief Initialize
222    */
223   //================================================================================
224
225   void CurveProjector::CurveSegment::Init(const gp_Pnt&       pf,
226                                           const gp_Pnt&       pl,
227                                           const double        uf,
228                                           const double        ul,
229                                           const double        tol,
230                                           Handle(Geom_Curve)& curve )
231   {
232     _pFirst  = pf;
233     _pLast   = pl;
234     _curve   = curve;
235     _length2 = pf.SquareDistance( pl );
236     _chord2  = Max( _line.     SquareDistance( curve->Value( uf + 0.25 * ( ul - uf ))),
237                     Max( _line.SquareDistance( curve->Value( uf + 0.5  * ( ul - uf ))),
238                          _line.SquareDistance( curve->Value( uf + 0.75 * ( ul - uf )))));
239     _chord2  = Max( tol, _chord2 );
240     _chord   = Sqrt( _chord2 );
241     _line.SetLocation( pf );
242     _line.SetDirection( gp_Vec( pf, pl ));
243
244     Bnd_Box bb;
245     BndLib_Add3dCurve::Add( GeomAdaptor_Curve( curve, uf, ul ), tol, bb );
246     Add( bb.CornerMin() );
247     Add( bb.CornerMax() );
248   }
249
250   //================================================================================
251   /*!
252    * \brief Return true if a point is close to the curve segment
253    *  \param [in] point - the point
254    *  \param [out] distance2 - distance to the curve
255    *  \param [out] u - parameter on the curve
256    *  \return bool - is the point is close to the curve segment
257    */
258   //================================================================================
259
260   bool CurveProjector::CurveSegment::IsOn( const gp_XYZ& point, double & distance2, double & u )
261   {
262     distance2 = _line.SquareDistance( point );
263     if ( distance2 > _chord2 )
264       return false;
265
266     // check if the point projection falls into the segment range
267     {
268       gp_Vec edge( _pFirst, _pLast );
269       gp_Vec n1p ( _pFirst, point  );
270       u = ( edge * n1p ) / _length2; // param [0,1] on the edge
271       if ( u < 0 )
272       {
273         if ( _pFirst.SquareDistance( point ) > _chord2 )
274           return false;
275       }
276       else if ( u > _chord )
277       {
278         if ( _pLast.SquareDistance( point ) > _chord2 )
279           return false;
280       }
281     }
282     gp_Pnt proj;
283     distance2 = ShapeAnalysis_Curve().Project( _curve, point, Precision::Confusion(),
284                                                proj, u, false );
285     distance2 *= distance2;
286     return true;
287   }
288
289   //================================================================================
290   /*!
291    * \brief Check if the segment is in contact with a box
292    */
293   //================================================================================
294
295   bool CurveProjector::CurveSegment::IsInContact( const Bnd_B3d& bb )
296   {
297     if ( bb.IsOut( _line.Position(), /*isRay=*/true, _chord ))
298       return false;
299
300     gp_Ax1 axRev = _line.Position().Reversed();
301     axRev.SetLocation( _pLast );
302     return !bb.IsOut( axRev, /*isRay=*/true, _chord );
303   }
304
305   //================================================================================
306   //================================================================================
307
308   int getSubmeshIDForCopiedMesh(const SMESHDS_Mesh* srcMeshDS, SMESH_Mesh* tgtMesh);
309
310   enum _ListenerDataType
311     {
312       WAIT_HYP_MODIF=1, // data indicating awaiting for valid parameters of src hyp
313       LISTEN_SRC_MESH, // data storing submesh depending on source mesh state
314       SRC_HYP // data storing ImportSource hyp
315     };
316   //================================================================================
317   /*!
318    * \brief _ListenerData holding ImportSource hyp holding in its turn
319    *  imported groups
320    */
321   struct _ListenerData : public SMESH_subMeshEventListenerData
322   {
323     const StdMeshers_ImportSource1D* _srcHyp;
324     _ListenerData(const StdMeshers_ImportSource1D* h, _ListenerDataType type=SRC_HYP):
325       SMESH_subMeshEventListenerData(/*isDeletable=*/true), _srcHyp(h)
326     {
327       myType = type;
328     }
329   };
330   //================================================================================
331   /*!
332    * \brief Comparator of sub-meshes
333    */
334   struct _SubLess
335   {
336     bool operator()(const SMESH_subMesh* sm1, const SMESH_subMesh* sm2 ) const
337     {
338       if ( sm1 == sm2 ) return false;
339       if ( !sm1 || !sm2 ) return sm1 < sm2;
340       const TopoDS_Shape& s1 = sm1->GetSubShape();
341       const TopoDS_Shape& s2 = sm2->GetSubShape();
342       TopAbs_ShapeEnum t1 = s1.IsNull() ? TopAbs_SHAPE : s1.ShapeType();
343       TopAbs_ShapeEnum t2 = s2.IsNull() ? TopAbs_SHAPE : s2.ShapeType();
344       if ( t1 == t2)
345         return (sm1 < sm2);
346       return t1 < t2; // to have: face < edge
347     }
348   };
349   //================================================================================
350   /*!
351    * \brief Container of data dedicated to one source mesh
352    */
353   struct _ImportData
354   {
355     const SMESH_Mesh* _srcMesh;
356     StdMeshers_Import_1D::TNodeNodeMap _n2n;
357     StdMeshers_Import_1D::TElemElemMap _e2e;
358
359     set< SMESH_subMesh*, _SubLess > _subM; // submeshes relating to this srcMesh
360     set< SMESH_subMesh*, _SubLess > _copyMeshSubM; // submeshes requesting mesh copying
361     set< SMESH_subMesh*, _SubLess > _copyGroupSubM; // submeshes requesting group copying
362     set< SMESH_subMesh*, _SubLess > _computedSubM;
363
364     SMESHDS_SubMesh*     _importMeshSubDS; // submesh storing a copy of _srcMesh
365     int                  _importMeshSubID; // id of _importMeshSubDS
366
367     _ImportData(const SMESH_Mesh* srcMesh=0):
368       _srcMesh(srcMesh), _importMeshSubDS(0),_importMeshSubID(-1) {}
369
370     void removeImportedMesh( SMESHDS_Mesh* meshDS )
371     {
372       if ( !_importMeshSubDS ) return;
373       SMDS_ElemIteratorPtr eIt = _importMeshSubDS->GetElements();
374       while ( eIt->more() )
375         meshDS->RemoveFreeElement( eIt->next(), 0, /*fromGroups=*/false );
376       SMDS_NodeIteratorPtr nIt = _importMeshSubDS->GetNodes();
377       while ( nIt->more() )
378         meshDS->RemoveFreeNode( nIt->next(), 0, /*fromGroups=*/false );
379       _importMeshSubDS->Clear();
380       _n2n.clear();
381       _e2e.clear();
382     }
383     void removeGroups( SMESH_subMesh* subM, const StdMeshers_ImportSource1D* srcHyp )
384     {
385       if ( !srcHyp ) return;
386       SMESH_Mesh*           tgtMesh = subM->GetFather();
387       const SMESHDS_Mesh* tgtMeshDS = tgtMesh->GetMeshDS();
388       const SMESHDS_Mesh* srcMeshDS = _srcMesh->GetMeshDS();
389       vector<SMESH_Group*>*  groups =
390         const_cast<StdMeshers_ImportSource1D*>(srcHyp)->GetResultGroups(*srcMeshDS,*tgtMeshDS);
391       if ( groups )
392       {
393         for ( unsigned i = 0; i < groups->size(); ++i )
394           tgtMesh->RemoveGroup( groups->at(i)->GetGroupDS()->GetID() );
395         groups->clear();
396       }
397     }
398     void trackHypParams( SMESH_subMesh* sm, const StdMeshers_ImportSource1D* srcHyp )
399     {
400       if ( !srcHyp ) return;
401       bool toCopyMesh, toCopyGroups;
402       srcHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
403
404       if ( toCopyMesh )_copyMeshSubM.insert( sm );
405       else             _copyMeshSubM.erase( sm );
406
407       if ( toCopyGroups ) _copyGroupSubM.insert( sm );
408       else                _copyGroupSubM.erase( sm );
409     }
410     void addComputed( SMESH_subMesh* sm )
411     {
412       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
413                                                                /*complexShapeFirst=*/true);
414       while ( smIt->more() )
415       {
416         sm = smIt->next();
417         switch ( sm->GetSubShape().ShapeType() )
418         {
419         case TopAbs_EDGE:
420           if ( SMESH_Algo::isDegenerated( TopoDS::Edge( sm->GetSubShape() )))
421             continue;
422         case TopAbs_FACE:
423           _subM.insert( sm );
424           if ( !sm->IsEmpty() )
425             _computedSubM.insert( sm );
426         case TopAbs_VERTEX:
427           break;
428         default:;
429         }
430       }
431     }
432   };
433   //================================================================================
434   /*!
435    * Listener notified on events relating to imported submesh
436    */
437   class _Listener : public SMESH_subMeshEventListener
438   {
439     typedef map< SMESH_Mesh*, list< _ImportData > > TMesh2ImpData;
440     TMesh2ImpData _tgtMesh2ImportData;
441
442     _Listener():SMESH_subMeshEventListener(/*isDeletable=*/false,
443                                            "StdMeshers_Import_1D::_Listener") {}
444
445   public:
446     // return pointer to a static listener
447     static _Listener* get() { static _Listener theListener; return &theListener; }
448
449     static _ImportData* getImportData(const SMESH_Mesh* srcMesh, SMESH_Mesh* tgtMesh);
450
451     static void storeImportSubmesh(SMESH_subMesh*                   importSub,
452                                    const SMESH_Mesh*                srcMesh,
453                                    const StdMeshers_ImportSource1D* srcHyp);
454
455     virtual void ProcessEvent(const int                       event,
456                               const int                       eventType,
457                               SMESH_subMesh*                  subMesh,
458                               SMESH_subMeshEventListenerData* data,
459                               const SMESH_Hypothesis*         hyp);
460     void removeSubmesh( SMESH_subMesh* sm, _ListenerData* data );
461     void clearSubmesh ( SMESH_subMesh* sm, _ListenerData* data, bool clearAllSub );
462     void clearN2N     ( SMESH_Mesh* tgtMesh );
463
464     // mark sm as missing src hyp with valid groups
465     static void waitHypModification(SMESH_subMesh* sm)
466     {
467       sm->SetEventListener
468         (get(), SMESH_subMeshEventListenerData::MakeData( sm, WAIT_HYP_MODIF ), sm);
469     }
470   };
471   //--------------------------------------------------------------------------------
472   /*!
473    * \brief Find or create ImportData for given meshes
474    */
475   _ImportData* _Listener::getImportData(const SMESH_Mesh* srcMesh,
476                                         SMESH_Mesh*       tgtMesh)
477   {
478     list< _ImportData >& dList = get()->_tgtMesh2ImportData[tgtMesh];
479     list< _ImportData >::iterator d = dList.begin();
480     for ( ; d != dList.end(); ++d )
481       if ( d->_srcMesh == srcMesh )
482         return &*d;
483     dList.push_back(_ImportData(srcMesh));
484     return &dList.back();
485   }
486
487   //--------------------------------------------------------------------------------
488   /*!
489    * \brief Remember an imported sub-mesh and set needed even listeners
490    *  \param importSub - submesh computed by Import algo
491    *  \param srcMesh - source mesh
492    *  \param srcHyp - ImportSource hypothesis
493    */
494   void _Listener::storeImportSubmesh(SMESH_subMesh*                   importSub,
495                                      const SMESH_Mesh*                srcMesh,
496                                      const StdMeshers_ImportSource1D* srcHyp)
497   {
498     // set listener to hear events of the submesh computed by "Import" algo
499     importSub->SetEventListener( get(), new _ListenerData(srcHyp), importSub );
500
501     // set listeners to hear events of the source mesh
502     SMESH_subMesh* smToNotify = importSub;
503     vector<SMESH_subMesh*> smToListen = srcHyp->GetSourceSubMeshes( srcMesh );
504     for ( size_t i = 0; i < smToListen.size(); ++i )
505     {
506       SMESH_subMeshEventListenerData* data = new _ListenerData(srcHyp, LISTEN_SRC_MESH);
507       data->mySubMeshes.push_back( smToNotify );
508       importSub->SetEventListener( get(), data, smToListen[i] );
509     }
510     // remember the submesh importSub and its sub-submeshes
511     _ImportData* iData = _Listener::getImportData( srcMesh, importSub->GetFather());
512     iData->trackHypParams( importSub, srcHyp );
513     iData->addComputed( importSub );
514     if ( !iData->_copyMeshSubM.empty() && iData->_importMeshSubID < 1 )
515     {
516       SMESH_Mesh* tgtMesh = importSub->GetFather();
517       iData->_importMeshSubID = getSubmeshIDForCopiedMesh( srcMesh->GetMeshDS(),tgtMesh);
518       iData->_importMeshSubDS = tgtMesh->GetMeshDS()->NewSubMesh( iData->_importMeshSubID );
519     }
520   }
521   //--------------------------------------------------------------------------------
522   /*!
523    * \brief Remove imported mesh and/or groups if needed
524    *  \param sm - submesh losing Import algo
525    *  \param data - data holding imported groups
526    */
527   void _Listener::removeSubmesh( SMESH_subMesh* sm, _ListenerData* data )
528   {
529     list< _ImportData > &  dList = _tgtMesh2ImportData[ sm->GetFather() ];
530     list< _ImportData >::iterator d = dList.begin();
531     for ( ; d != dList.end(); ++d )
532       if ( (*d)._subM.erase( sm ))
533       {
534         d->_computedSubM.erase( sm );
535         bool rmMesh   = d->_copyMeshSubM.erase( sm ) && d->_copyMeshSubM.empty();
536         bool rmGroups = (d->_copyGroupSubM.erase( sm ) && d->_copyGroupSubM.empty()) || rmMesh;
537         if ( rmMesh )
538           d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
539         if ( rmGroups && data && data->myType == SRC_HYP )
540           d->removeGroups( sm, data->_srcHyp );
541       }
542   }
543   //--------------------------------------------------------------------------------
544   /*!
545    * \brief Clear _ImportData::_n2n.
546    *        _n2n is useful within one mesh.Compute() only
547    */
548   void _Listener::clearN2N( SMESH_Mesh* tgtMesh )
549   {
550     list< _ImportData >& dList = get()->_tgtMesh2ImportData[tgtMesh];
551     list< _ImportData >::iterator d = dList.begin();
552     for ( ; d != dList.end(); ++d )
553       d->_n2n.clear();
554   }
555   //--------------------------------------------------------------------------------
556   /*!
557    * \brief Clear submeshes and remove imported mesh and/or groups if necessary
558    *  \param sm - cleared submesh
559    *  \param data - data holding imported groups
560    */
561   void _Listener::clearSubmesh(SMESH_subMesh* sm, _ListenerData* data, bool clearAllSub)
562   {
563     list< _ImportData > &  dList = _tgtMesh2ImportData[ sm->GetFather() ];
564     list< _ImportData >::iterator d = dList.begin();
565     for ( ; d != dList.end(); ++d )
566     {
567       if ( !d->_subM.count( sm )) continue;
568       if ( (*d)._computedSubM.erase( sm ) )
569       {
570         bool copyMesh = !d->_copyMeshSubM.empty();
571         if ( copyMesh || clearAllSub )
572         {
573           // remove imported mesh and groups
574           d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
575
576           if ( data && data->myType == SRC_HYP )
577             d->removeGroups( sm, data->_srcHyp );
578
579           // clear the rest submeshes
580           if ( !d->_computedSubM.empty() )
581           {
582             d->_computedSubM.clear();
583             set< SMESH_subMesh*, _SubLess>::iterator sub = d->_subM.begin();
584             for ( ; sub != d->_subM.end(); ++sub )
585             {
586               SMESH_subMesh* subM = *sub;
587               _ListenerData* hypData = (_ListenerData*) subM->GetEventListenerData( get() );
588               if ( hypData && hypData->myType == SRC_HYP )
589                 d->removeGroups( sm, hypData->_srcHyp );
590
591               subM->ComputeStateEngine( SMESH_subMesh::CLEAN );
592               if ( subM->GetSubShape().ShapeType() == TopAbs_FACE )
593                 subM->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
594             }
595           }
596         }
597         sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
598         if ( sm->GetSubShape().ShapeType() == TopAbs_FACE )
599           sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
600       }
601       if ( data && data->myType == SRC_HYP )
602         d->trackHypParams( sm, data->_srcHyp );
603       d->_n2n.clear();
604       d->_e2e.clear();
605     }
606   }
607   //--------------------------------------------------------------------------------
608   /*!
609    * \brief Remove imported mesh and/or groups
610    */
611   void _Listener::ProcessEvent(const int                       event,
612                                const int                       eventType,
613                                SMESH_subMesh*                  subMesh,
614                                SMESH_subMeshEventListenerData* data,
615                                const SMESH_Hypothesis*         /*hyp*/)
616   {
617     if ( data && data->myType == WAIT_HYP_MODIF )
618     {
619       // event of Import submesh
620       if ( SMESH_subMesh::MODIF_HYP  == event &&
621            SMESH_subMesh::ALGO_EVENT == eventType )
622       {
623         // re-call SetEventListener() to take into account valid parameters
624         // of ImportSource hypothesis
625         if ( SMESH_Algo* algo = subMesh->GetAlgo() )
626           algo->SetEventListener( subMesh );
627       }
628     }
629     else if ( data && data->myType == LISTEN_SRC_MESH )
630     {
631       // event of source mesh
632       if ( SMESH_subMesh::COMPUTE_EVENT == eventType )
633       {
634         switch ( event ) {
635         case SMESH_subMesh::CLEAN:
636           // source mesh cleaned -> clean target mesh
637           clearSubmesh( data->mySubMeshes.front(), (_ListenerData*) data, /*all=*/true );
638           break;
639         case SMESH_subMesh::SUBMESH_COMPUTED: {
640           // source mesh computed -> reset FAILED state of Import submeshes to
641           // READY_TO_COMPUTE
642           SMESH_Mesh* srcMesh = subMesh->GetFather();
643           if ( srcMesh->NbEdges() > 0 || srcMesh->NbFaces() > 0 )
644           {
645             SMESH_Mesh* m = data->mySubMeshes.front()->GetFather();
646             if ( SMESH_subMesh* sm1 = m->GetSubMeshContaining(1))
647             {
648               sm1->ComputeStateEngine(SMESH_subMesh::SUBMESH_COMPUTED );
649               sm1->ComputeSubMeshStateEngine( SMESH_subMesh::SUBMESH_COMPUTED );
650             }
651           }
652           break;
653         }
654         default:;
655         }
656       }
657       if ( !data->mySubMeshes.empty() )
658         clearN2N( data->mySubMeshes.front()->GetFather() );
659     }
660     else // event of Import submesh
661     {
662       // find out what happens: import hyp modified or removed
663       bool removeImport = false, modifHyp = false;
664       if ( SMESH_subMesh::ALGO_EVENT == eventType )
665         modifHyp = true;
666       if ( subMesh->GetAlgoState() != SMESH_subMesh::HYP_OK )
667       {
668         removeImport = true;
669       }
670       else if (( SMESH_subMesh::REMOVE_ALGO == event ||
671                  SMESH_subMesh::REMOVE_FATHER_ALGO == event ) &&
672                SMESH_subMesh::ALGO_EVENT == eventType )
673       {
674         SMESH_Algo* algo = subMesh->GetAlgo();
675         removeImport = ( strncmp( "Import", algo->GetName(), 6 ) != 0 );
676       }
677
678       if ( removeImport )
679       {
680         // treate removal of Import algo from subMesh
681         removeSubmesh( subMesh, (_ListenerData*) data );
682       }
683       else if ( modifHyp ||
684                 ( SMESH_subMesh::CLEAN         == event &&
685                   SMESH_subMesh::COMPUTE_EVENT == eventType))
686       {
687         // treate modification of ImportSource hypothesis
688         clearSubmesh( subMesh, (_ListenerData*) data, /*all=*/false );
689       }
690       else if ( SMESH_subMesh::CHECK_COMPUTE_STATE == event &&
691                 SMESH_subMesh::COMPUTE_EVENT       == eventType )
692       {
693         // check compute state of all submeshes impoting from same src mesh;
694         // this is to take into account 1D computed submeshes hidden by 2D import algo;
695         // else source mesh is not copied as _subM.size != _computedSubM.size()
696         list< _ImportData > &  dList = _tgtMesh2ImportData[ subMesh->GetFather() ];
697         list< _ImportData >::iterator d = dList.begin();
698         for ( ; d != dList.end(); ++d )
699           if ( d->_subM.count( subMesh ))
700           {
701             set<SMESH_subMesh*,_SubLess>::iterator smIt = d->_subM.begin();
702             for( ; smIt != d->_subM.end(); ++smIt )
703               if ( (*smIt)->IsMeshComputed() )
704                 d->_computedSubM.insert( *smIt);
705           }
706       }
707       // Clear _ImportData::_n2n if it's no more useful, i.e. when
708       // the event is not within mesh.Compute()
709       if ( SMESH_subMesh::ALGO_EVENT == eventType )
710         clearN2N( subMesh->GetFather() );
711     }
712   }
713
714   //================================================================================
715   /*!
716    * \brief Return an ID of submesh to store nodes and elements of a copied mesh
717    */
718   //================================================================================
719
720   int getSubmeshIDForCopiedMesh(const SMESHDS_Mesh* srcMeshDS,
721                                 SMESH_Mesh*         tgtMesh)
722   {
723     // To get SMESH_subMesh corresponding to srcMeshDS we need to have a shape
724     // for which SMESHDS_Mesh::IsGroupOfSubShapes() returns true.
725     // And this shape must be different from sub-shapes of the main shape.
726     // So we create a compound containing
727     // 1) some sub-shapes of SMESH_Mesh::PseudoShape() corresponding to
728     //    srcMeshDS->GetPersistentId()
729     // 2) the 1-st vertex of the main shape to assure
730     //    SMESHDS_Mesh::IsGroupOfSubShapes(shape)==true
731     TopoDS_Shape shapeForSrcMesh;
732     TopTools_IndexedMapOfShape pseudoSubShapes;
733     TopExp::MapShapes( SMESH_Mesh::PseudoShape(), pseudoSubShapes );
734
735     // index of pseudoSubShapes corresponding to srcMeshDS
736     int    subIndex = 1 + srcMeshDS->GetPersistentId() % pseudoSubShapes.Extent();
737     int nbSubShapes = 1 + srcMeshDS->GetPersistentId() / pseudoSubShapes.Extent();
738
739     // try to find already present shapeForSrcMesh
740     SMESHDS_Mesh* tgtMeshDS = tgtMesh->GetMeshDS();
741     for ( int i = tgtMeshDS->MaxShapeIndex(); i > 0 && shapeForSrcMesh.IsNull(); --i )
742     {
743       const TopoDS_Shape& s = tgtMeshDS->IndexToShape(i);
744       if ( s.ShapeType() != TopAbs_COMPOUND ) break;
745       TopoDS_Iterator sSubIt( s );
746       for ( int iSub = 0; iSub < nbSubShapes && sSubIt.More(); ++iSub, sSubIt.Next() )
747         if ( pseudoSubShapes( subIndex+iSub ).IsSame( sSubIt.Value()))
748           if ( iSub+1 == nbSubShapes )
749           {
750             shapeForSrcMesh = s;
751             break;
752           }
753     }
754     if ( shapeForSrcMesh.IsNull() )
755     {
756       // make a new shapeForSrcMesh
757       BRep_Builder aBuilder;
758       TopoDS_Compound comp;
759       aBuilder.MakeCompound( comp );
760       shapeForSrcMesh = comp;
761       for ( int iSub = 0; iSub < nbSubShapes; ++iSub )
762         if ( subIndex+iSub <= pseudoSubShapes.Extent() )
763           aBuilder.Add( comp, pseudoSubShapes( subIndex+iSub ));
764       TopExp_Explorer vExp( tgtMeshDS->ShapeToMesh(), TopAbs_VERTEX );
765       aBuilder.Add( comp, vExp.Current() );
766     }
767     SMESH_subMesh* sm = tgtMesh->GetSubMesh( shapeForSrcMesh );
768     SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
769     if ( !smDS )
770       smDS = tgtMeshDS->NewSubMesh( sm->GetId() );
771
772     // make ordinary submesh from a complex one
773     if ( smDS->IsComplexSubmesh() )
774     {
775       list< const SMESHDS_SubMesh* > subSM;
776       SMESHDS_SubMeshIteratorPtr smIt = smDS->GetSubMeshIterator();
777       while ( smIt->more() ) subSM.push_back( smIt->next() );
778       list< const SMESHDS_SubMesh* >::iterator sub = subSM.begin();
779       for ( ; sub != subSM.end(); ++sub)
780         smDS->RemoveSubMesh( *sub );
781     }
782     return sm->GetId();
783   }
784
785   //================================================================================
786   /*!
787    * \brief Return a submesh to store nodes and elements of a copied mesh
788    * and set event listeners in order to clear
789    * imported mesh and groups as soon as submesh state requires it
790    */
791   //================================================================================
792
793   SMESHDS_SubMesh* getSubmeshForCopiedMesh(const SMESH_Mesh*                    srcMesh,
794                                            SMESH_Mesh*                          tgtMesh,
795                                            const TopoDS_Shape&                  tgtShape,
796                                            StdMeshers_Import_1D::TNodeNodeMap*& n2n,
797                                            StdMeshers_Import_1D::TElemElemMap*& e2e,
798                                            bool &                               toCopyGroups)
799   {
800     StdMeshers_Import_1D::getMaps( srcMesh, tgtMesh, n2n,e2e );
801
802     _ImportData* iData = _Listener::getImportData(srcMesh,tgtMesh);
803
804     SMESH_subMesh* importedSM = tgtMesh->GetSubMesh( tgtShape );
805     iData->addComputed( importedSM );
806     if ( iData->_computedSubM.size() != iData->_subM.size() )
807       return 0; // not all submeshes computed yet
808
809     toCopyGroups = !iData->_copyGroupSubM.empty();
810
811     if ( !iData->_copyMeshSubM.empty())
812     {
813       // make submesh to store a copied mesh
814       int smID = getSubmeshIDForCopiedMesh( srcMesh->GetMeshDS(), tgtMesh );
815       SMESHDS_SubMesh* subDS = tgtMesh->GetMeshDS()->NewSubMesh( smID );
816
817       iData->_importMeshSubID = smID;
818       iData->_importMeshSubDS = subDS;
819       return subDS;
820     }
821     return 0;
822   }
823
824   //================================================================================
825   /*!
826    * \brief Return minimal square length of edges of 1D and 2D elements sharing the node
827    */
828   //================================================================================
829
830   double getMinEdgeLength2( const SMDS_MeshNode* n )
831   {
832     SMESH_NodeXYZ p = n;
833     double minLen2 = Precision::Infinite();
834     for ( SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator(); eIt->more();  )
835     {
836       const SMDS_MeshElement*      e = eIt->next();
837       const SMDSAbs_ElementType type = e->GetType();
838       if ( type != SMDSAbs_Edge && type != SMDSAbs_Face )
839         continue;
840       int i = e->GetNodeIndex( n );
841       int iNext = SMESH_MesherHelper::WrapIndex( i + 1, e->NbCornerNodes() );
842       minLen2 = Min( minLen2, p.SquareDistance( e->GetNode( iNext )));
843       if ( type != SMDSAbs_Face )
844         continue;
845       int iPrev = SMESH_MesherHelper::WrapIndex( i - 1, e->NbCornerNodes() );
846       minLen2 = Min( minLen2, p.SquareDistance( e->GetNode( iPrev )));
847     }
848     return minLen2;
849   }
850
851 } // namespace
852
853 //=============================================================================
854 /*!
855  * Creates StdMeshers_Import_1D
856  */
857 //=============================================================================
858
859 StdMeshers_Import_1D::StdMeshers_Import_1D(int hypId, SMESH_Gen * gen)
860   :SMESH_1D_Algo(hypId, gen), _sourceHyp(0)
861 {
862   _name = "Import_1D";
863   _shapeType = (1 << TopAbs_EDGE);
864
865   _compatibleHypothesis.push_back("ImportSource1D");
866 }
867
868 //=============================================================================
869 /*!
870  * Check presence of a hypothesis
871  */
872 //=============================================================================
873
874 bool StdMeshers_Import_1D::CheckHypothesis
875                          (SMESH_Mesh&                          aMesh,
876                           const TopoDS_Shape&                  aShape,
877                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
878 {
879   _sourceHyp = 0;
880
881   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
882   if ( hyps.size() == 0 )
883   {
884     aStatus = SMESH_Hypothesis::HYP_MISSING;
885     return false;  // can't work with no hypothesis
886   }
887
888   if ( hyps.size() > 1 )
889   {
890     aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
891     return false;
892   }
893
894   const SMESHDS_Hypothesis *theHyp = hyps.front();
895
896   string hypName = theHyp->GetName();
897
898   if (hypName == _compatibleHypothesis.front())
899   {
900     _sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
901     aStatus = _sourceHyp->GetGroups().empty() ? HYP_BAD_PARAMETER : HYP_OK;
902     if ( aStatus == HYP_BAD_PARAMETER )
903       _Listener::waitHypModification( aMesh.GetSubMesh( aShape ));
904     return aStatus == HYP_OK;
905   }
906
907   aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
908   return false;
909 }
910
911 //=============================================================================
912 /*!
913  * Import elements from the other mesh
914  */
915 //=============================================================================
916
917 bool StdMeshers_Import_1D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
918 {
919   if ( !_sourceHyp ) return false;
920
921   //MESSAGE("---------> StdMeshers_Import_1D::Compute");
922   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups(/*loaded=*/true);
923   if ( srcGroups.empty() )
924     return error("Invalid source groups");
925
926   SMESH_MesherHelper helper(theMesh);
927   helper.SetSubShape(theShape);
928   SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
929
930   const TopoDS_Edge& geomEdge = TopoDS::Edge( theShape );
931   const double edgeTol = BRep_Tool::Tolerance( geomEdge );
932   const int shapeID = tgtMesh->ShapeToIndex( geomEdge );
933
934
935   double geomTol = Precision::Confusion();
936   for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
937   {
938     const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
939     for ( SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements(); srcElems->more(); )
940     {
941       const SMDS_MeshElement* edge = srcElems->next();
942       geomTol = Sqrt( 0.5 * ( getMinEdgeLength2( edge->GetNode(0) ) +
943                               getMinEdgeLength2( edge->GetNode(1) ))) / 25;
944       iG = srcGroups.size();
945       break;
946     }
947   }
948   CurveProjector curveProjector( geomEdge, geomTol );
949
950   // get nodes on vertices
951   set<int> vertexIDs;
952   list < SMESH_TNodeXYZ > vertexNodes;
953   list < SMESH_TNodeXYZ >::iterator vNIt;
954   TopExp_Explorer vExp( theShape, TopAbs_VERTEX );
955   for ( ; vExp.More(); vExp.Next() )
956   {
957     const TopoDS_Vertex& v = TopoDS::Vertex( vExp.Current() );
958     if ( !vertexIDs.insert( tgtMesh->ShapeToIndex( v )).second )
959       continue; // closed edge
960     const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
961     if ( !n )
962     {
963       _gen->Compute(theMesh,v,/*anUpward=*/true);
964       n = SMESH_Algo::VertexNode( v, tgtMesh );
965       //MESSAGE("_gen->Compute " << n);
966       if ( !n ) return false; // very strange
967     }
968     vertexNodes.push_back( SMESH_TNodeXYZ( n ));
969     //MESSAGE("SMESH_Algo::VertexNode " << n->GetID() << " " << n->X() << " " << n->Y() << " " << n->Z() );
970   }
971
972   // import edges from groups
973   TNodeNodeMap* n2n;
974   TElemElemMap* e2e;
975   for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
976   {
977     const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
978
979     const int meshID = srcGroup->GetMesh()->GetPersistentId();
980     const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
981     if ( !srcMesh ) continue;
982     getMaps( srcMesh, &theMesh, n2n, e2e );
983
984     SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
985     vector<const SMDS_MeshNode*> newNodes;
986     while ( srcElems->more() ) // loop on group contents
987     {
988       const SMDS_MeshElement* edge = srcElems->next();
989       gp_XYZ middle = 0.5 * ( SMESH_NodeXYZ( edge->GetNode(0)) +
990                               SMESH_NodeXYZ( edge->GetNode(1)));
991       if ( curveProjector.IsOut( middle ))
992         continue;
993
994       // find or create nodes of a new edge
995       newNodes.resize( edge->NbNodes() );
996       newNodes.back() = 0;
997       int nbNodesOnVertex = 0;
998       SMDS_MeshElement::iterator node = edge->begin_nodes();
999       for ( size_t i = 0; i < newNodes.size(); ++i, ++node )
1000       {
1001         TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, nullptr )).first;
1002         if ( n2nIt->second )
1003         {
1004           int sId = n2nIt->second->getshapeId();
1005           if ( sId != shapeID )
1006           {
1007             if ( vertexIDs.count( sId ))
1008               ++nbNodesOnVertex;
1009             else
1010               break;
1011           }
1012         }
1013         else if ( !vertexNodes.empty() )
1014         {
1015           // find an existing vertex node
1016           double checktol = max(1.E-10, 10*edgeTol*edgeTol);
1017           for ( vNIt = vertexNodes.begin(); vNIt != vertexNodes.end(); ++vNIt)
1018             if ( vNIt->SquareDistance( *node ) < checktol)
1019             {
1020               (*n2nIt).second = vNIt->_node;
1021               vertexNodes.erase( vNIt );
1022               ++nbNodesOnVertex;
1023               break;
1024             }
1025         }
1026         if ( !n2nIt->second )
1027         {
1028           // find out if the node lies on theShape
1029           SMESH_NodeXYZ xyz = *node;
1030           double dist2, u;
1031           if ( curveProjector.IsOnCurve( xyz, dist2, u ))
1032           {
1033             // tolerance relative to the length of surrounding edges
1034             double mytol2 = getMinEdgeLength2( *node ) / 25 / 25;
1035             if ( dist2 < mytol2 )
1036             {
1037               SMDS_MeshNode* newNode = tgtMesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1038               n2nIt->second = newNode;
1039               tgtMesh->SetNodeOnEdge( newNode, shapeID, u );
1040             }
1041           }
1042         }
1043         if ( !(newNodes[i] = n2nIt->second ))
1044           break;
1045       }
1046       if ( !newNodes.back() )
1047       {
1048         //MESSAGE("not all nodes of edge lie on theShape");
1049         continue; // not all nodes of edge lie on theShape
1050       }
1051
1052       // make a new edge
1053       SMDS_MeshElement * newEdge;
1054       if ( newNodes.size() == 3 )
1055         newEdge = tgtMesh->AddEdge( newNodes[0], newNodes[1], newNodes[2] );
1056       else
1057         newEdge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
1058       tgtMesh->SetMeshElementOnShape( newEdge, shapeID );
1059       e2e->insert( make_pair( edge, newEdge ));
1060
1061       if ( nbNodesOnVertex >= 2 ) // EDGE is meshed by a sole segment
1062       {
1063         iG = srcGroups.size(); // stop looingp on groups
1064         break;
1065       }
1066     }  // loop on group contents
1067   } // loop on groups
1068
1069   if ( n2n->empty())
1070     return error("Empty source groups");
1071
1072   // check if the whole geom edge is covered by imported segments;
1073   // the check consist in passing by segments from one vetrex node to another
1074   bool isEdgeMeshed = false;
1075   if ( SMESHDS_SubMesh* tgtSM = tgtMesh->MeshElements( theShape ))
1076   {
1077     const TopoDS_Vertex& v = ( vExp.ReInit(), TopoDS::Vertex( vExp.Current() ));
1078     const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
1079     const SMDS_MeshElement* seg = 0;
1080     SMDS_ElemIteratorPtr segIt = n->GetInverseElementIterator(SMDSAbs_Edge);
1081     while ( segIt->more() && !seg )
1082       if ( !tgtSM->Contains( seg = segIt->next()))
1083         seg = 0;
1084     int nbPassedSegs = 0;
1085     while ( seg )
1086     {
1087       ++nbPassedSegs;
1088       const SMDS_MeshNode* n2 = seg->GetNode(0);
1089       n = ( n2 == n ? seg->GetNode(1) : n2 );
1090       if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
1091         break;
1092       const SMDS_MeshElement* seg2 = 0;
1093       segIt = n->GetInverseElementIterator(SMDSAbs_Edge);
1094       while ( segIt->more() && !seg2 )
1095         if ( seg == ( seg2 = segIt->next()))
1096           seg2 = 0;
1097       seg = seg2;
1098     }
1099     if (nbPassedSegs > 0 && tgtSM->NbElements() > nbPassedSegs )
1100       return error( "Source elements overlap one another");
1101
1102     isEdgeMeshed = ( tgtSM->NbElements() == nbPassedSegs &&
1103                      n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1104   }
1105   if ( !isEdgeMeshed )
1106     return error( "Source elements don't cover totally the geometrical edge" );
1107
1108   // copy meshes
1109   vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
1110   for ( size_t i = 0; i < srcMeshes.size(); ++i )
1111     importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
1112
1113   return true;
1114 }
1115
1116 //================================================================================
1117 /*!
1118  * \brief Copy mesh and groups
1119  */
1120 //================================================================================
1121
1122 void StdMeshers_Import_1D::importMesh(const SMESH_Mesh*          srcMesh,
1123                                       SMESH_Mesh &               tgtMesh,
1124                                       StdMeshers_ImportSource1D* srcHyp,
1125                                       const TopoDS_Shape&        tgtShape)
1126 {
1127   // get submesh to store the imported mesh
1128   TNodeNodeMap* n2n;
1129   TElemElemMap* e2e;
1130   bool toCopyGroups;
1131   SMESHDS_SubMesh* tgtSubMesh =
1132     getSubmeshForCopiedMesh( srcMesh, &tgtMesh, tgtShape, n2n, e2e, toCopyGroups );
1133   if ( !tgtSubMesh || tgtSubMesh->NbNodes() + tgtSubMesh->NbElements() > 0 )
1134     return; // not to copy srcMeshDS twice
1135
1136   SMESHDS_Mesh* tgtMeshDS = tgtMesh.GetMeshDS();
1137   SMESH_MeshEditor additor( &tgtMesh );
1138
1139   // 1. Copy mesh
1140
1141   SMESH_MeshEditor::ElemFeatures elemType;
1142   vector<const SMDS_MeshNode*> newNodes;
1143   const SMESHDS_Mesh* srcMeshDS = srcMesh->GetMeshDS();
1144   SMDS_ElemIteratorPtr eIt = srcMeshDS->elementsIterator();
1145   while ( eIt->more() )
1146   {
1147     const SMDS_MeshElement* elem = eIt->next();
1148     TElemElemMap::iterator e2eIt = e2e->insert( make_pair( elem, (SMDS_MeshElement*)0 )).first;
1149     if ( e2eIt->second ) continue; // already copied by Compute()
1150     newNodes.resize( elem->NbNodes() );
1151     SMDS_MeshElement::iterator node = elem->begin_nodes();
1152     for ( unsigned i = 0; i < newNodes.size(); ++i, ++node )
1153     {
1154       TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
1155       if ( !n2nIt->second )
1156       {
1157         (*n2nIt).second = tgtMeshDS->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
1158         tgtSubMesh->AddNode( n2nIt->second );
1159       }
1160       newNodes[i] = n2nIt->second;
1161     }
1162     const SMDS_MeshElement* newElem =
1163       tgtMeshDS->FindElement( newNodes, elem->GetType(), /*noMedium=*/false );
1164     if ( !newElem )
1165     {
1166       newElem = additor.AddElement( newNodes, elemType.Init( elem, /*basicOnly=*/false ));
1167       tgtSubMesh->AddElement( newElem );
1168     }
1169     if ( toCopyGroups )
1170       (*e2eIt).second = newElem;
1171   }
1172   // copy free nodes
1173   if ( srcMeshDS->NbNodes() > (int) n2n->size() )
1174   {
1175     SMDS_NodeIteratorPtr nIt = srcMeshDS->nodesIterator();
1176     while( nIt->more() )
1177     {
1178       const SMDS_MeshNode* node = nIt->next();
1179       if ( node->NbInverseElements() == 0 )
1180       {
1181         const SMDS_MeshNode* newNode = tgtMeshDS->AddNode( node->X(), node->Y(), node->Z());
1182         n2n->insert( make_pair( node, newNode ));
1183         tgtSubMesh->AddNode( newNode );
1184       }
1185     }
1186   }
1187
1188   // 2. Copy groups
1189
1190   vector<SMESH_Group*> resultGroups;
1191   if ( toCopyGroups )
1192   {
1193     // collect names of existing groups to assure uniqueness of group names within a type
1194     map< SMDSAbs_ElementType, set<string> > namesByType;
1195     SMESH_Mesh::GroupIteratorPtr groupIt = tgtMesh.GetGroups();
1196     while ( groupIt->more() )
1197     {
1198       SMESH_Group* tgtGroup = groupIt->next();
1199       namesByType[ tgtGroup->GetGroupDS()->GetType() ].insert( tgtGroup->GetName() );
1200     }
1201     if (srcMesh)
1202     {
1203       SMESH_Mesh::GroupIteratorPtr groupIt = srcMesh->GetGroups();
1204       while ( groupIt->more() )
1205       {
1206         SMESH_Group* srcGroup = groupIt->next();
1207         SMESHDS_GroupBase* srcGroupDS = srcGroup->GetGroupDS();
1208         string name = srcGroup->GetName();
1209         int nb = 1;
1210         while ( !namesByType[ srcGroupDS->GetType() ].insert( name ).second )
1211           name = SMESH_Comment(srcGroup->GetName()) << "_imported_" << nb++;
1212         SMESH_Group* newGroup = tgtMesh.AddGroup( srcGroupDS->GetType(), name.c_str() );
1213         SMESHDS_Group* newGroupDS = (SMESHDS_Group*)newGroup->GetGroupDS();
1214         resultGroups.push_back( newGroup );
1215
1216         eIt = srcGroupDS->GetElements();
1217         if ( srcGroupDS->GetType() == SMDSAbs_Node )
1218           while (eIt->more())
1219           {
1220             TNodeNodeMap::iterator n2nIt = n2n->find((const SMDS_MeshNode*) eIt->next() );
1221             if ( n2nIt != n2n->end() && n2nIt->second )
1222               newGroupDS->SMDSGroup().Add((*n2nIt).second );
1223           }
1224         else
1225           while (eIt->more())
1226           {
1227             TElemElemMap::iterator e2eIt = e2e->find( eIt->next() );
1228             if ( e2eIt != e2e->end() && e2eIt->second )
1229               newGroupDS->SMDSGroup().Add((*e2eIt).second );
1230           }
1231       }
1232     }
1233   }
1234   n2n->clear();
1235   e2e->clear();
1236
1237   // Remember created groups in order to remove them as soon as the srcHyp is
1238   // modified or something other similar happens. This imformation must be persistent,
1239   // for that store them in a hypothesis as it stores its values in the file anyway
1240   srcHyp->StoreResultGroups( resultGroups, *srcMeshDS, *tgtMeshDS );
1241 }
1242
1243 //=============================================================================
1244 /*!
1245  * \brief Set needed event listeners and create a submesh for a copied mesh
1246  *
1247  * This method is called only if a submesh has HYP_OK algo_state.
1248  */
1249 //=============================================================================
1250
1251 void StdMeshers_Import_1D::setEventListener(SMESH_subMesh*             subMesh,
1252                                             StdMeshers_ImportSource1D* sourceHyp)
1253 {
1254   if ( sourceHyp )
1255   {
1256     vector<SMESH_Mesh*> srcMeshes = sourceHyp->GetSourceMeshes();
1257     if ( srcMeshes.empty() )
1258       _Listener::waitHypModification( subMesh );
1259     for ( unsigned i = 0; i < srcMeshes.size(); ++i )
1260       // set a listener to remove the imported mesh and groups
1261       _Listener::storeImportSubmesh( subMesh, srcMeshes[i], sourceHyp );
1262   }
1263 }
1264 void StdMeshers_Import_1D::SetEventListener(SMESH_subMesh* subMesh)
1265 {
1266   if ( !_sourceHyp )
1267   {
1268     const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
1269     SMESH_Mesh*         tgtMesh  = subMesh->GetFather();
1270     Hypothesis_Status aStatus;
1271     CheckHypothesis( *tgtMesh, tgtShape, aStatus );
1272   }
1273   setEventListener( subMesh, _sourceHyp );
1274 }
1275
1276 void StdMeshers_Import_1D::SubmeshRestored(SMESH_subMesh* subMesh)
1277 {
1278   SetEventListener(subMesh);
1279 }
1280
1281 //=============================================================================
1282 /*!
1283  * Predict nb of mesh entities created by Compute()
1284  */
1285 //=============================================================================
1286
1287 bool StdMeshers_Import_1D::Evaluate(SMESH_Mesh &         theMesh,
1288                                     const TopoDS_Shape & theShape,
1289                                     MapShapeNbElems&     aResMap)
1290 {
1291   if ( !_sourceHyp ) return false;
1292
1293   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
1294   if ( srcGroups.empty() )
1295     return error("Invalid source groups");
1296
1297   vector<int> aVec(SMDSEntity_Last,0);
1298
1299   bool toCopyMesh, toCopyGroups;
1300   _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
1301   if ( toCopyMesh ) // the whole mesh is copied
1302   {
1303     vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
1304     for ( unsigned i = 0; i < srcMeshes.size(); ++i )
1305     {
1306       SMESH_subMesh* sm = getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
1307       if ( !sm || aResMap.count( sm )) continue; // already counted
1308       aVec.assign( SMDSEntity_Last, 0);
1309       const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
1310       for (int i = 0; i < SMDSEntity_Last; i++)
1311         aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
1312     }
1313   }
1314   else
1315   {
1316     SMESH_MesherHelper helper(theMesh);
1317
1318     const TopoDS_Edge& geomEdge = TopoDS::Edge( theShape );
1319     const double edgeTol = helper.MaxTolerance( geomEdge );
1320
1321     // take into account nodes on vertices
1322     TopExp_Explorer vExp( theShape, TopAbs_VERTEX );
1323     for ( ; vExp.More(); vExp.Next() )
1324       theMesh.GetSubMesh( vExp.Current())->Evaluate( aResMap );
1325
1326     // count edges imported from groups
1327     int nbEdges = 0, nbQuadEdges = 0;
1328     for ( size_t iG = 0; iG < srcGroups.size(); ++iG )
1329     {
1330       const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
1331       SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
1332       SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
1333       while ( srcElems->more() ) // loop on group contents
1334       {
1335         const SMDS_MeshElement* edge = srcElems->next();
1336         // find out if edge is located on geomEdge by projecting
1337         // a middle of edge to geomEdge
1338         SMESH_TNodeXYZ p1( edge->GetNode(0));
1339         SMESH_TNodeXYZ p2( edge->GetNode(1));
1340         gp_XYZ middle = ( p1 + p2 ) / 2.;
1341         tmpNode->setXYZ( middle.X(), middle.Y(), middle.Z());
1342         double u = 0;
1343         if ( helper.CheckNodeU( geomEdge, tmpNode, u, 10 * edgeTol, /*force=*/true ))
1344           ++( edge->IsQuadratic() ? nbQuadEdges : nbEdges);
1345       }
1346       helper.GetMeshDS()->RemoveNode(tmpNode);
1347     }
1348
1349     int nbNodes = nbEdges + 2 * nbQuadEdges - 1;
1350
1351     aVec[SMDSEntity_Node     ] = nbNodes;
1352     aVec[SMDSEntity_Edge     ] = nbEdges;
1353     aVec[SMDSEntity_Quad_Edge] = nbQuadEdges;
1354   }
1355
1356   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
1357   aResMap.insert( make_pair( sm, aVec ));
1358
1359   return true;
1360 }
1361
1362 //================================================================================
1363 /*!
1364  * \brief Return node-node and element-element maps for import of geiven source mesh
1365  */
1366 //================================================================================
1367
1368 void StdMeshers_Import_1D::getMaps(const SMESH_Mesh* srcMesh,
1369                                    SMESH_Mesh*       tgtMesh,
1370                                    TNodeNodeMap*&    n2n,
1371                                    TElemElemMap*&    e2e)
1372 {
1373   _ImportData* iData = _Listener::getImportData(srcMesh,tgtMesh);
1374   n2n = &iData->_n2n;
1375   e2e = &iData->_e2e;
1376   if ( iData->_copyMeshSubM.empty() )
1377   {
1378     // n2n->clear(); -- for sharing nodes on EDGEs
1379     e2e->clear();
1380   }
1381 }
1382
1383 //================================================================================
1384 /*!
1385  * \brief Return submesh corresponding to the copied mesh
1386  */
1387 //================================================================================
1388
1389 SMESH_subMesh* StdMeshers_Import_1D::getSubMeshOfCopiedMesh( SMESH_Mesh& tgtMesh,
1390                                                              SMESH_Mesh& srcMesh )
1391 {
1392   _ImportData* iData = _Listener::getImportData(&srcMesh,&tgtMesh);
1393   if ( iData->_copyMeshSubM.empty() ) return 0;
1394   SMESH_subMesh* sm = tgtMesh.GetSubMeshContaining( iData->_importMeshSubID );
1395   return sm;
1396 }
1397