1 // Copyright (C) 2007-2015 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 : implementaion of SMESH idl descriptions
24 // File : StdMeshers_Projection_3D.cxx
26 // Created : Fri Oct 20 11:37:07 2006
27 // Author : Edward AGAPOV (eap)
29 #include "StdMeshers_Projection_3D.hxx"
30 #include "StdMeshers_ProjectionSource3D.hxx"
32 #include "StdMeshers_ProjectionUtils.hxx"
34 #include "SMDS_PolyhedralVolumeOfNodes.hxx"
35 #include "SMDS_VolumeTool.hxx"
36 #include "SMESHDS_Hypothesis.hxx"
37 #include "SMESHDS_SubMesh.hxx"
38 #include "SMESH_Block.hxx"
39 #include "SMESH_Comment.hxx"
40 #include "SMESH_Gen.hxx"
41 #include "SMESH_Mesh.hxx"
42 #include "SMESH_MesherHelper.hxx"
43 #include "SMESH_Pattern.hxx"
44 #include "SMESH_subMesh.hxx"
45 #include "SMESH_subMeshEventListener.hxx"
47 #include "utilities.h"
50 #include <TopExp_Explorer.hxx>
53 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
54 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
55 #define SHOWYXZ(msg, xyz) \
56 //{gp_Pnt p(xyz); cout<<msg<< " ("<< p.X() << "; " <<p.Y() << "; " <<p.Z() << ") " <<endl; }
58 namespace TAssocTool = StdMeshers_ProjectionUtils;
61 //=======================================================================
62 //function : StdMeshers_Projection_3D
64 //=======================================================================
66 StdMeshers_Projection_3D::StdMeshers_Projection_3D(int hypId, int studyId, SMESH_Gen* gen)
67 :SMESH_3D_Algo(hypId, studyId, gen)
69 _name = "Projection_3D";
70 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID); // 1 bit per shape type
72 _compatibleHypothesis.push_back("ProjectionSource3D");
76 //================================================================================
80 //================================================================================
82 StdMeshers_Projection_3D::~StdMeshers_Projection_3D()
85 //=======================================================================
86 //function : CheckHypothesis
88 //=======================================================================
90 bool StdMeshers_Projection_3D::CheckHypothesis(SMESH_Mesh& aMesh,
91 const TopoDS_Shape& aShape,
92 SMESH_Hypothesis::Hypothesis_Status& aStatus)
94 // check aShape that must be a 6 faces block
96 if ( TAssocTool::Count( aShape, TopAbs_SHELL, 1 ) != 1 ||
97 TAssocTool::Count( aShape, TopAbs_FACE , 1 ) != 6 ||
98 TAssocTool::Count( aShape, TopAbs_EDGE , 1 ) != 12 ||
99 TAssocTool::Count( aShape, TopAbs_WIRE , 1 ) != 6 )
101 aStatus = HYP_BAD_GEOMETRY;
105 list <const SMESHDS_Hypothesis * >::const_iterator itl;
107 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
108 if ( hyps.size() == 0 )
110 aStatus = SMESH_Hypothesis::HYP_MISSING;
111 return false; // can't work with no hypothesis
114 if ( hyps.size() > 1 )
116 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
120 const SMESHDS_Hypothesis *theHyp = hyps.front();
122 string hypName = theHyp->GetName();
124 aStatus = SMESH_Hypothesis::HYP_OK;
126 if (hypName == "ProjectionSource3D")
128 _sourceHypo = static_cast<const StdMeshers_ProjectionSource3D *>(theHyp);
129 // Check hypo parameters
131 SMESH_Mesh* srcMesh = _sourceHypo->GetSourceMesh();
132 SMESH_Mesh* tgtMesh = & aMesh;
137 if ( _sourceHypo->HasVertexAssociation() )
140 TopoDS_Shape edge = TAssocTool::GetEdgeByVertices
141 ( srcMesh, _sourceHypo->GetSourceVertex(1), _sourceHypo->GetSourceVertex(2) );
142 if ( edge.IsNull() ||
143 !SMESH_MesherHelper::IsSubShape( edge, srcMesh ) ||
144 !SMESH_MesherHelper::IsSubShape( edge, _sourceHypo->GetSource3DShape() ))
146 SCRUTE((edge.IsNull()));
147 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, srcMesh )));
148 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, _sourceHypo->GetSource3DShape() )));
149 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
154 edge = TAssocTool::GetEdgeByVertices
155 ( tgtMesh, _sourceHypo->GetTargetVertex(1), _sourceHypo->GetTargetVertex(2) );
156 if ( edge.IsNull() ||
157 !SMESH_MesherHelper::IsSubShape( edge, tgtMesh ) ||
158 !SMESH_MesherHelper::IsSubShape( edge, aShape ))
160 SCRUTE((edge.IsNull()));
161 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, tgtMesh )));
162 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, aShape )));
163 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
167 // check a source shape
168 if ( !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh ) ||
169 ( srcMesh == tgtMesh && aShape == _sourceHypo->GetSource3DShape()))
171 SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh)));
172 SCRUTE((srcMesh == tgtMesh));
173 SCRUTE((aShape == _sourceHypo->GetSource3DShape()));
174 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
179 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
181 return ( aStatus == HYP_OK );
184 //=======================================================================
187 //=======================================================================
189 bool StdMeshers_Projection_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
194 SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
195 SMESH_Mesh * tgtMesh = & aMesh;
199 SMESHDS_Mesh * srcMeshDS = srcMesh->GetMeshDS();
200 SMESHDS_Mesh * tgtMeshDS = tgtMesh->GetMeshDS();
202 // get shell from shape3D
203 TopoDS_Shell srcShell, tgtShell;
204 TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
206 for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
207 srcShell = TopoDS::Shell( exp.Current() );
209 return error(COMPERR_BAD_SHAPE,
210 SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
212 exp.Init( aShape, TopAbs_SHELL );
213 for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
214 tgtShell = TopoDS::Shell( exp.Current() );
216 return error(COMPERR_BAD_SHAPE,
217 SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
219 // Check that shapes are blocks
220 if ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
221 SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
222 SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
223 return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
224 if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
225 SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
226 SMESH_MesherHelper::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
227 return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
229 // Assure that mesh on a source shape is computed
231 SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
232 //SMESH_subMesh* tgtSubMesh = tgtMesh->GetSubMesh( aShape );
235 if ( tgtMesh == srcMesh && !aShape.IsSame( _sourceHypo->GetSource3DShape() )) {
236 if ( !TAssocTool::MakeComputed( srcSubMesh ))
237 srcMeshError = TAssocTool::SourceNotComputedError( srcSubMesh, this );
240 if ( !srcSubMesh->IsMeshComputed() )
241 srcMeshError = TAssocTool::SourceNotComputedError();
244 // Find 2 pairs of corresponding vertices
246 TopoDS_Vertex tgtV000, tgtV100, srcV000, srcV100;
247 TAssocTool::TShapeShapeMap shape2ShapeMap;
249 if ( _sourceHypo->HasVertexAssociation() )
251 tgtV000 = _sourceHypo->GetTargetVertex(1);
252 tgtV100 = _sourceHypo->GetTargetVertex(2);
253 srcV000 = _sourceHypo->GetSourceVertex(1);
254 srcV100 = _sourceHypo->GetSourceVertex(2);
258 if ( !TAssocTool::FindSubShapeAssociation( tgtShell, tgtMesh, srcShell, srcMesh,
260 return error(COMPERR_BAD_SHAPE,"Topology of source and target shapes seems different" );
262 exp.Init( tgtShell, TopAbs_EDGE );
263 TopExp::Vertices( TopoDS::Edge( exp.Current() ), tgtV000, tgtV100 );
265 if ( !shape2ShapeMap.IsBound( tgtV000 ) || !shape2ShapeMap.IsBound( tgtV100 ))
266 return error("Association of sub-shapes failed" );
267 srcV000 = TopoDS::Vertex( shape2ShapeMap( tgtV000 ));
268 srcV100 = TopoDS::Vertex( shape2ShapeMap( tgtV100 ));
269 if ( !SMESH_MesherHelper::IsSubShape( srcV000, srcShell ) ||
270 !SMESH_MesherHelper::IsSubShape( srcV100, srcShell ))
271 return error("Incorrect association of sub-shapes" );
274 // Load 2 SMESH_Block's with src and tgt shells
276 SMESH_Block srcBlock, tgtBlock;
277 TopTools_IndexedMapOfOrientedShape scrShapes, tgtShapes;
278 if ( !tgtBlock.LoadBlockShapes( tgtShell, tgtV000, tgtV100, tgtShapes ))
279 return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. Not a block?");
281 if ( !srcBlock.LoadBlockShapes( srcShell, srcV000, srcV100, scrShapes ))
282 return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. Not a block?");
284 // Find matching nodes of src and tgt shells
286 TNodeNodeMap src2tgtNodeMap;
287 for ( int fId = SMESH_Block::ID_FirstF; fId < SMESH_Block::ID_Shell; ++fId )
289 // Corresponding sub-shapes
290 TopoDS_Face srcFace = TopoDS::Face( scrShapes( fId ));
291 TopoDS_Face tgtFace = TopoDS::Face( tgtShapes( fId ));
292 if ( _sourceHypo->HasVertexAssociation() ) { // associate face sub-shapes
293 shape2ShapeMap.Clear();
294 vector< int > edgeIdVec;
295 SMESH_Block::GetFaceEdgesIDs( fId, edgeIdVec );
296 for ( size_t i = 0; i < edgeIdVec.size(); ++i ) {
297 int eID = edgeIdVec[ i ];
298 shape2ShapeMap.Bind( scrShapes( eID ), tgtShapes( eID ));
300 vector< int > vertexIdVec;
301 SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
302 shape2ShapeMap.Bind( scrShapes( vertexIdVec[0]), tgtShapes( vertexIdVec[0]) );
303 shape2ShapeMap.Bind( scrShapes( vertexIdVec[1]), tgtShapes( vertexIdVec[1]) );
307 // Find matching nodes of tgt and src faces
308 TAssocTool::TNodeNodeMap faceMatchingNodes;
309 if ( ! TAssocTool::FindMatchingNodesOnFaces( srcFace, srcMesh, tgtFace, tgtMesh,
310 shape2ShapeMap, faceMatchingNodes ))
311 return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Mesh on faces #")
312 << srcMeshDS->ShapeToIndex( srcFace ) << " and "
313 << tgtMeshDS->ShapeToIndex( tgtFace ) << " seems different" );
315 // put found matching nodes of 2 faces to the global map
316 src2tgtNodeMap.insert( faceMatchingNodes.begin(), faceMatchingNodes.end() );
319 // ------------------
321 // ------------------
323 SMDS_VolumeTool volTool;
324 SMESH_MesherHelper helper( *tgtMesh );
325 helper.IsQuadraticSubMesh( aShape );
327 SMESHDS_SubMesh* srcSMDS = srcSubMesh->GetSubMeshDS();
328 SMDS_ElemIteratorPtr volIt = srcSMDS->GetElements();
329 while ( volIt->more() ) // loop on source volumes
331 const SMDS_MeshElement* srcVol = volIt->next();
332 if ( !srcVol || srcVol->GetType() != SMDSAbs_Volume )
334 int nbNodes = srcVol->NbNodes();
335 SMDS_VolumeTool::VolumeType volType = volTool.GetType( nbNodes );
336 if ( srcVol->IsQuadratic() )
337 nbNodes = volTool.NbCornerNodes( volType );
339 // Find or create a new tgt node for each node of a src volume
341 vector< const SMDS_MeshNode* > nodes( nbNodes );
342 for ( int i = 0; i < nbNodes; ++i )
344 const SMDS_MeshNode* srcNode = srcVol->GetNode( i );
345 const SMDS_MeshNode* tgtNode = 0;
346 TNodeNodeMap::iterator sN_tN = src2tgtNodeMap.find( srcNode );
347 if ( sN_tN != src2tgtNodeMap.end() ) // found
349 tgtNode = sN_tN->second;
351 else // Create a new tgt node
353 // compute normalized parameters of source node in srcBlock
354 gp_Pnt srcCoord = gpXYZ( srcNode );
356 if ( !srcBlock.ComputeParameters( srcCoord, srcParam ))
357 return error(SMESH_Comment("Can't compute normalized parameters ")
358 << "for source node " << srcNode->GetID());
359 // compute coordinates of target node by srcParam
361 if ( !tgtBlock.ShellPoint( srcParam, tgtXYZ ))
362 return error("Can't compute coordinates by normalized parameters");
364 SMDS_MeshNode* newNode = tgtMeshDS->AddNode( tgtXYZ.X(), tgtXYZ.Y(), tgtXYZ.Z() );
365 tgtMeshDS->SetNodeInVolume( newNode, helper.GetSubShapeID() );
367 src2tgtNodeMap.insert( make_pair( srcNode, tgtNode ));
369 nodes[ i ] = tgtNode;
372 // Create a new volume
374 SMDS_MeshVolume * tgtVol = 0;
375 int id = 0, force3d = false;
377 case SMDS_VolumeTool::TETRA :
378 case SMDS_VolumeTool::QUAD_TETRA:
379 tgtVol = helper.AddVolume( nodes[0],
382 nodes[3], id, force3d); break;
383 case SMDS_VolumeTool::PYRAM :
384 case SMDS_VolumeTool::QUAD_PYRAM:
385 tgtVol = helper.AddVolume( nodes[0],
389 nodes[4], id, force3d); break;
390 case SMDS_VolumeTool::PENTA :
391 case SMDS_VolumeTool::QUAD_PENTA:
392 tgtVol = helper.AddVolume( nodes[0],
397 nodes[5], id, force3d); break;
398 case SMDS_VolumeTool::HEXA :
399 case SMDS_VolumeTool::QUAD_HEXA :
400 tgtVol = helper.AddVolume( nodes[0],
407 nodes[7], id, force3d); break;
408 default: // polyhedron
409 const SMDS_VtkVolume * poly =
410 dynamic_cast<const SMDS_VtkVolume*>( srcVol );
412 RETURN_BAD_RESULT("Unexpected volume type");
413 if ( !poly->IsPoly())
414 RETURN_BAD_RESULT("Unexpected volume type");
415 tgtVol = tgtMeshDS->AddPolyhedralVolume( nodes, poly->GetQuantities() );
418 tgtMeshDS->SetMeshElementOnShape( tgtVol, helper.GetSubShapeID() );
420 } // loop on volumes of src shell
426 //=======================================================================
427 //function : Evaluate
429 //=======================================================================
431 bool StdMeshers_Projection_3D::Evaluate(SMESH_Mesh& aMesh,
432 const TopoDS_Shape& aShape,
433 MapShapeNbElems& aResMap)
438 SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
439 SMESH_Mesh * tgtMesh = & aMesh;
443 // get shell from shape3D
444 TopoDS_Shell srcShell, tgtShell;
445 TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
447 for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
448 srcShell = TopoDS::Shell( exp.Current() );
450 return error(COMPERR_BAD_SHAPE,
451 SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
453 exp.Init( aShape, TopAbs_SHELL );
454 for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
455 tgtShell = TopoDS::Shell( exp.Current() );
457 return error(COMPERR_BAD_SHAPE,
458 SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
460 // Check that shapes are blocks
461 if ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
462 SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
463 SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
464 return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
465 if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
466 SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
467 SMESH_MesherHelper::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
468 return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
470 // Assure that mesh on a source shape is computed
472 SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
474 if ( !srcSubMesh->IsMeshComputed() )
475 return error(COMPERR_BAD_INPUT_MESH,"Source mesh not computed");
478 std::vector<int> aVec(SMDSEntity_Last);
479 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
481 aVec[SMDSEntity_Node] = srcSubMesh->GetSubMeshDS()->NbNodes();
483 //bool quadratic = false;
484 SMDS_ElemIteratorPtr elemIt = srcSubMesh->GetSubMeshDS()->GetElements();
485 while ( elemIt->more() ) {
486 const SMDS_MeshElement* E = elemIt->next();
487 if( E->NbNodes()==4 ) {
488 aVec[SMDSEntity_Tetra]++;
490 else if( E->NbNodes()==5 ) {
491 aVec[SMDSEntity_Pyramid]++;
493 else if( E->NbNodes()==6 ) {
494 aVec[SMDSEntity_Penta]++;
496 else if( E->NbNodes()==8 ) {
497 aVec[SMDSEntity_Hexa]++;
499 else if( E->NbNodes()==10 && E->IsQuadratic() ) {
500 aVec[SMDSEntity_Quad_Tetra]++;
502 else if( E->NbNodes()==13 && E->IsQuadratic() ) {
503 aVec[SMDSEntity_Quad_Pyramid]++;
505 else if( E->NbNodes()==15 && E->IsQuadratic() ) {
506 aVec[SMDSEntity_Quad_Penta]++;
508 else if( E->NbNodes()==20 && E->IsQuadratic() ) {
509 aVec[SMDSEntity_Quad_Hexa]++;
512 aVec[SMDSEntity_Polyhedra]++;
516 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
517 aResMap.insert(std::make_pair(sm,aVec));
523 //=============================================================================
525 * \brief Sets a default event listener to submesh of the source shape
526 * \param subMesh - submesh where algo is set
528 * This method is called when a submesh gets HYP_OK algo_state.
529 * After being set, event listener is notified on each event of a submesh.
530 * Arranges that CLEAN event is translated from source submesh to
533 //=============================================================================
535 void StdMeshers_Projection_3D::SetEventListener(SMESH_subMesh* subMesh)
537 TAssocTool::SetEventListener( subMesh,
538 _sourceHypo->GetSource3DShape(),
539 _sourceHypo->GetSourceMesh() );
542 //================================================================================
544 * \brief Return true if the algorithm can mesh this shape
545 * \param [in] aShape - shape to check
546 * \param [in] toCheckAll - if true, this check returns OK if all shapes are OK,
547 * else, returns OK if at least one shape is OK
549 //================================================================================
551 bool StdMeshers_Projection_3D::IsApplicable(const TopoDS_Shape & aShape, bool toCheckAll)
553 TopExp_Explorer exp0( aShape, TopAbs_SOLID );
554 if ( !exp0.More() ) return false;
556 TopTools_IndexedMapOfOrientedShape blockShapes;
559 for ( ; exp0.More(); exp0.Next() )
561 int nbFoundShells = 0;
562 TopExp_Explorer exp1( exp0.Current(), TopAbs_SHELL );
563 for ( ; exp1.More(); exp1.Next(), ++nbFoundShells )
565 shell = TopoDS::Shell( exp1.Current() );
566 if ( nbFoundShells == 2 ) break;
568 if ( nbFoundShells != 1 ) {
569 if ( toCheckAll ) return false;
572 bool isBlock = SMESH_Block::FindBlockShapes( shell, v, v, blockShapes );
573 if ( toCheckAll && !isBlock ) return false;
574 if ( !toCheckAll && isBlock ) return true;