sonicDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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  sonicDyMFoam
26 
27 Group
28  grpCompressibleSolvers grpMovingMeshSolvers
29 
30 Description
31  Transient solver for trans-sonic/supersonic, turbulent flow of a
32  compressible gas, with optional mesh motion and mesh topology changes.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "dynamicFvMesh.H"
38 #include "psiThermo.H"
40 #include "pimpleControl.H"
41 #include "CorrectPhi.H"
42 #include "fvOptions.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #include "postProcess.H"
49 
50  #include "setRootCaseLists.H"
51  #include "createTime.H"
52  #include "createDynamicFvMesh.H"
53  #include "createDyMControls.H"
54  #include "createFields.H"
55  #include "createFieldRefs.H"
56  #include "createRhoUf.H"
57  #include "compressibleCourantNo.H"
58  #include "setInitialDeltaT.H"
59  #include "initContinuityErrs.H"
60 
61  turbulence->validate();
62 
63  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65  Info<< "\nStarting time loop\n" << endl;
66 
67  while (runTime.run())
68  {
69  #include "readDyMControls.H"
70 
71  {
72  // Store divrhoU from the previous mesh so that it can be mapped
73  // and used in correctPhi to ensure the corrected phi has the
74  // same divergence
75  volScalarField divrhoU
76  (
77  "divrhoU",
79  );
80 
81  #include "compressibleCourantNo.H"
82  #include "setDeltaT.H"
83 
84  runTime++;
85 
86  Info<< "Time = " << runTime.timeName() << nl << endl;
87 
88  // Store momentum to set rhoUf for introduced faces.
89  volVectorField rhoU("rhoU", rho*U);
90 
91  // Do any mesh changes
92  mesh.update();
93 
94  if (mesh.changing())
95  {
96  MRF.update();
97 
98  if (correctPhi)
99  {
100  // Calculate absolute flux from the mapped surface velocity
101  phi = mesh.Sf() & rhoUf;
102 
103  #include "correctPhi.H"
104 
105  // Make the fluxes relative to the mesh-motion
107  }
108  }
109 
110  if (checkMeshCourantNo)
111  {
112  #include "meshCourantNo.H"
113  }
114  }
115 
116  #include "rhoEqn.H"
117  Info<< "rhoEqn max/min : " << max(rho).value()
118  << " " << min(rho).value() << endl;
119 
120  // --- Pressure-velocity PIMPLE corrector loop
121  while (pimple.loop())
122  {
123  #include "UEqn.H"
124  #include "EEqn.H"
125 
126  // --- Pressure corrector loop
127  while (pimple.correct())
128  {
129  #include "pEqn.H"
130  }
131 
132  if (pimple.turbCorr())
133  {
134  turbulence->correct();
135  }
136  }
137 
138  runTime.write();
139 
140  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
141  << " ClockTime = " << runTime.elapsedClockTime() << " s"
142  << nl << endl;
143  }
144 
145  Info<< "End\n" << endl;
146 
147  return 0;
148 }
149 
150 
151 // ************************************************************************* //
pimpleNoLoopControl & pimple
surfaceScalarField & phi
IOMRFZoneList & MRF
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:265
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
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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:188
U
Definition: pEqn.H:72
messageStream Info
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
autoPtr< surfaceVectorField > rhoUf
Creates and initialises the velocity velocity field rhoUf.