Salome HOME
23076: [CEA 1499] Get in python all sub-shapes in error after Compute
[modules/smesh.git] / src / StdMeshers / StdMeshers_Projection_3D.cxx
1 // Copyright (C) 2007-2015  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 : 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 namespace TAssocTool = StdMeshers_ProjectionUtils;
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 ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
223        SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
224        SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
225     return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
226   if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
227        SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
228        SMESH_MesherHelper::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   string srcMeshError;
237   if ( tgtMesh == srcMesh  && !aShape.IsSame( _sourceHypo->GetSource3DShape() )) {
238     if ( !TAssocTool::MakeComputed( srcSubMesh ))
239       srcMeshError = TAssocTool::SourceNotComputedError( srcSubMesh, this );
240   }
241   else {
242     if ( !srcSubMesh->IsMeshComputed() )
243       srcMeshError = TAssocTool::SourceNotComputedError();
244   }
245
246   // Find 2 pairs of corresponding vertices
247
248   TopoDS_Vertex tgtV000, tgtV100, srcV000, srcV100;
249   TAssocTool::TShapeShapeMap shape2ShapeMap;
250
251   if ( _sourceHypo->HasVertexAssociation() )
252   {
253     tgtV000 = _sourceHypo->GetTargetVertex(1);
254     tgtV100 = _sourceHypo->GetTargetVertex(2);
255     srcV000 = _sourceHypo->GetSourceVertex(1);
256     srcV100 = _sourceHypo->GetSourceVertex(2);
257   }
258   else
259   {
260     if ( !TAssocTool::FindSubShapeAssociation( tgtShell, tgtMesh, srcShell, srcMesh,
261                                                shape2ShapeMap) )
262       return error(COMPERR_BAD_SHAPE,"Topology of source and target shapes seems different" );
263
264     exp.Init( tgtShell, TopAbs_EDGE );
265     TopExp::Vertices( TopoDS::Edge( exp.Current() ), tgtV000, tgtV100 );
266
267     if ( !shape2ShapeMap.IsBound( tgtV000 ) || !shape2ShapeMap.IsBound( tgtV100 ))
268       return error("Association of sub-shapes failed" );
269     srcV000 = TopoDS::Vertex( shape2ShapeMap( tgtV000 ));
270     srcV100 = TopoDS::Vertex( shape2ShapeMap( tgtV100 ));
271     if ( !SMESH_MesherHelper::IsSubShape( srcV000, srcShell ) ||
272          !SMESH_MesherHelper::IsSubShape( srcV100, srcShell ))
273       return error("Incorrect association of sub-shapes" );
274   }
275
276   // Load 2 SMESH_Block's with src and tgt shells
277
278   SMESH_Block srcBlock, tgtBlock;
279   TopTools_IndexedMapOfOrientedShape scrShapes, tgtShapes;
280   if ( !tgtBlock.LoadBlockShapes( tgtShell, tgtV000, tgtV100, tgtShapes ))
281     return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. Not a block?");
282
283   if ( !srcBlock.LoadBlockShapes( srcShell, srcV000, srcV100, scrShapes ))
284     return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. Not a block?");
285
286   // Find matching nodes of src and tgt shells
287
288   TNodeNodeMap src2tgtNodeMap;
289   for ( int fId = SMESH_Block::ID_FirstF; fId < SMESH_Block::ID_Shell; ++fId )
290   {
291     // Corresponding sub-shapes
292     TopoDS_Face srcFace = TopoDS::Face( scrShapes( fId ));
293     TopoDS_Face tgtFace = TopoDS::Face( tgtShapes( fId ));
294     if ( _sourceHypo->HasVertexAssociation() ) { // associate face sub-shapes
295       shape2ShapeMap.Clear();
296       vector< int > edgeIdVec;
297       SMESH_Block::GetFaceEdgesIDs( fId, edgeIdVec );
298       for ( int i = 0; i < edgeIdVec.size(); ++i ) {
299         int eID = edgeIdVec[ i ];
300         shape2ShapeMap.Bind( scrShapes( eID ), tgtShapes( eID ));
301         if ( i < 2 ) {
302           vector< int > vertexIdVec;
303           SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
304           shape2ShapeMap.Bind( scrShapes( vertexIdVec[0]), tgtShapes( vertexIdVec[0]) );
305           shape2ShapeMap.Bind( scrShapes( vertexIdVec[1]), tgtShapes( vertexIdVec[1]) );
306         }
307       }
308     }
309     // Find matching nodes of tgt and src faces
310     TAssocTool::TNodeNodeMap faceMatchingNodes;
311     if ( ! TAssocTool::FindMatchingNodesOnFaces( srcFace, srcMesh, tgtFace, tgtMesh, 
312                                                  shape2ShapeMap, faceMatchingNodes ))
313       return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Mesh on faces #")
314                    << srcMeshDS->ShapeToIndex( srcFace ) << " and "
315                    << tgtMeshDS->ShapeToIndex( tgtFace ) << " seems different" );
316
317     // put found matching nodes of 2 faces to the global map
318     src2tgtNodeMap.insert( faceMatchingNodes.begin(), faceMatchingNodes.end() );
319   }
320
321   // ------------------
322   // Make mesh
323   // ------------------
324
325   SMDS_VolumeTool volTool;
326   SMESH_MesherHelper helper( *tgtMesh );
327   helper.IsQuadraticSubMesh( aShape );
328
329   SMESHDS_SubMesh* srcSMDS = srcSubMesh->GetSubMeshDS();
330   SMDS_ElemIteratorPtr volIt = srcSMDS->GetElements();
331   while ( volIt->more() ) // loop on source volumes
332   {
333     const SMDS_MeshElement* srcVol = volIt->next();
334     if ( !srcVol || srcVol->GetType() != SMDSAbs_Volume )
335         continue;
336     int nbNodes = srcVol->NbNodes();
337     SMDS_VolumeTool::VolumeType  volType = volTool.GetType( nbNodes );
338     if ( srcVol->IsQuadratic() )
339       nbNodes = volTool.NbCornerNodes( volType );
340
341     // Find or create a new tgt node for each node of a src volume
342
343     vector< const SMDS_MeshNode* > nodes( nbNodes );
344     for ( int i = 0; i < nbNodes; ++i )
345     {
346       const SMDS_MeshNode* srcNode = srcVol->GetNode( i );
347       const SMDS_MeshNode* tgtNode = 0;
348       TNodeNodeMap::iterator sN_tN = src2tgtNodeMap.find( srcNode );
349       if ( sN_tN != src2tgtNodeMap.end() ) // found
350       {
351         tgtNode = sN_tN->second;
352       }
353       else // Create a new tgt node
354       {
355         // compute normalized parameters of source node in srcBlock
356         gp_Pnt srcCoord = gpXYZ( srcNode );
357         gp_XYZ srcParam;
358         if ( !srcBlock.ComputeParameters( srcCoord, srcParam ))
359           return error(SMESH_Comment("Can't compute normalized parameters ")
360                        << "for source node " << srcNode->GetID());
361         // compute coordinates of target node by srcParam
362         gp_XYZ tgtXYZ;
363         if ( !tgtBlock.ShellPoint( srcParam, tgtXYZ ))
364           return error("Can't compute coordinates by normalized parameters");
365         // add node
366         SMDS_MeshNode* newNode = tgtMeshDS->AddNode( tgtXYZ.X(), tgtXYZ.Y(), tgtXYZ.Z() );
367         tgtMeshDS->SetNodeInVolume( newNode, helper.GetSubShapeID() );
368         tgtNode = newNode;
369         src2tgtNodeMap.insert( make_pair( srcNode, tgtNode ));
370       }
371       nodes[ i ] = tgtNode;
372     }
373
374     // Create a new volume
375
376     SMDS_MeshVolume * tgtVol = 0;
377     int id = 0, force3d = false;
378     switch ( volType ) {
379     case SMDS_VolumeTool::TETRA     :
380     case SMDS_VolumeTool::QUAD_TETRA:
381       tgtVol = helper.AddVolume( nodes[0],
382                                  nodes[1],
383                                  nodes[2],
384                                  nodes[3], id, force3d); break;
385     case SMDS_VolumeTool::PYRAM     :
386     case SMDS_VolumeTool::QUAD_PYRAM:
387       tgtVol = helper.AddVolume( nodes[0],
388                                  nodes[1],
389                                  nodes[2],
390                                  nodes[3],
391                                  nodes[4], id, force3d); break;
392     case SMDS_VolumeTool::PENTA     :
393     case SMDS_VolumeTool::QUAD_PENTA:
394       tgtVol = helper.AddVolume( nodes[0],
395                                  nodes[1],
396                                  nodes[2],
397                                  nodes[3],
398                                  nodes[4],
399                                  nodes[5], id, force3d); break;
400     case SMDS_VolumeTool::HEXA      :
401     case SMDS_VolumeTool::QUAD_HEXA :
402       tgtVol = helper.AddVolume( nodes[0],
403                                  nodes[1],
404                                  nodes[2],
405                                  nodes[3],
406                                  nodes[4],
407                                  nodes[5],
408                                  nodes[6],
409                                  nodes[7], id, force3d); break;
410     default: // polyhedron
411       const SMDS_VtkVolume * poly =
412         dynamic_cast<const SMDS_VtkVolume*>( srcVol );
413       if ( !poly )
414         RETURN_BAD_RESULT("Unexpected volume type");
415       if ( !poly->IsPoly())
416         RETURN_BAD_RESULT("Unexpected volume type");
417       tgtVol = tgtMeshDS->AddPolyhedralVolume( nodes, poly->GetQuantities() );
418     }
419     if ( tgtVol ) {
420       tgtMeshDS->SetMeshElementOnShape( tgtVol, helper.GetSubShapeID() );
421     }
422   } // loop on volumes of src shell
423
424   return true;
425 }
426
427
428 //=======================================================================
429 //function : Evaluate
430 //purpose  : 
431 //=======================================================================
432
433 bool StdMeshers_Projection_3D::Evaluate(SMESH_Mesh& aMesh,
434                                         const TopoDS_Shape& aShape,
435                                         MapShapeNbElems& aResMap)
436 {
437   if ( !_sourceHypo )
438     return false;
439
440   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
441   SMESH_Mesh * tgtMesh = & aMesh;
442   if ( !srcMesh )
443     srcMesh = tgtMesh;
444
445   // get shell from shape3D
446   TopoDS_Shell srcShell, tgtShell;
447   TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
448   int nbShell;
449   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
450     srcShell = TopoDS::Shell( exp.Current() );
451   if ( nbShell != 1 )
452     return error(COMPERR_BAD_SHAPE,
453                  SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
454
455   exp.Init( aShape, TopAbs_SHELL );
456   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
457     tgtShell = TopoDS::Shell( exp.Current() );
458   if ( nbShell != 1 )
459     return error(COMPERR_BAD_SHAPE,
460                  SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
461
462   // Check that shapes are blocks
463   if ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
464        SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
465        SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
466     return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
467   if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
468        SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
469        SMESH_MesherHelper::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
470     return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
471
472   // Assure that mesh on a source shape is computed
473
474   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
475
476   if ( !srcSubMesh->IsMeshComputed() )
477     return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
478
479
480   std::vector<int> aVec(SMDSEntity_Last);
481   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
482
483   aVec[SMDSEntity_Node] = srcSubMesh->GetSubMeshDS()->NbNodes();
484
485   //bool quadratic = false;
486   SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
487   while ( elemIt->more() ) {
488     const SMDS_MeshElement* E  = elemIt->next();
489     if( E->NbNodes()==4 ) {
490       aVec[SMDSEntity_Tetra]++;
491     }
492     else if( E->NbNodes()==5 ) {
493       aVec[SMDSEntity_Pyramid]++;
494     }
495     else if( E->NbNodes()==6 ) {
496       aVec[SMDSEntity_Penta]++;
497     }
498     else if( E->NbNodes()==8 ) {
499       aVec[SMDSEntity_Hexa]++;
500     }
501     else if( E->NbNodes()==10 && E->IsQuadratic() ) {
502       aVec[SMDSEntity_Quad_Tetra]++;
503     }
504     else if( E->NbNodes()==13 && E->IsQuadratic() ) {
505       aVec[SMDSEntity_Quad_Pyramid]++;
506     }
507     else if( E->NbNodes()==15 && E->IsQuadratic() ) {
508       aVec[SMDSEntity_Quad_Penta]++;
509     }
510     else if( E->NbNodes()==20 && E->IsQuadratic() ) {
511       aVec[SMDSEntity_Quad_Hexa]++;
512     }
513     else {
514       aVec[SMDSEntity_Polyhedra]++;
515     }
516   }
517
518   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
519   aResMap.insert(std::make_pair(sm,aVec));
520
521   return true;
522 }
523
524
525 //=============================================================================
526 /*!
527  * \brief Sets a default event listener to submesh of the source shape
528   * \param subMesh - submesh where algo is set
529  *
530  * This method is called when a submesh gets HYP_OK algo_state.
531  * After being set, event listener is notified on each event of a submesh.
532  * Arranges that CLEAN event is translated from source submesh to
533  * the submesh
534  */
535 //=============================================================================
536
537 void StdMeshers_Projection_3D::SetEventListener(SMESH_subMesh* subMesh)
538 {
539   TAssocTool::SetEventListener( subMesh,
540                                 _sourceHypo->GetSource3DShape(),
541                                 _sourceHypo->GetSourceMesh() );
542 }
543   
544 //================================================================================
545 /*!
546  * \brief Return true if the algorithm can mesh this shape
547  *  \param [in] aShape - shape to check
548  *  \param [in] toCheckAll - if true, this check returns OK if all shapes are OK,
549  *              else, returns OK if at least one shape is OK
550  */
551 //================================================================================
552
553 bool StdMeshers_Projection_3D::IsApplicable(const TopoDS_Shape & aShape, bool toCheckAll)
554 {
555   TopExp_Explorer exp0( aShape, TopAbs_SOLID );
556   if ( !exp0.More() ) return false;
557
558   TopTools_IndexedMapOfOrientedShape blockShapes;
559   TopoDS_Vertex v;
560   TopoDS_Shell shell;
561   for ( ; exp0.More(); exp0.Next() )
562   {
563     int nbFoundShells = 0;
564     TopExp_Explorer exp1( exp0.Current(), TopAbs_SHELL );
565     for ( ; exp1.More(); exp1.Next(), ++nbFoundShells )
566     {
567       shell = TopoDS::Shell( exp1.Current() );
568       if ( nbFoundShells == 2 ) break;
569     }
570     if ( nbFoundShells != 1 ) {
571       if ( toCheckAll ) return false;
572       continue;
573     }   
574     bool isBlock = SMESH_Block::FindBlockShapes( shell, v, v, blockShapes );
575     if ( toCheckAll && !isBlock ) return false;
576     if ( !toCheckAll && isBlock ) return true;
577   }
578   return toCheckAll;
579 }