rhoPorousSimpleFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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  rhoPorousSimpleFoam
26 
27 Description
28  Steady-state solver for turbulent flow of compressible fluids, with
29  implicit or explicit porosity treatment and optional sources.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "fluidThermo.H"
37 #include "simpleControl.H"
38 #include "pressureReference.H"
39 #include "fvModels.H"
40 #include "fvConstraints.H"
41 #include "IOporosityModelList.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCaseLists.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControl.H"
53  #include "createFields.H"
54  #include "createZones.H"
55  #include "initContinuityErrs.H"
56 
57  turbulence->validate();
58 
59  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61  Info<< "\nStarting time loop\n" << endl;
62 
63  while (simple.loop(runTime))
64  {
65  Info<< "Time = " << runTime.userTimeName() << nl << endl;
66 
67  // Pressure-velocity SIMPLE corrector
68  {
69  #include "UEqn.H"
70  #include "EEqn.H"
71  #include "pEqn.H"
72  }
73 
74  turbulence->correct();
75  thermophysicalTransport->correct();
76 
77  runTime.write();
78  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
79  << " ClockTime = " << runTime.elapsedClockTime() << " s"
80  << nl << endl;
81  }
82 
83  Info<< "End\n" << endl;
84 
85  return 0;
86 }
87 
88 
89 // ************************************************************************* //
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
fluidReactionThermophysicalTransportModel & thermophysicalTransport
messageStream Info
Execute application functionObjects to post-process existing results.
simpleControl simple(mesh)