/*---------------------------------------------------------------------------*\
========= |
\\ / 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 .
Class
Foam::MomentumTransferPhaseSystem
Description
Class which models interfacial momentum transfer between a number of phases.
Drag, virtual mass, lift, wall lubrication and turbulent dispersion are all
modelled. The explicit contribution from the drag is omitted from the
transfer matrices, as this forms part of the solution of the pressure
equation.
SourceFiles
MomentumTransferPhaseSystem.C
\*---------------------------------------------------------------------------*/
#ifndef MomentumTransferPhaseSystem_H
#define MomentumTransferPhaseSystem_H
#include "phaseSystem.H"
#include "HashPtrTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class blendingMethod;
class blendedDragModel;
class blendedVirtualMassModel;
class blendedLiftModel;
class blendedWallLubricationModel;
class blendedTurbulentDispersionModel;
/*---------------------------------------------------------------------------*\
Class MomentumTransferPhaseSystem Declaration
\*---------------------------------------------------------------------------*/
template
class MomentumTransferPhaseSystem
:
public BasePhaseSystem
{
// Private typedefs
typedef HashPtrTable
<
volScalarField,
phaseInterfaceKey,
phaseInterfaceKey::hash
> KdTable;
typedef HashPtrTable
<
surfaceScalarField,
phaseInterfaceKey,
phaseInterfaceKey::hash
> KdfTable;
typedef HashPtrTable
<
volScalarField,
phaseInterfaceKey,
phaseInterfaceKey::hash
> VmTable;
typedef HashTable
<
autoPtr,
phaseInterfaceKey,
phaseInterfaceKey::hash
> dragModelTable;
typedef HashTable
<
autoPtr,
phaseInterfaceKey,
phaseInterfaceKey::hash
> virtualMassModelTable;
typedef HashTable
<
autoPtr,
phaseInterfaceKey,
phaseInterfaceKey::hash
> liftModelTable;
typedef HashTable
<
autoPtr,
phaseInterfaceKey,
phaseInterfaceKey::hash
> wallLubricationModelTable;
typedef HashTable
<
autoPtr,
phaseInterfaceKey,
phaseInterfaceKey::hash
> turbulentDispersionModelTable;
// Private Data
//- Drag coefficients
KdTable Kds_;
//- Face drag coefficients
KdfTable Kdfs_;
//- Virtual mass coefficients
VmTable Vms_;
// Sub Models
//- Drag models
dragModelTable dragModels_;
//- Virtual mass models
virtualMassModelTable virtualMassModels_;
//- Lift models
liftModelTable liftModels_;
//- Wall lubrication models
wallLubricationModelTable wallLubricationModels_;
//- Turbulent dispersion models
turbulentDispersionModelTable turbulentDispersionModels_;
protected:
// Protected Member Functions
//- Add momentum transfer terms which result from bulk mass transfers
void addDmdtUfs
(
const phaseSystem::dmdtfTable& dmdtfs,
phaseSystem::momentumTransferTable& eqns
);
void addTmpField
(
tmp& result,
const tmp& field
) const;
public:
// Constructors
//- Construct from fvMesh
MomentumTransferPhaseSystem(const fvMesh&);
//- Destructor
virtual ~MomentumTransferPhaseSystem();
// Member Functions
//- Return the momentum transfer matrices for the cell-based algorithm.
// This includes implicit and explicit forces that add into the cell
// UEqn in the normal way.
virtual autoPtr momentumTransfer();
//- As momentumTransfer, but for the face-based algorithm
virtual autoPtr momentumTransferf();
//- Return implicit force coefficients on the faces, for the face-based
// algorithm.
virtual PtrList KdVmfs() const;
//- Return the explicit force fluxes for the cell-based algorithm, that
// do not depend on phase mass/volume fluxes, and can therefore be
// evaluated outside the corrector loop. This includes things like
// lift, turbulent dispersion, and wall lubrication.
virtual PtrList Fs() const;
//- As Fs, but for the face-based algorithm
virtual PtrList Ffs() const;
//- Return the explicit drag force fluxes for the cell-based algorithm.
// These depend on phase mass/volume fluxes, and must therefore be
// evaluated inside the corrector loop.
virtual PtrList KdPhis() const;
//- As KdPhis, but for the face-based algorithm
virtual PtrList KdPhifs() const;
//- Return the implicit part of the drag force
virtual PtrList Kds() const;
//- Return the explicit part of the drag force for the cell-based
// algorithm. This is the cell-equivalent of KdPhis. These depend on
// phase velocities, and must therefore be evaluated inside the
// corrector loop.
virtual PtrList KdUs() const;
//- 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 alphaDByAf
(
const PtrList& rAUs,
const PtrList& rAUfs
) const;
//- Return the flux corrections for the cell-based algorithm. These
// depend on phase mass/volume fluxes, and must therefore be evaluated
// inside the corrector loop.
virtual PtrList ddtCorrs() const;
//- Set the cell and faces drag correction fields
virtual void dragCorrs
(
PtrList& dragCorrs,
PtrList& dragCorrf
) const;
//- Solve the drag system for the velocities and fluxes
virtual void partialElimination
(
const PtrList& rAUs,
const PtrList& KdUs,
const PtrList& alphafs,
const PtrList& rAUfs,
const PtrList& KdPhis
);
//- As partialElimination, but for the face-based algorithm. Only solves
// for the fluxes.
virtual void partialEliminationf
(
const PtrList& rAUfs,
const PtrList& alphafs,
const PtrList& KdPhifs
);
//- Read base phaseProperties dictionary
virtual bool read();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "MomentumTransferPhaseSystem.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //