Files
CFDEMcoupling-PFM/applications/solvers/cfdemSolverPimple/cfdemSolverPimple.C
2013-05-16 07:06:53 +02:00

119 lines
3.7 KiB
C

/*---------------------------------------------------------------------------*\
CFDEMcoupling - Open Source CFD-DEM coupling
CFDEMcoupling is part of the CFDEMproject
www.cfdem.com
Christoph Goniva, christoph.goniva@cfdem.com
Copyright (C) 2011 OpenFOAM Foundation
Copyright (C) 2009-2012 JKU, Linz
Copyright (C) 2012- DCS Computing GmbH,Linz
-------------------------------------------------------------------------------
License
This file is part of CFDEMcoupling.
CFDEMcoupling 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.
CFDEMcoupling 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 CFDEMcoupling. If not, see <http://www.gnu.org/licenses/>.
Application
cfdemSolverPimple
Description
Large time-step transient solver for incompressible, flow using the PIMPLE
(merged PISO-SIMPLE) algorithm. Turbulence modelling is generic,
i.e. laminar, RAS or LES may be selected.
The code is an evolution of the solver pimpleFoam in OpenFOAM(R) 2.1.x.,
where additional functionality for CFD-DEM coupling is added.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "singlePhaseTransportModel.H"
#include "turbulenceModel.H"
#include "pimpleControl.H"
#include "IObasicSourceList.H"
#include "cfdemCloud.H"
#include "implicitCouple.H"
#include "smoothingModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
#include "createFields.H"
#include "initContinuityErrs.H"
pimpleControl pimple(mesh);
// create cfdemCloud
#include "readGravitationalAcceleration.H"
cfdemCloud particleCloud(mesh);
#include "checkModelType.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
#include "readTimeControls.H"
#include "CourantNo.H"
#include "setDeltaT.H"
runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl;
// do particle stuff
Info << "- evolve()" << endl;
particleCloud.evolve(voidfraction,Us,U);
Info << "update Ksl.internalField()" << endl;
Ksl.oldTime().internalField() = particleCloud.momCoupleM(0).impMomSource();
particleCloud.smoothingM().smoothen(Ksl);
Ksl.correctBoundaryConditions();
// --- Pressure-velocity PIMPLE corrector loop
while (pimple.loop())
{
#include "UEqn.H"
// --- Pressure corrector loop
while (pimple.correct())
{
#include "pEqn.H"
}
if (pimple.turbCorr())
{
turbulence->correct();
}
}
runTime.write();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //