/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | www.openfoam.com \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2021 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 . Application overCompressibleInterDyMFoam Group grpMultiphaseSolvers Description Solver for two compressible, non-isothermal, immiscible fluids using VOF (i.e. volume of fluid) phase-fraction based interface capturing approach. This solver supports dynamic mesh motions including overset cases. The momentum and other fluid properties are of the "mixture" and a single momentum equation is solved. Either mixture or two-phase transport modelling may be selected. In the mixture approach, a single laminar, RAS or LES model is selected to model the momentum stress. In the Euler-Euler two-phase approach separate laminar, RAS or LES selected models are selected for each of the phases. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "dynamicFvMesh.H" #include "CMULES.H" #include "EulerDdtScheme.H" #include "CrankNicolsonDdtScheme.H" #include "subCycle.H" #include "compressibleInterPhaseTransportModel.H" #include "pimpleControl.H" #include "fvOptions.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[]) { argList::addNote ( "Solver for two compressible, non-isothermal, immiscible fluids" " using VOF phase-fraction based interface capturing approach.\n" "Supports dynamic mesh motions including overset cases." ); #include "postProcess.H" #include "addCheckCaseOptions.H" #include "setRootCaseLists.H" #include "createTime.H" #include "createDynamicFvMesh.H" pimpleControl pimple(mesh); #include "createTimeControls.H" #include "createDyMControls.H" #include "createFields.H" volScalarField& p = mixture.p(); volScalarField& T = mixture.T(); const volScalarField& psi1 = mixture.thermo1().psi(); const volScalarField& psi2 = mixture.thermo2().psi(); #include "correctPhi.H" #include "createUf.H" 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" faceMask = localMin(mesh).interpolate(cellMask.oldTime()); // Zero Uf on old faceMask (H-I) Uf *= faceMask; const surfaceVectorField Uint(fvc::interpolate(U)); // Update Uf and phi on new C-I faces Uf += (1-faceMask)*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 faceMask = localMin(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 "compressibleAlphaEqnSubCycle.H" const surfaceScalarField faceMask ( localMin(mesh).interpolate(cellMask) ); rhoPhi *= faceMask; turbulence.correctPhasePhi(); #include "UEqn.H" volScalarField divUp("divUp", fvc::div(fvc::absolute(phi, U), p)); #include "TEqn.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; } // ************************************************************************* //