polyCellSet: General cell set selection class

Description
    General cell set selection class for models that apply to sub-sets
    of the mesh.

    Currently supports cell selection from a set of points, a specified cellSet
    or cellZone or all of the cells.  The selection method can either be
    specified explicitly using the \c selectionMode entry or inferred from the
    presence of either a \c cellSet, \c cellZone or \c points entry.  The \c
    selectionMode entry is required to select \c all cells.

Usage
    Examples:
    \verbatim
        // Apply everywhere
        selectionMode   all;

        // Apply within a given cellSet
        selectionMode   cellSet; // Optional
        cellSet         rotor;

        // Apply within a given cellZone
        selectionMode   cellZone; // Optional
        cellSet         rotor;

        // Apply in cells containing a list of points
        selectionMode   points; // Optional
        points
        (
            (2.25 0.5 0)
            (2.75 0.5 0)
        );
    \endverbatim

Also used as the base-class for fvCellSet which additionally provides and
maintains the volume of the cell set.
This commit is contained in:
Henry Weller
2022-08-13 16:32:19 +01:00
parent 7fdde885fe
commit ac0eea9610
20 changed files with 577 additions and 351 deletions

View File

@ -24,13 +24,11 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "solidBodyMotionSolver.H" #include "solidBodyMotionSolver.H"
#include "addToRunTimeSelectionTable.H" #include "polyCellSet.H"
#include "transformField.H" #include "transformField.H"
#include "meshCellZones.H"
#include "cellSet.H"
#include "boolList.H"
#include "syncTools.H" #include "syncTools.H"
#include "polyTopoChangeMap.H" #include "polyTopoChangeMap.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
@ -61,52 +59,10 @@ Foam::solidBodyMotionSolver::solidBodyMotionSolver
moveAllCells_(false), moveAllCells_(false),
transform_(SBMFPtr_().transformation()) transform_(SBMFPtr_().transformation())
{ {
const word cellZoneName = const labelList cellIDs(polyCellSet(mesh, dict).cells());
coeffDict().lookupOrDefault<word>("cellZone", "none");
const word cellSetName = moveAllCells_ =
coeffDict().lookupOrDefault<word>("cellSet", "none"); returnReduce(cellIDs.size() == mesh.nCells(), andOp<bool>());
if ((cellZoneName != "none") && (cellSetName != "none"))
{
FatalIOErrorInFunction(coeffDict())
<< "Either cellZone OR cellSet can be supplied, but not both. "
<< "If neither is supplied, all cells will be included"
<< exit(FatalIOError);
}
labelList cellIDs;
if (cellZoneName != "none")
{
Info<< "Applying solid body motion to cellZone " << cellZoneName
<< endl;
label zoneID = mesh.cellZones().findZoneID(cellZoneName);
if (zoneID == -1)
{
FatalErrorInFunction
<< "Unable to find cellZone " << cellZoneName
<< ". Valid cellZones are:"
<< mesh.cellZones().names()
<< exit(FatalError);
}
cellIDs = mesh.cellZones()[zoneID];
}
if (cellSetName != "none")
{
Info<< "Applying solid body motion to cellSet " << cellSetName
<< endl;
cellSet set(mesh, cellSetName);
cellIDs = set.toc();
}
const label nCells = returnReduce(cellIDs.size(), sumOp<label>());
moveAllCells_ = nCells == 0;
if (moveAllCells_) if (moveAllCells_)
{ {
@ -120,15 +76,13 @@ Foam::solidBodyMotionSolver::solidBodyMotionSolver
forAll(cellIDs, i) forAll(cellIDs, i)
{ {
label celli = cellIDs[i]; const cell& c = mesh.cells()[cellIDs[i]];
const cell& c = mesh.cells()[celli];
forAll(c, j) forAll(c, j)
{ {
const face& f = mesh.faces()[c[j]]; const face& f = mesh.faces()[c[j]];
forAll(f, k) forAll(f, k)
{ {
label pointi = f[k]; movePts[f[k]] = true;
movePts[pointi] = true;
} }
} }
} }

View File

@ -26,116 +26,23 @@ License
#include "fvCellSet.H" #include "fvCellSet.H"
#include "volFields.H" #include "volFields.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
template<> const char* NamedEnum
<
fvCellSet::selectionModeType,
4
>::names[] =
{
"points",
"cellSet",
"cellZone",
"all"
};
const NamedEnum<fvCellSet::selectionModeType, 4>
fvCellSet::selectionModeTypeNames_;
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::fvCellSet::setCells()
{
Info<< incrIndent;
switch (selectionMode_)
{
case selectionModeType::points:
{
Info<< indent << "- selecting cells using points" << endl;
labelHashSet selectedCells;
forAll(points_, i)
{
label celli = mesh_.findCell(points_[i]);
if (celli >= 0)
{
selectedCells.insert(celli);
}
label globalCelli = returnReduce(celli, maxOp<label>());
if (globalCelli < 0)
{
WarningInFunction
<< "Unable to find owner cell for point " << points_[i]
<< endl;
}
}
cells_ = selectedCells.toc();
break;
}
case selectionModeType::cellSet:
{
Info<< indent
<< "- selecting cells using cellSet " << cellSetName_ << endl;
cellSet selectedCells(mesh_, cellSetName_);
cells_ = selectedCells.toc();
break;
}
case selectionModeType::cellZone:
{
Info<< indent
<< "- selecting cells using cellZone " << cellSetName_ << endl;
label zoneID = mesh_.cellZones().findZoneID(cellSetName_);
if (zoneID == -1)
{
FatalErrorInFunction
<< "Cannot find cellZone " << cellSetName_ << endl
<< "Valid cellZones are " << mesh_.cellZones().names()
<< exit(FatalError);
}
cells_ = mesh_.cellZones()[zoneID];
break;
}
case selectionModeType::all:
{
Info<< indent << "- selecting all cells" << endl;
cells_ = identity(mesh_.nCells());
break;
}
}
Info<< decrIndent;
}
void Foam::fvCellSet::setV() void Foam::fvCellSet::setV()
{ {
Info<< incrIndent; Info<< incrIndent;
const labelList& cells = this->cells();
V_ = 0; V_ = 0;
forAll(cells_, i) forAll(cells, i)
{ {
V_ += mesh_.V()[cells_[i]]; V_ += mesh_.V()[cells[i]];
} }
reduce(V_, sumOp<scalar>()); reduce(V_, sumOp<scalar>());
Info<< indent Info<< indent
<< "- selected " << returnReduce(cells_.size(), sumOp<label>()) << "- selected " << returnReduce(cells.size(), sumOp<label>())
<< " cell(s) with volume " << V_ << endl; << " cell(s) with volume " << V_ << endl;
Info<< decrIndent; Info<< decrIndent;
@ -150,12 +57,11 @@ Foam::fvCellSet::fvCellSet
const dictionary& dict const dictionary& dict
) )
: :
polyCellSet(mesh, dict),
mesh_(mesh), mesh_(mesh),
selectionMode_(selectionModeType::all),
cellSetName_(word::null),
V_(NaN) V_(NaN)
{ {
read(dict); setV();
} }
@ -169,94 +75,35 @@ Foam::fvCellSet::~fvCellSet()
void Foam::fvCellSet::movePoints() void Foam::fvCellSet::movePoints()
{ {
if (selectionMode_ == selectionModeType::points) polyCellSet::movePoints();
{
setCells();
}
setV(); setV();
} }
void Foam::fvCellSet::topoChange(const polyTopoChangeMap&) void Foam::fvCellSet::topoChange(const polyTopoChangeMap& map)
{ {
setCells(); polyCellSet::topoChange(map);
setV(); setV();
} }
void Foam::fvCellSet::mapMesh(const polyMeshMap&) void Foam::fvCellSet::mapMesh(const polyMeshMap& map)
{ {
setCells(); polyCellSet::mapMesh(map);
setV(); setV();
} }
void Foam::fvCellSet::distribute(const polyDistributionMap&) void Foam::fvCellSet::distribute(const polyDistributionMap& map)
{ {
setCells(); polyCellSet::distribute(map);
setV(); setV();
} }
bool Foam::fvCellSet::read(const dictionary& dict) bool Foam::fvCellSet::read(const dictionary& dict)
{ {
if (dict.found("selectionMode")) polyCellSet::read(dict);
{
selectionMode_ =
selectionModeTypeNames_.read(dict.lookup("selectionMode"));
}
else if (dict.found("points"))
{
selectionMode_ = selectionModeType::points;
}
else if (dict.found("cellSet"))
{
selectionMode_ = selectionModeType::cellSet;
}
else if (dict.lookup("cellZone"))
{
selectionMode_ = selectionModeType::cellZone;
}
else
{
FatalIOErrorInFunction(dict)
<< "selectionMode, points, cellSet or cellZone not specified. "
<< "Valid selectionMode types are" << selectionModeTypeNames_
<< exit(FatalIOError);
}
switch (selectionMode_)
{
case selectionModeType::points:
{
dict.lookup("points") >> points_;
break;
}
case selectionModeType::cellSet:
{
dict.lookup("cellSet") >> cellSetName_;
break;
}
case selectionModeType::cellZone:
{
dict.lookup("cellZone") >> cellSetName_;
break;
}
case selectionModeType::all:
{
break;
}
default:
{
FatalErrorInFunction
<< "Unknown selectionMode "
<< selectionModeTypeNames_[selectionMode_]
<< ". Valid selectionMode types are" << selectionModeTypeNames_
<< exit(FatalError);
}
}
setCells();
setV(); setV();
return true; return true;

View File

@ -22,7 +22,7 @@ License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class Class
Foam::fv::fvCellSet Foam::fvCellSet
Description Description
General cell set selection class for models that apply to sub-sets General cell set selection class for models that apply to sub-sets
@ -57,6 +57,9 @@ Usage
); );
\endverbatim \endverbatim
See also
Foam::polyCellSet
SourceFiles SourceFiles
fvCellSet.C fvCellSet.C
@ -65,7 +68,7 @@ SourceFiles
#ifndef fvCellSet_H #ifndef fvCellSet_H
#define fvCellSet_H #define fvCellSet_H
#include "cellSet.H" #include "polyCellSet.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -81,55 +84,19 @@ class polyDistributionMap;
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
class fvCellSet class fvCellSet
:
public polyCellSet
{ {
public:
// Public data
//- Enumeration for selection mode types
enum class selectionModeType
{
points,
cellSet,
cellZone,
all
};
//- Word list of selection mode type names
static const NamedEnum<selectionModeType, 4>
selectionModeTypeNames_;
private:
// Private data // Private data
const fvMesh& mesh_; const fvMesh& mesh_;
//- Cell selection mode
selectionModeType selectionMode_;
//- Name of cell set for "cellSet" and "cellZone" selectionMode
word cellSetName_;
//- List of points for "points" selectionMode
List<point> points_;
//- Set of cells to apply source to
mutable labelList cells_;
//- Sum of cell volumes //- Sum of cell volumes
mutable scalar V_; mutable scalar V_;
// Private functions // Private functions
//- Read the coefficients from the given dictionary
void readCoeffs(const dictionary& dict);
//- Set the cells
void setCells();
//- Set the sum of scalar volumes //- Set the sum of scalar volumes
void setV(); void setV();
@ -154,19 +121,9 @@ public:
// Access // Access
//- Return const access to the cell selection mode
inline const selectionModeType& selectionMode() const;
//- Return const access to the name of cell set for "cellSet"
// selectionMode
inline const word& cellSetName() const;
//- Return const access to the total cell volume //- Return const access to the total cell volume
inline scalar V() const; inline scalar V() const;
//- Return const access to the cell set
inline const labelList& cells() const;
// Mesh changes // Mesh changes

View File

@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org \\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
@ -25,29 +25,10 @@ License
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
inline const Foam::fvCellSet::selectionModeType&
Foam::fvCellSet::selectionMode() const
{
return selectionMode_;
}
inline const Foam::word& Foam::fvCellSet::cellSetName() const
{
return cellSetName_;
}
inline Foam::scalar Foam::fvCellSet::V() const inline Foam::scalar Foam::fvCellSet::V() const
{ {
return V_; return V_;
} }
inline const Foam::labelList& Foam::fvCellSet::cells() const
{
return cells_;
}
// ************************************************************************* // // ************************************************************************* //

View File

@ -138,8 +138,6 @@ void Foam::fv::actuationDiskSource::addSup
vectorField& Usource = eqn.source(); vectorField& Usource = eqn.source();
const vectorField& U = eqn.psi(); const vectorField& U = eqn.psi();
if (set_.V() > vSmall)
{
addActuationDiskAxialInertialResistance addActuationDiskAxialInertialResistance
( (
Usource, Usource,
@ -149,7 +147,6 @@ void Foam::fv::actuationDiskSource::addSup
geometricOneField(), geometricOneField(),
U U
); );
}
} }
@ -164,8 +161,6 @@ void Foam::fv::actuationDiskSource::addSup
vectorField& Usource = eqn.source(); vectorField& Usource = eqn.source();
const vectorField& U = eqn.psi(); const vectorField& U = eqn.psi();
if (set_.V() > vSmall)
{
addActuationDiskAxialInertialResistance addActuationDiskAxialInertialResistance
( (
Usource, Usource,
@ -175,7 +170,6 @@ void Foam::fv::actuationDiskSource::addSup
rho, rho,
U U
); );
}
} }
@ -191,8 +185,6 @@ void Foam::fv::actuationDiskSource::addSup
vectorField& Usource = eqn.source(); vectorField& Usource = eqn.source();
const vectorField& U = eqn.psi(); const vectorField& U = eqn.psi();
if (set_.V() > vSmall)
{
addActuationDiskAxialInertialResistance addActuationDiskAxialInertialResistance
( (
Usource, Usource,
@ -202,7 +194,6 @@ void Foam::fv::actuationDiskSource::addSup
rho, rho,
U U
); );
}
} }

