This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/immersedBoundary/pimpleDyMIbFoam/pimpleDyMIbFoam.C

138 lines
3.9 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend 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.
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.
You should have received a copy of the GNU General Public License
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
pimpleDyMFoam.C
Description
Transient solver for incompressible, flow of Newtonian fluids
with dynamic mesh using the PIMPLE (merged PISO-SIMPLE) algorithm.
Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
Consistent formulation without time-step and relaxation dependence by Jasak
Support for immersed boundary
Author
Hrvoje Jasak, Wikki Ltd. All rights reserved
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "singlePhaseTransportModel.H"
#include "turbulenceModel.H"
#include "dynamicFvMesh.H"
#include "pimpleControl.H"
#include "immersedBoundaryPolyPatch.H"
#include "immersedBoundaryFvPatch.H"
#include "emptyFvPatch.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
# include "setRootCase.H"
# include "createTime.H"
# include "createDynamicFvMesh.H"
pimpleControl pimple(mesh);
# include "initContinuityErrs.H"
# include "createIbMasks.H"
# include "createFields.H"
# include "createControls.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
# include "readControls.H"
# include "CourantNo.H"
# include "setDeltaT.H"
runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl;
bool meshChanged = mesh.update();
U.correctBoundaryConditions();
p.correctBoundaryConditions();
aU.correctBoundaryConditions();
phi.correctBoundaryConditions();
turbulence->correct();
# include "updateIbMasks.H"
# include "volContinuity.H"
if (runTime.outputTime())
{
volScalarField divMeshPhi("divMeshPhi", mag(fvc::surfaceIntegrate(mesh.phi())));
divMeshPhi.write();
}
if (checkMeshCourantNo)
{
# include "meshCourantNo.H"
}
// Fluxes will be corrected to absolute velocity
// HJ, 6/Feb/2009
# include "correctPhi.H"
// Make the fluxes relative to the mesh motion
fvc::makeRelative(phi, U);
// --- PIMPLE loop
while (pimple.loop())
{
# include "UEqn.H"
// --- PISO loop
while (pimple.correct())
{
# include "pEqn.H"
}
turbulence->correct();
}
runTime.write();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //