From 925813d32609dccbf255e85ce1c3dfcbc392f6ab Mon Sep 17 00:00:00 2001 From: mattijs <mattijs> Date: Mon, 15 Jan 2018 14:10:09 +0000 Subject: [PATCH] BUG: primitiveMesh: cellPoints, cellEdges inconsistent API. Fixes #703. --- src/OpenFOAM/meshes/primitiveMesh/primitiveMesh.H | 2 ++ .../meshes/primitiveMesh/primitiveMeshCellPoints.C | 13 +++++++------ .../meshes/primitiveMesh/primitiveMeshEdges.C | 13 +++++++------ .../polyTopoChange/polyTopoChange/combineFaces.C | 3 ++- .../meshRefinement/meshRefinementProblemCells.C | 3 ++- 5 files changed, 20 insertions(+), 14 deletions(-) diff --git a/src/OpenFOAM/meshes/primitiveMesh/primitiveMesh.H b/src/OpenFOAM/meshes/primitiveMesh/primitiveMesh.H index 28057107c4..23dd419349 100644 --- a/src/OpenFOAM/meshes/primitiveMesh/primitiveMesh.H +++ b/src/OpenFOAM/meshes/primitiveMesh/primitiveMesh.H @@ -822,6 +822,7 @@ public: const labelList& cellPoints ( const label celli, + labelHashSet&, DynamicList<label>& ) const; @@ -876,6 +877,7 @@ public: const labelList& cellEdges ( const label celli, + labelHashSet&, DynamicList<label>& ) const; diff --git a/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshCellPoints.C b/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshCellPoints.C index f7c0e4699a..0f6b1ce4d3 100644 --- a/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshCellPoints.C +++ b/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshCellPoints.C @@ -58,6 +58,7 @@ const Foam::labelListList& Foam::primitiveMesh::cellPoints() const const Foam::labelList& Foam::primitiveMesh::cellPoints ( const label celli, + labelHashSet& set, DynamicList<label>& storage ) const { @@ -70,7 +71,7 @@ const Foam::labelList& Foam::primitiveMesh::cellPoints const faceList& fcs = faces(); const labelList& cFaces = cells()[celli]; - labelSet_.clear(); + set.clear(); forAll(cFaces, i) { @@ -78,17 +79,17 @@ const Foam::labelList& Foam::primitiveMesh::cellPoints forAll(f, fp) { - labelSet_.insert(f[fp]); + set.insert(f[fp]); } } storage.clear(); - if (labelSet_.size() > storage.capacity()) + if (set.size() > storage.capacity()) { - storage.setCapacity(labelSet_.size()); + storage.setCapacity(set.size()); } - forAllConstIter(labelHashSet, labelSet_, iter) + forAllConstIter(labelHashSet, set, iter) { storage.append(iter.key()); } @@ -100,7 +101,7 @@ const Foam::labelList& Foam::primitiveMesh::cellPoints const Foam::labelList& Foam::primitiveMesh::cellPoints(const label celli) const { - return cellPoints(celli, labels_); + return cellPoints(celli, labelSet_, labels_); } diff --git a/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshEdges.C b/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshEdges.C index 3cf689595a..4790af73c1 100644 --- a/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshEdges.C +++ b/src/OpenFOAM/meshes/primitiveMesh/primitiveMeshEdges.C @@ -630,6 +630,7 @@ const Foam::labelList& Foam::primitiveMesh::faceEdges(const label facei) const const Foam::labelList& Foam::primitiveMesh::cellEdges ( const label celli, + labelHashSet& set, DynamicList<label>& storage ) const { @@ -641,7 +642,7 @@ const Foam::labelList& Foam::primitiveMesh::cellEdges { const labelList& cFaces = cells()[celli]; - labelSet_.clear(); + set.clear(); forAll(cFaces, i) { @@ -649,18 +650,18 @@ const Foam::labelList& Foam::primitiveMesh::cellEdges forAll(fe, feI) { - labelSet_.insert(fe[feI]); + set.insert(fe[feI]); } } storage.clear(); - if (labelSet_.size() > storage.capacity()) + if (set.size() > storage.capacity()) { - storage.setCapacity(labelSet_.size()); + storage.setCapacity(set.size()); } - forAllConstIter(labelHashSet, labelSet_, iter) + forAllConstIter(labelHashSet, set, iter) { storage.append(iter.key()); } @@ -672,7 +673,7 @@ const Foam::labelList& Foam::primitiveMesh::cellEdges const Foam::labelList& Foam::primitiveMesh::cellEdges(const label celli) const { - return cellEdges(celli, labels_); + return cellEdges(celli, labelSet_, labels_); } diff --git a/src/dynamicMesh/polyTopoChange/polyTopoChange/combineFaces.C b/src/dynamicMesh/polyTopoChange/polyTopoChange/combineFaces.C index 10e5ee81c6..6499498e3e 100644 --- a/src/dynamicMesh/polyTopoChange/polyTopoChange/combineFaces.C +++ b/src/dynamicMesh/polyTopoChange/polyTopoChange/combineFaces.C @@ -305,6 +305,7 @@ Foam::labelListList Foam::combineFaces::getMergeSets // Lists of faces that can be merged. DynamicList<labelList> allFaceSets(boundaryCells.size() / 10); // Storage for on-the-fly cell-edge addressing. + labelHashSet set; DynamicList<label> storage; // On all cells regionise the faces @@ -314,7 +315,7 @@ Foam::labelListList Foam::combineFaces::getMergeSets const cell& cFaces = mesh_.cells()[celli]; - const labelList& cEdges = mesh_.cellEdges(celli, storage); + const labelList& cEdges = mesh_.cellEdges(celli, set, storage); // Region per face Map<label> faceRegion(cFaces.size()); diff --git a/src/mesh/snappyHexMesh/meshRefinement/meshRefinementProblemCells.C b/src/mesh/snappyHexMesh/meshRefinement/meshRefinementProblemCells.C index d85faeb02c..a197d0c223 100644 --- a/src/mesh/snappyHexMesh/meshRefinement/meshRefinementProblemCells.C +++ b/src/mesh/snappyHexMesh/meshRefinement/meshRefinementProblemCells.C @@ -718,10 +718,11 @@ Foam::labelList Foam::meshRefinement::markFacesOnProblemCells // On-the-fly addressing storage. DynamicList<label> dynFEdges; DynamicList<label> dynCPoints; + labelHashSet pSet; forAll(cellLevel, celli) { - const labelList& cPoints = mesh_.cellPoints(celli, dynCPoints); + const labelList& cPoints = mesh_.cellPoints(celli, pSet, dynCPoints); // Get number of anchor points (pointLevel <= cellLevel) -- GitLab