rhoPimpleAdiabaticFoam.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) 2017 OpenCFD Ltd.
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  rhoPimpleAdiabaticFoam
28 
29 Description
30  Transient solver for laminar or turbulent flow of weakly compressible
31  fluids for low Mach number aeroacoustic applications.
32 
33  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
34  pseudo-transient simulations. The RCM interpolation is used as in
35 
36  \verbatim
37  Knacke, T. (2013).
38  Potential effects of Rhie & Chow type interpolations in airframe
39  noise simulations. In: Schram, C., Dénos, R., Lecomte E. (ed):
40  Accurate and efficient aeroacoustic prediction approaches for
41  airframe noise, VKI LS 2013-03.
42  \endverbatim
43 
44  Contact: info@upstream-cfd.com
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "fvCFD.H"
49 #include "fluidThermo.H"
51 #include "bound.H"
52 #include "pimpleControl.H"
53 #include "fvOptions.H"
54 
55 #include "ddtScheme.H"
56 #include "fvcCorrectAlpha.H"
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 int main(int argc, char *argv[])
61 {
62  argList::addNote
63  (
64  "Transient solver for laminar or turbulent flow"
65  " of weakly compressible fluids for low Mach number"
66  " aeroacoustic applications."
67  );
68 
69  #include "postProcess.H"
70 
71  #include "addCheckCaseOptions.H"
72  #include "setRootCaseLists.H"
73  #include "createTime.H"
74  #include "createMesh.H"
75  #include "createControl.H"
76  #include "createTimeControls.H"
77 
78  #include "createFields.H"
79  #include "createFvOptions.H"
80  #include "initContinuityErrs.H"
81 
82  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84  Info<< "\nStarting time loop\n" << endl;
85 
86  while (runTime.run())
87  {
88  #include "readTimeControls.H"
89  #include "compressibleCourantNo.H"
90  #include "setDeltaT.H"
91 
92  ++runTime;
93 
94  Info<< "Time = " << runTime.timeName() << nl << endl;
95 
96  if (pimple.nCorrPIMPLE() <= 1)
97  {
98  #include "rhoEqn.H"
99  }
100 
101  // --- Pressure-velocity PIMPLE corrector loop
102  while (pimple.loop())
103  {
104  U.storePrevIter();
105  rho.storePrevIter();
106  phi.storePrevIter();
107  phiByRho.storePrevIter();
108 
109  #include "UEqn.H"
110 
111  // --- Pressure corrector loop
112  while (pimple.correct())
113  {
114  #include "pEqn.H"
115  }
116 
117  #include "EEqn.H"
118 
119  if (pimple.turbCorr())
120  {
121  turbulence->correct();
122  }
123  }
124 
125  runTime.write();
126 
127  runTime.printExecutionTime(Info);
128  }
129 
130  Info<< "End\n" << endl;
131 
132  return 0;
133 }
134 
135 
136 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
fvOptions.H
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
fluidThermo.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rho
rho
Definition: readInitialConditions.H:88
phiByRho
phiByRho
Definition: pEqn.H:76
pimpleControl.H
createFvOptions.H
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
bound.H
Bound the given scalar field if it has gone unbounded.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
postProcess.H
Execute application functionObjects to post-process existing results.
ddtScheme.H
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
Required Variables.
createTime.H
fvCFD.H
fvcCorrectAlpha.H
Correct flux-U difference in the internal loop using relaxation factor.
turbulentFluidThermoModel.H