View File

@ -289,7 +289,7 @@ void Foam::fv::effectivenessHeatExchangerSource::addSup
} }
reduce(sumWeight, sumOp<scalar>()); reduce(sumWeight, sumOp<scalar>());
if (set_.V() > vSmall && mag(Qt) > vSmall) if (mag(Qt) > vSmall)
{ {
scalarField& heSource = eqn.source(); scalarField& heSource = eqn.source();

View File

@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org \\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
@ -80,8 +80,6 @@ void Foam::fv::radialActuationDiskSource::addSup
vectorField& Usource = eqn.source(); vectorField& Usource = eqn.source();
const vectorField& U = eqn.psi(); const vectorField& U = eqn.psi();
if (set_.V() > vSmall)
{
addRadialActuationDiskAxialInertialResistance addRadialActuationDiskAxialInertialResistance
( (
Usource, Usource,
@ -90,7 +88,6 @@ void Foam::fv::radialActuationDiskSource::addSup
geometricOneField(), geometricOneField(),
U U
); );
}
} }
@ -105,8 +102,6 @@ void Foam::fv::radialActuationDiskSource::addSup
vectorField& Usource = eqn.source(); vectorField& Usource = eqn.source();
const vectorField& U = eqn.psi(); const vectorField& U = eqn.psi();
if (set_.V() > vSmall)
{
addRadialActuationDiskAxialInertialResistance addRadialActuationDiskAxialInertialResistance
( (
Usource, Usource,
@ -115,7 +110,6 @@ void Foam::fv::radialActuationDiskSource::addSup
rho, rho,
U U
); );
}
} }

