From 042354a88ecbd64b6933c60ad441897e9bf4d84c Mon Sep 17 00:00:00 2001 From: Henry <Henry> Date: Mon, 29 Apr 2013 14:34:15 +0100 Subject: [PATCH] rhoPimpleDyMFoam: New version of rhoPimpleFoam for moving-mesh cases --- .../compressible/rhoPimpleFoam/Allwmake | 1 + .../rhoPimpleFoam/rhoPimpleDyMFoam/Make/files | 3 + .../rhoPimpleDyMFoam/Make/options | 27 ++++ .../rhoPimpleDyMFoam/correctPhi.H | 59 +++++++ .../rhoPimpleDyMFoam/createPcorrTypes.H | 13 ++ .../rhoPimpleFoam/rhoPimpleDyMFoam/pEqn.H | 114 ++++++++++++++ .../rhoPimpleDyMFoam/readControls.H | 7 + .../rhoPimpleDyMFoam/rhoPimpleDyMFoam.C | 146 ++++++++++++++++++ .../pimpleFoam/pimpleDyMFoam/Make/options | 19 ++- .../pimpleFoam/pimpleDyMFoam/UEqn.H | 23 --- .../pimpleFoam/pimpleDyMFoam/correctPhi.H | 20 +-- .../pimpleFoam/pimpleDyMFoam/createFields.H | 16 -- .../pimpleDyMFoam/createPcorrTypes.H | 13 ++ .../pimpleFoam/pimpleDyMFoam/pEqn.H | 6 +- .../pimpleFoam/pimpleDyMFoam/pimpleDyMFoam.C | 15 +- .../pimpleFoam/pimpleDyMFoam/readControls.H | 9 +- .../pimpleDyMFoam/movingCone/system/fvSchemes | 8 +- 17 files changed, 411 insertions(+), 88 deletions(-) create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/files create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/options create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/correctPhi.H create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/createPcorrTypes.H create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/pEqn.H create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/readControls.H create mode 100644 applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/rhoPimpleDyMFoam.C delete mode 100644 applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/UEqn.H create mode 100644 applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createPcorrTypes.H diff --git a/applications/solvers/compressible/rhoPimpleFoam/Allwmake b/applications/solvers/compressible/rhoPimpleFoam/Allwmake index ac06b7350a1..1b272c17f07 100755 --- a/applications/solvers/compressible/rhoPimpleFoam/Allwmake +++ b/applications/solvers/compressible/rhoPimpleFoam/Allwmake @@ -5,5 +5,6 @@ set -x wmake wmake rhoPimplecFoam wmake rhoLTSPimpleFoam +wmake rhoPimpleDyMFoam # ----------------------------------------------------------------- end-of-file diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/files b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/files new file mode 100644 index 00000000000..034b5c2b1b7 --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/files @@ -0,0 +1,3 @@ +rhoPimpleDyMFoam.C + +EXE = $(FOAM_APPBIN)/rhoPimpleDyMFoam diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/options b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/options new file mode 100644 index 00000000000..2093b620b5e --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/Make/options @@ -0,0 +1,27 @@ +EXE_INC = \ + -I.. \ + -I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \ + -I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel \ + -I$(LIB_SRC)/finiteVolume/cfdTools \ + -I$(LIB_SRC)/finiteVolume/lnInclude \ + -I$(LIB_SRC)/meshTools/lnInclude \ + -I$(LIB_SRC)/sampling/lnInclude \ + -I$(LIB_SRC)/fvOptions/lnInclude \ + -I$(LIB_SRC)/dynamicFvMesh/lnInclude \ + -I$(LIB_SRC)/dynamicMesh/lnInclude \ + -I$(LIB_SRC)/meshTools/lnInclude \ + +EXE_LIBS = \ + -lfluidThermophysicalModels \ + -lspecie \ + -lcompressibleTurbulenceModel \ + -lcompressibleRASModels \ + -lcompressibleLESModels \ + -lfiniteVolume \ + -lmeshTools \ + -lsampling \ + -lfvOptions \ + -ldynamicFvMesh \ + -ltopoChangerFvMesh \ + -ldynamicMesh \ + -lmeshTools diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/correctPhi.H b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/correctPhi.H new file mode 100644 index 00000000000..24316a16854 --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/correctPhi.H @@ -0,0 +1,59 @@ +{ + if (mesh.changing()) + { + forAll(U.boundaryField(), patchi) + { + if (U.boundaryField()[patchi].fixesValue()) + { + U.boundaryField()[patchi].initEvaluate(); + } + } + + forAll(U.boundaryField(), patchi) + { + if (U.boundaryField()[patchi].fixesValue()) + { + U.boundaryField()[patchi].evaluate(); + + phi.boundaryField()[patchi] = + rho.boundaryField()[patchi] + *( + U.boundaryField()[patchi] + & mesh.Sf().boundaryField()[patchi] + ); + } + } + } + + volScalarField pcorr + ( + IOobject + ( + "pcorr", + runTime.timeName(), + mesh, + IOobject::NO_READ, + IOobject::NO_WRITE + ), + mesh, + dimensionedScalar("pcorr", p.dimensions(), 0.0), + pcorrTypes + ); + + dimensionedScalar Dp("Dp", dimTime, 1.0); + + while (pimple.correctNonOrthogonal()) + { + fvScalarMatrix pcorrEqn + ( + fvm::laplacian(Dp, pcorr) == fvc::div(phi) - divrhoU + ); + + pcorrEqn.solve(); + + if (pimple.finalNonOrthogonalIter()) + { + phi -= pcorrEqn.flux(); + } + } +} diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/createPcorrTypes.H b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/createPcorrTypes.H new file mode 100644 index 00000000000..a602fd48430 --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/createPcorrTypes.H @@ -0,0 +1,13 @@ + wordList pcorrTypes + ( + p.boundaryField().size(), + zeroGradientFvPatchScalarField::typeName + ); + + for (label i=0; i<p.boundaryField().size(); i++) + { + if (p.boundaryField()[i].fixesValue()) + { + pcorrTypes[i] = fixedValueFvPatchScalarField::typeName; + } + } diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/pEqn.H b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/pEqn.H new file mode 100644 index 00000000000..7b4118658bd --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/pEqn.H @@ -0,0 +1,114 @@ +rho = thermo.rho(); +rho = max(rho, rhoMin); +rho = min(rho, rhoMax); +rho.relax(); + +volScalarField rAU(1.0/UEqn().A()); +volVectorField HbyA("HbyA", U); +HbyA = rAU*UEqn().H(); + +if (pimple.nCorrPISO() <= 1) +{ + UEqn.clear(); +} + +if (pimple.transonic()) +{ + surfaceScalarField phid + ( + "phid", + fvc::interpolate(psi) + *( + (fvc::interpolate(HbyA) & mesh.Sf()) + + fvc::ddtPhiCorr(rAU, rho, U, phiAbs) + ) + ); + + fvOptions.relativeFlux(fvc::interpolate(psi), phid); + + volScalarField Dp("Dp", rho*rAU); + + while (pimple.correctNonOrthogonal()) + { + fvScalarMatrix pEqn + ( + fvm::ddt(psi, p) + + fvm::div(phid, p) + - fvm::laplacian(Dp, p) + == + fvOptions(psi, p, rho.name()) + ); + + fvOptions.constrain(pEqn); + + pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter()))); + + if (pimple.finalNonOrthogonalIter()) + { + phi == pEqn.flux(); + } + } +} +else +{ + surfaceScalarField phiHbyA + ( + "phiHbyA", + fvc::interpolate(rho) + *( + (fvc::interpolate(HbyA) & mesh.Sf()) + - fvc::meshPhi(rho, U) + + fvc::ddtPhiCorr(rAU, rho, U, phiAbs) + ) + ); + + fvOptions.relativeFlux(fvc::interpolate(rho), phiHbyA); + + volScalarField Dp("Dp", rho*rAU); + + while (pimple.correctNonOrthogonal()) + { + // Pressure corrector + fvScalarMatrix pEqn + ( + fvm::ddt(psi, p) + + fvc::div(phiHbyA) + - fvm::laplacian(Dp, p) + == + fvOptions(psi, p, rho.name()) + ); + + fvOptions.constrain(pEqn); + + pEqn.solve(mesh.solver(p.select(pimple.finalInnerIter()))); + + if (pimple.finalNonOrthogonalIter()) + { + phi = phiHbyA + pEqn.flux(); + } + } +} + +#include "rhoEqn.H" +#include "compressibleContinuityErrs.H" + +// Explicitly relax pressure for momentum corrector +p.relax(); + +// Recalculate density from the relaxed pressure +rho = thermo.rho(); +rho = max(rho, rhoMin); +rho = min(rho, rhoMax); +rho.relax(); +Info<< "rho max/min : " << max(rho).value() + << " " << min(rho).value() << endl; + +U = HbyA - rAU*fvc::grad(p); +U.correctBoundaryConditions(); +fvOptions.correct(U); +K = 0.5*magSqr(U); + +if (thermo.dpdt()) +{ + dpdt = fvc::ddt(p); +} diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/readControls.H b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/readControls.H new file mode 100644 index 00000000000..cf509804604 --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/readControls.H @@ -0,0 +1,7 @@ + #include "readTimeControls.H" + + bool correctPhi = + pimple.dict().lookupOrDefault<Switch>("correctPhi", true); + + bool checkMeshCourantNo = + pimple.dict().lookupOrDefault<Switch>("checkMeshCourantNo", false); diff --git a/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/rhoPimpleDyMFoam.C b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/rhoPimpleDyMFoam.C new file mode 100644 index 00000000000..fb74f2b39a9 --- /dev/null +++ b/applications/solvers/compressible/rhoPimpleFoam/rhoPimpleDyMFoam/rhoPimpleDyMFoam.C @@ -0,0 +1,146 @@ +/*---------------------------------------------------------------------------*\ + ========= | + \\ / F ield | OpenFOAM: The Open Source CFD Toolbox + \\ / O peration | + \\ / A nd | Copyright (C) 2011-2013 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 + rhoPimpleFoam + +Description + Transient solver for laminar or turbulent flow of compressible fluids + for HVAC and similar applications. + + Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and + pseudo-transient simulations. + +\*---------------------------------------------------------------------------*/ + +#include "fvCFD.H" +#include "dynamicFvMesh.H" +#include "psiThermo.H" +#include "turbulenceModel.H" +#include "bound.H" +#include "pimpleControl.H" +#include "fvIOoptionList.H" + +// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + +int main(int argc, char *argv[]) +{ + #include "setRootCase.H" + #include "createTime.H" + #include "createDynamicFvMesh.H" + #include "initContinuityErrs.H" + + pimpleControl pimple(mesh); + + #include "readControls.H" + #include "createFields.H" + #include "createFvOptions.H" + #include "createPcorrTypes.H" + #include "CourantNo.H" + #include "setInitialDeltaT.H" + + // Create old-time absolute flux for ddtPhiCorr + surfaceScalarField phiAbs("phiAbs", phi); + + + // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // + + Info<< "\nStarting time loop\n" << endl; + + while (runTime.run()) + { + #include "readControls.H" + #include "compressibleCourantNo.H" + + // Make the fluxes absolute before mesh-motion + fvc::makeAbsolute(phi, rho, U); + + // Update absolute flux for ddtPhiCorr + phiAbs = phi; + + #include "setDeltaT.H" + + runTime++; + + Info<< "Time = " << runTime.timeName() << nl << endl; + + { + // Store divrhoU from the previous time-step/mesh for the correctPhi + volScalarField divrhoU(fvc::div(phi)); + + // Do any mesh changes + mesh.update(); + + if (mesh.changing() && correctPhi) + { + #include "correctPhi.H" + } + } + + // Make the fluxes relative to the mesh-motion + fvc::makeRelative(phi, rho, U); + + if (mesh.changing() && checkMeshCourantNo) + { + #include "meshCourantNo.H" + } + + if (pimple.nCorrPIMPLE() <= 1) + { + #include "rhoEqn.H" + Info<< "rhoEqn max/min : " << max(rho).value() + << " " << min(rho).value() << endl; + } + + // --- Pressure-velocity PIMPLE corrector loop + while (pimple.loop()) + { + #include "UEqn.H" + #include "EEqn.H" + + // --- Pressure corrector loop + while (pimple.correct()) + { + #include "pEqn.H" + } + + if (pimple.turbCorr()) + { + turbulence->correct(); + } + } + + runTime.write(); + + Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" + << " ClockTime = " << runTime.elapsedClockTime() << " s" + << nl << endl; + } + + Info<< "End\n" << endl; + + return 0; +} + + +// ************************************************************************* // diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/Make/options b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/Make/options index 5cc44fea11a..263fd7375f8 100644 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/Make/options +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/Make/options @@ -1,25 +1,24 @@ EXE_INC = \ -I.. \ - -I$(LIB_SRC)/dynamicFvMesh/lnInclude \ - -I$(LIB_SRC)/dynamicMesh/lnInclude \ - -I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/turbulenceModels/incompressible/turbulenceModel \ -I$(LIB_SRC)/transportModels \ -I$(LIB_SRC)/transportModels/incompressible/singlePhaseTransportModel \ -I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/fvOptions/lnInclude \ - -I$(LIB_SRC)/sampling/lnInclude - + -I$(LIB_SRC)/sampling/lnInclude \ + -I$(LIB_SRC)/dynamicFvMesh/lnInclude \ + -I$(LIB_SRC)/dynamicMesh/lnInclude \ + -I$(LIB_SRC)/meshTools/lnInclude \ EXE_LIBS = \ - -ldynamicFvMesh \ - -ltopoChangerFvMesh \ - -ldynamicMesh \ - -lmeshTools \ -lincompressibleTransportModels \ -lincompressibleTurbulenceModel \ -lincompressibleRASModels \ -lincompressibleLESModels \ -lfiniteVolume \ -lfvOptions \ - -lsampling + -lsampling \ + -ldynamicFvMesh \ + -ltopoChangerFvMesh \ + -ldynamicMesh \ + -lmeshTools diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/UEqn.H b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/UEqn.H deleted file mode 100644 index eaaa3e60c90..00000000000 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/UEqn.H +++ /dev/null @@ -1,23 +0,0 @@ -// Solve the Momentum equation - -tmp<fvVectorMatrix> UEqn -( - fvm::ddt(U) - + fvm::div(phi, U) - + turbulence->divDevReff(U) - == - fvOptions(U) -); - -UEqn().relax(); - -fvOptions.constrain(UEqn()); - -rAU = 1.0/UEqn().A(); - -if (pimple.momentumPredictor()) -{ - solve(UEqn() == -fvc::grad(p)); - - fvOptions.correct(U); -} diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/correctPhi.H b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/correctPhi.H index 9bed803d1e7..e1e897c18b2 100644 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/correctPhi.H +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/correctPhi.H @@ -22,20 +22,6 @@ } } - wordList pcorrTypes - ( - p.boundaryField().size(), - zeroGradientFvPatchScalarField::typeName - ); - - forAll(p.boundaryField(), patchI) - { - if (p.boundaryField()[patchI].fixesValue()) - { - pcorrTypes[patchI] = fixedValueFvPatchScalarField::typeName; - } - } - volScalarField pcorr ( IOobject @@ -51,11 +37,13 @@ pcorrTypes ); + dimensionedScalar Dp("Dp", dimTime, 1.0); + while (pimple.correctNonOrthogonal()) { fvScalarMatrix pcorrEqn ( - fvm::laplacian(rAU, pcorr) == fvc::div(phi) + fvm::laplacian(Dp, pcorr) == fvc::div(phi) ); pcorrEqn.setReference(pRefCell, pRefValue); @@ -68,6 +56,4 @@ } } -phi.oldTime() = phi; - #include "continuityErrs.H" diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createFields.H b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createFields.H index 16b3bd977d0..082ec86ea8e 100644 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createFields.H +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createFields.H @@ -40,19 +40,3 @@ ( incompressible::turbulenceModel::New(U, phi, laminarTransport) ); - - Info<< "Reading field rAU if present\n" << endl; - volScalarField rAU - ( - IOobject - ( - "rAU", - runTime.timeName(), - mesh, - IOobject::READ_IF_PRESENT, - IOobject::AUTO_WRITE - ), - mesh, - runTime.deltaT(), - zeroGradientFvPatchScalarField::typeName - ); diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createPcorrTypes.H b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createPcorrTypes.H new file mode 100644 index 00000000000..a602fd48430 --- /dev/null +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/createPcorrTypes.H @@ -0,0 +1,13 @@ + wordList pcorrTypes + ( + p.boundaryField().size(), + zeroGradientFvPatchScalarField::typeName + ); + + for (label i=0; i<p.boundaryField().size(); i++) + { + if (p.boundaryField()[i].fixesValue()) + { + pcorrTypes[i] = fixedValueFvPatchScalarField::typeName; + } + } diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pEqn.H b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pEqn.H index 4c48bb2a183..c9c0ce3c58e 100644 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pEqn.H +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pEqn.H @@ -10,13 +10,9 @@ surfaceScalarField phiHbyA ( "phiHbyA", (fvc::interpolate(HbyA) & mesh.Sf()) + + fvc::ddtPhiCorr(rAU, U, phiAbs) ); -if (ddtPhiCorr) -{ - phiHbyA += fvc::ddtPhiCorr(rAU, U, phi); -} - if (p.needReference()) { fvc::makeRelative(phiHbyA, U); diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pimpleDyMFoam.C b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pimpleDyMFoam.C index 3f68ebbcd34..883a20f6cff 100644 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pimpleDyMFoam.C +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/pimpleDyMFoam.C @@ -33,9 +33,9 @@ Description \*---------------------------------------------------------------------------*/ #include "fvCFD.H" +#include "dynamicFvMesh.H" #include "singlePhaseTransportModel.H" #include "turbulenceModel.H" -#include "dynamicFvMesh.H" #include "pimpleControl.H" #include "fvIOoptionList.H" @@ -44,15 +44,21 @@ Description int main(int argc, char *argv[]) { #include "setRootCase.H" - #include "createTime.H" #include "createDynamicFvMesh.H" #include "initContinuityErrs.H" + + pimpleControl pimple(mesh); + #include "createFields.H" #include "createFvOptions.H" #include "readTimeControls.H" + #include "createPcorrTypes.H" + #include "CourantNo.H" + #include "setInitialDeltaT.H" - pimpleControl pimple(mesh); + // Create old-time absolute flux for ddtPhiCorr + surfaceScalarField phiAbs("phiAbs", phi); // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // @@ -66,6 +72,9 @@ int main(int argc, char *argv[]) // Make the fluxes absolute fvc::makeAbsolute(phi, U); + // Update absolute flux for ddtPhiCorr + phiAbs = phi; + #include "setDeltaT.H" runTime++; diff --git a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/readControls.H b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/readControls.H index 2ad7f610f0f..1c53b5e147d 100644 --- a/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/readControls.H +++ b/applications/solvers/incompressible/pimpleFoam/pimpleDyMFoam/readControls.H @@ -1,12 +1,7 @@ #include "readTimeControls.H" - const dictionary& pimpleDict = pimple.dict(); - const bool correctPhi = - pimpleDict.lookupOrDefault("correctPhi", false); + pimple.dict().lookupOrDefault("correctPhi", false); const bool checkMeshCourantNo = - pimpleDict.lookupOrDefault("checkMeshCourantNo", false); - - const bool ddtPhiCorr = - pimpleDict.lookupOrDefault("ddtPhiCorr", true); + pimple.dict().lookupOrDefault("checkMeshCourantNo", false); diff --git a/tutorials/incompressible/pimpleDyMFoam/movingCone/system/fvSchemes b/tutorials/incompressible/pimpleDyMFoam/movingCone/system/fvSchemes index fe094084a1a..b7e6607f899 100644 --- a/tutorials/incompressible/pimpleDyMFoam/movingCone/system/fvSchemes +++ b/tutorials/incompressible/pimpleDyMFoam/movingCone/system/fvSchemes @@ -23,7 +23,6 @@ ddtSchemes gradSchemes { default Gauss linear; - grad(p) Gauss linear; } divSchemes @@ -35,18 +34,13 @@ divSchemes laplacianSchemes { - default none; - laplacian(nu,U) Gauss linear corrected; - laplacian(rAU,pcorr) Gauss linear corrected; - laplacian(rAU,p) Gauss linear corrected; + default Gauss linear corrected; laplacian(diffusivity,cellMotionU) Gauss linear uncorrected; - laplacian(nuEff,U) Gauss linear uncorrected; } interpolationSchemes { default linear; - interpolate(HbyA) linear; } snGradSchemes -- GitLab