Salome HOME
correct previous integration (Porting to Python 2.6)
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_3D.cxx
1 //  Copyright (C) 2007-2008  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 //  SMESH SMESH : implementaion of SMESH idl descriptions
23 // File      : StdMeshers_Projection_3D.cxx
24 // Module    : SMESH
25 // Created   : Fri Oct 20 11:37:07 2006
26 // Author    : Edward AGAPOV (eap)
27 //
28 #include "StdMeshers_Projection_3D.hxx"
29 #include "StdMeshers_ProjectionSource3D.hxx"
30
31 #include "StdMeshers_ProjectionUtils.hxx"
32
33 #include "SMDS_PolyhedralVolumeOfNodes.hxx"
34 #include "SMDS_VolumeTool.hxx"
35 #include "SMESHDS_Hypothesis.hxx"
36 #include "SMESHDS_SubMesh.hxx"
37 #include "SMESH_Block.hxx"
38 #include "SMESH_Comment.hxx"
39 #include "SMESH_Gen.hxx"
40 #include "SMESH_Mesh.hxx"
41 #include "SMESH_MesherHelper.hxx"
42 #include "SMESH_Pattern.hxx"
43 #include "SMESH_subMesh.hxx"
44 #include "SMESH_subMeshEventListener.hxx"
45
46 #include "utilities.h"
47
48 #include <TopExp.hxx>
49 #include <TopExp_Explorer.hxx>
50 #include <TopoDS.hxx>
51
52 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
53 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
54 #define SHOWYXZ(msg, xyz) // {\
55 // gp_Pnt p (xyz); \
56 // cout << msg << " ("<< p.X() << "; " <<p.Y() << "; " <<p.Z() << ") " <<endl;\
57 // }
58
59 typedef StdMeshers_ProjectionUtils TAssocTool;
60
61
62 //=======================================================================
63 //function : StdMeshers_Projection_3D
64 //purpose  : 
65 //=======================================================================
66
67 StdMeshers_Projection_3D::StdMeshers_Projection_3D(int hypId, int studyId, SMESH_Gen* gen)
68   :SMESH_3D_Algo(hypId, studyId, gen)
69 {
70   _name = "Projection_3D";
71   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);  // 1 bit per shape type
72
73   _compatibleHypothesis.push_back("ProjectionSource3D");
74   _sourceHypo = 0;
75 }
76
77 //================================================================================
78 /*!
79  * \brief Destructor
80  */
81 //================================================================================
82
83 StdMeshers_Projection_3D::~StdMeshers_Projection_3D()
84 {}
85
86 //=======================================================================
87 //function : CheckHypothesis
88 //purpose  : 
89 //=======================================================================
90
91 bool StdMeshers_Projection_3D::CheckHypothesis(SMESH_Mesh&                          aMesh,
92                                                const TopoDS_Shape&                  aShape,
93                                                SMESH_Hypothesis::Hypothesis_Status& aStatus)
94 {
95   // check aShape that must be a 6 faces block
96 /*  PAL16229
97   if ( TAssocTool::Count( aShape, TopAbs_SHELL, 1 ) != 1 ||
98        TAssocTool::Count( aShape, TopAbs_FACE , 1 ) != 6 ||
99        TAssocTool::Count( aShape, TopAbs_EDGE , 1 ) != 12 ||
100        TAssocTool::Count( aShape, TopAbs_WIRE , 1 ) != 6 )
101   {
102     aStatus = HYP_BAD_GEOMETRY;
103     return false;
104   }
105 */
106   list <const SMESHDS_Hypothesis * >::const_iterator itl;
107
108   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
109   if ( hyps.size() == 0 )
110   {
111     aStatus = SMESH_Hypothesis::HYP_MISSING;
112     return false;  // can't work with no hypothesis
113   }
114
115   if ( hyps.size() > 1 )
116   {
117     aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
118     return false;
119   }
120
121   const SMESHDS_Hypothesis *theHyp = hyps.front();
122
123   string hypName = theHyp->GetName();
124
125   aStatus = SMESH_Hypothesis::HYP_OK;
126
127   if (hypName == "ProjectionSource3D")
128   {
129     _sourceHypo = static_cast<const StdMeshers_ProjectionSource3D *>(theHyp);
130     // Check hypo parameters
131
132     SMESH_Mesh* srcMesh = _sourceHypo->GetSourceMesh();
133     SMESH_Mesh* tgtMesh = & aMesh;
134     if ( !srcMesh )
135       srcMesh = tgtMesh;
136
137     // check vertices
138     if ( _sourceHypo->HasVertexAssociation() )
139     {
140       // source vertices
141       TopoDS_Shape edge = TAssocTool::GetEdgeByVertices
142         ( srcMesh, _sourceHypo->GetSourceVertex(1), _sourceHypo->GetSourceVertex(2) );
143       if ( edge.IsNull() ||
144            !TAssocTool::IsSubShape( edge, srcMesh ) ||
145            !TAssocTool::IsSubShape( edge, _sourceHypo->GetSource3DShape() ))
146       {
147         SCRUTE((edge.IsNull()));
148         SCRUTE((TAssocTool::IsSubShape( edge, srcMesh )));
149         SCRUTE((TAssocTool::IsSubShape( edge, _sourceHypo->GetSource3DShape() )));
150         aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
151       }
152       else
153       {
154         // target vertices
155         edge = TAssocTool::GetEdgeByVertices
156           ( tgtMesh, _sourceHypo->GetTargetVertex(1), _sourceHypo->GetTargetVertex(2) );
157         if ( edge.IsNull() ||
158              !TAssocTool::IsSubShape( edge, tgtMesh ) ||
159              !TAssocTool::IsSubShape( edge, aShape ))
160         {
161           SCRUTE((edge.IsNull()));
162           SCRUTE((TAssocTool::IsSubShape( edge, tgtMesh )));
163           SCRUTE((TAssocTool::IsSubShape( edge, aShape )));
164           aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
165         }
166       }
167     }
168     // check a source shape
169     if ( !TAssocTool::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh ) ||
170          ( srcMesh == tgtMesh && aShape == _sourceHypo->GetSource3DShape()))
171     {
172       SCRUTE((TAssocTool::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh)));
173       SCRUTE((srcMesh == tgtMesh));
174       SCRUTE((aShape == _sourceHypo->GetSource3DShape()));
175       aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
176     }
177   }
178   else
179   {
180     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
181   }
182   return ( aStatus == HYP_OK );
183 }
184
185 //=======================================================================
186 //function : Compute
187 //purpose  : 
188 //=======================================================================
189
190 bool StdMeshers_Projection_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
191 {
192   if ( !_sourceHypo )
193     return false;
194
195   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
196   SMESH_Mesh * tgtMesh = & aMesh;
197   if ( !srcMesh )
198     srcMesh = tgtMesh;
199
200   SMESHDS_Mesh * srcMeshDS = srcMesh->GetMeshDS();
201   SMESHDS_Mesh * tgtMeshDS = tgtMesh->GetMeshDS();
202
203   // get shell from shape3D
204   TopoDS_Shell srcShell, tgtShell;
205   TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
206   int nbShell;
207   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
208     srcShell = TopoDS::Shell( exp.Current() );
209   if ( nbShell != 1 )
210     return error(COMPERR_BAD_SHAPE,
211                  SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
212
213   exp.Init( aShape, TopAbs_SHELL );
214   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
215     tgtShell = TopoDS::Shell( exp.Current() );
216   if ( nbShell != 1 )
217     return error(COMPERR_BAD_SHAPE,
218                  SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
219
220   // Check that shapes are blocks
221   if ( TAssocTool::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
222        TAssocTool::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
223        TAssocTool::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
224     return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
225   if ( TAssocTool::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
226        TAssocTool::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
227        TAssocTool::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
228     return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
229
230   // Assure that mesh on a source shape is computed
231
232   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
233   //SMESH_subMesh* tgtSubMesh = tgtMesh->GetSubMesh( aShape );
234
235   if ( tgtMesh == srcMesh && !aShape.IsSame( _sourceHypo->GetSource3DShape() )) {
236     if ( !TAssocTool::MakeComputed( srcSubMesh ))
237       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
238   }
239   else {
240     if ( !srcSubMesh->IsMeshComputed() )
241       return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
242   }
243
244   // Find 2 pairs of corresponding vertices
245
246   TopoDS_Vertex tgtV000, tgtV100, srcV000, srcV100;
247   TAssocTool::TShapeShapeMap shape2ShapeMap;
248
249   if ( _sourceHypo->HasVertexAssociation() )
250   {
251     tgtV000 = _sourceHypo->GetTargetVertex(1);
252     tgtV100 = _sourceHypo->GetTargetVertex(2);
253     srcV000 = _sourceHypo->GetSourceVertex(1);
254     srcV100 = _sourceHypo->GetSourceVertex(2);
255   }
256   else
257   {
258     if ( !TAssocTool::FindSubShapeAssociation( tgtShell, tgtMesh, srcShell, srcMesh,
259                                                shape2ShapeMap) )
260       return error(COMPERR_BAD_SHAPE,"Topology of source and target shapes seems different" );
261
262     exp.Init( tgtShell, TopAbs_EDGE );
263     TopExp::Vertices( TopoDS::Edge( exp.Current() ), tgtV000, tgtV100 );
264
265     if ( !shape2ShapeMap.IsBound( tgtV000 ) || !shape2ShapeMap.IsBound( tgtV100 ))
266       return error("Association of subshapes failed" );
267     srcV000 = TopoDS::Vertex( shape2ShapeMap( tgtV000 ));
268     srcV100 = TopoDS::Vertex( shape2ShapeMap( tgtV100 ));
269     if ( !TAssocTool::IsSubShape( srcV000, srcShell ) ||
270          !TAssocTool::IsSubShape( srcV100, srcShell ))
271       return error("Incorrect association of subshapes" );
272   }
273
274   // Load 2 SMESH_Block's with src and tgt shells
275
276   SMESH_Block srcBlock, tgtBlock;
277   TopTools_IndexedMapOfOrientedShape scrShapes, tgtShapes;
278   if ( !tgtBlock.LoadBlockShapes( tgtShell, tgtV000, tgtV100, tgtShapes ))
279     return error(COMPERR_BAD_SHAPE, "Can't detect block subshapes. Not a block?");
280
281   if ( !srcBlock.LoadBlockShapes( srcShell, srcV000, srcV100, scrShapes ))
282     return error(COMPERR_BAD_SHAPE, "Can't detect block subshapes. Not a block?");
283
284   // Find matching nodes of src and tgt shells
285
286   TNodeNodeMap src2tgtNodeMap;
287   for ( int fId = SMESH_Block::ID_FirstF; fId < SMESH_Block::ID_Shell; ++fId )
288   {
289     // Corresponding subshapes
290     TopoDS_Face srcFace = TopoDS::Face( scrShapes( fId ));
291     TopoDS_Face tgtFace = TopoDS::Face( tgtShapes( fId ));
292     if ( _sourceHypo->HasVertexAssociation() ) { // associate face subshapes
293       shape2ShapeMap.Clear();
294       vector< int > edgeIdVec;
295       SMESH_Block::GetFaceEdgesIDs( fId, edgeIdVec );
296       for ( int i = 0; i < edgeIdVec.size(); ++i ) {
297         int eID = edgeIdVec[ i ];
298         shape2ShapeMap.Bind( tgtShapes( eID ), scrShapes( eID ));
299         if ( i < 2 ) {
300           vector< int > vertexIdVec;
301           SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
302           shape2ShapeMap.Bind( tgtShapes( vertexIdVec[0] ), scrShapes( vertexIdVec[0] ));
303           shape2ShapeMap.Bind( tgtShapes( vertexIdVec[1] ), scrShapes( vertexIdVec[1] ));
304         }
305       }
306     }
307     // Find matching nodes of tgt and src faces
308     TNodeNodeMap faceMatchingNodes;
309     if ( ! TAssocTool::FindMatchingNodesOnFaces( srcFace, srcMesh, tgtFace, tgtMesh, 
310                                                  shape2ShapeMap, faceMatchingNodes ))
311       return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Mesh on faces #")
312                    << srcMeshDS->ShapeToIndex( srcFace ) << " and "
313                    << tgtMeshDS->ShapeToIndex( tgtFace ) << " seems different" );
314
315     // put found matching nodes of 2 faces to the global map
316     src2tgtNodeMap.insert( faceMatchingNodes.begin(), faceMatchingNodes.end() );
317   }
318
319   // ------------------
320   // Make mesh
321   // ------------------
322
323   SMDS_VolumeTool volTool;
324   SMESH_MesherHelper helper( *tgtMesh );
325   helper.IsQuadraticSubMesh( aShape );
326
327   SMESHDS_SubMesh* srcSMDS = srcSubMesh->GetSubMeshDS();
328   SMDS_ElemIteratorPtr volIt = srcSMDS->GetElements();
329   while ( volIt->more() ) // loop on source volumes
330   {
331     const SMDS_MeshElement* srcVol = volIt->next();
332     if ( !srcVol || srcVol->GetType() != SMDSAbs_Volume )
333         continue;
334     int nbNodes = srcVol->NbNodes();
335     SMDS_VolumeTool::VolumeType  volType = volTool.GetType( nbNodes );
336     if ( srcVol->IsQuadratic() )
337       nbNodes = volTool.NbCornerNodes( volType );
338
339     // Find or create a new tgt node for each node of a src volume
340
341     vector< const SMDS_MeshNode* > nodes( nbNodes );
342     for ( int i = 0; i < nbNodes; ++i )
343     {
344       const SMDS_MeshNode* srcNode = srcVol->GetNode( i );
345       const SMDS_MeshNode* tgtNode = 0;
346       TNodeNodeMap::iterator sN_tN = src2tgtNodeMap.find( srcNode );
347       if ( sN_tN != src2tgtNodeMap.end() ) // found
348       {
349         tgtNode = sN_tN->second;
350       }
351       else // Create a new tgt node
352       {
353         // compute normalized parameters of source node in srcBlock
354         gp_Pnt srcCoord = gpXYZ( srcNode );
355         gp_XYZ srcParam;
356         if ( !srcBlock.ComputeParameters( srcCoord, srcParam ))
357           return error(SMESH_Comment("Can't compute normalized parameters ")
358                        << "for source node " << srcNode->GetID());
359         // compute coordinates of target node by srcParam
360         gp_XYZ tgtXYZ;
361         if ( !tgtBlock.ShellPoint( srcParam, tgtXYZ ))
362           return error("Can't compute coordinates by normalized parameters");
363         // add node
364         SMDS_MeshNode* newNode = tgtMeshDS->AddNode( tgtXYZ.X(), tgtXYZ.Y(), tgtXYZ.Z() );
365         tgtMeshDS->SetNodeInVolume( newNode, helper.GetSubShapeID() );
366         tgtNode = newNode;
367         src2tgtNodeMap.insert( make_pair( srcNode, tgtNode ));
368       }
369       nodes[ i ] = tgtNode;
370     }
371
372     // Create a new volume
373
374     SMDS_MeshVolume * tgtVol = 0;
375     int id = 0, force3d = false;
376     switch ( volType ) {
377     case SMDS_VolumeTool::TETRA     :
378     case SMDS_VolumeTool::QUAD_TETRA:
379       tgtVol = helper.AddVolume( nodes[0],
380                                  nodes[1],
381                                  nodes[2],
382                                  nodes[3], id, force3d); break;
383     case SMDS_VolumeTool::PYRAM     :
384     case SMDS_VolumeTool::QUAD_PYRAM:
385       tgtVol = helper.AddVolume( nodes[0],
386                                  nodes[1],
387                                  nodes[2],
388                                  nodes[3],
389                                  nodes[4], id, force3d); break;
390     case SMDS_VolumeTool::PENTA     :
391     case SMDS_VolumeTool::QUAD_PENTA:
392       tgtVol = helper.AddVolume( nodes[0],
393                                  nodes[1],
394                                  nodes[2],
395                                  nodes[3],
396                                  nodes[4],
397                                  nodes[5], id, force3d); break;
398     case SMDS_VolumeTool::HEXA      :
399     case SMDS_VolumeTool::QUAD_HEXA :
400       tgtVol = helper.AddVolume( nodes[0],
401                                  nodes[1],
402                                  nodes[2],
403                                  nodes[3],
404                                  nodes[4],
405                                  nodes[5],
406                                  nodes[6],
407                                  nodes[7], id, force3d); break;
408     default: // polyhedron
409       const SMDS_PolyhedralVolumeOfNodes * poly =
410         dynamic_cast<const SMDS_PolyhedralVolumeOfNodes*>( srcVol );
411       if ( !poly )
412         RETURN_BAD_RESULT("Unexpected volume type");
413       tgtVol = tgtMeshDS->AddPolyhedralVolume( nodes, poly->GetQuanities() );
414     }
415     if ( tgtVol ) {
416       tgtMeshDS->SetMeshElementOnShape( tgtVol, helper.GetSubShapeID() );
417     }
418   } // loop on volumes of src shell
419
420   return true;
421 }
422
423
424 //=======================================================================
425 //function : Evaluate
426 //purpose  : 
427 //=======================================================================
428
429 bool StdMeshers_Projection_3D::Evaluate(SMESH_Mesh& aMesh,
430                                         const TopoDS_Shape& aShape,
431                                         MapShapeNbElems& aResMap)
432 {
433   if ( !_sourceHypo )
434     return false;
435
436   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
437   SMESH_Mesh * tgtMesh = & aMesh;
438   if ( !srcMesh )
439     srcMesh = tgtMesh;
440
441   // get shell from shape3D
442   TopoDS_Shell srcShell, tgtShell;
443   TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
444   int nbShell;
445   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
446     srcShell = TopoDS::Shell( exp.Current() );
447   if ( nbShell != 1 )
448     return error(COMPERR_BAD_SHAPE,
449                  SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
450
451   exp.Init( aShape, TopAbs_SHELL );
452   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
453     tgtShell = TopoDS::Shell( exp.Current() );
454   if ( nbShell != 1 )
455     return error(COMPERR_BAD_SHAPE,
456                  SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
457
458   // Check that shapes are blocks
459   if ( TAssocTool::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
460        TAssocTool::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
461        TAssocTool::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
462     return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
463   if ( TAssocTool::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
464        TAssocTool::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
465        TAssocTool::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
466     return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
467
468   // Assure that mesh on a source shape is computed
469
470   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
471
472   if ( !srcSubMesh->IsMeshComputed() )
473     return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
474
475
476   std::vector<int> aVec(SMDSEntity_Last);
477   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
478
479   aVec[SMDSEntity_Node] = srcSubMesh->GetSubMeshDS()->NbNodes();
480
481   //bool quadratic = false;
482   SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
483   while ( elemIt->more() ) {
484     const SMDS_MeshElement* E  = elemIt->next();
485     if( E->NbNodes()==4 ) {
486       aVec[SMDSEntity_Tetra]++;
487     }
488     else if( E->NbNodes()==5 ) {
489       aVec[SMDSEntity_Pyramid]++;
490     }
491     else if( E->NbNodes()==6 ) {
492       aVec[SMDSEntity_Penta]++;
493     }
494     else if( E->NbNodes()==8 ) {
495       aVec[SMDSEntity_Hexa]++;
496     }
497     else if( E->NbNodes()==10 && E->IsQuadratic() ) {
498       aVec[SMDSEntity_Quad_Tetra]++;
499     }
500     else if( E->NbNodes()==13 && E->IsQuadratic() ) {
501       aVec[SMDSEntity_Quad_Pyramid]++;
502     }
503     else if( E->NbNodes()==15 && E->IsQuadratic() ) {
504       aVec[SMDSEntity_Quad_Penta]++;
505     }
506     else if( E->NbNodes()==20 && E->IsQuadratic() ) {
507       aVec[SMDSEntity_Quad_Hexa]++;
508     }
509     else {
510       aVec[SMDSEntity_Polyhedra]++;
511     }
512   }
513
514   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
515   aResMap.insert(std::make_pair(sm,aVec));
516
517   return true;
518 }
519
520
521 //=============================================================================
522 /*!
523  * \brief Sets a default event listener to submesh of the source shape
524   * \param subMesh - submesh where algo is set
525  *
526  * This method is called when a submesh gets HYP_OK algo_state.
527  * After being set, event listener is notified on each event of a submesh.
528  * Arranges that CLEAN event is translated from source submesh to
529  * the submesh
530  */
531 //=============================================================================
532
533 void StdMeshers_Projection_3D::SetEventListener(SMESH_subMesh* subMesh)
534 {
535   TAssocTool::SetEventListener( subMesh,
536                                 _sourceHypo->GetSource3DShape(),
537                                 _sourceHypo->GetSourceMesh() );
538 }
539