Version: 6.3.1

src/StdMeshers/StdMeshers_HexaFromSkin_3D.cxx

Go to the documentation of this file.
00001 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
00002 //
00003 // This library is free software; you can redistribute it and/or
00004 // modify it under the terms of the GNU Lesser General Public
00005 // License as published by the Free Software Foundation; either
00006 // version 2.1 of the License.
00007 //
00008 // This library is distributed in the hope that it will be useful,
00009 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00010 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00011 // Lesser General Public License for more details.
00012 //
00013 // You should have received a copy of the GNU Lesser General Public
00014 // License along with this library; if not, write to the Free Software
00015 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00016 //
00017 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00018 //
00019 
00020 // File      : StdMeshers_HexaFromSkin_3D.cxx
00021 // Created   : Wed Jan 27 12:28:07 2010
00022 // Author    : Edward AGAPOV (eap)
00023 //
00024 #include "StdMeshers_HexaFromSkin_3D.hxx"
00025 
00026 #include "SMDS_VolumeOfNodes.hxx"
00027 #include "SMDS_VolumeTool.hxx"
00028 #include "SMESH_Block.hxx"
00029 #include "SMESH_MesherHelper.hxx"
00030 
00031 #include <gp_Ax2.hxx>
00032 
00033 //#include "utilities.h"
00034 #include <limits>
00035 
00036 // Define error message
00037 #ifdef _DEBUG_
00038 #define BAD_MESH_ERR \
00039   error(SMESH_Comment("Can't detect block-wise structure of the input 2D mesh.\n" \
00040                       __FILE__ ":" )<<__LINE__)
00041 #else
00042 #define BAD_MESH_ERR \
00043   error(SMESH_Comment("Can't detect block-wise structure of the input 2D mesh"))
00044 #endif
00045 
00046 // Debug output
00047 #define _DUMP_(msg) //cout << msg << endl
00048 
00049 
00050 
00051 namespace
00052 {
00053   enum EBoxSides 
00054     {
00055       B_BOTTOM=0, B_RIGHT, B_TOP, B_LEFT, B_FRONT, B_BACK, NB_BLOCK_SIDES
00056     };
00057 //   const char* SBoxSides[] = //!< names of block sides -- needed for DEBUG only
00058 //     {
00059 //       "BOTTOM", "RIGHT", "TOP", "LEFT", "FRONT", "BACK", "UNDEFINED"
00060 //     };
00061   enum EQuadEdge 
00062     {
00063       Q_BOTTOM = 0, Q_RIGHT, Q_TOP, Q_LEFT, NB_QUAD_SIDES
00064     };
00065 
00066 
00067   //================================================================================
00071   //================================================================================
00072 
00073   bool getEdgeEnds(EQuadEdge edge, bool& xMax1, bool& yMax1, bool& xMax2, bool& yMax2 )
00074   {
00075     xMax1=0, yMax1=0, xMax2=1, yMax2=1;
00076     switch( edge )
00077     {
00078     case Q_BOTTOM: yMax2 = 0; break;
00079     case Q_RIGHT:  xMax1 = 1; break;
00080     case Q_TOP:    yMax1 = 1; break;
00081     case Q_LEFT:   xMax2 = 0; break;
00082     default:
00083       return false;
00084     }
00085     return true;
00086   }
00087 
00088   //================================================================================
00092   struct _NodeDescriptor
00093   {
00094     int nbInverseFaces, nbNodesInInverseFaces;
00095 
00096     _NodeDescriptor(const SMDS_MeshNode* n): nbInverseFaces(0), nbNodesInInverseFaces(0)
00097     {
00098       if ( n )
00099       {
00100         set<const SMDS_MeshNode*> nodes;
00101         SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator(SMDSAbs_Face );
00102         while ( fIt->more() )
00103         {
00104           const SMDS_MeshElement* face = fIt->next();
00105           nodes.insert( face->begin_nodes(), face->end_nodes() );
00106           nbInverseFaces++;
00107         }
00108         nbNodesInInverseFaces = nodes.size();
00109       }
00110     }
00111     bool operator==(const _NodeDescriptor& other) const
00112     {
00113       return
00114         nbInverseFaces == other.nbInverseFaces &&
00115         nbNodesInInverseFaces == other.nbNodesInInverseFaces;
00116     }
00117   };
00118 
00119   //================================================================================
00125   //================================================================================
00126 
00127   bool isCornerNode( const SMDS_MeshNode* n )
00128   {
00129     return n && n->NbInverseElements( SMDSAbs_Face ) % 2;
00130   }
00131 
00132   //================================================================================
00136   //================================================================================
00137 
00138   bool isQuadrangle(const SMDS_MeshElement* e)
00139   {
00140     return ( e && e->NbCornerNodes() == 4 );
00141   }
00142 
00143   //================================================================================
00147   //================================================================================
00148 
00149   const SMDS_MeshNode* oppositeNode(const SMDS_MeshElement* quad, int iNode)
00150   {
00151     return quad->GetNode( (iNode+2) % 4 );
00152   }
00153 
00154   //================================================================================
00158   struct _Indexer
00159   {
00160     int _xSize, _ySize;
00161     _Indexer( int xSize=0, int ySize=0 ): _xSize(xSize), _ySize(ySize) {}
00162     int size() const { return _xSize * _ySize; }
00163     int operator()(int x, int y) const { return y * _xSize + x; }
00164   };
00165   //================================================================================
00169   class _OrientedIndexer : public _Indexer
00170   {
00171   public:
00172     enum OriFlags 
00173       {
00174         REV_X = 1, REV_Y = 2, SWAP_XY = 4, MAX_ORI = REV_X|REV_Y|SWAP_XY
00175       };
00176     _OrientedIndexer( const _Indexer& indexer, const int oriFlags ):
00177       _Indexer( indexer._xSize, indexer._ySize ),
00178       _xSize (indexer._xSize), _ySize(indexer._ySize),
00179       _xRevFun((oriFlags & REV_X) ? & reverse : & lazy),
00180       _yRevFun((oriFlags & REV_Y) ? & reverse : & lazy),
00181       _swapFun((oriFlags & SWAP_XY ) ? & swap : & lazy)
00182     {
00183       (*_swapFun)( _xSize, _ySize );
00184     }
00186     int operator()(int x, int y) const
00187     {
00188       (*_xRevFun)( x, const_cast<int&>( _xSize ));
00189       (*_yRevFun)( y, const_cast<int&>( _ySize ));
00190       (*_swapFun)( x, y );
00191       return _Indexer::operator()(x,y);
00192     }
00194     int corner(bool xMax, bool yMax) const
00195     {
00196       int x = xMax, y = yMax, size = 2;
00197       (*_xRevFun)( x, size );
00198       (*_yRevFun)( y, size );
00199       (*_swapFun)( x, y );
00200       return _Indexer::operator()(x ? _Indexer::_xSize-1 : 0 , y ? _Indexer::_ySize-1 : 0);
00201     }
00202     int xSize() const { return _xSize; }
00203     int ySize() const { return _ySize; }
00204   private:
00205     _Indexer _indexer;
00206     int _xSize, _ySize;
00207 
00208     typedef void (*TFun)(int& x, int& y);
00209     TFun _xRevFun, _yRevFun, _swapFun;
00210     
00211     static void lazy(int&, int&) {}
00212     static void reverse(int& x, int& size) { x = size - x - 1; }
00213     static void swap(int& x, int& y) { std::swap(x,y); }
00214   };
00215   //================================================================================
00219   struct _BlockSide
00220   {
00221     vector<const SMDS_MeshNode*> _grid;
00222     _Indexer                     _index;
00223     int                          _nbBlocksExpected;
00224     int                          _nbBlocksFound;
00225 
00226 #ifdef _DEBUG_ // want to get SIGSEGV in case of invalid index
00227 #define _grid_access_(pobj, i) pobj->_grid[ ((i) < pobj->_grid.size()) ? i : int(1e100)]
00228 #else
00229 #define _grid_access_(pobj, i) pobj->_grid[ i ]
00230 #endif
00231 
00232     const SMDS_MeshNode* getNode(int x, int y) const { return _grid_access_(this, _index( x,y ));}
00234     void setNode(int x, int y, const SMDS_MeshNode* n) { _grid_access_(this, _index( x,y )) = n; }
00236     SMESH_OrientedLink getEdge(EQuadEdge edge) const
00237     {
00238       bool x1, y1, x2, y2; getEdgeEnds( edge, x1, y1, x2, y2 );
00239       return SMESH_OrientedLink( getCornerNode ( x1, y1 ), getCornerNode( x2, y2 ));
00240     }
00242     const SMDS_MeshNode* getCornerNode(bool isXMax, bool isYMax) const
00243     {
00244       return getNode( isXMax ? _index._xSize-1 : 0 , isYMax ? _index._ySize-1 : 0 );
00245     }
00246     const SMDS_MeshElement* getCornerFace(const SMDS_MeshNode* cornerNode) const;
00248     bool isBound() const { return _nbBlocksExpected <= _nbBlocksFound; }
00250     gp_XYZ getXYZ(int x, int y) const { return SMESH_TNodeXYZ( getNode( x, y )); }
00252     gp_XYZ getGC() const
00253     {
00254       gp_XYZ xyz =
00255         getXYZ( 0, 0 ) +
00256         getXYZ( _index._xSize-1, 0 ) +
00257         getXYZ( 0, _index._ySize-1 ) +
00258         getXYZ( _index._xSize-1, _index._ySize-1 ) +
00259         getXYZ( _index._xSize/2, _index._ySize/2 );
00260       return xyz / 5;
00261     }
00263     int getNbFaces() const { return (_index._xSize-1) * (_index._ySize-1); }
00264   };
00265   //================================================================================
00269   struct _OrientedBlockSide
00270   {
00271     _BlockSide*       _side;
00272     _OrientedIndexer  _index;
00273 
00274     _OrientedBlockSide( _BlockSide* side=0, const int oriFlags=0 ):
00275       _side(side), _index(side ? side->_index : _Indexer(), oriFlags ) {}
00277     gp_XYZ xyz(int x, int y) const
00278     {
00279       return SMESH_TNodeXYZ( _grid_access_(_side, _index( x, y )) );
00280     }
00282     const SMDS_MeshNode* node(int x, int y) const
00283     {
00284       int i = _index( x, y );
00285       return ( i < 0 || i >= _side->_grid.size()) ? 0 : _side->_grid[i];
00286     }
00288     SMESH_OrientedLink edge(EQuadEdge edge) const
00289     {
00290       bool x1, y1, x2, y2; getEdgeEnds( edge, x1, y1, x2, y2 );
00291       return SMESH_OrientedLink( cornerNode ( x1, y1 ), cornerNode( x2, y2 ));
00292     }
00294     const SMDS_MeshNode* cornerNode(bool isXMax, bool isYMax) const
00295     {
00296       return _grid_access_(_side, _index.corner( isXMax, isYMax ));
00297     }
00299     int getHoriSize() const { return _index.xSize(); }
00300     int getVertSize() const  { return _index.ySize(); }
00302     operator bool() const { return _side; }
00304     const _BlockSide* operator->() const { return _side; }
00305     _BlockSide* operator->() { return _side; }
00306   };
00307   //================================================================================
00311   struct _Block
00312   {
00313     _OrientedBlockSide        _side[6]; // 6 sides of a sub-block
00314     set<const SMDS_MeshNode*> _corners;
00315 
00316     const _OrientedBlockSide& getSide(int i) const { return _side[i]; }
00317     bool setSide( int i, const _OrientedBlockSide& s)
00318     {
00319       if (( _side[i] = s ))
00320       {
00321         _corners.insert( s.cornerNode(0,0));
00322         _corners.insert( s.cornerNode(1,0));
00323         _corners.insert( s.cornerNode(0,1));
00324         _corners.insert( s.cornerNode(1,1));
00325       }
00326       return s;
00327     }
00328     void clear() { for (int i=0;i<6;++i) _side[i]=0; _corners.clear(); }
00329     bool hasSide( const _OrientedBlockSide& s) const
00330     {
00331       if ( s ) for (int i=0;i<6;++i) if ( _side[i] && _side[i]._side == s._side ) return true;
00332       return false;
00333     }
00334     int nbSides() const { int n=0; for (int i=0;i<6;++i) if ( _side[i] ) ++n; return n; }
00335     bool isValid() const;
00336   };
00337   //================================================================================
00341   class _Skin
00342   {
00343   public:
00344 
00345     int findBlocks(SMESH_Mesh& mesh);
00347     const _Block& getBlock(int i) const { return _blocks[i]; }
00349     const SMESH_Comment& error() const { return _error; }
00350 
00351   private:
00352     bool fillSide( _BlockSide&             side,
00353                    const SMDS_MeshElement* cornerQuad,
00354                    const SMDS_MeshNode*    cornerNode);
00355     bool fillRowsUntilCorner(const SMDS_MeshElement* quad,
00356                              const SMDS_MeshNode*    n1,
00357                              const SMDS_MeshNode*    n2,
00358                              vector<const SMDS_MeshNode*>& verRow1,
00359                              vector<const SMDS_MeshNode*>& verRow2,
00360                              bool alongN1N2 );
00361     _OrientedBlockSide findBlockSide( EBoxSides startBlockSide,
00362                                       EQuadEdge sharedSideEdge1,
00363                                       EQuadEdge sharedSideEdge2,
00364                                       bool      withGeometricAnalysis);
00366     void setSideBoundToBlock( _BlockSide& side )
00367     {
00368       if ( side._nbBlocksFound++, side.isBound() )
00369         for ( int e = 0; e < int(NB_QUAD_SIDES); ++e )
00370           _edge2sides[ side.getEdge( (EQuadEdge) e ) ].erase( &side );
00371     }
00373     int error(const SMESH_Comment& reason) { _error = reason; return 0; }
00374 
00375     SMESH_Comment      _error;
00376 
00377     list< _BlockSide > _allSides;
00378     vector< _Block >   _blocks;
00379 
00380     //map< const SMDS_MeshNode*, set< _BlockSide* > > _corner2sides;
00381     map< SMESH_OrientedLink, set< _BlockSide* > > _edge2sides;
00382   };
00383 
00384   //================================================================================
00388   //================================================================================
00389 
00390   int _Skin::findBlocks(SMESH_Mesh& mesh)
00391   {
00392     SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
00393 
00394     // Find a node at any block corner
00395 
00396     SMDS_NodeIteratorPtr nIt = meshDS->nodesIterator(/*idInceasingOrder=*/true);
00397     if ( !nIt->more() ) return error("Empty mesh");
00398 
00399     const SMDS_MeshNode* nCorner = 0;
00400     while ( nIt->more() )
00401     {
00402       nCorner = nIt->next();
00403       if ( isCornerNode( nCorner ))
00404         break;
00405       else
00406         nCorner = 0;
00407     }
00408     if ( !nCorner )
00409       return BAD_MESH_ERR;
00410 
00411     // --------------------------------------------------------------------
00412     // Find all block sides starting from mesh faces sharing the corner node
00413     // --------------------------------------------------------------------
00414 
00415     int nbFacesOnSides = 0;
00416     TIDSortedElemSet cornerFaces; // corner faces of found _BlockSide's
00417     list< const SMDS_MeshNode* > corners( 1, nCorner );
00418     list< const SMDS_MeshNode* >::iterator corner = corners.begin();
00419     while ( corner != corners.end() )
00420     {
00421       SMDS_ElemIteratorPtr faceIt = (*corner)->GetInverseElementIterator( SMDSAbs_Face );
00422       while ( faceIt->more() )
00423       {
00424         const SMDS_MeshElement* face = faceIt->next();
00425         if ( !cornerFaces.insert( face ).second )
00426           continue; // already loaded block side
00427 
00428         if ( !isQuadrangle( face ))
00429           return error("Non-quadrangle elements in the input mesh");
00430 
00431         if ( _allSides.empty() || !_allSides.back()._grid.empty() )
00432           _allSides.push_back( _BlockSide() );
00433 
00434         _BlockSide& side = _allSides.back();
00435         if ( !fillSide( side, face, *corner ) )
00436         {
00437           if ( !_error.empty() )
00438             return false;
00439         }
00440         else
00441         {
00442           for ( int isXMax = 0; isXMax < 2; ++isXMax )
00443             for ( int isYMax = 0; isYMax < 2; ++isYMax )
00444             {
00445               const SMDS_MeshNode* nCorner = side.getCornerNode(isXMax,isYMax );
00446               corners.push_back( nCorner );
00447               cornerFaces.insert( side.getCornerFace( nCorner ));
00448             }
00449           for ( int e = 0; e < int(NB_QUAD_SIDES); ++e )
00450             _edge2sides[ side.getEdge( (EQuadEdge) e ) ].insert( &side );
00451 
00452           nbFacesOnSides += side.getNbFaces();
00453         }
00454       }
00455       ++corner;
00456 
00457       // find block sides of other domains if any
00458       if ( corner == corners.end() && nbFacesOnSides < mesh.NbQuadrangles() )
00459       {
00460         while ( nIt->more() )
00461         {
00462           nCorner = nIt->next();
00463           if ( isCornerNode( nCorner ))
00464             corner = corners.insert( corner, nCorner );
00465         }
00466         nbFacesOnSides = mesh.NbQuadrangles();
00467       }
00468     }
00469     
00470     if ( _allSides.empty() )
00471       return BAD_MESH_ERR;
00472     if ( _allSides.back()._grid.empty() )
00473       _allSides.pop_back();
00474     _DUMP_("Nb detected sides "<< _allSides.size());
00475 
00476     // ---------------------------
00477     // Organize sides into blocks
00478     // ---------------------------
00479 
00480     // analyse sharing of sides by blocks and sort sides by nb of adjacent sides
00481     int nbBlockSides = 0; // total nb of block sides taking into account their sharing
00482     multimap<int, _BlockSide* > sortedSides;
00483     {
00484       list < _BlockSide >::iterator sideIt = _allSides.begin();
00485       for ( ; sideIt != _allSides.end(); ++sideIt )
00486       {
00487         _BlockSide& side = *sideIt;
00488         bool isSharedSide = true;
00489         int nbAdjacent = 0;
00490         for ( int e = 0; e < int(NB_QUAD_SIDES) && isSharedSide; ++e )
00491         {
00492           int nbAdj = _edge2sides[ side.getEdge( (EQuadEdge) e ) ].size();
00493           nbAdjacent += nbAdj;
00494           isSharedSide = ( nbAdj > 2 );
00495         }
00496         side._nbBlocksFound = 0;
00497         side._nbBlocksExpected = isSharedSide ? 2 : 1;
00498         nbBlockSides += side._nbBlocksExpected;
00499         sortedSides.insert( make_pair( nbAdjacent, & side ));
00500       }
00501     }
00502 
00503     // find sides of each block
00504     int nbBlocks = 0;
00505     while ( nbBlockSides >= 6 )
00506     {
00507       // get any side not bound to all blocks it belongs to
00508       multimap<int, _BlockSide*>::iterator i_side = sortedSides.begin();
00509       while ( i_side != sortedSides.end() && i_side->second->isBound())
00510         ++i_side;
00511 
00512       // start searching for block sides from the got side
00513       bool ok = true;
00514       if ( _blocks.empty() || _blocks.back()._side[B_FRONT] )
00515         _blocks.resize( _blocks.size() + 1 );
00516 
00517       _Block& block = _blocks.back();
00518       block.setSide( B_FRONT, i_side->second );
00519       setSideBoundToBlock( *i_side->second );
00520       nbBlockSides--;
00521 
00522       // edges of adjacent sides of B_FRONT corresponding to front's edges
00523       EQuadEdge edgeOfFront[4] = { Q_BOTTOM, Q_RIGHT, Q_TOP, Q_LEFT };
00524       EQuadEdge edgeOfAdj  [4] = { Q_BOTTOM, Q_LEFT, Q_BOTTOM, Q_LEFT };
00525       // first find all sides detectable w/o advanced analysis,
00526       // then repeat the search, which then may pass without advanced analysis
00527       for ( int advAnalys = 0; advAnalys < 2; ++advAnalys )
00528       {
00529         for ( int i = 0; (ok || !advAnalys) && i < NB_QUAD_SIDES; ++i )
00530           if ( !block._side[i] ) // try to find 4 sides adjacent to front side
00531             ok = block.setSide( i, findBlockSide( B_FRONT, edgeOfFront[i],edgeOfAdj[i],advAnalys));
00532         if ( ok || !advAnalys)
00533           if ( !block._side[B_BACK] && block._side[B_TOP] ) // try to find back side by top one
00534             ok = block.setSide( B_BACK, findBlockSide( B_TOP, Q_TOP, Q_TOP, advAnalys ));
00535         if ( !advAnalys ) ok = true;
00536       }
00537       ok = block.isValid();
00538       if ( ok )
00539       {
00540         // check if just found block is same as one of previously found
00541         bool isSame = false;
00542         for ( int i = 1; i < _blocks.size() && !isSame; ++i )
00543           isSame = ( block._corners == _blocks[i-1]._corners );
00544         ok = !isSame;
00545       }
00546 
00547       // count the found sides
00548       _DUMP_(endl << "** Block " << _blocks.size() << " valid: " << block.isValid());
00549       for (int i = 0; i < NB_BLOCK_SIDES; ++i )
00550       {
00551         _DUMP_("\tSide "<< SBoxSides[i] <<" "<< block._side[ i ]._side);
00552         if ( block._side[ i ] )
00553         {
00554           if ( ok && i != B_FRONT)
00555           {
00556             setSideBoundToBlock( *block._side[ i ]._side );
00557             nbBlockSides--;
00558           }
00559           _DUMP_("\t corners "<<
00560                  block._side[ i ].cornerNode(0,0)->GetID() << ", " <<
00561                  block._side[ i ].cornerNode(1,0)->GetID() << ", " <<
00562                  block._side[ i ].cornerNode(1,1)->GetID() << ", " <<
00563                  block._side[ i ].cornerNode(0,1)->GetID() << ", "<<endl);
00564         }
00565         else
00566         {
00567           _DUMP_("\t not found"<<endl);
00568         }
00569       }
00570       if ( !ok )
00571         block.clear();
00572       else
00573         nbBlocks++;
00574     }
00575     _DUMP_("Nb found blocks "<< nbBlocks <<endl);
00576 
00577     if ( nbBlocks == 0 && _error.empty() )
00578       return BAD_MESH_ERR;
00579 
00580     return nbBlocks;
00581   }
00582 
00583   //================================================================================
00587   //================================================================================
00588 
00589   bool _Skin::fillSide( _BlockSide&             side,
00590                         const SMDS_MeshElement* cornerQuad,
00591                         const SMDS_MeshNode*    nCorner)
00592   {
00593     // Find out size of block side mesured in nodes and by the way find two rows
00594     // of nodes in two directions.
00595 
00596     int x, y, nbX, nbY;
00597     const SMDS_MeshElement* firstQuad = cornerQuad;
00598     {
00599       // get a node on block edge
00600       int iCorner = firstQuad->GetNodeIndex( nCorner );
00601       const SMDS_MeshNode* nOnEdge = firstQuad->GetNode( (iCorner+1) % 4);
00602 
00603       // find out size of block side
00604       vector<const SMDS_MeshNode*> horRow1, horRow2, verRow1, verRow2;
00605       if ( !fillRowsUntilCorner( firstQuad, nCorner, nOnEdge, horRow1, horRow2, true ) ||
00606            !fillRowsUntilCorner( firstQuad, nCorner, nOnEdge, verRow1, verRow2, false ))
00607         return false;
00608       nbX = horRow1.size(), nbY = verRow1.size();
00609 
00610       // store found nodes
00611       side._index._xSize = horRow1.size();
00612       side._index._ySize = verRow1.size();
00613       side._grid.resize( side._index.size(), NULL );
00614 
00615       for ( x = 0; x < horRow1.size(); ++x )
00616       {
00617         side.setNode( x, 0, horRow1[x] );
00618         side.setNode( x, 1, horRow2[x] );
00619       }
00620       for ( y = 0; y < verRow1.size(); ++y )
00621       {
00622         side.setNode( 0, y, verRow1[y] );
00623         side.setNode( 1, y, verRow2[y] );
00624       }
00625     }
00626     // Find the rest nodes
00627 
00628     y = 1; // y of the row to fill
00629     TIDSortedElemSet emptySet, avoidSet;
00630     while ( ++y < nbY )
00631     {
00632       // get next firstQuad in the next row of quadrangles
00633       //
00634       //          n2up
00635       //     o---o               <- y row
00636       //     |   |
00637       //     o---o  o  o  o  o   <- found nodes
00638       //n1down    n2down       
00639       //
00640       int i1down, i2down, i2up;
00641       const SMDS_MeshNode* n1down = side.getNode( 0, y-1 );
00642       const SMDS_MeshNode* n2down = side.getNode( 1, y-1 );
00643       avoidSet.clear(); avoidSet.insert( firstQuad );
00644       firstQuad = SMESH_MeshEditor::FindFaceInSet( n1down, n2down, emptySet, avoidSet,
00645                                                    &i1down, &i2down);
00646       if ( !isQuadrangle( firstQuad ))
00647         return BAD_MESH_ERR;
00648 
00649       const SMDS_MeshNode* n2up = oppositeNode( firstQuad, i1down );
00650       avoidSet.clear(); avoidSet.insert( firstQuad );
00651 
00652       // find the rest nodes in the y-th row by faces in the row
00653 
00654       x = 1; 
00655       while ( ++x < nbX )
00656       {
00657         const SMDS_MeshElement* quad = SMESH_MeshEditor::FindFaceInSet( n2up, n2down, emptySet,
00658                                                                         avoidSet, &i2up, &i2down);
00659         if ( !isQuadrangle( quad ))
00660           return BAD_MESH_ERR;
00661 
00662         n2up   = oppositeNode( quad, i2down );
00663         n2down = oppositeNode( quad, i2up );
00664         avoidSet.clear(); avoidSet.insert( quad );
00665 
00666         side.setNode( x, y, n2up );
00667       }
00668     }
00669 
00670     // check side validity
00671     bool ok =
00672       side.getCornerFace( side.getCornerNode( 0, 0 )) &&
00673       side.getCornerFace( side.getCornerNode( 1, 0 )) &&
00674       side.getCornerFace( side.getCornerNode( 0, 1 )) &&
00675       side.getCornerFace( side.getCornerNode( 1, 1 ));
00676 
00677     return ok;
00678   }
00679 
00680   //================================================================================
00684   //================================================================================
00685 
00686   _OrientedBlockSide _Skin::findBlockSide( EBoxSides startBlockSide,
00687                                            EQuadEdge sharedSideEdge1,
00688                                            EQuadEdge sharedSideEdge2,
00689                                            bool      withGeometricAnalysis)
00690   {
00691     _Block& block = _blocks.back();
00692     _OrientedBlockSide& side1 = block._side[ startBlockSide ];
00693 
00694     // get corner nodes of the given block edge
00695     SMESH_OrientedLink edge = side1.edge( sharedSideEdge1 );
00696     const SMDS_MeshNode* n1 = edge.node1();
00697     const SMDS_MeshNode* n2 = edge.node2();
00698     if ( edge._reversed ) swap( n1, n2 );
00699 
00700     // find all sides sharing both nodes n1 and n2
00701     set< _BlockSide* > sidesOnEdge = _edge2sides[ edge ]; // copy a set
00702 
00703     // exclude loaded sides of block from sidesOnEdge
00704     for (int i = 0; i < NB_BLOCK_SIDES; ++i )
00705       if ( block._side[ i ] )
00706         sidesOnEdge.erase( block._side[ i ]._side );
00707 
00708     int nbSidesOnEdge = sidesOnEdge.size();
00709     _DUMP_("nbSidesOnEdge "<< nbSidesOnEdge << " " << n1->GetID() << "-" << n2->GetID() );
00710     if ( nbSidesOnEdge == 0 )
00711       return 0;
00712 
00713     _BlockSide* foundSide = 0;
00714     if ( nbSidesOnEdge == 1 )
00715     {
00716       foundSide = *sidesOnEdge.begin();
00717     }
00718     else
00719     {
00720       set< _BlockSide* >::iterator sideIt = sidesOnEdge.begin();
00721       int nbLoadedSides = block.nbSides();
00722       if ( nbLoadedSides > 1 )
00723       {
00724         // Find the side having more than 2 corners common with already loaded sides
00725         for (; !foundSide && sideIt != sidesOnEdge.end(); ++sideIt )
00726         {
00727           _BlockSide* sideI = *sideIt;
00728           int nbCommonCorners =
00729             block._corners.count( sideI->getCornerNode(0,0)) +
00730             block._corners.count( sideI->getCornerNode(1,0)) +
00731             block._corners.count( sideI->getCornerNode(0,1)) +
00732             block._corners.count( sideI->getCornerNode(1,1));
00733           if ( nbCommonCorners > 2 )
00734             foundSide = sideI;
00735         }
00736       }
00737 
00738       if ( !foundSide )
00739       {
00740         if ( !withGeometricAnalysis ) return 0;
00741 
00742         // Select one of found sides most close to startBlockSide
00743 
00744         gp_XYZ p1 ( n1->X(),n1->Y(),n1->Z()),  p2 (n2->X(),n2->Y(),n2->Z());
00745         gp_Vec p1p2( p1, p2 );
00746 
00747         const SMDS_MeshElement* face1 = side1->getCornerFace( n1 );
00748         gp_XYZ p1Op = SMESH_TNodeXYZ( oppositeNode( face1, face1->GetNodeIndex(n1)));
00749         gp_Vec side1Dir( p1, p1Op );
00750         gp_Ax2 pln( p1, p1p2, side1Dir ); // plane with normal p1p2 and X dir side1Dir
00751         _DUMP_("  Select adjacent for "<< side1._side << " - side dir ("
00752                << side1Dir.X() << ", " << side1Dir.Y() << ", " << side1Dir.Z() << ")" );
00753 
00754         map < double , _BlockSide* > angleOfSide;
00755         for (sideIt = sidesOnEdge.begin(); sideIt != sidesOnEdge.end(); ++sideIt )
00756         {
00757           _BlockSide* sideI = *sideIt;
00758           const SMDS_MeshElement* faceI = sideI->getCornerFace( n1 );
00759           gp_XYZ p1Op = SMESH_TNodeXYZ( oppositeNode( faceI, faceI->GetNodeIndex(n1)));
00760           gp_Vec sideIDir( p1, p1Op );
00761           // compute angle of (sideIDir projection to pln) and (X dir of pln)
00762           gp_Vec2d sideIDirProj( sideIDir * pln.XDirection(), sideIDir * pln.YDirection());
00763           double angle = sideIDirProj.Angle( gp::DX2d() );
00764           if ( angle < 0 ) angle += 2 * PI; // angle [0-2*PI]
00765           angleOfSide.insert( make_pair( angle, sideI ));
00766           _DUMP_("  "<< sideI << " - side dir ("
00767                  << sideIDir.X() << ", " << sideIDir.Y() << ", " << sideIDir.Z() << ")"
00768                  << " angle " << angle);
00769         }
00770         if ( nbLoadedSides == 1 )
00771         {
00772           double angF = angleOfSide.begin()->first, angL = angleOfSide.rbegin()->first;
00773           if ( angF > PI ) angF = 2*PI - angF;
00774           if ( angL > PI ) angL = 2*PI - angL;
00775           foundSide = angF < angL ? angleOfSide.begin()->second : angleOfSide.rbegin()->second;
00776         }
00777         else
00778         {
00779           gp_XYZ gc(0,0,0); // gravity center of already loaded block sides
00780           for (int i = 0; i < NB_BLOCK_SIDES; ++i )
00781             if ( block._side[ i ] )
00782               gc += block._side[ i ]._side->getGC();
00783           gc /= nbLoadedSides;
00784 
00785           gp_Vec gcDir( p1, gc );
00786           gp_Vec2d gcDirProj( gcDir * pln.XDirection(), gcDir * pln.YDirection());
00787           double gcAngle = gcDirProj.Angle( gp::DX2d() );
00788           foundSide = gcAngle < 0 ? angleOfSide.rbegin()->second : angleOfSide.begin()->second;
00789         }
00790       }
00791       _DUMP_("  selected "<< foundSide );
00792     }
00793 
00794     // Orient the found side correctly
00795 
00796     // corners of found side corresponding to nodes n1 and n2
00797     bool xMax1, yMax1, xMax2, yMax2;
00798     if ( !getEdgeEnds( sharedSideEdge2, xMax1, yMax1, xMax2, yMax2 ))
00799       return error(SMESH_Comment("Internal error at ")<<__FILE__<<":"<<__LINE__),
00800         _OrientedBlockSide(0);
00801 
00802     for ( int ori = 0; ori < _OrientedIndexer::MAX_ORI+1; ++ori )
00803     {
00804       _OrientedBlockSide orientedSide( foundSide, ori );
00805       const SMDS_MeshNode* n12 = orientedSide.cornerNode( xMax1, yMax1);
00806       const SMDS_MeshNode* n22 = orientedSide.cornerNode( xMax2, yMax2);
00807       if ( n1 == n12 && n2 == n22 )
00808         return orientedSide;
00809     }
00810     error(SMESH_Comment("Failed to orient a block side found by edge ")<<sharedSideEdge1
00811           << " of side " << startBlockSide
00812           << " of block " << _blocks.size());
00813     return 0;
00814   }
00815 
00816   //================================================================================
00822   //================================================================================
00823 
00824   bool _Skin::fillRowsUntilCorner(const SMDS_MeshElement*       quad,
00825                                   const SMDS_MeshNode*          n1,
00826                                   const SMDS_MeshNode*          n2,
00827                                   vector<const SMDS_MeshNode*>& row1,
00828                                   vector<const SMDS_MeshNode*>& row2,
00829                                   const bool                    alongN1N2 )
00830   {
00831     const SMDS_MeshNode* corner1 = n1;
00832 
00833     // Store nodes of quad in the rows and find new n1 and n2 to get
00834     // the next face so that new n2 is on block edge
00835     int i1 = quad->GetNodeIndex( n1 );
00836     int i2 = quad->GetNodeIndex( n2 );
00837     row1.clear(); row2.clear();
00838     row1.push_back( n1 );
00839     if ( alongN1N2 )
00840     {
00841       row1.push_back( n2 );
00842       row2.push_back( oppositeNode( quad, i2 ));
00843       row2.push_back( n1 = oppositeNode( quad, i1 ));
00844     }
00845     else
00846     {
00847       row2.push_back( n2 );
00848       row1.push_back( n2 = oppositeNode( quad, i2 ));
00849       row2.push_back( n1 = oppositeNode( quad, i1 ));
00850     }
00851 
00852     _NodeDescriptor nDesc( row1[1] );
00853     if ( nDesc == _NodeDescriptor( row1[0] ))
00854       return true; // row of 2 nodes
00855 
00856     // Find the rest nodes
00857     TIDSortedElemSet emptySet, avoidSet;
00858     //while ( !isCornerNode( n2 ))
00859     while ( nDesc == _NodeDescriptor( n2 ))
00860     {
00861       avoidSet.clear(); avoidSet.insert( quad );
00862       quad = SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2 );
00863       if ( !isQuadrangle( quad ))
00864         return BAD_MESH_ERR;
00865 
00866       row1.push_back( n2 = oppositeNode( quad, i1 ));
00867       row2.push_back( n1 = oppositeNode( quad, i2 ));
00868     }
00869     return n1 != corner1;
00870   }
00871 
00872   //================================================================================
00876   //================================================================================
00877 
00878   const SMDS_MeshElement* _BlockSide::getCornerFace(const SMDS_MeshNode* cornerNode) const
00879   {
00880     int x, y, isXMax, isYMax, found = 0;
00881     for ( isXMax = 0; isXMax < 2; ++isXMax )
00882     {
00883       for ( isYMax = 0; isYMax < 2; ++isYMax )
00884       {
00885         x = isXMax ? _index._xSize-1 : 0;
00886         y = isYMax ? _index._ySize-1 : 0;
00887         found = ( getNode(x,y) == cornerNode );
00888         if ( found ) break;
00889       }
00890       if ( found ) break;
00891     }
00892     if ( !found ) return 0;
00893     int dx = isXMax ? -1 : +1;
00894     int dy = isYMax ? -1 : +1;
00895     const SMDS_MeshNode* n1 = getNode(x,y);
00896     const SMDS_MeshNode* n2 = getNode(x+dx,y);
00897     const SMDS_MeshNode* n3 = getNode(x,y+dy);
00898     const SMDS_MeshNode* n4 = getNode(x+dx,y+dy);
00899     return SMDS_Mesh::FindFace(n1, n2, n3, n4 );
00900   }
00901 
00902   //================================================================================
00906   //================================================================================
00907 
00908   bool _Block::isValid() const
00909   {
00910     bool ok = ( nbSides() == 6 );
00911 
00912     // check only corners depending on side selection
00913     EBoxSides adjacent[4] = { B_BOTTOM, B_RIGHT, B_TOP, B_LEFT };
00914     EQuadEdge edgeAdj [4] = { Q_TOP,    Q_RIGHT, Q_TOP, Q_RIGHT };
00915     EQuadEdge edgeBack[4] = { Q_BOTTOM, Q_RIGHT, Q_TOP, Q_LEFT };
00916 
00917     for ( int i=0; ok && i < NB_QUAD_SIDES; ++i )
00918     { 
00919       SMESH_OrientedLink eBack = _side[ B_BACK      ].edge( edgeBack[i] );
00920       SMESH_OrientedLink eAdja = _side[ adjacent[i] ].edge( edgeAdj[i] );
00921       ok = ( eBack == eAdja );
00922     }
00923     return ok;
00924   }
00925 
00926 } // namespace
00927 
00928 //=======================================================================
00929 //function : StdMeshers_HexaFromSkin_3D
00930 //purpose  : 
00931 //=======================================================================
00932 
00933 StdMeshers_HexaFromSkin_3D::StdMeshers_HexaFromSkin_3D(int hypId, int studyId, SMESH_Gen* gen)
00934   :SMESH_3D_Algo(hypId, studyId, gen)
00935 {
00936   MESSAGE("StdMeshers_HexaFromSkin_3D::StdMeshers_HexaFromSkin_3D");
00937   _name = "HexaFromSkin_3D";
00938 }
00939 
00940 StdMeshers_HexaFromSkin_3D::~StdMeshers_HexaFromSkin_3D()
00941 {
00942   MESSAGE("StdMeshers_HexaFromSkin_3D::~StdMeshers_HexaFromSkin_3D");
00943 }
00944 
00945 //================================================================================
00949 //================================================================================
00950 
00951 bool StdMeshers_HexaFromSkin_3D::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
00952 {
00953   _Skin skin;
00954   int nbBlocks = skin.findBlocks(aMesh);
00955   if ( nbBlocks == 0 )
00956     return error( skin.error());
00957 
00958   vector< vector< const SMDS_MeshNode* > > columns;
00959   int x, xSize, y, ySize, z, zSize;
00960   _Indexer colIndex;
00961 
00962   for ( int i = 0; i < nbBlocks; ++i )
00963   {
00964     const _Block& block = skin.getBlock( i );
00965 
00966     // ------------------------------------------
00967     // Fill columns of nodes with existing nodes
00968     // ------------------------------------------
00969 
00970     xSize = block.getSide(B_BOTTOM).getHoriSize();
00971     ySize = block.getSide(B_BOTTOM).getVertSize();
00972     zSize = block.getSide(B_FRONT ).getVertSize();
00973     int X = xSize - 1, Y = ySize - 1, Z = zSize - 1;
00974     colIndex = _Indexer( xSize, ySize );
00975     columns.resize( colIndex.size() );
00976 
00977     // fill node columns by front and back box sides
00978     for ( x = 0; x < xSize; ++x ) {
00979       vector< const SMDS_MeshNode* >& column0 = columns[ colIndex( x, 0 )];
00980       vector< const SMDS_MeshNode* >& column1 = columns[ colIndex( x, Y )];
00981       column0.resize( zSize );
00982       column1.resize( zSize );
00983       for ( z = 0; z < zSize; ++z ) {
00984         column0[ z ] = block.getSide(B_FRONT).node( x, z );
00985         column1[ z ] = block.getSide(B_BACK) .node( x, z );
00986       }
00987     }
00988     // fill node columns by left and right box sides
00989     for ( y = 1; y < ySize-1; ++y ) {
00990       vector< const SMDS_MeshNode* >& column0 = columns[ colIndex( 0, y )];
00991       vector< const SMDS_MeshNode* >& column1 = columns[ colIndex( X, y )];
00992       column0.resize( zSize );
00993       column1.resize( zSize );
00994       for ( z = 0; z < zSize; ++z ) {
00995         column0[ z ] = block.getSide(B_LEFT) .node( y, z );
00996         column1[ z ] = block.getSide(B_RIGHT).node( y, z );
00997       }
00998     }
00999     // get nodes from top and bottom box sides
01000     for ( x = 1; x < xSize-1; ++x ) {
01001       for ( y = 1; y < ySize-1; ++y ) {
01002         vector< const SMDS_MeshNode* >& column = columns[ colIndex( x, y )];
01003         column.resize( zSize );
01004         column.front() = block.getSide(B_BOTTOM).node( x, y );
01005         column.back()  = block.getSide(B_TOP)   .node( x, y );
01006       }
01007     }
01008 
01009     // ----------------------------
01010     // Add internal nodes of a box
01011     // ----------------------------
01012     // projection points of internal nodes on box subshapes by which
01013     // coordinates of internal nodes are computed
01014     vector<gp_XYZ> pointOnShape( SMESH_Block::ID_Shell );
01015 
01016     // projections on vertices are constant
01017     pointOnShape[ SMESH_Block::ID_V000 ] = block.getSide(B_BOTTOM).xyz( 0, 0 );
01018     pointOnShape[ SMESH_Block::ID_V100 ] = block.getSide(B_BOTTOM).xyz( X, 0 );
01019     pointOnShape[ SMESH_Block::ID_V010 ] = block.getSide(B_BOTTOM).xyz( 0, Y );
01020     pointOnShape[ SMESH_Block::ID_V110 ] = block.getSide(B_BOTTOM).xyz( X, Y );
01021     pointOnShape[ SMESH_Block::ID_V001 ] = block.getSide(B_TOP).xyz( 0, 0 );
01022     pointOnShape[ SMESH_Block::ID_V101 ] = block.getSide(B_TOP).xyz( X, 0 );
01023     pointOnShape[ SMESH_Block::ID_V011 ] = block.getSide(B_TOP).xyz( 0, Y );
01024     pointOnShape[ SMESH_Block::ID_V111 ] = block.getSide(B_TOP).xyz( X, Y );
01025 
01026     for ( x = 1; x < xSize-1; ++x )
01027     {
01028       gp_XYZ params; // normalized parameters of internal node within a unit box
01029       params.SetCoord( 1, x / double(X) );
01030       for ( y = 1; y < ySize-1; ++y )
01031       {
01032         params.SetCoord( 2, y / double(Y) );
01033         // column to fill during z loop
01034         vector< const SMDS_MeshNode* >& column = columns[ colIndex( x, y )];
01035         // projections on horizontal edges
01036         pointOnShape[ SMESH_Block::ID_Ex00 ] = block.getSide(B_BOTTOM).xyz( x, 0 );
01037         pointOnShape[ SMESH_Block::ID_Ex10 ] = block.getSide(B_BOTTOM).xyz( x, Y );
01038         pointOnShape[ SMESH_Block::ID_E0y0 ] = block.getSide(B_BOTTOM).xyz( 0, y );
01039         pointOnShape[ SMESH_Block::ID_E1y0 ] = block.getSide(B_BOTTOM).xyz( X, y );
01040         pointOnShape[ SMESH_Block::ID_Ex01 ] = block.getSide(B_TOP).xyz( x, 0 );
01041         pointOnShape[ SMESH_Block::ID_Ex11 ] = block.getSide(B_TOP).xyz( x, Y );
01042         pointOnShape[ SMESH_Block::ID_E0y1 ] = block.getSide(B_TOP).xyz( 0, y );
01043         pointOnShape[ SMESH_Block::ID_E1y1 ] = block.getSide(B_TOP).xyz( X, y );
01044         // projections on horizontal sides
01045         pointOnShape[ SMESH_Block::ID_Fxy0 ] = block.getSide(B_BOTTOM).xyz( x, y );
01046         pointOnShape[ SMESH_Block::ID_Fxy1 ] = block.getSide(B_TOP)   .xyz( x, y );
01047         for ( z = 1; z < zSize-1; ++z ) // z loop
01048         {
01049           params.SetCoord( 3, z / double(Z) );
01050           // projections on vertical edges
01051           pointOnShape[ SMESH_Block::ID_E00z ] = block.getSide(B_FRONT).xyz( 0, z );    
01052           pointOnShape[ SMESH_Block::ID_E10z ] = block.getSide(B_FRONT).xyz( X, z );    
01053           pointOnShape[ SMESH_Block::ID_E01z ] = block.getSide(B_BACK).xyz( 0, z );    
01054           pointOnShape[ SMESH_Block::ID_E11z ] = block.getSide(B_BACK).xyz( X, z );
01055           // projections on vertical sides
01056           pointOnShape[ SMESH_Block::ID_Fx0z ] = block.getSide(B_FRONT).xyz( x, z );    
01057           pointOnShape[ SMESH_Block::ID_Fx1z ] = block.getSide(B_BACK) .xyz( x, z );    
01058           pointOnShape[ SMESH_Block::ID_F0yz ] = block.getSide(B_LEFT) .xyz( y, z );    
01059           pointOnShape[ SMESH_Block::ID_F1yz ] = block.getSide(B_RIGHT).xyz( y, z );
01060 
01061           // compute internal node coordinates
01062           gp_XYZ coords;
01063           SMESH_Block::ShellPoint( params, pointOnShape, coords );
01064           column[ z ] = aHelper->AddNode( coords.X(), coords.Y(), coords.Z() );
01065 
01066 #ifdef DEB_GRID
01067           // debug
01068           //cout << "----------------------------------------------------------------------"<<endl;
01069           //for ( int id = SMESH_Block::ID_V000; id < SMESH_Block::ID_Shell; ++id)
01070           //{
01071           //  gp_XYZ p = pointOnShape[ id ];
01072           //  SMESH_Block::DumpShapeID( id,cout)<<" ( "<<p.X()<<", "<<p.Y()<<", "<<p.Z()<<" )"<<endl;
01073           //}
01074           //cout << "Params: ( "<< params.X()<<", "<<params.Y()<<", "<<params.Z()<<" )"<<endl;
01075           //cout << "coords: ( "<< coords.X()<<", "<<coords.Y()<<", "<<coords.Z()<<" )"<<endl;
01076 #endif
01077         }
01078       }
01079     }
01080     // ----------------
01081     // Add hexahedrons
01082     // ----------------
01083 
01084     // find out orientation by a least distorted hexahedron (issue 0020855);
01085     // the last is defined by evaluating sum of face normals of 8 corner hexahedrons
01086     double badness = numeric_limits<double>::max();
01087     bool isForw = true;
01088     for ( int xMax = 0; xMax < 2; ++xMax )
01089       for ( int yMax = 0; yMax < 2; ++yMax )
01090         for ( int zMax = 0; zMax < 2; ++zMax )
01091         {
01092           x = xMax ? xSize-1 : 1;
01093           y = yMax ? ySize-1 : 1;
01094           z = zMax ? zSize-1 : 1;
01095           vector< const SMDS_MeshNode* >& col00 = columns[ colIndex( x-1, y-1 )];
01096           vector< const SMDS_MeshNode* >& col10 = columns[ colIndex( x  , y-1 )];
01097           vector< const SMDS_MeshNode* >& col01 = columns[ colIndex( x-1, y   )];
01098           vector< const SMDS_MeshNode* >& col11 = columns[ colIndex( x  , y )];
01099           
01100           const SMDS_MeshNode* n000 = col00[z-1];
01101           const SMDS_MeshNode* n100 = col10[z-1];
01102           const SMDS_MeshNode* n010 = col01[z-1];
01103           const SMDS_MeshNode* n110 = col11[z-1];
01104           const SMDS_MeshNode* n001 = col00[z];
01105           const SMDS_MeshNode* n101 = col10[z];
01106           const SMDS_MeshNode* n011 = col01[z];
01107           const SMDS_MeshNode* n111 = col11[z];
01108           SMDS_VolumeOfNodes probeVolume (n000,n010,n110,n100,
01109                                           n001,n011,n111,n101);
01110           SMDS_VolumeTool volTool( &probeVolume );
01111           double Nx=0.,Ny=0.,Nz=0.;
01112           for ( int iFace = 0; iFace < volTool.NbFaces(); ++iFace )
01113           {
01114             double nx,ny,nz;
01115             volTool.GetFaceNormal( iFace, nx,ny,nz );
01116             Nx += nx;
01117             Ny += ny;
01118             Nz += nz;
01119           }
01120           double quality = Nx*Nx + Ny*Ny + Nz*Nz;
01121           if ( quality < badness )
01122           {
01123             badness = quality;
01124             isForw = volTool.IsForward();
01125           }
01126         }
01127 
01128     // add elements
01129     for ( x = 0; x < xSize-1; ++x ) {
01130       for ( y = 0; y < ySize-1; ++y ) {
01131         vector< const SMDS_MeshNode* >& col00 = columns[ colIndex( x, y )];
01132         vector< const SMDS_MeshNode* >& col10 = columns[ colIndex( x+1, y )];
01133         vector< const SMDS_MeshNode* >& col01 = columns[ colIndex( x, y+1 )];
01134         vector< const SMDS_MeshNode* >& col11 = columns[ colIndex( x+1, y+1 )];
01135         // bottom face normal of a hexa mush point outside the volume
01136         if ( isForw )
01137           for ( z = 0; z < zSize-1; ++z )
01138             aHelper->AddVolume(col00[z],   col01[z],   col11[z],   col10[z],
01139                                col00[z+1], col01[z+1], col11[z+1], col10[z+1]);
01140         else
01141           for ( z = 0; z < zSize-1; ++z )
01142             aHelper->AddVolume(col00[z],   col10[z],   col11[z],   col01[z],
01143                                col00[z+1], col10[z+1], col11[z+1], col01[z+1]);
01144       }
01145     }
01146   } // loop on blocks
01147 
01148   return true;
01149 }
01150 
01151 //================================================================================
01155 //================================================================================
01156 
01157 bool StdMeshers_HexaFromSkin_3D::Evaluate(SMESH_Mesh &         aMesh,
01158                                           const TopoDS_Shape & aShape,
01159                                           MapShapeNbElems&     aResMap)
01160 {
01161   _Skin skin;
01162   int nbBlocks = skin.findBlocks(aMesh);
01163   if ( nbBlocks == 0 )
01164     return error( skin.error());
01165 
01166   bool secondOrder = aMesh.NbFaces( ORDER_QUADRATIC );
01167 
01168   int entity = secondOrder ? SMDSEntity_Quad_Hexa : SMDSEntity_Hexa;
01169   vector<int>& nbByType = aResMap[ aMesh.GetSubMesh( aShape )];
01170   if ( entity >= nbByType.size() )
01171     nbByType.resize( SMDSEntity_Last, 0 );
01172 
01173   for ( int i = 0; i < nbBlocks; ++i )
01174   {
01175     const _Block& block = skin.getBlock( i );
01176 
01177     int nbX = block.getSide(B_BOTTOM).getHoriSize();
01178     int nbY = block.getSide(B_BOTTOM).getVertSize();
01179     int nbZ = block.getSide(B_FRONT ).getVertSize();
01180 
01181     int nbHexa  = (nbX-1) * (nbY-1) * (nbZ-1);
01182     int nbNodes = (nbX-2) * (nbY-2) * (nbZ-2);
01183     if ( secondOrder )
01184       nbNodes +=
01185         (nbX-2) * (nbY-2) * (nbZ-1) +
01186         (nbX-2) * (nbY-1) * (nbZ-2) +
01187         (nbX-1) * (nbY-2) * (nbZ-2);
01188 
01189 
01190     nbByType[ entity ] += nbHexa;
01191     nbByType[ SMDSEntity_Node ] += nbNodes;
01192   }
01193 
01194   return true;
01195 }
01196 
01197 //================================================================================
01201 //================================================================================
01202 
01203 bool StdMeshers_HexaFromSkin_3D::CheckHypothesis(SMESH_Mesh&, const TopoDS_Shape&,
01204                                                  Hypothesis_Status& aStatus)
01205 {
01206   aStatus = SMESH_Hypothesis::HYP_OK;
01207   return true;
01208 }
01209 
01210 //================================================================================
01215 //================================================================================
01216 
01217 bool StdMeshers_HexaFromSkin_3D::Compute(SMESH_Mesh&, const TopoDS_Shape&)
01218 {
01219   return error("Algorithm can't work with geometrical shapes");
01220 }
01221 
Copyright © 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
Copyright © 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS