33template<
class CloudType>
43 C_(this->
coeffs().getScalar(
"C")),
44 gradInterForceInterpPtr_(nullptr)
48template<
class CloudType>
52 alphaName_(pf.alphaName_),
54 gradInterForceInterpPtr_(pf.gradInterForceInterpPtr_)
60template<
class CloudType>
67template<
class CloudType>
70 static word resultName(
"gradAlpha");
73 this->
mesh().template getObjectPtr<volVectorField>(resultName);
80 lookupObject<volScalarField>(alphaName_);
91 const volVectorField& gradInterForce = *resultPtr;
93 gradInterForceInterpPtr_.reset
95 interpolation<vector>::New
97 this->owner().solution().interpolationSchemes(),
104 gradInterForceInterpPtr_.clear();
114template<
class CloudType>
117 const typename CloudType::parcelType&
p,
118 const typename CloudType::parcelType::trackingData&
td,
125 forceSuSp value(Zero);
127 const interpolation<vector>& interp = gradInterForceInterpPtr_();
132 *interp.interpolate(
p.coordinates(),
p.currentTetIndices());
Vector force apply to particles to avoid the crossing of particles from one phase to the other....
virtual forceSuSp calcNonCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
virtual ~InterfaceForce()
Destructor.
virtual void cacheFields(const bool store)
Cache fields.
InterfaceForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict)
Construct from mesh.
Abstract base class for particle forces.
const fvMesh & mesh() const noexcept
Return the mesh database.
const CloudType & owner() const noexcept
Return const access to the cloud owner.
ParticleForce(CloudType &owner, const fvMesh &mesh, const dictionary &dict, const word &forceType, const bool readCoeffs)
Construct from mesh.
const dictionary & coeffs() const noexcept
Return the force coefficients dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Helper container for force Su and Sp terms.
const vector & Su() const
Return const access to the explicit contribution [kg.m/s2].
Mesh data needed to do the Finite Volume discretisation.
Abstract base class for volume field interpolation.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
Lookup type of boundary radiation properties.
bool store()
Register object with its registry and transfer ownership to the registry.
bool checkOut()
Remove object from registry, and remove all file watches.
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from Foam::string.
Calculate the gradient of the given field.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
scalarField Re(const UList< complex > &cmplx)
Extract real component.