Salome HOME
Merge from V6_3_BR 06/06/2011
[modules/smesh.git] / src / StdMeshers / StdMeshers_HexaFromSkin_3D.cxx
1 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // File      : StdMeshers_HexaFromSkin_3D.cxx
21 // Created   : Wed Jan 27 12:28:07 2010
22 // Author    : Edward AGAPOV (eap)
23 //
24 #include "StdMeshers_HexaFromSkin_3D.hxx"
25
26 #include "SMDS_VolumeOfNodes.hxx"
27 #include "SMDS_VolumeTool.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_MesherHelper.hxx"
30
31 #include <gp_Ax2.hxx>
32
33 //#include "utilities.h"
34 #include <limits>
35
36 // Define error message
37 #ifdef _DEBUG_
38 #define BAD_MESH_ERR \
39   error(SMESH_Comment("Can't detect block-wise structure of the input 2D mesh.\n" \
40                       __FILE__ ":" )<<__LINE__)
41 #else
42 #define BAD_MESH_ERR \
43   error(SMESH_Comment("Can't detect block-wise structure of the input 2D mesh"))
44 #endif
45
46 // Debug output
47 #define _DUMP_(msg) //cout << msg << endl
48
49
50
51 namespace
52 {
53   enum EBoxSides //!< sides of the block
54     {
55       B_BOTTOM=0, B_RIGHT, B_TOP, B_LEFT, B_FRONT, B_BACK, NB_BLOCK_SIDES
56     };
57 //   const char* SBoxSides[] = //!< names of block sides -- needed for DEBUG only
58 //     {
59 //       "BOTTOM", "RIGHT", "TOP", "LEFT", "FRONT", "BACK", "UNDEFINED"
60 //     };
61   enum EQuadEdge //!< edges of quadrangle side
62     {
63       Q_BOTTOM = 0, Q_RIGHT, Q_TOP, Q_LEFT, NB_QUAD_SIDES
64     };
65
66
67   //================================================================================
68   /*!
69    * \brief return logical coordinates (i.e. min or max) of ends of edge
70    */
71   //================================================================================
72
73   bool getEdgeEnds(EQuadEdge edge, bool& xMax1, bool& yMax1, bool& xMax2, bool& yMax2 )
74   {
75     xMax1=0, yMax1=0, xMax2=1, yMax2=1;
76     switch( edge )
77     {
78     case Q_BOTTOM: yMax2 = 0; break;
79     case Q_RIGHT:  xMax1 = 1; break;
80     case Q_TOP:    yMax1 = 1; break;
81     case Q_LEFT:   xMax2 = 0; break;
82     default:
83       return false;
84     }
85     return true;
86   }
87
88   //================================================================================
89   /*!
90    * \brief Description of node used to detect corner nodes
91    */
92   struct _NodeDescriptor
93   {
94     int nbInverseFaces, nbNodesInInverseFaces;
95
96     _NodeDescriptor(const SMDS_MeshNode* n): nbInverseFaces(0), nbNodesInInverseFaces(0)
97     {
98       if ( n )
99       {
100         set<const SMDS_MeshNode*> nodes;
101         SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator(SMDSAbs_Face );
102         while ( fIt->more() )
103         {
104           const SMDS_MeshElement* face = fIt->next();
105           nodes.insert( face->begin_nodes(), face->end_nodes() );
106           nbInverseFaces++;
107         }
108         nbNodesInInverseFaces = nodes.size();
109       }
110     }
111     bool operator==(const _NodeDescriptor& other) const
112     {
113       return
114         nbInverseFaces == other.nbInverseFaces &&
115         nbNodesInInverseFaces == other.nbNodesInInverseFaces;
116     }
117   };
118
119   //================================================================================
120   /*!
121    * \brief return true if a node is at block corner
122    *
123    * This check is valid for simple cases only
124    */
125   //================================================================================
126
127   bool isCornerNode( const SMDS_MeshNode* n )
128   {
129     return n && n->NbInverseElements( SMDSAbs_Face ) % 2;
130   }
131
132   //================================================================================
133   /*!
134    * \brief check element type
135    */
136   //================================================================================
137
138   bool isQuadrangle(const SMDS_MeshElement* e)
139   {
140     return ( e && e->NbCornerNodes() == 4 );
141   }
142
143   //================================================================================
144   /*!
145    * \brief return opposite node of a quadrangle face
146    */
147   //================================================================================
148
149   const SMDS_MeshNode* oppositeNode(const SMDS_MeshElement* quad, int iNode)
150   {
151     return quad->GetNode( (iNode+2) % 4 );
152   }
153
154   //================================================================================
155   /*!
156    * \brief Convertor of a pair of integers to a sole index
157    */
158   struct _Indexer
159   {
160     int _xSize, _ySize;
161     _Indexer( int xSize=0, int ySize=0 ): _xSize(xSize), _ySize(ySize) {}
162     int size() const { return _xSize * _ySize; }
163     int operator()(int x, int y) const { return y * _xSize + x; }
164   };
165   //================================================================================
166   /*!
167    * \brief Oriented convertor of a pair of integers to a sole index 
168    */
169   class _OrientedIndexer : public _Indexer
170   {
171   public:
172     enum OriFlags //!< types of block side orientation
173       {
174         REV_X = 1, REV_Y = 2, SWAP_XY = 4, MAX_ORI = REV_X|REV_Y|SWAP_XY
175       };
176     _OrientedIndexer( const _Indexer& indexer, const int oriFlags ):
177       _Indexer( indexer._xSize, indexer._ySize ),
178       _xSize (indexer._xSize), _ySize(indexer._ySize),
179       _xRevFun((oriFlags & REV_X) ? & reverse : & lazy),
180       _yRevFun((oriFlags & REV_Y) ? & reverse : & lazy),
181       _swapFun((oriFlags & SWAP_XY ) ? & swap : & lazy)
182     {
183       (*_swapFun)( _xSize, _ySize );
184     }
185     //!< Return index by XY
186     int operator()(int x, int y) const
187     {
188       (*_xRevFun)( x, const_cast<int&>( _xSize ));
189       (*_yRevFun)( y, const_cast<int&>( _ySize ));
190       (*_swapFun)( x, y );
191       return _Indexer::operator()(x,y);
192     }
193     //!< Return index for a corner
194     int corner(bool xMax, bool yMax) const
195     {
196       int x = xMax, y = yMax, size = 2;
197       (*_xRevFun)( x, size );
198       (*_yRevFun)( y, size );
199       (*_swapFun)( x, y );
200       return _Indexer::operator()(x ? _Indexer::_xSize-1 : 0 , y ? _Indexer::_ySize-1 : 0);
201     }
202     int xSize() const { return _xSize; }
203     int ySize() const { return _ySize; }
204   private:
205     _Indexer _indexer;
206     int _xSize, _ySize;
207
208     typedef void (*TFun)(int& x, int& y);
209     TFun _xRevFun, _yRevFun, _swapFun;
210     
211     static void lazy(int&, int&) {}
212     static void reverse(int& x, int& size) { x = size - x - 1; }
213     static void swap(int& x, int& y) { std::swap(x,y); }
214   };
215   //================================================================================
216   /*!
217    * \brief Structure corresponding to the meshed side of block
218    */
219   struct _BlockSide
220   {
221     vector<const SMDS_MeshNode*> _grid;
222     _Indexer                     _index;
223     int                          _nbBlocksExpected;
224     int                          _nbBlocksFound;
225
226 #ifdef _DEBUG_ // want to get SIGSEGV in case of invalid index
227 #define _grid_access_(pobj, i) pobj->_grid[ ((i) < pobj->_grid.size()) ? i : int(1e100)]
228 #else
229 #define _grid_access_(pobj, i) pobj->_grid[ i ]
230 #endif
231     //!< Return node at XY
232     const SMDS_MeshNode* getNode(int x, int y) const { return _grid_access_(this, _index( x,y ));}
233     //!< Set node at XY
234     void setNode(int x, int y, const SMDS_MeshNode* n) { _grid_access_(this, _index( x,y )) = n; }
235     //!< Return an edge
236     SMESH_OrientedLink getEdge(EQuadEdge edge) const
237     {
238       bool x1, y1, x2, y2; getEdgeEnds( edge, x1, y1, x2, y2 );
239       return SMESH_OrientedLink( getCornerNode ( x1, y1 ), getCornerNode( x2, y2 ));
240     }
241     //!< Return a corner node
242     const SMDS_MeshNode* getCornerNode(bool isXMax, bool isYMax) const
243     {
244       return getNode( isXMax ? _index._xSize-1 : 0 , isYMax ? _index._ySize-1 : 0 );
245     }
246     const SMDS_MeshElement* getCornerFace(const SMDS_MeshNode* cornerNode) const;
247     //!< True if all blocks this side belongs to have beem found
248     bool isBound() const { return _nbBlocksExpected <= _nbBlocksFound; }
249     //!< Return coordinates of node at XY
250     gp_XYZ getXYZ(int x, int y) const { return SMESH_TNodeXYZ( getNode( x, y )); }
251     //!< Return gravity center of the four corners and the middle node
252     gp_XYZ getGC() const
253     {
254       gp_XYZ xyz =
255         getXYZ( 0, 0 ) +
256         getXYZ( _index._xSize-1, 0 ) +
257         getXYZ( 0, _index._ySize-1 ) +
258         getXYZ( _index._xSize-1, _index._ySize-1 ) +
259         getXYZ( _index._xSize/2, _index._ySize/2 );
260       return xyz / 5;
261     }
262     //!< Return number of mesh faces
263     int getNbFaces() const { return (_index._xSize-1) * (_index._ySize-1); }
264   };
265   //================================================================================
266   /*!
267    * \brief _BlockSide with changed orientation
268    */
269   struct _OrientedBlockSide
270   {
271     _BlockSide*       _side;
272     _OrientedIndexer  _index;
273
274     _OrientedBlockSide( _BlockSide* side=0, const int oriFlags=0 ):
275       _side(side), _index(side ? side->_index : _Indexer(), oriFlags ) {}
276     //!< return coordinates by XY
277     gp_XYZ xyz(int x, int y) const
278     {
279       return SMESH_TNodeXYZ( _grid_access_(_side, _index( x, y )) );
280     }
281     //!< safely return a node by XY
282     const SMDS_MeshNode* node(int x, int y) const
283     {
284       int i = _index( x, y );
285       return ( i < 0 || i >= _side->_grid.size()) ? 0 : _side->_grid[i];
286     }
287     //!< Return an edge
288     SMESH_OrientedLink edge(EQuadEdge edge) const
289     {
290       bool x1, y1, x2, y2; getEdgeEnds( edge, x1, y1, x2, y2 );
291       return SMESH_OrientedLink( cornerNode ( x1, y1 ), cornerNode( x2, y2 ));
292     }
293     //!< Return a corner node
294     const SMDS_MeshNode* cornerNode(bool isXMax, bool isYMax) const
295     {
296       return _grid_access_(_side, _index.corner( isXMax, isYMax ));
297     }
298     //!< return its size in nodes
299     int getHoriSize() const { return _index.xSize(); }
300     int getVertSize() const  { return _index.ySize(); }
301     //!< True if _side has been initialized
302     operator bool() const { return _side; }
303     //! Direct access to _side
304     const _BlockSide* operator->() const { return _side; }
305     _BlockSide* operator->() { return _side; }
306   };
307   //================================================================================
308   /*!
309    * \brief Meshed skin of block
310    */
311   struct _Block
312   {
313     _OrientedBlockSide        _side[6]; // 6 sides of a sub-block
314     set<const SMDS_MeshNode*> _corners;
315
316     const _OrientedBlockSide& getSide(int i) const { return _side[i]; }
317     bool setSide( int i, const _OrientedBlockSide& s)
318     {
319       if (( _side[i] = s ))
320       {
321         _corners.insert( s.cornerNode(0,0));
322         _corners.insert( s.cornerNode(1,0));
323         _corners.insert( s.cornerNode(0,1));
324         _corners.insert( s.cornerNode(1,1));
325       }
326       return s;
327     }
328     void clear() { for (int i=0;i<6;++i) _side[i]=0; _corners.clear(); }
329     bool hasSide( const _OrientedBlockSide& s) const
330     {
331       if ( s ) for (int i=0;i<6;++i) if ( _side[i] && _side[i]._side == s._side ) return true;
332       return false;
333     }
334     int nbSides() const { int n=0; for (int i=0;i<6;++i) if ( _side[i] ) ++n; return n; }
335     bool isValid() const;
336   };
337   //================================================================================
338   /*!
339    * \brief Skin mesh possibly containing several meshed blocks
340    */
341   class _Skin
342   {
343   public:
344
345     int findBlocks(SMESH_Mesh& mesh);
346     //!< return i-th block
347     const _Block& getBlock(int i) const { return _blocks[i]; }
348     //!< return error description
349     const SMESH_Comment& error() const { return _error; }
350
351   private:
352     bool fillSide( _BlockSide&             side,
353                    const SMDS_MeshElement* cornerQuad,
354                    const SMDS_MeshNode*    cornerNode);
355     bool fillRowsUntilCorner(const SMDS_MeshElement* quad,
356                              const SMDS_MeshNode*    n1,
357                              const SMDS_MeshNode*    n2,
358                              vector<const SMDS_MeshNode*>& verRow1,
359                              vector<const SMDS_MeshNode*>& verRow2,
360                              bool alongN1N2 );
361     _OrientedBlockSide findBlockSide( EBoxSides startBlockSide,
362                                       EQuadEdge sharedSideEdge1,
363                                       EQuadEdge sharedSideEdge2,
364                                       bool      withGeometricAnalysis);
365     //!< update own data and data of the side bound to block
366     void setSideBoundToBlock( _BlockSide& side )
367     {
368       if ( side._nbBlocksFound++, side.isBound() )
369         for ( int e = 0; e < int(NB_QUAD_SIDES); ++e )
370           _edge2sides[ side.getEdge( (EQuadEdge) e ) ].erase( &side );
371     }
372     //!< store reason of error
373     int error(const SMESH_Comment& reason) { _error = reason; return 0; }
374
375     SMESH_Comment      _error;
376
377     list< _BlockSide > _allSides;
378     vector< _Block >   _blocks;
379
380     //map< const SMDS_MeshNode*, set< _BlockSide* > > _corner2sides;
381     map< SMESH_OrientedLink, set< _BlockSide* > > _edge2sides;
382   };
383
384   //================================================================================
385   /*!
386    * \brief Find and return number of submeshes corresponding to blocks
387    */
388   //================================================================================
389
390   int _Skin::findBlocks(SMESH_Mesh& mesh)
391   {
392     SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
393
394     // Find a node at any block corner
395
396     SMDS_NodeIteratorPtr nIt = meshDS->nodesIterator(/*idInceasingOrder=*/true);
397     if ( !nIt->more() ) return error("Empty mesh");
398
399     const SMDS_MeshNode* nCorner = 0;
400     while ( nIt->more() )
401     {
402       nCorner = nIt->next();
403       if ( isCornerNode( nCorner ))
404         break;
405       else
406         nCorner = 0;
407     }
408     if ( !nCorner )
409       return BAD_MESH_ERR;
410
411     // --------------------------------------------------------------------
412     // Find all block sides starting from mesh faces sharing the corner node
413     // --------------------------------------------------------------------
414
415     int nbFacesOnSides = 0;
416     TIDSortedElemSet cornerFaces; // corner faces of found _BlockSide's
417     list< const SMDS_MeshNode* > corners( 1, nCorner );
418     list< const SMDS_MeshNode* >::iterator corner = corners.begin();
419     while ( corner != corners.end() )
420     {
421       SMDS_ElemIteratorPtr faceIt = (*corner)->GetInverseElementIterator( SMDSAbs_Face );
422       while ( faceIt->more() )
423       {
424         const SMDS_MeshElement* face = faceIt->next();
425         if ( !cornerFaces.insert( face ).second )
426           continue; // already loaded block side
427
428         if ( !isQuadrangle( face ))
429           return error("Non-quadrangle elements in the input mesh");
430
431         if ( _allSides.empty() || !_allSides.back()._grid.empty() )
432           _allSides.push_back( _BlockSide() );
433
434         _BlockSide& side = _allSides.back();
435         if ( !fillSide( side, face, *corner ) )
436         {
437           if ( !_error.empty() )
438             return false;
439         }
440         else
441         {
442           for ( int isXMax = 0; isXMax < 2; ++isXMax )
443             for ( int isYMax = 0; isYMax < 2; ++isYMax )
444             {
445               const SMDS_MeshNode* nCorner = side.getCornerNode(isXMax,isYMax );
446               corners.push_back( nCorner );
447               cornerFaces.insert( side.getCornerFace( nCorner ));
448             }
449           for ( int e = 0; e < int(NB_QUAD_SIDES); ++e )
450             _edge2sides[ side.getEdge( (EQuadEdge) e ) ].insert( &side );
451
452           nbFacesOnSides += side.getNbFaces();
453         }
454       }
455       ++corner;
456
457       // find block sides of other domains if any
458       if ( corner == corners.end() && nbFacesOnSides < mesh.NbQuadrangles() )
459       {
460         while ( nIt->more() )
461         {
462           nCorner = nIt->next();
463           if ( isCornerNode( nCorner ))
464             corner = corners.insert( corner, nCorner );
465         }
466         nbFacesOnSides = mesh.NbQuadrangles();
467       }
468     }
469     
470     if ( _allSides.empty() )
471       return BAD_MESH_ERR;
472     if ( _allSides.back()._grid.empty() )
473       _allSides.pop_back();
474     _DUMP_("Nb detected sides "<< _allSides.size());
475
476     // ---------------------------
477     // Organize sides into blocks
478     // ---------------------------
479
480     // analyse sharing of sides by blocks and sort sides by nb of adjacent sides
481     int nbBlockSides = 0; // total nb of block sides taking into account their sharing
482     multimap<int, _BlockSide* > sortedSides;
483     {
484       list < _BlockSide >::iterator sideIt = _allSides.begin();
485       for ( ; sideIt != _allSides.end(); ++sideIt )
486       {
487         _BlockSide& side = *sideIt;
488         bool isSharedSide = true;
489         int nbAdjacent = 0;
490         for ( int e = 0; e < int(NB_QUAD_SIDES) && isSharedSide; ++e )
491         {
492           int nbAdj = _edge2sides[ side.getEdge( (EQuadEdge) e ) ].size();
493           nbAdjacent += nbAdj;
494           isSharedSide = ( nbAdj > 2 );
495         }
496         side._nbBlocksFound = 0;
497         side._nbBlocksExpected = isSharedSide ? 2 : 1;
498         nbBlockSides += side._nbBlocksExpected;
499         sortedSides.insert( make_pair( nbAdjacent, & side ));
500       }
501     }
502
503     // find sides of each block
504     int nbBlocks = 0;
505     while ( nbBlockSides >= 6 )
506     {
507       // get any side not bound to all blocks it belongs to
508       multimap<int, _BlockSide*>::iterator i_side = sortedSides.begin();
509       while ( i_side != sortedSides.end() && i_side->second->isBound())
510         ++i_side;
511
512       // start searching for block sides from the got side
513       bool ok = true;
514       if ( _blocks.empty() || _blocks.back()._side[B_FRONT] )
515         _blocks.resize( _blocks.size() + 1 );
516
517       _Block& block = _blocks.back();
518       block.setSide( B_FRONT, i_side->second );
519       setSideBoundToBlock( *i_side->second );
520       nbBlockSides--;
521
522       // edges of adjacent sides of B_FRONT corresponding to front's edges
523       EQuadEdge edgeOfFront[4] = { Q_BOTTOM, Q_RIGHT, Q_TOP, Q_LEFT };
524       EQuadEdge edgeOfAdj  [4] = { Q_BOTTOM, Q_LEFT, Q_BOTTOM, Q_LEFT };
525       // first find all sides detectable w/o advanced analysis,
526       // then repeat the search, which then may pass without advanced analysis
527       for ( int advAnalys = 0; advAnalys < 2; ++advAnalys )
528       {
529         for ( int i = 0; (ok || !advAnalys) && i < NB_QUAD_SIDES; ++i )
530           if ( !block._side[i] ) // try to find 4 sides adjacent to front side
531             ok = block.setSide( i, findBlockSide( B_FRONT, edgeOfFront[i],edgeOfAdj[i],advAnalys));
532         if ( ok || !advAnalys)
533           if ( !block._side[B_BACK] && block._side[B_TOP] ) // try to find back side by top one
534             ok = block.setSide( B_BACK, findBlockSide( B_TOP, Q_TOP, Q_TOP, advAnalys ));
535         if ( !advAnalys ) ok = true;
536       }
537       ok = block.isValid();
538       if ( ok )
539       {
540         // check if just found block is same as one of previously found
541         bool isSame = false;
542         for ( int i = 1; i < _blocks.size() && !isSame; ++i )
543           isSame = ( block._corners == _blocks[i-1]._corners );
544         ok = !isSame;
545       }
546
547       // count the found sides
548       _DUMP_(endl << "** Block " << _blocks.size() << " valid: " << block.isValid());
549       for (int i = 0; i < NB_BLOCK_SIDES; ++i )
550       {
551         _DUMP_("\tSide "<< SBoxSides[i] <<" "<< block._side[ i ]._side);
552         if ( block._side[ i ] )
553         {
554           if ( ok && i != B_FRONT)
555           {
556             setSideBoundToBlock( *block._side[ i ]._side );
557             nbBlockSides--;
558           }
559           _DUMP_("\t corners "<<
560                  block._side[ i ].cornerNode(0,0)->GetID() << ", " <<
561                  block._side[ i ].cornerNode(1,0)->GetID() << ", " <<
562                  block._side[ i ].cornerNode(1,1)->GetID() << ", " <<
563                  block._side[ i ].cornerNode(0,1)->GetID() << ", "<<endl);
564         }
565         else
566         {
567           _DUMP_("\t not found"<<endl);
568         }
569       }
570       if ( !ok )
571         block.clear();
572       else
573         nbBlocks++;
574     }
575     _DUMP_("Nb found blocks "<< nbBlocks <<endl);
576
577     if ( nbBlocks == 0 && _error.empty() )
578       return BAD_MESH_ERR;
579
580     return nbBlocks;
581   }
582
583   //================================================================================
584   /*!
585    * \brief Fill block side data starting from its corner quadrangle
586    */
587   //================================================================================
588
589   bool _Skin::fillSide( _BlockSide&             side,
590                         const SMDS_MeshElement* cornerQuad,
591                         const SMDS_MeshNode*    nCorner)
592   {
593     // Find out size of block side mesured in nodes and by the way find two rows
594     // of nodes in two directions.
595
596     int x, y, nbX, nbY;
597     const SMDS_MeshElement* firstQuad = cornerQuad;
598     {
599       // get a node on block edge
600       int iCorner = firstQuad->GetNodeIndex( nCorner );
601       const SMDS_MeshNode* nOnEdge = firstQuad->GetNode( (iCorner+1) % 4);
602
603       // find out size of block side
604       vector<const SMDS_MeshNode*> horRow1, horRow2, verRow1, verRow2;
605       if ( !fillRowsUntilCorner( firstQuad, nCorner, nOnEdge, horRow1, horRow2, true ) ||
606            !fillRowsUntilCorner( firstQuad, nCorner, nOnEdge, verRow1, verRow2, false ))
607         return false;
608       nbX = horRow1.size(), nbY = verRow1.size();
609
610       // store found nodes
611       side._index._xSize = horRow1.size();
612       side._index._ySize = verRow1.size();
613       side._grid.resize( side._index.size(), NULL );
614
615       for ( x = 0; x < horRow1.size(); ++x )
616       {
617         side.setNode( x, 0, horRow1[x] );
618         side.setNode( x, 1, horRow2[x] );
619       }
620       for ( y = 0; y < verRow1.size(); ++y )
621       {
622         side.setNode( 0, y, verRow1[y] );
623         side.setNode( 1, y, verRow2[y] );
624       }
625     }
626     // Find the rest nodes
627
628     y = 1; // y of the row to fill
629     TIDSortedElemSet emptySet, avoidSet;
630     while ( ++y < nbY )
631     {
632       // get next firstQuad in the next row of quadrangles
633       //
634       //          n2up
635       //     o---o               <- y row
636       //     |   |
637       //     o---o  o  o  o  o   <- found nodes
638       //n1down    n2down       
639       //
640       int i1down, i2down, i2up;
641       const SMDS_MeshNode* n1down = side.getNode( 0, y-1 );
642       const SMDS_MeshNode* n2down = side.getNode( 1, y-1 );
643       avoidSet.clear(); avoidSet.insert( firstQuad );
644       firstQuad = SMESH_MeshEditor::FindFaceInSet( n1down, n2down, emptySet, avoidSet,
645                                                    &i1down, &i2down);
646       if ( !isQuadrangle( firstQuad ))
647         return BAD_MESH_ERR;
648
649       const SMDS_MeshNode* n2up = oppositeNode( firstQuad, i1down );
650       avoidSet.clear(); avoidSet.insert( firstQuad );
651
652       // find the rest nodes in the y-th row by faces in the row
653
654       x = 1; 
655       while ( ++x < nbX )
656       {
657         const SMDS_MeshElement* quad = SMESH_MeshEditor::FindFaceInSet( n2up, n2down, emptySet,
658                                                                         avoidSet, &i2up, &i2down);
659         if ( !isQuadrangle( quad ))
660           return BAD_MESH_ERR;
661
662         n2up   = oppositeNode( quad, i2down );
663         n2down = oppositeNode( quad, i2up );
664         avoidSet.clear(); avoidSet.insert( quad );
665
666         side.setNode( x, y, n2up );
667       }
668     }
669
670     // check side validity
671     bool ok =
672       side.getCornerFace( side.getCornerNode( 0, 0 )) &&
673       side.getCornerFace( side.getCornerNode( 1, 0 )) &&
674       side.getCornerFace( side.getCornerNode( 0, 1 )) &&
675       side.getCornerFace( side.getCornerNode( 1, 1 ));
676
677     return ok;
678   }
679
680   //================================================================================
681   /*!
682    * \brief Try to find a block side adjacent to the given side by given edge
683    */
684   //================================================================================
685
686   _OrientedBlockSide _Skin::findBlockSide( EBoxSides startBlockSide,
687                                            EQuadEdge sharedSideEdge1,
688                                            EQuadEdge sharedSideEdge2,
689                                            bool      withGeometricAnalysis)
690   {
691     _Block& block = _blocks.back();
692     _OrientedBlockSide& side1 = block._side[ startBlockSide ];
693
694     // get corner nodes of the given block edge
695     SMESH_OrientedLink edge = side1.edge( sharedSideEdge1 );
696     const SMDS_MeshNode* n1 = edge.node1();
697     const SMDS_MeshNode* n2 = edge.node2();
698     if ( edge._reversed ) swap( n1, n2 );
699
700     // find all sides sharing both nodes n1 and n2
701     set< _BlockSide* > sidesOnEdge = _edge2sides[ edge ]; // copy a set
702
703     // exclude loaded sides of block from sidesOnEdge
704     for (int i = 0; i < NB_BLOCK_SIDES; ++i )
705       if ( block._side[ i ] )
706         sidesOnEdge.erase( block._side[ i ]._side );
707
708     int nbSidesOnEdge = sidesOnEdge.size();
709     _DUMP_("nbSidesOnEdge "<< nbSidesOnEdge << " " << n1->GetID() << "-" << n2->GetID() );
710     if ( nbSidesOnEdge == 0 )
711       return 0;
712
713     _BlockSide* foundSide = 0;
714     if ( nbSidesOnEdge == 1 )
715     {
716       foundSide = *sidesOnEdge.begin();
717     }
718     else
719     {
720       set< _BlockSide* >::iterator sideIt = sidesOnEdge.begin();
721       int nbLoadedSides = block.nbSides();
722       if ( nbLoadedSides > 1 )
723       {
724         // Find the side having more than 2 corners common with already loaded sides
725         for (; !foundSide && sideIt != sidesOnEdge.end(); ++sideIt )
726         {
727           _BlockSide* sideI = *sideIt;
728           int nbCommonCorners =
729             block._corners.count( sideI->getCornerNode(0,0)) +
730             block._corners.count( sideI->getCornerNode(1,0)) +
731             block._corners.count( sideI->getCornerNode(0,1)) +
732             block._corners.count( sideI->getCornerNode(1,1));
733           if ( nbCommonCorners > 2 )
734             foundSide = sideI;
735         }
736       }
737
738       if ( !foundSide )
739       {
740         if ( !withGeometricAnalysis ) return 0;
741
742         // Select one of found sides most close to startBlockSide
743
744         gp_XYZ p1 ( n1->X(),n1->Y(),n1->Z()),  p2 (n2->X(),n2->Y(),n2->Z());
745         gp_Vec p1p2( p1, p2 );
746
747         const SMDS_MeshElement* face1 = side1->getCornerFace( n1 );
748         gp_XYZ p1Op = SMESH_TNodeXYZ( oppositeNode( face1, face1->GetNodeIndex(n1)));
749         gp_Vec side1Dir( p1, p1Op );
750         gp_Ax2 pln( p1, p1p2, side1Dir ); // plane with normal p1p2 and X dir side1Dir
751         _DUMP_("  Select adjacent for "<< side1._side << " - side dir ("
752                << side1Dir.X() << ", " << side1Dir.Y() << ", " << side1Dir.Z() << ")" );
753
754         map < double , _BlockSide* > angleOfSide;
755         for (sideIt = sidesOnEdge.begin(); sideIt != sidesOnEdge.end(); ++sideIt )
756         {
757           _BlockSide* sideI = *sideIt;
758           const SMDS_MeshElement* faceI = sideI->getCornerFace( n1 );
759           gp_XYZ p1Op = SMESH_TNodeXYZ( oppositeNode( faceI, faceI->GetNodeIndex(n1)));
760           gp_Vec sideIDir( p1, p1Op );
761           // compute angle of (sideIDir projection to pln) and (X dir of pln)
762           gp_Vec2d sideIDirProj( sideIDir * pln.XDirection(), sideIDir * pln.YDirection());
763           double angle = sideIDirProj.Angle( gp::DX2d() );
764           if ( angle < 0 ) angle += 2 * PI; // angle [0-2*PI]
765           angleOfSide.insert( make_pair( angle, sideI ));
766           _DUMP_("  "<< sideI << " - side dir ("
767                  << sideIDir.X() << ", " << sideIDir.Y() << ", " << sideIDir.Z() << ")"
768                  << " angle " << angle);
769         }
770         if ( nbLoadedSides == 1 )
771         {
772           double angF = angleOfSide.begin()->first, angL = angleOfSide.rbegin()->first;
773           if ( angF > PI ) angF = 2*PI - angF;
774           if ( angL > PI ) angL = 2*PI - angL;
775           foundSide = angF < angL ? angleOfSide.begin()->second : angleOfSide.rbegin()->second;
776         }
777         else
778         {
779           gp_XYZ gc(0,0,0); // gravity center of already loaded block sides
780           for (int i = 0; i < NB_BLOCK_SIDES; ++i )
781             if ( block._side[ i ] )
782               gc += block._side[ i ]._side->getGC();
783           gc /= nbLoadedSides;
784
785           gp_Vec gcDir( p1, gc );
786           gp_Vec2d gcDirProj( gcDir * pln.XDirection(), gcDir * pln.YDirection());
787           double gcAngle = gcDirProj.Angle( gp::DX2d() );
788           foundSide = gcAngle < 0 ? angleOfSide.rbegin()->second : angleOfSide.begin()->second;
789         }
790       }
791       _DUMP_("  selected "<< foundSide );
792     }
793
794     // Orient the found side correctly
795
796     // corners of found side corresponding to nodes n1 and n2
797     bool xMax1, yMax1, xMax2, yMax2;
798     if ( !getEdgeEnds( sharedSideEdge2, xMax1, yMax1, xMax2, yMax2 ))
799       return error(SMESH_Comment("Internal error at ")<<__FILE__<<":"<<__LINE__),
800         _OrientedBlockSide(0);
801
802     for ( int ori = 0; ori < _OrientedIndexer::MAX_ORI+1; ++ori )
803     {
804       _OrientedBlockSide orientedSide( foundSide, ori );
805       const SMDS_MeshNode* n12 = orientedSide.cornerNode( xMax1, yMax1);
806       const SMDS_MeshNode* n22 = orientedSide.cornerNode( xMax2, yMax2);
807       if ( n1 == n12 && n2 == n22 )
808         return orientedSide;
809     }
810     error(SMESH_Comment("Failed to orient a block side found by edge ")<<sharedSideEdge1
811           << " of side " << startBlockSide
812           << " of block " << _blocks.size());
813     return 0;
814   }
815
816   //================================================================================
817   /*!
818    * \brief: Fill rows (which are actually columns,if !alongN1N2) of nodes starting
819    * from the given quadrangle until another block corner encounters.
820    *  n1 and n2 are at bottom of quad, n1 is at block corner.
821    */
822   //================================================================================
823
824   bool _Skin::fillRowsUntilCorner(const SMDS_MeshElement*       quad,
825                                   const SMDS_MeshNode*          n1,
826                                   const SMDS_MeshNode*          n2,
827                                   vector<const SMDS_MeshNode*>& row1,
828                                   vector<const SMDS_MeshNode*>& row2,
829                                   const bool                    alongN1N2 )
830   {
831     const SMDS_MeshNode* corner1 = n1;
832
833     // Store nodes of quad in the rows and find new n1 and n2 to get
834     // the next face so that new n2 is on block edge
835     int i1 = quad->GetNodeIndex( n1 );
836     int i2 = quad->GetNodeIndex( n2 );
837     row1.clear(); row2.clear();
838     row1.push_back( n1 );
839     if ( alongN1N2 )
840     {
841       row1.push_back( n2 );
842       row2.push_back( oppositeNode( quad, i2 ));
843       row2.push_back( n1 = oppositeNode( quad, i1 ));
844     }
845     else
846     {
847       row2.push_back( n2 );
848       row1.push_back( n2 = oppositeNode( quad, i2 ));
849       row2.push_back( n1 = oppositeNode( quad, i1 ));
850     }
851
852     _NodeDescriptor nDesc( row1[1] );
853     if ( nDesc == _NodeDescriptor( row1[0] ))
854       return true; // row of 2 nodes
855
856     // Find the rest nodes
857     TIDSortedElemSet emptySet, avoidSet;
858     //while ( !isCornerNode( n2 ))
859     while ( nDesc == _NodeDescriptor( n2 ))
860     {
861       avoidSet.clear(); avoidSet.insert( quad );
862       quad = SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2 );
863       if ( !isQuadrangle( quad ))
864         return BAD_MESH_ERR;
865
866       row1.push_back( n2 = oppositeNode( quad, i1 ));
867       row2.push_back( n1 = oppositeNode( quad, i2 ));
868     }
869     return n1 != corner1;
870   }
871
872   //================================================================================
873   /*!
874    * \brief Return a corner face by a corner node
875    */
876   //================================================================================
877
878   const SMDS_MeshElement* _BlockSide::getCornerFace(const SMDS_MeshNode* cornerNode) const
879   {
880     int x, y, isXMax, isYMax, found = 0;
881     for ( isXMax = 0; isXMax < 2; ++isXMax )
882     {
883       for ( isYMax = 0; isYMax < 2; ++isYMax )
884       {
885         x = isXMax ? _index._xSize-1 : 0;
886         y = isYMax ? _index._ySize-1 : 0;
887         found = ( getNode(x,y) == cornerNode );
888         if ( found ) break;
889       }
890       if ( found ) break;
891     }
892     if ( !found ) return 0;
893     int dx = isXMax ? -1 : +1;
894     int dy = isYMax ? -1 : +1;
895     const SMDS_MeshNode* n1 = getNode(x,y);
896     const SMDS_MeshNode* n2 = getNode(x+dx,y);
897     const SMDS_MeshNode* n3 = getNode(x,y+dy);
898     const SMDS_MeshNode* n4 = getNode(x+dx,y+dy);
899     return SMDS_Mesh::FindFace(n1, n2, n3, n4 );
900   }
901
902   //================================================================================
903   /*!
904    * \brief Checks own validity
905    */
906   //================================================================================
907
908   bool _Block::isValid() const
909   {
910     bool ok = ( nbSides() == 6 );
911
912     // check only corners depending on side selection
913     EBoxSides adjacent[4] = { B_BOTTOM, B_RIGHT, B_TOP, B_LEFT };
914     EQuadEdge edgeAdj [4] = { Q_TOP,    Q_RIGHT, Q_TOP, Q_RIGHT };
915     EQuadEdge edgeBack[4] = { Q_BOTTOM, Q_RIGHT, Q_TOP, Q_LEFT };
916
917     for ( int i=0; ok && i < NB_QUAD_SIDES; ++i )
918     { 
919       SMESH_OrientedLink eBack = _side[ B_BACK      ].edge( edgeBack[i] );
920       SMESH_OrientedLink eAdja = _side[ adjacent[i] ].edge( edgeAdj[i] );
921       ok = ( eBack == eAdja );
922     }
923     return ok;
924   }
925
926 } // namespace
927
928 //=======================================================================
929 //function : StdMeshers_HexaFromSkin_3D
930 //purpose  : 
931 //=======================================================================
932
933 StdMeshers_HexaFromSkin_3D::StdMeshers_HexaFromSkin_3D(int hypId, int studyId, SMESH_Gen* gen)
934   :SMESH_3D_Algo(hypId, studyId, gen)
935 {
936   MESSAGE("StdMeshers_HexaFromSkin_3D::StdMeshers_HexaFromSkin_3D");
937   _name = "HexaFromSkin_3D";
938 }
939
940 StdMeshers_HexaFromSkin_3D::~StdMeshers_HexaFromSkin_3D()
941 {
942   MESSAGE("StdMeshers_HexaFromSkin_3D::~StdMeshers_HexaFromSkin_3D");
943 }
944
945 //================================================================================
946 /*!
947  * \brief Main method, which generates hexaheda
948  */
949 //================================================================================
950
951 bool StdMeshers_HexaFromSkin_3D::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
952 {
953   _Skin skin;
954   int nbBlocks = skin.findBlocks(aMesh);
955   if ( nbBlocks == 0 )
956     return error( skin.error());
957
958   vector< vector< const SMDS_MeshNode* > > columns;
959   int x, xSize, y, ySize, z, zSize;
960   _Indexer colIndex;
961
962   for ( int i = 0; i < nbBlocks; ++i )
963   {
964     const _Block& block = skin.getBlock( i );
965
966     // ------------------------------------------
967     // Fill columns of nodes with existing nodes
968     // ------------------------------------------
969
970     xSize = block.getSide(B_BOTTOM).getHoriSize();
971     ySize = block.getSide(B_BOTTOM).getVertSize();
972     zSize = block.getSide(B_FRONT ).getVertSize();
973     int X = xSize - 1, Y = ySize - 1, Z = zSize - 1;
974     colIndex = _Indexer( xSize, ySize );
975     columns.resize( colIndex.size() );
976
977     // fill node columns by front and back box sides
978     for ( x = 0; x < xSize; ++x ) {
979       vector< const SMDS_MeshNode* >& column0 = columns[ colIndex( x, 0 )];
980       vector< const SMDS_MeshNode* >& column1 = columns[ colIndex( x, Y )];
981       column0.resize( zSize );
982       column1.resize( zSize );
983       for ( z = 0; z < zSize; ++z ) {
984         column0[ z ] = block.getSide(B_FRONT).node( x, z );
985         column1[ z ] = block.getSide(B_BACK) .node( x, z );
986       }
987     }
988     // fill node columns by left and right box sides
989     for ( y = 1; y < ySize-1; ++y ) {
990       vector< const SMDS_MeshNode* >& column0 = columns[ colIndex( 0, y )];
991       vector< const SMDS_MeshNode* >& column1 = columns[ colIndex( X, y )];
992       column0.resize( zSize );
993       column1.resize( zSize );
994       for ( z = 0; z < zSize; ++z ) {
995         column0[ z ] = block.getSide(B_LEFT) .node( y, z );
996         column1[ z ] = block.getSide(B_RIGHT).node( y, z );
997       }
998     }
999     // get nodes from top and bottom box sides
1000     for ( x = 1; x < xSize-1; ++x ) {
1001       for ( y = 1; y < ySize-1; ++y ) {
1002         vector< const SMDS_MeshNode* >& column = columns[ colIndex( x, y )];
1003         column.resize( zSize );
1004         column.front() = block.getSide(B_BOTTOM).node( x, y );
1005         column.back()  = block.getSide(B_TOP)   .node( x, y );
1006       }
1007     }
1008
1009     // ----------------------------
1010     // Add internal nodes of a box
1011     // ----------------------------
1012     // projection points of internal nodes on box subshapes by which
1013     // coordinates of internal nodes are computed
1014     vector<gp_XYZ> pointOnShape( SMESH_Block::ID_Shell );
1015
1016     // projections on vertices are constant
1017     pointOnShape[ SMESH_Block::ID_V000 ] = block.getSide(B_BOTTOM).xyz( 0, 0 );
1018     pointOnShape[ SMESH_Block::ID_V100 ] = block.getSide(B_BOTTOM).xyz( X, 0 );
1019     pointOnShape[ SMESH_Block::ID_V010 ] = block.getSide(B_BOTTOM).xyz( 0, Y );
1020     pointOnShape[ SMESH_Block::ID_V110 ] = block.getSide(B_BOTTOM).xyz( X, Y );
1021     pointOnShape[ SMESH_Block::ID_V001 ] = block.getSide(B_TOP).xyz( 0, 0 );
1022     pointOnShape[ SMESH_Block::ID_V101 ] = block.getSide(B_TOP).xyz( X, 0 );
1023     pointOnShape[ SMESH_Block::ID_V011 ] = block.getSide(B_TOP).xyz( 0, Y );
1024     pointOnShape[ SMESH_Block::ID_V111 ] = block.getSide(B_TOP).xyz( X, Y );
1025
1026     for ( x = 1; x < xSize-1; ++x )
1027     {
1028       gp_XYZ params; // normalized parameters of internal node within a unit box
1029       params.SetCoord( 1, x / double(X) );
1030       for ( y = 1; y < ySize-1; ++y )
1031       {
1032         params.SetCoord( 2, y / double(Y) );
1033         // column to fill during z loop
1034         vector< const SMDS_MeshNode* >& column = columns[ colIndex( x, y )];
1035         // projections on horizontal edges
1036         pointOnShape[ SMESH_Block::ID_Ex00 ] = block.getSide(B_BOTTOM).xyz( x, 0 );
1037         pointOnShape[ SMESH_Block::ID_Ex10 ] = block.getSide(B_BOTTOM).xyz( x, Y );
1038         pointOnShape[ SMESH_Block::ID_E0y0 ] = block.getSide(B_BOTTOM).xyz( 0, y );
1039         pointOnShape[ SMESH_Block::ID_E1y0 ] = block.getSide(B_BOTTOM).xyz( X, y );
1040         pointOnShape[ SMESH_Block::ID_Ex01 ] = block.getSide(B_TOP).xyz( x, 0 );
1041         pointOnShape[ SMESH_Block::ID_Ex11 ] = block.getSide(B_TOP).xyz( x, Y );
1042         pointOnShape[ SMESH_Block::ID_E0y1 ] = block.getSide(B_TOP).xyz( 0, y );
1043         pointOnShape[ SMESH_Block::ID_E1y1 ] = block.getSide(B_TOP).xyz( X, y );
1044         // projections on horizontal sides
1045         pointOnShape[ SMESH_Block::ID_Fxy0 ] = block.getSide(B_BOTTOM).xyz( x, y );
1046         pointOnShape[ SMESH_Block::ID_Fxy1 ] = block.getSide(B_TOP)   .xyz( x, y );
1047         for ( z = 1; z < zSize-1; ++z ) // z loop
1048         {
1049           params.SetCoord( 3, z / double(Z) );
1050           // projections on vertical edges
1051           pointOnShape[ SMESH_Block::ID_E00z ] = block.getSide(B_FRONT).xyz( 0, z );    
1052           pointOnShape[ SMESH_Block::ID_E10z ] = block.getSide(B_FRONT).xyz( X, z );    
1053           pointOnShape[ SMESH_Block::ID_E01z ] = block.getSide(B_BACK).xyz( 0, z );    
1054           pointOnShape[ SMESH_Block::ID_E11z ] = block.getSide(B_BACK).xyz( X, z );
1055           // projections on vertical sides
1056           pointOnShape[ SMESH_Block::ID_Fx0z ] = block.getSide(B_FRONT).xyz( x, z );    
1057           pointOnShape[ SMESH_Block::ID_Fx1z ] = block.getSide(B_BACK) .xyz( x, z );    
1058           pointOnShape[ SMESH_Block::ID_F0yz ] = block.getSide(B_LEFT) .xyz( y, z );    
1059           pointOnShape[ SMESH_Block::ID_F1yz ] = block.getSide(B_RIGHT).xyz( y, z );
1060
1061           // compute internal node coordinates
1062           gp_XYZ coords;
1063           SMESH_Block::ShellPoint( params, pointOnShape, coords );
1064           column[ z ] = aHelper->AddNode( coords.X(), coords.Y(), coords.Z() );
1065
1066 #ifdef DEB_GRID
1067           // debug
1068           //cout << "----------------------------------------------------------------------"<<endl;
1069           //for ( int id = SMESH_Block::ID_V000; id < SMESH_Block::ID_Shell; ++id)
1070           //{
1071           //  gp_XYZ p = pointOnShape[ id ];
1072           //  SMESH_Block::DumpShapeID( id,cout)<<" ( "<<p.X()<<", "<<p.Y()<<", "<<p.Z()<<" )"<<endl;
1073           //}
1074           //cout << "Params: ( "<< params.X()<<", "<<params.Y()<<", "<<params.Z()<<" )"<<endl;
1075           //cout << "coords: ( "<< coords.X()<<", "<<coords.Y()<<", "<<coords.Z()<<" )"<<endl;
1076 #endif
1077         }
1078       }
1079     }
1080     // ----------------
1081     // Add hexahedrons
1082     // ----------------
1083
1084     // find out orientation by a least distorted hexahedron (issue 0020855);
1085     // the last is defined by evaluating sum of face normals of 8 corner hexahedrons
1086     double badness = numeric_limits<double>::max();
1087     bool isForw = true;
1088     for ( int xMax = 0; xMax < 2; ++xMax )
1089       for ( int yMax = 0; yMax < 2; ++yMax )
1090         for ( int zMax = 0; zMax < 2; ++zMax )
1091         {
1092           x = xMax ? xSize-1 : 1;
1093           y = yMax ? ySize-1 : 1;
1094           z = zMax ? zSize-1 : 1;
1095           vector< const SMDS_MeshNode* >& col00 = columns[ colIndex( x-1, y-1 )];
1096           vector< const SMDS_MeshNode* >& col10 = columns[ colIndex( x  , y-1 )];
1097           vector< const SMDS_MeshNode* >& col01 = columns[ colIndex( x-1, y   )];
1098           vector< const SMDS_MeshNode* >& col11 = columns[ colIndex( x  , y )];
1099           
1100           const SMDS_MeshNode* n000 = col00[z-1];
1101           const SMDS_MeshNode* n100 = col10[z-1];
1102           const SMDS_MeshNode* n010 = col01[z-1];
1103           const SMDS_MeshNode* n110 = col11[z-1];
1104           const SMDS_MeshNode* n001 = col00[z];
1105           const SMDS_MeshNode* n101 = col10[z];
1106           const SMDS_MeshNode* n011 = col01[z];
1107           const SMDS_MeshNode* n111 = col11[z];
1108           SMDS_VolumeOfNodes probeVolume (n000,n010,n110,n100,
1109                                           n001,n011,n111,n101);
1110           SMDS_VolumeTool volTool( &probeVolume );
1111           double Nx=0.,Ny=0.,Nz=0.;
1112           for ( int iFace = 0; iFace < volTool.NbFaces(); ++iFace )
1113           {
1114             double nx,ny,nz;
1115             volTool.GetFaceNormal( iFace, nx,ny,nz );
1116             Nx += nx;
1117             Ny += ny;
1118             Nz += nz;
1119           }
1120           double quality = Nx*Nx + Ny*Ny + Nz*Nz;
1121           if ( quality < badness )
1122           {
1123             badness = quality;
1124             isForw = volTool.IsForward();
1125           }
1126         }
1127
1128     // add elements
1129     for ( x = 0; x < xSize-1; ++x ) {
1130       for ( y = 0; y < ySize-1; ++y ) {
1131         vector< const SMDS_MeshNode* >& col00 = columns[ colIndex( x, y )];
1132         vector< const SMDS_MeshNode* >& col10 = columns[ colIndex( x+1, y )];
1133         vector< const SMDS_MeshNode* >& col01 = columns[ colIndex( x, y+1 )];
1134         vector< const SMDS_MeshNode* >& col11 = columns[ colIndex( x+1, y+1 )];
1135         // bottom face normal of a hexa mush point outside the volume
1136         if ( isForw )
1137           for ( z = 0; z < zSize-1; ++z )
1138             aHelper->AddVolume(col00[z],   col01[z],   col11[z],   col10[z],
1139                                col00[z+1], col01[z+1], col11[z+1], col10[z+1]);
1140         else
1141           for ( z = 0; z < zSize-1; ++z )
1142             aHelper->AddVolume(col00[z],   col10[z],   col11[z],   col01[z],
1143                                col00[z+1], col10[z+1], col11[z+1], col01[z+1]);
1144       }
1145     }
1146   } // loop on blocks
1147
1148   return true;
1149 }
1150
1151 //================================================================================
1152 /*!
1153  * \brief Evaluate nb of hexa
1154  */
1155 //================================================================================
1156
1157 bool StdMeshers_HexaFromSkin_3D::Evaluate(SMESH_Mesh &         aMesh,
1158                                           const TopoDS_Shape & aShape,
1159                                           MapShapeNbElems&     aResMap)
1160 {
1161   _Skin skin;
1162   int nbBlocks = skin.findBlocks(aMesh);
1163   if ( nbBlocks == 0 )
1164     return error( skin.error());
1165
1166   bool secondOrder = aMesh.NbFaces( ORDER_QUADRATIC );
1167
1168   int entity = secondOrder ? SMDSEntity_Quad_Hexa : SMDSEntity_Hexa;
1169   vector<int>& nbByType = aResMap[ aMesh.GetSubMesh( aShape )];
1170   if ( entity >= nbByType.size() )
1171     nbByType.resize( SMDSEntity_Last, 0 );
1172
1173   for ( int i = 0; i < nbBlocks; ++i )
1174   {
1175     const _Block& block = skin.getBlock( i );
1176
1177     int nbX = block.getSide(B_BOTTOM).getHoriSize();
1178     int nbY = block.getSide(B_BOTTOM).getVertSize();
1179     int nbZ = block.getSide(B_FRONT ).getVertSize();
1180
1181     int nbHexa  = (nbX-1) * (nbY-1) * (nbZ-1);
1182     int nbNodes = (nbX-2) * (nbY-2) * (nbZ-2);
1183     if ( secondOrder )
1184       nbNodes +=
1185         (nbX-2) * (nbY-2) * (nbZ-1) +
1186         (nbX-2) * (nbY-1) * (nbZ-2) +
1187         (nbX-1) * (nbY-2) * (nbZ-2);
1188
1189
1190     nbByType[ entity ] += nbHexa;
1191     nbByType[ SMDSEntity_Node ] += nbNodes;
1192   }
1193
1194   return true;
1195 }
1196
1197 //================================================================================
1198 /*!
1199  * \brief Abstract method must be defined but does nothing
1200  */
1201 //================================================================================
1202
1203 bool StdMeshers_HexaFromSkin_3D::CheckHypothesis(SMESH_Mesh&, const TopoDS_Shape&,
1204                                                  Hypothesis_Status& aStatus)
1205 {
1206   aStatus = SMESH_Hypothesis::HYP_OK;
1207   return true;
1208 }
1209
1210 //================================================================================
1211 /*!
1212  * \brief Abstract method must be defined but just reports an error as this
1213  *  algo is not intended to work with shapes
1214  */
1215 //================================================================================
1216
1217 bool StdMeshers_HexaFromSkin_3D::Compute(SMESH_Mesh&, const TopoDS_Shape&)
1218 {
1219   return error("Algorithm can't work with geometrical shapes");
1220 }
1221