48void Foam::porosityModels::fixedCoeff::apply
66 const label celli =
cells[i];
68 const tensor Cd =
rho*(alphaZones[j] + betaZones[j]*
mag(
U[celli]));
69 const scalar isoCd =
tr(Cd);
71 Udiag[celli] += V[celli]*isoCd;
72 Usource[celli] -= V[celli]*((Cd -
I*isoCd) &
U[celli]);
78void Foam::porosityModels::fixedCoeff::apply
85 forAll(cellZoneIDs_, zoneI)
94 const label celli =
cells[i];
95 const label j = fieldIndex(i);
107Foam::porosityModels::fixedCoeff::fixedCoeff
110 const word& modelType,
113 const wordRe& cellZoneName
116 porosityModel(
name, modelType,
mesh,
dict, cellZoneName),
119 alpha_(cellZoneIDs_.size()),
120 beta_(cellZoneIDs_.size())
135 alphaCoeff.xx() = alphaXYZ_.value().x();
136 alphaCoeff.yy() = alphaXYZ_.value().y();
137 alphaCoeff.zz() = alphaXYZ_.value().z();
141 betaCoeff.xx() = betaXYZ_.value().x();
142 betaCoeff.yy() = betaXYZ_.value().y();
143 betaCoeff.zz() = betaXYZ_.value().z();
147 forAll(cellZoneIDs_, zonei)
149 alpha_[zonei].resize(1);
150 beta_[zonei].resize(1);
152 alpha_[zonei] = csys().transform(alphaCoeff);
153 beta_[zonei] = csys().transform(betaCoeff);
158 forAll(cellZoneIDs_, zonei)
163 mesh_.cellZones()[cellZoneIDs_[zonei]]
166 alpha_[zonei] =
csys().transform(cc, alphaCoeff);
167 beta_[zonei] =
csys().transform(cc, betaCoeff);
184 const scalar rhoRef = coeffs_.get<scalar>(
"rhoRef");
186 apply(Udiag, Usource, V,
U, rhoRef);
188 force = Udiag*
U - Usource;
254 dict_.writeEntry(name_,
os);
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const Cmpt & yy() const noexcept
const Cmpt & xx() const noexcept
const Cmpt & zz() const noexcept
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Top level model for porosity models.
const fvMesh & mesh_
Reference to the mesh database.
dictionary coeffs_
Model coefficients dictionary.
porosityModel(const porosityModel &)=delete
No copy construct.
const dictionary dict_
Dictionary used for model construction.
const coordinateSystem & csys() const
Local coordinate system.
const dictionary & dict() const
Return dictionary used for model construction.
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
const word & name() const
Return const access to the porosity model name.
labelList cellZoneIDs_
Cell zone IDs.
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu)
Return the force over the cell zone(s).
label fieldIndex(const label index) const
Return label index.
Fixed coefficient form of porosity model.
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
bool writeData(Ostream &os) const
Write.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
Tensor of scalars, i.e. Tensor<scalar>.
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
OBJstream os(runTime.globalPath()/outputName)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
List< label > labelList
A List of labels.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static const Identity< scalar > I
UIndirectList< point > pointUIndList
UIndirectList of point.
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fvMatrix< vector > fvVectorMatrix
const dimensionSet dimForce
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.