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-2018 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"
50 #include "localMin.H"
51 #include "interpolationCellPoint.H"
52 #include "transform.H"
53 #include "fvMeshSubset.H"
54 #include "oversetAdjustPhi.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 int main(int argc, char *argv[])
59 {
60  argList::addNote
61  (
62  "Transient solver for incompressible, turbulent flow"
63  " on a moving mesh."
64  );
65 
66  #include "postProcess.H"
67 
68  #include "setRootCaseLists.H"
69  #include "createTime.H"
70  #include "createDynamicFvMesh.H"
71  #include "initContinuityErrs.H"
72 
73  pimpleControl pimple(mesh);
74 
75  #include "createFields.H"
76  #include "createUf.H"
77  #include "createMRF.H"
78  #include "createFvOptions.H"
79  #include "createControls.H"
80  #include "CourantNo.H"
81  #include "setInitialDeltaT.H"
82 
83  turbulence->validate();
84 
85  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87  Info<< "\nStarting time loop\n" << endl;
88 
89  while (runTime.run())
90  {
91  #include "readControls.H"
92  #include "CourantNo.H"
93 
94  #include "setDeltaT.H"
95 
96  ++runTime;
97 
98  Info<< "Time = " << runTime.timeName() << nl << endl;
99 
100  bool changed = mesh.update();
101 
102  if (changed)
103  {
104  #include "setCellMask.H"
105  #include "setInterpolatedCells.H"
106 
107  surfaceScalarField faceMaskOld
108  (
109  localMin<scalar>(mesh).interpolate(cellMask.oldTime())
110  );
111 
112  // Zero Uf on old faceMask (H-I)
113  Uf *= faceMaskOld;
114  // Update Uf and phi on new C-I faces
115  Uf += (1-faceMaskOld)*fvc::interpolate(U);
116  phi = mesh.Sf() & Uf;
117 
118  // Zero phi on current H-I
120  (
121  localMin<scalar>(mesh).interpolate(cellMask)
122  );
123  phi *= faceMask;
124  }
125 
126 
127  if (mesh.changing() && correctPhi)
128  {
129  // Calculate absolute flux from the mapped surface velocity
130  #include "correctPhi.H"
131  }
132 
133  // Make the flux relative to the mesh motion
135 
136  if (mesh.changing() && checkMeshCourantNo)
137  {
138  #include "meshCourantNo.H"
139  }
140 
141  // --- Pressure-velocity PIMPLE corrector loop
142  while (pimple.loop())
143  {
144  #include "UEqn.H"
145 
146  // --- Pressure corrector loop
147  while (pimple.correct())
148  {
149  #include "pEqn.H"
150  }
151 
152  if (pimple.turbCorr())
153  {
154  laminarTransport.correct();
155  turbulence->correct();
156  }
157  }
158 
159  runTime.write();
160 
161  runTime.printExecutionTime(Info);
162  }
163 
164  Info<< "End\n" << endl;
165 
166  return 0;
167 }
168 
169 
170 // ************************************************************************* //
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
singlePhaseTransportModel.H
fvMeshSubset.H
turbulentTransportModel.H
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
faceMask
surfaceScalarField faceMask(localMin< scalar >(mesh).interpolate(cellMask))
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
createUf.H
Creates and initialises the velocity velocity field Uf.
setCellMask.H
Sets blocked cells mask field.
setInterpolatedCells.H
Sets blocked cells mask field.
correctPhi
correctPhi
Definition: readDyMControls.H:3
interpolationCellPoint.H
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
createFvOptions.H
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
U
U
Definition: pEqn.H:72
localMin.H
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.
createTime.H
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
transform.H
3D tensor transformation operations.
cellCellStencilObject.H
zeroGradientFvPatchFields.H
oversetAdjustPhi.H
Adjust the balance of fluxes on the faces between interpolated and calculated to obey continuity.
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
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.