138 lines
3.8 KiB
C
138 lines
3.8 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
|
|
\\/ 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 3 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, see <http://www.gnu.org/licenses/>.
|
|
|
|
Application
|
|
dnsFoam
|
|
|
|
Description
|
|
Direct numerical simulation solver for boxes of isotropic turbulence.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
#include "Kmesh.H"
|
|
#include "UOprocess.H"
|
|
#include "fft.H"
|
|
#include "calcEk.H"
|
|
#include "graph.H"
|
|
#include "pisoControl.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
#include "postProcess.H"
|
|
|
|
#include "setRootCase.H"
|
|
#include "createTime.H"
|
|
#include "createMeshNoClear.H"
|
|
#include "createControl.H"
|
|
#include "createFields.H"
|
|
#include "initContinuityErrs.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Info<< nl << "Starting time loop" << endl;
|
|
|
|
while (runTime.loop())
|
|
{
|
|
Info<< "Time = " << runTime.timeName() << nl << endl;
|
|
|
|
force.primitiveFieldRef() = ReImSum
|
|
(
|
|
fft::reverseTransform
|
|
(
|
|
K/(mag(K) + 1.0e-6) ^ forceGen.newField(), K.nn()
|
|
)
|
|
);
|
|
|
|
#include "globalProperties.H"
|
|
|
|
fvVectorMatrix UEqn
|
|
(
|
|
fvm::ddt(U)
|
|
+ fvm::div(phi, U)
|
|
- fvm::laplacian(nu, U)
|
|
==
|
|
force
|
|
);
|
|
|
|
solve(UEqn == -fvc::grad(p));
|
|
|
|
|
|
// --- PISO loop
|
|
while (piso.correct())
|
|
{
|
|
volScalarField rAU(1.0/UEqn.A());
|
|
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
|
|
volVectorField HbyA(constrainHbyA(rAU*UEqn.H(), U, p));
|
|
surfaceScalarField phiHbyA
|
|
(
|
|
"phiHbyA",
|
|
fvc::flux(HbyA)
|
|
+ rAUf*fvc::ddtCorr(U, phi)
|
|
);
|
|
|
|
// Update the pressure BCs to ensure flux consistency
|
|
constrainPressure(p, U, phiHbyA, rAUf);
|
|
|
|
fvScalarMatrix pEqn
|
|
(
|
|
fvm::laplacian(rAUf, p) == fvc::div(phiHbyA)
|
|
);
|
|
|
|
pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
|
|
|
|
phi = phiHbyA - pEqn.flux();
|
|
|
|
#include "continuityErrs.H"
|
|
|
|
U = HbyA - rAU*fvc::grad(p);
|
|
U.correctBoundaryConditions();
|
|
}
|
|
|
|
runTime.write();
|
|
|
|
if (runTime.writeTime())
|
|
{
|
|
calcEk(U, K).write
|
|
(
|
|
runTime.path()/"graphs"/runTime.timeName(),
|
|
"Ek",
|
|
runTime.graphFormat()
|
|
);
|
|
}
|
|
|
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
|
<< nl << endl;
|
|
}
|
|
|
|
Info<< "End\n" << endl;
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|