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