Salome HOME
d8bce83cd0f1c38568cee8cfb8d72067e6736b6d
[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); cout<<msg<< " ("<< p.X() << "; " <<p.Y() << "; " <<p.Z() << ") " <<endl; }
57
58 namespace TAssocTool = StdMeshers_ProjectionUtils;
59
60
61 //=======================================================================
62 //function : StdMeshers_Projection_3D
63 //purpose  : 
64 //=======================================================================
65
66 StdMeshers_Projection_3D::StdMeshers_Projection_3D(int hypId, int studyId, SMESH_Gen* gen)
67   :SMESH_3D_Algo(hypId, studyId, gen)
68 {
69   _name = "Projection_3D";
70   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);  // 1 bit per shape type
71
72   _compatibleHypothesis.push_back("ProjectionSource3D");
73   _sourceHypo = 0;
74 }
75
76 //================================================================================
77 /*!
78  * \brief Destructor
79  */
80 //================================================================================
81
82 StdMeshers_Projection_3D::~StdMeshers_Projection_3D()
83 {}
84
85 //=======================================================================
86 //function : CheckHypothesis
87 //purpose  : 
88 //=======================================================================
89
90 bool StdMeshers_Projection_3D::CheckHypothesis(SMESH_Mesh&                          aMesh,
91                                                const TopoDS_Shape&                  aShape,
92                                                SMESH_Hypothesis::Hypothesis_Status& aStatus)
93 {
94   // check aShape that must be a 6 faces block
95 /*  PAL16229
96   if ( TAssocTool::Count( aShape, TopAbs_SHELL, 1 ) != 1 ||
97        TAssocTool::Count( aShape, TopAbs_FACE , 1 ) != 6 ||
98        TAssocTool::Count( aShape, TopAbs_EDGE , 1 ) != 12 ||
99        TAssocTool::Count( aShape, TopAbs_WIRE , 1 ) != 6 )
100   {
101     aStatus = HYP_BAD_GEOMETRY;
102     return false;
103   }
104 */
105   list <const SMESHDS_Hypothesis * >::const_iterator itl;
106
107   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
108   if ( hyps.size() == 0 )
109   {
110     aStatus = SMESH_Hypothesis::HYP_MISSING;
111     return false;  // can't work with no hypothesis
112   }
113
114   if ( hyps.size() > 1 )
115   {
116     aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
117     return false;
118   }
119
120   const SMESHDS_Hypothesis *theHyp = hyps.front();
121
122   string hypName = theHyp->GetName();
123
124   aStatus = SMESH_Hypothesis::HYP_OK;
125
126   if (hypName == "ProjectionSource3D")
127   {
128     _sourceHypo = static_cast<const StdMeshers_ProjectionSource3D *>(theHyp);
129     // Check hypo parameters
130
131     SMESH_Mesh* srcMesh = _sourceHypo->GetSourceMesh();
132     SMESH_Mesh* tgtMesh = & aMesh;
133     if ( !srcMesh )
134       srcMesh = tgtMesh;
135
136     // check vertices
137     if ( _sourceHypo->HasVertexAssociation() )
138     {
139       // source vertices
140       TopoDS_Shape edge = TAssocTool::GetEdgeByVertices
141         ( srcMesh, _sourceHypo->GetSourceVertex(1), _sourceHypo->GetSourceVertex(2) );
142       if ( edge.IsNull() ||
143            !SMESH_MesherHelper::IsSubShape( edge, srcMesh ) ||
144            !SMESH_MesherHelper::IsSubShape( edge, _sourceHypo->GetSource3DShape() ))
145       {
146         SCRUTE((edge.IsNull()));
147         SCRUTE((SMESH_MesherHelper::IsSubShape( edge, srcMesh )));
148         SCRUTE((SMESH_MesherHelper::IsSubShape( edge, _sourceHypo->GetSource3DShape() )));
149         aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
150       }
151       else
152       {
153         // target vertices
154         edge = TAssocTool::GetEdgeByVertices
155           ( tgtMesh, _sourceHypo->GetTargetVertex(1), _sourceHypo->GetTargetVertex(2) );
156         if ( edge.IsNull() ||
157              !SMESH_MesherHelper::IsSubShape( edge, tgtMesh ) ||
158              !SMESH_MesherHelper::IsSubShape( edge, aShape ))
159         {
160           SCRUTE((edge.IsNull()));
161           SCRUTE((SMESH_MesherHelper::IsSubShape( edge, tgtMesh )));
162           SCRUTE((SMESH_MesherHelper::IsSubShape( edge, aShape )));
163           aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
164         }
165       }
166     }
167     // check a source shape
168     if ( !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh ) ||
169          ( srcMesh == tgtMesh && aShape == _sourceHypo->GetSource3DShape()))
170     {
171       SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh)));
172       SCRUTE((srcMesh == tgtMesh));
173       SCRUTE((aShape == _sourceHypo->GetSource3DShape()));
174       aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
175     }
176   }
177   else
178   {
179     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
180   }
181   return ( aStatus == HYP_OK );
182 }
183
184 //=======================================================================
185 //function : Compute
186 //purpose  : 
187 //=======================================================================
188
189 bool StdMeshers_Projection_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
190 {
191   if ( !_sourceHypo )
192     return false;
193
194   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
195   SMESH_Mesh * tgtMesh = & aMesh;
196   if ( !srcMesh )
197     srcMesh = tgtMesh;
198
199   SMESHDS_Mesh * srcMeshDS = srcMesh->GetMeshDS();
200   SMESHDS_Mesh * tgtMeshDS = tgtMesh->GetMeshDS();
201
202   // get shell from shape3D
203   TopoDS_Shell srcShell, tgtShell;
204   TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
205   int nbShell;
206   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
207     srcShell = TopoDS::Shell( exp.Current() );
208   if ( nbShell != 1 )
209     return error(COMPERR_BAD_SHAPE,
210                  SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
211
212   exp.Init( aShape, TopAbs_SHELL );
213   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
214     tgtShell = TopoDS::Shell( exp.Current() );
215   if ( nbShell != 1 )
216     return error(COMPERR_BAD_SHAPE,
217                  SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
218
219   // Check that shapes are blocks
220   if ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
221        SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
222        SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
223     return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
224   if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
225        SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
226        SMESH_MesherHelper::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
227     return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
228
229   // Assure that mesh on a source shape is computed
230
231   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
232   //SMESH_subMesh* tgtSubMesh = tgtMesh->GetSubMesh( aShape );
233
234   string srcMeshError;
235   if ( tgtMesh == srcMesh  && !aShape.IsSame( _sourceHypo->GetSource3DShape() )) {
236     if ( !TAssocTool::MakeComputed( srcSubMesh ))
237       srcMeshError = TAssocTool::SourceNotComputedError( srcSubMesh, this );
238   }
239   else {
240     if ( !srcSubMesh->IsMeshComputed() )
241       srcMeshError = TAssocTool::SourceNotComputedError();
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 sub-shapes failed" );
267     srcV000 = TopoDS::Vertex( shape2ShapeMap( tgtV000 ));
268     srcV100 = TopoDS::Vertex( shape2ShapeMap( tgtV100 ));
269     if ( !SMESH_MesherHelper::IsSubShape( srcV000, srcShell ) ||
270          !SMESH_MesherHelper::IsSubShape( srcV100, srcShell ))
271       return error("Incorrect association of sub-shapes" );
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 sub-shapes. Not a block?");
280
281   if ( !srcBlock.LoadBlockShapes( srcShell, srcV000, srcV100, scrShapes ))
282     return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. 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 sub-shapes
290     TopoDS_Face srcFace = TopoDS::Face( scrShapes( fId ));
291     TopoDS_Face tgtFace = TopoDS::Face( tgtShapes( fId ));
292     if ( _sourceHypo->HasVertexAssociation() ) { // associate face sub-shapes
293       shape2ShapeMap.Clear();
294       vector< int > edgeIdVec;
295       SMESH_Block::GetFaceEdgesIDs( fId, edgeIdVec );
296       for ( size_t i = 0; i < edgeIdVec.size(); ++i ) {
297         int eID = edgeIdVec[ i ];
298         shape2ShapeMap.Bind( scrShapes( eID ), tgtShapes( eID ));
299         if ( i < 2 ) {
300           vector< int > vertexIdVec;
301           SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
302           shape2ShapeMap.Bind( scrShapes( vertexIdVec[0]), tgtShapes( vertexIdVec[0]) );
303           shape2ShapeMap.Bind( scrShapes( vertexIdVec[1]), tgtShapes( vertexIdVec[1]) );
304         }
305       }
306     }
307     // Find matching nodes of tgt and src faces
308     TAssocTool::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_VtkVolume * poly =
410         dynamic_cast<const SMDS_VtkVolume*>( srcVol );
411       if ( !poly )
412         RETURN_BAD_RESULT("Unexpected volume type");
413       if ( !poly->IsPoly())
414         RETURN_BAD_RESULT("Unexpected volume type");
415       tgtVol = tgtMeshDS->AddPolyhedralVolume( nodes, poly->GetQuantities() );
416     }
417     if ( tgtVol ) {
418       tgtMeshDS->SetMeshElementOnShape( tgtVol, helper.GetSubShapeID() );
419     }
420   } // loop on volumes of src shell
421
422   return true;
423 }
424
425
426 //=======================================================================
427 //function : Evaluate
428 //purpose  : 
429 //=======================================================================
430
431 bool StdMeshers_Projection_3D::Evaluate(SMESH_Mesh& aMesh,
432                                         const TopoDS_Shape& aShape,
433                                         MapShapeNbElems& aResMap)
434 {
435   if ( !_sourceHypo )
436     return false;
437
438   SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
439   SMESH_Mesh * tgtMesh = & aMesh;
440   if ( !srcMesh )
441     srcMesh = tgtMesh;
442
443   // get shell from shape3D
444   TopoDS_Shell srcShell, tgtShell;
445   TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
446   int nbShell;
447   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
448     srcShell = TopoDS::Shell( exp.Current() );
449   if ( nbShell != 1 )
450     return error(COMPERR_BAD_SHAPE,
451                  SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
452
453   exp.Init( aShape, TopAbs_SHELL );
454   for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
455     tgtShell = TopoDS::Shell( exp.Current() );
456   if ( nbShell != 1 )
457     return error(COMPERR_BAD_SHAPE,
458                  SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
459
460   // Check that shapes are blocks
461   if ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
462        SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
463        SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
464     return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
465   if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
466        SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
467        SMESH_MesherHelper::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
468     return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
469
470   // Assure that mesh on a source shape is computed
471
472   SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
473
474   if ( !srcSubMesh->IsMeshComputed() )
475     return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
476
477
478   std::vector<int> aVec(SMDSEntity_Last);
479   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
480
481   aVec[SMDSEntity_Node] = srcSubMesh->GetSubMeshDS()->NbNodes();
482
483   //bool quadratic = false;
484   SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
485   while ( elemIt->more() ) {
486     const SMDS_MeshElement* E  = elemIt->next();
487     if( E->NbNodes()==4 ) {
488       aVec[SMDSEntity_Tetra]++;
489     }
490     else if( E->NbNodes()==5 ) {
491       aVec[SMDSEntity_Pyramid]++;
492     }
493     else if( E->NbNodes()==6 ) {
494       aVec[SMDSEntity_Penta]++;
495     }
496     else if( E->NbNodes()==8 ) {
497       aVec[SMDSEntity_Hexa]++;
498     }
499     else if( E->NbNodes()==10 && E->IsQuadratic() ) {
500       aVec[SMDSEntity_Quad_Tetra]++;
501     }
502     else if( E->NbNodes()==13 && E->IsQuadratic() ) {
503       aVec[SMDSEntity_Quad_Pyramid]++;
504     }
505     else if( E->NbNodes()==15 && E->IsQuadratic() ) {
506       aVec[SMDSEntity_Quad_Penta]++;
507     }
508     else if( E->NbNodes()==20 && E->IsQuadratic() ) {
509       aVec[SMDSEntity_Quad_Hexa]++;
510     }
511     else {
512       aVec[SMDSEntity_Polyhedra]++;
513     }
514   }
515
516   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
517   aResMap.insert(std::make_pair(sm,aVec));
518
519   return true;
520 }
521
522
523 //=============================================================================
524 /*!
525  * \brief Sets a default event listener to submesh of the source shape
526   * \param subMesh - submesh where algo is set
527  *
528  * This method is called when a submesh gets HYP_OK algo_state.
529  * After being set, event listener is notified on each event of a submesh.
530  * Arranges that CLEAN event is translated from source submesh to
531  * the submesh
532  */
533 //=============================================================================
534
535 void StdMeshers_Projection_3D::SetEventListener(SMESH_subMesh* subMesh)
536 {
537   TAssocTool::SetEventListener( subMesh,
538                                 _sourceHypo->GetSource3DShape(),
539                                 _sourceHypo->GetSourceMesh() );
540 }
541   
542 //================================================================================
543 /*!
544  * \brief Return true if the algorithm can mesh this shape
545  *  \param [in] aShape - shape to check
546  *  \param [in] toCheckAll - if true, this check returns OK if all shapes are OK,
547  *              else, returns OK if at least one shape is OK
548  */
549 //================================================================================
550
551 bool StdMeshers_Projection_3D::IsApplicable(const TopoDS_Shape & aShape, bool toCheckAll)
552 {
553   TopExp_Explorer exp0( aShape, TopAbs_SOLID );
554   if ( !exp0.More() ) return false;
555
556   TopTools_IndexedMapOfOrientedShape blockShapes;
557   TopoDS_Vertex v;
558   TopoDS_Shell shell;
559   for ( ; exp0.More(); exp0.Next() )
560   {
561     int nbFoundShells = 0;
562     TopExp_Explorer exp1( exp0.Current(), TopAbs_SHELL );
563     for ( ; exp1.More(); exp1.Next(), ++nbFoundShells )
564     {
565       shell = TopoDS::Shell( exp1.Current() );
566       if ( nbFoundShells == 2 ) break;
567     }
568     if ( nbFoundShells != 1 ) {
569       if ( toCheckAll ) return false;
570       continue;
571     }   
572     bool isBlock = SMESH_Block::FindBlockShapes( shell, v, v, blockShapes );
573     if ( toCheckAll && !isBlock ) return false;
574     if ( !toCheckAll && isBlock ) return true;
575   }
576   return toCheckAll;
577 }