40template<
class CloudType>
66 alphaMin_(this->
coeffDict().getScalar(
"alphaMin")),
67 rhoMin_(this->
coeffDict().getScalar(
"rhoMin"))
69 alpha_ = this->owner().theta();
74template<
class CloudType>
82 phiCorrect_(cm.phiCorrect_()),
83 uCorrect_(cm.uCorrect_()),
84 applyLimiting_(cm.applyLimiting_),
85 applyGravity_(cm.applyGravity_),
86 alphaMin_(cm.alphaMin_),
95template<
class CloudType>
102template<
class CloudType>
116 const auto& rhoAverage =
121 const auto& uAverage =
126 const auto& uSqrAverage =
132 mesh.setFluxRequired(alpha_.name());
138 alpha_ =
max(this->owner().theta(), alphaMin_);
139 alpha_.correctBoundaryConditions();
151 rho.primitiveFieldRef() =
max(rhoAverage.primitiveField(), rhoMin_);
152 rho.correctBoundaryConditions();
165 auto& tauPrime = ttauPrime.ref();
167 tauPrime.primitiveFieldRef() =
168 this->particleStressModel_->dTaudTheta
170 alpha_.primitiveField(),
171 rho.primitiveField(),
175 tauPrime.correctBoundaryConditions();
212 alphaEqn +=
fvm::div(phiGByA(), alpha_);
243 U.correctBoundaryConditions();
253 phiCorrect_.ref() -= phiGByA();
256 forAll(phiCorrect_(), facei)
259 const scalar phiCurr =
phi[facei];
260 scalar& phiCorr = phiCorrect_.ref()[facei];
264 if (phiCurr*phiCorr < 0)
270 else if (phiCorr > 0)
272 phiCorr =
max(phiCorr - phiCurr, 0);
276 phiCorr =
min(phiCorr - phiCurr, 0);
282 phiCorrect_.ref() += phiGByA();
293 uCorrect_->correctBoundaryConditions();
310template<
class CloudType>
320 const label celli =
p.cell();
321 const label facei =
p.tetFace();
324 const vector U = uCorrect_()[celli];
328 const scalar nMag =
mag(nHat);
333 const label patchi =
mesh.boundaryMesh().whichPatch(facei);
336 phi = phiCorrect_()[facei];
341 phiCorrect_().boundaryField()[patchi]
343 mesh.boundaryMesh()[patchi].whichFace(facei)
348 const scalar t =
p.coordinates()[0];
352 return U + (1.0 - t)*nHat*(
phi/nMag - (
U & nHat));
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
Base class for lagrangian averaging methods.
virtual tmp< Field< Type > > primitiveField() const =0
Return an internal field of the average.
const CloudType & owner() const
Return const access to the owner cloud.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Base class for packing models.
autoPtr< ParticleStressModel > particleStressModel_
Protected data.
PackingModel(CloudType &owner)
Construct null from owner.
Implicit model for applying an inter-particle stress to the particles.
Implicit(const dictionary &dict, CloudType &owner)
Construct from components.
virtual ~Implicit()
Destructor.
virtual void cacheFields(const bool store)
Calculate the inter particles stresses.
virtual vector velocityCorrection(typename CloudType::parcelType &p, const scalar deltaT) const
Calculate the velocity correction.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux() const
Return the face-flux field from the matrix.
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
Template invariant parts for fvPatchField.
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Lookup type of boundary radiation properties.
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
Reconstruct volField from a face flux field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
const dimensionSet dimPressure
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
DSMCCloud< dsmcParcel > CloudType
const dimensionSet dimless
Dimensionless.
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< volScalarField > trho
#define forAll(list, i)
Loop across all elements in list.