1 // Copyright (C) 2007-2012 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.
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 : implementaion of SMESH idl descriptions
24 // File : StdMeshers_RadialPrism_3D.cxx
26 // Created : Fri Oct 20 11:37:07 2006
27 // Author : Edward AGAPOV (eap)
30 #include "StdMeshers_RadialPrism_3D.hxx"
32 #include <Basics_OCCTVersion.hxx>
34 #include "StdMeshers_ProjectionUtils.hxx"
35 #include "StdMeshers_NumberOfLayers.hxx"
36 #include "StdMeshers_LayerDistribution.hxx"
37 #include "StdMeshers_Prism_3D.hxx"
38 #include "StdMeshers_Regular_1D.hxx"
40 #include "SMDS_MeshNode.hxx"
41 #include "SMESHDS_SubMesh.hxx"
42 #include "SMESH_Gen.hxx"
43 #include "SMESH_Mesh.hxx"
44 #include "SMESH_MesherHelper.hxx"
45 #include "SMESH_subMesh.hxx"
47 #include "utilities.h"
49 #include <BRepAdaptor_Curve.hxx>
50 #include <BRepBuilderAPI_MakeEdge.hxx>
51 #if OCC_VERSION_LARGE > 0x06050400
52 #include <BRepClass3d.hxx>
54 #include <BRepTools.hxx>
56 #include <BRep_Tool.hxx>
57 #include <TopExp_Explorer.hxx>
59 #include <TopoDS_Shell.hxx>
60 #include <TopoDS_Solid.hxx>
61 #include <TopTools_MapOfShape.hxx>
67 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
68 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
70 typedef StdMeshers_ProjectionUtils TAssocTool;
72 //=======================================================================
73 //function : StdMeshers_RadialPrism_3D
75 //=======================================================================
77 StdMeshers_RadialPrism_3D::StdMeshers_RadialPrism_3D(int hypId, int studyId, SMESH_Gen* gen)
78 :SMESH_3D_Algo(hypId, studyId, gen)
80 _name = "RadialPrism_3D";
81 _shapeType = (1 << TopAbs_SOLID); // 1 bit per shape type
83 _compatibleHypothesis.push_back("LayerDistribution");
84 _compatibleHypothesis.push_back("NumberOfLayers");
86 myDistributionHypo = 0;
89 //================================================================================
93 //================================================================================
95 StdMeshers_RadialPrism_3D::~StdMeshers_RadialPrism_3D()
98 //=======================================================================
99 //function : CheckHypothesis
101 //=======================================================================
103 bool StdMeshers_RadialPrism_3D::CheckHypothesis(SMESH_Mesh& aMesh,
104 const TopoDS_Shape& aShape,
105 SMESH_Hypothesis::Hypothesis_Status& aStatus)
107 // check aShape that must have 2 shells
109 if ( TAssocTool::Count( aShape, TopAbs_SOLID, 0 ) != 1 ||
110 TAssocTool::Count( aShape, TopAbs_SHELL, 0 ) != 2 )
112 aStatus = HYP_BAD_GEOMETRY;
117 myDistributionHypo = 0;
119 list <const SMESHDS_Hypothesis * >::const_iterator itl;
121 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
122 if ( hyps.size() == 0 )
124 aStatus = SMESH_Hypothesis::HYP_MISSING;
125 return false; // can't work with no hypothesis
128 if ( hyps.size() > 1 )
130 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
134 const SMESHDS_Hypothesis *theHyp = hyps.front();
136 string hypName = theHyp->GetName();
138 if (hypName == "NumberOfLayers")
140 myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
141 aStatus = SMESH_Hypothesis::HYP_OK;
144 if (hypName == "LayerDistribution")
146 myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
147 aStatus = SMESH_Hypothesis::HYP_OK;
150 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
154 //=======================================================================
157 //=======================================================================
159 bool StdMeshers_RadialPrism_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
162 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
164 myHelper = new SMESH_MesherHelper( aMesh );
165 myHelper->IsQuadraticSubMesh( aShape );
166 // to delete helper at exit from Compute()
167 std::auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
170 TopoDS_Solid solid = TopoDS::Solid( aShape );
171 #if OCC_VERSION_LARGE > 0x06050400
172 TopoDS_Shell outerShell = BRepClass3d::OuterShell( solid );
174 TopoDS_Shell outerShell = BRepTools::OuterShell( solid );
176 TopoDS_Shape innerShell;
178 for ( TopoDS_Iterator It (solid); It.More(); It.Next(), ++nbShells )
179 if ( !outerShell.IsSame( It.Value() ))
180 innerShell = It.Value();
182 return error(COMPERR_BAD_SHAPE, SMESH_Comment("Must be 2 shells but not ")<<nbShells);
184 // ----------------------------------
185 // Associate sub-shapes of the shells
186 // ----------------------------------
188 TAssocTool::TShapeShapeMap shape2ShapeMap;
189 if ( !TAssocTool::FindSubShapeAssociation( innerShell, &aMesh,
192 return error(COMPERR_BAD_SHAPE,"Topology of inner and outer shells seems different" );
194 // ------------------
196 // ------------------
198 TNode2ColumnMap node2columnMap;
199 myLayerPositions.clear();
201 for ( exp.Init( outerShell, TopAbs_FACE ); exp.More(); exp.Next() )
203 // Corresponding sub-shapes
204 TopoDS_Face outFace = TopoDS::Face( exp.Current() );
206 if ( !shape2ShapeMap.IsBound( outFace, /*isOut=*/true )) {
207 return error(SMESH_Comment("Corresponding inner face not found for face #" )
208 << meshDS->ShapeToIndex( outFace ));
210 inFace = TopoDS::Face( shape2ShapeMap( outFace, /*isOut=*/true ));
213 // Find matching nodes of in and out faces
214 TNodeNodeMap nodeIn2OutMap;
215 if ( ! TAssocTool::FindMatchingNodesOnFaces( inFace, &aMesh, outFace, &aMesh,
216 shape2ShapeMap, nodeIn2OutMap ))
217 return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Mesh on faces #")
218 << meshDS->ShapeToIndex( outFace ) << " and "
219 << meshDS->ShapeToIndex( inFace ) << " seems different" );
223 SMDS_ElemIteratorPtr faceIt = meshDS->MeshElements( inFace )->GetElements();
224 while ( faceIt->more() ) // loop on faces on inFace
226 const SMDS_MeshElement* face = faceIt->next();
227 if ( !face || face->GetType() != SMDSAbs_Face )
229 int nbNodes = face->NbNodes();
230 if ( face->IsQuadratic() )
233 // find node columns for each node
234 vector< const TNodeColumn* > columns( nbNodes );
235 for ( int i = 0; i < nbNodes; ++i )
237 const SMDS_MeshNode* nIn = face->GetNode( i );
238 TNode2ColumnMap::iterator n_col = node2columnMap.find( nIn );
239 if ( n_col != node2columnMap.end() ) {
240 columns[ i ] = & n_col->second;
243 TNodeNodeMap::iterator nInOut = nodeIn2OutMap.find( nIn );
244 if ( nInOut == nodeIn2OutMap.end() )
245 RETURN_BAD_RESULT("No matching node for "<< nIn->GetID() <<
246 " in face "<< face->GetID());
247 columns[ i ] = makeNodeColumn( node2columnMap, nIn, nInOut->second );
251 StdMeshers_Prism_3D::AddPrisms( columns, myHelper );
253 } // loop on faces of out shell
258 //================================================================================
260 * \brief Create a column of nodes from outNode to inNode
261 * \param n2ColMap - map of node columns to add a created column
262 * \param outNode - botton node of a column
263 * \param inNode - top node of a column
264 * \retval const TNodeColumn* - a new column pointer
266 //================================================================================
268 TNodeColumn* StdMeshers_RadialPrism_3D::makeNodeColumn( TNode2ColumnMap& n2ColMap,
269 const SMDS_MeshNode* outNode,
270 const SMDS_MeshNode* inNode)
272 SMESHDS_Mesh * meshDS = myHelper->GetMeshDS();
273 int shapeID = myHelper->GetSubShapeID();
275 if ( myLayerPositions.empty() ) {
276 gp_Pnt pIn = gpXYZ( inNode ), pOut = gpXYZ( outNode );
277 computeLayerPositions( pIn, pOut );
279 int nbSegments = myLayerPositions.size() + 1;
281 TNode2ColumnMap::iterator n_col =
282 n2ColMap.insert( make_pair( outNode, TNodeColumn() )).first;
283 TNodeColumn & column = n_col->second;
284 column.resize( nbSegments + 1 );
285 column.front() = outNode;
286 column.back() = inNode;
288 gp_XYZ p1 = gpXYZ( outNode );
289 gp_XYZ p2 = gpXYZ( inNode );
290 for ( int z = 1; z < nbSegments; ++z )
292 double r = myLayerPositions[ z - 1 ];
293 gp_XYZ p = ( 1 - r ) * p1 + r * p2;
294 SMDS_MeshNode* n = meshDS->AddNode( p.X(), p.Y(), p.Z() );
295 meshDS->SetNodeInVolume( n, shapeID );
302 //================================================================================
303 //================================================================================
305 * \brief Class computing layers distribution using data of
306 * StdMeshers_LayerDistribution hypothesis
308 //================================================================================
309 //================================================================================
311 class TNodeDistributor: public StdMeshers_Regular_1D
313 list <const SMESHDS_Hypothesis *> myUsedHyps;
315 // -----------------------------------------------------------------------------
316 static TNodeDistributor* GetDistributor(SMESH_Mesh& aMesh)
318 const int myID = -1000;
319 map < int, SMESH_1D_Algo * > & algoMap = aMesh.GetGen()->_map1D_Algo;
320 map < int, SMESH_1D_Algo * >::iterator id_algo = algoMap.find( myID );
321 if ( id_algo == algoMap.end() )
322 return new TNodeDistributor( myID, 0, aMesh.GetGen() );
323 return static_cast< TNodeDistributor* >( id_algo->second );
325 // -----------------------------------------------------------------------------
326 bool Compute( vector< double > & positions,
330 const StdMeshers_LayerDistribution* hyp)
332 double len = pIn.Distance( pOut );
333 if ( len <= DBL_MIN ) return error("Too close points of inner and outer shells");
335 if ( !hyp || !hyp->GetLayerDistribution() )
336 return error( "Invalid LayerDistribution hypothesis");
338 myUsedHyps.push_back( hyp->GetLayerDistribution() );
340 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pIn, pOut );
341 SMESH_Hypothesis::Hypothesis_Status aStatus;
342 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
343 return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
344 "with LayerDistribution hypothesis");
346 BRepAdaptor_Curve C3D(edge);
347 double f = C3D.FirstParameter(), l = C3D.LastParameter();
348 list< double > params;
349 if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false ))
350 return error("StdMeshers_Regular_1D failed to compute layers distribution");
353 positions.reserve( params.size() );
354 for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
355 positions.push_back( *itU / len );
359 // -----------------------------------------------------------------------------
360 TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
361 : StdMeshers_Regular_1D( hypId, studyId, gen)
364 // -----------------------------------------------------------------------------
365 virtual const list <const SMESHDS_Hypothesis *> &
366 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
370 // -----------------------------------------------------------------------------
373 //================================================================================
375 * \brief Compute positions of nodes between the internal and the external surfaces
376 * \retval bool - is a success
378 //================================================================================
380 bool StdMeshers_RadialPrism_3D::computeLayerPositions(const gp_Pnt& pIn,
385 int nbSegments = myNbLayerHypo->GetNumberOfLayers();
386 myLayerPositions.resize( nbSegments - 1 );
387 for ( int z = 1; z < nbSegments; ++z )
388 myLayerPositions[ z - 1 ] = double( z )/ double( nbSegments );
391 if ( myDistributionHypo ) {
392 SMESH_Mesh * mesh = myHelper->GetMesh();
393 if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( myLayerPositions, pIn, pOut,
394 *mesh, myDistributionHypo ))
396 error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
400 RETURN_BAD_RESULT("Bad hypothesis");
404 //=======================================================================
405 //function : Evaluate
407 //=======================================================================
409 bool StdMeshers_RadialPrism_3D::Evaluate(SMESH_Mesh& aMesh,
410 const TopoDS_Shape& aShape,
411 MapShapeNbElems& aResMap)
414 TopoDS_Solid solid = TopoDS::Solid( aShape );
415 #if OCC_VERSION_LARGE > 0x06050400
416 TopoDS_Shell outerShell = BRepClass3d::OuterShell( solid );
418 TopoDS_Shell outerShell = BRepTools::OuterShell( solid );
420 TopoDS_Shape innerShell;
422 for ( TopoDS_Iterator It (solid); It.More(); It.Next(), ++nbShells )
423 if ( !outerShell.IsSame( It.Value() ))
424 innerShell = It.Value();
425 if ( nbShells != 2 ) {
426 std::vector<int> aResVec(SMDSEntity_Last);
427 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
428 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
429 aResMap.insert(std::make_pair(sm,aResVec));
430 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
431 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
435 // Associate sub-shapes of the shells
436 TAssocTool::TShapeShapeMap shape2ShapeMap;
437 if ( !TAssocTool::FindSubShapeAssociation( outerShell, &aMesh,
440 std::vector<int> aResVec(SMDSEntity_Last);
441 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
442 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
443 aResMap.insert(std::make_pair(sm,aResVec));
444 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
445 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
449 // get info for outer shell
450 int nb0d_Out=0, nb2d_3_Out=0, nb2d_4_Out=0;
451 //TopTools_SequenceOfShape FacesOut;
452 for (TopExp_Explorer exp(outerShell, TopAbs_FACE); exp.More(); exp.Next()) {
453 //FacesOut.Append(exp.Current());
454 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
455 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
456 std::vector<int> aVec = (*anIt).second;
457 nb0d_Out += aVec[SMDSEntity_Node];
458 nb2d_3_Out += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
459 nb2d_4_Out += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
462 TopTools_MapOfShape tmpMap;
463 for (TopExp_Explorer exp(outerShell, TopAbs_EDGE); exp.More(); exp.Next()) {
464 if( tmpMap.Contains( exp.Current() ) )
466 tmpMap.Add( exp.Current() );
467 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
468 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
469 std::vector<int> aVec = (*anIt).second;
470 nb0d_Out += aVec[SMDSEntity_Node];
471 nb1d_Out += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
474 for (TopExp_Explorer exp(outerShell, TopAbs_VERTEX); exp.More(); exp.Next()) {
475 if( tmpMap.Contains( exp.Current() ) )
477 tmpMap.Add( exp.Current() );
481 // get info for inner shell
482 int nb0d_In=0, nb2d_3_In=0, nb2d_4_In=0;
483 //TopTools_SequenceOfShape FacesIn;
484 for (TopExp_Explorer exp(innerShell, TopAbs_FACE); exp.More(); exp.Next()) {
485 //FacesIn.Append(exp.Current());
486 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
487 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
488 std::vector<int> aVec = (*anIt).second;
489 nb0d_In += aVec[SMDSEntity_Node];
490 nb2d_3_In += Max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
491 nb2d_4_In += Max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
495 bool IsQuadratic = false;
497 for (TopExp_Explorer exp(innerShell, TopAbs_EDGE); exp.More(); exp.Next()) {
498 if( tmpMap.Contains( exp.Current() ) )
500 tmpMap.Add( exp.Current() );
501 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
502 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
503 std::vector<int> aVec = (*anIt).second;
504 nb0d_In += aVec[SMDSEntity_Node];
505 nb1d_In += Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
507 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
512 for (TopExp_Explorer exp(innerShell, TopAbs_VERTEX); exp.More(); exp.Next()) {
513 if( tmpMap.Contains( exp.Current() ) )
515 tmpMap.Add( exp.Current() );
519 bool IsOK = (nb0d_Out==nb0d_In) && (nb1d_Out==nb1d_In) &&
520 (nb2d_3_Out==nb2d_3_In) && (nb2d_4_Out==nb2d_4_In);
522 std::vector<int> aResVec(SMDSEntity_Last);
523 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
524 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
525 aResMap.insert(std::make_pair(sm,aResVec));
526 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
527 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
532 if( myNbLayerHypo ) {
533 nbLayers = myNbLayerHypo->GetNumberOfLayers();
535 if ( myDistributionHypo ) {
536 if ( !myDistributionHypo->GetLayerDistribution() ) {
537 std::vector<int> aResVec(SMDSEntity_Last);
538 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
539 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
540 aResMap.insert(std::make_pair(sm,aResVec));
541 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
542 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
545 TopExp_Explorer exp(outerShell, TopAbs_VERTEX);
546 TopoDS_Vertex Vout = TopoDS::Vertex(exp.Current());
547 TopoDS_Vertex Vin = TopoDS::Vertex( shape2ShapeMap(Vout) );
548 if ( myLayerPositions.empty() ) {
549 gp_Pnt pIn = BRep_Tool::Pnt(Vin);
550 gp_Pnt pOut = BRep_Tool::Pnt(Vout);
551 computeLayerPositions( pIn, pOut );
553 nbLayers = myLayerPositions.size() + 1;
556 std::vector<int> aResVec(SMDSEntity_Last);
557 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
559 aResVec[SMDSEntity_Quad_Penta] = nb2d_3_Out * nbLayers;
560 aResVec[SMDSEntity_Quad_Hexa] = nb2d_4_Out * nbLayers;
561 int nb1d = ( nb2d_3_Out*3 + nb2d_4_Out*4 ) / 2;
562 aResVec[SMDSEntity_Node] = nb0d_Out * ( 2*nbLayers - 1 ) - nb1d * nbLayers;
565 aResVec[SMDSEntity_Node] = nb0d_Out * ( nbLayers - 1 );
566 aResVec[SMDSEntity_Penta] = nb2d_3_Out * nbLayers;
567 aResVec[SMDSEntity_Hexa] = nb2d_4_Out * nbLayers;
569 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
570 aResMap.insert(std::make_pair(sm,aResVec));