1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by the
13 Free Software Foundation, either version 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
28 Incompressible laminar CFD code for simulation of a single bubble rising
29 in a stil liquid. Interface between fluid phases is tracked using moving
32 \*---------------------------------------------------------------------------*/
35 #include "dynamicFvMesh.H"
36 #include "freeSurface.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43 # include "setRootCase.H"
45 # include "createTime.H"
47 # include "createDynamicFvMesh.H"
49 # include "createFields.H"
51 # include "initContinuityErrs.H"
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 Info << "\nStarting time loop\n" << endl;
57 for (runTime++; !runTime.end(); runTime++)
59 Info << "Time = " << runTime.value() << endl << endl;
61 # include "readPISOControls.H"
63 interface.moveMeshPointsForOldFreeSurfDisplacement();
65 interface.updateDisplacementDirections();
67 interface.predictPoints();
69 Info<< "\nMax surface Courant Number = "
70 << interface.maxCourantNumber() << endl << endl;
72 for (int corr=0; corr<nOuterCorr; corr++)
74 // Update interface bc
75 interface.updateBoundaryConditions();
77 // Make the fluxes relative
78 phi -= fvc::meshPhi(rho, U);
80 # include "CourantNo.H"
85 + fvm::div(fvc::interpolate(rho)*phi, U, "div(phi,U)")
86 - fvm::laplacian(mu, U)
89 solve(UEqn == -fvc::grad(p));
92 for (int i=0; i<nCorr; i++)
94 volScalarField AU = UEqn.A();
98 phi = (fvc::interpolate(U) & mesh.Sf());
100 # include "scalePhi.H"
102 // Non-orthogonal pressure corrector loop
103 for (label nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
107 fvm::laplacian(1.0/AU, p)
111 # include "setReference.H"
115 if (nonOrth == nNonOrthCorr)
121 # include "continuityErrs.H"
123 // Momentum corrector
124 U -= fvc::grad(p)/AU;
125 U.correctBoundaryConditions();
128 interface.correctPoints();
130 # include "freeSurfaceContinuityErrs.H"
133 # include "volContinuity.H"
135 Info << "Total surface tension force: "
136 << interface.totalSurfaceTensionForce() << endl;
139 interface.totalViscousForce()
140 + interface.totalPressureForce();
142 Info << "Total force: " << totalForce << endl;
146 Info << "ExecutionTime = "
147 << scalar(runTime.elapsedCpuTime())
148 << " s\n" << endl << endl;
151 Info << "End\n" << endl;
156 // ************************************************************************* //