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