Files
OpenFOAM-12/applications/modules/multiphaseEuler/phaseSystem/phaseModels/StationaryPhaseModel/StationaryPhaseModel.H
Henry Weller 5e03874bbb multiphaseEuler: Updated to us the new phaseSolidThermophysicalTransportModel class
for thermophysical transport within stationary solid phases.  This provides a
consistent interface to heat transport within solids for single and now
multiphase solvers so that for example the wallHeatFlux functionObject can now
be used with multiphaseEuler, see tutorials/multiphaseEuler/boilingBed.
Also this development supports anisotropic thermal conductivity within the
stationary solid regions which was not possible previously.

The tutorials/multiphaseEuler/bed and tutorials/multiphaseEuler/boilingBed
tutorial cases have been updated for phaseSolidThermophysicalTransportModel by
changing the thermo type in physicalProperties.solid to heSolidThermo.  This
change will need to be made to all multiphaseEuler cases involving stationary
phases.
2023-10-11 14:53:09 +01:00

187 lines
5.9 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2015-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::StationaryPhaseModel
Description
Class which represents a stationary (and therefore probably solid) phase.
Generates, but does not store, zero velocity and flux field and turbulent
quantities. Throws an error when non-const access is requested to the
motion fields or when the momentum equation is requested. Usage must
must protect against such calls.
See also
MovingPhaseModel
SourceFiles
StationaryPhaseModel.C
\*---------------------------------------------------------------------------*/
#ifndef StationaryPhaseModel_H
#define StationaryPhaseModel_H
#include "phaseModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class StationaryPhaseModel Declaration
\*---------------------------------------------------------------------------*/
template<class BasePhaseModel>
class StationaryPhaseModel
:
public BasePhaseModel
{
public:
// Constructors
StationaryPhaseModel
(
const phaseSystem& fluid,
const word& phaseName,
const bool referencePhase,
const label index
);
//- Destructor
virtual ~StationaryPhaseModel();
// Member Functions
// Momentum
//- Return whether the phase is stationary
virtual bool stationary() const;
//- Return the momentum equation
virtual tmp<fvVectorMatrix> UEqn();
//- Return the momentum equation for the face-based algorithm
virtual tmp<fvVectorMatrix> UfEqn();
//- Return the velocity
virtual tmp<volVectorField> U() const;
//- Access the velocity
virtual volVectorField& URef();
//- Access the velocity
virtual const volVectorField& URef() const;
//- Return the volumetric flux
virtual tmp<surfaceScalarField> phi() const;
//- Access the volumetric flux
virtual surfaceScalarField& phiRef();
//- Access the volumetric flux
virtual const surfaceScalarField& phiRef() const;
//- Return the face velocity
// Required for moving mesh cases
virtual const autoPtr<surfaceVectorField>& Uf() const;
//- Access the face velocity
// Required for moving mesh cases
virtual surfaceVectorField& UfRef();
//- Access the face velocity
// Required for moving mesh cases
virtual const surfaceVectorField& UfRef() const;
//- Return the volumetric flux of the phase
virtual tmp<surfaceScalarField> alphaPhi() const;
//- Access the volumetric flux of the phase
virtual surfaceScalarField& alphaPhiRef();
//- Access the volumetric flux of the phase
virtual const surfaceScalarField& alphaPhiRef() const;
//- Return the mass flux of the phase
virtual tmp<surfaceScalarField> alphaRhoPhi() const;
//- Access the mass flux of the phase
virtual surfaceScalarField& alphaRhoPhiRef();
//- Access the mass flux of the phase
virtual const surfaceScalarField& alphaRhoPhiRef() const;
//- Return the velocity transport matrix
virtual tmp<fvVectorMatrix> UgradU() const;
//- Return the substantive acceleration matrix
virtual tmp<fvVectorMatrix> DUDt() const;
//- Return the continuity error
virtual tmp<volScalarField> continuityError() const;
//- Return the phase kinetic energy
virtual tmp<volScalarField> K() const;
// Compressibility (variable density)
//- Return the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual const autoPtr<volScalarField>& divU() const;
//- Set the phase dilatation rate (d(alpha)/dt + div(alpha*phi))
virtual void divU(tmp<volScalarField> divU);
// Momentum transport
//- Return the turbulent kinetic energy
virtual tmp<volScalarField> k() const;
//- Return the face-phase-pressure'
// (derivative of phase-pressure w.r.t. phase-fraction)
virtual tmp<surfaceScalarField> pPrimef() const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "StationaryPhaseModel.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //