mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
PDRFoam: Added support for fvOptions
This commit is contained in:
@ -1,7 +1,7 @@
|
||||
{
|
||||
volScalarField& hea = thermo.he();
|
||||
|
||||
solve
|
||||
fvScalarMatrix EaEqn
|
||||
(
|
||||
betav*fvm::ddt(rho, hea) + mvConvection->fvmDiv(phi, hea)
|
||||
+ betav*fvc::ddt(rho, K) + fvc::div(phi, K)
|
||||
@ -16,7 +16,16 @@
|
||||
: -betav*dpdt
|
||||
)
|
||||
- fvm::laplacian(Db, hea)
|
||||
+ betav*fvOptions(rho, hea)
|
||||
);
|
||||
|
||||
EaEqn.relax();
|
||||
|
||||
fvOptions.constrain(EaEqn);
|
||||
|
||||
EaEqn.solve();
|
||||
|
||||
fvOptions.correct(hea);
|
||||
|
||||
thermo.correct();
|
||||
}
|
||||
|
||||
@ -2,7 +2,7 @@ if (ign.ignited())
|
||||
{
|
||||
volScalarField& heau = thermo.heu();
|
||||
|
||||
solve
|
||||
fvScalarMatrix heauEqn
|
||||
(
|
||||
betav*fvm::ddt(rho, heau) + mvConvection->fvmDiv(phi, heau)
|
||||
+ (betav*fvc::ddt(rho, K) + fvc::div(phi, K))*rho/thermo.rhou()
|
||||
@ -23,5 +23,13 @@ if (ign.ignited())
|
||||
// A possible solution would be to solve for ftu as well as ft.
|
||||
//- fvm::div(muEff*fvc::grad(b)/(b + 0.001), heau)
|
||||
//+ fvm::Sp(fvc::div(muEff*fvc::grad(b)/(b + 0.001)), heau)
|
||||
==
|
||||
betav*fvOptions(rho, heau)
|
||||
);
|
||||
|
||||
fvOptions.constrain(heauEqn);
|
||||
|
||||
heauEqn.solve();
|
||||
|
||||
fvOptions.correct(heau);
|
||||
}
|
||||
|
||||
@ -77,6 +77,7 @@ Description
|
||||
#include "Switch.H"
|
||||
#include "bound.H"
|
||||
#include "pimpleControl.H"
|
||||
#include "fvIOoptionList.H"
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -93,11 +94,13 @@ int main(int argc, char *argv[])
|
||||
#include "readGravitationalAcceleration.H"
|
||||
#include "createFields.H"
|
||||
#include "createMRF.H"
|
||||
#include "createFvOptions.H"
|
||||
#include "initContinuityErrs.H"
|
||||
#include "createTimeControls.H"
|
||||
#include "compressibleCourantNo.H"
|
||||
#include "setInitialDeltaT.H"
|
||||
|
||||
turbulence->validate();
|
||||
scalar StCoNum = 0.0;
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -87,6 +87,7 @@ int main(int argc, char *argv[])
|
||||
#include "compressibleCourantNo.H"
|
||||
#include "setInitialDeltaT.H"
|
||||
|
||||
turbulence->validate();
|
||||
scalar StCoNum = 0.0;
|
||||
|
||||
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
||||
|
||||
@ -7,13 +7,17 @@
|
||||
+ turbulence->divDevRhoReff(U)
|
||||
==
|
||||
betav*rho*g
|
||||
+ betav*fvOptions(rho, U)
|
||||
);
|
||||
|
||||
fvOptions.constrain(UEqn);
|
||||
|
||||
volSymmTensorField invA(inv(I*UEqn.A() + drag->Dcu()));
|
||||
|
||||
if (pimple.momentumPredictor())
|
||||
{
|
||||
U = invA & (UEqn.H() - betav*fvc::grad(p));
|
||||
U.correctBoundaryConditions();
|
||||
fvOptions.correct(U);
|
||||
K = 0.5*magSqr(U);
|
||||
}
|
||||
|
||||
@ -73,6 +73,8 @@ if (ign.ignited())
|
||||
+ fvm::div(phiSt, b)
|
||||
- fvm::Sp(fvc::div(phiSt), b)
|
||||
- fvm::laplacian(Db, b)
|
||||
==
|
||||
betav*fvOptions(rho, b)
|
||||
);
|
||||
|
||||
|
||||
@ -82,8 +84,14 @@ if (ign.ignited())
|
||||
|
||||
// Solve for b
|
||||
// ~~~~~~~~~~~
|
||||
bEqn.relax();
|
||||
|
||||
fvOptions.constrain(bEqn);
|
||||
|
||||
bEqn.solve();
|
||||
|
||||
fvOptions.correct(b);
|
||||
|
||||
Info<< "min(b) = " << min(b).value() << endl;
|
||||
|
||||
if (composition.contains("ft"))
|
||||
|
||||
@ -25,6 +25,8 @@ if (pimple.transonic())
|
||||
betav*fvm::ddt(psi, p)
|
||||
+ fvm::div(phid, p)
|
||||
- fvm::laplacian(rho*invA, p)
|
||||
==
|
||||
betav*fvOptions(psi, p, rho.name())
|
||||
);
|
||||
|
||||
pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
|
||||
@ -53,6 +55,8 @@ else
|
||||
betav*fvm::ddt(psi, p)
|
||||
+ fvc::div(phiHbyA)
|
||||
- fvm::laplacian(rho*invA, p)
|
||||
==
|
||||
betav*fvOptions(psi, p, rho.name())
|
||||
);
|
||||
|
||||
pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
|
||||
@ -69,6 +73,7 @@ else
|
||||
|
||||
U = HbyA - (invA & (betav*fvc::grad(p)));
|
||||
U.correctBoundaryConditions();
|
||||
fvOptions.correct(U);
|
||||
K = 0.5*magSqr(U);
|
||||
|
||||
if (thermo.dpdt())
|
||||
|
||||
Reference in New Issue
Block a user