Optimize performance of mesh deletion upon hyp modification

This commit is contained in:
eap 2016-11-09 18:59:43 +03:00
parent d2a8dd635e
commit 6b1de62331
15 changed files with 217 additions and 215 deletions

View File

@ -46,7 +46,7 @@
different dimensions. Compatible algos have equal geometries in "input" and "output".
need-hyp - (optional) Boolean. Does the algo require a hypothesis or not. Default is "false".
need-geom - (optional) [true, fasle, never]. Can the algo work w/o geometry or not.
Default is "true" "never" means that the algo can't work with geometry.
Default is "true". "never" means that the algo can't work with geometry.
support-submeshes - (optional) Boolean. Does an multi-dimensional algo support sub-meshes.
Default is "false".

View File

@ -130,7 +130,6 @@ SMDS_Mesh::SMDS_Mesh():
myNodeIDFactory(new SMDS_MeshNodeIDFactory()),
myElementIDFactory(new SMDS_MeshElementIDFactory()),
myModified(false), myModifTime(0), myCompactTime(0),
myNodeMin(0), myNodeMax(0),
myHasConstructionEdges(false), myHasConstructionFaces(false),
myHasInverseElements(true),
xmin(0), xmax(0), ymin(0), ymax(0), zmin(0), zmax(0)
@ -3424,7 +3423,7 @@ void SMDS_Mesh::RemoveFreeElement(const SMDS_MeshElement * elem)
*/
bool SMDS_Mesh::Contains (const SMDS_MeshElement* elem) const
{
// we should not imply on validity of *elem, so iterate on containers
// we should not rely on validity of *elem, so iterate on containers
// of all types in the hope of finding <elem> somewhere there
SMDS_NodeIteratorPtr itn = nodesIterator();
while (itn->more())
@ -3444,7 +3443,7 @@ bool SMDS_Mesh::Contains (const SMDS_MeshElement* elem) const
int SMDS_Mesh::MaxNodeID() const
{
return myNodeMax;
return myNodeIDFactory->GetMaxID();
}
//=======================================================================
@ -3454,7 +3453,7 @@ int SMDS_Mesh::MaxNodeID() const
int SMDS_Mesh::MinNodeID() const
{
return myNodeMin;
return myNodeIDFactory->GetMinID();
}
//=======================================================================
@ -4632,44 +4631,6 @@ SMDS_MeshVolume* SMDS_Mesh::AddVolumeWithID(const SMDS_MeshNode * n1,
return volvtk;
}
void SMDS_Mesh::updateNodeMinMax()
{
myNodeMin = 0;
if (myNodes.size() == 0)
{
myNodeMax=0;
return;
}
while ( !myNodes[myNodeMin] && myNodeMin < (int)myNodes.size() )
myNodeMin++;
myNodeMax=myNodes.size()-1;
while (!myNodes[myNodeMax] && (myNodeMin>=0))
myNodeMin--;
}
void SMDS_Mesh::incrementNodesCapacity(int nbNodes)
{
// int val = myCellIdSmdsToVtk.size();
// MESSAGE(" ------------------- resize myCellIdSmdsToVtk " << val << " --> " << val + nbNodes);
// myCellIdSmdsToVtk.resize(val + nbNodes, -1); // fill new elements with -1
int val = myNodes.size();
myNodes.resize(val +nbNodes, 0);
}
void SMDS_Mesh::incrementCellsCapacity(int nbCells)
{
int val = myCellIdVtkToSmds.size();
myCellIdVtkToSmds.resize(val + nbCells, -1); // fill new elements with -1
val = myCells.size();
myNodes.resize(val +nbCells, 0);
}
void SMDS_Mesh::adjustStructure()
{
myGrid->GetPoints()->GetData()->SetNumberOfTuples(myNodeIDFactory->GetMaxID());
}
void SMDS_Mesh::dumpGrid(string ficdump)
{
// vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
@ -4730,28 +4691,28 @@ int SMDS_Mesh::fromVtkToSmds(int vtkid)
throw SALOME_Exception(LOCALIZED ("vtk id out of bounds"));
}
void SMDS_Mesh::updateBoundingBox()
{
xmin = 0; xmax = 0;
ymin = 0; ymax = 0;
zmin = 0; zmax = 0;
vtkPoints *points = myGrid->GetPoints();
int myNodesSize = this->myNodes.size();
for (int i = 0; i < myNodesSize; i++)
{
if (SMDS_MeshNode *n = myNodes[i])
{
double coords[3];
points->GetPoint(n->myVtkID, coords);
if (coords[0] < xmin) xmin = coords[0];
else if (coords[0] > xmax) xmax = coords[0];
if (coords[1] < ymin) ymin = coords[1];
else if (coords[1] > ymax) ymax = coords[1];
if (coords[2] < zmin) zmin = coords[2];
else if (coords[2] > zmax) zmax = coords[2];
}
}
}
// void SMDS_Mesh::updateBoundingBox()
// {
// xmin = 0; xmax = 0;
// ymin = 0; ymax = 0;
// zmin = 0; zmax = 0;
// vtkPoints *points = myGrid->GetPoints();
// int myNodesSize = this->myNodes.size();
// for (int i = 0; i < myNodesSize; i++)
// {
// if (SMDS_MeshNode *n = myNodes[i])
// {
// double coords[3];
// points->GetPoint(n->myVtkID, coords);
// if (coords[0] < xmin) xmin = coords[0];
// else if (coords[0] > xmax) xmax = coords[0];
// if (coords[1] < ymin) ymin = coords[1];
// else if (coords[1] > ymax) ymax = coords[1];
// if (coords[2] < zmin) zmin = coords[2];
// else if (coords[2] > zmax) zmax = coords[2];
// }
// }
// }
double SMDS_Mesh::getMaxDim()
{

View File

@ -739,14 +739,10 @@ public:
typedef std::vector<SMDS_MeshNode *> SetOfNodes;
typedef std::vector<SMDS_MeshCell *> SetOfCells;
void updateNodeMinMax();
void updateBoundingBox();
//void updateBoundingBox();
double getMaxDim();
int fromVtkToSmds(int vtkid);
void incrementNodesCapacity(int nbNodes);
void incrementCellsCapacity(int nbCells);
void adjustStructure();
void dumpGrid(std::string ficdump="dumpGrid");
static int chunkSize;
@ -809,10 +805,10 @@ protected:
int myMeshId;
//! actual nodes coordinates, cells definition and reverse connectivity are stored in a vtkUnstructuredGrid
SMDS_UnstructuredGrid* myGrid;
SMDS_UnstructuredGrid* myGrid;
//! Small objects like SMDS_MeshNode are allocated by chunks to limit memory costs of new
ObjectPool<SMDS_MeshNode>* myNodePool;
ObjectPool<SMDS_MeshNode>* myNodePool;
//! Small objects like SMDS_VtkVolume are allocated by chunks to limit memory costs of new
ObjectPool<SMDS_VtkVolume>* myVolumePool;
@ -820,32 +816,27 @@ protected:
ObjectPool<SMDS_VtkEdge>* myEdgePool;
ObjectPool<SMDS_BallElement>* myBallPool;
//! SMDS_MeshNodes refer to vtk nodes (vtk id = index in myNodes),store reference to this mesh, and sub-shape
SetOfNodes myNodes;
//! SMDS_MeshCells refer to vtk cells (vtk id != index in myCells),store reference to this mesh, and sub-shape
SetOfCells myCells;
//! SMDS_MeshNodes refer to vtk nodes (vtk id != index in myNodes),store reference to this mesh, and sub-shape
SetOfNodes myNodes;
SetOfCells myCells;
//! a buffer to speed up elements addition by excluding some memory allocation
std::vector<vtkIdType> myNodeIds;
std::vector<vtkIdType> myNodeIds;
//! for cells only: index = ID in vtkUnstructuredGrid, value = ID for SMDS users
std::vector<int> myCellIdVtkToSmds;
std::vector<int> myCellIdVtkToSmds;
SMDS_Mesh * myParent;
std::list<SMDS_Mesh *> myChildren;
SMDS_MeshNodeIDFactory *myNodeIDFactory;
SMDS_MeshElementIDFactory *myElementIDFactory;
SMDS_MeshInfo myInfo;
SMDS_Mesh * myParent;
std::list<SMDS_Mesh *> myChildren;
SMDS_MeshNodeIDFactory * myNodeIDFactory;
SMDS_MeshElementIDFactory * myElementIDFactory;
SMDS_MeshInfo myInfo;
//! any add, remove or change of node or cell
bool myModified;
//! use a counter to keep track of modifications
unsigned long myModifTime, myCompactTime;
int myNodeMin;
int myNodeMax;
bool myHasConstructionEdges;
bool myHasConstructionFaces;
bool myHasInverseElements;

View File

@ -55,13 +55,13 @@ void SMDS_MeshElement::init(int id, ShortType meshId, LongType shapeId )
void SMDS_MeshElement::Print(ostream & OS) const
{
OS << "dump of mesh element" << endl;
OS << "dump of mesh element" << endl;
}
ostream & operator <<(ostream & OS, const SMDS_MeshElement * ME)
{
ME->Print(OS);
return OS;
ME->Print(OS);
return OS;
}
///////////////////////////////////////////////////////////////////////////////
@ -70,7 +70,7 @@ ostream & operator <<(ostream & OS, const SMDS_MeshElement * ME)
///////////////////////////////////////////////////////////////////////////////
SMDS_ElemIteratorPtr SMDS_MeshElement::nodesIterator() const
{
return elementsIterator(SMDSAbs_Node);
return elementsIterator(SMDSAbs_Node);
}
///////////////////////////////////////////////////////////////////////////////
@ -79,7 +79,7 @@ SMDS_ElemIteratorPtr SMDS_MeshElement::nodesIterator() const
///////////////////////////////////////////////////////////////////////////////
SMDS_ElemIteratorPtr SMDS_MeshElement::edgesIterator() const
{
return elementsIterator(SMDSAbs_Edge);
return elementsIterator(SMDSAbs_Edge);
}
///////////////////////////////////////////////////////////////////////////////
@ -88,7 +88,7 @@ SMDS_ElemIteratorPtr SMDS_MeshElement::edgesIterator() const
///////////////////////////////////////////////////////////////////////////////
SMDS_ElemIteratorPtr SMDS_MeshElement::facesIterator() const
{
return elementsIterator(SMDSAbs_Face);
return elementsIterator(SMDSAbs_Face);
}
///////////////////////////////////////////////////////////////////////////////
@ -96,14 +96,14 @@ SMDS_ElemIteratorPtr SMDS_MeshElement::facesIterator() const
///////////////////////////////////////////////////////////////////////////////
int SMDS_MeshElement::NbNodes() const
{
int nbnodes=0;
SMDS_ElemIteratorPtr it=nodesIterator();
while(it->more())
{
it->next();
nbnodes++;
}
return nbnodes;
int nbnodes=0;
SMDS_ElemIteratorPtr it=nodesIterator();
while(it->more())
{
it->next();
nbnodes++;
}
return nbnodes;
}
///////////////////////////////////////////////////////////////////////////////
@ -111,14 +111,14 @@ int SMDS_MeshElement::NbNodes() const
///////////////////////////////////////////////////////////////////////////////
int SMDS_MeshElement::NbEdges() const
{
int nbedges=0;
SMDS_ElemIteratorPtr it=edgesIterator();
while(it->more())
{
it->next();
nbedges++;
}
return nbedges;
int nbedges=0;
SMDS_ElemIteratorPtr it=edgesIterator();
while(it->more())
{
it->next();
nbedges++;
}
return nbedges;
}
///////////////////////////////////////////////////////////////////////////////
@ -126,14 +126,14 @@ int SMDS_MeshElement::NbEdges() const
///////////////////////////////////////////////////////////////////////////////
int SMDS_MeshElement::NbFaces() const
{
int nbfaces=0;
SMDS_ElemIteratorPtr it=facesIterator();
while(it->more())
{
it->next();
nbfaces++;
}
return nbfaces;
int nbfaces=0;
SMDS_ElemIteratorPtr it=facesIterator();
while(it->more())
{
it->next();
nbfaces++;
}
return nbfaces;
}
///////////////////////////////////////////////////////////////////////////////
@ -145,7 +145,7 @@ class SMDS_MeshElement_MyIterator:public SMDS_ElemIterator
{
const SMDS_MeshElement * myElement;
bool myMore;
public:
public:
SMDS_MeshElement_MyIterator(const SMDS_MeshElement * element):
myElement(element),myMore(true) {}
@ -228,28 +228,28 @@ SMDS_NodeIteratorPtr SMDS_MeshElement::nodeIterator() const
bool operator<(const SMDS_MeshElement& e1, const SMDS_MeshElement& e2)
{
if(e1.GetType()!=e2.GetType()) return false;
switch(e1.GetType())
{
case SMDSAbs_Node:
return static_cast<const SMDS_MeshNode &>(e1) <
static_cast<const SMDS_MeshNode &>(e2);
if(e1.GetType()!=e2.GetType()) return false;
switch(e1.GetType())
{
case SMDSAbs_Node:
return static_cast<const SMDS_MeshNode &>(e1) <
static_cast<const SMDS_MeshNode &>(e2);
case SMDSAbs_Edge:
return static_cast<const SMDS_MeshEdge &>(e1) <
static_cast<const SMDS_MeshEdge &>(e2);
case SMDSAbs_Edge:
return static_cast<const SMDS_MeshEdge &>(e1) <
static_cast<const SMDS_MeshEdge &>(e2);
case SMDSAbs_Face:
return static_cast<const SMDS_MeshFace &>(e1) <
static_cast<const SMDS_MeshFace &>(e2);
case SMDSAbs_Face:
return static_cast<const SMDS_MeshFace &>(e1) <
static_cast<const SMDS_MeshFace &>(e2);
case SMDSAbs_Volume:
return static_cast<const SMDS_MeshVolume &>(e1) <
static_cast<const SMDS_MeshVolume &>(e2);
case SMDSAbs_Volume:
return static_cast<const SMDS_MeshVolume &>(e1) <
static_cast<const SMDS_MeshVolume &>(e2);
default : MESSAGE("Internal Error");
}
return false;
default : MESSAGE("Internal Error");
}
return false;
}
bool SMDS_MeshElement::IsValidIndex(const int ind) const
@ -291,11 +291,11 @@ int SMDS_MeshElement::NbCornerNodes() const
}
//================================================================================
/*!
* \brief Check if a node belongs to the element
* \param node - the node to check
* \retval int - node index within the element, -1 if not found
*/
/*!
* \brief Check if a node belongs to the element
* \param node - the node to check
* \retval int - node index within the element, -1 if not found
*/
//================================================================================
int SMDS_MeshElement::GetNodeIndex( const SMDS_MeshNode* node ) const

View File

@ -51,13 +51,12 @@ public:
virtual void Clear();
protected:
void updateMinMax() const;
virtual void updateMinMax() const;
void updateMinMax(int id) const
{
if (id > myMax) myMax = id;
if (id < myMin) myMin = id;
}
};
#endif

View File

@ -50,8 +50,8 @@ int SMDS_MeshIDFactory::GetFreeID()
else
{
set<int>::iterator i = myPoolOfID.begin();
newid = *i;//myPoolOfID.top();
myPoolOfID.erase( i );//myPoolOfID.pop();
newid = *i;
myPoolOfID.erase( i );
}
return newid;
}
@ -77,11 +77,13 @@ void SMDS_MeshIDFactory::ReleaseID(int ID, int vtkId)
while ( i != myPoolOfID.begin() && myMaxID == *i ) {
--myMaxID; --i;
}
if ( myMaxID == *i )
if ( myMaxID == *i ) {
--myMaxID; // begin of myPoolOfID reached
else
++i;
myPoolOfID.erase( i, myPoolOfID.end() );
myPoolOfID.clear();
}
else if ( myMaxID < ID-1 ) {
myPoolOfID.erase( ++i, myPoolOfID.end() );
}
}
}
}

