Loading...
Searching...
No Matches
liquidFilmFoam.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2017 Wikki Ltd
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 liquidFilmFoam
28
29Group
30 grpFiniteAreaSolvers
31
32Description
33 Transient solver for incompressible, laminar flow of Newtonian fluids in
34 liquid film formulation.
35
36Author
37 Zeljko Tukovic, FMENA
38 Hrvoje Jasak, Wikki Ltd.
39
40\*---------------------------------------------------------------------------*/
41
42#include "fvCFD.H"
43#include "faCFD.H"
44#include "loopControl.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48int main(int argc, char *argv[])
49{
50 argList::addNote
51 (
52 "Transient solver for incompressible, laminar flow"
53 " of Newtonian fluids in liquid film formulation."
54 );
55
56 #include "setRootCaseLists.H"
57 #include "createTime.H"
58 #include "createMesh.H"
59 #include "createFaMesh.H"
60 #include "readGravitationalAcceleration.H"
62 #include "createFaFields.H"
63 #include "createFvFields.H"
64 #include "createTimeControls.H"
65
66 Info<< "\nStarting time loop\n" << endl;
67
68 while (runTime.run())
69 {
70 #include "readSolutionControls.H"
71 #include "readTimeControls.H"
72 #include "surfaceCourantNo.H"
73 #include "capillaryCourantNo.H"
74 #include "setDeltaT.H"
75
76 ++runTime;
77
78 Info<< "Time = " << runTime.timeName() << nl << endl;
79
80 while (iters.loop())
81 {
82 phi2s = fac::interpolate(h)*phis;
83
84 #include "calcFrictionFactor.H"
85
86 faVectorMatrix UsEqn
87 (
88 fam::ddt(h, Us)
89 + fam::div(phi2s, Us)
90 + fam::Sp
91 (
92 0.0125
93 *frictionFactor.internalField()
94 *mag(Us.internalField()),
95 Us
96 )
97 ==
98 Gs*h
99 - fam::Sp(Sd, Us)
100 );
101
102 UsEqn.relax();
103 solve(UsEqn == - fac::grad(ps*h)/rhol + ps*fac::grad(h)/rhol);
104
105 areaScalarField UsA(UsEqn.A());
106
107 Us = UsEqn.H()/UsA;
108 Us.correctBoundaryConditions();
109
110 phis =
111 (fac::interpolate(Us) & aMesh.Le())
112 - fac::interpolate(1.0/(rhol*UsA))*fac::lnGrad(ps*h)*aMesh.magLe()
113 + fac::interpolate(ps/(rhol*UsA))*fac::lnGrad(h)*aMesh.magLe();
114
115 faScalarMatrix hEqn
116 (
117 fam::ddt(h)
118 + fam::div(phis, h)
119 ==
120 Sm
121 - fam::Sp
122 (
123 Sd/(h + dimensionedScalar("small", dimLength, SMALL)),
124 h
125 )
126 );
127
128 hEqn.relax();
129 hEqn.solve();
130
131 phi2s = hEqn.flux();
132
133 // Bound h
134 h.primitiveFieldRef() = max
135 (
136 max
137 (
138 h.primitiveField(),
139 fac::average(max(h, h0))().primitiveField()
140 *pos(h0.value() - h.primitiveField())
141 ),
142 h0.value()
143 );
144
145 ps = rhol*Gn*h - sigma*fac::laplacian(h);
146 ps.correctBoundaryConditions();
147
148 Us -= (1.0/(rhol*UsA))*fac::grad(ps*h)
149 - (ps/(rhol*UsA))*fac::grad(h);
150 Us.correctBoundaryConditions();
151 }
152
153 if (runTime.writeTime())
154 {
155 vsm.mapToVolume(h, H.boundaryFieldRef());
156 vsm.mapToVolume(Us, U.boundaryFieldRef());
157
158 runTime.write();
159 }
160
161 runTime.printExecutionTime(Info);
162 }
163
164 Info<< "End\n" << endl;
165
166 return 0;
167}
168
169
170// ************************************************************************* //
Template specialisation for scalar faMatrix.
U
Definition pEqn.H:72
engineTime & runTime
Required Variables.
const volSurfaceMapping vsm(aMesh)
volScalarField H(IOobject("H", runTime.timeName(), mesh.thisDb(), IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimLength, Zero))
Required Classes.
Read the control parameters used by setDeltaT.
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
faMatrix< vector > faVectorMatrix
Definition faMatrices.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
scalar h0
loopControl iters(runTime, aMesh.solutionDict(), "solution")
Read the control parameters used by setDeltaT.
volScalarField & h
CEqn solve()
dimensionedScalar rhol("rhol", dimDensity, transportProperties)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
edgeScalarField phis(IOobject("phis", runTime.timeName(), aMesh.thisDb(), IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
Calculates and outputs the mean and maximum Courant Numbers for the Finite Area method.