sonicDyMFoam.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  sonicDyMFoam
28 
29 Group
30  grpCompressibleSolvers grpMovingMeshSolvers
31 
32 Description
33  Transient solver for trans-sonic/supersonic, turbulent flow of a
34  compressible gas, with optional mesh motion and mesh topology changes.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
39 #include "dynamicFvMesh.H"
40 #include "psiThermo.H"
42 #include "pimpleControl.H"
43 #include "CorrectPhi.H"
44 #include "fvOptions.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
50  argList::addNote
51  (
52  "Transient solver for trans-sonic/supersonic, turbulent flow"
53  " of a compressible gas.\n"
54  "With optional mesh motion and mesh topology changes."
55  );
56 
57  #include "postProcess.H"
58 
59  #include "setRootCaseLists.H"
60  #include "createTime.H"
61  #include "createDynamicFvMesh.H"
62  #include "createDyMControls.H"
63  #include "createFields.H"
64  #include "createFieldRefs.H"
65  #include "createRhoUf.H"
66  #include "compressibleCourantNo.H"
67  #include "setInitialDeltaT.H"
68  #include "initContinuityErrs.H"
69 
70  turbulence->validate();
71 
72  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74  Info<< "\nStarting time loop\n" << endl;
75 
76  while (runTime.run())
77  {
78  #include "readDyMControls.H"
79 
80  {
81  // Store divrhoU from the previous mesh so that it can be mapped
82  // and used in correctPhi to ensure the corrected phi has the
83  // same divergence
84  volScalarField divrhoU
85  (
86  "divrhoU",
88  );
89 
90  #include "compressibleCourantNo.H"
91  #include "setDeltaT.H"
92 
93  ++runTime;
94 
95  Info<< "Time = " << runTime.timeName() << nl << endl;
96 
97  // Store momentum to set rhoUf for introduced faces.
98  volVectorField rhoU("rhoU", rho*U);
99 
100  // Do any mesh changes
101  mesh.update();
102 
103  if (mesh.changing())
104  {
105  MRF.update();
106 
107  if (correctPhi)
108  {
109  // Calculate absolute flux from the mapped surface velocity
110  phi = mesh.Sf() & rhoUf;
111 
112  #include "correctPhi.H"
113 
114  // Make the fluxes relative to the mesh-motion
116  }
117  }
118 
119  if (checkMeshCourantNo)
120  {
121  #include "meshCourantNo.H"
122  }
123  }
124 
125  #include "rhoEqn.H"
126  Info<< "rho min/max : " << min(rho).value() << " " << max(rho).value()
127  << endl;
128 
129  // --- Pressure-velocity PIMPLE corrector loop
130  while (pimple.loop())
131  {
132  #include "UEqn.H"
133  #include "EEqn.H"
134 
135  // --- Pressure corrector loop
136  while (pimple.correct())
137  {
138  #include "pEqn.H"
139  }
140 
141  if (pimple.turbCorr())
142  {
143  turbulence->correct();
144  }
145  }
146 
147  runTime.write();
148 
149  runTime.printExecutionTime(Info);
150  }
151 
152  Info<< "End\n" << endl;
153 
154  return 0;
155 }
156 
157 
158 // ************************************************************************* //
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
dynamicFvMesh & mesh
IOMRFZoneList & MRF
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
pimpleControl & pimple
Info<< "Reading field U\"<< 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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
messageStream Info
Information stream (stdout output on master, null elsewhere)
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
autoPtr< surfaceVectorField > rhoUf
Creates and initialises the velocity velocity field rhoUf.