2010-05-12 13:27:55 +00:00
|
|
|
/*---------------------------------------------------------------------------*\
|
|
|
|
========= |
|
|
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
|
|
\\ / O peration |
|
|
|
|
\\ / A nd | Copyright held by original author
|
|
|
|
\\/ 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 2 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, write to the Free Software Foundation,
|
|
|
|
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
|
|
|
|
|
|
|
|
Application
|
|
|
|
interDyMFoam
|
|
|
|
|
|
|
|
Description
|
|
|
|
Solver for 2 incompressible, isothermal immiscible fluids using a VOF
|
|
|
|
(volume of fluid) phase-fraction based interface capturing approach,
|
|
|
|
with optional mesh motion and mesh topology changes including adaptive
|
|
|
|
re-meshing.
|
|
|
|
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
|
|
|
|
#include "fvCFD.H"
|
|
|
|
#include "dynamicFvMesh.H"
|
|
|
|
#include "MULES.H"
|
|
|
|
#include "subCycle.H"
|
|
|
|
#include "interfaceProperties.H"
|
|
|
|
#include "twoPhaseMixture.H"
|
2010-09-22 18:13:13 +00:00
|
|
|
#include "turbulenceModel.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
int main(int argc, char *argv[])
|
|
|
|
{
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "setRootCase.H"
|
|
|
|
# include "createTime.H"
|
|
|
|
# include "createDynamicFvMesh.H"
|
|
|
|
# include "readGravitationalAcceleration.H"
|
|
|
|
# include "readPISOControls.H"
|
|
|
|
# include "initContinuityErrs.H"
|
|
|
|
# include "createFields.H"
|
|
|
|
# include "readTimeControls.H"
|
|
|
|
# include "correctPhi.H"
|
|
|
|
# include "CourantNo.H"
|
|
|
|
# include "setInitialDeltaT.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
|
|
|
|
while (runTime.run())
|
|
|
|
{
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "readControls.H"
|
|
|
|
# include "CourantNo.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// Make the fluxes absolute
|
|
|
|
fvc::makeAbsolute(phi, U);
|
|
|
|
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "setDeltaT.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
runTime++;
|
|
|
|
|
|
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
|
2011-02-27 02:17:43 +00:00
|
|
|
bool meshChanged = mesh.update();
|
|
|
|
reduce(meshChanged, orOp<bool>());
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2011-02-27 02:17:43 +00:00
|
|
|
# include "volContinuity.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
volScalarField gh("gh", g & mesh.C());
|
|
|
|
surfaceScalarField ghf("ghf", g & mesh.Cf());
|
|
|
|
|
2011-02-27 02:17:43 +00:00
|
|
|
if (correctPhi && meshChanged)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "correctPhi.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
// Make the fluxes relative to the mesh motion
|
|
|
|
fvc::makeRelative(phi, U);
|
|
|
|
|
2011-02-27 02:17:43 +00:00
|
|
|
if (checkMeshCourantNo)
|
2010-05-12 13:27:55 +00:00
|
|
|
{
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "meshCourantNo.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
twoPhaseProperties.correct();
|
|
|
|
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "alphaEqnSubCycle.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "UEqn.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
|
|
|
|
// --- PISO loop
|
|
|
|
for (int corr=0; corr<nCorr; corr++)
|
|
|
|
{
|
2010-11-17 12:50:34 +00:00
|
|
|
# include "pEqn.H"
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
p = pd + rho*gh;
|
|
|
|
|
|
|
|
if (pd.needReference())
|
|
|
|
{
|
|
|
|
p += dimensionedScalar
|
|
|
|
(
|
|
|
|
"p",
|
|
|
|
p.dimensions(),
|
|
|
|
pRefValue - getRefCellValue(p, pdRefCell)
|
|
|
|
);
|
|
|
|
}
|
|
|
|
|
|
|
|
turbulence->correct();
|
|
|
|
|
|
|
|
runTime.write();
|
|
|
|
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
|
|
<< nl << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
|
2010-09-22 18:13:13 +00:00
|
|
|
return 0;
|
2010-05-12 13:27:55 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
// ************************************************************************* //
|