View File

@ -82,13 +82,10 @@ SMDS_MeshNode::~SMDS_MeshNode()
//purpose :
//=======================================================================
void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * parent)
void SMDS_MeshNode::RemoveInverseElement(const SMDS_MeshElement * elem)
{
//MESSAGE("RemoveInverseElement " << myID << " " << parent->GetID());
const SMDS_MeshCell* cell = dynamic_cast<const SMDS_MeshCell*>(parent);
MYASSERT(cell);
if ( SMDS_Mesh::_meshList[myMeshId]->getGrid()->HasLinks() )
SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, cell->getVtkId());
SMDS_Mesh::_meshList[myMeshId]->getGrid()->RemoveReferenceToCell(myVtkID, elem->getVtkId());
}
//=======================================================================
@ -313,12 +310,10 @@ vtkIdType SMDS_MeshNode::GetVtkType() const
//=======================================================================
void SMDS_MeshNode::AddInverseElement(const SMDS_MeshElement* ME)
{
const SMDS_MeshCell *cell = dynamic_cast<const SMDS_MeshCell*> (ME);
assert(cell);
SMDS_UnstructuredGrid* grid = SMDS_Mesh::_meshList[myMeshId]->getGrid();
vtkCellLinks *Links = grid->GetLinks();
Links->ResizeCellList(myVtkID, 1);
Links->AddCellReference(cell->getVtkId(), myVtkID);
Links->AddCellReference(ME->getVtkId(), myVtkID);
}
//=======================================================================

