1 // Copyright (C) 2007-2016 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 : SMESH_Gen.cxx
25 // Author : Paul RASCLE, EDF
31 #include "SMESH_Gen.hxx"
33 #include "SMDS_Mesh.hxx"
34 #include "SMDS_MeshElement.hxx"
35 #include "SMDS_MeshNode.hxx"
36 #include "SMESHDS_Document.hxx"
37 #include "SMESH_HypoFilter.hxx"
38 #include "SMESH_Mesh.hxx"
39 #include "SMESH_MesherHelper.hxx"
40 #include "SMESH_subMesh.hxx"
42 #include "utilities.h"
44 #include "Utils_ExceptHandlers.hxx"
46 #include <TopoDS_Iterator.hxx>
57 //#include <vtkDebugLeaks.h>
60 //=============================================================================
64 //=============================================================================
66 SMESH_Gen::SMESH_Gen()
70 _segmentation = _nbSegments = 10;
71 SMDS_Mesh::_meshList.clear();
72 _compute_canceled = false;
73 //vtkDebugLeaks::SetExitError(0);
76 //=============================================================================
80 //=============================================================================
82 SMESH_Gen::~SMESH_Gen()
84 std::map < int, StudyContextStruct * >::iterator i_sc = _mapStudyContext.begin();
85 for ( ; i_sc != _mapStudyContext.end(); ++i_sc )
87 delete i_sc->second->myDocument;
92 //=============================================================================
94 * Creates a mesh in a study.
95 * if (theIsEmbeddedMode) { mesh modification commands are not logged }
97 //=============================================================================
99 SMESH_Mesh* SMESH_Gen::CreateMesh(int theStudyId, bool theIsEmbeddedMode)
100 throw(SALOME_Exception)
102 Unexpect aCatch(SalomeException);
104 // Get studyContext, create it if it does'nt exist, with a SMESHDS_Document
105 StudyContextStruct *aStudyContext = GetStudyContext(theStudyId);
107 // create a new SMESH_mesh object
108 SMESH_Mesh *aMesh = new SMESH_Mesh(_localId++,
112 aStudyContext->myDocument);
113 aStudyContext->mapMesh[_localId-1] = aMesh;
118 //=============================================================================
122 //=============================================================================
124 bool SMESH_Gen::Compute(SMESH_Mesh & aMesh,
125 const TopoDS_Shape & aShape,
126 const bool aShapeOnly /*=false*/,
127 const bool anUpward /*=false*/,
128 const ::MeshDimension aDim /*=::MeshDim_3D*/,
129 TSetOfInt* aShapesId /*=0*/)
135 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
137 const bool includeSelf = true;
138 const bool complexShapeFirst = true;
139 const int globalAlgoDim = 100;
141 SMESH_subMeshIteratorPtr smIt;
143 // Fix of Issue 22150. Due to !BLSURF->OnlyUnaryInput(), BLSURF computes edges
144 // that must be computed by Projection 1D-2D when Projection asks to compute
146 SMESH_subMesh::compute_event computeEvent =
147 aShapeOnly ? SMESH_subMesh::COMPUTE_SUBMESH : SMESH_subMesh::COMPUTE;
149 if ( anUpward ) // is called from the below code in this method
151 // ===============================================
152 // Mesh all the sub-shapes starting from vertices
153 // ===============================================
155 smIt = sm->getDependsOnIterator(includeSelf, !complexShapeFirst);
156 while ( smIt->more() )
158 SMESH_subMesh* smToCompute = smIt->next();
160 // do not mesh vertices of a pseudo shape
161 const TopoDS_Shape& shape = smToCompute->GetSubShape();
162 const TopAbs_ShapeEnum shapeType = shape.ShapeType();
163 if ( !aMesh.HasShapeToMesh() && shapeType == TopAbs_VERTEX )
166 // check for preview dimension limitations
167 if ( aShapesId && GetShapeDim( shapeType ) > (int)aDim )
169 // clear compute state not to show previous compute errors
170 // if preview invoked less dimension less than previous
171 smToCompute->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
175 if (smToCompute->GetComputeState() == SMESH_subMesh::READY_TO_COMPUTE)
177 if (_compute_canceled)
179 setCurrentSubMesh( smToCompute );
180 smToCompute->ComputeStateEngine( computeEvent );
181 setCurrentSubMesh( NULL );
184 // we check all the sub-meshes here and detect if any of them failed to compute
185 if (smToCompute->GetComputeState() == SMESH_subMesh::FAILED_TO_COMPUTE &&
186 ( shapeType != TopAbs_EDGE || !SMESH_Algo::isDegenerated( TopoDS::Edge( shape ))))
188 else if ( aShapesId )
189 aShapesId->insert( smToCompute->GetId() );
191 //aMesh.GetMeshDS()->Modified();
196 // ================================================================
197 // Apply algos that do NOT require discreteized boundaries
198 // ("all-dimensional") and do NOT support sub-meshes, starting from
199 // the most complex shapes and collect sub-meshes with algos that
200 // DO support sub-meshes
201 // ================================================================
203 list< SMESH_subMesh* > smWithAlgoSupportingSubmeshes[4]; // for each dim
205 // map to sort sm with same dim algos according to dim of
206 // the shape the algo assigned to (issue 0021217).
207 // Other issues influenced the algo applying order:
208 // 21406, 21556, 21893, 20206
209 multimap< int, SMESH_subMesh* > shDim2sm;
210 multimap< int, SMESH_subMesh* >::reverse_iterator shDim2smIt;
211 TopoDS_Shape algoShape;
212 int prevShapeDim = -1, aShapeDim;
214 smIt = sm->getDependsOnIterator(includeSelf, complexShapeFirst);
215 while ( smIt->more() )
217 SMESH_subMesh* smToCompute = smIt->next();
218 if ( smToCompute->GetComputeState() != SMESH_subMesh::READY_TO_COMPUTE )
221 const TopoDS_Shape& aSubShape = smToCompute->GetSubShape();
222 aShapeDim = GetShapeDim( aSubShape );
223 if ( aShapeDim < 1 ) break;
225 // check for preview dimension limitations
226 if ( aShapesId && aShapeDim > (int)aDim )
229 SMESH_Algo* algo = GetAlgo( smToCompute, &algoShape );
230 if ( algo && !algo->NeedDiscreteBoundary() )
232 if ( algo->SupportSubmeshes() )
234 // reload sub-meshes from shDim2sm into smWithAlgoSupportingSubmeshes
235 // so that more local algos to go first
236 if ( prevShapeDim != aShapeDim )
238 prevShapeDim = aShapeDim;
239 for ( shDim2smIt = shDim2sm.rbegin(); shDim2smIt != shDim2sm.rend(); ++shDim2smIt )
240 if ( shDim2smIt->first == globalAlgoDim )
241 smWithAlgoSupportingSubmeshes[ aShapeDim ].push_back( shDim2smIt->second );
243 smWithAlgoSupportingSubmeshes[ aShapeDim ].push_front( shDim2smIt->second );
246 // add smToCompute to shDim2sm map
247 if ( algoShape.IsSame( aMesh.GetShapeToMesh() ))
249 aShapeDim = globalAlgoDim; // to compute last
253 aShapeDim = GetShapeDim( algoShape );
254 if ( algoShape.ShapeType() == TopAbs_COMPOUND )
256 TopoDS_Iterator it( algoShape );
257 aShapeDim += GetShapeDim( it.Value() );
260 shDim2sm.insert( make_pair( aShapeDim, smToCompute ));
262 else // Compute w/o support of sub-meshes
264 if (_compute_canceled)
266 setCurrentSubMesh( smToCompute );
267 smToCompute->ComputeStateEngine( computeEvent );
268 setCurrentSubMesh( NULL );
270 aShapesId->insert( smToCompute->GetId() );
274 // reload sub-meshes from shDim2sm into smWithAlgoSupportingSubmeshes
275 for ( shDim2smIt = shDim2sm.rbegin(); shDim2smIt != shDim2sm.rend(); ++shDim2smIt )
276 if ( shDim2smIt->first == globalAlgoDim )
277 smWithAlgoSupportingSubmeshes[3].push_back( shDim2smIt->second );
279 smWithAlgoSupportingSubmeshes[0].push_front( shDim2smIt->second );
281 // ======================================================
282 // Apply all-dimensional algorithms supporing sub-meshes
283 // ======================================================
285 std::vector< SMESH_subMesh* > smVec;
286 for ( aShapeDim = 0; aShapeDim < 4; ++aShapeDim )
288 // ------------------------------------------------
289 // sort list of sub-meshes according to mesh order
290 // ------------------------------------------------
291 smVec.assign( smWithAlgoSupportingSubmeshes[ aShapeDim ].begin(),
292 smWithAlgoSupportingSubmeshes[ aShapeDim ].end() );
293 aMesh.SortByMeshOrder( smVec );
295 // ------------------------------------------------------------
296 // compute sub-meshes with local uni-dimensional algos under
297 // sub-meshes with all-dimensional algos
298 // ------------------------------------------------------------
299 // start from lower shapes
300 for ( size_t i = 0; i < smVec.size(); ++i )
304 // get a shape the algo is assigned to
305 if ( !GetAlgo( sm, & algoShape ))
306 continue; // strange...
308 // look for more local algos
309 smIt = sm->getDependsOnIterator(!includeSelf, !complexShapeFirst);
310 while ( smIt->more() )
312 SMESH_subMesh* smToCompute = smIt->next();
314 const TopoDS_Shape& aSubShape = smToCompute->GetSubShape();
315 const int aShapeDim = GetShapeDim( aSubShape );
316 //if ( aSubShape.ShapeType() == TopAbs_VERTEX ) continue;
317 if ( aShapeDim < 1 ) continue;
319 // check for preview dimension limitations
320 if ( aShapesId && GetShapeDim( aSubShape.ShapeType() ) > (int)aDim )
323 SMESH_HypoFilter filter( SMESH_HypoFilter::IsAlgo() );
325 .And( SMESH_HypoFilter::IsApplicableTo( aSubShape ))
326 .And( SMESH_HypoFilter::IsMoreLocalThan( algoShape, aMesh ));
328 if ( SMESH_Algo* subAlgo = (SMESH_Algo*) aMesh.GetHypothesis( smToCompute, filter, true))
330 if ( ! subAlgo->NeedDiscreteBoundary() ) continue;
331 SMESH_Hypothesis::Hypothesis_Status status;
332 if ( subAlgo->CheckHypothesis( aMesh, aSubShape, status ))
333 // mesh a lower smToCompute starting from vertices
334 Compute( aMesh, aSubShape, aShapeOnly, /*anUpward=*/true, aDim, aShapesId );
338 // --------------------------------
339 // apply the all-dimensional algos
340 // --------------------------------
341 for ( size_t i = 0; i < smVec.size(); ++i )
344 if ( sm->GetComputeState() == SMESH_subMesh::READY_TO_COMPUTE)
346 const TopAbs_ShapeEnum shapeType = sm->GetSubShape().ShapeType();
347 // check for preview dimension limitations
348 if ( aShapesId && GetShapeDim( shapeType ) > (int)aDim )
351 if (_compute_canceled)
353 setCurrentSubMesh( sm );
354 sm->ComputeStateEngine( computeEvent );
355 setCurrentSubMesh( NULL );
357 aShapesId->insert( sm->GetId() );
360 } // loop on shape dimensions
362 // -----------------------------------------------
363 // mesh the rest sub-shapes starting from vertices
364 // -----------------------------------------------
365 ret = Compute( aMesh, aShape, aShapeOnly, /*anUpward=*/true, aDim, aShapesId );
370 SMESHDS_Mesh *myMesh = aMesh.GetMeshDS();
371 //MESSAGE("*** compactMesh after compute");
372 myMesh->compactMesh();
374 // fix quadratic mesh by bending iternal links near concave boundary
375 if ( aShape.IsSame( aMesh.GetShapeToMesh() ) &&
376 !aShapesId && // not preview
377 ret ) // everything is OK
379 SMESH_MesherHelper aHelper( aMesh );
380 if ( aHelper.IsQuadraticMesh() != SMESH_MesherHelper::LINEAR )
382 aHelper.FixQuadraticElements( sm->GetComputeError() );
388 //=============================================================================
390 * Prepare Compute a mesh
392 //=============================================================================
393 void SMESH_Gen::PrepareCompute(SMESH_Mesh & aMesh,
394 const TopoDS_Shape & aShape)
396 _compute_canceled = false;
397 resetCurrentSubMesh();
400 //=============================================================================
402 * Cancel Compute a mesh
404 //=============================================================================
405 void SMESH_Gen::CancelCompute(SMESH_Mesh & aMesh,
406 const TopoDS_Shape & aShape)
408 _compute_canceled = true;
409 if ( const SMESH_subMesh* sm = GetCurrentSubMesh() )
411 const_cast< SMESH_subMesh* >( sm )->ComputeStateEngine( SMESH_subMesh::COMPUTE_CANCELED );
413 resetCurrentSubMesh();
416 //================================================================================
418 * \brief Returns a sub-mesh being currently computed
420 //================================================================================
422 const SMESH_subMesh* SMESH_Gen::GetCurrentSubMesh() const
424 return _sm_current.empty() ? 0 : _sm_current.back();
427 //================================================================================
429 * \brief Sets a sub-mesh being currently computed.
431 * An algorithm can call Compute() for a sub-shape, hence we keep a stack of sub-meshes
433 //================================================================================
435 void SMESH_Gen::setCurrentSubMesh(SMESH_subMesh* sm)
438 _sm_current.push_back( sm );
440 else if ( !_sm_current.empty() )
441 _sm_current.pop_back();
444 void SMESH_Gen::resetCurrentSubMesh()
449 //=============================================================================
453 //=============================================================================
455 bool SMESH_Gen::Evaluate(SMESH_Mesh & aMesh,
456 const TopoDS_Shape & aShape,
457 MapShapeNbElems& aResMap,
459 TSetOfInt* aShapesId)
463 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
465 const bool includeSelf = true;
466 const bool complexShapeFirst = true;
467 SMESH_subMeshIteratorPtr smIt;
469 if ( anUpward ) { // is called from below code here
470 // -----------------------------------------------
471 // mesh all the sub-shapes starting from vertices
472 // -----------------------------------------------
473 smIt = sm->getDependsOnIterator(includeSelf, !complexShapeFirst);
474 while ( smIt->more() ) {
475 SMESH_subMesh* smToCompute = smIt->next();
477 // do not mesh vertices of a pseudo shape
478 const TopAbs_ShapeEnum shapeType = smToCompute->GetSubShape().ShapeType();
479 //if ( !aMesh.HasShapeToMesh() && shapeType == TopAbs_VERTEX )
481 if ( !aMesh.HasShapeToMesh() ) {
482 if( shapeType == TopAbs_VERTEX || shapeType == TopAbs_WIRE ||
483 shapeType == TopAbs_SHELL )
487 smToCompute->Evaluate(aResMap);
489 aShapesId->insert( smToCompute->GetId() );
494 // -----------------------------------------------------------------
495 // apply algos that DO NOT require Discreteized boundaries and DO NOT
496 // support sub-meshes, starting from the most complex shapes
497 // and collect sub-meshes with algos that DO support sub-meshes
498 // -----------------------------------------------------------------
499 list< SMESH_subMesh* > smWithAlgoSupportingSubmeshes;
500 smIt = sm->getDependsOnIterator(includeSelf, complexShapeFirst);
501 while ( smIt->more() ) {
502 SMESH_subMesh* smToCompute = smIt->next();
503 const TopoDS_Shape& aSubShape = smToCompute->GetSubShape();
504 const int aShapeDim = GetShapeDim( aSubShape );
505 if ( aShapeDim < 1 ) break;
507 SMESH_Algo* algo = GetAlgo( smToCompute );
508 if ( algo && !algo->NeedDiscreteBoundary() ) {
509 if ( algo->SupportSubmeshes() ) {
510 smWithAlgoSupportingSubmeshes.push_front( smToCompute );
513 smToCompute->Evaluate(aResMap);
515 aShapesId->insert( smToCompute->GetId() );
520 // ------------------------------------------------------------
521 // sort list of meshes according to mesh order
522 // ------------------------------------------------------------
523 std::vector< SMESH_subMesh* > smVec( smWithAlgoSupportingSubmeshes.begin(),
524 smWithAlgoSupportingSubmeshes.end() );
525 aMesh.SortByMeshOrder( smVec );
527 // ------------------------------------------------------------
528 // compute sub-meshes under shapes with algos that DO NOT require
529 // Discreteized boundaries and DO support sub-meshes
530 // ------------------------------------------------------------
531 // start from lower shapes
532 for ( size_t i = 0; i < smVec.size(); ++i )
536 // get a shape the algo is assigned to
537 TopoDS_Shape algoShape;
538 if ( !GetAlgo( sm, & algoShape ))
539 continue; // strange...
541 // look for more local algos
542 smIt = sm->getDependsOnIterator(!includeSelf, !complexShapeFirst);
543 while ( smIt->more() ) {
544 SMESH_subMesh* smToCompute = smIt->next();
546 const TopoDS_Shape& aSubShape = smToCompute->GetSubShape();
547 const int aShapeDim = GetShapeDim( aSubShape );
548 if ( aShapeDim < 1 ) continue;
550 SMESH_HypoFilter filter( SMESH_HypoFilter::IsAlgo() );
552 .And( SMESH_HypoFilter::IsApplicableTo( aSubShape ))
553 .And( SMESH_HypoFilter::IsMoreLocalThan( algoShape, aMesh ));
555 if ( SMESH_Algo* subAlgo = (SMESH_Algo*) aMesh.GetHypothesis( smToCompute, filter, true ))
557 if ( ! subAlgo->NeedDiscreteBoundary() ) continue;
558 SMESH_Hypothesis::Hypothesis_Status status;
559 if ( subAlgo->CheckHypothesis( aMesh, aSubShape, status ))
560 // mesh a lower smToCompute starting from vertices
561 Evaluate( aMesh, aSubShape, aResMap, /*anUpward=*/true, aShapesId );
565 // ----------------------------------------------------------
566 // apply the algos that do not require Discreteized boundaries
567 // ----------------------------------------------------------
568 for ( size_t i = 0; i < smVec.size(); ++i )
571 sm->Evaluate(aResMap);
573 aShapesId->insert( sm->GetId() );
576 // -----------------------------------------------
577 // mesh the rest sub-shapes starting from vertices
578 // -----------------------------------------------
579 ret = Evaluate( aMesh, aShape, aResMap, /*anUpward=*/true, aShapesId );
582 MESSAGE( "VSR - SMESH_Gen::Evaluate() finished, OK = " << ret);
587 //=======================================================================
588 //function : checkConformIgnoredAlgos
590 //=======================================================================
592 static bool checkConformIgnoredAlgos(SMESH_Mesh& aMesh,
593 SMESH_subMesh* aSubMesh,
594 const SMESH_Algo* aGlobIgnoAlgo,
595 const SMESH_Algo* aLocIgnoAlgo,
597 set<SMESH_subMesh*>& aCheckedMap,
598 list< SMESH_Gen::TAlgoStateError > & theErrors)
601 if ( aSubMesh->GetSubShape().ShapeType() == TopAbs_VERTEX)
607 const list<const SMESHDS_Hypothesis*>& listHyp =
608 aMesh.GetMeshDS()->GetHypothesis( aSubMesh->GetSubShape() );
609 list<const SMESHDS_Hypothesis*>::const_iterator it=listHyp.begin();
610 for ( ; it != listHyp.end(); it++)
612 const SMESHDS_Hypothesis * aHyp = *it;
613 if (aHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO)
616 const SMESH_Algo* algo = dynamic_cast<const SMESH_Algo*> (aHyp);
619 if ( aLocIgnoAlgo ) // algo is hidden by a local algo of upper dim
621 theErrors.push_back( SMESH_Gen::TAlgoStateError() );
622 theErrors.back().Set( SMESH_Hypothesis::HYP_HIDDEN_ALGO, algo, false );
623 INFOS( "Local <" << algo->GetName() << "> is hidden by local <"
624 << aLocIgnoAlgo->GetName() << ">");
628 bool isGlobal = (aMesh.IsMainShape( aSubMesh->GetSubShape() ));
629 int dim = algo->GetDim();
630 int aMaxGlobIgnoDim = ( aGlobIgnoAlgo ? aGlobIgnoAlgo->GetDim() : -1 );
631 bool isNeededDim = ( aGlobIgnoAlgo ? aGlobIgnoAlgo->NeedLowerHyps( dim ) : false );
633 if (( dim < aMaxGlobIgnoDim && !isNeededDim ) &&
634 ( isGlobal || !aGlobIgnoAlgo->SupportSubmeshes() ))
636 // algo is hidden by a global algo
637 theErrors.push_back( SMESH_Gen::TAlgoStateError() );
638 theErrors.back().Set( SMESH_Hypothesis::HYP_HIDDEN_ALGO, algo, true );
639 INFOS( ( isGlobal ? "Global" : "Local" )
640 << " <" << algo->GetName() << "> is hidden by global <"
641 << aGlobIgnoAlgo->GetName() << ">");
643 else if ( !algo->NeedDiscreteBoundary() && !isGlobal)
645 // local algo is not hidden and hides algos on sub-shapes
646 if (checkConform && !aSubMesh->IsConform( algo ))
649 checkConform = false; // no more check conformity
650 INFOS( "ERROR: Local <" << algo->GetName() <<
651 "> would produce not conform mesh: "
652 "<Not Conform Mesh Allowed> hypotesis is missing");
653 theErrors.push_back( SMESH_Gen::TAlgoStateError() );
654 theErrors.back().Set( SMESH_Hypothesis::HYP_NOTCONFORM, algo, false );
657 // sub-algos will be hidden by a local <algo> if <algo> does not support sub-meshes
658 if ( algo->SupportSubmeshes() )
660 SMESH_subMeshIteratorPtr revItSub =
661 aSubMesh->getDependsOnIterator( /*includeSelf=*/false, /*complexShapeFirst=*/true);
662 bool checkConform2 = false;
663 while ( revItSub->more() )
665 SMESH_subMesh* sm = revItSub->next();
666 checkConformIgnoredAlgos (aMesh, sm, aGlobIgnoAlgo,
667 algo, checkConform2, aCheckedMap, theErrors);
668 aCheckedMap.insert( sm );
677 //=======================================================================
678 //function : checkMissing
679 //purpose : notify on missing hypothesis
680 // Return false if algo or hipothesis is missing
681 //=======================================================================
683 static bool checkMissing(SMESH_Gen* aGen,
685 SMESH_subMesh* aSubMesh,
686 const int aTopAlgoDim,
688 const bool checkNoAlgo,
689 set<SMESH_subMesh*>& aCheckedMap,
690 list< SMESH_Gen::TAlgoStateError > & theErrors)
692 switch ( aSubMesh->GetSubShape().ShapeType() )
696 case TopAbs_SOLID: break; // check this sub-mesh, it can be meshed
698 return true; // not meshable sub-mesh
700 if ( aCheckedMap.count( aSubMesh ))
704 SMESH_Algo* algo = 0;
706 switch (aSubMesh->GetAlgoState())
708 case SMESH_subMesh::NO_ALGO: {
711 // should there be any algo?
712 int shapeDim = SMESH_Gen::GetShapeDim( aSubMesh->GetSubShape() );
713 if (aTopAlgoDim > shapeDim)
715 MESSAGE( "ERROR: " << shapeDim << "D algorithm is missing" );
717 theErrors.push_back( SMESH_Gen::TAlgoStateError() );
718 theErrors.back().Set( SMESH_Hypothesis::HYP_MISSING, shapeDim, true );
723 case SMESH_subMesh::MISSING_HYP: {
724 // notify if an algo missing hyp is attached to aSubMesh
725 algo = aSubMesh->GetAlgo();
727 bool IsGlobalHypothesis = aGen->IsGlobalHypothesis( algo, aMesh );
728 if (!IsGlobalHypothesis || !globalChecked[ algo->GetDim() ])
730 TAlgoStateErrorName errName = SMESH_Hypothesis::HYP_MISSING;
731 SMESH_Hypothesis::Hypothesis_Status status;
732 algo->CheckHypothesis( aMesh, aSubMesh->GetSubShape(), status );
733 if ( status == SMESH_Hypothesis::HYP_BAD_PARAMETER ) {
734 MESSAGE( "ERROR: hypothesis of " << (IsGlobalHypothesis ? "Global " : "Local ")
735 << "<" << algo->GetName() << "> has a bad parameter value");
737 } else if ( status == SMESH_Hypothesis::HYP_BAD_GEOMETRY ) {
738 MESSAGE( "ERROR: " << (IsGlobalHypothesis ? "Global " : "Local ")
739 << "<" << algo->GetName() << "> assigned to mismatching geometry");
742 MESSAGE( "ERROR: " << (IsGlobalHypothesis ? "Global " : "Local ")
743 << "<" << algo->GetName() << "> misses some hypothesis");
745 if (IsGlobalHypothesis)
746 globalChecked[ algo->GetDim() ] = true;
747 theErrors.push_back( SMESH_Gen::TAlgoStateError() );
748 theErrors.back().Set( errName, algo, IsGlobalHypothesis );
753 case SMESH_subMesh::HYP_OK:
754 algo = aSubMesh->GetAlgo();
756 if (!algo->NeedDiscreteBoundary())
758 SMESH_subMeshIteratorPtr itsub = aSubMesh->getDependsOnIterator( /*includeSelf=*/false,
759 /*complexShapeFirst=*/false);
760 while ( itsub->more() )
761 aCheckedMap.insert( itsub->next() );
767 // do not check under algo that hides sub-algos or
768 // re-start checking NO_ALGO state
770 bool isTopLocalAlgo =
771 ( aTopAlgoDim <= algo->GetDim() && !aGen->IsGlobalHypothesis( algo, aMesh ));
772 if (!algo->NeedDiscreteBoundary() || isTopLocalAlgo)
774 bool checkNoAlgo2 = ( algo->NeedDiscreteBoundary() );
775 SMESH_subMeshIteratorPtr itsub = aSubMesh->getDependsOnIterator( /*includeSelf=*/false,
776 /*complexShapeFirst=*/true);
777 while ( itsub->more() )
779 // sub-meshes should not be checked further more
780 SMESH_subMesh* sm = itsub->next();
784 //check algo on sub-meshes
785 int aTopAlgoDim2 = algo->GetDim();
786 if (!checkMissing (aGen, aMesh, sm, aTopAlgoDim2,
787 globalChecked, checkNoAlgo2, aCheckedMap, theErrors))
790 if (sm->GetAlgoState() == SMESH_subMesh::NO_ALGO )
791 checkNoAlgo2 = false;
794 aCheckedMap.insert( sm );
800 //=======================================================================
801 //function : CheckAlgoState
802 //purpose : notify on bad state of attached algos, return false
803 // if Compute() would fail because of some algo bad state
804 //=======================================================================
806 bool SMESH_Gen::CheckAlgoState(SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)
808 list< TAlgoStateError > errors;
809 return GetAlgoState( aMesh, aShape, errors );
812 //=======================================================================
813 //function : GetAlgoState
814 //purpose : notify on bad state of attached algos, return false
815 // if Compute() would fail because of some algo bad state
816 // theErrors list contains problems description
817 //=======================================================================
819 bool SMESH_Gen::GetAlgoState(SMESH_Mesh& theMesh,
820 const TopoDS_Shape& theShape,
821 list< TAlgoStateError > & theErrors)
824 bool hasAlgo = false;
826 SMESH_subMesh* sm = theMesh.GetSubMesh(theShape);
827 const SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
828 TopoDS_Shape mainShape = meshDS->ShapeToMesh();
834 const SMESH_Algo* aGlobAlgoArr[] = {0,0,0,0};
836 const list<const SMESHDS_Hypothesis*>& listHyp = meshDS->GetHypothesis( mainShape );
837 list<const SMESHDS_Hypothesis*>::const_iterator it=listHyp.begin();
838 for ( ; it != listHyp.end(); it++)
840 const SMESHDS_Hypothesis * aHyp = *it;
841 if (aHyp->GetType() == SMESHDS_Hypothesis::PARAM_ALGO)
844 const SMESH_Algo* algo = dynamic_cast<const SMESH_Algo*> (aHyp);
847 int dim = algo->GetDim();
848 aGlobAlgoArr[ dim ] = algo;
853 // --------------------------------------------------------
854 // info on algos that will be ignored because of ones that
855 // don't NeedDiscreteBoundary() attached to super-shapes,
856 // check that a conform mesh will be produced
857 // --------------------------------------------------------
860 // find a global algo possibly hiding sub-algos
862 const SMESH_Algo* aGlobIgnoAlgo = 0;
863 for (dim = 3; dim > 0; dim--)
865 if (aGlobAlgoArr[ dim ] &&
866 !aGlobAlgoArr[ dim ]->NeedDiscreteBoundary() /*&&
867 !aGlobAlgoArr[ dim ]->SupportSubmeshes()*/ )
869 aGlobIgnoAlgo = aGlobAlgoArr[ dim ];
874 set<SMESH_subMesh*> aCheckedSubs;
875 bool checkConform = ( !theMesh.IsNotConformAllowed() );
877 // loop on theShape and its sub-shapes
878 SMESH_subMeshIteratorPtr revItSub = sm->getDependsOnIterator( /*includeSelf=*/true,
879 /*complexShapeFirst=*/true);
880 while ( revItSub->more() )
882 SMESH_subMesh* smToCheck = revItSub->next();
883 if ( smToCheck->GetSubShape().ShapeType() == TopAbs_VERTEX)
886 if ( aCheckedSubs.insert( smToCheck ).second ) // not yet checked
887 if (!checkConformIgnoredAlgos (theMesh, smToCheck, aGlobIgnoAlgo,
888 0, checkConform, aCheckedSubs, theErrors))
891 if ( smToCheck->GetAlgoState() != SMESH_subMesh::NO_ALGO )
895 // ----------------------------------------------------------------
896 // info on missing hypothesis and find out if all needed algos are
898 // ----------------------------------------------------------------
900 // find max dim of global algo
902 for (dim = 3; dim > 0; dim--)
904 if (aGlobAlgoArr[ dim ])
910 bool checkNoAlgo = theMesh.HasShapeToMesh() ? bool( aTopAlgoDim ) : false;
911 bool globalChecked[] = { false, false, false, false };
913 // loop on theShape and its sub-shapes
914 aCheckedSubs.clear();
915 revItSub = sm->getDependsOnIterator( /*includeSelf=*/true, /*complexShapeFirst=*/true);
916 while ( revItSub->more() )
918 SMESH_subMesh* smToCheck = revItSub->next();
919 if ( smToCheck->GetSubShape().ShapeType() == TopAbs_VERTEX)
922 if (!checkMissing (this, theMesh, smToCheck, aTopAlgoDim,
923 globalChecked, checkNoAlgo, aCheckedSubs, theErrors))
926 if (smToCheck->GetAlgoState() == SMESH_subMesh::NO_ALGO )
933 theErrors.push_back( TAlgoStateError() );
934 theErrors.back().Set( SMESH_Hypothesis::HYP_MISSING, theMesh.HasShapeToMesh() ? 1 : 3, true );
940 //=======================================================================
941 //function : IsGlobalHypothesis
942 //purpose : check if theAlgo is attached to the main shape
943 //=======================================================================
945 bool SMESH_Gen::IsGlobalHypothesis(const SMESH_Hypothesis* theHyp, SMESH_Mesh& aMesh)
947 SMESH_HypoFilter filter( SMESH_HypoFilter::Is( theHyp ));
948 return aMesh.GetHypothesis( aMesh.GetMeshDS()->ShapeToMesh(), filter, false );
951 //================================================================================
953 * \brief Return paths to xml files of plugins
955 //================================================================================
957 std::vector< std::string > SMESH_Gen::GetPluginXMLPaths()
959 // Get paths to xml files of plugins
960 vector< string > xmlPaths;
962 if ( const char* meshersList = getenv("SMESH_MeshersList") )
964 string meshers = meshersList, plugin;
965 string::size_type from = 0, pos;
966 while ( from < meshers.size() )
968 // cut off plugin name
969 pos = meshers.find( ':', from );
970 if ( pos != string::npos )
971 plugin = meshers.substr( from, pos-from );
973 plugin = meshers.substr( from ), pos = meshers.size();
976 // get PLUGIN_ROOT_DIR path
977 string rootDirVar, pluginSubDir = plugin;
978 if ( plugin == "StdMeshers" )
979 rootDirVar = "SMESH", pluginSubDir = "smesh";
981 for ( pos = 0; pos < plugin.size(); ++pos )
982 rootDirVar += toupper( plugin[pos] );
983 rootDirVar += "_ROOT_DIR";
985 const char* rootDir = getenv( rootDirVar.c_str() );
986 if ( !rootDir || strlen(rootDir) == 0 )
988 rootDirVar = plugin + "_ROOT_DIR"; // HexoticPLUGIN_ROOT_DIR
989 rootDir = getenv( rootDirVar.c_str() );
990 if ( !rootDir || strlen(rootDir) == 0 ) continue;
993 // get a separator from rootDir
994 for ( pos = strlen( rootDir )-1; pos >= 0 && sep.empty(); --pos )
995 if ( rootDir[pos] == '/' || rootDir[pos] == '\\' )
1001 if (sep.empty() ) sep = "\\";
1003 if (sep.empty() ) sep = "/";
1006 // get a path to resource file
1007 string xmlPath = rootDir;
1008 if ( xmlPath[ xmlPath.size()-1 ] != sep[0] )
1010 xmlPath += "share" + sep + "salome" + sep + "resources" + sep;
1011 for ( pos = 0; pos < pluginSubDir.size(); ++pos )
1012 xmlPath += tolower( pluginSubDir[pos] );
1013 xmlPath += sep + plugin + ".xml";
1016 fileOK = (GetFileAttributes(xmlPath.c_str()) != INVALID_FILE_ATTRIBUTES);
1018 fileOK = (access(xmlPath.c_str(), F_OK) == 0);
1021 xmlPaths.push_back( xmlPath );
1028 //=============================================================================
1030 * Finds algo to mesh a shape. Optionally returns a shape the found algo is bound to
1032 //=============================================================================
1034 SMESH_Algo *SMESH_Gen::GetAlgo(SMESH_Mesh & aMesh,
1035 const TopoDS_Shape & aShape,
1036 TopoDS_Shape* assignedTo)
1038 return GetAlgo( aMesh.GetSubMesh( aShape ), assignedTo );
1041 //=============================================================================
1043 * Finds algo to mesh a sub-mesh. Optionally returns a shape the found algo is bound to
1045 //=============================================================================
1047 SMESH_Algo *SMESH_Gen::GetAlgo(SMESH_subMesh * aSubMesh,
1048 TopoDS_Shape* assignedTo)
1050 if ( !aSubMesh ) return 0;
1052 const TopoDS_Shape & aShape = aSubMesh->GetSubShape();
1053 SMESH_Mesh& aMesh = *aSubMesh->GetFather();
1055 SMESH_HypoFilter filter( SMESH_HypoFilter::IsAlgo() );
1056 filter.And( filter.IsApplicableTo( aShape ));
1058 typedef SMESH_Algo::Features AlgoData;
1060 TopoDS_Shape assignedToShape;
1062 (SMESH_Algo*) aMesh.GetHypothesis( aSubMesh, filter, true, &assignedToShape );
1065 aShape.ShapeType() == TopAbs_FACE &&
1066 !aShape.IsSame( assignedToShape ) &&
1067 SMESH_MesherHelper::NbAncestors( aShape, aMesh, TopAbs_SOLID ) > 1 )
1069 // Issue 0021559. If there is another 2D algo with different types of output
1070 // elements that can be used to mesh aShape, and 3D algos on adjacent SOLIDs
1071 // have different types of input elements, we choose a most appropriate 2D algo.
1073 // try to find a concurrent 2D algo
1074 filter.AndNot( filter.Is( algo ));
1075 TopoDS_Shape assignedToShape2;
1077 (SMESH_Algo*) aMesh.GetHypothesis( aSubMesh, filter, true, &assignedToShape2 );
1078 if ( algo2 && // algo found
1079 !assignedToShape2.IsSame( aMesh.GetShapeToMesh() ) && // algo is local
1080 ( SMESH_MesherHelper::GetGroupType( assignedToShape2 ) == // algo of the same level
1081 SMESH_MesherHelper::GetGroupType( assignedToShape )) &&
1082 aMesh.IsOrderOK( aMesh.GetSubMesh( assignedToShape2 ), // no forced order
1083 aMesh.GetSubMesh( assignedToShape )))
1085 // get algos on the adjacent SOLIDs
1086 filter.Init( filter.IsAlgo() ).And( filter.HasDim( 3 ));
1087 vector< SMESH_Algo* > algos3D;
1088 PShapeIteratorPtr solidIt = SMESH_MesherHelper::GetAncestors( aShape, aMesh,
1090 while ( const TopoDS_Shape* solid = solidIt->next() )
1091 if ( SMESH_Algo* algo3D = (SMESH_Algo*) aMesh.GetHypothesis( *solid, filter, true ))
1093 algos3D.push_back( algo3D );
1094 filter.AndNot( filter.HasName( algo3D->GetName() ));
1096 // check compatibility of algos
1097 if ( algos3D.size() > 1 )
1099 const AlgoData& algoData = algo->SMESH_Algo::GetFeatures();
1100 const AlgoData& algoData2 = algo2->SMESH_Algo::GetFeatures();
1101 const AlgoData& algoData3d0 = algos3D[0]->SMESH_Algo::GetFeatures();
1102 const AlgoData& algoData3d1 = algos3D[1]->SMESH_Algo::GetFeatures();
1103 if (( algoData2.IsCompatible( algoData3d0 ) &&
1104 algoData2.IsCompatible( algoData3d1 ))
1106 !(algoData.IsCompatible( algoData3d0 ) &&
1107 algoData.IsCompatible( algoData3d1 )))
1113 if ( assignedTo && algo )
1114 * assignedTo = assignedToShape;
1119 //=============================================================================
1121 * Returns StudyContextStruct for a study
1123 //=============================================================================
1125 StudyContextStruct *SMESH_Gen::GetStudyContext(int studyId)
1127 // Get studyContext, create it if it does'nt exist, with a SMESHDS_Document
1129 if (_mapStudyContext.find(studyId) == _mapStudyContext.end())
1131 _mapStudyContext[studyId] = new StudyContextStruct;
1132 _mapStudyContext[studyId]->myDocument = new SMESHDS_Document(studyId);
1134 StudyContextStruct *myStudyContext = _mapStudyContext[studyId];
1135 return myStudyContext;
1138 //================================================================================
1140 * \brief Return shape dimension by TopAbs_ShapeEnum
1142 //================================================================================
1144 int SMESH_Gen::GetShapeDim(const TopAbs_ShapeEnum & aShapeType)
1146 static vector<int> dim;
1149 dim.resize( TopAbs_SHAPE, -1 );
1150 dim[ TopAbs_COMPOUND ] = MeshDim_3D;
1151 dim[ TopAbs_COMPSOLID ] = MeshDim_3D;
1152 dim[ TopAbs_SOLID ] = MeshDim_3D;
1153 dim[ TopAbs_SHELL ] = MeshDim_2D;
1154 dim[ TopAbs_FACE ] = MeshDim_2D;
1155 dim[ TopAbs_WIRE ] = MeshDim_1D;
1156 dim[ TopAbs_EDGE ] = MeshDim_1D;
1157 dim[ TopAbs_VERTEX ] = MeshDim_0D;
1159 return dim[ aShapeType ];
1162 //=============================================================================
1164 * Genarate a new id unique withing this Gen
1166 //=============================================================================
1168 int SMESH_Gen::GetANewId()