mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
STYLE: use meshCells() instead of cutCells() for consistency
- other sampled surface types use meshCells() throughout. Only cuttingPlane was different.
This commit is contained in:
@ -55,7 +55,7 @@ void Foam::cuttingPlane::calcCutCells
|
|||||||
listSize = cellIdLabels.size();
|
listSize = cellIdLabels.size();
|
||||||
}
|
}
|
||||||
|
|
||||||
cutCells_.setSize(listSize);
|
meshCells_.setSize(listSize);
|
||||||
label cutcelli(0);
|
label cutcelli(0);
|
||||||
|
|
||||||
// Find the cut cells by detecting any cell that uses points with
|
// Find the cut cells by detecting any cell that uses points with
|
||||||
@ -70,9 +70,7 @@ void Foam::cuttingPlane::calcCutCells
|
|||||||
}
|
}
|
||||||
|
|
||||||
const labelList& cEdges = cellEdges[celli];
|
const labelList& cEdges = cellEdges[celli];
|
||||||
|
|
||||||
label nCutEdges = 0;
|
label nCutEdges = 0;
|
||||||
|
|
||||||
forAll(cEdges, i)
|
forAll(cEdges, i)
|
||||||
{
|
{
|
||||||
const edge& e = edges[cEdges[i]];
|
const edge& e = edges[cEdges[i]];
|
||||||
@ -87,8 +85,7 @@ void Foam::cuttingPlane::calcCutCells
|
|||||||
|
|
||||||
if (nCutEdges > 2)
|
if (nCutEdges > 2)
|
||||||
{
|
{
|
||||||
cutCells_[cutcelli++] = celli;
|
meshCells_[cutcelli++] = celli;
|
||||||
|
|
||||||
break;
|
break;
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -96,7 +93,7 @@ void Foam::cuttingPlane::calcCutCells
|
|||||||
}
|
}
|
||||||
|
|
||||||
// Set correct list size
|
// Set correct list size
|
||||||
cutCells_.setSize(cutcelli);
|
meshCells_.setSize(cutcelli);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
@ -114,7 +111,7 @@ void Foam::cuttingPlane::intersectEdges
|
|||||||
// Per edge -1 or the label of the intersection point
|
// Per edge -1 or the label of the intersection point
|
||||||
edgePoint.setSize(edges.size());
|
edgePoint.setSize(edges.size());
|
||||||
|
|
||||||
DynamicList<point> dynCuttingPoints(4*cutCells_.size());
|
DynamicList<point> dynCuttingPoints(4*meshCells_.size());
|
||||||
|
|
||||||
forAll(edges, edgeI)
|
forAll(edges, edgeI)
|
||||||
{
|
{
|
||||||
@ -258,15 +255,15 @@ void Foam::cuttingPlane::walkCellCuts
|
|||||||
const pointField& cutPoints = this->points();
|
const pointField& cutPoints = this->points();
|
||||||
|
|
||||||
// use dynamic lists to handle triangulation and/or missed cuts
|
// use dynamic lists to handle triangulation and/or missed cuts
|
||||||
DynamicList<face> dynCutFaces(cutCells_.size());
|
DynamicList<face> dynCutFaces(meshCells_.size());
|
||||||
DynamicList<label> dynCutCells(cutCells_.size());
|
DynamicList<label> dynCutCells(meshCells_.size());
|
||||||
|
|
||||||
// scratch space for calculating the face vertices
|
// scratch space for calculating the face vertices
|
||||||
DynamicList<label> faceVerts(10);
|
DynamicList<label> faceVerts(10);
|
||||||
|
|
||||||
forAll(cutCells_, i)
|
forAll(meshCells_, i)
|
||||||
{
|
{
|
||||||
label celli = cutCells_[i];
|
label celli = meshCells_[i];
|
||||||
|
|
||||||
// Find the starting edge to walk from.
|
// Find the starting edge to walk from.
|
||||||
const labelList& cEdges = mesh.cellEdges()[celli];
|
const labelList& cEdges = mesh.cellEdges()[celli];
|
||||||
@ -330,20 +327,18 @@ void Foam::cuttingPlane::walkCellCuts
|
|||||||
}
|
}
|
||||||
|
|
||||||
this->storedFaces().transfer(dynCutFaces);
|
this->storedFaces().transfer(dynCutFaces);
|
||||||
cutCells_.transfer(dynCutCells);
|
meshCells_.transfer(dynCutCells);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||||
|
|
||||||
// Construct without cutting
|
|
||||||
Foam::cuttingPlane::cuttingPlane(const plane& pln)
|
Foam::cuttingPlane::cuttingPlane(const plane& pln)
|
||||||
:
|
:
|
||||||
plane(pln)
|
plane(pln)
|
||||||
{}
|
{}
|
||||||
|
|
||||||
|
|
||||||
// Construct from plane and mesh reference, restricted to a list of cells
|
|
||||||
Foam::cuttingPlane::cuttingPlane
|
Foam::cuttingPlane::cuttingPlane
|
||||||
(
|
(
|
||||||
const plane& pln,
|
const plane& pln,
|
||||||
@ -369,7 +364,7 @@ void Foam::cuttingPlane::reCut
|
|||||||
)
|
)
|
||||||
{
|
{
|
||||||
MeshStorage::clear();
|
MeshStorage::clear();
|
||||||
cutCells_.clear();
|
meshCells_.clear();
|
||||||
|
|
||||||
const scalarField dotProducts((mesh.points() - refPoint()) & normal());
|
const scalarField dotProducts((mesh.points() - refPoint()) & normal());
|
||||||
|
|
||||||
@ -391,7 +386,7 @@ void Foam::cuttingPlane::remapFaces
|
|||||||
const labelUList& faceMap
|
const labelUList& faceMap
|
||||||
)
|
)
|
||||||
{
|
{
|
||||||
// recalculate the cells cut
|
// Recalculate the cells cut
|
||||||
if (notNull(faceMap) && faceMap.size())
|
if (notNull(faceMap) && faceMap.size())
|
||||||
{
|
{
|
||||||
MeshStorage::remapFaces(faceMap);
|
MeshStorage::remapFaces(faceMap);
|
||||||
@ -399,9 +394,9 @@ void Foam::cuttingPlane::remapFaces
|
|||||||
List<label> newCutCells(faceMap.size());
|
List<label> newCutCells(faceMap.size());
|
||||||
forAll(faceMap, facei)
|
forAll(faceMap, facei)
|
||||||
{
|
{
|
||||||
newCutCells[facei] = cutCells_[faceMap[facei]];
|
newCutCells[facei] = meshCells_[faceMap[facei]];
|
||||||
}
|
}
|
||||||
cutCells_.transfer(newCutCells);
|
meshCells_.transfer(newCutCells);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
@ -420,7 +415,7 @@ void Foam::cuttingPlane::operator=(const cuttingPlane& rhs)
|
|||||||
|
|
||||||
static_cast<MeshStorage&>(*this) = rhs;
|
static_cast<MeshStorage&>(*this) = rhs;
|
||||||
static_cast<plane&>(*this) = rhs;
|
static_cast<plane&>(*this) = rhs;
|
||||||
cutCells_ = rhs.cutCells();
|
meshCells_ = rhs.meshCells();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
@ -70,7 +70,7 @@ class cuttingPlane
|
|||||||
// Private data
|
// Private data
|
||||||
|
|
||||||
//- List of cells cut by the plane
|
//- List of cells cut by the plane
|
||||||
labelList cutCells_;
|
labelList meshCells_;
|
||||||
|
|
||||||
|
|
||||||
// Private Member Functions
|
// Private Member Functions
|
||||||
@ -121,7 +121,7 @@ protected:
|
|||||||
|
|
||||||
// Protected Member Functions
|
// Protected Member Functions
|
||||||
|
|
||||||
//- Recut mesh with existing planeDesc, restricted to a list of cells
|
//- Recut mesh with existing plane, restricted to a list of cells
|
||||||
void reCut
|
void reCut
|
||||||
(
|
(
|
||||||
const primitiveMesh&,
|
const primitiveMesh&,
|
||||||
@ -157,15 +157,15 @@ public:
|
|||||||
}
|
}
|
||||||
|
|
||||||
//- Return List of cells cut by the plane
|
//- Return List of cells cut by the plane
|
||||||
const labelList& cutCells() const
|
const labelList& meshCells() const
|
||||||
{
|
{
|
||||||
return cutCells_;
|
return meshCells_;
|
||||||
}
|
}
|
||||||
|
|
||||||
//- Return true or false to question: have any cells been cut?
|
//- Return true or false to question: have any cells been cut?
|
||||||
bool cut() const
|
bool cut() const
|
||||||
{
|
{
|
||||||
return cutCells_.size();
|
return meshCells_.size();
|
||||||
}
|
}
|
||||||
|
|
||||||
//- Sample the cell field
|
//- Sample the cell field
|
||||||
|
|||||||
@ -28,8 +28,6 @@ Description
|
|||||||
|
|
||||||
#include "cuttingPlane.H"
|
#include "cuttingPlane.H"
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
||||||
|
|
||||||
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
|
||||||
|
|
||||||
template<class Type>
|
template<class Type>
|
||||||
@ -38,7 +36,7 @@ Foam::tmp<Foam::Field<Type>> Foam::cuttingPlane::sample
|
|||||||
const Field<Type>& fld
|
const Field<Type>& fld
|
||||||
) const
|
) const
|
||||||
{
|
{
|
||||||
return tmp<Field<Type>>(new Field<Type>(fld, cutCells()));
|
return tmp<Field<Type>>(new Field<Type>(fld, meshCells()));
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|||||||
@ -169,7 +169,7 @@ public:
|
|||||||
//- For every face original cell in mesh
|
//- For every face original cell in mesh
|
||||||
const labelList& meshCells() const
|
const labelList& meshCells() const
|
||||||
{
|
{
|
||||||
return cuttingPlane::cutCells();
|
return cuttingPlane::meshCells();
|
||||||
}
|
}
|
||||||
|
|
||||||
//- Sample field on surface
|
//- Sample field on surface
|
||||||
|
|||||||
Reference in New Issue
Block a user