/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | foam-extend: Open Source CFD \\ / O peration | Version: 3.2 \\ / 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 viscoelasticFluidFoam Description Transient solver for incompressible, laminar flow of viscoelastic fluids. Author Jovani L. Favero and Hrvoje Jasak. All rights reserved \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "viscoelasticModel.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { # include "setRootCase.H" # include "createTime.H" # include "createMesh.H" # include "createFields.H" # include "initContinuityErrs.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (runTime.run()) { # include "readPISOControls.H" # include "readTimeControls.H" # include "CourantNo.H" # include "setDeltaT.H" runTime++; Info<< "Time = " << runTime.timeName() << nl << endl; // Pressure-velocity SIMPLE corrector loop for (int corr = 0; corr < nCorr; corr++) { // Momentum predictor tmp UEqn ( fvm::ddt(U) + fvm::div(phi, U) - visco.divTau(U) ); UEqn().relax(); solve(UEqn() == -fvc::grad(p)); p.boundaryField().updateCoeffs(); volScalarField rUA = 1.0/UEqn().A(); U = rUA*UEqn().H(); UEqn.clear(); phi = fvc::interpolate(U) & mesh.Sf(); adjustPhi(phi, U, p); // Store pressure for under-relaxation p.storePrevIter(); // Non-orthogonal pressure corrector loop for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++) { fvScalarMatrix pEqn ( fvm::laplacian(rUA, p) == fvc::div(phi) ); pEqn.setReference(pRefCell, pRefValue); pEqn.solve(); if (nonOrth == nNonOrthCorr) { phi -= pEqn.flux(); } } # include "continuityErrs.H" // Explicitly relax pressure for momentum corrector p.relax(); // Momentum corrector U -= rUA*fvc::grad(p); U.correctBoundaryConditions(); visco.correct(); } runTime.write(); Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s\n\n" << endl; } Info<< "End\n" << endl; return(0); } // ************************************************************************* //