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