Go to the documentation of this file.
34 namespace porosityModels
53 forAll(cellZoneIDs_, zoneI)
63 const label j = fieldIndex(i);
64 const tensor Cd =
rho*(alphaZones[j] + betaZones[j]*
mag(
U[cellI]));
65 const scalar isoCd =
tr(Cd);
67 Udiag[cellI] += V[cellI]*isoCd;
68 Usource[cellI] -= V[cellI]*((Cd -
I*isoCd) &
U[cellI]);
82 forAll(cellZoneIDs_, zoneI)
92 const label j = fieldIndex(i);
107 const word& modelType,
110 const word& cellZoneName
116 alpha_(cellZoneIDs_.size()),
117 beta_(cellZoneIDs_.size())
119 adjustNegativeResistance(alphaXYZ_);
120 adjustNegativeResistance(betaXYZ_);
122 calcTransformModelData();
136 if (coordSys_.R().uniform())
138 forAll (cellZoneIDs_, zoneI)
140 alpha_[zoneI].setSize(1);
141 beta_[zoneI].setSize(1);
144 alpha_[zoneI][0].
xx() = alphaXYZ_.value().x();
145 alpha_[zoneI][0].yy() = alphaXYZ_.value().y();
146 alpha_[zoneI][0].zz() = alphaXYZ_.value().z();
147 alpha_[zoneI][0] = coordSys_.R().transformTensor(alpha_[zoneI][0]);
150 beta_[zoneI][0].
xx() = betaXYZ_.value().x();
151 beta_[zoneI][0].yy() = betaXYZ_.value().y();
152 beta_[zoneI][0].zz() = betaXYZ_.value().z();
153 beta_[zoneI][0] = coordSys_.R().transformTensor(beta_[zoneI][0]);
158 forAll(cellZoneIDs_, zoneI)
162 alpha_[zoneI].setSize(
cells.size());
163 beta_[zoneI].setSize(
cells.size());
168 alpha_[zoneI][i].
xx() = alphaXYZ_.value().x();
169 alpha_[zoneI][i].yy() = alphaXYZ_.value().y();
170 alpha_[zoneI][i].zz() = alphaXYZ_.value().z();
173 beta_[zoneI][i].
xx() = betaXYZ_.value().x();
174 beta_[zoneI][i].yy() = betaXYZ_.value().y();
175 beta_[zoneI][i].zz() = betaXYZ_.value().z();
180 alpha_[zoneI] =
R.transformTensor(alpha_[zoneI],
cells);
181 beta_[zoneI] =
R.transformTensor(beta_[zoneI],
cells);
198 scalar rhoRef =
readScalar(coeffs_.lookup(
"rhoRef"));
200 apply(Udiag, Usource, V,
U, rhoRef);
202 force = Udiag*
U - Usource;
219 coeffs_.lookup(
"rhoRef") >>
rho;
222 apply(Udiag, Usource, V,
U,
rho);
241 coeffs_.lookup(
"rhoRef") >>
rho;
244 apply(Udiag, Usource, V,
U,
rho);
259 coeffs_.lookup(
"rhoRef") >>
rho;
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Templated 3D tensor derived from VectorSpace adding construction from 9 components,...
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
A class for handling words, derived from string.
const dimensionedScalar mu
Atomic mass unit.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
virtual void calcForce(const volVectorField &U, const volScalarField &rho, const volScalarField &mu, vectorField &force) const
Calculate the porosity force.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
fixedCoeff(const fixedCoeff &)
Disallow default bitwise copy construct.
defineTypeNameAndDebug(DarcyForchheimer, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
const dimensionSet dimForce
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
#define R(A, B, C, D, E, F, K, M)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool writeData(Ostream &os) const
Write.
Pre-declare SubField and related Field type.
virtual Ostream & write(const token &)=0
Write next token to stream.
Abstract base class for coordinate rotation.
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
static const sphericalTensor I(1)
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
A list of keyword definitions, which are a keyword followed by any number of values (e....
void apply(scalarField &Udiag, vectorField &Usource, const scalarField &V, const vectorField &U, const scalar rho) const
Apply.
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
Ostream & indent(Ostream &os)
Indent stream.
virtual ~fixedCoeff()
Destructor.
Fixed coefficient form of porosity model.
Top level model for porosity models.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Generic GeometricField class.
addToRunTimeSelectionTable(porosityModel, DarcyForchheimer, mesh)
word name(const complex &)
Return a string representation of a complex.