1 // Copyright (C) 2007-2020 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // SMESH SMESH : implementation of SMESH idl descriptions
24 // File : StdMeshers_Hexa_3D.cxx
25 // Moved here from SMESH_Hexa_3D.cxx
26 // Author : Paul RASCLE, EDF
29 #include "StdMeshers_Hexa_3D.hxx"
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESH_Comment.hxx"
33 #include "SMESH_Gen.hxx"
34 #include "SMESH_Mesh.hxx"
35 #include "SMESH_MesherHelper.hxx"
36 #include "SMESH_subMesh.hxx"
37 #include "StdMeshers_BlockRenumber.hxx"
38 #include "StdMeshers_CompositeHexa_3D.hxx"
39 #include "StdMeshers_FaceSide.hxx"
40 #include "StdMeshers_HexaFromSkin_3D.hxx"
41 #include "StdMeshers_Penta_3D.hxx"
42 #include "StdMeshers_Prism_3D.hxx"
43 #include "StdMeshers_Quadrangle_2D.hxx"
44 #include "StdMeshers_ViscousLayers.hxx"
46 #include <BRep_Tool.hxx>
47 #include <Bnd_B3d.hxx>
49 #include <TopExp_Explorer.hxx>
50 #include <TopTools_MapOfShape.hxx>
51 #include <TopTools_SequenceOfShape.hxx>
54 #include "utilities.h"
55 #include "Utils_ExceptHandlers.hxx"
59 typedef SMESH_Comment TComm;
63 static SMESH_ComputeErrorPtr ComputePentahedralMesh(SMESH_Mesh &,
65 SMESH_ProxyMesh* proxyMesh=0);
67 static bool EvaluatePentahedralMesh(SMESH_Mesh &, const TopoDS_Shape &,
70 //=============================================================================
74 //=============================================================================
76 StdMeshers_Hexa_3D::StdMeshers_Hexa_3D(int hypId, SMESH_Gen * gen)
77 :SMESH_3D_Algo(hypId, gen)
80 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID); // 1 bit /shape type
81 _requireShape = false;
82 _compatibleHypothesis.push_back("ViscousLayers");
83 _compatibleHypothesis.push_back("BlockRenumber");
84 _quadAlgo = new StdMeshers_Quadrangle_2D( gen->GetANewId(), _gen );
87 //=============================================================================
91 //=============================================================================
93 StdMeshers_Hexa_3D::~StdMeshers_Hexa_3D()
99 //=============================================================================
101 * Retrieves defined hypotheses
103 //=============================================================================
105 bool StdMeshers_Hexa_3D::CheckHypothesis
107 const TopoDS_Shape& aShape,
108 SMESH_Hypothesis::Hypothesis_Status& aStatus)
110 // check nb of faces in the shape
112 aStatus = SMESH_Hypothesis::HYP_BAD_GEOMETRY;
114 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next())
121 _viscousLayersHyp = nullptr;
122 _blockRenumberHyp = nullptr;
124 const list<const SMESHDS_Hypothesis*>& hyps =
125 GetUsedHypothesis(aMesh, aShape, /*ignoreAuxiliary=*/false);
126 list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
127 if ( h == hyps.end())
129 aStatus = SMESH_Hypothesis::HYP_OK;
133 // only StdMeshers_ViscousLayers can be used
135 for ( ; h != hyps.end(); ++h )
137 if ( !_viscousLayersHyp &&
138 (_viscousLayersHyp = dynamic_cast< const StdMeshers_ViscousLayers*> ( *h )))
140 if ( !_blockRenumberHyp &&
141 (_blockRenumberHyp = dynamic_cast< const StdMeshers_BlockRenumber*> ( *h )))
145 if ((int) hyps.size() != (bool)_viscousLayersHyp + (bool)_blockRenumberHyp )
146 aStatus = HYP_INCOMPATIBLE;
149 if ( _viscousLayersHyp )
150 if ( !error( _viscousLayersHyp->CheckHypothesis( aMesh, aShape, aStatus )))
151 aStatus = HYP_BAD_PARAMETER;
153 if ( _blockRenumberHyp && aStatus == HYP_OK )
154 error( _blockRenumberHyp->CheckHypothesis( aMesh, aShape ));
157 return aStatus == HYP_OK;
162 //=============================================================================
164 typedef boost::shared_ptr< FaceQuadStruct > FaceQuadStructPtr;
165 typedef std::vector<gp_XYZ> TXYZColumn;
167 // symbolic names of box sides
168 enum EBoxSides{ B_BOTTOM=0, B_RIGHT, B_TOP, B_LEFT, B_FRONT, B_BACK, B_NB_SIDES };
170 // symbolic names of sides of quadrangle
171 enum EQuadSides{ Q_BOTTOM=0, Q_RIGHT, Q_TOP, Q_LEFT, Q_NB_SIDES };
173 enum EAxes{ COO_X=1, COO_Y, COO_Z };
175 //=============================================================================
177 * \brief Container of nodes of structured mesh on a qudrangular geom FACE
182 FaceQuadStructPtr _quad;
184 // map of (node parameter on EDGE) to (column (vector) of nodes)
185 TParam2ColumnMap _u2nodesMap;
187 // node column's taken from _u2nodesMap taking into account sub-shape orientation
188 vector<TNodeColumn> _columns;
190 // columns of normalized parameters of nodes within the unitary cube
191 vector<TXYZColumn> _ijkColumns;
193 // geometry of a cube side
196 const SMDS_MeshNode* GetNode(int iCol, int iRow) const
198 return _columns[iCol][iRow];
200 gp_XYZ GetXYZ(int iCol, int iRow) const
202 return SMESH_TNodeXYZ( GetNode( iCol, iRow ));
204 gp_XYZ& GetIJK(int iCol, int iRow)
206 return _ijkColumns[iCol][iRow];
210 //================================================================================
212 * \brief Converter of a pair of integers to a sole index
217 _Indexer( int xSize, int ySize ): _xSize(xSize), _ySize(ySize) {}
218 int size() const { return _xSize * _ySize; }
219 int operator()(const int x, const int y) const { return y * _xSize + x; }
222 //================================================================================
224 * \brief Finds FaceQuadStruct having a side equal to a given one and rearranges
225 * the found FaceQuadStruct::side to have the given side at a Q_BOTTOM place
227 FaceQuadStructPtr getQuadWithBottom( StdMeshers_FaceSidePtr side,
228 FaceQuadStructPtr quad[ 6 ])
230 FaceQuadStructPtr foundQuad;
231 for ( int i = 1; i < 6; ++i )
233 if ( !quad[i] ) continue;
234 for ( size_t iS = 0; iS < quad[i]->side.size(); ++iS )
236 const StdMeshers_FaceSidePtr side2 = quad[i]->side[iS];
237 if (( side->FirstVertex().IsSame( side2->FirstVertex() ) ||
238 side->FirstVertex().IsSame( side2->LastVertex() ))
240 ( side->LastVertex().IsSame( side2->FirstVertex() ) ||
241 side->LastVertex().IsSame( side2->LastVertex() ))
244 if ( iS != Q_BOTTOM )
246 vector< FaceQuadStruct::Side > newSides;
247 for ( size_t j = iS; j < quad[i]->side.size(); ++j )
248 newSides.push_back( quad[i]->side[j] );
249 for ( size_t j = 0; j < iS; ++j )
250 newSides.push_back( quad[i]->side[j] );
251 quad[i]->side.swap( newSides );
253 foundQuad.swap(quad[i]);
261 //================================================================================
263 * \brief Put quads to aCubeSide in the order of enum EBoxSides
265 //================================================================================
267 bool arrangeQuads( FaceQuadStructPtr quad[ 6 ], _FaceGrid aCubeSide[ 6 ], bool reverseBottom )
269 swap( aCubeSide[B_BOTTOM]._quad, quad[0] );
271 swap( aCubeSide[B_BOTTOM]._quad->side[ Q_RIGHT],// direct the bottom normal inside cube
272 aCubeSide[B_BOTTOM]._quad->side[ Q_LEFT ] );
274 aCubeSide[B_FRONT]._quad = getQuadWithBottom( aCubeSide[B_BOTTOM]._quad->side[Q_BOTTOM], quad );
275 aCubeSide[B_RIGHT]._quad = getQuadWithBottom( aCubeSide[B_BOTTOM]._quad->side[Q_RIGHT ], quad );
276 aCubeSide[B_BACK ]._quad = getQuadWithBottom( aCubeSide[B_BOTTOM]._quad->side[Q_TOP ], quad );
277 aCubeSide[B_LEFT ]._quad = getQuadWithBottom( aCubeSide[B_BOTTOM]._quad->side[Q_LEFT ], quad );
278 if ( aCubeSide[B_FRONT ]._quad )
279 aCubeSide[B_TOP]._quad = getQuadWithBottom( aCubeSide[B_FRONT ]._quad->side[Q_TOP ], quad );
281 for ( int i = 1; i < 6; ++i )
282 if ( !aCubeSide[i]._quad )
287 //================================================================================
289 * \brief Rearrange block sides according to StdMeshers_BlockRenumber hypothesis
291 //================================================================================
293 bool arrangeForRenumber( _FaceGrid blockSide[ 6 ],
295 TopoDS_Vertex& v001 )
297 std::swap( blockSide[B_BOTTOM]._quad->side[ Q_RIGHT],// restore after arrangeQuads()
298 blockSide[B_BOTTOM]._quad->side[ Q_LEFT ] );
301 TopTools_MapOfShape cornerVertices;
302 cornerVertices.Add( blockSide[B_BOTTOM]._quad->side[Q_BOTTOM].grid->LastVertex() );
303 cornerVertices.Add( blockSide[B_BOTTOM]._quad->side[Q_BOTTOM].grid->FirstVertex() );
304 cornerVertices.Add( blockSide[B_BOTTOM]._quad->side[Q_TOP ].grid->LastVertex() );
305 cornerVertices.Add( blockSide[B_BOTTOM]._quad->side[Q_TOP ].grid->FirstVertex() );
306 cornerVertices.Add( blockSide[B_TOP ]._quad->side[Q_BOTTOM].grid->FirstVertex() );
307 cornerVertices.Add( blockSide[B_TOP ]._quad->side[Q_BOTTOM].grid->LastVertex() );
308 cornerVertices.Add( blockSide[B_TOP ]._quad->side[Q_TOP ].grid->FirstVertex() );
309 cornerVertices.Add( blockSide[B_TOP ]._quad->side[Q_TOP ].grid->LastVertex() );
313 // block CS is not defined;
314 // renumber only if the block has an edge parallel to an axis of global CS
316 v000 = StdMeshers_RenumberHelper::GetVertex000( cornerVertices );
320 for ( auto it = cornerVertices.cbegin(); it != cornerVertices.cend(); ++it )
321 bbox.Add( BRep_Tool::Pnt( TopoDS::Vertex( *it )));
322 double tol = 1e-5 * Sqrt( bbox.SquareExtent() );
324 // get block edges starting at v000
326 std::vector< StdMeshers_FaceSidePtr > edgesAtV000;
327 std::vector< gp_Vec > edgeDir;
328 std::vector< int > iParallel; // 0 - none, 1 - X, 2 - Y, 3 - Z
329 TopTools_MapOfShape lastVertices;
330 for ( int iQ = 0; iQ < 6; ++iQ )
332 FaceQuadStructPtr quad = blockSide[iQ]._quad;
333 for ( size_t iS = 0; iS < quad->side.size() && edgesAtV000.size() < 3; ++iS )
335 StdMeshers_FaceSidePtr edge = quad->side[iS];
336 TopoDS_Vertex v1 = edge->FirstVertex(), v2 = edge->LastVertex();
337 if (( v1.IsSame( v000 ) && !lastVertices.Contains( v2 )) ||
338 ( v2.IsSame( v000 ) && !lastVertices.Contains( v1 )))
340 bool reverse = v2.IsSame( v000 );
343 lastVertices.Add( v2 );
345 edgesAtV000.push_back( edge );
347 gp_Pnt pf = BRep_Tool::Pnt( v1 );
348 gp_Pnt pl = BRep_Tool::Pnt( v2 );
349 gp_Vec vec( pf, pl );
350 edgeDir.push_back( vec );
352 iParallel.push_back( 0 );
353 if ( !v001.IsNull() )
355 if ( v001.IsSame( v2 ))
356 iParallel.back() = 3;
360 bool isStraight = true;
361 for ( int iE = 0; iE < edge->NbEdges() && isStraight; ++iE )
362 isStraight = SMESH_Algo::IsStraight( edge->Edge( iE ));
364 // is parallel to a GCS axis?
367 int nbDiff = (( Abs( vec.X() ) > tol ) +
368 ( Abs( vec.Y() ) > tol ) +
369 ( Abs( vec.Z() ) > tol ) );
371 iParallel.back() = ( Abs( vec.X() ) > tol ) ? 1 : ( Abs( vec.Y() ) > tol ) ? 2 : 3;
375 edgeDir.back() = gp_Vec( pf, edge->Value3d( reverse ? 0.99 : 0.01 ));
381 if ( std::accumulate( iParallel.begin(), iParallel.end(), 0 ) == 0 )
384 // find edge OZ and edge OX
385 StdMeshers_FaceSidePtr edgeOZ, edgeOX;
386 auto iZIt = std::find( iParallel.begin(), iParallel.end(), 3 );
387 if ( iZIt != iParallel.end() )
389 int i = std::distance( iParallel.begin(), iZIt );
390 edgeOZ = edgesAtV000[ i ];
391 int iE1 = SMESH_MesherHelper::WrapIndex( i + 1, edgesAtV000.size() );
392 int iE2 = SMESH_MesherHelper::WrapIndex( i + 2, edgesAtV000.size() );
393 if (( edgeDir[ iE1 ] ^ edgeDir[ iE2 ] ) * edgeDir[ i ] < 0 )
394 std::swap( iE1, iE2 );
395 edgeOX = edgesAtV000[ iE1 ];
399 for ( size_t i = 0; i < edgesAtV000.size(); ++i )
401 if ( !iParallel[ i ] )
403 int iE1 = SMESH_MesherHelper::WrapIndex( i + 1, edgesAtV000.size() );
404 int iE2 = SMESH_MesherHelper::WrapIndex( i + 2, edgesAtV000.size() );
405 if (( edgeDir[ iE1 ] ^ edgeDir[ iE2 ] ) * edgeDir[ i ] < 0 )
406 std::swap( iE1, iE2 );
407 edgeOZ = edgesAtV000[ iParallel[i] == 1 ? iE2 : iE1 ];
408 edgeOX = edgesAtV000[ iParallel[i] == 1 ? i : iE1 ];
413 if ( !edgeOZ || !edgeOX )
416 TopoDS_Vertex v100 = edgeOX->LastVertex();
417 if ( v100.IsSame( v000 ))
418 v100 = edgeOX->FirstVertex();
420 // Find the left quad, one including v000 but not v100
422 for ( int iQ = 0; iQ < 6; ++iQ )
424 FaceQuadStructPtr quad = blockSide[iQ]._quad;
425 bool hasV000 = false, hasV100 = false;
426 for ( size_t iS = 0; iS < quad->side.size(); ++iS )
428 StdMeshers_FaceSidePtr edge = quad->side[iS];
429 if ( edge->FirstVertex().IsSame( v000 ) || edge->LastVertex().IsSame( v000 ))
431 if ( edge->FirstVertex().IsSame( v100 ) || edge->LastVertex().IsSame( v100 ))
434 if ( hasV000 && !hasV100 )
436 // orient the left quad
437 for ( int i = 0; i < 4; ++i )
439 if ( quad->side[Q_BOTTOM].grid->Edge(0).IsSame( edgeOZ->Edge(0) ))
441 quad->shift( 1, true );
444 FaceQuadStructPtr quads[ 6 ];
445 quads[0].swap( blockSide[iQ]._quad );
446 for ( int i = 1, j = 0; i < 6; ++i, ++j )
447 if ( blockSide[ j ]._quad )
448 quads[ i ].swap( blockSide[ j ]._quad );
452 return arrangeQuads( quads, blockSide, false/* true*/ );
458 //================================================================================
460 * \brief Returns true if the 1st base node of sideGrid1 belongs to sideGrid2
462 //================================================================================
464 bool beginsAtSide( const _FaceGrid& sideGrid1,
465 const _FaceGrid& sideGrid2,
466 SMESH_ProxyMesh::Ptr proxymesh )
468 const TNodeColumn& col0 = sideGrid2._u2nodesMap.begin()->second;
469 const TNodeColumn& col1 = sideGrid2._u2nodesMap.rbegin()->second;
470 const SMDS_MeshNode* n00 = col0.front();
471 const SMDS_MeshNode* n01 = col0.back();
472 const SMDS_MeshNode* n10 = col1.front();
473 const SMDS_MeshNode* n11 = col1.back();
474 const SMDS_MeshNode* n = (sideGrid1._u2nodesMap.begin()->second)[0];
477 n00 = proxymesh->GetProxyNode( n00 );
478 n10 = proxymesh->GetProxyNode( n10 );
479 n01 = proxymesh->GetProxyNode( n01 );
480 n11 = proxymesh->GetProxyNode( n11 );
481 n = proxymesh->GetProxyNode( n );
483 return ( n == n00 || n == n01 || n == n10 || n == n11 );
486 //================================================================================
488 * \brief Fill in _FaceGrid::_ijkColumns
489 * \param [in,out] fg - a _FaceGrid
490 * \param [in] i1 - coordinate index along _columns
491 * \param [in] i2 - coordinate index along _columns[i]
492 * \param [in] v3 - value of the constant parameter
494 //================================================================================
496 void computeIJK( _FaceGrid& fg, int i1, int i2, double v3 )
498 gp_XYZ ijk( v3, v3, v3 );
499 const size_t nbCol = fg._columns.size();
500 const size_t nbRow = fg._columns[0].size();
502 fg._ijkColumns.resize( nbCol );
503 for ( size_t i = 0; i < nbCol; ++i )
504 fg._ijkColumns[ i ].resize( nbRow, ijk );
506 vector< double > len( nbRow );
508 for ( size_t i = 0; i < nbCol; ++i )
510 gp_Pnt pPrev = fg.GetXYZ( i, 0 );
511 for ( size_t j = 1; j < nbRow; ++j )
513 gp_Pnt p = fg.GetXYZ( i, j );
514 len[ j ] = len[ j-1 ] + p.Distance( pPrev );
517 for ( size_t j = 0; j < nbRow; ++j )
518 fg.GetIJK( i, j ).SetCoord( i2, len[ j ]/len.back() );
522 for ( size_t j = 0; j < nbRow; ++j )
524 gp_Pnt pPrev = fg.GetXYZ( 0, j );
525 for ( size_t i = 1; i < nbCol; ++i )
527 gp_Pnt p = fg.GetXYZ( i, j );
528 len[ i ] = len[ i-1 ] + p.Distance( pPrev );
531 for ( size_t i = 0; i < nbCol; ++i )
532 fg.GetIJK( i, j ).SetCoord( i1, len[ i ]/len.back() );
537 //=============================================================================
539 * Generates hexahedron mesh on hexaedron like form using algorithm from
540 * "Application de l'interpolation transfinie � la cr�ation de maillages
541 * C0 ou G1 continus sur des triangles, quadrangles, tetraedres, pentaedres
542 * et hexaedres d�form�s."
543 * Alain PERONNET - 8 janvier 1999
545 //=============================================================================
547 bool StdMeshers_Hexa_3D::Compute(SMESH_Mesh & aMesh,
548 const TopoDS_Shape & aShape)
550 // PAL14921. Enable catching std::bad_alloc and Standard_OutOfMemory outside
551 //Unexpect aCatch(SalomeException);
552 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
554 // Shape verification
555 // ----------------------
557 // shape must be a solid (or a shell) with 6 faces
558 TopExp_Explorer exp(aShape,TopAbs_SHELL);
560 return error(COMPERR_BAD_SHAPE, "No SHELL in the geometry");
561 if ( exp.Next(), exp.More() )
562 return error(COMPERR_BAD_SHAPE, "More than one SHELL in the geometry");
564 TopTools_IndexedMapOfShape FF, EE;
565 TopExp::MapShapes( aShape, TopAbs_FACE, FF);
566 if ( FF.Extent() != 6)
568 static StdMeshers_CompositeHexa_3D compositeHexa(_gen->GetANewId(), _gen);
569 compositeHexa.SetHypothesis( _blockRenumberHyp );
570 if ( !compositeHexa.Compute( aMesh, aShape ))
571 return error( compositeHexa.GetComputeError() );
573 return _blockRenumberHyp ? error( _blockRenumberHyp->CheckHypothesis( aMesh, aShape )) : true;
576 // Find sides of a cube
577 // ---------------------
579 // tool creating quadratic elements if needed
580 SMESH_MesherHelper helper (aMesh);
581 _quadraticMesh = helper.IsQuadraticSubMesh(aShape);
583 TopExp::MapShapes( aShape, TopAbs_EDGE, EE );
584 SMESH_MesherHelper* faceHelper = ( EE.Size() == 12 ) ? 0 : &helper;
586 FaceQuadStructPtr quad[ 6 ];
587 for ( int i = 0; i < 6; ++i )
590 faceHelper->SetSubShape( FF( i+1 ));
591 if ( !( quad[i] = FaceQuadStructPtr( _quadAlgo->CheckNbEdges( aMesh, FF( i+1 ),
592 /*considerMesh=*/true,
594 return error( _quadAlgo->GetComputeError() );
595 if ( quad[i]->side.size() != 4 )
596 return error( COMPERR_BAD_SHAPE, "Not a quadrangular box side" );
599 // put quads in a proper order
600 _FaceGrid aCubeSide[ 6 ];
601 if ( !arrangeQuads( quad, aCubeSide, true ))
602 return error( COMPERR_BAD_SHAPE );
605 // Make viscous layers
606 // --------------------
608 SMESH_ProxyMesh::Ptr proxymesh;
609 if ( _viscousLayersHyp )
611 proxymesh = _viscousLayersHyp->Compute( aMesh, aShape, /*makeN2NMap=*/ true );
616 // Check if there are triangles on cube sides
617 // -------------------------------------------
619 if ( aMesh.NbTriangles() > 0 )
621 for ( int i = 0; i < 6; ++i )
623 const TopoDS_Face& sideF = aCubeSide[i]._quad->face;
624 const SMESHDS_SubMesh* smDS =
625 proxymesh ? proxymesh->GetSubMesh( sideF ) : meshDS->MeshElements( sideF );
626 if ( !SMESH_MesherHelper::IsSameElemGeometry( smDS, SMDSGeom_QUADRANGLE,
627 /*nullSubMeshRes=*/false ))
629 SMESH_ComputeErrorPtr err = ComputePentahedralMesh(aMesh, aShape, proxymesh.get());
635 // Arrange sides according to _blockRenumberHyp
636 bool toRenumber = _blockRenumberHyp;
639 TopoDS_Vertex v000, v001;
640 _blockRenumberHyp->IsSolidIncluded( aMesh, aShape, v000, v001 );
642 toRenumber = arrangeForRenumber( aCubeSide, v000, v001 );
647 meshDS->CompactMesh(); // remove numbering holes
651 // Check presence of regular grid mesh on FACEs of the cube
652 // ------------------------------------------------------------
654 for ( int i = 0; i < 6; ++i )
656 const TopoDS_Face& F = aCubeSide[i]._quad->face;
657 StdMeshers_FaceSidePtr baseQuadSide = aCubeSide[i]._quad->side[ Q_BOTTOM ];
658 list<TopoDS_Edge> baseEdges( baseQuadSide->Edges().begin(), baseQuadSide->Edges().end() );
660 // assure correctness of node positions on baseE:
661 // helper.GetNodeU() will fix positions if they are wrong
662 helper.ToFixNodeParameters( true );
663 for ( int iE = 0; iE < baseQuadSide->NbEdges(); ++iE )
665 const TopoDS_Edge& baseE = baseQuadSide->Edge( iE );
666 if ( SMESHDS_SubMesh* smDS = meshDS->MeshElements( baseE ))
669 helper.SetSubShape( baseE );
670 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
671 while ( eIt->more() )
673 const SMDS_MeshElement* e = eIt->next();
674 // expect problems on a composite side
675 try { helper.GetNodeU( baseE, e->GetNode(0), e->GetNode(1), &ok); }
677 try { helper.GetNodeU( baseE, e->GetNode(1), e->GetNode(0), &ok); }
685 helper.LoadNodeColumns( aCubeSide[i]._u2nodesMap, F, baseEdges, meshDS, proxymesh.get());
688 // check if the loaded grid corresponds to nb of quadrangles on the FACE
689 const SMESHDS_SubMesh* faceSubMesh =
690 proxymesh ? proxymesh->GetSubMesh( F ) : meshDS->MeshElements( F );
691 const int nbQuads = faceSubMesh->NbElements();
692 const int nbHor = aCubeSide[i]._u2nodesMap.size() - 1;
693 const int nbVer = aCubeSide[i]._u2nodesMap.begin()->second.size() - 1;
694 ok = ( nbQuads == nbHor * nbVer );
698 SMESH_ComputeErrorPtr err = ComputePentahedralMesh(aMesh, aShape, proxymesh.get());
703 // Orient loaded grids of cube sides along axis of the unitary cube coord system
705 isReverse[B_BOTTOM] = beginsAtSide( aCubeSide[B_BOTTOM], aCubeSide[B_RIGHT ], proxymesh );
706 isReverse[B_TOP ] = beginsAtSide( aCubeSide[B_TOP ], aCubeSide[B_RIGHT ], proxymesh );
707 isReverse[B_FRONT ] = beginsAtSide( aCubeSide[B_FRONT ], aCubeSide[B_RIGHT ], proxymesh );
708 isReverse[B_BACK ] = beginsAtSide( aCubeSide[B_BACK ], aCubeSide[B_RIGHT ], proxymesh );
709 isReverse[B_LEFT ] = beginsAtSide( aCubeSide[B_LEFT ], aCubeSide[B_BACK ], proxymesh );
710 isReverse[B_RIGHT ] = beginsAtSide( aCubeSide[B_RIGHT ], aCubeSide[B_BACK ], proxymesh );
711 for ( int i = 0; i < 6; ++i )
713 aCubeSide[i]._columns.resize( aCubeSide[i]._u2nodesMap.size() );
715 size_t iFwd = 0, iRev = aCubeSide[i]._columns.size()-1;
716 size_t* pi = isReverse[i] ? &iRev : &iFwd;
717 TParam2ColumnMap::iterator u2nn = aCubeSide[i]._u2nodesMap.begin();
718 for ( ; iFwd < aCubeSide[i]._columns.size(); --iRev, ++iFwd, ++u2nn )
719 aCubeSide[i]._columns[ *pi ].swap( u2nn->second );
721 aCubeSide[i]._u2nodesMap.clear();
725 for ( int i = 0; i < 6; ++i )
726 for ( size_t j = 0; j < aCubeSide[i]._columns.size(); ++j)
727 for ( size_t k = 0; k < aCubeSide[i]._columns[j].size(); ++k)
729 const SMDS_MeshNode* & n = aCubeSide[i]._columns[j][k];
730 n = proxymesh->GetProxyNode( n );
733 // 4) Create internal nodes of the cube
734 // -------------------------------------
736 helper.SetSubShape( aShape );
737 helper.SetElementsOnShape(true);
739 // shortcuts to sides
740 _FaceGrid* fBottom = & aCubeSide[ B_BOTTOM ];
741 _FaceGrid* fRight = & aCubeSide[ B_RIGHT ];
742 _FaceGrid* fTop = & aCubeSide[ B_TOP ];
743 _FaceGrid* fLeft = & aCubeSide[ B_LEFT ];
744 _FaceGrid* fFront = & aCubeSide[ B_FRONT ];
745 _FaceGrid* fBack = & aCubeSide[ B_BACK ];
747 // cube size measured in nb of nodes
748 size_t x, xSize = fBottom->_columns.size() , X = xSize - 1;
749 size_t y, ySize = fLeft->_columns.size() , Y = ySize - 1;
750 size_t z, zSize = fLeft->_columns[0].size(), Z = zSize - 1;
752 // check sharing of FACEs (IPAL54417)
753 if ( fFront ->_columns.size() != xSize ||
754 fBack ->_columns.size() != xSize ||
755 fTop ->_columns.size() != xSize ||
757 fRight ->_columns.size() != ySize ||
758 fTop ->_columns[0].size() != ySize ||
759 fBottom->_columns[0].size() != ySize ||
761 fRight ->_columns[0].size() != zSize ||
762 fFront ->_columns[0].size() != zSize ||
763 fBack ->_columns[0].size() != zSize )
764 return error( COMPERR_BAD_SHAPE, "Not sewed faces" );
766 // columns of internal nodes "rising" from nodes of fBottom
767 _Indexer colIndex( xSize, ySize );
768 vector< vector< const SMDS_MeshNode* > > columns( colIndex.size() );
770 // fill node columns by front and back box sides
771 for ( x = 0; x < xSize; ++x ) {
772 vector< const SMDS_MeshNode* >& column0 = columns[ colIndex( x, 0 )];
773 vector< const SMDS_MeshNode* >& column1 = columns[ colIndex( x, Y )];
774 column0.resize( zSize );
775 column1.resize( zSize );
776 for ( z = 0; z < zSize; ++z ) {
777 column0[ z ] = fFront->GetNode( x, z );
778 column1[ z ] = fBack ->GetNode( x, z );
781 // fill node columns by left and right box sides
782 for ( y = 1; y < ySize-1; ++y ) {
783 vector< const SMDS_MeshNode* >& column0 = columns[ colIndex( 0, y )];
784 vector< const SMDS_MeshNode* >& column1 = columns[ colIndex( X, y )];
785 column0.resize( zSize );
786 column1.resize( zSize );
787 for ( z = 0; z < zSize; ++z ) {
788 column0[ z ] = fLeft ->GetNode( y, z );
789 column1[ z ] = fRight->GetNode( y, z );
792 // get nodes from top and bottom box sides
793 for ( x = 1; x < xSize-1; ++x ) {
794 for ( y = 1; y < ySize-1; ++y ) {
795 vector< const SMDS_MeshNode* >& column = columns[ colIndex( x, y )];
796 column.resize( zSize );
797 column.front() = fBottom->GetNode( x, y );
798 column.back() = fTop ->GetNode( x, y );
802 // compute normalized parameters of nodes on sides (PAL23189)
803 computeIJK( *fBottom, COO_X, COO_Y, /*z=*/0. );
804 computeIJK( *fRight, COO_Y, COO_Z, /*x=*/1. );
805 computeIJK( *fTop, COO_X, COO_Y, /*z=*/1. );
806 computeIJK( *fLeft, COO_Y, COO_Z, /*x=*/0. );
807 computeIJK( *fFront, COO_X, COO_Z, /*y=*/0. );
808 computeIJK( *fBack, COO_X, COO_Z, /*y=*/1. );
810 StdMeshers_RenumberHelper renumHelper( aMesh, _blockRenumberHyp );
812 // projection points of the internal node on cube sub-shapes by which
813 // coordinates of the internal node are computed
814 vector<gp_XYZ> pointsOnShapes( SMESH_Block::ID_Shell );
816 // projections on vertices are constant
817 pointsOnShapes[ SMESH_Block::ID_V000 ] = fBottom->GetXYZ( 0, 0 );
818 pointsOnShapes[ SMESH_Block::ID_V100 ] = fBottom->GetXYZ( X, 0 );
819 pointsOnShapes[ SMESH_Block::ID_V010 ] = fBottom->GetXYZ( 0, Y );
820 pointsOnShapes[ SMESH_Block::ID_V110 ] = fBottom->GetXYZ( X, Y );
821 pointsOnShapes[ SMESH_Block::ID_V001 ] = fTop->GetXYZ( 0, 0 );
822 pointsOnShapes[ SMESH_Block::ID_V101 ] = fTop->GetXYZ( X, 0 );
823 pointsOnShapes[ SMESH_Block::ID_V011 ] = fTop->GetXYZ( 0, Y );
824 pointsOnShapes[ SMESH_Block::ID_V111 ] = fTop->GetXYZ( X, Y );
826 gp_XYZ params; // normalized parameters of an internal node within the unit box
829 for ( y = 0; y < ySize; ++y )
831 vector< const SMDS_MeshNode* >& column0y = columns[ colIndex( 0, y )];
832 for ( z = 0; z < zSize; ++z )
833 renumHelper.AddReplacingNode( column0y[ z ] );
836 for ( x = 1; x < xSize-1; ++x )
840 vector< const SMDS_MeshNode* >& columnX0 = columns[ colIndex( x, 0 )];
841 for ( z = 0; z < zSize; ++z )
842 renumHelper.AddReplacingNode( columnX0[ z ] );
845 const double rX = x / double(X);
846 for ( y = 1; y < ySize-1; ++y )
848 const double rY = y / double(Y);
850 // a column to fill in during z loop
851 vector< const SMDS_MeshNode* >& column = columns[ colIndex( x, y )];
852 // projection points on horizontal edges
853 pointsOnShapes[ SMESH_Block::ID_Ex00 ] = fBottom->GetXYZ( x, 0 );
854 pointsOnShapes[ SMESH_Block::ID_Ex10 ] = fBottom->GetXYZ( x, Y );
855 pointsOnShapes[ SMESH_Block::ID_E0y0 ] = fBottom->GetXYZ( 0, y );
856 pointsOnShapes[ SMESH_Block::ID_E1y0 ] = fBottom->GetXYZ( X, y );
857 pointsOnShapes[ SMESH_Block::ID_Ex01 ] = fTop->GetXYZ( x, 0 );
858 pointsOnShapes[ SMESH_Block::ID_Ex11 ] = fTop->GetXYZ( x, Y );
859 pointsOnShapes[ SMESH_Block::ID_E0y1 ] = fTop->GetXYZ( 0, y );
860 pointsOnShapes[ SMESH_Block::ID_E1y1 ] = fTop->GetXYZ( X, y );
861 // projection points on horizontal faces
862 pointsOnShapes[ SMESH_Block::ID_Fxy0 ] = fBottom->GetXYZ( x, y );
863 pointsOnShapes[ SMESH_Block::ID_Fxy1 ] = fTop ->GetXYZ( x, y );
866 renumHelper.AddReplacingNode( column[ 0 ] );
868 for ( z = 1; z < zSize-1; ++z ) // z loop
870 const double rZ = z / double(Z);
872 const gp_XYZ& pBo = fBottom->GetIJK( x, y );
873 const gp_XYZ& pTo = fTop ->GetIJK( x, y );
874 const gp_XYZ& pFr = fFront ->GetIJK( x, z );
875 const gp_XYZ& pBa = fBack ->GetIJK( x, z );
876 const gp_XYZ& pLe = fLeft ->GetIJK( y, z );
877 const gp_XYZ& pRi = fRight ->GetIJK( y, z );
878 params.SetCoord( 1, 0.5 * ( pBo.X() * ( 1. - rZ ) + pTo.X() * rZ +
879 pFr.X() * ( 1. - rY ) + pBa.X() * rY ));
880 params.SetCoord( 2, 0.5 * ( pBo.Y() * ( 1. - rZ ) + pTo.Y() * rZ +
881 pLe.Y() * ( 1. - rX ) + pRi.Y() * rX ));
882 params.SetCoord( 3, 0.5 * ( pFr.Z() * ( 1. - rY ) + pBa.Z() * rY +
883 pLe.Z() * ( 1. - rX ) + pRi.Z() * rX ));
885 // projection points on vertical edges
886 pointsOnShapes[ SMESH_Block::ID_E00z ] = fFront->GetXYZ( 0, z );
887 pointsOnShapes[ SMESH_Block::ID_E10z ] = fFront->GetXYZ( X, z );
888 pointsOnShapes[ SMESH_Block::ID_E01z ] = fBack->GetXYZ( 0, z );
889 pointsOnShapes[ SMESH_Block::ID_E11z ] = fBack->GetXYZ( X, z );
890 // projection points on vertical faces
891 pointsOnShapes[ SMESH_Block::ID_Fx0z ] = fFront->GetXYZ( x, z );
892 pointsOnShapes[ SMESH_Block::ID_Fx1z ] = fBack ->GetXYZ( x, z );
893 pointsOnShapes[ SMESH_Block::ID_F0yz ] = fLeft ->GetXYZ( y, z );
894 pointsOnShapes[ SMESH_Block::ID_F1yz ] = fRight->GetXYZ( y, z );
896 // compute internal node coordinates
898 SMESH_Block::ShellPoint( params, pointsOnShapes, coords );
899 column[ z ] = helper.AddNode( coords.X(), coords.Y(), coords.Z() );
903 renumHelper.AddReplacingNode( column[ Z ] );
908 vector< const SMDS_MeshNode* >& columnX0 = columns[ colIndex( x, Y )];
909 for ( z = 0; z < zSize; ++z )
910 renumHelper.AddReplacingNode( columnX0[ z ] );
915 for ( y = 0; y < ySize; ++y )
917 vector< const SMDS_MeshNode* >& columnXy = columns[ colIndex( X, y )];
918 for ( z = 0; z < zSize; ++z )
919 renumHelper.AddReplacingNode( columnXy[ z ] );
922 // side data no more needed, free memory
923 for ( int i = 0; i < 6; ++i )
924 SMESHUtils::FreeVector( aCubeSide[i]._columns );
926 // 5) Create hexahedrons
927 // ---------------------
929 for ( x = 0; x < xSize-1; ++x ) {
930 for ( y = 0; y < ySize-1; ++y ) {
931 vector< const SMDS_MeshNode* >& col00 = columns[ colIndex( x, y )];
932 vector< const SMDS_MeshNode* >& col10 = columns[ colIndex( x+1, y )];
933 vector< const SMDS_MeshNode* >& col01 = columns[ colIndex( x, y+1 )];
934 vector< const SMDS_MeshNode* >& col11 = columns[ colIndex( x+1, y+1 )];
935 for ( z = 0; z < zSize-1; ++z )
937 // bottom face normal of a hexa mush point outside the volume
939 helper.AddVolume(col00[z], col01[z], col01[z+1], col00[z+1],
940 col10[z], col11[z], col11[z+1], col10[z+1]);
942 helper.AddVolume(col00[z], col01[z], col11[z], col10[z],
943 col00[z+1], col01[z+1], col11[z+1], col10[z+1]);
949 renumHelper.DoReplaceNodes();
952 if ( _blockRenumberHyp )
954 return error( _blockRenumberHyp->CheckHypothesis( aMesh, aShape ));
960 //=============================================================================
964 //=============================================================================
966 bool StdMeshers_Hexa_3D::Evaluate(SMESH_Mesh & aMesh,
967 const TopoDS_Shape & aShape,
968 MapShapeNbElems& aResMap)
970 vector < SMESH_subMesh * >meshFaces;
971 TopTools_SequenceOfShape aFaces;
972 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
973 aFaces.Append(exp.Current());
974 SMESH_subMesh *aSubMesh = aMesh.GetSubMeshContaining(exp.Current());
976 meshFaces.push_back(aSubMesh);
978 if (meshFaces.size() != 6) {
979 //return error(COMPERR_BAD_SHAPE, TComm(meshFaces.size())<<" instead of 6 faces in a block");
980 static StdMeshers_CompositeHexa_3D compositeHexa(-10, aMesh.GetGen());
981 return compositeHexa.Evaluate(aMesh, aShape, aResMap);
986 //TopoDS_Shape aFace = meshFaces[i]->GetSubShape();
987 TopoDS_Shape aFace = aFaces.Value(i+1);
988 SMESH_Algo *algo = _gen->GetAlgo(aMesh, aFace);
990 std::vector<int> aResVec(SMDSEntity_Last);
991 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
992 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
993 aResMap.insert(std::make_pair(sm,aResVec));
994 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
995 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
998 string algoName = algo->GetName();
999 bool isAllQuad = false;
1000 if (algoName == "Quadrangle_2D") {
1001 MapShapeNbElemsItr anIt = aResMap.find(meshFaces[i]);
1002 if( anIt == aResMap.end() ) continue;
1003 std::vector<int> aVec = (*anIt).second;
1004 int nbtri = Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
1008 if ( ! isAllQuad ) {
1009 return EvaluatePentahedralMesh(aMesh, aShape, aResMap);
1013 // find number of 1d elems for 1 face
1015 TopTools_MapOfShape Edges1;
1016 bool IsQuadratic = false;
1017 bool IsFirst = true;
1018 for (TopExp_Explorer exp(aFaces.Value(1), TopAbs_EDGE); exp.More(); exp.Next()) {
1019 Edges1.Add(exp.Current());
1020 SMESH_subMesh *sm = aMesh.GetSubMesh(exp.Current());
1022 MapShapeNbElemsItr anIt = aResMap.find(sm);
1023 if( anIt == aResMap.end() ) continue;
1024 std::vector<int> aVec = (*anIt).second;
1025 nb1d += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
1027 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1032 // find face opposite to 1 face
1034 for(i=2; i<=6; i++) {
1035 bool IsOpposite = true;
1036 for(TopExp_Explorer exp(aFaces.Value(i), TopAbs_EDGE); exp.More(); exp.Next()) {
1037 if( Edges1.Contains(exp.Current()) ) {
1047 // find number of 2d elems on side faces
1049 for(i=2; i<=6; i++) {
1050 if( i == OppNum ) continue;
1051 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[i-1] );
1052 if( anIt == aResMap.end() ) continue;
1053 std::vector<int> aVec = (*anIt).second;
1054 nb2d += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1057 MapShapeNbElemsItr anIt = aResMap.find( meshFaces[0] );
1058 std::vector<int> aVec = (*anIt).second;
1059 int nb2d_face0 = Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
1060 int nb0d_face0 = aVec[SMDSEntity_Node];
1062 std::vector<int> aResVec(SMDSEntity_Last);
1063 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
1065 aResVec[SMDSEntity_Quad_Hexa] = nb2d_face0 * ( nb2d/nb1d );
1066 int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
1067 aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
1070 aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
1071 aResVec[SMDSEntity_Hexa] = nb2d_face0 * ( nb2d/nb1d );
1073 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1074 aResMap.insert(std::make_pair(sm,aResVec));
1079 //================================================================================
1081 * \brief Computes hexahedral mesh from 2D mesh of block
1083 //================================================================================
1085 bool StdMeshers_Hexa_3D::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
1087 static StdMeshers_HexaFromSkin_3D * algo = 0;
1089 SMESH_Gen* gen = aMesh.GetGen();
1090 algo = new StdMeshers_HexaFromSkin_3D( gen->GetANewId(), gen );
1092 algo->InitComputeError();
1093 algo->Compute( aMesh, aHelper );
1094 return error( algo->GetComputeError());
1097 //================================================================================
1099 * \brief Return true if the algorithm can mesh this shape
1100 * \param [in] aShape - shape to check
1101 * \param [in] toCheckAll - if true, this check returns OK if all shapes are OK,
1102 * else, returns OK if at least one shape is OK
1104 //================================================================================
1106 bool StdMeshers_Hexa_3D::IsApplicable( const TopoDS_Shape & aShape, bool toCheckAll )
1108 TopExp_Explorer exp0( aShape, TopAbs_SOLID );
1109 if ( !exp0.More() ) return false;
1111 for ( ; exp0.More(); exp0.Next() )
1113 int nbFoundShells = 0;
1114 TopExp_Explorer exp1( exp0.Current(), TopAbs_SHELL );
1115 for ( ; exp1.More(); exp1.Next(), ++nbFoundShells)
1116 if ( nbFoundShells == 2 ) break;
1117 if ( nbFoundShells != 1 ) {
1118 if ( toCheckAll ) return false;
1121 exp1.Init( exp0.Current(), TopAbs_FACE );
1122 int nbEdges = SMESH_MesherHelper::Count( exp1.Current(), TopAbs_EDGE, /*ignoreSame=*/true );
1123 bool ok = ( nbEdges > 3 );
1124 if ( toCheckAll && !ok ) return false;
1125 if ( !toCheckAll && ok ) return true;
1130 //=======================================================================
1131 //function : ComputePentahedralMesh
1133 //=======================================================================
1135 SMESH_ComputeErrorPtr ComputePentahedralMesh(SMESH_Mesh & aMesh,
1136 const TopoDS_Shape & aShape,
1137 SMESH_ProxyMesh* proxyMesh)
1139 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New();
1142 err->myName = COMPERR_BAD_INPUT_MESH;
1143 err->myComment = "Can't build pentahedral mesh on viscous layers";
1147 StdMeshers_Penta_3D anAlgo;
1149 bOK=anAlgo.Compute(aMesh, aShape);
1151 err = anAlgo.GetComputeError();
1153 if ( !bOK && anAlgo.ErrorStatus() == 5 )
1155 static StdMeshers_Prism_3D * aPrism3D = 0;
1157 SMESH_Gen* gen = aMesh.GetGen();
1158 aPrism3D = new StdMeshers_Prism_3D( gen->GetANewId(), gen );
1160 SMESH_Hypothesis::Hypothesis_Status aStatus;
1161 if ( aPrism3D->CheckHypothesis( aMesh, aShape, aStatus ) ) {
1162 aPrism3D->InitComputeError();
1163 bOK = aPrism3D->Compute( aMesh, aShape );
1164 err = aPrism3D->GetComputeError();
1171 //=======================================================================
1172 //function : EvaluatePentahedralMesh
1174 //=======================================================================
1176 bool EvaluatePentahedralMesh(SMESH_Mesh & aMesh,
1177 const TopoDS_Shape & aShape,
1178 MapShapeNbElems& aResMap)
1180 StdMeshers_Penta_3D anAlgo;
1181 bool bOK = anAlgo.Evaluate(aMesh, aShape, aResMap);
1183 //err = anAlgo.GetComputeError();
1184 //if ( !bOK && anAlgo.ErrorStatus() == 5 )
1186 static StdMeshers_Prism_3D * aPrism3D = 0;
1188 SMESH_Gen* gen = aMesh.GetGen();
1189 aPrism3D = new StdMeshers_Prism_3D( gen->GetANewId(), gen );
1191 SMESH_Hypothesis::Hypothesis_Status aStatus;
1192 if ( aPrism3D->CheckHypothesis( aMesh, aShape, aStatus ) ) {
1193 return aPrism3D->Evaluate(aMesh, aShape, aResMap);