Skip to content
Snippets Groups Projects
Commit 89955507 authored by Henry Weller's avatar Henry Weller
Browse files

PDRFoam: Added support for fvOptions

parent b8ffb82b
Branches
Tags
No related merge requests found
{ {
volScalarField& hea = thermo.he(); volScalarField& hea = thermo.he();
solve fvScalarMatrix EaEqn
( (
betav*fvm::ddt(rho, hea) + mvConvection->fvmDiv(phi, hea) betav*fvm::ddt(rho, hea) + mvConvection->fvmDiv(phi, hea)
+ betav*fvc::ddt(rho, K) + fvc::div(phi, K) + betav*fvc::ddt(rho, K) + fvc::div(phi, K)
...@@ -16,7 +16,16 @@ ...@@ -16,7 +16,16 @@
: -betav*dpdt : -betav*dpdt
) )
- fvm::laplacian(Db, hea) - fvm::laplacian(Db, hea)
+ betav*fvOptions(rho, hea)
); );
EaEqn.relax();
fvOptions.constrain(EaEqn);
EaEqn.solve();
fvOptions.correct(hea);
thermo.correct(); thermo.correct();
} }
...@@ -2,7 +2,7 @@ if (ign.ignited()) ...@@ -2,7 +2,7 @@ if (ign.ignited())
{ {
volScalarField& heau = thermo.heu(); volScalarField& heau = thermo.heu();
solve fvScalarMatrix heauEqn
( (
betav*fvm::ddt(rho, heau) + mvConvection->fvmDiv(phi, heau) betav*fvm::ddt(rho, heau) + mvConvection->fvmDiv(phi, heau)
+ (betav*fvc::ddt(rho, K) + fvc::div(phi, K))*rho/thermo.rhou() + (betav*fvc::ddt(rho, K) + fvc::div(phi, K))*rho/thermo.rhou()
...@@ -23,5 +23,13 @@ if (ign.ignited()) ...@@ -23,5 +23,13 @@ if (ign.ignited())
// A possible solution would be to solve for ftu as well as ft. // A possible solution would be to solve for ftu as well as ft.
//- fvm::div(muEff*fvc::grad(b)/(b + 0.001), heau) //- fvm::div(muEff*fvc::grad(b)/(b + 0.001), heau)
//+ fvm::Sp(fvc::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 ...@@ -77,6 +77,7 @@ Description
#include "Switch.H" #include "Switch.H"
#include "bound.H" #include "bound.H"
#include "pimpleControl.H" #include "pimpleControl.H"
#include "fvIOoptionList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -93,11 +94,13 @@ int main(int argc, char *argv[]) ...@@ -93,11 +94,13 @@ int main(int argc, char *argv[])
#include "readGravitationalAcceleration.H" #include "readGravitationalAcceleration.H"
#include "createFields.H" #include "createFields.H"
#include "createMRF.H" #include "createMRF.H"
#include "createFvOptions.H"
#include "initContinuityErrs.H" #include "initContinuityErrs.H"
#include "createTimeControls.H" #include "createTimeControls.H"
#include "compressibleCourantNo.H" #include "compressibleCourantNo.H"
#include "setInitialDeltaT.H" #include "setInitialDeltaT.H"
turbulence->validate();
scalar StCoNum = 0.0; scalar StCoNum = 0.0;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
......
...@@ -87,6 +87,7 @@ int main(int argc, char *argv[]) ...@@ -87,6 +87,7 @@ int main(int argc, char *argv[])
#include "compressibleCourantNo.H" #include "compressibleCourantNo.H"
#include "setInitialDeltaT.H" #include "setInitialDeltaT.H"
turbulence->validate();
scalar StCoNum = 0.0; scalar StCoNum = 0.0;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
......
...@@ -7,13 +7,17 @@ ...@@ -7,13 +7,17 @@
+ turbulence->divDevRhoReff(U) + turbulence->divDevRhoReff(U)
== ==
betav*rho*g betav*rho*g
+ betav*fvOptions(rho, U)
); );
fvOptions.constrain(UEqn);
volSymmTensorField invA(inv(I*UEqn.A() + drag->Dcu())); volSymmTensorField invA(inv(I*UEqn.A() + drag->Dcu()));
if (pimple.momentumPredictor()) if (pimple.momentumPredictor())
{ {
U = invA & (UEqn.H() - betav*fvc::grad(p)); U = invA & (UEqn.H() - betav*fvc::grad(p));
U.correctBoundaryConditions(); U.correctBoundaryConditions();
fvOptions.correct(U);
K = 0.5*magSqr(U); K = 0.5*magSqr(U);
} }
...@@ -73,6 +73,8 @@ if (ign.ignited()) ...@@ -73,6 +73,8 @@ if (ign.ignited())
+ fvm::div(phiSt, b) + fvm::div(phiSt, b)
- fvm::Sp(fvc::div(phiSt), b) - fvm::Sp(fvc::div(phiSt), b)
- fvm::laplacian(Db, b) - fvm::laplacian(Db, b)
==
betav*fvOptions(rho, b)
); );
...@@ -82,8 +84,14 @@ if (ign.ignited()) ...@@ -82,8 +84,14 @@ if (ign.ignited())
// Solve for b // Solve for b
// ~~~~~~~~~~~ // ~~~~~~~~~~~
bEqn.relax();
fvOptions.constrain(bEqn);
bEqn.solve(); bEqn.solve();
fvOptions.correct(b);
Info<< "min(b) = " << min(b).value() << endl; Info<< "min(b) = " << min(b).value() << endl;
if (composition.contains("ft")) if (composition.contains("ft"))
......
...@@ -25,6 +25,8 @@ if (pimple.transonic()) ...@@ -25,6 +25,8 @@ if (pimple.transonic())
betav*fvm::ddt(psi, p) betav*fvm::ddt(psi, p)
+ fvm::div(phid, p) + fvm::div(phid, p)
- fvm::laplacian(rho*invA, p) - fvm::laplacian(rho*invA, p)
==
betav*fvOptions(psi, p, rho.name())
); );
pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter()))); pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
...@@ -53,6 +55,8 @@ else ...@@ -53,6 +55,8 @@ else
betav*fvm::ddt(psi, p) betav*fvm::ddt(psi, p)
+ fvc::div(phiHbyA) + fvc::div(phiHbyA)
- fvm::laplacian(rho*invA, p) - fvm::laplacian(rho*invA, p)
==
betav*fvOptions(psi, p, rho.name())
); );
pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter()))); pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter())));
...@@ -69,6 +73,7 @@ else ...@@ -69,6 +73,7 @@ else
U = HbyA - (invA & (betav*fvc::grad(p))); U = HbyA - (invA & (betav*fvc::grad(p)));
U.correctBoundaryConditions(); U.correctBoundaryConditions();
fvOptions.correct(U);
K = 0.5*magSqr(U); K = 0.5*magSqr(U);
if (thermo.dpdt()) if (thermo.dpdt())
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment