/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 4.1 \\ / 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 . Application icoDyMFoamEngine Description Transient solver for incompressible, laminar flow of Newtonian fluids with dynamic mesh for in-cylinder internal combustion engine simulations Author Hrvoje Jasak, Wikki Ltd. All rights reserved \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "dynamicFvMesh.H" #include "engineTime.H" #include "pisoControl.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { # include "setRootCase.H" # include "createEngineTime.H" # include "createDynamicFvMesh.H" pisoControl piso(mesh); # include "initContinuityErrs.H" # include "createFields.H" # include "createControls.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (runTime.run()) { # include "readControls.H" # include "checkTotalVolume.H" # include "CourantNo.H" // Make the fluxes absolute fvc::makeAbsolute(phi, U); # include "setDeltaT.H" runTime++; Info<< "Time = " << runTime.timeName() << nl << endl; bool meshChanged = mesh.update(); reduce(meshChanged, orOp()); if (meshChanged) { # include "checkTotalVolume.H" # include "correctPhi.H" # include "CourantNo.H" } // Make the fluxes relative to the mesh motion fvc::makeRelative(phi, U); if (checkMeshCourantNo) { # include "meshCourantNo.H" } # include "UEqn.H" // --- PISO loop while (piso.correct()) { rUA = 1.0/UEqn.A(); U = rUA*UEqn.H(); phi = (fvc::interpolate(U) & mesh.Sf()); //+ fvc::ddtPhiCorr(rUA, U, phi); adjustPhi(phi, U, p); // Non-orthogonal pressure corrector loop while (piso.correctNonOrthogonal()) { fvScalarMatrix pEqn ( fvm::laplacian(rUA, p) == fvc::div(phi) ); pEqn.setReference(pRefCell, pRefValue); pEqn.solve ( mesh.solutionDict().solver(p.select(piso.finalInnerIter())) ); if (piso.finalNonOrthogonalIter()) { phi -= pEqn.flux(); } } # include "continuityErrs.H" // Make the fluxes relative to the mesh motion fvc::makeRelative(phi, U); U -= rUA*fvc::grad(p); U.correctBoundaryConditions(); } runTime.write(); Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" << " ClockTime = " << runTime.elapsedClockTime() << " s" << nl << endl; } Info<< "End\n" << endl; return(0); } // ************************************************************************* //