/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2020 OpenCFD Ltd. ------------------------------------------------------------------------------- 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 . \*---------------------------------------------------------------------------*/ { volScalarField& rDeltaT = trDeltaT.ref(); const dictionary& pimpleDict = pimple.dict(); // Maximum flow Courant number scalar maxCo(pimpleDict.get("maxCo")); // Maximum time scale scalar maxDeltaT(pimpleDict.getOrDefault("maxDeltaT", GREAT)); // Smoothing parameter (0-1) when smoothing iterations > 0 scalar rDeltaTSmoothingCoeff ( pimpleDict.getOrDefault("rDeltaTSmoothingCoeff", 0.1) ); // Damping coefficient (1-0) scalar rDeltaTDampingCoeff ( pimpleDict.getOrDefault("rDeltaTDampingCoeff", 0.2) ); // Maximum change in cell temperature per iteration // (relative to previous value) scalar alphaTemp(pimpleDict.getOrDefault("alphaTemp", 0.05)); Info<< "Time scales min/max:" << endl; // Cache old reciprocal time scale field volScalarField rDeltaT0("rDeltaT0", rDeltaT); // Flow time scale { rDeltaT.ref() = ( fvc::surfaceSum(mag(phi))()() /((2*maxCo)*mesh.V()*rho()) ); // Limit the largest time scale rDeltaT.max(1/maxDeltaT); Info<< " Flow = " << gMin(1/rDeltaT.primitiveField()) << ", " << gMax(1/rDeltaT.primitiveField()) << endl; } // Reaction source time scale { volScalarField::Internal rDeltaTT ( mag ( parcels.hsTrans()/(mesh.V()*runTime.deltaT()) + Qdot ) /( alphaTemp *rho() *thermo.Cp()()() *T() ) ); Info<< " Temperature = " << gMin(1/(rDeltaTT.field() + VSMALL)) << ", " << gMax(1/(rDeltaTT.field() + VSMALL)) << endl; rDeltaT.ref() = max ( rDeltaT(), rDeltaTT ); } // Update the boundary values of the reciprocal time-step rDeltaT.correctBoundaryConditions(); // Spatially smooth the time scale field if (rDeltaTSmoothingCoeff < 1.0) { fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff); } // Limit rate of change of time scale // - reduce as much as required // - only increase at a fraction of old time scale if ( rDeltaTDampingCoeff < 1.0 && runTime.timeIndex() > runTime.startTimeIndex() + 1 ) { rDeltaT = max ( rDeltaT, (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0 ); } Info<< " Overall = " << gMin(1/rDeltaT.primitiveField()) << ", " << gMax(1/rDeltaT.primitiveField()) << endl; } // ************************************************************************* //