41 template<
class BasePhaseModel>
45 word phiName(IOobject::groupName(
"phi", this->
name()));
47 typeIOobject<surfaceScalarField> phiHeader
50 U.mesh().time().name(),
55 if (phiHeader.headerOk())
57 Info<<
"Reading face flux field " << phiName <<
endl;
59 return tmp<surfaceScalarField>
66 U.mesh().time().name(),
77 Info<<
"Calculating face flux field " << phiName <<
endl;
81 U.boundaryField().size(),
82 calculatedFvPatchScalarField::typeName
87 if (!
U.boundaryField()[
patchi].assignable())
89 phiTypes[
patchi] = fixedValueFvPatchScalarField::typeName;
93 return tmp<surfaceScalarField>
100 U.mesh().time().name(),
115 template<
class BasePhaseModel>
119 const word& phaseName,
120 const bool referencePhase,
124 BasePhaseModel(fluid, phaseName, referencePhase, index),
130 fluid.mesh().time().
name(),
143 fluid.mesh().time().
name(),
156 fluid.mesh().time().
name(),
183 fluid.mesh().time().
name(),
193 if (fluid.
mesh().
dynamic() || this->fluid().MRF().size())
215 template<
class BasePhaseModel>
222 template<
class BasePhaseModel>
233 template<
class BasePhaseModel>
240 template<
class BasePhaseModel>
243 BasePhaseModel::correctKinematics();
247 K_.ref() = 0.5*
magSqr(this->
U());
252 template<
class BasePhaseModel>
255 BasePhaseModel::predictMomentumTransport();
256 momentumTransport_->predict();
260 template<
class BasePhaseModel>
263 BasePhaseModel::correctMomentumTransport();
264 momentumTransport_->correct();
268 template<
class BasePhaseModel>
287 template<
class BasePhaseModel>
294 template<
class BasePhaseModel>
305 +
fvm::SuSp(-this->continuityError(), U_)
307 + momentumTransport_->divDevTau(U_)
312 template<
class BasePhaseModel>
326 + momentumTransport_->divDevTau(U_)
331 template<
class BasePhaseModel>
339 template<
class BasePhaseModel>
347 template<
class BasePhaseModel>
355 template<
class BasePhaseModel>
363 template<
class BasePhaseModel>
371 template<
class BasePhaseModel>
379 template<
class BasePhaseModel>
387 template<
class BasePhaseModel>
398 <<
"Uf has not been allocated."
406 template<
class BasePhaseModel>
417 <<
"Uf has not been allocated."
425 template<
class BasePhaseModel>
433 template<
class BasePhaseModel>
441 template<
class BasePhaseModel>
449 template<
class BasePhaseModel>
457 template<
class BasePhaseModel>
465 template<
class BasePhaseModel>
473 template<
class BasePhaseModel>
482 + this->fluid().MRF().DDt(U_);
486 template<
class BasePhaseModel>
494 template<
class BasePhaseModel>
498 return continuityError_;
502 template<
class BasePhaseModel>
520 template<
class BasePhaseModel>
528 template<
class BasePhaseModel>
544 template<
class BasePhaseModel>
548 return momentumTransport_->k();
552 template<
class BasePhaseModel>
556 return momentumTransport_->pPrimef();
#define forAll(list, i)
Loop across all elements in list.
Generic GeometricField class.
static const GeometricField< Type, PatchField, GeoMesh > & null()
Return a null geometric field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
writeOption writeOpt() const
static word groupName(Name name, const word &group)
virtual void correctUf()
Correct the face velocity for moving meshes.
virtual void correctKinematics()
Correct the kinematics.
virtual tmp< fvVectorMatrix > UEqn()
Return the momentum equation.
volVectorField U_
Velocity field.
virtual void correct()
Correct the phase properties other than the thermo.
virtual tmp< volScalarField > k() const
Return the turbulent kinetic energy.
virtual tmp< fvVectorMatrix > DUDt() const
Return the substantive acceleration matrix.
virtual ~MovingPhaseModel()
Destructor.
virtual tmp< volVectorField > U() const
Return the velocity.
virtual void correctMomentumTransport()
Correct the momentumTransport.
virtual surfaceScalarField & alphaRhoPhiRef()
Access the mass flux of the phase.
virtual const autoPtr< volScalarField > & divU() const
Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual surfaceScalarField & phiRef()
Access the volumetric flux.
virtual tmp< surfaceScalarField > alphaPhi() const
Return the volumetric flux of the phase.
virtual tmp< volScalarField > K() const
Return the phase kinetic energy.
surfaceScalarField phi_
Flux.
MovingPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual const autoPtr< surfaceVectorField > & Uf() const
Return the face velocity.
virtual surfaceVectorField & UfRef()
Access the face velocity.
virtual surfaceScalarField & alphaPhiRef()
Access the volumetric flux of the phase.
virtual tmp< surfaceScalarField > alphaRhoPhi() const
Return the mass flux of the phase.
virtual tmp< surfaceScalarField > phi() const
Return the volumetric flux.
virtual tmp< fvVectorMatrix > UfEqn()
Return the momentum equation for the face-based algorithm.
virtual volVectorField & URef()
Access the velocity.
virtual bool stationary() const
Return whether the phase is stationary.
virtual void predictMomentumTransport()
Predict the momentumTransport.
virtual tmp< fvVectorMatrix > UgradU() const
Return the velocity transport matrix.
virtual tmp< volScalarField > continuityError() const
Return the continuity error.
virtual void correctContinuityError(const volScalarField &source)
Correct the continuity error.
autoPtr< surfaceVectorField > Uf_
Face velocity field.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Dimension set for the base types.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
bool dynamic() const
Is this mesh dynamic?
const polyMesh & mesh() const
Return reference to polyMesh.
Abstract base class for turbulence models (RAS, LES and laminar).
Class to represent a system of phases and model interfacial transfers between them.
const fvMesh & mesh() const
Return the mesh.
A class for managing temporary objects.
T * ptr() const
Return tmp pointer for reuse.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the face-flux of the given 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 implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
tmp< VolField< Type > > DDt(const surfaceScalarField &phi, const VolField< Type > &psi)
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
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< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< word > wordList
A List of words.
VolField< vector > volVectorField
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimTime
const dimensionSet dimDensity
VolField< scalar > volScalarField
SurfaceField< vector > surfaceVectorField
dimensioned< scalar > magSqr(const dimensioned< Type > &)