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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Application
27 SRFPimpleFoam
28
29Group
30 grpIncompressibleSolvers
31
32Description
33 Large time-step transient solver for incompressible flow
34 in a single rotating frame.
35
36 Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
37
38\*---------------------------------------------------------------------------*/
39
40#include "fvCFD.H"
43#include "pimpleControl.H"
44#include "SRFModel.H"
45#include "fvOptions.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49int main(int argc, char *argv[])
50{
51 argList::addNote
52 (
53 "Large time-step transient solver for incompressible flow"
54 " in a single rotating frame."
55 );
56
57 #include "postProcess.H"
58
59 #include "addCheckCaseOptions.H"
60 #include "setRootCaseLists.H"
61 #include "createTime.H"
62 #include "createMesh.H"
63 #include "createControl.H"
64 #include "createTimeControls.H"
65 #include "createFields.H"
66 #include "initContinuityErrs.H"
67
68 turbulence->validate();
69
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71
72 Info<< "\nStarting time loop\n" << endl;
73
74 while (runTime.run())
75 {
76 #include "readTimeControls.H"
77 #include "CourantNo.H"
78 #include "setDeltaT.H"
79
80 ++runTime;
81
82 Info<< "Time = " << runTime.timeName() << nl << endl;
83
84 // --- Pressure-velocity PIMPLE corrector loop
85 while (pimple.loop())
86 {
87 #include "UrelEqn.H"
88
89 // --- Pressure corrector loop
90 while (pimple.correct())
91 {
92 #include "pEqn.H"
93 }
94
95 // Update the absolute velocity
96 U = Urel + SRF->U();
97
98 if (pimple.turbCorr())
99 {
100 laminarTransport.correct();
101 turbulence->correct();
102 }
103 }
104
105 runTime.write();
106
107 runTime.printExecutionTime(Info);
108 }
109
110 Info<< "End\n" << endl;
111
112 return 0;
113}
114
115
116// ************************************************************************* //
Required Classes.
pimpleControl & pimple
U
Definition: pEqn.H:72
engineTime & runTime
Required Variables.
Read the control parameters used by setDeltaT.
compressible::turbulenceModel & turbulence
Urel
Definition: pEqn.H:56
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.
singlePhaseTransportModel laminarTransport(U, phi)
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))