The new flexible and extensible modular solvers structure already provides most
of the simulation functionality needed for single phase, multiphase,
multicomponent etc. fluid flow problems as well as a very effective method of
combining these with solid heat transfer, solid stress, surface film to solve
complex multi-region, multi-physics problems and are now the primary mechanism
for the further development of OpenFOAM simulation capability in future. To
emphasis this for both users and developers the applications/solvers directory
has been separated into applications/modules containing all the solver modules:
├── modules
│ ├── compressibleMultiphaseVoF
│ ├── compressibleVoF
│ ├── film
│ ├── fluid
│ ├── fluidSolver
│ ├── functions
│ ├── incompressibleDenseParticleFluid
│ ├── incompressibleDriftFlux
│ ├── incompressibleFluid
│ ├── incompressibleMultiphaseVoF
│ ├── incompressibleVoF
│ ├── isothermalFilm
│ ├── isothermalFluid
│ ├── movingMesh
│ ├── multicomponentFluid
│ ├── multiphaseEuler
│ ├── multiphaseVoFSolver
│ ├── shockFluid
│ ├── solid
│ ├── solidDisplacement
│ ├── twoPhaseSolver
│ ├── twoPhaseVoFSolver
│ ├── VoFSolver
│ └── XiFluid
applications/solvers containing the foamRun and foamMultiRun solver applications
which instantiate and execute the chosen solver modules and also standalone
solver applications for special initialisation and test activities:
├── solvers
│ ├── boundaryFoam
│ ├── chemFoam
│ ├── foamMultiRun
│ ├── foamRun
│ └── potentialFoam
and applications/legacy containing legacy solver applications which are not
currently being actively developed but the functionality of which will be merged
into the solver modules or form the basis of new solver modules as the need
arises:
├── legacy
│ ├── basic
│ │ ├── financialFoam
│ │ └── laplacianFoam
│ ├── combustion
│ │ └── PDRFoam
│ ├── compressible
│ │ └── rhoPorousSimpleFoam
│ ├── electromagnetics
│ │ ├── electrostaticFoam
│ │ ├── magneticFoam
│ │ └── mhdFoam
│ ├── incompressible
│ │ ├── adjointShapeOptimisationFoam
│ │ ├── dnsFoam
│ │ ├── icoFoam
│ │ ├── porousSimpleFoam
│ │ └── shallowWaterFoam
│ └── lagrangian
│ ├── dsmcFoam
│ ├── mdEquilibrationFoam
│ └── mdFoam
Correspondingly the tutorials directory structure has been reorganised with the
modular solver directories at the top level with names that make it easier for
users to find example cases relating to their particular requirements and a
legacy sub-directory containing cases corresponding to the legacy solver
applications listed above:
├── compressibleMultiphaseVoF
│ └── damBreak4phaseLaminar
├── compressibleVoF
│ ├── ballValve
│ ├── climbingRod
│ ├── damBreak
│ ├── depthCharge2D
│ ├── depthCharge3D
│ ├── sloshingTank2D
│ └── throttle
├── film
│ └── rivuletPanel
├── fluid
│ ├── aerofoilNACA0012
│ ├── aerofoilNACA0012Steady
│ ├── angledDuct
│ ├── angledDuctExplicitFixedCoeff
│ ├── angledDuctLTS
│ ├── annularThermalMixer
│ ├── BernardCells
│ ├── blockedChannel
│ ├── buoyantCavity
│ ├── cavity
│ ├── decompressionTank
│ ├── externalCoupledCavity
│ ├── forwardStep
│ ├── helmholtzResonance
│ ├── hotRadiationRoom
│ ├── hotRadiationRoomFvDOM
│ ├── hotRoom
│ ├── hotRoomBoussinesq
│ ├── hotRoomBoussinesqSteady
│ ├── hotRoomComfort
│ ├── iglooWithFridges
│ ├── mixerVessel2DMRF
│ ├── nacaAirfoil
│ ├── pitzDaily
│ ├── prism
│ ├── shockTube
│ ├── squareBend
│ ├── squareBendLiq
│ └── squareBendLiqSteady
├── incompressibleDenseParticleFluid
│ ├── column
│ ├── cyclone
│ ├── Goldschmidt
│ ├── GoldschmidtMPPIC
│ └── injectionChannel
├── incompressibleDriftFlux
│ ├── dahl
│ ├── mixerVessel2DMRF
│ └── tank3D
├── incompressibleFluid
│ ├── airFoil2D
│ ├── ballValve
│ ├── blockedChannel
│ ├── cavity
│ ├── cavityCoupledU
│ ├── channel395
│ ├── drivaerFastback
│ ├── ductSecondaryFlow
│ ├── elipsekkLOmega
│ ├── flowWithOpenBoundary
│ ├── hopperParticles
│ ├── impeller
│ ├── mixerSRF
│ ├── mixerVessel2D
│ ├── mixerVessel2DMRF
│ ├── mixerVesselHorizontal2DParticles
│ ├── motorBike
│ ├── motorBikeSteady
│ ├── movingCone
│ ├── offsetCylinder
│ ├── oscillatingInlet
│ ├── pipeCyclic
│ ├── pitzDaily
│ ├── pitzDailyLES
│ ├── pitzDailyLESDevelopedInlet
│ ├── pitzDailyLTS
│ ├── pitzDailyPulse
│ ├── pitzDailyScalarTransport
│ ├── pitzDailySteady
│ ├── pitzDailySteadyExperimentalInlet
│ ├── pitzDailySteadyMappedToPart
│ ├── pitzDailySteadyMappedToRefined
│ ├── planarContraction
│ ├── planarCouette
│ ├── planarPoiseuille
│ ├── porousBlockage
│ ├── propeller
│ ├── roomResidenceTime
│ ├── rotor2DRotating
│ ├── rotor2DSRF
│ ├── rotorDisk
│ ├── T3A
│ ├── TJunction
│ ├── TJunctionFan
│ ├── turbineSiting
│ ├── waveSubSurface
│ ├── windAroundBuildings
│ └── wingMotion
├── incompressibleMultiphaseVoF
│ ├── damBreak4phase
│ ├── damBreak4phaseFineLaminar
│ ├── damBreak4phaseLaminar
│ └── mixerVessel2DMRF
├── incompressibleVoF
│ ├── angledDuct
│ ├── capillaryRise
│ ├── cavitatingBullet
│ ├── climbingRod
│ ├── containerDischarge2D
│ ├── damBreak
│ ├── damBreakLaminar
│ ├── damBreakPorousBaffle
│ ├── damBreakWithObstacle
│ ├── DTCHull
│ ├── DTCHullMoving
│ ├── DTCHullWave
│ ├── floatingObject
│ ├── floatingObjectWaves
│ ├── forcedUpstreamWave
│ ├── mixerVessel
│ ├── mixerVessel2DMRF
│ ├── mixerVesselHorizontal2D
│ ├── nozzleFlow2D
│ ├── planingHullW3
│ ├── propeller
│ ├── sloshingCylinder
│ ├── sloshingTank2D
│ ├── sloshingTank2D3DoF
│ ├── sloshingTank3D
│ ├── sloshingTank3D3DoF
│ ├── sloshingTank3D6DoF
│ ├── testTubeMixer
│ ├── waterChannel
│ ├── wave
│ ├── wave3D
│ └── weirOverflow
├── isothermalFilm
│ └── rivuletPanel
├── isothermalFluid
│ ├── potentialFreeSurfaceMovingOscillatingBox
│ └── potentialFreeSurfaceOscillatingBox
├── legacy
│ ├── basic
│ │ ├── financialFoam
│ │ │ └── europeanCall
│ │ └── laplacianFoam
│ │ └── flange
│ ├── combustion
│ │ └── PDRFoam
│ │ └── flamePropagationWithObstacles
│ ├── compressible
│ │ └── rhoPorousSimpleFoam
│ │ ├── angledDuctExplicit
│ │ └── angledDuctImplicit
│ ├── electromagnetics
│ │ ├── electrostaticFoam
│ │ │ └── chargedWire
│ │ └── mhdFoam
│ │ └── hartmann
│ ├── incompressible
│ │ ├── adjointShapeOptimisationFoam
│ │ │ └── pitzDaily
│ │ ├── dnsFoam
│ │ │ └── boxTurb16
│ │ ├── icoFoam
│ │ │ ├── cavity
│ │ │ └── elbow
│ │ ├── porousSimpleFoam
│ │ │ ├── angledDuctExplicit
│ │ │ └── angledDuctImplicit
│ │ └── shallowWaterFoam
│ │ └── squareBump
│ ├── lagrangian
│ │ ├── dsmcFoam
│ │ │ ├── freeSpacePeriodic
│ │ │ ├── freeSpaceStream
│ │ │ ├── supersonicCorner
│ │ │ └── wedge15Ma5
│ │ ├── mdEquilibrationFoam
│ │ │ ├── periodicCubeArgon
│ │ │ └── periodicCubeWater
│ │ └── mdFoam
│ │ └── nanoNozzle
├── mesh
│ ├── blockMesh
│ │ ├── pipe
│ │ ├── sphere
│ │ ├── sphere7
│ │ └── sphere7ProjectedEdges
│ ├── refineMesh
│ │ └── refineFieldDirs
│ └── snappyHexMesh
│ ├── flange
│ └── pipe
├── movingMesh
│ └── SnakeRiverCanyon
├── multicomponentFluid
│ ├── aachenBomb
│ ├── counterFlowFlame2D
│ ├── counterFlowFlame2D_GRI
│ ├── counterFlowFlame2D_GRI_TDAC
│ ├── counterFlowFlame2DLTS
│ ├── counterFlowFlame2DLTS_GRI_TDAC
│ ├── DLR_A_LTS
│ ├── filter
│ ├── lockExchange
│ ├── membrane
│ ├── nc7h16
│ ├── parcelInBox
│ ├── SandiaD_LTS
│ ├── simplifiedSiwek
│ ├── smallPoolFire2D
│ ├── smallPoolFire3D
│ ├── verticalChannel
│ ├── verticalChannelLTS
│ └── verticalChannelSteady
├── multiphaseEuler
│ ├── bed
│ ├── bubbleColumn
│ ├── bubbleColumnEvaporating
│ ├── bubbleColumnEvaporatingDissolving
│ ├── bubbleColumnEvaporatingReacting
│ ├── bubbleColumnIATE
│ ├── bubbleColumnLaminar
│ ├── bubbleColumnLES
│ ├── bubblePipe
│ ├── damBreak4phase
│ ├── fluidisedBed
│ ├── fluidisedBedLaminar
│ ├── Grossetete
│ ├── hydrofoil
│ ├── injection
│ ├── LBend
│ ├── mixerVessel2D
│ ├── mixerVessel2DMRF
│ ├── pipeBend
│ ├── steamInjection
│ ├── titaniaSynthesis
│ ├── titaniaSynthesisSurface
│ ├── wallBoilingIATE
│ ├── wallBoilingPolydisperse
│ └── wallBoilingPolydisperseTwoGroups
├── multiRegion
│ ├── CHT
│ │ ├── circuitBoardCooling
│ │ ├── coolingCylinder2D
│ │ ├── coolingSphere
│ │ ├── heatedDuct
│ │ ├── heatExchanger
│ │ ├── multiphaseCoolingCylinder2D
│ │ ├── reverseBurner
│ │ ├── shellAndTubeHeatExchanger
│ │ ├── VoFcoolingCylinder2D
│ │ └── wallBoiling
│ └── film
│ ├── cylinder
│ ├── cylinderDripping
│ ├── cylinderVoF
│ ├── hotBoxes
│ ├── rivuletBox
│ ├── rivuletPanel
│ ├── splashPanel
│ └── VoFToFilm
├── potentialFoam
│ ├── cylinder
│ └── pitzDaily
├── resources
│ ├── blockMesh
│ ├── geometry
│ └── thermoData
├── shockFluid
│ ├── biconic25-55Run35
│ ├── forwardStep
│ ├── LadenburgJet60psi
│ ├── movingCone
│ ├── obliqueShock
│ ├── shockTube
│ └── wedge15Ma5
├── solidDisplacement
│ ├── beamEndLoad
│ └── plateHole
└── XiFluid
├── kivaTest
└── moriyoshiHomogeneous
696 lines
21 KiB
C++
696 lines
21 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::phaseSystem
|
|
|
|
Description
|
|
Class to represent a system of phases and model interfacial transfers
|
|
between them.
|
|
|
|
SourceFiles
|
|
phaseSystem.C
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#ifndef phaseSystem_H
|
|
#define phaseSystem_H
|
|
|
|
#include "IOdictionary.H"
|
|
|
|
#include "phaseModel.H"
|
|
#include "phaseInterface.H"
|
|
#include "phaseInterfaceKey.H"
|
|
#include "HashPtrTable.H"
|
|
#include "PtrListDictionary.H"
|
|
#include "hashedWordList.H"
|
|
|
|
#include "pimpleNoLoopControl.H"
|
|
|
|
#include "IOMRFZoneList.H"
|
|
#include "fvModels.H"
|
|
#include "fvConstraints.H"
|
|
|
|
#include "volFields.H"
|
|
#include "surfaceFields.H"
|
|
#include "fvMatricesFwd.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
namespace Foam
|
|
{
|
|
|
|
class interfaceSurfaceTensionModel;
|
|
class pressureReference;
|
|
class nonOrthogonalSolutionControl;
|
|
|
|
/*---------------------------------------------------------------------------*\
|
|
Class phaseSystem Declaration
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
class phaseSystem
|
|
:
|
|
public IOdictionary
|
|
{
|
|
public:
|
|
|
|
// Public Typedefs
|
|
|
|
typedef HashPtrTable<fvVectorMatrix> momentumTransferTable;
|
|
|
|
typedef HashPtrTable<fvScalarMatrix> heatTransferTable;
|
|
|
|
typedef HashPtrTable<fvScalarMatrix> specieTransferTable;
|
|
|
|
typedef PtrListDictionary<phaseModel> phaseModelList;
|
|
|
|
typedef UPtrList<phaseModel> phaseModelPartialList;
|
|
|
|
typedef
|
|
HashPtrTable
|
|
<
|
|
volScalarField,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>
|
|
dmdtfTable;
|
|
|
|
typedef
|
|
HashPtrTable
|
|
<
|
|
HashPtrTable<volScalarField>,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>
|
|
dmidtfTable;
|
|
|
|
|
|
|
|
protected:
|
|
|
|
// Protected typedefs
|
|
|
|
typedef
|
|
HashTable
|
|
<
|
|
autoPtr<interfaceSurfaceTensionModel>,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>
|
|
interfaceSurfaceTensionModelTable;
|
|
|
|
typedef HashTable<scalar, phaseInterfaceKey, phaseInterfaceKey::hash>
|
|
cAlphaTable;
|
|
|
|
|
|
// Protected data
|
|
|
|
//- Reference to the mesh
|
|
const fvMesh& mesh_;
|
|
|
|
//- Reference to pimpleNoLoopControl
|
|
const pimpleNoLoopControl& pimple_;
|
|
|
|
//- Optional MRF zones
|
|
IOMRFZoneList MRF_;
|
|
|
|
//- Name of optional reference phase which is not solved for
|
|
// but obtained from the sum of the other phases
|
|
word referencePhaseName_;
|
|
|
|
//- Phase models
|
|
phaseModelList phaseModels_;
|
|
|
|
//- Moving phase models
|
|
phaseModelPartialList movingPhaseModels_;
|
|
|
|
//- Stationary phase models
|
|
phaseModelPartialList stationaryPhaseModels_;
|
|
|
|
//- Anisothermal phase models
|
|
phaseModelPartialList anisothermalPhaseModels_;
|
|
|
|
//- Multi-component phase models
|
|
phaseModelPartialList multicomponentPhaseModels_;
|
|
|
|
//- Total volumetric flux
|
|
surfaceScalarField phi_;
|
|
|
|
//- Rate of change of pressure
|
|
volScalarField dpdt_;
|
|
|
|
//- Interface compression coefficients
|
|
cAlphaTable cAlphas_;
|
|
|
|
//- Stabilisation for normalisation of the interface normal
|
|
const dimensionedScalar deltaN_;
|
|
|
|
|
|
// Sub Models
|
|
|
|
//- Surface tension models
|
|
interfaceSurfaceTensionModelTable interfaceSurfaceTensionModels_;
|
|
|
|
|
|
//- Flag to indicate that returned lists of fields are "complete"; i.e.,
|
|
// that an absence of force is returned as a constructed list of zeros,
|
|
// rather than a null pointer
|
|
static const bool fillFields_ = false;
|
|
|
|
|
|
// Protected member functions
|
|
|
|
//- Calculate and return the mixture flux
|
|
tmp<surfaceScalarField> calcPhi
|
|
(
|
|
const phaseModelList& phaseModels
|
|
) const;
|
|
|
|
//- Return the sum of the phase fractions of the moving phases
|
|
tmp<volScalarField> sumAlphaMoving() const;
|
|
|
|
//- Re-normalise the velocity of the phases
|
|
// around the specified mixture mean
|
|
void setMixtureU(const volVectorField& Um);
|
|
|
|
|
|
// Functions required for interface compression
|
|
|
|
//- Normal to interface between two phases
|
|
// Used for interface compression
|
|
tmp<surfaceVectorField> nHatfv
|
|
(
|
|
const volScalarField& alpha1,
|
|
const volScalarField& alpha2
|
|
) const;
|
|
|
|
//- Normal to interface between two phases dotted with face areas
|
|
// Used for interface compression
|
|
tmp<surfaceScalarField> nHatf
|
|
(
|
|
const volScalarField& alpha1,
|
|
const volScalarField& alpha2
|
|
) const;
|
|
|
|
//- Curvature of interface between two phases
|
|
// Used for interface compression
|
|
tmp<volScalarField> K
|
|
(
|
|
const phaseModel& alpha1,
|
|
const phaseModel& alpha2
|
|
) const;
|
|
|
|
|
|
// Sub-model construction
|
|
|
|
//- Return the dictionary containing interfacial model or value
|
|
// settings for the given name. Performs backwards compatibility
|
|
// conversions.
|
|
template<class Type>
|
|
dictionary interfacialDict(const word& name) const;
|
|
|
|
|
|
public:
|
|
|
|
//- Runtime type information
|
|
TypeName("phaseSystem");
|
|
|
|
//- Default name of the phase properties dictionary
|
|
static const word propertiesName;
|
|
|
|
|
|
// Declare runtime construction
|
|
|
|
declareRunTimeSelectionTable
|
|
(
|
|
autoPtr,
|
|
phaseSystem,
|
|
dictionary,
|
|
(
|
|
const fvMesh& mesh
|
|
),
|
|
(mesh)
|
|
);
|
|
|
|
|
|
// Constructors
|
|
|
|
//- Construct from fvMesh
|
|
phaseSystem(const fvMesh& mesh);
|
|
|
|
|
|
// Selectors
|
|
|
|
static autoPtr<phaseSystem> New
|
|
(
|
|
const fvMesh& mesh
|
|
);
|
|
|
|
|
|
//- Destructor
|
|
virtual ~phaseSystem();
|
|
|
|
|
|
// Member Functions
|
|
|
|
// Access
|
|
|
|
//- Return the mesh
|
|
inline const fvMesh& mesh() const;
|
|
|
|
//- Return pimpleNoLoopControl
|
|
inline const pimpleNoLoopControl& pimple() const;
|
|
|
|
//- Return the phase models
|
|
inline const phaseModelList& phases() const;
|
|
|
|
//- Access the phase models
|
|
inline phaseModelList& phases();
|
|
|
|
//- Return the models for phases that are moving
|
|
inline const phaseModelPartialList& movingPhases() const;
|
|
|
|
//- Access the models for phases that are moving
|
|
inline phaseModelPartialList& movingPhases();
|
|
|
|
//- Return the models for phases that are stationary
|
|
inline const phaseModelPartialList& stationaryPhases() const;
|
|
|
|
//- Access the models for phases that are stationary
|
|
inline phaseModelPartialList& stationaryPhases();
|
|
|
|
//- Return the models for phases that have variable temperature
|
|
inline const phaseModelPartialList& anisothermalPhases() const;
|
|
|
|
//- Access the models for phases that have variable temperature
|
|
inline phaseModelPartialList& anisothermalPhases();
|
|
|
|
//- Return the models for phases that have multiple species
|
|
inline const phaseModelPartialList& multicomponentPhases() const;
|
|
|
|
//- Access the models for phases that have multiple species
|
|
inline phaseModelPartialList& multicomponentPhases();
|
|
|
|
//- Return the phase not given as an argument in a two-phase system
|
|
// An error is generated if the system is not two-phase
|
|
inline const phaseModel& otherPhase(const phaseModel& phase) const;
|
|
|
|
//- Return the mixture flux
|
|
inline const surfaceScalarField& phi() const;
|
|
|
|
//- Access the mixture flux
|
|
inline surfaceScalarField& phi();
|
|
|
|
//- Return the rate of change of the pressure
|
|
inline const volScalarField& dpdt() const;
|
|
|
|
//- Access the rate of change of the pressure
|
|
inline volScalarField& dpdt();
|
|
|
|
//- Return MRF zones
|
|
inline const IOMRFZoneList& MRF() const;
|
|
|
|
//- Access the fvModels
|
|
inline Foam::fvModels& fvModels(fvMesh& mesh);
|
|
|
|
//- Access the fvModels
|
|
inline const Foam::fvModels& fvModels() const;
|
|
|
|
//- Access the fvConstraints
|
|
inline Foam::fvConstraints& fvConstraints(fvMesh& mesh);
|
|
|
|
//- Access the fvConstraints
|
|
inline const Foam::fvConstraints& fvConstraints() const;
|
|
|
|
|
|
// Sub-model construction
|
|
|
|
//- Return the model name. This is the same as the model's typename
|
|
// but without "Model" on the end.
|
|
template<class ModelType>
|
|
word modelName() const;
|
|
|
|
//- Generate interfacial-model lists
|
|
template<class ModelType, class ... InterfaceTypes>
|
|
void generateInterfacialModels
|
|
(
|
|
const dictionary& dict,
|
|
const phaseInterface& interface,
|
|
PtrList<phaseInterface>& interfaces,
|
|
PtrList<ModelType>& models
|
|
) const;
|
|
|
|
//- Generate interfacial-model tables
|
|
template<class ModelType>
|
|
void generateInterfacialModels
|
|
(
|
|
const dictionary& dict,
|
|
HashTable
|
|
<
|
|
autoPtr<ModelType>,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>& models
|
|
) const;
|
|
|
|
//- Generate interfacial-model tables
|
|
template<class ModelType>
|
|
void generateInterfacialModels
|
|
(
|
|
HashTable
|
|
<
|
|
autoPtr<ModelType>,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>& models
|
|
) const;
|
|
|
|
//- Generate interfacial-model tables
|
|
template<class ValueType>
|
|
void generateInterfacialValues
|
|
(
|
|
const dictionary& dict,
|
|
HashTable
|
|
<
|
|
ValueType,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>& values
|
|
) const;
|
|
|
|
//- Generate interfacial-model tables
|
|
template<class ValueType>
|
|
void generateInterfacialValues
|
|
(
|
|
const word& valueName,
|
|
HashTable
|
|
<
|
|
ValueType,
|
|
phaseInterfaceKey,
|
|
phaseInterfaceKey::hash
|
|
>& values
|
|
) const;
|
|
|
|
//- Return the dictionary from which to construct a low-level
|
|
// sub-model. Checks that there is just one sub-dictionary then
|
|
// returns it.
|
|
template<class ModelType>
|
|
static const dictionary& modelSubDict(const dictionary& dict);
|
|
|
|
//- Check that mass transfer is supported across the given interface
|
|
template<class ModelType>
|
|
void validateMassTransfer(const phaseInterface& interface) const;
|
|
|
|
|
|
// Sub-model lookup
|
|
|
|
//- Check availability of a sub model for a given interface
|
|
template<class ModelType>
|
|
bool foundInterfacialModel(const phaseInterface& interface) const;
|
|
|
|
//- Return a sub model for an interface
|
|
template<class ModelType>
|
|
const ModelType& lookupInterfacialModel
|
|
(
|
|
const phaseInterface& interface
|
|
) const;
|
|
|
|
|
|
// Field construction
|
|
|
|
//- Fill up gaps in a phase-indexed list of fields with zeros
|
|
template
|
|
<
|
|
class Type,
|
|
template<class> class PatchField,
|
|
class GeoMesh
|
|
>
|
|
void fillFields
|
|
(
|
|
const word& name,
|
|
const dimensionSet& dims,
|
|
PtrList<GeometricField<Type, PatchField, GeoMesh>>& fieldList
|
|
) const;
|
|
|
|
//- Fill up gaps in a phase-indexed table of fields with zeros
|
|
template
|
|
<
|
|
class Type,
|
|
template<class> class PatchField,
|
|
class GeoMesh
|
|
>
|
|
void fillFields
|
|
(
|
|
const word& name,
|
|
const dimensionSet& dims,
|
|
HashPtrTable<GeometricField<Type, PatchField, GeoMesh>>&
|
|
fieldTable
|
|
) const;
|
|
|
|
|
|
// Properties
|
|
|
|
//- Return the mixture density
|
|
tmp<volScalarField> rho() const;
|
|
|
|
//- Return the mixture velocity
|
|
tmp<volVectorField> U() const;
|
|
|
|
//- Return the surface tension coefficient for an interface
|
|
tmp<volScalarField> sigma(const phaseInterfaceKey& key) const;
|
|
|
|
//- Return the surface tension coefficient for an interface on a
|
|
// patch
|
|
tmp<scalarField> sigma
|
|
(
|
|
const phaseInterfaceKey& key,
|
|
const label patchi
|
|
) const;
|
|
|
|
//- Indicator of the proximity of the interface
|
|
// Field values are 1 near and 0 away for the interface.
|
|
tmp<volScalarField> nearInterface() const;
|
|
|
|
//- Stabilisation for normalisation of the interface normal
|
|
inline const dimensionedScalar& deltaN() const;
|
|
|
|
//- Return the mass transfer rate for an interface
|
|
virtual tmp<volScalarField> dmdtf
|
|
(
|
|
const phaseInterfaceKey& key
|
|
) const;
|
|
|
|
//- Return the mass transfer rates for each phase
|
|
virtual PtrList<volScalarField> dmdts() const;
|
|
|
|
//- Return the mass transfer pressure implicit coefficients
|
|
// for each phase
|
|
virtual PtrList<volScalarField> d2mdtdps() const;
|
|
|
|
//- Return incompressibility
|
|
bool incompressible() const;
|
|
|
|
|
|
// Transfers
|
|
|
|
//- Return the momentum transfer matrices for the cell-based
|
|
// algorithm
|
|
virtual autoPtr<momentumTransferTable> momentumTransfer() = 0;
|
|
|
|
//- Return the momentum transfer matrices for the face-based
|
|
// algorithm
|
|
virtual autoPtr<momentumTransferTable> momentumTransferf() = 0;
|
|
|
|
//- Return the implicit force coefficients for the face-based
|
|
// algorithm
|
|
virtual PtrList<surfaceScalarField> KdVmfs() const = 0;
|
|
|
|
//- Return the force fluxes for the cell-based algorithm
|
|
virtual PtrList<surfaceScalarField> Fs() const = 0;
|
|
|
|
//- Return the force fluxes for the face-based algorithm
|
|
virtual PtrList<surfaceScalarField> Ffs() const = 0;
|
|
|
|
//- Return the force fluxes for the cell-based algorithm
|
|
virtual PtrList<surfaceScalarField> KdPhis() const = 0;
|
|
|
|
//- Return the force fluxes for the face-based algorithm
|
|
virtual PtrList<surfaceScalarField> KdPhifs() const = 0;
|
|
|
|
//- Return the implicit part of the drag force
|
|
virtual PtrList<volScalarField> Kds() const = 0;
|
|
|
|
//- Return the explicit part of the drag force
|
|
virtual PtrList<volVectorField> KdUs() const = 0;
|
|
|
|
//- Returns true if the phase pressure is treated implicitly
|
|
// in the phase fraction equation
|
|
virtual bool implicitPhasePressure(const phaseModel& phase) const;
|
|
|
|
//- Returns true if the phase pressure is treated implicitly
|
|
// in the phase fraction equation for any phase
|
|
virtual bool implicitPhasePressure() const;
|
|
|
|
//- Return the phase diffusivity
|
|
// divided by the momentum central coefficient
|
|
virtual tmp<surfaceScalarField> alphaDByAf
|
|
(
|
|
const PtrList<volScalarField>& rAUs,
|
|
const PtrList<surfaceScalarField>& rAUfs
|
|
) const = 0;
|
|
|
|
//- Return the flux corrections for the cell-based algorithm
|
|
virtual PtrList<surfaceScalarField> ddtCorrs() const = 0;
|
|
|
|
//- Set the cell and faces drag correction fields
|
|
virtual void dragCorrs
|
|
(
|
|
PtrList<volVectorField>& dragCorrs,
|
|
PtrList<surfaceScalarField>& dragCorrf
|
|
) const = 0;
|
|
|
|
//- Solve the drag system for the new velocities and fluxes
|
|
virtual void partialElimination
|
|
(
|
|
const PtrList<volScalarField>& rAUs,
|
|
const PtrList<volVectorField>& KdUs,
|
|
const PtrList<surfaceScalarField>& alphafs,
|
|
const PtrList<surfaceScalarField>& rAUfs,
|
|
const PtrList<surfaceScalarField>& KdPhis
|
|
) = 0;
|
|
|
|
//- Solve the drag system for the new fluxes
|
|
virtual void partialEliminationf
|
|
(
|
|
const PtrList<surfaceScalarField>& rAUfs,
|
|
const PtrList<surfaceScalarField>& alphafs,
|
|
const PtrList<surfaceScalarField>& KdPhifs
|
|
) = 0;
|
|
|
|
//- Re-normalise the flux of the phases
|
|
// around the specified mixture mean
|
|
void setMixturePhi
|
|
(
|
|
const PtrList<surfaceScalarField>& alphafs,
|
|
const surfaceScalarField& phim
|
|
);
|
|
|
|
//- Return the heat transfer matrices
|
|
virtual autoPtr<heatTransferTable> heatTransfer() const = 0;
|
|
|
|
//- Return the specie transfer matrices
|
|
virtual autoPtr<specieTransferTable> specieTransfer() const = 0;
|
|
|
|
//- Return the surface tension force
|
|
tmp<surfaceScalarField> surfaceTension
|
|
(
|
|
const phaseModel& phase
|
|
) const;
|
|
|
|
|
|
// Evolution
|
|
|
|
//- Solve for the phase fractions
|
|
virtual void solve
|
|
(
|
|
const PtrList<volScalarField>& rAUs,
|
|
const PtrList<surfaceScalarField>& rAUfs
|
|
);
|
|
|
|
//- Correct the fluid properties other than those listed below
|
|
virtual void correct();
|
|
|
|
//- Correct the continuity errors
|
|
virtual void correctContinuityError();
|
|
|
|
//- Correct the kinematics
|
|
virtual void correctKinematics();
|
|
|
|
//- Correct the thermodynamics
|
|
virtual void correctThermo();
|
|
|
|
//- Correct the reactions
|
|
virtual void correctReactions();
|
|
|
|
//- Correct the species mass fractions
|
|
virtual void correctSpecies();
|
|
|
|
//- Predict the momentumTransport
|
|
virtual void predictMomentumTransport();
|
|
|
|
//- Predict the energy transport e.g. alphat
|
|
virtual void predictThermophysicalTransport();
|
|
|
|
//- Correct the momentumTransport
|
|
virtual void correctMomentumTransport();
|
|
|
|
//- Correct the energy transport e.g. alphat
|
|
virtual void correctThermophysicalTransport();
|
|
|
|
//- Update the fluid properties for mesh changes
|
|
virtual void meshUpdate();
|
|
|
|
//- Correct fixed-flux BCs to be consistent with the velocity BCs
|
|
void correctBoundaryFlux();
|
|
|
|
void correctPhi
|
|
(
|
|
const volScalarField& p_rgh,
|
|
const autoPtr<volScalarField>& divU,
|
|
const pressureReference& pressureReference,
|
|
nonOrthogonalSolutionControl& pimple
|
|
);
|
|
|
|
|
|
// IO
|
|
|
|
//- Read base phaseProperties dictionary
|
|
virtual bool read();
|
|
};
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
tmp<volScalarField> byDt(const volScalarField& vf);
|
|
tmp<surfaceScalarField> byDt(const surfaceScalarField& sf);
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
} // End namespace Foam
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
#include "phaseSystemI.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
#ifdef NoRepository
|
|
#include "phaseSystemTemplates.C"
|
|
#endif
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
#endif
|
|
|
|
// ************************************************************************* //
|