View File

@ -182,6 +182,8 @@ $(cellZoneSources)/setToCellZone/setToCellZone.C
pointZoneSources = sets/pointZoneSources pointZoneSources = sets/pointZoneSources
$(pointZoneSources)/setToPointZone/setToPointZone.C $(pointZoneSources)/setToPointZone/setToPointZone.C
sets/polyCellSet/polyCellSet.C
momentOfInertia/momentOfInertia.C momentOfInertia/momentOfInertia.C
surfaceSets/surfaceSets.C surfaceSets/surfaceSets.C

View File

@ -0,0 +1,242 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "polyCellSet.H"
#include "polyMesh.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
template<> const char* NamedEnum
<
polyCellSet::selectionModeType,
4
>::names[] =
{
"points",
"cellSet",
"cellZone",
"all"
};
const NamedEnum<polyCellSet::selectionModeType, 4>
polyCellSet::selectionModeTypeNames_;
}
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::polyCellSet::setCells()
{
Info<< incrIndent;
switch (selectionMode_)
{
case selectionModeType::points:
{
Info<< indent << "- selecting cells using points" << endl;
labelHashSet selectedCells;
forAll(points_, i)
{
label celli = mesh_.findCell(points_[i]);
if (celli >= 0)
{
selectedCells.insert(celli);
}
label globalCelli = returnReduce(celli, maxOp<label>());
if (globalCelli < 0)
{
WarningInFunction
<< "Unable to find owner cell for point " << points_[i]
<< endl;
}
}
cells_ = selectedCells.toc();
break;
}
case selectionModeType::cellSet:
{
Info<< indent
<< "- selecting cells using cellSet " << cellSetName_ << endl;
cells_ = cellSet(mesh_, cellSetName_).toc();
break;
}
case selectionModeType::cellZone:
{
Info<< indent
<< "- selecting cells using cellZone " << cellSetName_ << endl;
label zoneID = mesh_.cellZones().findZoneID(cellSetName_);
if (zoneID == -1)
{
FatalErrorInFunction
<< "Cannot find cellZone " << cellSetName_ << endl
<< "Valid cellZones are " << mesh_.cellZones().names()
<< exit(FatalError);
}
cells_ = mesh_.cellZones()[zoneID];
break;
}
case selectionModeType::all:
{
Info<< indent << "- selecting all cells" << endl;
cells_ = identity(mesh_.nCells());
break;
}
}
Info<< decrIndent;
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::polyCellSet::polyCellSet
(
const polyMesh& mesh,
const dictionary& dict
)
:
mesh_(mesh),
selectionMode_(selectionModeType::all),
cellSetName_(word::null)
{
read(dict);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::polyCellSet::~polyCellSet()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::polyCellSet::movePoints()
{
if (selectionMode_ == selectionModeType::points)
{
setCells();
}
}
void Foam::polyCellSet::topoChange(const polyTopoChangeMap&)
{
setCells();
}
void Foam::polyCellSet::mapMesh(const polyMeshMap&)
{
setCells();
}
void Foam::polyCellSet::distribute(const polyDistributionMap&)
{
setCells();
}
bool Foam::polyCellSet::read(const dictionary& dict)
{
if (dict.found("selectionMode"))
{
selectionMode_ =
selectionModeTypeNames_.read(dict.lookup("selectionMode"));
}
else if (dict.found("points"))
{
selectionMode_ = selectionModeType::points;
}
else if (dict.found("cellSet"))
{
selectionMode_ = selectionModeType::cellSet;
}
else if (dict.found("cellZone"))
{
selectionMode_ = selectionModeType::cellZone;
}
else
{
FatalIOErrorInFunction(dict)
<< "selectionMode, cellZone, cellSet or points not specified."
<< nl << "Valid selectionModes:"
<< selectionModeTypeNames_.sortedToc()
<< exit(FatalIOError);
}
switch (selectionMode_)
{
case selectionModeType::points:
{
dict.lookup("points") >> points_;
break;
}
case selectionModeType::cellSet:
{
dict.lookup("cellSet") >> cellSetName_;
break;
}
case selectionModeType::cellZone:
{
dict.lookup("cellZone") >> cellSetName_;
break;
}
case selectionModeType::all:
{
break;
}
default:
{
FatalErrorInFunction
<< "Unknown selectionMode "
<< selectionModeTypeNames_[selectionMode_]
<< nl << "Valid selectionModes:"
<< selectionModeTypeNames_.toc()
<< exit(FatalError);
}
}
setCells();
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,193 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::polyCellSet
Description
General cell set selection class for models that apply to sub-sets
of the mesh.
Currently supports cell selection from a set of points, a specified cellSet
or cellZone or all of the cells. The selection method can either be
specified explicitly using the \c selectionMode entry or inferred from the
presence of either a \c cellSet, \c cellZone or \c points entry. The \c
selectionMode entry is required to select \c all cells.
Usage
Examples:
\verbatim
// Apply everywhere
selectionMode all;
// Apply within a given cellSet
selectionMode cellSet; // Optional
cellSet rotor;
// Apply within a given cellZone
selectionMode cellZone; // Optional
cellSet rotor;
// Apply in cells containing a list of points
selectionMode points; // Optional
points
(
(2.25 0.5 0)
(2.75 0.5 0)
);
\endverbatim
SourceFiles
polyCellSet.C
\*---------------------------------------------------------------------------*/
#ifndef polyCellSet_H
#define polyCellSet_H
#include "cellSet.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class polyMesh;
class polyMeshMap;
class polyDistributionMap;
/*---------------------------------------------------------------------------*\
Class polyCellSet Declaration
\*---------------------------------------------------------------------------*/
class polyCellSet
{
public:
// Public data
//- Enumeration for selection mode types
enum class selectionModeType
{
points,
cellSet,
cellZone,
all
};
//- Word list of selection mode type names
static const NamedEnum<selectionModeType, 4>
selectionModeTypeNames_;
private:
// Private data
const polyMesh& mesh_;
//- Cell selection mode
selectionModeType selectionMode_;
//- Name of cell set for "cellSet" and "cellZone" selectionMode
word cellSetName_;
//- List of points for "points" selectionMode
List<point> points_;
//- Set of cells to apply source to
mutable labelList cells_;
// Private functions
//- Set the cells
void setCells();
public:
// Constructors
//- Construct from mesh and dictionary
polyCellSet
(
const polyMesh& mesh,
const dictionary& dict
);
//- Destructor
~polyCellSet();
// Member Functions
// Access
//- Return const access to the cell selection mode
inline const selectionModeType& selectionMode() const;
//- Return const access to the name of cell set for "cellSet"
// selectionMode
inline const word& cellSetName() const;
//- Return const access to the cell set
inline const labelList& cells() const;
// Mesh changes
//- Update for mesh motion
void movePoints();
//- Update topology using the given map
void topoChange(const polyTopoChangeMap&);
//- Update from another mesh using the given map
void mapMesh(const polyMeshMap&);
//- Redistribute or update using the given distribution map
void distribute(const polyDistributionMap&);
// IO
//- Read coefficients dictionary
bool read(const dictionary& dict);
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#include "polyCellSetI.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,47 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
inline const Foam::polyCellSet::selectionModeType&
Foam::polyCellSet::selectionMode() const
{
return selectionMode_;
}
inline const Foam::word& Foam::polyCellSet::cellSetName() const
{
return cellSetName_;
}
inline const Foam::labelList& Foam::polyCellSet::cells() const
{
return cells_;
}
// ************************************************************************* //

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction rotatingMotion; solidBodyMotionFunction rotatingMotion;
origin (0 0 0); origin (0 0 0);

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction SDA; solidBodyMotionFunction SDA;
CofG (0 0 0); CofG (0 0 0);

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction multiMotion; solidBodyMotionFunction multiMotion;
oscillation oscillation

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction SDA; solidBodyMotionFunction SDA;
CofG (0 0 0); CofG (0 0 0);

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction SDA; solidBodyMotionFunction SDA;
CofG (0 0 0); CofG (0 0 0);

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction SDA; solidBodyMotionFunction SDA;
CofG (0 0 0); CofG (0 0 0);

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction SDA; solidBodyMotionFunction SDA;
CofG (0 0 0); CofG (0 0 0);

View File

@ -22,6 +22,8 @@ mover
motionSolver solidBody; motionSolver solidBody;
selectionMode all;
solidBodyMotionFunction sixDoFMotion; solidBodyMotionFunction sixDoFMotion;
CofG (0 0 0); CofG (0 0 0);

View File

@ -24,6 +24,8 @@ mover
solidBodyMotionFunction multiMotion; solidBodyMotionFunction multiMotion;
selectionMode all;
// Table rotating in z axis // Table rotating in z axis
rotatingTable rotatingTable
{ {