overPimpleDyMFoam.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-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2022 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  overPimpleDyMFoam
29 
30 Group
31  grpIncompressibleSolvers grpMovingMeshSolvers
32 
33 Description
34  Transient solver for incompressible flow of Newtonian fluids
35  on a moving mesh using the PIMPLE (merged PISO-SIMPLE) algorithm.
36 
37  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "dynamicFvMesh.H"
45 #include "pimpleControl.H"
46 #include "fvOptions.H"
47 
48 #include "cellCellStencilObject.H"
49 #include "localMin.H"
50 #include "oversetAdjustPhi.H"
51 #include "oversetPatchPhiErr.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  argList::addNote
58  (
59  "Transient solver for incompressible, turbulent flow"
60  " on a moving mesh."
61  );
62 
63  #include "postProcess.H"
64 
65  #include "setRootCaseLists.H"
66  #include "createTime.H"
67  #include "createDynamicFvMesh.H"
68  #include "createDyMControls.H"
69  #include "initContinuityErrs.H"
70 
71  #include "createFields.H"
72  #include "createUf.H"
73  #include "createMRF.H"
74  #include "createFvOptions.H"
75  #include "createControls.H"
76  #include "CourantNo.H"
77  #include "setInitialDeltaT.H"
78 
79  turbulence->validate();
80 
81  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83  Info<< "\nStarting time loop\n" << endl;
84 
85  while (runTime.run())
86  {
87  #include "readDyMControls.H"
88  #include "readOversetDyMControls.H"
89 
90  #include "CourantNo.H"
91 
92  #include "setDeltaT.H"
93 
94  ++runTime;
95 
96  Info<< "Time = " << runTime.timeName() << nl << endl;
97 
98  mesh.update();
99 
100  if (mesh.changing())
101  {
102  #include "setCellMask.H"
103  #include "setInterpolatedCells.H"
104  #include "correctPhiFaceMask.H"
105 
107 
108  if (checkMeshCourantNo)
109  {
110  #include "meshCourantNo.H"
111  }
112  }
113 
114  // --- Pressure-velocity PIMPLE corrector loop
115  while (pimple.loop())
116  {
117  #include "UEqn.H"
118 
119  // --- Pressure corrector loop
120  while (pimple.correct())
121  {
122  #include "pEqn.H"
123  }
124 
125  if (pimple.turbCorr())
126  {
127  laminarTransport.correct();
128  turbulence->correct();
129  }
130  }
131 
132  runTime.write();
133 
134  runTime.printExecutionTime(Info);
135  }
136 
137  Info<< "End\n" << endl;
138 
139  return 0;
140 }
141 
142 
143 // ************************************************************************* //
Sets blocked cells mask field.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
checkMeshCourantNo
Sets blocked cells mask field.
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
oversetPatchPhiErr
Creates and initialises the velocity velocity field Uf.
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
Correct phi on new faces C-I faces.
Calculates and outputs the mean and maximum Courant Numbers.
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
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity...