FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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, flow of Newtonian fluids
29  on a moving mesh using the PIMPLE (merged PISO-SIMPLE) algorithm.
30 
31  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
32 
33 Usage
34  - PimpleDyMFoam [OPTION]
35 
36  @param -case <dir> \n
37  Specify the case directory
38 
39  @param -parallel \n
40  Run the case in parallel
41 
42  @param -help \n
43  Display short usage message
44 
45  @param -doc \n
46  Display Doxygen documentation page
47 
48  @param -srcDoc \n
49  Display source code
50 
51 \*---------------------------------------------------------------------------*/
52 
53 #include <finiteVolume/fvCFD.H>
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 int main(int argc, char *argv[])
61 {
62  #include <OpenFOAM/setRootCase.H>
63 
64  #include <OpenFOAM/createTime.H>
68  #include "createFields.H"
70 
71  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73  Info<< "\nStarting time loop\n" << endl;
74 
75  while (runTime.run())
76  {
77  #include "readControls.H"
78  #include <finiteVolume/CourantNo.H>
79 
80  // Make the fluxes absolute
82 
83  #include <finiteVolume/setDeltaT.H>
84 
85  runTime++;
86 
87  Info<< "Time = " << runTime.timeName() << nl << endl;
88 
89  mesh.update();
90 
91  if (mesh.changing() && correctPhi)
92  {
93  #include "correctPhi.H"
94  }
95 
96  // Make the fluxes relative to the mesh motion
98 
100  {
102  }
103 
104  // --- PIMPLE loop
105  for (int ocorr=0; ocorr<nOuterCorr; ocorr++)
106  {
107  if (nOuterCorr != 1)
108  {
109  p.storePrevIter();
110  }
111 
112  #include "UEqn.H"
113 
114  // --- PISO loop
115  for (int corr=0; corr<nCorr; corr++)
116  {
117  rAU = 1.0/UEqn.A();
118 
119  U = rAU*UEqn.H();
120  phi = (fvc::interpolate(U) & mesh.Sf());
121 
122  if (p.needReference())
123  {
125  adjustPhi(phi, U, p);
127  }
128 
129  for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
130  {
131  fvScalarMatrix pEqn
132  (
133  fvm::laplacian(rAU, p) == fvc::div(phi)
134  );
135 
136  pEqn.setReference(pRefCell, pRefValue);
137 
138  if
139  (
140  ocorr == nOuterCorr-1
141  && corr == nCorr-1
142  && nonOrth == nNonOrthCorr)
143  {
144  pEqn.solve(mesh.solver(p.name() + "Final"));
145  }
146  else
147  {
148  pEqn.solve(mesh.solver(p.name()));
149  }
150 
151  if (nonOrth == nNonOrthCorr)
152  {
153  phi -= pEqn.flux();
154  }
155  }
156 
158 
159  // Explicitly relax pressure for momentum corrector
160  if (ocorr != nOuterCorr-1)
161  {
162  p.relax();
163  }
164 
165  // Make the fluxes relative to the mesh motion
167 
168  U -= rAU*fvc::grad(p);
169  U.correctBoundaryConditions();
170  }
171 
172  turbulence->correct();
173  }
174 
175  runTime.write();
176 
177  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
178  << " ClockTime = " << runTime.elapsedClockTime() << " s"
179  << nl << endl;
180  }
181 
182  Info<< "End\n" << endl;
183 
184  return 0;
185 }
186 
187 
188 // ************************ vim: set sw=4 sts=4 et: ************************ //