overRhoPimpleDyMFoam.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-2015 OpenFOAM Foundation
9  Copyright (C) 2016-2017 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  overRhoPimpleDyMFoam
29 
30 Group
31  grpCompressibleSolvers grpMovingMeshSolvers
32 
33 Description
34  Transient solver for laminar or turbulent flow of compressible fluids
35  for HVAC and similar applications.
36 
37  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
38  pseudo-transient simulations.
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #include "fvCFD.H"
43 #include "dynamicFvMesh.H"
44 #include "fluidThermo.H"
46 #include "bound.H"
47 #include "pimpleControl.H"
48 #include "pressureControl.H"
49 #include "CorrectPhi.H"
50 #include "fvOptions.H"
51 #include "localEulerDdtScheme.H"
52 #include "fvcSmooth.H"
53 #include "cellCellStencilObject.H"
54 #include "localMin.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 int main(int argc, char *argv[])
59 {
60  argList::addNote
61  (
62  "Transient solver for compressible turbulent flow.\n"
63  "With optional mesh motion and mesh topology changes."
64  );
65 
66  #include "setRootCaseLists.H"
67  #include "createTime.H"
68  #include "createDynamicFvMesh.H"
69  #include "createDyMControls.H"
70  #include "createRDeltaT.H"
71  #include "initContinuityErrs.H"
72  #include "createFields.H"
73  #include "createMRF.H"
74  #include "createFvOptions.H"
75  #include "createRhoUfIfPresent.H"
76  #include "createControls.H"
77 
78  turbulence->validate();
79 
80  if (!LTS)
81  {
82  #include "compressibleCourantNo.H"
83  #include "setInitialDeltaT.H"
84  }
85 
86  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87 
88  Info<< "\nStarting time loop\n" << endl;
89 
90  while (runTime.run())
91  {
92  #include "readControls.H"
93  #include "readDyMControls.H"
94 
95 
96  // Store divrhoU from the previous mesh so that it can be mapped
97  // and used in correctPhi to ensure the corrected phi has the
98  // same divergence
99  autoPtr<volScalarField> divrhoU;
100  if (correctPhi)
101  {
102  divrhoU.reset
103  (
104  new volScalarField
105  (
106  "divrhoU",
108  )
109  );
110  }
111 
112  if (LTS)
113  {
114  #include "setRDeltaT.H"
115  }
116  else
117  {
118  #include "compressibleCourantNo.H"
119  #include "setDeltaT.H"
120  }
121 
122  ++runTime;
123 
124  Info<< "Time = " << runTime.timeName() << nl << endl;
125 
126  // --- Pressure-velocity PIMPLE corrector loop
127  while (pimple.loop())
128  {
129  if (pimple.firstIter() || moveMeshOuterCorrectors)
130  {
131 
132  // Do any mesh changes
133  mesh.update();
134 
135  if (mesh.changing())
136  {
137  MRF.update();
138 
139  #include "setCellMask.H"
140 
141  const surfaceScalarField faceMaskOld
142  (
143  localMin<scalar>(mesh).interpolate(cellMask.oldTime())
144  );
145 
146  // Zero Uf on old faceMask (H-I)
147  rhoUf() *= faceMaskOld;
148 
150 
151  // Update Uf and phi on new C-I faces
152  rhoUf() += (1-faceMaskOld)*rhoUfint;
153 
154  // Update Uf boundary
155  forAll(rhoUf().boundaryField(), patchI)
156  {
157  rhoUf().boundaryFieldRef()[patchI] =
158  rhoUfint.boundaryField()[patchI];
159  }
160 
161  // Calculate absolute flux from the mapped surface velocity
162  phi = mesh.Sf() & rhoUf();
163 
164  if (correctPhi)
165  {
166  #include "correctPhi.H"
167  }
168 
169  // Zero phi on current H-I
171  (
172  localMin<scalar>(mesh).interpolate(cellMask)
173  );
174 
175  phi *= faceMask;
176  U *= cellMask;
177 
178  // Make the fluxes relative to the mesh-motion
180 
181  }
182 
183  if (checkMeshCourantNo)
184  {
185  #include "meshCourantNo.H"
186  }
187  }
188 
189  if (pimple.firstIter() && !pimple.SIMPLErho())
190  {
191  #include "rhoEqn.H"
192  }
193 
194  #include "UEqn.H"
195  #include "EEqn.H"
196 
197  // --- Pressure corrector loop
198  while (pimple.correct())
199  {
200  #include "pEqn.H"
201  }
202 
203  if (pimple.turbCorr())
204  {
205  turbulence->correct();
206  }
207  }
208 
209  rho = thermo.rho();
210 
211  runTime.write();
212 
213  runTime.printExecutionTime(Info);
214  }
215 
216  Info<< "End\n" << endl;
217 
218  return 0;
219 }
220 
221 
222 // ************************************************************************* //
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
faceMask
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
createRDeltaT.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
setCellMask.H
Sets blocked cells mask field.
correctPhi
correctPhi
Definition: readDyMControls.H:3
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
createFvOptions.H
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
localEulerDdtScheme.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
bound.H
Bound the given scalar field if it has gone unbounded.
createRhoUfIfPresent.H
Creates and initialises the velocity field rhoUf if required.
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
readDyMControls.H
LTS
bool LTS
Definition: createRDeltaT.H:1
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
pressureControl.H
U
U
Definition: pEqn.H:72
localMin.H
rhoUf
autoPtr< surfaceVectorField > rhoUf
Definition: createRhoUfIfPresent.H:33
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
CorrectPhi.H
createTime.H
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:59
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition: readDyMControls.H:15
cellCellStencilObject.H
createDynamicFvMesh.H
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
turbulentFluidThermoModel.H
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...