/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011 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 . Application interFoam Description Solver for 2 incompressible, isothermal immiscible fluids using a VOF (volume of fluid) phase-fraction based interface capturing approach. The momentum and other fluid properties are of the "mixture" and a single momentum equation is solved. Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected. For a two-fluid approach see twoPhaseEulerFoam. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "MULES.H" #include "subCycle.H" #include "interfaceProperties.H" #include "twoPhaseMixture.H" #include "turbulenceModel.H" #include "interpolationTable.H" #include "pimpleControl.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { // VVPF b: const dimensionedScalar gunits("gunits", dimensionSet(0,1,-2,0,0,0,0), 9.81); const dimensionedScalar gunits2("gunits", dimensionSet(0,1,-2,0,0,0,0), 1.5); #include "setRootCase.H" #include "createTime.H" #include "createMesh.H" pimpleControl pimple(mesh); #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 "readTimeControls.H" #include "CourantNo.H" #include "alphaCourantNo.H" #include "setDeltaT.H" runTime++; Info<< "Time = " << runTime.timeName() << nl << endl; //--------------------------- //simplest possible solution //just check whether the center of gravity of the droplet is in the left or rigth half of the channel, //assuming the mesh is set up so that the center of origin of the x-axis is in the middle of the channel double alphaSum = 0; forAll(mesh.cells(),i) //forAll(mesh.C(),i) works as well { double currentPosition = mesh.C()[i][0]; if (currentPosition < 0.) { alphaSum -= alpha1[i]; //*mesh.V()[i]; multiplication with cell voulme only necessary for non-uniform mesh } else if (currentPosition > 0.) { alphaSum += alpha1[i]; } } Info << "\nalphaSum: " << alphaSum << "\n"; if (alphaSum < 0.) { g=gunits*vector(0,-1,0) + gunits2*vector(1,0,0); } else if (alphaSum > 0.) { g=gunits*vector(0,-1,0) + gunits2*vector(-1,0,0); } else { g=gunits*vector(0,-1,0); } // VVPF c: Info << " AnalyzeGravity " << " " << g.component(vector::X).value() << " " << g.component(vector::Y).value() << " " << g.component(vector::Z).value() << endl; // VVPF d: // Comes from the file "createFields.H". You dont have to delete the lines // in "createFields.H". The field gh is just initilaized there. Info<< "Calculating field g.h\n" << endl; volScalarField gh("gh", g & mesh.C()); surfaceScalarField ghf("ghf", g & mesh.Cf()); //--------------------------- twoPhaseProperties.correct(); #include "alphaEqnSubCycle.H" // --- Pressure-velocity PIMPLE corrector loop while (pimple.loop()) { #include "UEqn.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; } // ************************************************************************* //