34 namespace porosityModels
44 void Foam::porosityModels::fixedCoeff::apply
60 const scalar isoCd =
tr(Cd);
62 Udiag[celli] += V[celli]*isoCd;
63 Usource[celli] -= V[celli]*((Cd -
I*isoCd) &
U[celli]);
68 void Foam::porosityModels::fixedCoeff::apply
80 const label j = fieldIndex(i);
82 const tensor beta = beta_[j];
97 const word& cellZoneName
103 rhoRefFound_(coeffDict.
found(
"rhoRef")),
104 rhoRef_(coeffDict.lookupOrDefault<scalar>(
"rhoRef", 1.0))
123 if (coordSys_.R().uniform())
129 alpha_[0].xx() = alphaXYZ_.value().x();
130 alpha_[0].yy() = alphaXYZ_.value().y();
131 alpha_[0].zz() = alphaXYZ_.value().z();
132 alpha_[0] = coordSys_.R().transform(
Zero, alpha_[0]);
135 beta_[0].xx() = betaXYZ_.value().x();
136 beta_[0].yy() = betaXYZ_.value().y();
137 beta_[0].zz() = betaXYZ_.value().z();
138 beta_[0] = coordSys_.R().transform(
Zero, beta_[0]);
144 alpha_.setSize(
cells.size());
145 beta_.setSize(
cells.size());
150 alpha_[i].xx() = alphaXYZ_.value().x();
151 alpha_[i].yy() = alphaXYZ_.value().y();
152 alpha_[i].zz() = alphaXYZ_.value().z();
155 beta_[i].xx() = betaXYZ_.value().x();
156 beta_[i].yy() = betaXYZ_.value().y();
157 beta_[i].zz() = betaXYZ_.value().z();
165 alpha_ =
R.transform(alpha_);
166 beta_ =
R.transform(beta_);
189 apply(Udiag, Usource, V,
U, rhoRef_);
191 force = Udiag*
U - Usource;
211 apply(Udiag, Usource, V,
U, rhoRef_);
229 apply(AU,
U, rhoRef_);
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A List with indirect addressing.
Abstract base class for coordinate rotation.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
const cellZoneList & cellZones() const
Return cell zones.
Top level model for porosity models.
const fvMesh & mesh_
Reference to the mesh database.
word zoneName_
Name of cellZone.
void adjustNegativeResistance(dimensionedVector &resist)
Adjust negative resistance values to be multiplier of max value.
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.
fixedCoeff(const word &name, const fvMesh &mesh, const dictionary &dict, const dictionary &coeffDict, const word &cellZoneName)
virtual ~fixedCoeff()
Destructor.
virtual void correct(fvVectorMatrix &UEqn) const
Add resistance.
virtual void calcTransformModelData()
Transform the model data wrt mesh changes.
A class for handling words, derived from string.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
volScalarField scalarField(fieldObject, mesh)
volVectorField vectorField(fieldObject, mesh)
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar mu
Atomic mass unit.
addToRunTimeSelectionTable(porosityModel, powerLawLopesdaCosta, mesh)
defineTypeNameAndDebug(powerLawLopesdaCosta, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Tensor< scalar > tensor
Tensor of scalars.
const dimensionSet dimless
const dimensionSet dimLength
static const Identity< scalar > I
const dimensionSet dimForce
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
Field< vector > vectorField
Specialisation of Field<T> for vector.
static scalar R(const scalar a, const scalar x)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.