/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright held by original author \\/ 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 2 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, write to the Free Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA Application interDyMFoam 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 "MULES.H" #include "subCycle.H" #include "interfaceProperties.H" #include "twoPhaseMixture.H" #include "incompressible/RASModel/RASModel.H" #include "EulerDdtScheme.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { #include "setRootCase.H" #include "createTime.H" #include "createDynamicFvMesh.H" #include "readEnvironmentalProperties.H" #include "readPISOControls.H" #include "initContinuityErrs.H" #include "createFields.H" #include "readTimeControls.H" #include "correctPhi.H" #include "CourantNo.H" #include "setInitialDeltaT.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (runTime.run()) { #include "readControls.H" #include "CourantNo.H" // Make the fluxes absolute fvc::makeAbsolute(phi, U); #include "setDeltaT.H" runTime++; Info<< "Time = " << runTime.timeName() << nl << endl; scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime(); // Do any mesh changes mesh.update(); if (mesh.changing()) { Info<< "Execution time for mesh.update() = " << runTime.elapsedCpuTime() - timeBeforeMeshUpdate << " s" << endl; } volScalarField gh("gh", g & mesh.C()); surfaceScalarField ghf("ghf", g & mesh.Cf()); if (mesh.changing() && correctPhi) { #include "correctPhi.H" } // Make the fluxes relative to the mesh motion fvc::makeRelative(phi, U); if (mesh.changing() && checkMeshCourantNo) { #include "meshCourantNo.H" } twoPhaseProperties.correct(); #include "gammaEqnSubCycle.H" #include "UEqn.H" // --- PISO loop for (int corr=0; corrcorrect(); runTime.write(); Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" << " ClockTime = " << runTime.elapsedClockTime() << " s" << nl << endl; } Info<< "End\n" << endl; return(0); } // ************************************************************************* //