2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
2013-12-11 16:09:41 +00:00
|
|
|
\\ / F ield | foam-extend: Open Source CFD
|
2015-05-17 13:32:07 +00:00
|
|
|
\\ / O peration | Version: 3.2
|
|
|
|
\\ / A nd | Web: http://www.foam-extend.org
|
|
|
|
\\/ M anipulation | For copyright notice see file Copyright
|
2010-05-12 13:27:55 +00:00
|
|
|
-------------------------------------------------------------------------------
|
|
|
|
License
|
2013-12-11 16:09:41 +00:00
|
|
|
This file is part of foam-extend.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend is free software: you can redistribute it and/or modify it
|
2010-05-12 13:27:55 +00:00
|
|
|
under the terms of the GNU General Public License as published by the
|
2013-12-11 16:09:41 +00:00
|
|
|
Free Software Foundation, either version 3 of the License, or (at your
|
2010-05-12 13:27:55 +00:00
|
|
|
option) any later version.
|
|
|
|
|
2013-12-11 16:09:41 +00:00
|
|
|
foam-extend 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.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
2013-12-11 16:09:41 +00:00
|
|
|
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
Application
|
|
|
|
boussinesqBuoyantFoam
|
|
|
|
|
|
|
|
Description
|
|
|
|
Transient solver for buoyancy-driven laminar flow of incompressible
|
|
|
|
Newtonian fluids using the Boussinesq Model
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
2016-05-27 15:46:05 +00:00
|
|
|
#include "pisoControl.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
|
|
|
|
|
|
|
# include "setRootCase.H"
|
|
|
|
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createMesh.H"
|
2016-05-27 15:46:05 +00:00
|
|
|
|
|
|
|
pisoControl piso(mesh);
|
|
|
|
|
2010-05-12 13:27:55 +00:00
|
|
|
# include "readTransportProperties.H"
|
2010-09-22 18:13:13 +00:00
|
|
|
# include "readGravitationalAcceleration.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
# include "createFields.H"
|
|
|
|
# include "initContinuityErrs.H"
|
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
|
2016-05-27 15:46:05 +00:00
|
|
|
while (runTime.loop())
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
|
|
|
|
# include "CourantNo.H"
|
|
|
|
|
|
|
|
fvVectorMatrix UEqn
|
|
|
|
(
|
|
|
|
fvm::ddt(U)
|
|
|
|
+ fvm::div(phi, U)
|
|
|
|
- fvm::laplacian(nu, U)
|
|
|
|
==
|
|
|
|
-beta*(T - T0)*g
|
|
|
|
);
|
|
|
|
|
|
|
|
solve(UEqn == -fvc::grad(p));
|
|
|
|
|
|
|
|
// --- PISO loop
|
|
|
|
|
2016-05-27 15:46:05 +00:00
|
|
|
while (piso.correct())
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
volScalarField rUA = 1.0/UEqn.A();
|
|
|
|
|
|
|
|
U = rUA*UEqn.H();
|
2010-08-26 14:22:03 +00:00
|
|
|
phi = (fvc::interpolate(U) & mesh.Sf())
|
2010-05-12 13:27:55 +00:00
|
|
|
+ fvc::ddtPhiCorr(rUA, U, phi);
|
|
|
|
|
|
|
|
adjustPhi(phi, U, p);
|
|
|
|
|
2016-05-27 15:46:05 +00:00
|
|
|
while (piso.correctNonOrthogonal())
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
fvScalarMatrix pEqn
|
|
|
|
(
|
|
|
|
fvm::laplacian(rUA, p) == fvc::div(phi)
|
|
|
|
);
|
|
|
|
|
|
|
|
pEqn.setReference(pRefCell, pRefValue);
|
|
|
|
pEqn.solve();
|
|
|
|
|
2016-05-27 15:46:05 +00:00
|
|
|
if (piso.finalNonOrthogonalIter())
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
|
|
|
phi -= pEqn.flux();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
# include "continuityErrs.H"
|
|
|
|
|
|
|
|
U -= rUA*fvc::grad(p);
|
|
|
|
U.correctBoundaryConditions();
|
|
|
|
}
|
|
|
|
|
|
|
|
// Solve energy equation
|
|
|
|
solve
|
|
|
|
(
|
|
|
|
fvm::ddt(T)
|
|
|
|
+ fvm::div(phi, T)
|
|
|
|
- fvm::laplacian(DT, T)
|
|
|
|
);
|
|
|
|
|
|
|
|
rho = rho0*(scalar(1) - beta*(T - T0));
|
|
|
|
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
|
|
<< nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
|
|
|
return(0);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|