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