Files
openfoam/applications/solvers/multiphase/interFoam/overInterDyMFoam/overInterDyMFoam.C
sergio 5daa38d5b4 ENH:
Update of overRhoPimpleDyMFoam and overInterDyMFoam solvers.
Adding corresponding tutorials with best possible settings
The main effort was put on reducing pressure spikes as the
stencil change with hole cells on the background mesh.
2018-10-17 10:10:06 -07:00

240 lines
6.8 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2014 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2016-2017 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 <http://www.gnu.org/licenses/>.
Application
overInterDyMFoam
Group
grpMultiphaseSolvers grpMovingMeshSolvers
Description
Solver for 2 incompressible, isothermal immiscible fluids using a VOF
(volume of fluid) phase-fraction based interface capturing approach,
with optional mesh motion and mesh topology changes including adaptive
re-meshing.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "dynamicFvMesh.H"
#include "CMULES.H"
#include "EulerDdtScheme.H"
#include "localEulerDdtScheme.H"
#include "CrankNicolsonDdtScheme.H"
#include "subCycle.H"
#include "immiscibleIncompressibleTwoPhaseMixture.H"
#include "turbulentTransportModel.H"
#include "pimpleControl.H"
#include "fvOptions.H"
#include "CorrectPhi.H"
#include "fvcSmooth.H"
#include "cellCellStencilObject.H"
#include "localMin.H"
#include "interpolationCellPoint.H"
#include "transform.H"
#include "fvMeshSubset.H"
#include "oversetAdjustPhi.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "postProcess.H"
#include "setRootCase.H"
#include "createTime.H"
#include "createDynamicFvMesh.H"
#include "initContinuityErrs.H"
pimpleControl pimple(mesh);
#include "createTimeControls.H"
#include "createDyMControls.H"
#include "createFields.H"
#include "createAlphaFluxes.H"
#include "createFvOptions.H"
volScalarField rAU
(
IOobject
(
"rAU",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
);
#include "correctPhi.H"
#include "createUf.H"
turbulence->validate();
if (!LTS)
{
#include "CourantNo.H"
#include "setInitialDeltaT.H"
}
#include "setCellMask.H"
#include "setInterpolatedCells.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
#include "readControls.H"
if (LTS)
{
#include "setRDeltaT.H"
}
else
{
#include "CourantNo.H"
#include "alphaCourantNo.H"
#include "setDeltaT.H"
}
++runTime;
Info<< "Time = " << runTime.timeName() << nl << endl;
// --- Pressure-velocity PIMPLE corrector loop
while (pimple.loop())
{
if (pimple.firstIter() || moveMeshOuterCorrectors)
{
scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
mesh.update();
if (mesh.changing())
{
Info<< "Execution time for mesh.update() = "
<< runTime.elapsedCpuTime() - timeBeforeMeshUpdate
<< " s" << endl;
// Do not apply previous time-step mesh compression flux
// if the mesh topology changed
if (mesh.topoChanging())
{
talphaPhi1Corr0.clear();
}
gh = (g & mesh.C()) - ghRef;
ghf = (g & mesh.Cf()) - ghRef;
// Update cellMask field for blocking out hole cells
#include "setCellMask.H"
#include "setInterpolatedCells.H"
const surfaceScalarField faceMaskOld
(
localMin<scalar>(mesh).interpolate(cellMask.oldTime())
);
// Zero Uf on old faceMask (H-I)
Uf *= faceMaskOld;
const surfaceVectorField Uint(fvc::interpolate(U));
// Update Uf and phi on new C-I faces
Uf += (1-faceMaskOld)*Uint;
// Update Uf boundary
forAll(Uf.boundaryField(), patchI)
{
Uf.boundaryFieldRef()[patchI] =
Uint.boundaryField()[patchI];
}
phi = mesh.Sf() & Uf;
// Correct phi on individual regions
if (correctPhi)
{
#include "correctPhi.H"
}
mixture.correct();
// Zero phi on current H-I
const surfaceScalarField faceMask
(
localMin<scalar>(mesh).interpolate(cellMask)
);
phi *= faceMask;
U *= cellMask;
// Make the flux relative to the mesh motion
fvc::makeRelative(phi, U);
}
if (mesh.changing() && checkMeshCourantNo)
{
#include "meshCourantNo.H"
}
}
#include "alphaControls.H"
#include "alphaEqnSubCycle.H"
const surfaceScalarField faceMask
(
localMin<scalar>(mesh).interpolate(cellMask)
);
rhoPhi *= faceMask;
mixture.correct();
#include "UEqn.H"
// --- Pressure corrector loop
while (pimple.correct())
{
#include "pEqn.H"
}
if (pimple.turbCorr())
{
turbulence->correct();
}
}
runTime.write();
runTime.printExecutionTime(Info);
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //