Salome HOME
PR: relax constraints on node distances on StdMeshers_import_1D
[modules/smesh.git] / src / StdMeshers / StdMeshers_Import_1D.cxx
1 // Copyright (C) 2007-2013  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.
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 : implementaion 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_MesherHelper.hxx"
40 #include "SMESH_subMesh.hxx"
41 #include "SMESH_subMeshEventListener.hxx"
42
43 #include "Utils_SALOME_Exception.hxx"
44 #include "utilities.h"
45
46 #include <BRep_Builder.hxx>
47 #include <BRep_Tool.hxx>
48 #include <TopExp.hxx>
49 #include <TopExp_Explorer.hxx>
50 #include <TopoDS.hxx>
51 #include <TopoDS_Compound.hxx>
52 #include <TopoDS_Edge.hxx>
53 #include <TopoDS_Vertex.hxx>
54
55 using namespace std;
56
57 //=============================================================================
58 /*!
59  * Creates StdMeshers_Import_1D
60  */
61 //=============================================================================
62
63 StdMeshers_Import_1D::StdMeshers_Import_1D(int hypId, int studyId, SMESH_Gen * gen)
64   :SMESH_1D_Algo(hypId, studyId, gen), _sourceHyp(0)
65 {
66   MESSAGE("StdMeshers_Import_1D::StdMeshers_Import_1D");
67   _name = "Import_1D";
68   _shapeType = (1 << TopAbs_EDGE);
69
70   _compatibleHypothesis.push_back("ImportSource1D");
71 }
72
73 //=============================================================================
74 /*!
75  * Check presence of a hypothesis
76  */
77 //=============================================================================
78
79 bool StdMeshers_Import_1D::CheckHypothesis
80                          (SMESH_Mesh&                          aMesh,
81                           const TopoDS_Shape&                  aShape,
82                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
83 {
84   _sourceHyp = 0;
85
86   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
87   if ( hyps.size() == 0 )
88   {
89     aStatus = SMESH_Hypothesis::HYP_MISSING;
90     return false;  // can't work with no hypothesis
91   }
92
93   if ( hyps.size() > 1 )
94   {
95     aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
96     return false;
97   }
98
99   const SMESHDS_Hypothesis *theHyp = hyps.front();
100
101   string hypName = theHyp->GetName();
102
103   if (hypName == _compatibleHypothesis.front())
104   {
105     _sourceHyp = (StdMeshers_ImportSource1D *)theHyp;
106     aStatus = SMESH_Hypothesis::HYP_OK;
107     return true;
108   }
109
110   aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
111   return true;
112 }
113
114 //================================================================================
115 namespace // INTERNAL STUFF
116 //================================================================================
117 {
118   int getSubmeshIDForCopiedMesh(const SMESHDS_Mesh* srcMeshDS, SMESH_Mesh* tgtMesh);
119
120   enum _ListenerDataType
121     {
122       WAIT_HYP_MODIF=1, // data indicating awaiting for valid parameters of src hyp
123       LISTEN_SRC_MESH, // data storing submesh depending on source mesh state
124       SRC_HYP // data storing ImportSource hyp
125     };
126   //================================================================================
127   /*!
128    * \brief _ListenerData holding ImportSource hyp holding in its turn
129    *  imported groups
130    */
131   struct _ListenerData : public SMESH_subMeshEventListenerData
132   {
133     const StdMeshers_ImportSource1D* _srcHyp;
134     _ListenerData(const StdMeshers_ImportSource1D* h, _ListenerDataType type=SRC_HYP):
135       SMESH_subMeshEventListenerData(/*isDeletable=*/true), _srcHyp(h)
136     {
137       myType = type;
138     }
139   };
140   //================================================================================
141   /*!
142    * \brief Comparator of sub-meshes
143    */
144   struct _SubLess
145   {
146     bool operator()(const SMESH_subMesh* sm1, const SMESH_subMesh* sm2 ) const
147     {
148       if ( sm1 == sm2 ) return false;
149       if ( !sm1 || !sm2 ) return sm1 < sm2;
150       const TopoDS_Shape& s1 = sm1->GetSubShape();
151       const TopoDS_Shape& s2 = sm2->GetSubShape();
152       TopAbs_ShapeEnum t1 = s1.IsNull() ? TopAbs_SHAPE : s1.ShapeType();
153       TopAbs_ShapeEnum t2 = s2.IsNull() ? TopAbs_SHAPE : s2.ShapeType();
154       if ( t1 == t2)
155         return (sm1 < sm2);
156       return t1 < t2; // to have: face < edge
157     }
158   };
159   //================================================================================
160   /*!
161    * \brief Container of data dedicated to one source mesh
162    */
163   struct _ImportData
164   {
165     const SMESH_Mesh* _srcMesh;
166     StdMeshers_Import_1D::TNodeNodeMap _n2n;
167     StdMeshers_Import_1D::TElemElemMap _e2e;
168
169     set< SMESH_subMesh*, _SubLess > _subM; // submeshes relating to this srcMesh
170     set< SMESH_subMesh*, _SubLess > _copyMeshSubM; // submeshes requesting mesh copying
171     set< SMESH_subMesh*, _SubLess > _copyGroupSubM; // submeshes requesting group copying
172     set< SMESH_subMesh*, _SubLess > _computedSubM;
173
174     SMESHDS_SubMesh*     _importMeshSubDS; // submesh storing a copy of _srcMesh
175     int                  _importMeshSubID; // id of _importMeshSubDS
176
177     _ImportData(const SMESH_Mesh* srcMesh=0):
178       _srcMesh(srcMesh), _importMeshSubDS(0),_importMeshSubID(-1) {}
179
180     void removeImportedMesh( SMESHDS_Mesh* meshDS )
181     {
182       if ( !_importMeshSubDS ) return;
183       SMDS_ElemIteratorPtr eIt = _importMeshSubDS->GetElements();
184       while ( eIt->more() )
185         meshDS->RemoveFreeElement( eIt->next(), _importMeshSubDS, /*fromGroups=*/false );
186       SMDS_NodeIteratorPtr nIt = _importMeshSubDS->GetNodes();
187       while ( nIt->more() )
188         meshDS->RemoveFreeNode( nIt->next(), _importMeshSubDS, /*fromGroups=*/false );
189       _n2n.clear();
190       _e2e.clear();
191     }
192     void removeGroups( SMESH_subMesh* subM, const StdMeshers_ImportSource1D* srcHyp )
193     {
194       if ( !srcHyp ) return;
195       SMESH_Mesh*           tgtMesh = subM->GetFather();
196       const SMESHDS_Mesh* tgtMeshDS = tgtMesh->GetMeshDS();
197       const SMESHDS_Mesh* srcMeshDS = _srcMesh->GetMeshDS();
198       vector<SMESH_Group*>*  groups =
199         const_cast<StdMeshers_ImportSource1D*>(srcHyp)->GetResultGroups(*srcMeshDS,*tgtMeshDS);
200       if ( groups )
201       {
202         for ( unsigned i = 0; i < groups->size(); ++i )
203           tgtMesh->RemoveGroup( groups->at(i)->GetGroupDS()->GetID() );
204         groups->clear();
205       }
206     }
207     void trackHypParams( SMESH_subMesh* sm, const StdMeshers_ImportSource1D* srcHyp )
208     {
209       if ( !srcHyp ) return;
210       bool toCopyMesh, toCopyGroups;
211       srcHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
212
213       if ( toCopyMesh )_copyMeshSubM.insert( sm );
214       else             _copyMeshSubM.erase( sm );
215
216       if ( toCopyGroups ) _copyGroupSubM.insert( sm );
217       else                _copyGroupSubM.erase( sm );
218     }
219     void addComputed( SMESH_subMesh* sm )
220     {
221       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
222                                                                /*complexShapeFirst=*/true);
223       while ( smIt->more() )
224       {
225         sm = smIt->next();
226         switch ( sm->GetSubShape().ShapeType() )
227         {
228         case TopAbs_EDGE:
229         case TopAbs_FACE:
230           _subM.insert( sm );
231           if ( !sm->IsEmpty() )
232             _computedSubM.insert( sm );
233         case TopAbs_VERTEX:
234           break;
235         default:;
236         }
237       }
238     }
239   };
240   //================================================================================
241   /*!
242    * Listener notified on events relating to imported submesh
243    */
244   class _Listener : public SMESH_subMeshEventListener
245   {
246     typedef map< SMESH_Mesh*, list< _ImportData > > TMesh2ImpData;
247     TMesh2ImpData _tgtMesh2ImportData;
248
249     _Listener():SMESH_subMeshEventListener(/*isDeletable=*/false,
250                                            "StdMeshers_Import_1D::_Listener") {}
251
252   public:
253     // return poiter to a static listener
254     static _Listener* get() { static _Listener theListener; return &theListener; }
255
256     static _ImportData* getImportData(const SMESH_Mesh* srcMesh, SMESH_Mesh* tgtMesh);
257
258     static void storeImportSubmesh(SMESH_subMesh*                   importSub,
259                                    const SMESH_Mesh*                srcMesh,
260                                    const StdMeshers_ImportSource1D* srcHyp);
261
262     virtual void ProcessEvent(const int                       event,
263                               const int                       eventType,
264                               SMESH_subMesh*                  subMesh,
265                               SMESH_subMeshEventListenerData* data,
266                               const SMESH_Hypothesis*         hyp);
267     void removeSubmesh( SMESH_subMesh* sm, _ListenerData* data );
268     void clearSubmesh ( SMESH_subMesh* sm, _ListenerData* data, bool clearAllSub );
269
270     // mark sm as missing src hyp with valid groups
271     static void waitHypModification(SMESH_subMesh* sm)
272     {
273       sm->SetEventListener
274         (get(), SMESH_subMeshEventListenerData::MakeData( sm, WAIT_HYP_MODIF ), sm);
275     }
276   };
277   //--------------------------------------------------------------------------------
278   /*!
279    * \brief Find or create ImportData for given meshes
280    */
281   _ImportData* _Listener::getImportData(const SMESH_Mesh* srcMesh,
282                                         SMESH_Mesh*       tgtMesh)
283   {
284     list< _ImportData >& dList = get()->_tgtMesh2ImportData[tgtMesh];
285     list< _ImportData >::iterator d = dList.begin();
286     for ( ; d != dList.end(); ++d )
287       if ( d->_srcMesh == srcMesh )
288         return &*d;
289     dList.push_back(_ImportData(srcMesh));
290     return &dList.back();
291   }
292
293   //--------------------------------------------------------------------------------
294   /*!
295    * \brief Remember an imported sub-mesh and set needed even listeners
296    *  \param importSub - submesh computed by Import algo
297    *  \param srcMesh - source mesh
298    *  \param srcHyp - ImportSource hypothesis
299    */
300   void _Listener::storeImportSubmesh(SMESH_subMesh*                   importSub,
301                                      const SMESH_Mesh*                srcMesh,
302                                      const StdMeshers_ImportSource1D* srcHyp)
303   {
304     // set listener to hear events of the submesh computed by "Import" algo
305     importSub->SetEventListener( get(), new _ListenerData(srcHyp), importSub );
306
307     // set listeners to hear events of the source mesh
308     SMESH_subMesh* smToNotify = importSub;
309     vector<SMESH_subMesh*> smToListen = srcHyp->GetSourceSubMeshes( srcMesh );
310     for ( size_t i = 0; i < smToListen.size(); ++i )
311     {
312       SMESH_subMeshEventListenerData* data = new _ListenerData(srcHyp, LISTEN_SRC_MESH);
313       data->mySubMeshes.push_back( smToNotify );
314       importSub->SetEventListener( get(), data, smToListen[i] );
315     }
316     // remember the submesh importSub and its sub-submeshes
317     _ImportData* iData = _Listener::getImportData( srcMesh, importSub->GetFather());
318     iData->trackHypParams( importSub, srcHyp );
319     iData->addComputed( importSub );
320     if ( !iData->_copyMeshSubM.empty() && iData->_importMeshSubID < 1 )
321     {
322       SMESH_Mesh* tgtMesh = importSub->GetFather();
323       iData->_importMeshSubID = getSubmeshIDForCopiedMesh( srcMesh->GetMeshDS(),tgtMesh);
324       iData->_importMeshSubDS = tgtMesh->GetMeshDS()->NewSubMesh( iData->_importMeshSubID );
325     }
326   }
327   //--------------------------------------------------------------------------------
328   /*!
329    * \brief Remove imported mesh and/or groups if needed
330    *  \param sm - submesh loosing Import algo
331    *  \param data - data holding imported groups
332    */
333   void _Listener::removeSubmesh( SMESH_subMesh* sm, _ListenerData* data )
334   {
335     list< _ImportData > &  dList = _tgtMesh2ImportData[ sm->GetFather() ];
336     list< _ImportData >::iterator d = dList.begin();
337     for ( ; d != dList.end(); ++d )
338       if ( (*d)._subM.erase( sm ))
339       {
340         d->_computedSubM.erase( sm );
341         bool rmMesh   = d->_copyMeshSubM.erase( sm ) && d->_copyMeshSubM.empty();
342         bool rmGroups = (d->_copyGroupSubM.erase( sm ) && d->_copyGroupSubM.empty()) || rmMesh;
343         if ( rmMesh )
344           d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
345         if ( rmGroups && data )
346           d->removeGroups( sm, data->_srcHyp );
347       }
348   }
349   //--------------------------------------------------------------------------------
350   /*!
351    * \brief Clear submeshes and remove imported mesh and/or groups if necessary
352    *  \param sm - cleared submesh
353    *  \param data - data holding imported groups
354    */
355   void _Listener::clearSubmesh(SMESH_subMesh* sm, _ListenerData* data, bool clearAllSub)
356   {
357     list< _ImportData > &  dList = _tgtMesh2ImportData[ sm->GetFather() ];
358     list< _ImportData >::iterator d = dList.begin();
359     for ( ; d != dList.end(); ++d )
360     {
361       if ( !d->_subM.count( sm )) continue;
362       if ( (*d)._computedSubM.erase( sm ) )
363       {
364         bool copyMesh = !d->_copyMeshSubM.empty();
365         if ( copyMesh || clearAllSub )
366         {
367           // remove imported mesh and groups
368           d->removeImportedMesh( sm->GetFather()->GetMeshDS() );
369
370           if ( data )
371             d->removeGroups( sm, data->_srcHyp );
372
373           // clear the rest submeshes
374           if ( !d->_computedSubM.empty() )
375           {
376             d->_computedSubM.clear();
377             set< SMESH_subMesh*, _SubLess>::iterator sub = d->_subM.begin();
378             for ( ; sub != d->_subM.end(); ++sub )
379             {
380               SMESH_subMesh* subM = *sub;
381               _ListenerData* hypData = (_ListenerData*) subM->GetEventListenerData( get() );
382               if ( hypData )
383                 d->removeGroups( sm, hypData->_srcHyp );
384
385               subM->ComputeStateEngine( SMESH_subMesh::CLEAN );
386               if ( subM->GetSubShape().ShapeType() == TopAbs_FACE )
387                 subM->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
388             }
389           }
390         }
391         sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
392         if ( sm->GetSubShape().ShapeType() == TopAbs_FACE )
393           sm->ComputeSubMeshStateEngine( SMESH_subMesh::CLEAN );
394       }
395       if ( data )
396         d->trackHypParams( sm, data->_srcHyp );
397       d->_n2n.clear();
398       d->_e2e.clear();
399     }
400   }
401   //--------------------------------------------------------------------------------
402   /*!
403    * \brief Remove imported mesh and/or groups
404    */
405   void _Listener::ProcessEvent(const int                       event,
406                                const int                       eventType,
407                                SMESH_subMesh*                  subMesh,
408                                SMESH_subMeshEventListenerData* data,
409                                const SMESH_Hypothesis*         /*hyp*/)
410   {
411     if ( data && data->myType == WAIT_HYP_MODIF )
412     {
413       // event of Import submesh
414       if ( SMESH_subMesh::MODIF_HYP  == event &&
415            SMESH_subMesh::ALGO_EVENT == eventType )
416       {
417         // re-call SetEventListener() to take into account valid parameters
418         // of ImportSource hypothesis
419         if ( SMESH_Algo* algo = subMesh->GetAlgo() )
420           algo->SetEventListener( subMesh );
421       }
422     }
423     else if ( data && data->myType == LISTEN_SRC_MESH )
424     {
425       // event of source mesh
426       if ( SMESH_subMesh::COMPUTE_EVENT == eventType )
427       {
428         switch ( event ) {
429         case SMESH_subMesh::CLEAN:
430           // source mesh cleaned -> clean target mesh
431           clearSubmesh( data->mySubMeshes.front(), (_ListenerData*) data, /*all=*/true );
432           break;
433         case SMESH_subMesh::SUBMESH_COMPUTED: {
434           // source mesh computed -> reset FAILED state of Import submeshes to
435           // READY_TO_COMPUTE
436           SMESH_Mesh* srcMesh = subMesh->GetFather();
437           if ( srcMesh->NbEdges() > 0 || srcMesh->NbFaces() > 0 )
438           {
439             SMESH_Mesh* m = data->mySubMeshes.front()->GetFather();
440             if ( SMESH_subMesh* sm1 = m->GetSubMeshContaining(1))
441             {
442               sm1->ComputeStateEngine(SMESH_subMesh::SUBMESH_COMPUTED );
443               sm1->ComputeSubMeshStateEngine( SMESH_subMesh::SUBMESH_COMPUTED );
444             }
445           }
446           break;
447         }
448         default:;
449         }
450       }
451     }
452     else // event of Import submesh
453     {
454       // find out what happens: import hyp modified or removed
455       bool removeImport = false, modifHyp = false;
456       if ( SMESH_subMesh::ALGO_EVENT == eventType )
457         modifHyp = true;
458       if ( subMesh->GetAlgoState() != SMESH_subMesh::HYP_OK )
459       {
460         removeImport = true;
461       }
462       else if (( SMESH_subMesh::REMOVE_ALGO == event ||
463                  SMESH_subMesh::REMOVE_FATHER_ALGO == event ) &&
464                SMESH_subMesh::ALGO_EVENT == eventType )
465       {
466         SMESH_Algo* algo = subMesh->GetAlgo();
467         removeImport = ( strncmp( "Import", algo->GetName(), 6 ) != 0 );
468       }
469
470       if ( removeImport )
471       {
472         // treate removal of Import algo from subMesh
473         removeSubmesh( subMesh, (_ListenerData*) data );
474       }
475       else if ( modifHyp ||
476                 ( SMESH_subMesh::CLEAN         == event &&
477                   SMESH_subMesh::COMPUTE_EVENT == eventType))
478       {
479         // treate modification of ImportSource hypothesis
480         clearSubmesh( subMesh, (_ListenerData*) data, /*all=*/false );
481       }
482       else if ( SMESH_subMesh::CHECK_COMPUTE_STATE == event &&
483                 SMESH_subMesh::COMPUTE_EVENT       == eventType )
484       {
485         // check compute state of all submeshes impoting from same src mesh;
486         // this is to take into account 1D computed submeshes hidden by 2D import algo;
487         // else source mesh is not copied as _subM.size != _computedSubM.size()
488         list< _ImportData > &  dList = _tgtMesh2ImportData[ subMesh->GetFather() ];
489         list< _ImportData >::iterator d = dList.begin();
490         for ( ; d != dList.end(); ++d )
491           if ( d->_subM.count( subMesh ))
492           {
493             set<SMESH_subMesh*,_SubLess>::iterator smIt = d->_subM.begin();
494             for( ; smIt != d->_subM.end(); ++smIt )
495               if ( (*smIt)->IsMeshComputed() )
496                 d->_computedSubM.insert( *smIt);
497           }
498       }
499     }
500   }
501
502   //================================================================================
503   /*!
504    * \brief Return an ID of submesh to store nodes and elements of a copied mesh
505    */
506   //================================================================================
507
508   int getSubmeshIDForCopiedMesh(const SMESHDS_Mesh* srcMeshDS,
509                                 SMESH_Mesh*         tgtMesh)
510   {
511     // To get SMESH_subMesh corresponding to srcMeshDS we need to have a shape
512     // for which SMESHDS_Mesh::IsGroupOfSubShapes() returns true.
513     // And this shape must be different from sub-shapes of the main shape.
514     // So we create a compound containing
515     // 1) some sub-shapes of SMESH_Mesh::PseudoShape() corresponding to
516     //    srcMeshDS->GetPersistentId()
517     // 2) the 1-st vertex of the main shape to assure
518     //    SMESHDS_Mesh::IsGroupOfSubShapes(shape)==true
519     TopoDS_Shape shapeForSrcMesh;
520     TopTools_IndexedMapOfShape pseudoSubShapes;
521     TopExp::MapShapes( SMESH_Mesh::PseudoShape(), pseudoSubShapes );
522
523     // index of pseudoSubShapes corresponding to srcMeshDS
524     int    subIndex = 1 + srcMeshDS->GetPersistentId() % pseudoSubShapes.Extent();
525     int nbSubShapes = 1 + srcMeshDS->GetPersistentId() / pseudoSubShapes.Extent();
526
527     // try to find already present shapeForSrcMesh
528     SMESHDS_Mesh* tgtMeshDS = tgtMesh->GetMeshDS();
529     for ( int i = tgtMeshDS->MaxShapeIndex(); i > 0 && shapeForSrcMesh.IsNull(); --i )
530     {
531       const TopoDS_Shape& s = tgtMeshDS->IndexToShape(i);
532       if ( s.ShapeType() != TopAbs_COMPOUND ) break;
533       TopoDS_Iterator sSubIt( s );
534       for ( int iSub = 0; iSub < nbSubShapes && sSubIt.More(); ++iSub, sSubIt.Next() )
535         if ( pseudoSubShapes( subIndex+iSub ).IsSame( sSubIt.Value()))
536           if ( iSub+1 == nbSubShapes )
537           {
538             shapeForSrcMesh = s;
539             break;
540           }
541     }
542     if ( shapeForSrcMesh.IsNull() )
543     {
544       // make a new shapeForSrcMesh
545       BRep_Builder aBuilder;
546       TopoDS_Compound comp;
547       aBuilder.MakeCompound( comp );
548       shapeForSrcMesh = comp;
549       for ( int iSub = 0; iSub < nbSubShapes; ++iSub )
550         aBuilder.Add( comp, pseudoSubShapes( subIndex+iSub ));
551       TopExp_Explorer vExp( tgtMeshDS->ShapeToMesh(), TopAbs_VERTEX );
552       aBuilder.Add( comp, vExp.Current() );
553     }
554     SMESH_subMesh* sm = tgtMesh->GetSubMesh( shapeForSrcMesh );
555     SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
556     if ( !smDS )
557       smDS = tgtMeshDS->NewSubMesh( sm->GetId() );
558
559     // make ordinary submesh from a complex one
560     if ( smDS->IsComplexSubmesh() )
561     {
562       list< const SMESHDS_SubMesh* > subSM;
563       SMESHDS_SubMeshIteratorPtr smIt = smDS->GetSubMeshIterator();
564       while ( smIt->more() ) subSM.push_back( smIt->next() );
565       list< const SMESHDS_SubMesh* >::iterator sub = subSM.begin();
566       for ( ; sub != subSM.end(); ++sub)
567         smDS->RemoveSubMesh( *sub );
568     }
569     return sm->GetId();
570   }
571
572   //================================================================================
573   /*!
574    * \brief Return a submesh to store nodes and elements of a copied mesh
575    * and set event listeners in order to clear
576    * imported mesh and groups as soon as submesh state requires it
577    */
578   //================================================================================
579
580   SMESHDS_SubMesh* getSubmeshForCopiedMesh(const SMESH_Mesh*                    srcMesh,
581                                            SMESH_Mesh*                          tgtMesh,
582                                            const TopoDS_Shape&                  tgtShape,
583                                            StdMeshers_Import_1D::TNodeNodeMap*& n2n,
584                                            StdMeshers_Import_1D::TElemElemMap*& e2e,
585                                            bool &                               toCopyGroups)
586   {
587     StdMeshers_Import_1D::getMaps( srcMesh, tgtMesh, n2n,e2e );
588
589     _ImportData* iData = _Listener::getImportData(srcMesh,tgtMesh);
590
591     SMESH_subMesh* importedSM = tgtMesh->GetSubMesh( tgtShape );
592     iData->addComputed( importedSM );
593     if ( iData->_computedSubM.size() != iData->_subM.size() )
594       return 0; // not all submeshes computed yet
595
596     toCopyGroups = !iData->_copyGroupSubM.empty();
597
598     if ( !iData->_copyMeshSubM.empty())
599     {
600       // make submesh to store a copied mesh
601       int smID = getSubmeshIDForCopiedMesh( srcMesh->GetMeshDS(), tgtMesh );
602       SMESHDS_SubMesh* subDS = tgtMesh->GetMeshDS()->NewSubMesh( smID );
603
604       iData->_importMeshSubID = smID;
605       iData->_importMeshSubDS = subDS;
606       return subDS;
607     }
608     return 0;
609   }
610
611 } // namespace
612
613
614 //=============================================================================
615 /*!
616  * Import elements from the other mesh
617  */
618 //=============================================================================
619
620 bool StdMeshers_Import_1D::Compute(SMESH_Mesh & theMesh, const TopoDS_Shape & theShape)
621 {
622   if ( !_sourceHyp ) return false;
623
624   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
625   if ( srcGroups.empty() )
626     return error("Invalid source groups");
627
628   SMESH_MesherHelper helper(theMesh);
629   helper.SetSubShape(theShape);
630   SMESHDS_Mesh* tgtMesh = theMesh.GetMeshDS();
631
632   const TopoDS_Edge& geomEdge = TopoDS::Edge( theShape );
633   const double edgeTol = BRep_Tool::Tolerance( geomEdge );
634   const int shapeID = tgtMesh->ShapeToIndex( geomEdge );
635
636   set<int> subShapeIDs;
637   subShapeIDs.insert( shapeID );
638
639   // get nodes on vertices
640   list < SMESH_TNodeXYZ > vertexNodes;
641   list < SMESH_TNodeXYZ >::iterator vNIt;
642   TopExp_Explorer vExp( theShape, TopAbs_VERTEX );
643   for ( ; vExp.More(); vExp.Next() )
644   {
645     const TopoDS_Vertex& v = TopoDS::Vertex( vExp.Current() );
646     if ( !subShapeIDs.insert( tgtMesh->ShapeToIndex( v )).second )
647       continue; // closed edge
648     const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
649     if ( !n )
650     {
651       _gen->Compute(theMesh,v,/*anUpward=*/true);
652       n = SMESH_Algo::VertexNode( v, tgtMesh );
653       if ( !n ) return false; // very strange
654     }
655     vertexNodes.push_back( SMESH_TNodeXYZ( n ));
656   }
657
658   // import edges from groups
659   TNodeNodeMap* n2n;
660   TElemElemMap* e2e;
661   for ( int iG = 0; iG < srcGroups.size(); ++iG )
662   {
663     const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
664
665     const int meshID = srcGroup->GetMesh()->GetPersistentId();
666     const SMESH_Mesh* srcMesh = GetMeshByPersistentID( meshID );
667     if ( !srcMesh ) continue;
668     getMaps( srcMesh, &theMesh, n2n, e2e );
669
670     SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
671     vector<const SMDS_MeshNode*> newNodes;
672     SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
673     double u = 0;
674     while ( srcElems->more() ) // loop on group contents
675     {
676       const SMDS_MeshElement* edge = srcElems->next();
677       // find or create nodes of a new edge
678       newNodes.resize( edge->NbNodes() );
679       newNodes.back() = 0;
680       SMDS_MeshElement::iterator node = edge->begin_nodes();
681       SMESH_TNodeXYZ a(edge->GetNode(0));
682       // --- define a tolerance relative to the length of an edge
683       double mytol = a.Distance(edge->GetNode(edge->NbNodes()-1))/25;
684       //MESSAGE("mytol = " << mytol);
685       for ( unsigned i = 0; i < newNodes.size(); ++i, ++node )
686       {
687         TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
688         if ( n2nIt->second )
689         {
690           if ( !subShapeIDs.count( n2nIt->second->getshapeId() ))
691             break;
692         }
693         else
694         {
695           // find an existing vertex node
696           double checktol = max(1.E-10, 10*edgeTol*edgeTol);
697           for ( vNIt = vertexNodes.begin(); vNIt != vertexNodes.end(); ++vNIt)
698             if ( vNIt->SquareDistance( *node ) < checktol)
699             {
700               //MESSAGE("SquareDistance " << vNIt->SquareDistance( *node ) << " checktol " << checktol);
701               (*n2nIt).second = vNIt->_node;
702               vertexNodes.erase( vNIt );
703               break;
704             }
705         }
706         if ( !n2nIt->second )
707         {
708           // find out if node lies on theShape
709           tmpNode->setXYZ( (*node)->X(), (*node)->Y(), (*node)->Z());
710           if ( helper.CheckNodeU( geomEdge, tmpNode, u, mytol, /*force=*/true ))
711           {
712             SMDS_MeshNode* newNode = tgtMesh->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
713             n2nIt->second = newNode;
714             tgtMesh->SetNodeOnEdge( newNode, shapeID, u );
715             //MESSAGE("u=" << u);
716           }
717         }
718         if ( !(newNodes[i] = n2nIt->second ))
719           break;
720       }
721       if ( !newNodes.back() )
722       {
723         //MESSAGE("not all nodes of edge lie on theShape");
724         continue; // not all nodes of edge lie on theShape
725       }
726
727       // make a new edge
728       SMDS_MeshElement * newEdge;
729       if ( newNodes.size() == 3 )
730         newEdge = tgtMesh->AddEdge( newNodes[0], newNodes[1], newNodes[2] );
731       else
732         newEdge = tgtMesh->AddEdge( newNodes[0], newNodes[1]);
733       //MESSAGE("add Edge");
734       tgtMesh->SetMeshElementOnShape( newEdge, shapeID );
735       e2e->insert( make_pair( edge, newEdge ));
736     }
737     helper.GetMeshDS()->RemoveNode(tmpNode);
738   }
739   if ( n2n->empty())
740     return error("Empty source groups");
741
742   // check if the whole geom edge is covered by imported segments;
743   // the check consist in passing by segments from one vetrex node to another
744   bool isEdgeMeshed = false;
745   if ( SMESHDS_SubMesh* tgtSM = tgtMesh->MeshElements( theShape ))
746   {
747     const TopoDS_Vertex& v = ( vExp.ReInit(), TopoDS::Vertex( vExp.Current() ));
748     const SMDS_MeshNode* n = SMESH_Algo::VertexNode( v, tgtMesh );
749     const SMDS_MeshElement* seg = 0;
750     SMDS_ElemIteratorPtr segIt = n->GetInverseElementIterator(SMDSAbs_Edge);
751     while ( segIt->more() && !seg )
752       if ( !tgtSM->Contains( seg = segIt->next()))
753         seg = 0;
754     int nbPassedSegs = 0;
755     while ( seg )
756     {
757       ++nbPassedSegs;
758       const SMDS_MeshNode* n2 = seg->GetNode(0);
759       n = ( n2 == n ? seg->GetNode(1) : n2 );
760       if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
761         break;
762       const SMDS_MeshElement* seg2 = 0;
763       segIt = n->GetInverseElementIterator(SMDSAbs_Edge);
764       while ( segIt->more() && !seg2 )
765         if ( seg == ( seg2 = segIt->next()))
766           seg2 = 0;
767       seg = seg2;
768     }
769     if (nbPassedSegs > 0 && tgtSM->NbElements() > nbPassedSegs )
770       return error( "Source elements overlap one another");
771
772     isEdgeMeshed = ( tgtSM->NbElements() == nbPassedSegs &&
773                      n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
774   }
775   if ( !isEdgeMeshed )
776     return error( "Source elements don't cover totally the geometrical edge" );
777
778   // copy meshes
779   vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
780   for ( unsigned i = 0; i < srcMeshes.size(); ++i )
781     importMesh( srcMeshes[i], theMesh, _sourceHyp, theShape );
782
783   return true;
784 }
785
786 //================================================================================
787 /*!
788  * \brief Copy mesh and groups
789  */
790 //================================================================================
791
792 void StdMeshers_Import_1D::importMesh(const SMESH_Mesh*          srcMesh,
793                                       SMESH_Mesh &               tgtMesh,
794                                       StdMeshers_ImportSource1D* srcHyp,
795                                       const TopoDS_Shape&        tgtShape)
796 {
797   // get submesh to store the imported mesh
798   TNodeNodeMap* n2n;
799   TElemElemMap* e2e;
800   bool toCopyGroups;
801   SMESHDS_SubMesh* tgtSubMesh =
802     getSubmeshForCopiedMesh( srcMesh, &tgtMesh, tgtShape, n2n, e2e, toCopyGroups );
803   if ( !tgtSubMesh || tgtSubMesh->NbNodes() + tgtSubMesh->NbElements() > 0 )
804     return; // not to copy srcMeshDS twice
805
806   SMESHDS_Mesh* tgtMeshDS = tgtMesh.GetMeshDS();
807   SMESH_MeshEditor additor( &tgtMesh );
808
809   // 1. Copy mesh
810
811   vector<const SMDS_MeshNode*> newNodes;
812   const SMESHDS_Mesh* srcMeshDS = srcMesh->GetMeshDS();
813   SMDS_ElemIteratorPtr eIt = srcMeshDS->elementsIterator();
814   while ( eIt->more() )
815   {
816     const SMDS_MeshElement* elem = eIt->next();
817     TElemElemMap::iterator e2eIt = e2e->insert( make_pair( elem, (SMDS_MeshElement*)0 )).first;
818     if ( e2eIt->second ) continue; // already copied by Compute()
819     newNodes.resize( elem->NbNodes() );
820     SMDS_MeshElement::iterator node = elem->begin_nodes();
821     for ( unsigned i = 0; i < newNodes.size(); ++i, ++node )
822     {
823       TNodeNodeMap::iterator n2nIt = n2n->insert( make_pair( *node, (SMDS_MeshNode*)0 )).first;
824       if ( !n2nIt->second )
825       {
826         (*n2nIt).second = tgtMeshDS->AddNode( (*node)->X(), (*node)->Y(), (*node)->Z());
827         tgtSubMesh->AddNode( n2nIt->second );
828       }
829       newNodes[i] = n2nIt->second;
830     }
831     const SMDS_MeshElement* newElem =
832       tgtMeshDS->FindElement( newNodes, elem->GetType(), /*noMedium=*/false );
833     if ( !newElem )
834     {
835       newElem = additor.AddElement( newNodes, elem->GetType(), elem->IsPoly());
836       tgtSubMesh->AddElement( newElem );
837     }
838     if ( toCopyGroups )
839       (*e2eIt).second = newElem;
840   }
841   // copy free nodes
842   if ( srcMeshDS->NbNodes() > n2n->size() )
843   {
844     SMDS_NodeIteratorPtr nIt = srcMeshDS->nodesIterator();
845     while( nIt->more() )
846     {
847       const SMDS_MeshNode* node = nIt->next();
848       if ( node->NbInverseElements() == 0 )
849       {
850         const SMDS_MeshNode* newNode = tgtMeshDS->AddNode( node->X(), node->Y(), node->Z());
851         n2n->insert( make_pair( node, newNode ));
852         tgtSubMesh->AddNode( newNode );
853       }
854     }
855   }
856
857   // 2. Copy groups
858
859   vector<SMESH_Group*> resultGroups;
860   if ( toCopyGroups )
861   {
862     // collect names of existing groups to assure uniqueness of group names within a type
863     map< SMDSAbs_ElementType, set<string> > namesByType;
864     SMESH_Mesh::GroupIteratorPtr groupIt = tgtMesh.GetGroups();
865     while ( groupIt->more() )
866     {
867       SMESH_Group* tgtGroup = groupIt->next();
868       namesByType[ tgtGroup->GetGroupDS()->GetType() ].insert( tgtGroup->GetName() );
869     }
870     if (srcMesh)
871     {
872       SMESH_Mesh::GroupIteratorPtr groupIt = srcMesh->GetGroups();
873       while ( groupIt->more() )
874       {
875         SMESH_Group* srcGroup = groupIt->next();
876         SMESHDS_GroupBase* srcGroupDS = srcGroup->GetGroupDS();
877         string name = srcGroup->GetName();
878         int nb = 1;
879         while ( !namesByType[ srcGroupDS->GetType() ].insert( name ).second )
880           name = SMESH_Comment(srcGroup->GetName()) << "_imported_" << nb++;
881         SMESH_Group* newGroup = tgtMesh.AddGroup( srcGroupDS->GetType(), name.c_str(), nb );
882         SMESHDS_Group* newGroupDS = (SMESHDS_Group*)newGroup->GetGroupDS();
883         resultGroups.push_back( newGroup );
884
885         eIt = srcGroupDS->GetElements();
886         if ( srcGroupDS->GetType() == SMDSAbs_Node )
887           while (eIt->more())
888           {
889             TNodeNodeMap::iterator n2nIt = n2n->find((const SMDS_MeshNode*) eIt->next() );
890             if ( n2nIt != n2n->end() && n2nIt->second )
891               newGroupDS->SMDSGroup().Add((*n2nIt).second );
892           }
893         else
894           while (eIt->more())
895           {
896             TElemElemMap::iterator e2eIt = e2e->find( eIt->next() );
897             if ( e2eIt != e2e->end() && e2eIt->second )
898               newGroupDS->SMDSGroup().Add((*e2eIt).second );
899           }
900       }
901     }
902   }
903   n2n->clear();
904   e2e->clear();
905
906   // Remember created groups in order to remove them as soon as the srcHyp is
907   // modified or something other similar happens. This imformation must be persistent,
908   // for that store them in a hypothesis as it stores its values in the file anyway
909   srcHyp->StoreResultGroups( resultGroups, *srcMeshDS, *tgtMeshDS );
910 }
911
912 //=============================================================================
913 /*!
914  * \brief Set needed event listeners and create a submesh for a copied mesh
915  *
916  * This method is called only if a submesh has HYP_OK algo_state.
917  */
918 //=============================================================================
919
920 void StdMeshers_Import_1D::setEventListener(SMESH_subMesh*             subMesh,
921                                             StdMeshers_ImportSource1D* sourceHyp)
922 {
923   if ( sourceHyp )
924   {
925     vector<SMESH_Mesh*> srcMeshes = sourceHyp->GetSourceMeshes();
926     if ( srcMeshes.empty() )
927       _Listener::waitHypModification( subMesh );
928     for ( unsigned i = 0; i < srcMeshes.size(); ++i )
929       // set a listener to remove the imported mesh and groups
930       _Listener::storeImportSubmesh( subMesh, srcMeshes[i], sourceHyp );
931   }
932 }
933 void StdMeshers_Import_1D::SetEventListener(SMESH_subMesh* subMesh)
934 {
935   if ( !_sourceHyp )
936   {
937     const TopoDS_Shape& tgtShape = subMesh->GetSubShape();
938     SMESH_Mesh*         tgtMesh  = subMesh->GetFather();
939     Hypothesis_Status aStatus;
940     CheckHypothesis( *tgtMesh, tgtShape, aStatus );
941   }
942   setEventListener( subMesh, _sourceHyp );
943 }
944
945 void StdMeshers_Import_1D::SubmeshRestored(SMESH_subMesh* subMesh)
946 {
947   SetEventListener(subMesh);
948 }
949
950 //=============================================================================
951 /*!
952  * Predict nb of mesh entities created by Compute()
953  */
954 //=============================================================================
955
956 bool StdMeshers_Import_1D::Evaluate(SMESH_Mesh &         theMesh,
957                                     const TopoDS_Shape & theShape,
958                                     MapShapeNbElems&     aResMap)
959 {
960   if ( !_sourceHyp ) return false;
961
962   const vector<SMESH_Group*>& srcGroups = _sourceHyp->GetGroups();
963   if ( srcGroups.empty() )
964     return error("Invalid source groups");
965
966   vector<int> aVec(SMDSEntity_Last,0);
967
968   bool toCopyMesh, toCopyGroups;
969   _sourceHyp->GetCopySourceMesh(toCopyMesh, toCopyGroups);
970   if ( toCopyMesh ) // the whole mesh is copied
971   {
972     vector<SMESH_Mesh*> srcMeshes = _sourceHyp->GetSourceMeshes();
973     for ( unsigned i = 0; i < srcMeshes.size(); ++i )
974     {
975       SMESH_subMesh* sm = getSubMeshOfCopiedMesh( theMesh, *srcMeshes[i]);
976       if ( !sm || aResMap.count( sm )) continue; // already counted
977       aVec.assign( SMDSEntity_Last, 0);
978       const SMDS_MeshInfo& aMeshInfo = srcMeshes[i]->GetMeshDS()->GetMeshInfo();
979       for (int i = 0; i < SMDSEntity_Last; i++)
980         aVec[i] = aMeshInfo.NbEntities((SMDSAbs_EntityType)i);
981     }
982   }
983   else
984   {
985     SMESH_MesherHelper helper(theMesh);
986
987     const TopoDS_Edge& geomEdge = TopoDS::Edge( theShape );
988     const double edgeTol = helper.MaxTolerance( geomEdge );
989
990     // take into account nodes on vertices
991     TopExp_Explorer vExp( theShape, TopAbs_VERTEX );
992     for ( ; vExp.More(); vExp.Next() )
993       theMesh.GetSubMesh( vExp.Current())->Evaluate( aResMap );
994
995     // count edges imported from groups
996     int nbEdges = 0, nbQuadEdges = 0;
997     for ( int iG = 0; iG < srcGroups.size(); ++iG )
998     {
999       const SMESHDS_GroupBase* srcGroup = srcGroups[iG]->GetGroupDS();
1000       SMDS_ElemIteratorPtr srcElems = srcGroup->GetElements();
1001       SMDS_MeshNode *tmpNode = helper.AddNode(0,0,0);
1002       while ( srcElems->more() ) // loop on group contents
1003       {
1004         const SMDS_MeshElement* edge = srcElems->next();
1005         // find out if edge is located on geomEdge by projecting
1006         // a middle of edge to geomEdge
1007         SMESH_TNodeXYZ p1( edge->GetNode(0));
1008         SMESH_TNodeXYZ p2( edge->GetNode(1));
1009         gp_XYZ middle = ( p1 + p2 ) / 2.;
1010         tmpNode->setXYZ( middle.X(), middle.Y(), middle.Z());
1011         double u = 0;
1012         if ( helper.CheckNodeU( geomEdge, tmpNode, u, 10 * edgeTol, /*force=*/true ))
1013           ++( edge->IsQuadratic() ? nbQuadEdges : nbEdges);
1014       }
1015       helper.GetMeshDS()->RemoveNode(tmpNode);
1016     }
1017
1018     int nbNodes = nbEdges + 2 * nbQuadEdges - 1;
1019
1020     aVec[SMDSEntity_Node     ] = nbNodes;
1021     aVec[SMDSEntity_Edge     ] = nbEdges;
1022     aVec[SMDSEntity_Quad_Edge] = nbQuadEdges;
1023   }
1024
1025   SMESH_subMesh * sm = theMesh.GetSubMesh(theShape);
1026   aResMap.insert(make_pair(sm,aVec));
1027
1028   return true;
1029 }
1030
1031 //================================================================================
1032 /*!
1033  * \brief Return node-node and element-element maps for import of geiven source mesh
1034  */
1035 //================================================================================
1036
1037 void StdMeshers_Import_1D::getMaps(const SMESH_Mesh* srcMesh,
1038                                    SMESH_Mesh*       tgtMesh,
1039                                    TNodeNodeMap*&    n2n,
1040                                    TElemElemMap*&    e2e)
1041 {
1042   _ImportData* iData = _Listener::getImportData(srcMesh,tgtMesh);
1043   n2n = &iData->_n2n;
1044   e2e = &iData->_e2e;
1045   if ( iData->_copyMeshSubM.empty() )
1046   {
1047     n2n->clear();
1048     e2e->clear();
1049   }
1050 }
1051
1052 //================================================================================
1053 /*!
1054  * \brief Return submesh corresponding to the copied mesh
1055  */
1056 //================================================================================
1057
1058 SMESH_subMesh* StdMeshers_Import_1D::getSubMeshOfCopiedMesh( SMESH_Mesh& tgtMesh,
1059                                                              SMESH_Mesh& srcMesh )
1060 {
1061   _ImportData* iData = _Listener::getImportData(&srcMesh,&tgtMesh);
1062   if ( iData->_copyMeshSubM.empty() ) return 0;
1063   SMESH_subMesh* sm = tgtMesh.GetSubMeshContaining( iData->_importMeshSubID );
1064   return sm;
1065 }
1066