View File

@ -126,9 +126,21 @@ int SMDS_MeshNodeIDFactory::GetMinID() const
void SMDS_MeshNodeIDFactory::updateMinMax() const
{
myMesh->updateNodeMinMax();
myMin = myMesh->MinNodeID();
myMax = myMesh->MaxNodeID();
myMin = INT_MAX;
myMax = 0;
for (size_t i = 0; i < myMesh->myNodes.size(); i++)
{
if (myMesh->myNodes[i])
{
int id = myMesh->myNodes[i]->GetID();
if (id > myMax)
myMax = id;
if (id < myMin)
myMin = id;
}
}
if (myMin == INT_MAX)
myMin = 0;
}
SMDS_ElemIteratorPtr SMDS_MeshNodeIDFactory::elementsIterator() const

View File

@ -49,7 +49,7 @@ public:
virtual void emptyPool(int maxId);
protected:
void updateMinMax() const;
virtual void updateMinMax() const;
void updateMinMax(int id) const
{
if (id > myMax)

View File

@ -139,10 +139,12 @@ vtkPoints* SMDS_UnstructuredGrid::GetPoints()
int SMDS_UnstructuredGrid::InsertNextLinkedCell(int type, int npts, vtkIdType *pts)
{
if ( !this->Links )
BuildLinks();
if ( !this->Links ) // don't create Links until they are needed
{
return this->InsertNextCell(type, npts, pts);
}
if (type != VTK_POLYHEDRON)
if ( type != VTK_POLYHEDRON )
return vtkUnstructuredGrid::InsertNextLinkedCell(type, npts, pts);
// --- type = VTK_POLYHEDRON

View File

@ -1208,48 +1208,56 @@ void SMESH_Mesh::NotifySubMeshesHypothesisModification(const SMESH_Hypothesis* h
SMESH_Algo *algo;
const SMESH_HypoFilter* compatibleHypoKind;
list <const SMESHDS_Hypothesis * > usedHyps;
// keep sub-meshes not to miss ones whose state can change due to notifying others
vector< SMESH_subMesh* > smToNotify;
bool allMeshedEdgesNotified = true;
SMESH_subMeshIteratorPtr smIt( _subMeshHolder->GetIterator() );
while ( smIt->more() )
{
SMESH_subMesh* aSubMesh = smIt->next();
bool toNotify = false;
// if aSubMesh meshing depends on hyp,
// we call aSubMesh->AlgoStateEngine( MODIF_HYP, hyp ) that causes either
// 1) clearing already computed aSubMesh or
// 2) changing algo_state from MISSING_HYP to HYP_OK when parameters of hyp becomes valid,
// other possible changes are not interesting. (IPAL0052457 - assigning hyp performance pb)
if ( aSubMesh->GetComputeState() != SMESH_subMesh::COMPUTE_OK &&
aSubMesh->GetComputeState() != SMESH_subMesh::FAILED_TO_COMPUTE &&
aSubMesh->GetAlgoState() != SMESH_subMesh::MISSING_HYP &&
!hyp->DataDependOnParams() )
continue;
const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
if (( aSubMesh->IsApplicableHypotesis( hyp )) &&
( algo = aSubMesh->GetAlgo() ) &&
( compatibleHypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() )) &&
( compatibleHypoKind->IsOk( hyp, aSubShape )))
if ( aSubMesh->GetComputeState() == SMESH_subMesh::COMPUTE_OK ||
aSubMesh->GetComputeState() == SMESH_subMesh::FAILED_TO_COMPUTE ||
aSubMesh->GetAlgoState() == SMESH_subMesh::MISSING_HYP ||
hyp->DataDependOnParams() )
{
// check if hyp is used by algo
usedHyps.clear();
if ( GetHypotheses( aSubMesh, *compatibleHypoKind, usedHyps, true ) &&
find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() )
const TopoDS_Shape & aSubShape = aSubMesh->GetSubShape();
if (( aSubMesh->IsApplicableHypotesis( hyp )) &&
( algo = aSubMesh->GetAlgo() ) &&
( compatibleHypoKind = algo->GetCompatibleHypoFilter( !hyp->IsAuxiliary() )) &&
( compatibleHypoKind->IsOk( hyp, aSubShape )))
{
smToNotify.push_back( aSubMesh );
// check if hyp is used by algo
usedHyps.clear();
toNotify = ( GetHypotheses( aSubMesh, *compatibleHypoKind, usedHyps, true ) &&
std::find( usedHyps.begin(), usedHyps.end(), hyp ) != usedHyps.end() );
}
}
if ( toNotify )
smToNotify.push_back( aSubMesh );
if ( !aSubMesh->IsEmpty() &&
aSubMesh->GetSubShape().ShapeType() == TopAbs_EDGE &&
!toNotify )
allMeshedEdgesNotified = false;
}
for ( size_t i = 0; i < smToNotify.size(); ++i )
{
smToNotify[i]->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
const_cast< SMESH_Hypothesis*>( hyp ));
}
// if all meshed EDGEs will be notified then the notification is equivalent
// to the whole mesh clearing
if ( allMeshedEdgesNotified )
Clear();
else
// notify in reverse order to avoid filling of the pool of IDs
for ( int i = smToNotify.size()-1; i >= 0; --i )
smToNotify[i]->AlgoStateEngine(SMESH_subMesh::MODIF_HYP,
const_cast< SMESH_Hypothesis*>( hyp ));
HasModificationsToDiscard(); // to reset _isModified flag if mesh becomes empty
GetMeshDS()->Modified();

View File

@ -1272,24 +1272,34 @@ void SMESH_subMesh::DumpAlgoState(bool isMain)
static void cleanSubMesh( SMESH_subMesh * subMesh )
{
if (subMesh) {
if (SMESHDS_SubMesh * subMeshDS = subMesh->GetSubMeshDS()) {
if (SMESHDS_SubMesh * subMeshDS = subMesh->GetSubMeshDS())
{
SMESHDS_Mesh * meshDS = subMesh->GetFather()->GetMeshDS();
SMDS_ElemIteratorPtr ite = subMeshDS->GetElements();
while (ite->more()) {
const SMDS_MeshElement * elt = ite->next();
//MESSAGE( " RM elt: "<<elt->GetID()<<" ( "<<elt->NbNodes()<<" )" );
//meshDS->RemoveElement(elt);
meshDS->RemoveFreeElement(elt, 0);
int nbElems = subMeshDS->NbElements();
if ( nbElems > 0 )
{
// start from elem with max ID to avoid filling the pool of IDs
const SMDS_MeshElement * lastElem = subMeshDS->GetElement( nbElems-1 );
bool rev = ( lastElem->GetID() == meshDS->MaxElementID() );
SMDS_ElemIteratorPtr ite = subMeshDS->GetElements( rev );
while (ite->more()) {
const SMDS_MeshElement * elt = ite->next();
meshDS->RemoveFreeElement(elt, 0);
}
}
SMDS_NodeIteratorPtr itn = subMeshDS->GetNodes();
while (itn->more()) {
const SMDS_MeshNode * node = itn->next();
//MESSAGE( " RM node: "<<node->GetID());
if ( node->NbInverseElements() == 0 )
meshDS->RemoveFreeNode(node, 0);
else // for StdMeshers_CompositeSegment_1D: node in one submesh, edge in another
meshDS->RemoveNode(node);
int nbNodes = subMeshDS->NbNodes();
if ( nbNodes > 0 )
{
const SMDS_MeshNode * lastNode = subMeshDS->GetNode( nbNodes-1 );
bool rev = ( lastNode->GetID() == meshDS->MaxNodeID() );
SMDS_NodeIteratorPtr itn = subMeshDS->GetNodes( rev );
while (itn->more()) {
const SMDS_MeshNode * node = itn->next();
if ( node->NbInverseElements() == 0 )
meshDS->RemoveFreeNode(node, 0);
else // for StdMeshers_CompositeSegment_1D: node in one submesh, edge in another
meshDS->RemoveNode(node);
}
}
subMeshDS->Clear();
}

View File

@ -361,13 +361,23 @@ public:
//purpose :
//=======================================================================
SMDS_ElemIteratorPtr SMESHDS_SubMesh::GetElements() const
SMDS_ElemIteratorPtr SMESHDS_SubMesh::GetElements( bool reverse ) const
{
if ( IsComplexSubmesh() )
return SMDS_ElemIteratorPtr( new MyElemIterator( mySubMeshes ));
return SMDS_ElemIteratorPtr
(new MySetIterator<const SMDS_MeshElement*, std::vector<const SMDS_MeshElement*> >(myElements));
if ( reverse )
{
typedef
SMDS_SetIterator< const SMDS_MeshElement*,
std::vector< const SMDS_MeshElement* >::const_reverse_iterator > RIter;
return SMDS_ElemIteratorPtr( new RIter( myElements.rbegin(), myElements.rend() ));
}
typedef
SMDS_SetIterator< const SMDS_MeshElement*,
std::vector< const SMDS_MeshElement* >::const_iterator > FIter;
return SMDS_ElemIteratorPtr( new FIter( myElements.begin(), myElements.end() ));
}
//=======================================================================
@ -375,13 +385,23 @@ SMDS_ElemIteratorPtr SMESHDS_SubMesh::GetElements() const
//purpose :
//=======================================================================
SMDS_NodeIteratorPtr SMESHDS_SubMesh::GetNodes() const
SMDS_NodeIteratorPtr SMESHDS_SubMesh::GetNodes( bool reverse ) const
{
if ( IsComplexSubmesh() )
return SMDS_NodeIteratorPtr( new MyNodeIterator( mySubMeshes ));
return SMDS_NodeIteratorPtr
(new MySetIterator<const SMDS_MeshNode*, std::vector<const SMDS_MeshNode*> >(myNodes));
if ( reverse )
{
typedef
SMDS_SetIterator< const SMDS_MeshNode*,
std::vector< const SMDS_MeshNode* >::const_reverse_iterator > RIter;
return SMDS_NodeIteratorPtr( new RIter( myNodes.rbegin(), myNodes.rend() ));
}
typedef
SMDS_SetIterator< const SMDS_MeshNode*,
std::vector< const SMDS_MeshNode* >::const_iterator > FIter;
return SMDS_NodeIteratorPtr( new FIter( myNodes.begin(), myNodes.end() ));
}
//=======================================================================
@ -478,7 +498,7 @@ void SMESHDS_SubMesh::RemoveAllSubmeshes()
//=======================================================================
//function : ContainsSubMesh
//purpose :
//purpose :
//=======================================================================
bool SMESHDS_SubMesh::ContainsSubMesh( const SMESHDS_SubMesh* theSubMesh ) const
@ -488,7 +508,7 @@ bool SMESHDS_SubMesh::ContainsSubMesh( const SMESHDS_SubMesh* theSubMesh ) const
//=======================================================================
//function : GetSubMeshIterator
//purpose :
//purpose :
//=======================================================================
SMESHDS_SubMeshIteratorPtr SMESHDS_SubMesh::GetSubMeshIterator() const

View File

@ -65,9 +65,9 @@ class SMESHDS_EXPORT SMESHDS_SubMesh
// for both types
virtual int NbElements() const;
virtual SMDS_ElemIteratorPtr GetElements() const;
virtual SMDS_ElemIteratorPtr GetElements(bool reverse=false) const;
virtual int NbNodes() const;
virtual SMDS_NodeIteratorPtr GetNodes() const;
virtual SMDS_NodeIteratorPtr GetNodes(bool reverse=false) const;
virtual bool Contains(const SMDS_MeshElement * ME) const; // check if elem or node is in
virtual bool IsQuadratic() const;

View File

@ -94,7 +94,9 @@ void SMESHGUI_GenericHypothesisCreator::create( SMESH::SMESH_Hypothesis_ptr init
void SMESHGUI_GenericHypothesisCreator::create( bool isAlgo,
const QString& theHypName,
QWidget* theParent, QObject* obj, const QString& slot )
QWidget* theParent,
QObject* obj,
const QString& slot )
{
myIsCreate = true;