38 template<
class CloudType>
50 this->owner().
name() +
":alpha",
51 this->owner().db().time().
name(),
58 zeroGradientFvPatchScalarField::typeName
62 applyLimiting_(this->coeffDict().lookup(
"applyLimiting")),
63 applyGravity_(this->coeffDict().lookup(
"applyGravity")),
64 alphaMin_(this->coeffDict().template lookup<scalar>(
"alphaMin")),
65 rhoMin_(this->coeffDict().template lookup<scalar>(
"rhoMin"))
67 alpha_ = this->
owner().alpha();
71 template<
class CloudType>
81 cm.phiCorrect_.
valid()
82 ? cm.phiCorrect_().
clone()
88 ? cm.uCorrect_().
clone()
91 applyLimiting_(cm.applyLimiting_),
92 applyGravity_(cm.applyGravity_),
93 alphaMin_(cm.alphaMin_),
100 template<
class CloudType>
107 template<
class CloudType>
144 alpha_ =
max(this->owner().
alpha(), alphaMin_);
145 alpha_.correctBoundaryConditions();
153 this->owner().db().time().
name(),
163 rho.correctBoundaryConditions();
175 this->owner().db().time().
name(),
186 this->particleStressModel_->dTaudTheta
188 alpha_.primitiveField(),
189 rho.primitiveField(),
229 alphaEqn +=
fvm::div(phiGByA(), alpha_);
253 this->owner().db().time().
name(),
263 U.correctBoundaryConditions();
273 phiCorrect_.ref() -= deltaT*phiGByA();
276 forAll(phiCorrect_(), facei)
279 const scalar phiCurr = phi[facei];
280 scalar& phiCorr = phiCorrect_.ref()[facei];
284 if (phiCurr*phiCorr < 0)
290 else if (phiCorr > 0)
292 phiCorr =
max(phiCorr - phiCurr, 0);
296 phiCorr =
min(phiCorr - phiCurr, 0);
302 phiCorrect_.ref() += deltaT*phiGByA();
312 uCorrect_->correctBoundaryConditions();
328 template<
class CloudType>
338 const label celli =
p.cell();
339 const label facei =
p.tetFace();
342 const vector U = uCorrect_()[celli];
346 const scalar nMag =
mag(nHat);
354 phi = phiCorrect_()[facei];
359 phiCorrect_().boundaryField()[
patchi]
366 const scalar t =
p.coordinates()[0];
370 return U + (1.0 - t)*nHat*(phi/nMag - (
U & nHat));
#define forAll(list, i)
Loop across all elements in list.
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.
Templated base class for dsmc cloud.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef()
Return a reference to the primitive field.
void correctBoundaryConditions()
Correct boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Base class for packing models.
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 any number of values (e....
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Mesh data needed to do the Finite Volume discretisation.
const fvSchemes & schemes() const
Return the fvSchemes.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const polyMesh & mesh() const
Return reference to polyMesh.
void setFluxRequired(const word &name) const
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
virtual const faceList & faces() const
Return raw faces.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
virtual const pointField & points() const
Return raw points.
virtual void cacheFields(const bool store)
Cache dependent sub-model fields.
A class for managing temporary objects.
A class for handling words, derived from string.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Reconstruct volField from a face flux field.
Calculate the matrix for the second-order temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
bool valid(const PtrList< ModelType > &l)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &vf)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dimensionSet dimPressure
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimless
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimDensity
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const word cloudName(propsDict.lookup("cloudName"))