OpenFOAM: API Guide
v2112
The open source CFD toolbox
multiphaseInterFoam.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) 2011-2017 OpenFOAM Foundation
9
-------------------------------------------------------------------------------
10
License
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
26
Application
27
multiphaseInterFoam
28
29
Group
30
grpMultiphaseSolvers
31
32
Description
33
Solver for N incompressible fluids which captures the interfaces and
34
includes surface-tension and contact-angle effects for each phase, with
35
optional mesh motion and mesh topology changes.
36
37
Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
38
39
\*---------------------------------------------------------------------------*/
40
41
#include "
fvCFD.H
"
42
#include "
dynamicFvMesh.H
"
43
#include "
multiphaseMixture.H
"
44
#include "
turbulentTransportModel.H
"
45
#include "
pimpleControl.H
"
46
#include "
fvOptions.H
"
47
#include "
CorrectPhi.H
"
48
49
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51
int
main(
int
argc,
char
*argv[])
52
{
53
argList::addNote
54
(
55
"Solver for N incompressible fluids which captures the interfaces and"
56
" includes surface-tension and contact-angle effects for each phase.\n"
57
"With optional mesh motion and mesh topology changes."
58
);
59
60
#include "
postProcess.H
"
61
62
#include "
addCheckCaseOptions.H
"
63
#include "
setRootCaseLists.H
"
64
#include "
createTime.H
"
65
#include "
createDynamicFvMesh.H
"
66
#include "initContinuityErrs.H"
67
#include "createDyMControls.H"
68
#include "createFields.H"
69
#include "
initCorrectPhi.H
"
70
#include "
createUfIfPresent.H
"
71
72
turbulence
->validate();
73
74
#include "CourantNo.H"
75
#include "setInitialDeltaT.H"
76
77
const
surfaceScalarField
&
rhoPhi
(
mixture
.rhoPhi());
78
79
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80
81
Info
<<
"\nStarting time loop\n"
<<
endl
;
82
83
while
(
runTime
.run())
84
{
85
#include "
readDyMControls.H
"
86
#include "CourantNo.H"
87
#include "alphaCourantNo.H"
88
#include "setDeltaT.H"
89
90
++
runTime
;
91
92
Info
<<
"Time = "
<<
runTime
.timeName() <<
nl
<<
endl
;
93
94
// --- Pressure-velocity PIMPLE corrector loop
95
while
(
pimple
.loop())
96
{
97
if
(
pimple
.firstIter() ||
moveMeshOuterCorrectors
)
98
{
99
scalar timeBeforeMeshUpdate =
runTime
.elapsedCpuTime();
100
101
mesh
.update();
102
103
if
(
mesh
.changing())
104
{
105
Info
<<
"Execution time for mesh.update() = "
106
<<
runTime
.elapsedCpuTime() - timeBeforeMeshUpdate
107
<<
" s"
<<
endl
;
108
109
gh
= (
g
&
mesh
.C()) - ghRef;
110
ghf
= (
g
&
mesh
.Cf()) - ghRef;
111
112
MRF
.update();
113
114
if
(
correctPhi
)
115
{
116
// Calculate absolute flux
117
// from the mapped surface velocity
118
phi
=
mesh
.Sf() &
Uf
();
119
120
#include "correctPhi.H"
121
122
// Make the flux relative to the mesh motion
123
fvc::makeRelative(
phi
,
U
);
124
125
mixture
.correct();
126
}
127
128
if
(
checkMeshCourantNo
)
129
{
130
#include "
meshCourantNo.H
"
131
}
132
}
133
}
134
135
mixture
.solve();
136
rho
=
mixture
.rho();
137
138
#include "UEqn.H"
139
140
// --- Pressure corrector loop
141
while
(
pimple
.correct())
142
{
143
#include "pEqn.H"
144
}
145
146
if
(
pimple
.turbCorr())
147
{
148
turbulence
->correct();
149
}
150
}
151
152
runTime
.write();
153
154
runTime
.printExecutionTime(Info);
155
}
156
157
Info
<<
"End\n"
<<
endl
;
158
159
return
0;
160
}
161
162
163
// ************************************************************************* //
CorrectPhi.H
addCheckCaseOptions.H
Required Classes.
rhoPhi
rhoPhi
Definition:
rhoEqn.H:10
g
const uniformDimensionedVectorField & g
Definition:
createFluidFields.H:26
phi
surfaceScalarField & phi
Definition:
setRegionFluidFields.H:8
ghf
const surfaceScalarField & ghf
Definition:
setRegionFluidFields.H:18
MRF
IOMRFZoneList & MRF
Definition:
setRegionFluidFields.H:22
gh
const volScalarField & gh
Definition:
setRegionFluidFields.H:17
pimple
pimpleControl & pimple
Definition:
setRegionFluidFields.H:56
U
U
Definition:
pEqn.H:72
createDynamicFvMesh.H
mesh
dynamicFvMesh & mesh
Definition:
createDynamicFvMesh.H:6
runTime
engineTime & runTime
Definition:
createEngineTime.H:13
createTime.H
createUfIfPresent.H
Creates and initialises the velocity field Uf if required.
Uf
autoPtr< surfaceVectorField > Uf
Definition:
createUfIfPresent.H:33
dynamicFvMesh.H
turbulence
compressible::turbulenceModel & turbulence
Definition:
setRegionFluidFields.H:30
fvCFD.H
fvOptions.H
initCorrectPhi.H
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
multiphaseMixture.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition:
surfaceFieldsFwd.H:66
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition:
Ostream.H:372
Foam::nl
constexpr char nl
The newline '\n' character (0x0a)
Definition:
Ostream.H:53
pimpleControl.H
postProcess.H
Execute application functionObjects to post-process existing results.
readDyMControls.H
checkMeshCourantNo
checkMeshCourantNo
Definition:
readDyMControls.H:9
correctPhi
correctPhi
Definition:
readDyMControls.H:3
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition:
readDyMControls.H:15
rho
rho
Definition:
readInitialConditions.H:88
setRootCaseLists.H
mixture
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture
Definition:
createFields.H:39
turbulentTransportModel.H
applications
solvers
multiphase
multiphaseInterFoam
multiphaseInterFoam.C
Generated by
1.9.5
OPENFOAM® is a registered
trademark
of OpenCFD Ltd.