SRFPimpleFoam.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 | Copyright (C) 2011-2015 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  SRFPimpleFoam
26 
27 Group
28  grpIncompressibleSolvers
29 
30 Description
31  Large time-step transient solver for incompressible, flow in a single
32  rotating frame using the PIMPLE (merged PISO-SIMPLE) algorithm.
33 
34  Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
41 #include "pimpleControl.H"
42 #include "SRFModel.H"
43 #include "fvOptions.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  #include "setRootCase.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52 
53  pimpleControl pimple(mesh);
54 
55  #include "createTimeControls.H"
56  #include "createFields.H"
57  #include "createFvOptions.H"
58  #include "initContinuityErrs.H"
59 
60  turbulence->validate();
61 
62  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64  Info<< "\nStarting time loop\n" << endl;
65 
66  while (runTime.run())
67  {
68  #include "readTimeControls.H"
69  #include "CourantNo.H"
70  #include "setDeltaT.H"
71 
72  runTime++;
73 
74  Info<< "Time = " << runTime.timeName() << nl << endl;
75 
76  // --- Pressure-velocity PIMPLE corrector loop
77  while (pimple.loop())
78  {
79  #include "UrelEqn.H"
80 
81  // --- Pressure corrector loop
82  while (pimple.correct())
83  {
84  #include "pEqn.H"
85  }
86 
87  // Update the absolute velocity
88  U = Urel + SRF->U();
89 
90  if (pimple.turbCorr())
91  {
92  laminarTransport.correct();
93  turbulence->correct();
94  }
95  }
96 
97  runTime.write();
98 
99  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
100  << " ClockTime = " << runTime.elapsedClockTime() << " s"
101  << nl << endl;
102  }
103 
104  Info<< "End\n" << endl;
105 
106  return 0;
107 }
108 
109 
110 // ************************************************************************* //
SRF
Info<< "Reading field p\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Urel\n"<< endl;volVectorField Urel(IOobject("Urel", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading/calculating face flux field phi\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), linearInterpolate(Urel) &mesh.Sf());label pRefCell=0;scalar pRefValue=0.0;setRefCell(p, simple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());Info<< "Creating SRF model\n"<< endl;autoPtr< SRF::SRFModel > SRF(SRF::SRFModel::New(Urel))
fvOptions.H
singlePhaseTransportModel.H
turbulentTransportModel.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
pimpleControl.H
U
U
Definition: pEqn.H:46
createFvOptions.H
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Urel
Urel
Definition: pEqn.H:45
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
SRFModel.H
setRootCase.H
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
createTimeControls.H
Read the control parameters used by setDeltaT.
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
createTime.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1