/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation ------------------------------------------------------------------------------- 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 . Application dnsFoam Group grpDNSSolvers 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[]) { argList::addNote ( "Direct numerical simulation for boxes of isotropic turbulence." ); #include "postProcess.H" #include "setRootCaseLists.H" #include "createTime.H" #include "createMeshNoClear.H" #include "createControl.H" #include "createFields.H" #include "initContinuityErrs.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // label ntot = 1; forAll(K.nn(), idim) { ntot *= K.nn()[idim]; } const scalar recRootN = 1.0/Foam::sqrt(scalar(ntot)); 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() )*recRootN ); #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() ); } runTime.printExecutionTime(Info); } Info<< "End\n" << endl; return 0; } // ************************************************************************* //