multiphaseInterDyMFoam.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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  multiphaseInterFoam
26 
27 Description
28  Solver for n incompressible fluids which captures the interfaces and
29  includes surface-tension and contact-angle effects for each phase, with
30  optional mesh motion and mesh topology changes.
31 
32  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "dynamicFvMesh.H"
38 #include "multiphaseMixture.H"
40 #include "pimpleControl.H"
41 #include "fvOptions.H"
42 #include "CorrectPhi.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "postProcess.H"
49 
50  #include "setRootCase.H"
51  #include "createTime.H"
52  #include "createDynamicFvMesh.H"
53  #include "initContinuityErrs.H"
54  #include "createControl.H"
55  #include "createTimeControls.H"
56  #include "createDyMControls.H"
57  #include "createFields.H"
58  #include "createFvOptions.H"
59 
61  (
62  IOobject
63  (
64  "rAU",
65  runTime.timeName(),
66  mesh,
67  IOobject::READ_IF_PRESENT,
68  IOobject::AUTO_WRITE
69  ),
70  mesh,
71  dimensionedScalar("rAUf", dimTime/rho.dimensions(), 1.0)
72  );
73 
74  #include "correctPhi.H"
75  #include "createUf.H"
76  #include "CourantNo.H"
77  #include "setInitialDeltaT.H"
78 
79  const surfaceScalarField& rhoPhi(mixture.rhoPhi());
80 
81  turbulence->validate();
82 
83  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85  Info<< "\nStarting time loop\n" << endl;
86 
87  while (runTime.run())
88  {
89  #include "readControls.H"
90  #include "CourantNo.H"
91  #include "alphaCourantNo.H"
92 
93  #include "setDeltaT.H"
94 
95  runTime++;
96 
97  Info<< "Time = " << runTime.timeName() << nl << endl;
98 
99  // --- Pressure-velocity PIMPLE corrector loop
100  while (pimple.loop())
101  {
102  if (pimple.firstIter() || moveMeshOuterCorrectors)
103  {
104  scalar timeBeforeMeshUpdate = runTime.elapsedCpuTime();
105 
106  mesh.update();
107 
108  if (mesh.changing())
109  {
110  Info<< "Execution time for mesh.update() = "
111  << runTime.elapsedCpuTime() - timeBeforeMeshUpdate
112  << " s" << endl;
113 
114  gh = (g & mesh.C()) - ghRef;
115  ghf = (g & mesh.Cf()) - ghRef;
116  }
117 
118  if (mesh.changing() && correctPhi)
119  {
120  // Calculate absolute flux from the mapped surface velocity
121  phi = mesh.Sf() & Uf;
122 
123  #include "correctPhi.H"
124 
125  // Make the flux relative to the mesh motion
127 
128  mixture.correct();
129  }
130 
131  if (mesh.changing() && checkMeshCourantNo)
132  {
133  #include "meshCourantNo.H"
134  }
135  }
136 
137  mixture.solve();
138  rho = mixture.rho();
139 
140  #include "UEqn.H"
141 
142  // --- Pressure corrector loop
143  while (pimple.correct())
144  {
145  #include "pEqn.H"
146  }
147 
148  if (pimple.turbCorr())
149  {
150  turbulence->correct();
151  }
152  }
153 
154  runTime.write();
155 
156  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
157  << " ClockTime = " << runTime.elapsedClockTime() << " s"
158  << nl << endl;
159  }
160 
161  Info<< "End\n" << endl;
162 
163  return 0;
164 }
165 
166 
167 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
surfaceScalarField & phi
U
Definition: pEqn.H:83
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Uf
Definition: pEqn.H:67
const surfaceScalarField & ghf
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
static const char nl
Definition: Ostream.H:262
const dimensionedVector & g
Info<< "Reading field p_rgh\n"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\n"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
moveMeshOuterCorrectors
Definition: readControls.H:8
Calculates and outputs the mean and maximum Courant Numbers.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
rhoPhi
Definition: rhoEqn.H:10
const volScalarField & gh
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
checkMeshCourantNo
Definition: readControls.H:5
correctPhi
Definition: readControls.H:3
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
Read the control parameters used by setDeltaT.
volScalarField rAU(1.0/UEqn.A())