overBuoyantPimpleDyMFoam.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) 2019-2022 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  overBuoyantPimpleDymFoam
28 
29 Group
30  grpHeatTransferSolvers
31 
32 Description
33  Transient solver for buoyant, turbulent flow of compressible fluids for
34  ventilation and heat-transfer with overset feature
35 
36  Turbulence is modelled using a run-time selectable compressible RAS or
37  LES model.
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "dynamicFvMesh.H"
43 #include "rhoThermo.H"
45 #include "radiationModel.H"
46 #include "fvOptions.H"
47 #include "pimpleControl.H"
48 #include "pressureControl.H"
49 
50 #include "CorrectPhi.H"
51 #include "cellCellStencilObject.H"
52 #include "localMin.H"
53 #include "oversetAdjustPhi.H"
54 
55 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56 
57 int main(int argc, char *argv[])
58 {
59  argList::addNote
60  (
61  "Transient solver for buoyant, turbulent fluid flow"
62  " of compressible fluids, including radiation."
63  );
64 
65  #include "postProcess.H"
66 
67  #include "addCheckCaseOptions.H"
68  #include "setRootCaseLists.H"
69  #include "createTime.H"
70  #include "createDynamicFvMesh.H"
71  #include "createDyMControls.H"
72  #include "createFields.H"
73  #include "createFieldRefs.H"
74  #include "initContinuityErrs.H"
75 
76  #include "createRhoUfIfPresent.H"
77  #include "createControls.H"
78 
79  #include "compressibleCourantNo.H"
80  #include "setInitialDeltaT.H"
81 
82  turbulence->validate();
83 
84  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86  Info<< "\nStarting time loop\n" << endl;
87 
88  while (runTime.run())
89  {
90  #include "readDyMControls.H"
91 
92  #include "compressibleCourantNo.H"
93  #include "setDeltaT.H"
94 
95  // Store divrhoU from the previous mesh so that it can be mapped
96  // and used in correctPhi to ensure the corrected phi has the
97  // same divergence
98  autoPtr<volScalarField> divrhoU;
99  if (correctPhi)
100  {
101  divrhoU.reset
102  (
103  new volScalarField
104  (
105  "divrhoU",
107  )
108  );
109  }
110 
111 
112  ++runTime;
113 
114  Info<< "Time = " << runTime.timeName() << nl << endl;
115 
116  // --- Pressure-velocity PIMPLE corrector loop
117  while (pimple.loop())
118  {
119  if (pimple.firstIter() || moveMeshOuterCorrectors)
120  {
121  // Do any mesh changes
122  mesh.update();
123 
124  if (mesh.changing())
125  {
126  MRF.update();
127 
128  #include "setCellMask.H"
129  #include "setInterpolatedCells.H"
130  #include "correctRhoPhiFaceMask.H"
131 
132  if (correctPhi)
133  {
134  #include "correctPhi.H"
135  }
136 
137  // Make the fluxes relative to the mesh-motion
139  }
140 
141  if (checkMeshCourantNo)
142  {
143  #include "meshCourantNo.H"
144  }
145  }
146 
147  if (pimple.firstIter())
148  {
149  #include "rhoEqn.H"
150  }
151 
152  #include "UEqn.H"
153  #include "EEqn.H"
154 
155  // --- Pressure corrector loop
156  while (pimple.correct())
157  {
158  #include "pEqn.H"
159  }
160 
161  if (pimple.turbCorr())
162  {
163  turbulence->correct();
164  }
165  }
166 
167  rho = thermo.rho();
168 
169  runTime.write();
170 
171  runTime.printExecutionTime(Info);
172  }
173 
174  Info<< "End\n" << endl;
175 
176  return 0;
177 }
178 
179 
180 // ************************************************************************* //
Sets blocked cells mask field.
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
Sets blocked cells mask field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
IOMRFZoneList & MRF
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
moveMeshOuterCorrectors
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)
Correct phi on new faces C-I faces.
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
Required Classes.
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity...