1 // Copyright (C) 2007-2021 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_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_VolumeTool.hxx"
35 #include "SMESHDS_Hypothesis.hxx"
36 #include "SMESHDS_Mesh.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;
62 //=======================================================================
63 //function : StdMeshers_Projection_3D
65 //=======================================================================
67 StdMeshers_Projection_3D::StdMeshers_Projection_3D(int hypId, SMESH_Gen* gen)
68 :SMESH_3D_Algo(hypId, gen)
70 _name = "Projection_3D";
71 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID); // 1 bit per shape type
73 _compatibleHypothesis.push_back("ProjectionSource3D");
77 //================================================================================
81 //================================================================================
83 StdMeshers_Projection_3D::~StdMeshers_Projection_3D()
86 //=======================================================================
87 //function : CheckHypothesis
89 //=======================================================================
91 bool StdMeshers_Projection_3D::CheckHypothesis(SMESH_Mesh& aMesh,
92 const TopoDS_Shape& aShape,
93 SMESH_Hypothesis::Hypothesis_Status& aStatus)
95 // check aShape that must be a 6 faces block
97 if ( TAssocTool::Count( aShape, TopAbs_SHELL, 1 ) != 1 ||
98 TAssocTool::Count( aShape, TopAbs_FACE , 1 ) != 6 ||
99 TAssocTool::Count( aShape, TopAbs_EDGE , 1 ) != 12 ||
100 TAssocTool::Count( aShape, TopAbs_WIRE , 1 ) != 6 )
102 aStatus = HYP_BAD_GEOMETRY;
106 list <const SMESHDS_Hypothesis * >::const_iterator itl;
108 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
109 if ( hyps.size() == 0 )
111 aStatus = SMESH_Hypothesis::HYP_MISSING;
112 return false; // can't work with no hypothesis
115 if ( hyps.size() > 1 )
117 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
121 const SMESHDS_Hypothesis *theHyp = hyps.front();
123 string hypName = theHyp->GetName();
125 aStatus = SMESH_Hypothesis::HYP_OK;
127 if (hypName == "ProjectionSource3D")
129 _sourceHypo = static_cast<const StdMeshers_ProjectionSource3D *>(theHyp);
130 // Check hypo parameters
132 SMESH_Mesh* srcMesh = _sourceHypo->GetSourceMesh();
133 SMESH_Mesh* tgtMesh = & aMesh;
138 if ( _sourceHypo->HasVertexAssociation() )
141 TopoDS_Shape edge = TAssocTool::GetEdgeByVertices
142 ( srcMesh, _sourceHypo->GetSourceVertex(1), _sourceHypo->GetSourceVertex(2) );
143 if ( edge.IsNull() ||
144 !SMESH_MesherHelper::IsSubShape( edge, srcMesh ) ||
145 !SMESH_MesherHelper::IsSubShape( edge, _sourceHypo->GetSource3DShape() ))
147 SCRUTE((edge.IsNull()));
148 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, srcMesh )));
149 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, _sourceHypo->GetSource3DShape() )));
150 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
155 edge = TAssocTool::GetEdgeByVertices
156 ( tgtMesh, _sourceHypo->GetTargetVertex(1), _sourceHypo->GetTargetVertex(2) );
157 if ( edge.IsNull() ||
158 !SMESH_MesherHelper::IsSubShape( edge, tgtMesh ) ||
159 !SMESH_MesherHelper::IsSubShape( edge, aShape ))
161 SCRUTE((edge.IsNull()));
162 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, tgtMesh )));
163 SCRUTE((SMESH_MesherHelper::IsSubShape( edge, aShape )));
164 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
168 // check a source shape
169 if ( !SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh ) ||
170 ( srcMesh == tgtMesh && aShape == _sourceHypo->GetSource3DShape()))
172 SCRUTE((SMESH_MesherHelper::IsSubShape( _sourceHypo->GetSource3DShape(), srcMesh)));
173 SCRUTE((srcMesh == tgtMesh));
174 SCRUTE((aShape == _sourceHypo->GetSource3DShape()));
175 aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
180 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
182 return ( aStatus == HYP_OK );
185 //=======================================================================
188 //=======================================================================
190 bool StdMeshers_Projection_3D::Compute(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
195 SMESH_Mesh * srcMesh = _sourceHypo->GetSourceMesh();
196 SMESH_Mesh * tgtMesh = & aMesh;
200 SMESHDS_Mesh * srcMeshDS = srcMesh->GetMeshDS();
201 SMESHDS_Mesh * tgtMeshDS = tgtMesh->GetMeshDS();
203 // get shell from shape3D
204 TopoDS_Shell srcShell, tgtShell;
205 TopExp_Explorer exp( _sourceHypo->GetSource3DShape(), TopAbs_SHELL );
207 for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
208 srcShell = TopoDS::Shell( exp.Current() );
210 return error(COMPERR_BAD_SHAPE,
211 SMESH_Comment("Source shape must have 1 shell but not ") << nbShell);
213 exp.Init( aShape, TopAbs_SHELL );
214 for ( nbShell = 0; exp.More(); exp.Next(), ++nbShell )
215 tgtShell = TopoDS::Shell( exp.Current() );
217 return error(COMPERR_BAD_SHAPE,
218 SMESH_Comment("Target shape must have 1 shell but not ") << nbShell);
220 // Check that shapes are blocks
221 if ( SMESH_MesherHelper::Count( tgtShell, TopAbs_FACE , 1 ) != 6 ||
222 SMESH_MesherHelper::Count( tgtShell, TopAbs_EDGE , 1 ) != 12 ||
223 SMESH_MesherHelper::Count( tgtShell, TopAbs_WIRE , 1 ) != 6 )
224 return error(COMPERR_BAD_SHAPE, "Target shape is not a block");
225 if ( SMESH_MesherHelper::Count( srcShell, TopAbs_FACE , 1 ) != 6 ||
226 SMESH_MesherHelper::Count( srcShell, TopAbs_EDGE , 1 ) != 12 ||
227 SMESH_MesherHelper::Count( srcShell, TopAbs_WIRE , 1 ) != 6 )
228 return error(COMPERR_BAD_SHAPE, "Source shape is not a block");
230 // Assure that mesh on a source shape is computed
232 SMESH_subMesh* srcSubMesh = srcMesh->GetSubMesh( _sourceHypo->GetSource3DShape() );
233 //SMESH_subMesh* tgtSubMesh = tgtMesh->GetSubMesh( aShape );
236 if ( tgtMesh == srcMesh && !aShape.IsSame( _sourceHypo->GetSource3DShape() )) {
237 if ( !TAssocTool::MakeComputed( srcSubMesh ))
238 srcMeshError = TAssocTool::SourceNotComputedError( srcSubMesh, this );
241 if ( !srcSubMesh->IsMeshComputed() )
242 srcMeshError = TAssocTool::SourceNotComputedError();
245 // Find 2 pairs of corresponding vertices
247 TopoDS_Vertex tgtV000, tgtV100, srcV000, srcV100;
248 TAssocTool::TShapeShapeMap shape2ShapeMap;
250 if ( _sourceHypo->HasVertexAssociation() )
252 tgtV000 = _sourceHypo->GetTargetVertex(1);
253 tgtV100 = _sourceHypo->GetTargetVertex(2);
254 srcV000 = _sourceHypo->GetSourceVertex(1);
255 srcV100 = _sourceHypo->GetSourceVertex(2);
259 if ( !TAssocTool::FindSubShapeAssociation( tgtShell, tgtMesh, srcShell, srcMesh,
261 return error(COMPERR_BAD_SHAPE,"Topology of source and target shapes seems different" );
263 exp.Init( tgtShell, TopAbs_EDGE );
264 TopExp::Vertices( TopoDS::Edge( exp.Current() ), tgtV000, tgtV100 );
266 if ( !shape2ShapeMap.IsBound( tgtV000 ) || !shape2ShapeMap.IsBound( tgtV100 ))
267 return error("Association of sub-shapes failed" );
268 srcV000 = TopoDS::Vertex( shape2ShapeMap( tgtV000 ));
269 srcV100 = TopoDS::Vertex( shape2ShapeMap( tgtV100 ));
270 if ( !SMESH_MesherHelper::IsSubShape( srcV000, srcShell ) ||
271 !SMESH_MesherHelper::IsSubShape( srcV100, srcShell ))
272 return error("Incorrect association of sub-shapes" );
275 // Load 2 SMESH_Block's with src and tgt shells
277 SMESH_Block srcBlock, tgtBlock;
278 TopTools_IndexedMapOfOrientedShape scrShapes, tgtShapes;
279 if ( !tgtBlock.LoadBlockShapes( tgtShell, tgtV000, tgtV100, tgtShapes ))
280 return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. Not a block?");
282 if ( !srcBlock.LoadBlockShapes( srcShell, srcV000, srcV100, scrShapes ))
283 return error(COMPERR_BAD_SHAPE, "Can't detect block sub-shapes. Not a block?");
285 // Find matching nodes of src and tgt shells
287 TNodeNodeMap src2tgtNodeMap;
288 for ( int fId = SMESH_Block::ID_FirstF; fId < SMESH_Block::ID_Shell; ++fId )
290 // Corresponding sub-shapes
291 TopoDS_Face srcFace = TopoDS::Face( scrShapes( fId ));
292 TopoDS_Face tgtFace = TopoDS::Face( tgtShapes( fId ));
293 if ( _sourceHypo->HasVertexAssociation() ) { // associate face sub-shapes
294 shape2ShapeMap.Clear();
295 vector< int > edgeIdVec;
296 SMESH_Block::GetFaceEdgesIDs( fId, edgeIdVec );
297 for ( size_t i = 0; i < edgeIdVec.size(); ++i ) {
298 int eID = edgeIdVec[ i ];
299 shape2ShapeMap.Bind( scrShapes( eID ), tgtShapes( eID ));
301 vector< int > vertexIdVec;
302 SMESH_Block::GetEdgeVertexIDs( eID, vertexIdVec );
303 shape2ShapeMap.Bind( scrShapes( vertexIdVec[0]), tgtShapes( vertexIdVec[0]) );
304 shape2ShapeMap.Bind( scrShapes( vertexIdVec[1]), tgtShapes( vertexIdVec[1]) );
308 // Find matching nodes of tgt and src faces
309 TAssocTool::TNodeNodeMap faceMatchingNodes;
310 if ( ! TAssocTool::FindMatchingNodesOnFaces( srcFace, srcMesh, tgtFace, tgtMesh,
311 shape2ShapeMap, faceMatchingNodes ))
312 return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Mesh on faces #")
313 << srcMeshDS->ShapeToIndex( srcFace ) << " and "
314 << tgtMeshDS->ShapeToIndex( tgtFace ) << " seems different" );
316 // put found matching nodes of 2 faces to the global map
317 src2tgtNodeMap.insert( faceMatchingNodes.begin(), faceMatchingNodes.end() );
320 // ------------------
322 // ------------------
324 SMDS_VolumeTool volTool;
325 SMESH_MesherHelper helper( *tgtMesh );
326 helper.IsQuadraticSubMesh( aShape );
328 SMESHDS_SubMesh* srcSMDS = srcSubMesh->GetSubMeshDS();
329 SMDS_ElemIteratorPtr volIt = srcSMDS->GetElements();
330 while ( volIt->more() ) // loop on source volumes
332 const SMDS_MeshElement* srcVol = volIt->next();
333 if ( !srcVol || srcVol->GetType() != SMDSAbs_Volume )
335 int nbNodes = srcVol->NbNodes();
336 SMDS_VolumeTool::VolumeType volType = volTool.GetType( nbNodes );
337 if ( srcVol->IsQuadratic() )
338 nbNodes = volTool.NbCornerNodes( volType );
340 // Find or create a new tgt node for each node of a src volume
342 vector< const SMDS_MeshNode* > nodes( nbNodes );
343 for ( int i = 0; i < nbNodes; ++i )
345 const SMDS_MeshNode* srcNode = srcVol->GetNode( i );
346 const SMDS_MeshNode* tgtNode = 0;
347 TNodeNodeMap::iterator sN_tN = src2tgtNodeMap.find( srcNode );
348 if ( sN_tN != src2tgtNodeMap.end() ) // found
350 tgtNode = sN_tN->second;
352 else // Create a new tgt node
354 // compute normalized parameters of source node in srcBlock
355 gp_Pnt srcCoord = gpXYZ( srcNode );
357 if ( !srcBlock.ComputeParameters( srcCoord, srcParam ))
358 return error(SMESH_Comment("Can't compute normalized parameters ")
359 << "for source node " << srcNode->GetID());
360 // compute coordinates of target node by srcParam
362 if ( !tgtBlock.ShellPoint( srcParam, tgtXYZ ))
363 return error("Can't compute coordinates by normalized parameters");
365 SMDS_MeshNode* newNode = tgtMeshDS->AddNode( tgtXYZ.X(), tgtXYZ.Y(), tgtXYZ.Z() );
366 tgtMeshDS->SetNodeInVolume( newNode, helper.GetSubShapeID() );
368 src2tgtNodeMap.insert( make_pair( srcNode, tgtNode ));
370 nodes[ i ] = tgtNode;
373 // Create a new volume
375 SMDS_MeshVolume * tgtVol = 0;
376 int id = 0, force3d = false;
378 case SMDS_VolumeTool::TETRA :
379 case SMDS_VolumeTool::QUAD_TETRA:
380 tgtVol = helper.AddVolume( nodes[0],
383 nodes[3], id, force3d); break;
384 case SMDS_VolumeTool::PYRAM :
385 case SMDS_VolumeTool::QUAD_PYRAM:
386 tgtVol = helper.AddVolume( nodes[0],
390 nodes[4], id, force3d); break;
391 case SMDS_VolumeTool::PENTA :
392 case SMDS_VolumeTool::QUAD_PENTA:
393 tgtVol = helper.AddVolume( nodes[0],
398 nodes[5], id, force3d); break;
399 case SMDS_VolumeTool::HEXA :
400 case SMDS_VolumeTool::QUAD_HEXA :
401 tgtVol = helper.AddVolume( nodes[0],
408 nodes[7], id, force3d); break;
409 default: // polyhedron
410 const SMDS_MeshVolume * poly = tgtMeshDS->DownCast<SMDS_MeshVolume>( 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<smIdType> 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;