ENH: driftFluxFoam: made viscosity and relative velocity modelling run-time selectable
This commit is contained in:
parent
02e20513d9
commit
6625be1dc4
9
applications/solvers/multiphase/driftFluxFoam/Allwclean
Executable file
9
applications/solvers/multiphase/driftFluxFoam/Allwclean
Executable file
@ -0,0 +1,9 @@
|
||||
#!/bin/sh
|
||||
cd ${0%/*} || exit 1 # run from this directory
|
||||
set -x
|
||||
|
||||
wclean libso viscosityModels
|
||||
wclean libso relativeVelocityModels
|
||||
wclean
|
||||
|
||||
# ----------------------------------------------------------------- end-of-file
|
9
applications/solvers/multiphase/driftFluxFoam/Allwmake
Executable file
9
applications/solvers/multiphase/driftFluxFoam/Allwmake
Executable file
@ -0,0 +1,9 @@
|
||||
#!/bin/sh
|
||||
cd ${0%/*} || exit 1 # run from this directory
|
||||
set -x
|
||||
|
||||
wmake libso viscosityModels
|
||||
wmake libso relativeVelocityModels
|
||||
wmake
|
||||
|
||||
# ----------------------------------------------------------------- end-of-file
|
@ -2,10 +2,18 @@ EXE_INC = \
|
||||
-I$(LIB_SRC)/finiteVolume/lnInclude \
|
||||
-I$(LIB_SRC)/meshTools/lnInclude \
|
||||
-I$(LIB_SRC)/sampling/lnInclude \
|
||||
-I$(LIB_SRC)/fvOptions/lnInclude
|
||||
-I$(LIB_SRC)/fvOptions/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels \
|
||||
-I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels/incompressible/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude \
|
||||
-I./relativeVelocityModels/lnInclude
|
||||
|
||||
EXE_LIBS = \
|
||||
-lfiniteVolume \
|
||||
-lmeshTools \
|
||||
-lsampling \
|
||||
-lfvOptions
|
||||
-lfvOptions \
|
||||
-lincompressibleTransportModels \
|
||||
-ldriftFluxTransportModels \
|
||||
-ldriftFluxRelativeVelocityModels
|
||||
|
@ -4,11 +4,7 @@
|
||||
(
|
||||
fvm::ddt(rho, U)
|
||||
+ fvm::div(rhoPhi, U)
|
||||
+ fvc::div
|
||||
(
|
||||
(alpha1/(scalar(1.001) - alpha1))*((rho2*rho1)/rho)*Vdj*Vdj,
|
||||
"div(phiVdj,Vdj)"
|
||||
)
|
||||
+ fvc::div(uRelModel.tau(), "div(phiUkm,Ukm)")
|
||||
- fvm::laplacian(muEff, U)
|
||||
- fvc::div(muEff*dev2(T(fvc::grad(U))))
|
||||
==
|
||||
|
@ -13,7 +13,7 @@
|
||||
|
||||
surfaceScalarField phir
|
||||
(
|
||||
rho2*(mesh.Sf() & fvc::interpolate(Vdj/rho))
|
||||
mesh.Sf() & fvc::interpolate(uRelModel.Udm())
|
||||
);
|
||||
|
||||
if (nAlphaSubCycles > 1)
|
||||
|
@ -1,20 +0,0 @@
|
||||
if (VdjModel == "general")
|
||||
{
|
||||
Vdj = V0*
|
||||
(
|
||||
exp(-a*max(alpha1 - alphaMin, scalar(0)))
|
||||
- exp(-a1*max(alpha1 - alphaMin, scalar(0)))
|
||||
);
|
||||
}
|
||||
else if (VdjModel == "simple")
|
||||
{
|
||||
Vdj = V0*pow(10.0, -a*max(alpha1, scalar(0)));
|
||||
}
|
||||
else
|
||||
{
|
||||
FatalErrorIn(args.executable())
|
||||
<< "Unknown VdjModel : " << VdjModel
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
Vdj.correctBoundaryConditions();
|
@ -1,39 +0,0 @@
|
||||
{
|
||||
mul = muc +
|
||||
plasticViscosity
|
||||
(
|
||||
plasticViscosityCoeff,
|
||||
plasticViscosityExponent,
|
||||
alpha1
|
||||
);
|
||||
|
||||
if (BinghamPlastic)
|
||||
{
|
||||
volScalarField tauy = yieldStress
|
||||
(
|
||||
yieldStressCoeff,
|
||||
yieldStressExponent,
|
||||
yieldStressOffset,
|
||||
alpha1
|
||||
);
|
||||
|
||||
mul =
|
||||
tauy/
|
||||
(
|
||||
mag(fvc::grad(U))
|
||||
+ 1.0e-4*
|
||||
(
|
||||
tauy
|
||||
+ dimensionedScalar
|
||||
(
|
||||
"deltaTauy",
|
||||
tauy.dimensions(),
|
||||
1.0e-15
|
||||
)
|
||||
)/mul
|
||||
)
|
||||
+ mul;
|
||||
}
|
||||
|
||||
mul = min(mul, muMax);
|
||||
}
|
@ -12,22 +12,6 @@
|
||||
mesh
|
||||
);
|
||||
|
||||
Info<< "Reading field alpha1\n" << endl;
|
||||
volScalarField alpha1
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"alpha1",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::MUST_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mesh
|
||||
);
|
||||
|
||||
volScalarField alpha2("alpha2", scalar(1) - alpha1);
|
||||
|
||||
Info<< "Reading field U\n" << endl;
|
||||
volVectorField U
|
||||
(
|
||||
@ -45,7 +29,17 @@
|
||||
#include "createPhi.H"
|
||||
|
||||
|
||||
// Transport
|
||||
// ~~~~~~~~~
|
||||
|
||||
Info<< "Reading transportProperties\n" << endl;
|
||||
incompressibleTwoPhaseMixture twoPhaseProperties(U, phi);
|
||||
|
||||
volScalarField& alpha1(twoPhaseProperties.alpha1());
|
||||
volScalarField& alpha2(twoPhaseProperties.alpha2());
|
||||
|
||||
const dimensionedScalar& rho1 = twoPhaseProperties.rho1();
|
||||
const dimensionedScalar& rho2 = twoPhaseProperties.rho2();
|
||||
|
||||
IOdictionary transportProperties
|
||||
(
|
||||
@ -59,40 +53,7 @@
|
||||
)
|
||||
);
|
||||
|
||||
|
||||
dimensionedScalar rho1(transportProperties.lookup("rho1"));
|
||||
dimensionedScalar rho2(transportProperties.lookup("rho2"));
|
||||
|
||||
dimensionedScalar muc(transportProperties.lookup("muc"));
|
||||
dimensionedScalar muMax(transportProperties.lookup("muMax"));
|
||||
|
||||
dimensionedScalar plasticViscosityCoeff
|
||||
(
|
||||
transportProperties.lookup("plasticViscosityCoeff")
|
||||
);
|
||||
|
||||
dimensionedScalar plasticViscosityExponent
|
||||
(
|
||||
transportProperties.lookup("plasticViscosityExponent")
|
||||
);
|
||||
|
||||
dimensionedScalar yieldStressCoeff
|
||||
(
|
||||
transportProperties.lookup("yieldStressCoeff")
|
||||
);
|
||||
|
||||
dimensionedScalar yieldStressExponent
|
||||
(
|
||||
transportProperties.lookup("yieldStressExponent")
|
||||
);
|
||||
|
||||
dimensionedScalar yieldStressOffset
|
||||
(
|
||||
transportProperties.lookup("yieldStressOffset")
|
||||
);
|
||||
|
||||
Switch BinghamPlastic(transportProperties.lookup("BinghamPlastic"));
|
||||
|
||||
// Mixture density
|
||||
volScalarField rho
|
||||
(
|
||||
IOobject
|
||||
@ -121,63 +82,24 @@
|
||||
fvc::interpolate(rho)*phi
|
||||
);
|
||||
|
||||
Info<< "Calculating field mul\n" << endl;
|
||||
volScalarField mul
|
||||
|
||||
// Relative Velocity
|
||||
// ~~~~~~~~~~~~~~~~~
|
||||
|
||||
autoPtr<relativeVelocityModel> uRelModelPtr
|
||||
(
|
||||
IOobject
|
||||
relativeVelocityModel::New
|
||||
(
|
||||
"mul",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
muc
|
||||
+ plasticViscosity
|
||||
(
|
||||
plasticViscosityCoeff,
|
||||
plasticViscosityExponent,
|
||||
alpha1
|
||||
transportProperties,
|
||||
twoPhaseProperties
|
||||
)
|
||||
);
|
||||
|
||||
|
||||
Info<< "Initialising field Vdj\n" << endl;
|
||||
volVectorField Vdj
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"Vdj",
|
||||
runTime.timeName(),
|
||||
mesh,
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mesh,
|
||||
dimensionedVector("0.0", U.dimensions(), vector::zero),
|
||||
U.boundaryField().types()
|
||||
);
|
||||
relativeVelocityModel& uRelModel(uRelModelPtr());
|
||||
|
||||
|
||||
Info<< "Selecting Drift-Flux model " << endl;
|
||||
|
||||
const word VdjModel(transportProperties.lookup("VdjModel"));
|
||||
|
||||
Info<< tab << VdjModel << " selected\n" << endl;
|
||||
|
||||
const dictionary& VdjModelCoeffs
|
||||
(
|
||||
transportProperties.subDict(VdjModel + "Coeffs")
|
||||
);
|
||||
|
||||
dimensionedVector V0(VdjModelCoeffs.lookup("V0"));
|
||||
|
||||
dimensionedScalar a(VdjModelCoeffs.lookup("a"));
|
||||
|
||||
dimensionedScalar a1(VdjModelCoeffs.lookup("a1"));
|
||||
|
||||
dimensionedScalar alphaMin(VdjModelCoeffs.lookup("alphaMin"));
|
||||
|
||||
// Turbulence
|
||||
// ~~~~~~~~~~
|
||||
|
||||
IOdictionary RASProperties
|
||||
(
|
||||
@ -284,7 +206,6 @@
|
||||
<< "wallFunctionCoeffs" << wallFunctionDict << endl;
|
||||
}
|
||||
|
||||
|
||||
nearWallDist y(mesh);
|
||||
|
||||
Info<< "Reading field k\n" << endl;
|
||||
@ -329,7 +250,6 @@
|
||||
Cmu*rho*sqr(k)/epsilon
|
||||
);
|
||||
|
||||
|
||||
Info<< "Calculating field muEff\n" << endl;
|
||||
volScalarField muEff
|
||||
(
|
||||
@ -341,10 +261,13 @@
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
mut + mul
|
||||
mut + twoPhaseProperties.mu()
|
||||
);
|
||||
|
||||
|
||||
// Pressure
|
||||
// ~~~~~~~~
|
||||
|
||||
Info<< "Calculating field (g.h)f\n" << endl;
|
||||
volScalarField gh("gh", g & mesh.C());
|
||||
surfaceScalarField ghf("gh", g & mesh.Cf());
|
||||
@ -385,4 +308,7 @@
|
||||
}
|
||||
|
||||
|
||||
// MULES Correction
|
||||
// ~~~~~~~~~~~~~~~~
|
||||
|
||||
tmp<surfaceScalarField> tphiAlphaCorr0;
|
||||
|
@ -36,12 +36,12 @@ Description
|
||||
#include "fvCFD.H"
|
||||
#include "CMULES.H"
|
||||
#include "subCycle.H"
|
||||
#include "incompressibleTwoPhaseMixture.H"
|
||||
#include "relativeVelocityModel.H"
|
||||
#include "nearWallDist.H"
|
||||
#include "wallFvPatch.H"
|
||||
#include "bound.H"
|
||||
#include "Switch.H"
|
||||
#include "plasticViscosity.H"
|
||||
#include "yieldStress.H"
|
||||
#include "pimpleControl.H"
|
||||
#include "fvIOoptionList.H"
|
||||
#include "fixedFluxPressureFvPatchScalarField.H"
|
||||
@ -82,9 +82,11 @@ int main(int argc, char *argv[])
|
||||
{
|
||||
#include "alphaControls.H"
|
||||
|
||||
#include "calcVdj.H"
|
||||
uRelModel.update();
|
||||
|
||||
#include "alphaEqnSubCycle.H"
|
||||
#include "correctViscosity.H"
|
||||
|
||||
twoPhaseProperties.correct();
|
||||
|
||||
#include "UEqn.H"
|
||||
|
||||
|
@ -21,6 +21,8 @@ if (turbulence)
|
||||
Cmu*k/sigmak*(g & fvc::grad(rho))/(epsilon + epsilonMin)
|
||||
);
|
||||
|
||||
volScalarField muc(twoPhaseProperties.nuModel2().nu()*rho2);
|
||||
|
||||
#include "wallFunctions.H"
|
||||
|
||||
// Dissipation equation
|
||||
@ -75,4 +77,4 @@ if (turbulence)
|
||||
#include "wallViscosity.H"
|
||||
}
|
||||
|
||||
muEff = mut + mul;
|
||||
muEff = mut + twoPhaseProperties.mu();
|
||||
|
@ -1,21 +0,0 @@
|
||||
volScalarField plasticViscosity
|
||||
(
|
||||
const dimensionedScalar& plasticViscosityCoeff,
|
||||
const dimensionedScalar& plasticViscosityExponent,
|
||||
const volScalarField& Alpha
|
||||
)
|
||||
{
|
||||
tmp<volScalarField> tfld
|
||||
(
|
||||
plasticViscosityCoeff*
|
||||
(
|
||||
pow
|
||||
(
|
||||
10.0,
|
||||
plasticViscosityExponent*Alpha + SMALL
|
||||
) - scalar(1)
|
||||
)
|
||||
);
|
||||
|
||||
return tfld();
|
||||
}
|
@ -0,0 +1,5 @@
|
||||
relativeVelocityModel/relativeVelocityModel.C
|
||||
simple/simple.C
|
||||
general/general.C
|
||||
|
||||
LIB = $(FOAM_LIBBIN)/libdriftFluxRelativeVelocityModels
|
@ -0,0 +1,10 @@
|
||||
EXE_INC = \
|
||||
-I$(LIB_SRC)/finiteVolume/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels \
|
||||
-I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels/incompressible/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels/interfaceProperties/lnInclude
|
||||
|
||||
LIB_LIBS = \
|
||||
-lfiniteVolume \
|
||||
-lincompressibleTransportModels
|
@ -0,0 +1,78 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "general.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace relativeVelocityModels
|
||||
{
|
||||
defineTypeNameAndDebug(general, 0);
|
||||
addToRunTimeSelectionTable(relativeVelocityModel, general, dictionary);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::relativeVelocityModels::general::general
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
)
|
||||
:
|
||||
relativeVelocityModel(dict, mixture),
|
||||
a_(dict.lookup("a")),
|
||||
a1_(dict.lookup("a1")),
|
||||
V0_(dict.lookup("V0")),
|
||||
residualAlpha_(dict.lookup("residualAlpha"))
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::relativeVelocityModels::general::~general()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::tmp<Foam::volVectorField>
|
||||
Foam::relativeVelocityModels::general::Ur() const
|
||||
{
|
||||
return
|
||||
V0_
|
||||
*(
|
||||
exp(-a_*max(alphaD_ - residualAlpha_, scalar(0)))
|
||||
- exp(-a1_*max(alphaD_ - residualAlpha_, scalar(0)))
|
||||
)
|
||||
/max(alphaC_, residualAlpha_);
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,106 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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::general
|
||||
|
||||
Description
|
||||
General relative velocity model
|
||||
|
||||
SourceFiles
|
||||
general.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef general_H
|
||||
#define general_H
|
||||
|
||||
#include "relativeVelocityModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace relativeVelocityModels
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class general Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class general
|
||||
:
|
||||
public relativeVelocityModel
|
||||
{
|
||||
// Private data
|
||||
|
||||
//- a coefficient
|
||||
dimensionedScalar a_;
|
||||
|
||||
//- a1 coefficient
|
||||
dimensionedScalar a1_;
|
||||
|
||||
//- Drift velocity
|
||||
dimensionedVector V0_;
|
||||
|
||||
//- Residual phase fraction
|
||||
dimensionedScalar residualAlpha_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("general");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
general
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
~general();
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Relative velocity
|
||||
virtual tmp<volVectorField> Ur() const;
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace relativeVelocityModels
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,175 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "relativeVelocityModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
defineTypeNameAndDebug(relativeVelocityModel, 0);
|
||||
defineRunTimeSelectionTable(relativeVelocityModel, dictionary);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::relativeVelocityModel::relativeVelocityModel
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
)
|
||||
:
|
||||
mixture_(mixture),
|
||||
|
||||
continuousPhaseName_(dict.lookup("continuousPhase")),
|
||||
|
||||
alphaC_
|
||||
(
|
||||
mixture.phase1Name() == continuousPhaseName_
|
||||
? mixture.alpha1()
|
||||
: mixture.alpha2()
|
||||
),
|
||||
|
||||
alphaD_
|
||||
(
|
||||
mixture.phase1Name() == continuousPhaseName_
|
||||
? mixture.alpha2()
|
||||
: mixture.alpha1()
|
||||
),
|
||||
|
||||
rhoC_
|
||||
(
|
||||
mixture.phase1Name() == continuousPhaseName_
|
||||
? mixture.rho1()
|
||||
: mixture.rho2()
|
||||
),
|
||||
|
||||
rhoD_
|
||||
(
|
||||
mixture.phase1Name() == continuousPhaseName_
|
||||
? mixture.rho2()
|
||||
: mixture.rho1()
|
||||
),
|
||||
|
||||
Udm_
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"Udm",
|
||||
alphaC_.time().timeName(),
|
||||
alphaC_.mesh()
|
||||
),
|
||||
alphaC_.mesh(),
|
||||
dimensionedVector("Udm", dimVelocity, vector::zero),
|
||||
mixture.U().boundaryField().types()
|
||||
),
|
||||
|
||||
tau_
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
"Udm",
|
||||
alphaC_.time().timeName(),
|
||||
alphaC_.mesh()
|
||||
),
|
||||
alphaC_.mesh(),
|
||||
dimensionedSymmTensor
|
||||
(
|
||||
"Udm",
|
||||
sqr(dimVelocity)*dimDensity,
|
||||
symmTensor::zero
|
||||
)
|
||||
)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::autoPtr<Foam::relativeVelocityModel> Foam::relativeVelocityModel::New
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
)
|
||||
{
|
||||
word modelType(dict.lookup(typeName));
|
||||
|
||||
Info<< "Selecting relative velocity model " << modelType << endl;
|
||||
|
||||
dictionaryConstructorTable::iterator cstrIter =
|
||||
dictionaryConstructorTablePtr_->find(modelType);
|
||||
|
||||
if (cstrIter == dictionaryConstructorTablePtr_->end())
|
||||
{
|
||||
FatalErrorIn
|
||||
(
|
||||
"relativeVelocityModel::New"
|
||||
"("
|
||||
"const dictionary&"
|
||||
")"
|
||||
) << "Unknown time scale model type " << modelType
|
||||
<< ", constructor not in hash table" << nl << nl
|
||||
<< " Valid time scale model types are:" << nl
|
||||
<< dictionaryConstructorTablePtr_->sortedToc()
|
||||
<< abort(FatalError);
|
||||
}
|
||||
|
||||
return
|
||||
autoPtr<relativeVelocityModel>
|
||||
(
|
||||
cstrIter()
|
||||
(
|
||||
dict.subDict(modelType + "Coeffs"),
|
||||
mixture
|
||||
)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::relativeVelocityModel::~relativeVelocityModel()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
||||
|
||||
void Foam::relativeVelocityModel::update()
|
||||
{
|
||||
tmp<volVectorField> URel(Ur());
|
||||
|
||||
tmp<volScalarField> betaC(alphaC_*rhoC_);
|
||||
tmp<volScalarField> betaD(alphaD_*rhoD_);
|
||||
tmp<volScalarField> rhoM(betaC() + betaD());
|
||||
|
||||
tmp<volVectorField> Udm = URel()*betaC()/rhoM;
|
||||
tmp<volVectorField> Ucm = Udm() - URel;
|
||||
|
||||
Udm_ = Udm();
|
||||
tau_ = betaD*sqr(Udm) + betaC*sqr(Ucm);
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,158 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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::relativeVelocityModel
|
||||
|
||||
Description
|
||||
|
||||
SourceFiles
|
||||
relativeVelocityModel.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef relativeVelocityModel_H
|
||||
#define relativeVelocityModel_H
|
||||
|
||||
#include "fvCFD.H"
|
||||
#include "dictionary.H"
|
||||
#include "incompressibleTwoPhaseMixture.H"
|
||||
#include "runTimeSelectionTables.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class relativeVelocityModel Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class relativeVelocityModel
|
||||
{
|
||||
// Private Member Functions
|
||||
|
||||
//- Disallow default bitwise copy construct
|
||||
relativeVelocityModel(const relativeVelocityModel&);
|
||||
|
||||
//- Disallow default bitwise assignment
|
||||
void operator=(const relativeVelocityModel&);
|
||||
|
||||
protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Mixture properties
|
||||
const incompressibleTwoPhaseMixture& mixture_;
|
||||
|
||||
//- Name of the continuous phase
|
||||
const word continuousPhaseName_;
|
||||
|
||||
//- Continuous phase fraction
|
||||
const volScalarField& alphaC_;
|
||||
|
||||
//- Dispersed phase fraction
|
||||
const volScalarField& alphaD_;
|
||||
|
||||
//- Continuous density
|
||||
const dimensionedScalar& rhoC_;
|
||||
|
||||
//- Dispersed density
|
||||
const dimensionedScalar& rhoD_;
|
||||
|
||||
//- Dispersed diffusion velocity
|
||||
volVectorField Udm_;
|
||||
|
||||
//- Stress
|
||||
volSymmTensorField tau_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("relativeVelocityModel");
|
||||
|
||||
//- Declare runtime constructor selection table
|
||||
declareRunTimeSelectionTable
|
||||
(
|
||||
autoPtr,
|
||||
relativeVelocityModel,
|
||||
dictionary,
|
||||
(const dictionary& dict, const incompressibleTwoPhaseMixture& mixture),
|
||||
(dict, mixture)
|
||||
);
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
relativeVelocityModel
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
);
|
||||
|
||||
|
||||
// Selector
|
||||
static autoPtr<relativeVelocityModel> New
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
virtual ~relativeVelocityModel();
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Calculate the relative velocity of the dispersed phase
|
||||
virtual tmp<volVectorField> Ur() const = 0;
|
||||
|
||||
//- Return the diffusion velocity of the dispersed phase
|
||||
const volVectorField& Udm() const
|
||||
{
|
||||
return Udm_;
|
||||
}
|
||||
|
||||
//- Return the stress tensor due to the phase transport
|
||||
const volSymmTensorField& tau() const
|
||||
{
|
||||
return tau_;
|
||||
}
|
||||
|
||||
//- Update the stored diffusion velocity and stress
|
||||
void update();
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,74 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "simple.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace relativeVelocityModels
|
||||
{
|
||||
defineTypeNameAndDebug(simple, 0);
|
||||
addToRunTimeSelectionTable(relativeVelocityModel, simple, dictionary);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::relativeVelocityModels::simple::simple
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
)
|
||||
:
|
||||
relativeVelocityModel(dict, mixture),
|
||||
a_(dict.lookup("a")),
|
||||
V0_(dict.lookup("V0")),
|
||||
residualAlpha_(dict.lookup("residualAlpha"))
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::relativeVelocityModels::simple::~simple()
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::tmp<Foam::volVectorField>
|
||||
Foam::relativeVelocityModels::simple::Ur() const
|
||||
{
|
||||
return
|
||||
V0_
|
||||
*pow(scalar(10), -a_*max(alphaD_, scalar(0)))
|
||||
/max(alphaC_, residualAlpha_);
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,103 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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::simple
|
||||
|
||||
Description
|
||||
Simple relative velocity model
|
||||
|
||||
SourceFiles
|
||||
simple.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef simple_H
|
||||
#define simple_H
|
||||
|
||||
#include "relativeVelocityModel.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace relativeVelocityModels
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class simple Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class simple
|
||||
:
|
||||
public relativeVelocityModel
|
||||
{
|
||||
// Private data
|
||||
|
||||
//- a coefficient
|
||||
dimensionedScalar a_;
|
||||
|
||||
//- Drift velocity
|
||||
dimensionedVector V0_;
|
||||
|
||||
//- Residual phase fraction
|
||||
dimensionedScalar residualAlpha_;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("simple");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
simple
|
||||
(
|
||||
const dictionary& dict,
|
||||
const incompressibleTwoPhaseMixture& mixture
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
~simple();
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Relative velocity
|
||||
virtual tmp<volVectorField> Ur() const;
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace relativeVelocityModels
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,131 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "BinghamPlastic.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "surfaceFields.H"
|
||||
#include "fvc.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace viscosityModels
|
||||
{
|
||||
defineTypeNameAndDebug(BinghamPlastic, 0);
|
||||
|
||||
addToRunTimeSelectionTable
|
||||
(
|
||||
viscosityModel,
|
||||
BinghamPlastic,
|
||||
dictionary
|
||||
);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
Foam::tmp<Foam::volScalarField>
|
||||
Foam::viscosityModels::BinghamPlastic::correctionNu
|
||||
(
|
||||
const dimensionedScalar& rhoc,
|
||||
const dimensionedScalar& rhop,
|
||||
const volScalarField& nuc
|
||||
) const
|
||||
{
|
||||
volScalarField
|
||||
tauy
|
||||
(
|
||||
yieldStressCoeff_
|
||||
*(
|
||||
pow
|
||||
(
|
||||
scalar(10),
|
||||
yieldStressExponent_
|
||||
*(max(alpha_, scalar(0)) + yieldStressOffset_)
|
||||
)
|
||||
- pow
|
||||
(
|
||||
scalar(10),
|
||||
yieldStressExponent_*yieldStressOffset_
|
||||
)
|
||||
)
|
||||
);
|
||||
|
||||
volScalarField
|
||||
nup
|
||||
(
|
||||
plastic::correctionNu(rhoc, rhop, nuc)
|
||||
);
|
||||
|
||||
dimensionedScalar tauySmall("tauySmall", tauy.dimensions(), SMALL);
|
||||
|
||||
return
|
||||
tauy
|
||||
/(
|
||||
mag(fvc::grad(U_))
|
||||
+ 1.0e-4*(tauy + tauySmall)/(nup + (rhoc/rhop)*nuc)
|
||||
)
|
||||
+ nup;
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::viscosityModels::BinghamPlastic::BinghamPlastic
|
||||
(
|
||||
const word& name,
|
||||
const dictionary& viscosityProperties,
|
||||
const volVectorField& U,
|
||||
const surfaceScalarField& phi
|
||||
)
|
||||
:
|
||||
plastic(name, viscosityProperties, U, phi, typeName),
|
||||
yieldStressCoeff_(plasticCoeffs_.lookup("yieldStressCoeff")),
|
||||
yieldStressExponent_(plasticCoeffs_.lookup("yieldStressExponent")),
|
||||
yieldStressOffset_(plasticCoeffs_.lookup("yieldStressOffset")),
|
||||
U_(U)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
||||
|
||||
bool Foam::viscosityModels::BinghamPlastic::read
|
||||
(
|
||||
const dictionary& viscosityProperties
|
||||
)
|
||||
{
|
||||
plastic::read(viscosityProperties);
|
||||
|
||||
plasticCoeffs_.lookup("yieldStressCoeff") >> yieldStressCoeff_;
|
||||
plasticCoeffs_.lookup("yieldStressExponent") >> yieldStressExponent_;
|
||||
plasticCoeffs_.lookup("yieldStressOffset") >> yieldStressOffset_;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,122 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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::viscosityModels::BinghamPlastic
|
||||
|
||||
Description
|
||||
Viscosity correction model for Bingham plastics.
|
||||
|
||||
SourceFiles
|
||||
BinghamPlastic.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef BinghamPlastic_H
|
||||
#define BinghamPlastic_H
|
||||
|
||||
#include "plastic.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace viscosityModels
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class BinghamPlastic Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class BinghamPlastic
|
||||
:
|
||||
public plastic
|
||||
{
|
||||
protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Yield stress coefficient
|
||||
dimensionedScalar yieldStressCoeff_;
|
||||
|
||||
//- Yield stress exponent
|
||||
dimensionedScalar yieldStressExponent_;
|
||||
|
||||
//- Yield stress offset
|
||||
dimensionedScalar yieldStressOffset_;
|
||||
|
||||
//- Velocity
|
||||
const volVectorField& U_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Calculate and return the laminar viscosity correction
|
||||
virtual tmp<volScalarField> correctionNu
|
||||
(
|
||||
const dimensionedScalar& rhoc,
|
||||
const dimensionedScalar& rhop,
|
||||
const volScalarField& nuc
|
||||
) const;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("BinghamPlastic");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
BinghamPlastic
|
||||
(
|
||||
const word& name,
|
||||
const dictionary& viscosityProperties,
|
||||
const volVectorField& U,
|
||||
const surfaceScalarField& phi
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
~BinghamPlastic()
|
||||
{}
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Read transportProperties dictionary
|
||||
bool read(const dictionary& viscosityProperties);
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace viscosityModels
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,4 @@
|
||||
plastic/plastic.C
|
||||
BinghamPlastic/BinghamPlastic.C
|
||||
|
||||
LIB = $(FOAM_LIBBIN)/libdriftFluxTransportModels
|
@ -0,0 +1,10 @@
|
||||
EXE_INC = \
|
||||
-I$(LIB_SRC)/finiteVolume/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels \
|
||||
-I$(LIB_SRC)/transportModels/twoPhaseMixture/lnInclude \
|
||||
-I$(LIB_SRC)/transportModels/incompressible/lnInclude
|
||||
|
||||
LIB_LIBS = \
|
||||
-ltwoPhaseMixture \
|
||||
-lincompressibleTransportModels \
|
||||
-lfiniteVolume
|
@ -0,0 +1,200 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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/>.
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#include "plastic.H"
|
||||
#include "addToRunTimeSelectionTable.H"
|
||||
#include "surfaceFields.H"
|
||||
#include "incompressibleTwoPhaseMixture.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
namespace viscosityModels
|
||||
{
|
||||
defineTypeNameAndDebug(plastic, 0);
|
||||
|
||||
addToRunTimeSelectionTable
|
||||
(
|
||||
viscosityModel,
|
||||
plastic,
|
||||
dictionary
|
||||
);
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
|
||||
|
||||
Foam::tmp<Foam::volScalarField>
|
||||
Foam::viscosityModels::plastic::calcNu() const
|
||||
{
|
||||
const incompressibleTwoPhaseMixture& twoPhaseProperties =
|
||||
alpha_.mesh().lookupObject<incompressibleTwoPhaseMixture>
|
||||
(
|
||||
"transportProperties"
|
||||
);
|
||||
|
||||
bool isThisIsPhase1(&twoPhaseProperties.nuModel1() == this);
|
||||
|
||||
dimensionedScalar
|
||||
rhoc
|
||||
(
|
||||
isThisIsPhase1
|
||||
? twoPhaseProperties.rho2()
|
||||
: twoPhaseProperties.rho1()
|
||||
);
|
||||
|
||||
dimensionedScalar
|
||||
rhop
|
||||
(
|
||||
isThisIsPhase1
|
||||
? twoPhaseProperties.rho1()
|
||||
: twoPhaseProperties.rho2()
|
||||
);
|
||||
|
||||
volScalarField
|
||||
nuc
|
||||
(
|
||||
(
|
||||
isThisIsPhase1
|
||||
? twoPhaseProperties.nuModel2()
|
||||
: twoPhaseProperties.nuModel1()
|
||||
).nu()
|
||||
);
|
||||
|
||||
volScalarField
|
||||
nup
|
||||
(
|
||||
correctionNu(rhoc, rhop, nuc)
|
||||
);
|
||||
|
||||
return
|
||||
max
|
||||
(
|
||||
nuMin_,
|
||||
min
|
||||
(
|
||||
nuMax_,
|
||||
(
|
||||
nup + (rhoc/rhop)*nuc*alpha_
|
||||
)
|
||||
)
|
||||
)
|
||||
/max(alpha_, SMALL);
|
||||
}
|
||||
|
||||
|
||||
Foam::tmp<Foam::volScalarField>
|
||||
Foam::viscosityModels::plastic::correctionNu
|
||||
(
|
||||
const dimensionedScalar& rhoc,
|
||||
const dimensionedScalar& rhop,
|
||||
const volScalarField& nuc
|
||||
) const
|
||||
{
|
||||
return
|
||||
plasticViscosityCoeff_
|
||||
*(
|
||||
pow
|
||||
(
|
||||
scalar(10),
|
||||
plasticViscosityExponent_*alpha_
|
||||
) - scalar(1)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
||||
|
||||
Foam::viscosityModels::plastic::plastic
|
||||
(
|
||||
const word& name,
|
||||
const dictionary& viscosityProperties,
|
||||
const volVectorField& U,
|
||||
const surfaceScalarField& phi,
|
||||
const word modelName
|
||||
)
|
||||
:
|
||||
viscosityModel(name, viscosityProperties, U, phi),
|
||||
plasticCoeffs_(viscosityProperties.subDict(modelName + "Coeffs")),
|
||||
plasticViscosityCoeff_
|
||||
(
|
||||
plasticCoeffs_.lookup("plasticViscosityCoeff")
|
||||
),
|
||||
plasticViscosityExponent_
|
||||
(
|
||||
plasticCoeffs_.lookup("plasticViscosityExponent")
|
||||
),
|
||||
nuMin_(plasticCoeffs_.lookup("nuMin")),
|
||||
nuMax_(plasticCoeffs_.lookup("nuMax")),
|
||||
alpha_
|
||||
(
|
||||
U.mesh().lookupObject<volScalarField>
|
||||
(
|
||||
IOobject::groupName
|
||||
(
|
||||
viscosityProperties.lookupOrDefault<word>("alpha", "alpha"),
|
||||
viscosityProperties.dictName()
|
||||
)
|
||||
)
|
||||
),
|
||||
nu_
|
||||
(
|
||||
IOobject
|
||||
(
|
||||
name,
|
||||
U_.time().timeName(),
|
||||
U_.db(),
|
||||
IOobject::NO_READ,
|
||||
IOobject::AUTO_WRITE
|
||||
),
|
||||
U_.mesh(),
|
||||
dimensionedScalar("nu", dimViscosity, 0)
|
||||
)
|
||||
{}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
||||
|
||||
bool Foam::viscosityModels::plastic::read
|
||||
(
|
||||
const dictionary& viscosityProperties
|
||||
)
|
||||
{
|
||||
viscosityModel::read(viscosityProperties);
|
||||
|
||||
plasticCoeffs_ = viscosityProperties.subDict(typeName + "Coeffs");
|
||||
|
||||
plasticCoeffs_.lookup("k") >> plasticViscosityCoeff_;
|
||||
plasticCoeffs_.lookup("n") >> plasticViscosityExponent_;
|
||||
plasticCoeffs_.lookup("nuMin") >> nuMin_;
|
||||
plasticCoeffs_.lookup("nuMax") >> nuMax_;
|
||||
|
||||
return true;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
@ -0,0 +1,158 @@
|
||||
/*---------------------------------------------------------------------------*\
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2014 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::viscosityModels::plastic
|
||||
|
||||
Description
|
||||
Viscosity correction model for a generic power-law plastic.
|
||||
|
||||
SourceFiles
|
||||
plastic.C
|
||||
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
#ifndef plastic_H
|
||||
#define plastic_H
|
||||
|
||||
#include "viscosityModel.H"
|
||||
#include "dimensionedScalar.H"
|
||||
#include "volFields.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
|
||||
class incompressibleTwoPhaseMixture;
|
||||
|
||||
namespace viscosityModels
|
||||
{
|
||||
|
||||
/*---------------------------------------------------------------------------*\
|
||||
Class plastic Declaration
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
||||
class plastic
|
||||
:
|
||||
public viscosityModel
|
||||
{
|
||||
protected:
|
||||
|
||||
// Protected data
|
||||
|
||||
//- Dictionary
|
||||
dictionary plasticCoeffs_;
|
||||
|
||||
//- Plastic viscosity coefficient
|
||||
dimensionedScalar plasticViscosityCoeff_;
|
||||
|
||||
//- Plastic viscosity exponent
|
||||
dimensionedScalar plasticViscosityExponent_;
|
||||
|
||||
//- Minimum viscosity
|
||||
dimensionedScalar nuMin_;
|
||||
|
||||
//- Maximum viscosity
|
||||
dimensionedScalar nuMax_;
|
||||
|
||||
//- Plastic phase fraction
|
||||
const volScalarField& alpha_;
|
||||
|
||||
//- Viscosity
|
||||
volScalarField nu_;
|
||||
|
||||
|
||||
// Protected Member Functions
|
||||
|
||||
//- Calculate and return the laminar viscosity
|
||||
virtual tmp<volScalarField> calcNu() const;
|
||||
|
||||
//- Calculate and return the laminar viscosity correction
|
||||
virtual tmp<volScalarField> correctionNu
|
||||
(
|
||||
const dimensionedScalar& rhoc,
|
||||
const dimensionedScalar& rhop,
|
||||
const volScalarField& nuc
|
||||
) const;
|
||||
|
||||
|
||||
public:
|
||||
|
||||
//- Runtime type information
|
||||
TypeName("plastic");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
plastic
|
||||
(
|
||||
const word& name,
|
||||
const dictionary& viscosityProperties,
|
||||
const volVectorField& U,
|
||||
const surfaceScalarField& phi,
|
||||
const word modelName=typeName
|
||||
);
|
||||
|
||||
|
||||
//- Destructor
|
||||
~plastic()
|
||||
{}
|
||||
|
||||
|
||||
// Member Functions
|
||||
|
||||
//- Return the laminar viscosity
|
||||
tmp<volScalarField> nu() const
|
||||
{
|
||||
return nu_;
|
||||
}
|
||||
|
||||
//- Return the laminar viscosity for patch
|
||||
tmp<scalarField> nu(const label patchi) const
|
||||
{
|
||||
return nu_.boundaryField()[patchi];
|
||||
}
|
||||
|
||||
//- Correct the laminar viscosity
|
||||
void correct()
|
||||
{
|
||||
nu_ = calcNu();
|
||||
}
|
||||
|
||||
//- Read transportProperties dictionary
|
||||
bool read(const dictionary& viscosityProperties);
|
||||
};
|
||||
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
} // End namespace viscosityModels
|
||||
} // End namespace Foam
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
#endif
|
||||
|
||||
// ************************************************************************* //
|
@ -4,7 +4,6 @@
|
||||
const scalar Cmu25 = ::pow(Cmu.value(), 0.25);
|
||||
const scalar Cmu75 = ::pow(Cmu.value(), 0.75);
|
||||
const scalar kappa_ = kappa.value();
|
||||
const scalar muc_ = muc.value();
|
||||
|
||||
const fvPatchList& patches = mesh.boundary();
|
||||
|
||||
@ -34,6 +33,7 @@
|
||||
if (isA<wallFvPatch>(curPatch))
|
||||
{
|
||||
const scalarField& mutw = mut.boundaryField()[patchi];
|
||||
const scalarField& mucw = muc.boundaryField()[patchi];
|
||||
|
||||
scalarField magFaceGradU
|
||||
(
|
||||
@ -55,7 +55,7 @@
|
||||
/(kappa_*y[patchi][facei]);
|
||||
|
||||
G[faceCelli] +=
|
||||
(mutw[facei] + muc_)
|
||||
(mutw[facei] + mucw[facei])
|
||||
*magFaceGradU[facei]
|
||||
*Cmu25*::sqrt(k[faceCelli])
|
||||
/(kappa_*y[patchi][facei]);
|
||||
|
@ -2,8 +2,6 @@
|
||||
const scalar Cmu25 = ::pow(Cmu.value(), 0.25);
|
||||
const scalar kappa_ = kappa.value();
|
||||
const scalar E_ = E.value();
|
||||
const scalar muc_ = muc.value();
|
||||
const scalar nuc_ = muc_/rho2.value();
|
||||
|
||||
const fvPatchList& patches = mesh.boundary();
|
||||
|
||||
@ -14,6 +12,7 @@
|
||||
if (isA<wallFvPatch>(curPatch))
|
||||
{
|
||||
scalarField& mutw = mut.boundaryField()[patchi];
|
||||
const scalarField& mucw = muc.boundaryField()[patchi];
|
||||
|
||||
forAll(curPatch, facei)
|
||||
{
|
||||
@ -21,12 +20,12 @@
|
||||
|
||||
scalar yPlus =
|
||||
Cmu25*y[patchi][facei]*::sqrt(k[faceCelli])
|
||||
/nuc_;
|
||||
/(mucw[facei]/rho2.value());
|
||||
|
||||
if (yPlus > 11.6)
|
||||
{
|
||||
mutw[facei] =
|
||||
muc_*(yPlus*kappa_/::log(E_*yPlus) - 1);
|
||||
mucw[facei]*(yPlus*kappa_/::log(E_*yPlus) - 1);
|
||||
}
|
||||
else
|
||||
{
|
||||
|
@ -1,27 +0,0 @@
|
||||
volScalarField yieldStress
|
||||
(
|
||||
const dimensionedScalar& yieldStressCoeff,
|
||||
const dimensionedScalar& yieldStressExponent,
|
||||
const dimensionedScalar& yieldStressOffset,
|
||||
const volScalarField& alpha1
|
||||
)
|
||||
{
|
||||
tmp<volScalarField> tfld
|
||||
(
|
||||
yieldStressCoeff*
|
||||
(
|
||||
pow
|
||||
(
|
||||
10.0,
|
||||
yieldStressExponent*(max(alpha1, scalar(0)) + yieldStressOffset)
|
||||
)
|
||||
- pow
|
||||
(
|
||||
10.0,
|
||||
yieldStressExponent*yieldStressOffset
|
||||
)
|
||||
)
|
||||
);
|
||||
|
||||
return tfld();
|
||||
}
|
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -28,6 +28,13 @@ License
|
||||
#include "surfaceFields.H"
|
||||
#include "fvc.H"
|
||||
|
||||
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
|
||||
|
||||
namespace Foam
|
||||
{
|
||||
defineTypeNameAndDebug(incompressibleTwoPhaseMixture, 0);
|
||||
}
|
||||
|
||||
|
||||
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
|
||||
|
||||
|
@ -2,7 +2,7 @@
|
||||
========= |
|
||||
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
||||
\\ / O peration |
|
||||
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
|
||||
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
|
||||
\\/ M anipulation |
|
||||
-------------------------------------------------------------------------------
|
||||
License
|
||||
@ -80,6 +80,9 @@ protected:
|
||||
|
||||
public:
|
||||
|
||||
TypeName("incompressibleTwoPhaseMixture");
|
||||
|
||||
|
||||
// Constructors
|
||||
|
||||
//- Construct from components
|
||||
@ -121,6 +124,12 @@ public:
|
||||
return rho2_;
|
||||
};
|
||||
|
||||
//- Return const-access to the mixture velocity
|
||||
const volVectorField& U() const
|
||||
{
|
||||
return U_;
|
||||
}
|
||||
|
||||
//- Return the dynamic laminar viscosity
|
||||
tmp<volScalarField> mu() const;
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -10,7 +10,7 @@ FoamFile
|
||||
version 2.0;
|
||||
format ascii;
|
||||
class volScalarField;
|
||||
object alpha1;
|
||||
object alpha.sludge;
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -15,41 +15,45 @@ FoamFile
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
muc muc [ 1 -1 -1 0 0 0 0 ] 0.00178;
|
||||
phases (sludge water);
|
||||
|
||||
plasticViscosityCoeff plasticViscosityCoeff [ 1 -1 -1 0 0 0 0 ] 0.00023143;
|
||||
|
||||
plasticViscosityExponent plasticViscosityExponent [ 0 0 0 0 0 0 0 ] 179.26;
|
||||
|
||||
BinghamPlastic on;
|
||||
|
||||
yieldStressCoeff yieldStressCoeff [ 1 -1 -2 0 0 0 0 ] 0.00042189;
|
||||
|
||||
yieldStressExponent yieldStressExponent [ 0 0 0 0 0 0 0 ] 1050.8;
|
||||
|
||||
yieldStressOffset yieldStressOffset [ 0 0 0 0 0 0 0 ] 0;
|
||||
|
||||
muMax muMax [ 1 -1 -1 0 0 0 0 ] 10.0;
|
||||
|
||||
rho1 rho1 [ 1 -3 0 0 0 0 0 ] 1996;
|
||||
rho2 rho2 [ 1 -3 0 0 0 0 0 ] 996;
|
||||
|
||||
VdjModel simple;
|
||||
|
||||
simpleCoeffs
|
||||
sludge
|
||||
{
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.002198 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 285.84;
|
||||
a1 a1 [ 0 0 0 0 0 0 0 ] 0;
|
||||
alphaMin alphaMin [ 0 0 0 0 0 0 0 ] 0;
|
||||
transportModel BinghamPlastic;
|
||||
|
||||
"(plastic|BinghamPlastic)Coeffs"
|
||||
{
|
||||
plasticViscosityCoeff plasticViscosityCoeff [ 0 2 -1 0 0 0 0 ] 1.1595e-07;
|
||||
plasticViscosityExponent plasticViscosityExponent [ 0 0 0 0 0 0 0 ] 179.26;
|
||||
|
||||
yieldStressCoeff yieldStressCoeff [ 0 2 -2 0 0 0 0 ] 2.1137e-07;
|
||||
yieldStressExponent yieldStressExponent [ 0 0 0 0 0 0 0 ] 1050.8;
|
||||
yieldStressOffset yieldStressOffset [ 0 0 0 0 0 0 0 ] 0;
|
||||
|
||||
nuMin nuMin [ 0 2 -1 0 0 0 0 ] 1e-10;
|
||||
nuMax nuMax [ 0 2 -1 0 0 0 0 ] 5e-3;
|
||||
}
|
||||
|
||||
rho rho [ 1 -3 0 0 0 0 0 ] 1996;
|
||||
}
|
||||
|
||||
generalCoeffs
|
||||
water
|
||||
{
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.0018 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 1e-05;
|
||||
transportModel Newtonian;
|
||||
|
||||
nu nu [ 0 2 -1 0 0 0 0 ] 1.7871e-06;
|
||||
rho rho [ 1 -3 0 0 0 0 0 ] 996;
|
||||
}
|
||||
|
||||
relativeVelocityModel simple;
|
||||
|
||||
"(simple|general)Coeffs"
|
||||
{
|
||||
continuousPhase water;
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.002198 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 285.84;
|
||||
a1 a1 [ 0 0 0 0 0 0 0 ] 0.1;
|
||||
alphaMin alphaMin [ 0 0 0 0 0 0 0 ] 2e-05;
|
||||
residualAlpha residualAlpha [ 0 0 0 0 0 0 0 ] 0;
|
||||
}
|
||||
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -30,9 +30,9 @@ divSchemes
|
||||
default none;
|
||||
|
||||
div(rhoPhi,U) Gauss linearUpwind grad(U);
|
||||
div(phiVdj,Vdj) Gauss linear;
|
||||
div(phi,alpha) Gauss vanLeer;
|
||||
div(phirb,alpha) Gauss linear;
|
||||
div(phiUkm,Ukm) Gauss linear;
|
||||
"div\(phi,alpha.*\)" Gauss vanLeer;
|
||||
"div\(phirb,alpha.*\)" Gauss linear;
|
||||
div(rhoPhi,k) Gauss limitedLinear 1;
|
||||
div(rhoPhi,epsilon) Gauss limitedLinear 1;
|
||||
|
||||
@ -58,7 +58,7 @@ fluxRequired
|
||||
{
|
||||
default no;
|
||||
p_rgh;
|
||||
alpha1;
|
||||
"alpha.*";
|
||||
}
|
||||
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -17,7 +17,7 @@ FoamFile
|
||||
|
||||
solvers
|
||||
{
|
||||
"alpha1.*"
|
||||
"alpha.*"
|
||||
{
|
||||
nAlphaCorr 2;
|
||||
nAlphaSubCycles 1;
|
||||
@ -33,7 +33,7 @@ solvers
|
||||
minIter 1;
|
||||
}
|
||||
|
||||
alpha1Diffusion
|
||||
"alpha.*Diffusion"
|
||||
{
|
||||
solver PCG;
|
||||
preconditioner DIC;
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -10,7 +10,7 @@ FoamFile
|
||||
version 2.0;
|
||||
format ascii;
|
||||
class volScalarField;
|
||||
object alpha1;
|
||||
object alpha.sludge;
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -15,41 +15,45 @@ FoamFile
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
muc muc [ 1 -1 -1 0 0 0 0 ] 0.00178;
|
||||
phases (sludge water);
|
||||
|
||||
plasticViscosityCoeff plasticViscosityCoeff [ 1 -1 -1 0 0 0 0 ] 0.00023143;
|
||||
|
||||
plasticViscosityExponent plasticViscosityExponent [ 0 0 0 0 0 0 0 ] 179.26;
|
||||
|
||||
BinghamPlastic on;
|
||||
|
||||
yieldStressCoeff yieldStressCoeff [ 1 -1 -2 0 0 0 0 ] 0.00042189;
|
||||
|
||||
yieldStressExponent yieldStressExponent [ 0 0 0 0 0 0 0 ] 1050.8;
|
||||
|
||||
yieldStressOffset yieldStressOffset [ 0 0 0 0 0 0 0 ] 0;
|
||||
|
||||
muMax muMax [ 1 -1 -1 0 0 0 0 ] 10.0;
|
||||
|
||||
rho1 rho1 [ 1 -3 0 0 0 0 0 ] 1996;
|
||||
rho2 rho2 [ 1 -3 0 0 0 0 0 ] 996;
|
||||
|
||||
VdjModel simple;
|
||||
|
||||
simpleCoeffs
|
||||
sludge
|
||||
{
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.002198 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 285.84;
|
||||
a1 a1 [ 0 0 0 0 0 0 0 ] 0;
|
||||
alphaMin alphaMin [ 0 0 0 0 0 0 0 ] 0;
|
||||
transportModel BinghamPlastic;
|
||||
|
||||
"(plastic|BinghamPlastic)Coeffs"
|
||||
{
|
||||
plasticViscosityCoeff plasticViscosityCoeff [ 0 2 -1 0 0 0 0 ] 1.1595e-07;
|
||||
plasticViscosityExponent plasticViscosityExponent [ 0 0 0 0 0 0 0 ] 179.26;
|
||||
|
||||
yieldStressCoeff yieldStressCoeff [ 0 2 -2 0 0 0 0 ] 2.1137e-07;
|
||||
yieldStressExponent yieldStressExponent [ 0 0 0 0 0 0 0 ] 1050.8;
|
||||
yieldStressOffset yieldStressOffset [ 0 0 0 0 0 0 0 ] 0;
|
||||
|
||||
nuMin nuMin [ 0 2 -1 0 0 0 0 ] 1e-10;
|
||||
nuMax nuMax [ 0 2 -1 0 0 0 0 ] 5e-3;
|
||||
}
|
||||
|
||||
rho rho [ 1 -3 0 0 0 0 0 ] 1996;
|
||||
}
|
||||
|
||||
generalCoeffs
|
||||
water
|
||||
{
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.0018 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 1e-05;
|
||||
transportModel Newtonian;
|
||||
|
||||
nu nu [ 0 2 -1 0 0 0 0 ] 1.7871e-06;
|
||||
rho rho [ 1 -3 0 0 0 0 0 ] 996;
|
||||
}
|
||||
|
||||
relativeVelocityModel simple;
|
||||
|
||||
"(simple|general)Coeffs"
|
||||
{
|
||||
continuousPhase water;
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.002198 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 285.84;
|
||||
a1 a1 [ 0 0 0 0 0 0 0 ] 0.1;
|
||||
alphaMin alphaMin [ 0 0 0 0 0 0 0 ] 2e-05;
|
||||
residualAlpha residualAlpha [ 0 0 0 0 0 0 0 ] 0;
|
||||
}
|
||||
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -30,9 +30,9 @@ divSchemes
|
||||
default none;
|
||||
|
||||
div(rhoPhi,U) Gauss linearUpwind grad(U);
|
||||
div(phiVdj,Vdj) Gauss linear;
|
||||
div(phi,alpha) Gauss vanLeer;
|
||||
div(phirb,alpha) Gauss linear;
|
||||
div(phiUkm,Ukm) Gauss linear;
|
||||
"div\(phi,alpha.*\)" Gauss vanLeer;
|
||||
"div\(phirb,alpha.*\)" Gauss linear;
|
||||
div(rhoPhi,k) Gauss limitedLinear 1;
|
||||
div(rhoPhi,epsilon) Gauss limitedLinear 1;
|
||||
|
||||
@ -58,7 +58,7 @@ fluxRequired
|
||||
{
|
||||
default no;
|
||||
p_rgh;
|
||||
alpha1;
|
||||
"alpha.*";
|
||||
}
|
||||
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -17,7 +17,7 @@ FoamFile
|
||||
|
||||
solvers
|
||||
{
|
||||
"alpha1.*"
|
||||
"alpha.*"
|
||||
{
|
||||
nAlphaCorr 2;
|
||||
nAlphaSubCycles 1;
|
||||
@ -33,7 +33,7 @@ solvers
|
||||
minIter 1;
|
||||
}
|
||||
|
||||
alpha1Diffusion
|
||||
"alpha.*Diffusion"
|
||||
{
|
||||
solver PCG;
|
||||
preconditioner DIC;
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -10,7 +10,7 @@ FoamFile
|
||||
version 2.0;
|
||||
format ascii;
|
||||
class volScalarField;
|
||||
object alpha1;
|
||||
object alpha.sludge;
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -15,41 +15,45 @@ FoamFile
|
||||
}
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
muc muc [ 1 -1 -1 0 0 0 0 ] 0.00178;
|
||||
phases (sludge water);
|
||||
|
||||
plasticViscosityCoeff plasticViscosityCoeff [ 1 -1 -1 0 0 0 0 ] 0.00023143;
|
||||
|
||||
plasticViscosityExponent plasticViscosityExponent [ 0 0 0 0 0 0 0 ] 0.17926;
|
||||
|
||||
BinghamPlastic on;
|
||||
|
||||
yieldStressCoeff yieldStressCoeff [ 1 -1 -2 0 0 0 0 ] 5.5469e-07;
|
||||
|
||||
yieldStressExponent yieldStressExponent [ 0 0 0 0 0 0 0 ] 95.25;
|
||||
|
||||
yieldStressOffset yieldStressOffset [ 0 0 0 0 0 0 0 ] 0;
|
||||
|
||||
muMax muMax [ 1 -1 -1 0 0 0 0 ] 10.0;
|
||||
|
||||
rho1 rho1 [ 1 -3 0 0 0 0 0 ] 1042;
|
||||
rho2 rho2 [ 1 -3 0 0 0 0 0 ] 1000;
|
||||
|
||||
VdjModel simple;
|
||||
|
||||
simpleCoeffs
|
||||
sludge
|
||||
{
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.002198 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 8.84;
|
||||
a1 a1 [ 0 0 0 0 0 0 0 ] 0;
|
||||
alphaMin alphaMin [ 0 0 0 0 0 0 0 ] 0;
|
||||
transportModel BinghamPlastic;
|
||||
|
||||
"(plastic|BinghamPlastic)Coeffs"
|
||||
{
|
||||
plasticViscosityCoeff plasticViscosityCoeff [ 0 2 -1 0 0 0 0 ] 2.2210e-07;
|
||||
plasticViscosityExponent plasticViscosityExponent [ 0 0 0 0 0 0 0 ] 0.17926;
|
||||
|
||||
yieldStressCoeff yieldStressCoeff [ 0 2 -2 0 0 0 0 ] 5.3233e-10;
|
||||
yieldStressExponent yieldStressExponent [ 0 0 0 0 0 0 0 ] 95.25;
|
||||
yieldStressOffset yieldStressOffset [ 0 0 0 0 0 0 0 ] 0;
|
||||
|
||||
nuMin nuMin [ 0 2 -1 0 0 0 0 ] 1e-10;
|
||||
nuMax nuMax [ 0 2 -1 0 0 0 0 ] 5e-3;
|
||||
}
|
||||
|
||||
rho rho [ 1 -3 0 0 0 0 0 ] 1042;
|
||||
}
|
||||
|
||||
generalCoeffs
|
||||
water
|
||||
{
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.0018 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 1e-05;
|
||||
transportModel Newtonian;
|
||||
|
||||
nu nu [ 0 2 -1 0 0 0 0 ] 1.78e-06;
|
||||
rho rho [ 1 -3 0 0 0 0 0 ] 1000;
|
||||
}
|
||||
|
||||
relativeVelocityModel simple;
|
||||
|
||||
"(simple|general)Coeffs"
|
||||
{
|
||||
continuousPhase water;
|
||||
V0 V0 [ 0 1 -1 0 0 0 0 ] ( 0 -0.002198 0 );
|
||||
a a [ 0 0 0 0 0 0 0 ] 8.84;
|
||||
a1 a1 [ 0 0 0 0 0 0 0 ] 0.1;
|
||||
alphaMin alphaMin [ 0 0 0 0 0 0 0 ] 2e-05;
|
||||
residualAlpha residualAlpha [ 0 0 0 0 0 0 0 ] 0;
|
||||
}
|
||||
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -30,9 +30,9 @@ divSchemes
|
||||
default none;
|
||||
|
||||
div(rhoPhi,U) Gauss linearUpwind grad(U);
|
||||
div(phiVdj,Vdj) Gauss linear;
|
||||
div(phi,alpha) Gauss vanLeer;
|
||||
div(phirb,alpha) Gauss linear;
|
||||
div(phiUkm,Ukm) Gauss linear;
|
||||
"div\(phi,alpha.*\)" Gauss vanLeer;
|
||||
"div\(phirb,alpha.*\)" Gauss linear;
|
||||
div(rhoPhi,k) Gauss limitedLinear 1;
|
||||
div(rhoPhi,epsilon) Gauss limitedLinear 1;
|
||||
|
||||
@ -58,7 +58,7 @@ fluxRequired
|
||||
{
|
||||
default no;
|
||||
p_rgh;
|
||||
alpha1;
|
||||
"alpha.*";
|
||||
}
|
||||
|
||||
|
||||
|
@ -1,7 +1,7 @@
|
||||
/*--------------------------------*- C++ -*----------------------------------*\
|
||||
| ========= | |
|
||||
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
|
||||
| \\ / O peration | Version: dev |
|
||||
| \\ / O peration | Version: 2.3.x |
|
||||
| \\ / A nd | Web: www.OpenFOAM.org |
|
||||
| \\/ M anipulation | |
|
||||
\*---------------------------------------------------------------------------*/
|
||||
@ -17,7 +17,7 @@ FoamFile
|
||||
|
||||
solvers
|
||||
{
|
||||
"alpha1.*"
|
||||
"alpha.*"
|
||||
{
|
||||
nAlphaCorr 2;
|
||||
nAlphaSubCycles 1;
|
||||
@ -33,7 +33,7 @@ solvers
|
||||
minIter 1;
|
||||
}
|
||||
|
||||
alpha1Diffusion
|
||||
"alpha.*Diffusion"
|
||||
{
|
||||
solver PCG;
|
||||
preconditioner DIC;
|
||||
|
Loading…
Reference in New Issue
Block a user