Salome HOME
Fix polyhedron on penta
[modules/smesh.git] / src / StdMeshers / StdMeshers_PolyhedronPerSolid_3D.cxx
1 // Copyright (C) 2007-2020  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 // File      : StdMeshers_PolyhedronPerSolid_3D.cxx
24 // Module    : SMESH
25 // Created   : Fri Oct 20 11:37:07 2006
26 // Author    : Edward AGAPOV (eap)
27 //
28 #include "StdMeshers_PolyhedronPerSolid_3D.hxx"
29
30 #include "SMESHDS_Mesh.hxx"
31 #include "SMESH_ControlsDef.hxx"
32 #include "SMESH_Gen.hxx"
33 #include "SMESH_Mesh.hxx"
34 #include "SMESH_MeshAlgos.hxx"
35 #include "SMESH_MeshEditor.hxx"
36 #include "SMESH_MesherHelper.hxx"
37 #include "SMESH_ProxyMesh.hxx"
38 #include "SMESH_subMesh.hxx"
39 #include "StdMeshers_PolygonPerFace_2D.hxx"
40 #include "StdMeshers_Regular_1D.hxx"
41 #include "StdMeshers_ViscousLayers.hxx"
42
43 #include <TopExp_Explorer.hxx>
44
45 #include <vector>
46 #include <set>
47
48 namespace
49 {
50   struct _EdgeMesher : public StdMeshers_Regular_1D
51   {
52     _EdgeMesher( int hypId, SMESH_Gen* gen )
53       : StdMeshers_Regular_1D( hypId, gen )
54     {
55       _hypType = NB_SEGMENTS;
56       _ivalue[ NB_SEGMENTS_IND ] = 1;
57     }
58   };
59
60   //=======================================================================
61   //function : addHexa
62   //purpose  :
63   //=======================================================================
64
65   const SMDS_MeshElement* addHexa( std::vector< const SMDS_MeshElement* >& faces,
66                                    const std::vector< int > &              quantities,
67                                    SMESH_MesherHelper &                    helper )
68   {
69     const SMDS_MeshElement* newHexa = 0;
70
71     // check nb of nodes in faces
72     for ( size_t i = 0; i < quantities.size(); ++i )
73       if ( quantities[ i ] != 4 )
74         return newHexa;
75
76     // look for a top face
77     const SMDS_MeshElement* topFace = 0;
78     const SMDS_MeshElement* botFace = faces[0];
79     std::vector< const SMDS_MeshNode* > nodes( 16 ); // last 8 is a working buffer
80     nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
81     for ( size_t iF = 1; iF < faces.size() &&  !topFace; ++iF )
82     {
83       bool hasCommonNode = false;
84       for ( int iN = 0; iN < quantities[ 0 ] &&  !hasCommonNode; ++iN )
85         hasCommonNode = ( faces[ iF ]->GetNodeIndex( nodes[ iN ]) >= 0 );
86
87       if ( !hasCommonNode )
88         topFace = faces[ iF ];
89     }
90
91     nodes.resize( 8 ); // set top nodes after hexa nodes - [8-11]
92     nodes.insert( nodes.end(), topFace->begin_nodes(), topFace->end_nodes() );
93     nodes.resize( 12 );
94     nodes.insert( nodes.end(), nodes.begin() + 8, nodes.begin() + 12 );
95
96     // find corresponding nodes of top and bottom by finding a side face including 2 node of each
97     SMESHDS_Mesh* mesh = helper.GetMeshDS();
98     const SMDS_MeshElement* sideFace = 0;
99     size_t i;
100     for ( i = 8; i < nodes.size()-1 &&  !sideFace; ++i )
101     {
102       sideFace = mesh->FindFace( nodes[0], nodes[1], nodes[ i ], nodes[ i + 1 ]);
103     }
104     if ( !sideFace )
105       return newHexa;
106
107     --i; // restore after ++i in the loop
108     bool botOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 0 ], nodes[ 1 ] );
109     bool topOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ i ], nodes[ i + 1 ] );
110     if ( botOriRight == topOriRight )
111     {
112       nodes[ 4 ] = nodes[ i + 1 ];
113       nodes[ 5 ] = nodes[ i + 0 ];
114       nodes[ 6 ] = nodes[ i + 3 ];
115       nodes[ 7 ] = nodes[ i + 2 ];
116     }
117     else
118     {
119       nodes[ 4 ] = nodes[ i + 0 ];
120       nodes[ 5 ] = nodes[ i + 1 ];
121       nodes[ 6 ] = nodes[ i + 2 ];
122       nodes[ 7 ] = nodes[ i + 3 ];
123     }
124
125     newHexa = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ], nodes[ 3 ],
126                                 nodes[ 4 ], nodes[ 5 ], nodes[ 6 ], nodes[ 7 ]);
127
128     return newHexa;
129   }
130
131   //=======================================================================
132   //function : addTetra
133   //purpose  :
134   //=======================================================================
135
136   const SMDS_MeshElement* addTetra( std::vector< const SMDS_MeshElement* >& faces,
137                                     const std::vector< int > &              quantities,
138                                     SMESH_MesherHelper &                    helper )
139   {
140     const SMDS_MeshElement* newTetra = 0;
141
142     // check nb of nodes in faces
143     for ( size_t i = 0; i < quantities.size(); ++i )
144       if ( quantities[ i ] != 3 )
145         return newTetra;
146
147     const SMDS_MeshElement* botFace = faces[0];
148
149     std::vector< const SMDS_MeshNode* > nodes( 6 );
150     nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
151     nodes.resize( 3 );
152
153     const SMDS_MeshNode* topNode = 0;
154     for ( size_t i = 0; i < 3 &&  !topNode; ++i )
155     {
156       topNode = faces[ 1 ]->GetNode( i );
157       if ( botFace->GetNodeIndex( topNode ) >= 0 )
158         topNode = 0;
159     }
160     if ( !topNode )
161       return newTetra;
162
163     nodes.push_back( topNode );
164
165     newTetra = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ], nodes[ 3 ]);
166
167     return newTetra;
168   }
169
170   //=======================================================================
171   //function : addPenta
172   //purpose  :
173   //=======================================================================
174
175   const SMDS_MeshElement* addPenta( std::vector< const SMDS_MeshElement* >& faces,
176                                     const std::vector< int > &              quantities,
177                                     SMESH_MesherHelper &                    helper )
178   {
179     const SMDS_MeshElement* newPenta = 0;
180
181     // check nb of nodes in faces and find triangle faces
182     int trias[2] = { -1, -1 };
183     for ( size_t i = 0; i < quantities.size(); ++i )
184       if ( quantities[ i ] != 4 )
185       {
186         if ( quantities[ i ] != 3 )
187           return newPenta;
188         int iTria = ( trias[0] != -1 );
189         if ( trias[ iTria ] != -1 )
190           return newPenta;
191         trias[ iTria ] = i;
192       }
193     if ( trias[1] == -1 )
194       return newPenta;
195
196     int iSide = trias[0] + 1;
197     if ( iSide == trias[1] )
198       ++iSide;
199     if (iSide == 5)
200       // use first side (otherwise, out of bounds)
201       iSide = 0;
202
203     const SMDS_MeshElement* botFace  = faces[ trias[0]];
204     const SMDS_MeshElement* topFace  = faces[ trias[1]];
205     const SMDS_MeshElement* sideFace = faces[ iSide ];
206     const SMDS_MeshNode* nodes[ 6 ] = { 0,0,0,0,0,0 };
207     for ( int i = 0 ; i < 3; ++i )
208     {
209       const SMDS_MeshNode* botNode = botFace->GetNode( i );
210       if ( sideFace->GetNodeIndex( botNode ) < 0 )
211         nodes[2] = botNode;
212       else
213         nodes[ bool( nodes[0] )] = botNode;
214
215       const SMDS_MeshNode* topNode = topFace->GetNode( i );
216       if ( sideFace->GetNodeIndex( topNode ) < 0 )
217         nodes[5] = topNode;
218       else
219         nodes[ 3 + bool( nodes[3]) ] = topNode;
220     }
221     bool botOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 0 ], nodes[ 1 ]);
222     bool topOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 3 ], nodes[ 4 ]);
223     if ( botOriRight == topOriRight )
224       std::swap( nodes[ 3 ], nodes[ 4 ]);
225
226     newPenta = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ],
227                                  nodes[ 3 ], nodes[ 4 ], nodes[ 5 ]);
228
229     return newPenta;
230   }
231
232   //=======================================================================
233   //function : addPyra
234   //purpose  :
235   //=======================================================================
236
237   const SMDS_MeshElement* addPyra( std::vector< const SMDS_MeshElement* >& faces,
238                                    const std::vector< int > &              quantities,
239                                    SMESH_MesherHelper &                    helper )
240   {
241     const SMDS_MeshElement* newPyra = 0;
242
243     // check nb of nodes in faces
244     int iBot = -1;
245     for ( size_t i = 0; i < quantities.size(); ++i )
246       if ( quantities[ i ] != 3 )
247       {
248         if ( quantities[ i ] != 4 || iBot != -1 )
249           return newPyra;
250         iBot = i;
251       }
252
253     const SMDS_MeshElement* botFace = faces[ iBot ];
254
255     std::vector< const SMDS_MeshNode* > nodes( 8 );
256     nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
257     nodes.resize( 4 );
258
259     const SMDS_MeshNode* topNode = 0;
260     for ( size_t i = 0; i < 4 &&  !topNode; ++i )
261     {
262       topNode = faces[ 1 ]->GetNode( i );
263       if ( botFace->GetNodeIndex( topNode ) >= 0 )
264         topNode = 0;
265     }
266     if ( !topNode )
267       return newPyra;
268
269     nodes.push_back( topNode );
270
271     newPyra = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ], nodes[ 3 ], nodes[4] );
272
273     return newPyra;
274   }
275
276   //=======================================================================
277   //function : addHPrism
278   //purpose  : add hexagonal prism
279   //=======================================================================
280
281   const SMDS_MeshElement* addHPrism( std::vector< const SMDS_MeshElement* >& faces,
282                                      const std::vector< int > &              quantities,
283                                      SMESH_MesherHelper &                    helper )
284   {
285     const SMDS_MeshElement* newHexPrism = 0;
286
287     // check nb of nodes in faces and find hexagons
288     int hexa[2] = { -1, -1 };
289     for ( size_t i = 0; i < quantities.size(); ++i )
290       if ( quantities[ i ] != 4 )
291       {
292         if ( quantities[ i ] != 6 )
293           return newHexPrism;
294         int iHex = ( hexa[0] != -1 );
295         if ( hexa[ iHex ] != -1 )
296           return newHexPrism;
297         hexa[ iHex ] = i;
298       }
299     if ( hexa[1] == -1 )
300       return newHexPrism;
301
302     int iSide = hexa[0] + 1;
303     if ( iSide == hexa[1] )
304       ++iSide;
305
306     const SMDS_MeshElement* botFace = faces[ hexa[ 0 ]];
307     const SMDS_MeshElement* topFace = faces[ hexa[ 1 ]];
308     std::vector< const SMDS_MeshNode* > nodes( 24 ); // last 12 is a working buffer
309
310     nodes.assign( botFace->begin_nodes(), botFace->end_nodes() );
311     nodes.resize( 12 ); // set top nodes after hexa nodes - [12-17]
312     nodes.insert( nodes.end(), topFace->begin_nodes(), topFace->end_nodes() );
313     nodes.resize( 18 );
314     nodes.insert( nodes.end(), nodes.begin() + 12, nodes.begin() + 18 );
315
316     // find corresponding nodes of top and bottom by finding a side face including 2 node of each
317     SMESHDS_Mesh* mesh = helper.GetMeshDS();
318     const SMDS_MeshElement* sideFace = 0;
319     size_t i;
320     for ( i = 12; i < nodes.size()-1 &&  !sideFace; ++i )
321     {
322       sideFace = mesh->FindFace( nodes[0], nodes[1], nodes[ i ], nodes[ i + 1 ]);
323     }
324     if ( !sideFace )
325       return newHexPrism;
326
327     --i; // restore after ++i in the loop
328     bool botOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ 0 ], nodes[ 1 ] );
329     bool topOriRight = SMESH_MeshAlgos::IsRightOrder( sideFace, nodes[ i ], nodes[ i + 1 ] );
330     if ( botOriRight == topOriRight )
331     {
332       nodes[ 6  ] = nodes[ i + 1 ];
333       nodes[ 7  ] = nodes[ i + 0 ];
334       nodes[ 8  ] = nodes[ i + 5 ];
335       nodes[ 9  ] = nodes[ i + 4 ];
336       nodes[ 10 ] = nodes[ i + 3 ];
337       nodes[ 11 ] = nodes[ i + 2 ];
338     }
339     else
340     {
341       nodes[ 6  ] = nodes[ i + 0 ];
342       nodes[ 7  ] = nodes[ i + 1 ];
343       nodes[ 8  ] = nodes[ i + 2 ];
344       nodes[ 9  ] = nodes[ i + 3 ];
345       nodes[ 10 ] = nodes[ i + 4 ];
346       nodes[ 11 ] = nodes[ i + 5 ];
347     }
348
349     newHexPrism = helper.AddVolume( nodes[ 0 ], nodes[ 1 ], nodes[ 2 ],
350                                     nodes[ 3 ], nodes[ 4 ], nodes[ 5 ],
351                                     nodes[ 6 ], nodes[ 7 ], nodes[ 8 ],
352                                     nodes[ 9 ], nodes[10 ], nodes[11 ]);
353
354     return newHexPrism;
355   }
356
357   //=======================================================================
358   //function : addPoly
359   //purpose  :
360   //=======================================================================
361
362   const SMDS_MeshElement* addPoly( std::vector< const SMDS_MeshElement* >& faces,
363                                    const std::vector< int > &              quantities,
364                                    SMESH_MesherHelper &                    helper )
365   {
366     const SMDS_MeshElement* newPoly = 0;
367
368     std::vector< const SMDS_MeshNode* > nodes;
369     for ( size_t iF = 0; iF < faces.size(); ++iF )
370       nodes.insert( nodes.end(), faces[iF]->begin_nodes(), faces[iF]->end_nodes() );
371
372     newPoly = helper.AddPolyhedralVolume( nodes, quantities );
373
374     return newPoly;
375   }
376
377 } // namespace
378
379 //=======================================================================
380 //function : StdMeshers_PolyhedronPerSolid_3D
381 //purpose  :
382 //=======================================================================
383
384 StdMeshers_PolyhedronPerSolid_3D::StdMeshers_PolyhedronPerSolid_3D(int        hypId,
385                                                                    SMESH_Gen* gen)
386   :SMESH_3D_Algo(hypId, gen),
387    myEdgeMesher( new _EdgeMesher( gen->GetANewId(), gen )),
388    myFaceMesher( new StdMeshers_PolygonPerFace_2D( gen->GetANewId(), gen ))
389 {
390   _name = "PolyhedronPerSolid_3D";
391   _requireDiscreteBoundary = false;
392   _supportSubmeshes = true;
393   _compatibleHypothesis.push_back("ViscousLayers");
394   _neededLowerHyps[0] = _neededLowerHyps[1] = _neededLowerHyps[2] = true;
395 }
396
397 //=======================================================================
398 //function : ~StdMeshers_PolyhedronPerSolid_3D
399 //purpose  :
400 //=======================================================================
401
402 StdMeshers_PolyhedronPerSolid_3D::~StdMeshers_PolyhedronPerSolid_3D()
403 {
404   delete myEdgeMesher;
405   delete myFaceMesher;
406 }
407
408 //=======================================================================
409 //function : CheckHypothesis
410 //purpose  :
411 //=======================================================================
412
413 bool StdMeshers_PolyhedronPerSolid_3D::CheckHypothesis(SMESH_Mesh&         theMesh,
414                                                        const TopoDS_Shape& theShape,
415                                                        Hypothesis_Status&  theStatus)
416 {
417   myViscousLayersHyp = NULL;
418
419   const std::list<const SMESHDS_Hypothesis*>& hyps =
420     GetUsedHypothesis( theMesh, theShape, /*ignoreAuxiliary=*/false);
421   std::list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
422   if ( h == hyps.end())
423   {
424     theStatus = SMESH_Hypothesis::HYP_OK;
425     return true;
426   }
427
428   // only StdMeshers_ViscousLayers can be used
429   theStatus = HYP_OK;
430   for ( ; h != hyps.end(); ++h )
431   {
432     if ( !(myViscousLayersHyp = dynamic_cast< const StdMeshers_ViscousLayers*> ( *h )))
433       break;
434   }
435   if ( !myViscousLayersHyp )
436     theStatus = HYP_INCOMPATIBLE;
437   else
438     error( myViscousLayersHyp->CheckHypothesis( theMesh, theShape, theStatus ));
439
440   return theStatus == HYP_OK;
441 }
442
443 //=======================================================================
444 //function : Compute
445 //purpose  :
446 //=======================================================================
447
448 bool StdMeshers_PolyhedronPerSolid_3D::Compute(SMESH_Mesh&         theMesh,
449                                                const TopoDS_Shape& theShape)
450 {
451   const SMDS_MeshElement* newVolume = 0;
452
453   SMESH_subMesh* sm = theMesh.GetSubMesh( theShape );
454   SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator( /*includeSelf=*/true,
455                                                             /*complexFirst=*/false);
456   while ( smIt->more() )
457   {
458     sm = smIt->next();
459     if ( !sm->IsEmpty() )
460       continue;
461
462     const TopoDS_Shape & shape = sm->GetSubShape();
463     switch ( shape.ShapeType() )
464     {
465     case TopAbs_VERTEX:
466       sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
467       break;
468
469     case TopAbs_EDGE:
470       sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
471       if ( sm->IsEmpty() )
472         myEdgeMesher->Compute( theMesh, shape );
473       break;
474
475     case TopAbs_FACE:
476       sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
477       if ( sm->IsEmpty() && !myFaceMesher->Compute( theMesh, shape ))
478       {
479         sm->GetComputeError() = myFaceMesher->GetComputeError();
480         sm->GetComputeError()->myAlgo = myFaceMesher;
481         return false;
482       }
483       break;
484
485     case TopAbs_SOLID:
486     {
487       SMESH_MesherHelper helper( theMesh );
488       helper.SetElementsOnShape( true );
489       _quadraticMesh = helper.IsQuadraticSubMesh( shape );
490
491       SMESH_ProxyMesh::Ptr proxymesh( new SMESH_ProxyMesh( theMesh ));
492       if ( myViscousLayersHyp )
493       {
494         proxymesh = myViscousLayersHyp->Compute( theMesh, theShape );
495         if ( !proxymesh )
496           return false;
497       }
498
499       std::vector< const SMDS_MeshElement* > faces;
500       faces.reserve( 20 );
501
502       for ( TopExp_Explorer faceEx( shape, TopAbs_FACE ); faceEx.More(); faceEx.Next() )
503       {
504         const SMESHDS_SubMesh* smDS = proxymesh->GetSubMesh( faceEx.Current() );
505         for ( SMDS_ElemIteratorPtr faceIt = smDS->GetElements(); faceIt->more(); )
506           faces.push_back( faceIt->next() );
507       }
508
509       bool useMediumNodes = false;
510       if ( !_quadraticMesh && theMesh.GetMeshDS()->GetMeshInfo().NbFaces( ORDER_QUADRATIC ))
511         for ( size_t i = 0; i < faces.size() &&  !useMediumNodes ; ++i )
512           useMediumNodes = faces[ i ]->IsQuadratic();
513
514       std::vector< int > quantities( faces.size() );
515       std::set< const SMDS_MeshNode* > nodes;
516       for ( size_t i = 0; i < faces.size(); ++i )
517       {
518         quantities[ i ] = useMediumNodes ? faces[ i ]->NbNodes() : faces[ i ]->NbCornerNodes();
519         for ( int iN = 0; iN < quantities[ i ]; ++iN )
520           nodes.insert( faces[ i ]->GetNode( iN ));
521       }
522
523       const size_t nbNodes = nodes.size(), nbFaces = faces.size();
524       if (      nbNodes == 8  && nbFaces == 6 ) newVolume = addHexa  ( faces, quantities, helper );
525       else if ( nbNodes == 4  && nbFaces == 4 ) newVolume = addTetra ( faces, quantities, helper );
526       else if ( nbNodes == 6  && nbFaces == 5 ) newVolume = addPenta ( faces, quantities, helper );
527       else if ( nbNodes == 5  && nbFaces == 5 ) newVolume = addPyra  ( faces, quantities, helper );
528       else if ( nbNodes == 12 && nbFaces == 8 ) newVolume = addHPrism( faces, quantities, helper );
529       if ( !newVolume )
530         newVolume = addPoly ( faces, quantities, helper );
531
532       if ( newVolume )
533       {
534         SMESH::Controls::BadOrientedVolume checker;
535         checker.SetMesh( theMesh.GetMeshDS() );
536         if ( checker.IsSatisfy( newVolume->GetID() ))
537         {
538           SMESH_MeshEditor editor( &theMesh );
539           editor.Reorient( newVolume );
540         }
541       }
542     }
543     default:;
544
545     } // switch ( shape.ShapeType() )
546   } // loop on sub-meshes
547
548   return newVolume;
549 }
550
551 //=======================================================================
552 //function : Evaluate
553 //purpose  :
554 //=======================================================================
555
556 bool StdMeshers_PolyhedronPerSolid_3D::Evaluate(SMESH_Mesh&         theMesh,
557                                                 const TopoDS_Shape& theShape,
558                                                 MapShapeNbElems&    theResMap)
559 {
560   _quadraticMesh = false;
561
562   SMESH_subMesh* sm = theMesh.GetSubMesh( theShape );
563   SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator( /*includeSelf=*/true,
564                                                             /*complexFirst=*/false);
565   while ( smIt->more() )
566   {
567     sm = smIt->next();
568
569     MapShapeNbElems::iterator sm2vec = theResMap.find( sm );
570     if ( sm2vec != theResMap.end() && !sm2vec->second.empty() )
571       continue;
572
573     const TopoDS_Shape & shape = sm->GetSubShape();
574     switch ( shape.ShapeType() )
575     {
576     case TopAbs_EDGE:
577       myEdgeMesher->Evaluate( theMesh, shape, theResMap );
578       break;
579
580     case TopAbs_FACE:
581     {
582       myFaceMesher->Evaluate( theMesh, shape, theResMap );
583       std::vector<int> & quantities = theResMap[ sm ];
584       _quadraticMesh = ( !quantities.empty() &&
585                          ( quantities[ SMDSEntity_Quad_Triangle   ] +
586                            quantities[ SMDSEntity_Quad_Quadrangle ] +
587                            quantities[ SMDSEntity_Quad_Polygon    ]));
588       break;
589     }
590
591     case TopAbs_SOLID:
592     {
593       std::vector<int> & quantities = theResMap[ sm ];
594       quantities.resize( SMDSEntity_Last, 0 );
595
596       SMESH_MesherHelper helper( theMesh );
597       const int nbNodes = helper.Count( shape, TopAbs_VERTEX, /*ignoreSame=*/true );
598       const int nbFaces = helper.Count( shape, TopAbs_FACE,   /*ignoreSame=*/false );
599
600       if (      nbNodes == 8 && nbFaces == 6 )
601         quantities[ _quadraticMesh ? SMDSEntity_Quad_Hexa : SMDSEntity_Hexa ] = 1;
602       else if ( nbNodes == 4 && nbFaces == 4 )
603         quantities[ _quadraticMesh ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra ] = 1;
604       else if ( nbNodes == 6 && nbFaces == 5 )
605         quantities[ _quadraticMesh ? SMDSEntity_Quad_Penta : SMDSEntity_Penta ] = 1;
606       else if ( nbNodes == 5 && nbFaces == 5 )
607         quantities[ _quadraticMesh ? SMDSEntity_Quad_Pyramid : SMDSEntity_Pyramid ] = 1;
608       else if ( nbNodes == 12 && nbFaces == 8 )
609         quantities[ /*_quadraticMesh ? SMDSEntity_Quad_Pyramid :*/ SMDSEntity_Hexagonal_Prism ] = 1;
610       else
611         quantities[ /*_quadraticMesh ? SMDSEntity_Quad_Polyhedra : */SMDSEntity_Polyhedra ] = 1;
612
613       return true;
614     }
615     default:;
616
617     } // switch ( shape.ShapeType() )
618   } // loop on sub-meshes
619
620   return false;
621 }