1 // Copyright (C) 2007-2008 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
22 // SMESH SMESH : implementaion of SMESH idl descriptions
23 // File : StdMeshers_RadialPrism_3D.cxx
25 // Created : Fri Oct 20 11:37:07 2006
26 // Author : Edward AGAPOV (eap)
28 #include "StdMeshers_RadialPrism_3D.hxx"
30 #include "StdMeshers_ProjectionUtils.hxx"
31 #include "StdMeshers_NumberOfLayers.hxx"
32 #include "StdMeshers_LayerDistribution.hxx"
33 #include "StdMeshers_Prism_3D.hxx"
34 #include "StdMeshers_Regular_1D.hxx"
36 #include "SMDS_MeshNode.hxx"
37 #include "SMESHDS_SubMesh.hxx"
38 #include "SMESH_Gen.hxx"
39 #include "SMESH_Mesh.hxx"
40 #include "SMESH_MesherHelper.hxx"
41 #include "SMESH_subMesh.hxx"
43 #include "utilities.h"
45 #include <BRepAdaptor_Curve.hxx>
46 #include <BRepBuilderAPI_MakeEdge.hxx>
47 #include <BRepTools.hxx>
48 #include <BRep_Tool.hxx>
49 #include <TopExp_Explorer.hxx>
51 #include <TopoDS_Shell.hxx>
52 #include <TopoDS_Solid.hxx>
53 #include <TopTools_MapOfShape.hxx>
60 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
61 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
63 typedef StdMeshers_ProjectionUtils TAssocTool;
65 //=======================================================================
66 //function : StdMeshers_RadialPrism_3D
68 //=======================================================================
70 StdMeshers_RadialPrism_3D::StdMeshers_RadialPrism_3D(int hypId, int studyId, SMESH_Gen* gen)
71 :SMESH_3D_Algo(hypId, studyId, gen)
73 _name = "RadialPrism_3D";
74 _shapeType = (1 << TopAbs_SOLID); // 1 bit per shape type
76 _compatibleHypothesis.push_back("LayerDistribution");
77 _compatibleHypothesis.push_back("NumberOfLayers");
79 myDistributionHypo = 0;
82 //================================================================================
86 //================================================================================
88 StdMeshers_RadialPrism_3D::~StdMeshers_RadialPrism_3D()
91 //=======================================================================
92 //function : CheckHypothesis
94 //=======================================================================
96 bool StdMeshers_RadialPrism_3D::CheckHypothesis(SMESH_Mesh& aMesh,
97 const TopoDS_Shape& aShape,
98 SMESH_Hypothesis::Hypothesis_Status& aStatus)
100 // check aShape that must have 2 shells
102 if ( TAssocTool::Count( aShape, TopAbs_SOLID, 0 ) != 1 ||
103 TAssocTool::Count( aShape, TopAbs_SHELL, 0 ) != 2 )
105 aStatus = HYP_BAD_GEOMETRY;
110 myDistributionHypo = 0;
112 list <const SMESHDS_Hypothesis * >::const_iterator itl;
114 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
115 if ( hyps.size() == 0 )
117 aStatus = SMESH_Hypothesis::HYP_MISSING;
118 return false; // can't work with no hypothesis
121 if ( hyps.size() > 1 )
123 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
127 const SMESHDS_Hypothesis *theHyp = hyps.front();
129 string hypName = theHyp->GetName();
131 if (hypName == "NumberOfLayers")
133 myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
134 aStatus = SMESH_Hypothesis::HYP_OK;
137 if (hypName == "LayerDistribution")
139 myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
140 aStatus = SMESH_Hypothesis::HYP_OK;
143 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
147 //=======================================================================
150 //=======================================================================
152 bool StdMeshers_RadialPrism_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
155 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
157 myHelper = new SMESH_MesherHelper( aMesh );
158 myHelper->IsQuadraticSubMesh( aShape );
159 // to delete helper at exit from Compute()
160 std::auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
163 TopoDS_Solid solid = TopoDS::Solid( aShape );
164 TopoDS_Shell outerShell = BRepTools::OuterShell( solid );
165 TopoDS_Shape innerShell;
167 for ( TopoDS_Iterator It (solid); It.More(); It.Next(), ++nbShells )
168 if ( !outerShell.IsSame( It.Value() ))
169 innerShell = It.Value();
171 return error(COMPERR_BAD_SHAPE, SMESH_Comment("Must be 2 shells but not ")<<nbShells);
173 // ----------------------------------
174 // Associate subshapes of the shells
175 // ----------------------------------
177 TAssocTool::TShapeShapeMap shape2ShapeMap;
178 if ( !TAssocTool::FindSubShapeAssociation( outerShell, &aMesh,
181 return error(COMPERR_BAD_SHAPE,"Topology of inner and outer shells seems different" );
183 // ------------------
185 // ------------------
187 TNode2ColumnMap node2columnMap;
188 myLayerPositions.clear();
190 for ( exp.Init( outerShell, TopAbs_FACE ); exp.More(); exp.Next() )
192 // Corresponding subshapes
193 TopoDS_Face outFace = TopoDS::Face( exp.Current() );
195 if ( !shape2ShapeMap.IsBound( outFace )) {
196 return error(SMESH_Comment("Corresponding inner face not found for face #" )
197 << meshDS->ShapeToIndex( outFace ));
199 inFace = TopoDS::Face( shape2ShapeMap( outFace ));
202 // Find matching nodes of in and out faces
203 TNodeNodeMap nodeIn2OutMap;
204 if ( ! TAssocTool::FindMatchingNodesOnFaces( inFace, &aMesh, outFace, &aMesh,
205 shape2ShapeMap, nodeIn2OutMap ))
206 return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Mesh on faces #")
207 << meshDS->ShapeToIndex( outFace ) << " and "
208 << meshDS->ShapeToIndex( inFace ) << " seems different" );
212 SMDS_ElemIteratorPtr faceIt = meshDS->MeshElements( inFace )->GetElements();
213 while ( faceIt->more() ) // loop on faces on inFace
215 const SMDS_MeshElement* face = faceIt->next();
216 if ( !face || face->GetType() != SMDSAbs_Face )
218 int nbNodes = face->NbNodes();
219 if ( face->IsQuadratic() )
222 // find node columns for each node
223 vector< const TNodeColumn* > columns( nbNodes );
224 for ( int i = 0; i < nbNodes; ++i )
226 const SMDS_MeshNode* nIn = face->GetNode( i );
227 TNode2ColumnMap::iterator n_col = node2columnMap.find( nIn );
228 if ( n_col != node2columnMap.end() ) {
229 columns[ i ] = & n_col->second;
232 TNodeNodeMap::iterator nInOut = nodeIn2OutMap.find( nIn );
233 if ( nInOut == nodeIn2OutMap.end() )
234 RETURN_BAD_RESULT("No matching node for "<< nIn->GetID() <<
235 " in face "<< face->GetID());
236 columns[ i ] = makeNodeColumn( node2columnMap, nIn, nInOut->second );
240 StdMeshers_Prism_3D::AddPrisms( columns, myHelper );
242 } // loop on faces of out shell
247 //================================================================================
249 * \brief Create a column of nodes from outNode to inNode
250 * \param n2ColMap - map of node columns to add a created column
251 * \param outNode - botton node of a column
252 * \param inNode - top node of a column
253 * \retval const TNodeColumn* - a new column pointer
255 //================================================================================
257 TNodeColumn* StdMeshers_RadialPrism_3D::makeNodeColumn( TNode2ColumnMap& n2ColMap,
258 const SMDS_MeshNode* outNode,
259 const SMDS_MeshNode* inNode)
261 SMESHDS_Mesh * meshDS = myHelper->GetMeshDS();
262 int shapeID = myHelper->GetSubShapeID();
264 if ( myLayerPositions.empty() ) {
265 gp_Pnt pIn = gpXYZ( inNode ), pOut = gpXYZ( outNode );
266 computeLayerPositions( pIn, pOut );
268 int nbSegments = myLayerPositions.size() + 1;
270 TNode2ColumnMap::iterator n_col =
271 n2ColMap.insert( make_pair( outNode, TNodeColumn() )).first;
272 TNodeColumn & column = n_col->second;
273 column.resize( nbSegments + 1 );
274 column.front() = outNode;
275 column.back() = inNode;
277 gp_XYZ p1 = gpXYZ( outNode );
278 gp_XYZ p2 = gpXYZ( inNode );
279 for ( int z = 1; z < nbSegments; ++z )
281 double r = myLayerPositions[ z - 1 ];
282 gp_XYZ p = ( 1 - r ) * p1 + r * p2;
283 SMDS_MeshNode* n = meshDS->AddNode( p.X(), p.Y(), p.Z() );
284 meshDS->SetNodeInVolume( n, shapeID );
291 //================================================================================
292 //================================================================================
294 * \brief Class computing layers distribution using data of
295 * StdMeshers_LayerDistribution hypothesis
297 //================================================================================
298 //================================================================================
300 class TNodeDistributor: public StdMeshers_Regular_1D
302 list <const SMESHDS_Hypothesis *> myUsedHyps;
304 // -----------------------------------------------------------------------------
305 static TNodeDistributor* GetDistributor(SMESH_Mesh& aMesh)
307 const int myID = -1000;
308 map < int, SMESH_1D_Algo * > & algoMap = aMesh.GetGen()->_map1D_Algo;
309 map < int, SMESH_1D_Algo * >::iterator id_algo = algoMap.find( myID );
310 if ( id_algo == algoMap.end() )
311 return new TNodeDistributor( myID, 0, aMesh.GetGen() );
312 return static_cast< TNodeDistributor* >( id_algo->second );
314 // -----------------------------------------------------------------------------
315 bool Compute( vector< double > & positions,
319 const StdMeshers_LayerDistribution* hyp)
321 double len = pIn.Distance( pOut );
322 if ( len <= DBL_MIN ) return error("Too close points of inner and outer shells");
324 if ( !hyp || !hyp->GetLayerDistribution() )
325 return error( "Invalid LayerDistribution hypothesis");
327 myUsedHyps.push_back( hyp->GetLayerDistribution() );
329 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pIn, pOut );
330 SMESH_Hypothesis::Hypothesis_Status aStatus;
331 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
332 return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
333 "with LayerDistribution hypothesis");
335 BRepAdaptor_Curve C3D(edge);
336 double f = C3D.FirstParameter(), l = C3D.LastParameter();
337 list< double > params;
338 if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false ))
339 return error("StdMeshers_Regular_1D failed to compute layers distribution");
342 positions.reserve( params.size() );
343 for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
344 positions.push_back( *itU / len );
348 // -----------------------------------------------------------------------------
349 TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
350 : StdMeshers_Regular_1D( hypId, studyId, gen)
353 // -----------------------------------------------------------------------------
354 virtual const list <const SMESHDS_Hypothesis *> &
355 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
359 // -----------------------------------------------------------------------------
362 //================================================================================
364 * \brief Compute positions of nodes between the internal and the external surfaces
365 * \retval bool - is a success
367 //================================================================================
369 bool StdMeshers_RadialPrism_3D::computeLayerPositions(const gp_Pnt& pIn,
374 int nbSegments = myNbLayerHypo->GetNumberOfLayers();
375 myLayerPositions.resize( nbSegments - 1 );
376 for ( int z = 1; z < nbSegments; ++z )
377 myLayerPositions[ z - 1 ] = double( z )/ double( nbSegments );
380 if ( myDistributionHypo ) {
381 SMESH_Mesh * mesh = myHelper->GetMesh();
382 if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( myLayerPositions, pIn, pOut,
383 *mesh, myDistributionHypo ))
385 error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
389 RETURN_BAD_RESULT("Bad hypothesis");
393 //=======================================================================
394 //function : Evaluate
396 //=======================================================================
398 bool StdMeshers_RadialPrism_3D::Evaluate(SMESH_Mesh& aMesh,
399 const TopoDS_Shape& aShape,
400 MapShapeNbElems& aResMap)
403 TopoDS_Solid solid = TopoDS::Solid( aShape );
404 TopoDS_Shell outerShell = BRepTools::OuterShell( solid );
405 TopoDS_Shape innerShell;
407 for ( TopoDS_Iterator It (solid); It.More(); It.Next(), ++nbShells )
408 if ( !outerShell.IsSame( It.Value() ))
409 innerShell = It.Value();
410 if ( nbShells != 2 ) {
411 std::vector<int> aResVec(17);
412 for(int i=0; i<17; i++) aResVec[i] = 0;
413 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
414 aResMap.insert(std::make_pair(sm,aResVec));
415 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
416 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
420 // Associate subshapes of the shells
421 TAssocTool::TShapeShapeMap shape2ShapeMap;
422 if ( !TAssocTool::FindSubShapeAssociation( outerShell, &aMesh,
425 std::vector<int> aResVec(17);
426 for(int i=0; i<17; i++) aResVec[i] = 0;
427 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
428 aResMap.insert(std::make_pair(sm,aResVec));
429 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
430 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
434 // get info for outer shell
435 int nb0d_Out=0, nb2d_3_Out=0, nb2d_4_Out=0;
436 //TopTools_SequenceOfShape FacesOut;
437 for (TopExp_Explorer exp(outerShell, TopAbs_FACE); exp.More(); exp.Next()) {
438 //FacesOut.Append(exp.Current());
439 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
440 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
441 std::vector<int> aVec = (*anIt).second;
443 nb2d_3_Out += Max(aVec[3],aVec[4]);
444 nb2d_4_Out += Max(aVec[5],aVec[6]);
447 TopTools_MapOfShape tmpMap;
448 for (TopExp_Explorer exp(outerShell, TopAbs_EDGE); exp.More(); exp.Next()) {
449 if( tmpMap.Contains( exp.Current() ) )
451 tmpMap.Add( exp.Current() );
452 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
453 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
454 std::vector<int> aVec = (*anIt).second;
456 nb1d_Out += Max(aVec[1],aVec[2]);
459 for (TopExp_Explorer exp(outerShell, TopAbs_VERTEX); exp.More(); exp.Next()) {
460 if( tmpMap.Contains( exp.Current() ) )
462 tmpMap.Add( exp.Current() );
466 // get info for inner shell
467 int nb0d_In=0, nb2d_3_In=0, nb2d_4_In=0;
468 //TopTools_SequenceOfShape FacesIn;
469 for (TopExp_Explorer exp(innerShell, TopAbs_FACE); exp.More(); exp.Next()) {
470 //FacesIn.Append(exp.Current());
471 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
472 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
473 std::vector<int> aVec = (*anIt).second;
475 nb2d_3_In += Max(aVec[3],aVec[4]);
476 nb2d_4_In += Max(aVec[5],aVec[6]);
480 bool IsQuadratic = false;
482 for (TopExp_Explorer exp(innerShell, TopAbs_EDGE); exp.More(); exp.Next()) {
483 if( tmpMap.Contains( exp.Current() ) )
485 tmpMap.Add( exp.Current() );
486 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
487 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
488 std::vector<int> aVec = (*anIt).second;
490 nb1d_In += Max(aVec[1],aVec[2]);
492 IsQuadratic = (aVec[2] > aVec[1]);
497 for (TopExp_Explorer exp(innerShell, TopAbs_VERTEX); exp.More(); exp.Next()) {
498 if( tmpMap.Contains( exp.Current() ) )
500 tmpMap.Add( exp.Current() );
504 bool IsOK = (nb0d_Out==nb0d_In) && (nb1d_Out==nb1d_In) &&
505 (nb2d_3_Out==nb2d_3_In) && (nb2d_4_Out==nb2d_4_In);
507 std::vector<int> aResVec(17);
508 for(int i=0; i<17; i++) aResVec[i] = 0;
509 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
510 aResMap.insert(std::make_pair(sm,aResVec));
511 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
512 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
517 if( myNbLayerHypo ) {
518 nbLayers = myNbLayerHypo->GetNumberOfLayers();
520 if ( myDistributionHypo ) {
521 if ( !myDistributionHypo->GetLayerDistribution() ) {
522 std::vector<int> aResVec(17);
523 for(int i=0; i<17; 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));
530 TopExp_Explorer exp(outerShell, TopAbs_VERTEX);
531 TopoDS_Vertex Vout = TopoDS::Vertex(exp.Current());
532 TopoDS_Vertex Vin = TopoDS::Vertex( shape2ShapeMap(Vout) );
533 if ( myLayerPositions.empty() ) {
534 gp_Pnt pIn = BRep_Tool::Pnt(Vin);
535 gp_Pnt pOut = BRep_Tool::Pnt(Vout);
536 computeLayerPositions( pIn, pOut );
538 nbLayers = myLayerPositions.size() + 1;
541 std::vector<int> aResVec(17);
542 for(int i=0; i<17; i++) aResVec[i] = 0;
544 aResVec[13] = nb2d_3_Out * nbLayers;
545 aResVec[15] = nb2d_4_Out * nbLayers;
546 int nb1d = ( nb2d_3_Out*3 + nb2d_4_Out*4 ) / 2;
547 aResVec[0] = nb0d_Out * ( 2*nbLayers - 1 ) - nb1d * nbLayers;
550 aResVec[0] = nb0d_Out * ( nbLayers - 1 );
551 aResVec[12] = nb2d_3_Out * nbLayers;
552 aResVec[14] = nb2d_4_Out * nbLayers;
554 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
555 aResMap.insert(std::make_pair(sm,aResVec));