pimpleDyMFoam.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  pimpleDyMFoam.C
26 
27 Description
28  Transient solver for incompressible, turbulent flow of Newtonian fluids
29  on a moving mesh.
30 
31  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "dynamicFvMesh.H"
39 #include "pimpleControl.H"
40 #include "CorrectPhi.H"
41 #include "fvOptions.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCase.H"
50  #include "createTime.H"
51  #include "createDynamicFvMesh.H"
52  #include "initContinuityErrs.H"
53  #include "createControls.H"
54  #include "createFields.H"
55  #include "createUf.H"
56  #include "createFvOptions.H"
57  #include "CourantNo.H"
58  #include "setInitialDeltaT.H"
59 
60  turbulence->validate();
61 
62  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64  Info<< "\nStarting time loop\n" << endl;
65 
66  while (runTime.run())
67  {
68  #include "readControls.H"
69  #include "CourantNo.H"
70 
71  #include "setDeltaT.H"
72 
73  runTime++;
74 
75  Info<< "Time = " << runTime.timeName() << nl << endl;
76 
77  mesh.update();
78 
79  // Calculate absolute flux from the mapped surface velocity
80  phi = mesh.Sf() & Uf;
81 
82  if (mesh.changing() && correctPhi)
83  {
84  #include "correctPhi.H"
85  }
86 
87  // Make the flux relative to the mesh motion
89 
90  if (mesh.changing() && checkMeshCourantNo)
91  {
92  #include "meshCourantNo.H"
93  }
94 
95  // --- Pressure-velocity PIMPLE corrector loop
96  while (pimple.loop())
97  {
98  #include "UEqn.H"
99 
100  // --- Pressure corrector loop
101  while (pimple.correct())
102  {
103  #include "pEqn.H"
104  }
105 
106  if (pimple.turbCorr())
107  {
108  laminarTransport.correct();
109  turbulence->correct();
110  }
111  }
112 
113  runTime.write();
114 
115  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
116  << " ClockTime = " << runTime.elapsedClockTime() << " s"
117  << nl << endl;
118  }
119 
120  Info<< "End\n" << endl;
121 
122  return 0;
123 }
124 
125 
126 // ************************************************************************* //
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
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
Creates and initialises the velocity velocity field Uf.
static const char nl
Definition: Ostream.H:262
Calculates and outputs the mean and maximum Courant Numbers.
messageStream Info
checkMeshCourantNo
Definition: readControls.H:5
correctPhi
Definition: readControls.H:3
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75