buoyantBoussinesqSimpleFoam.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  buoyantBoussinesqSimpleFoam
26 
27 Group
28  grpHeatTransferSolvers
29 
30 Description
31  Steady-state solver for buoyant, turbulent flow of incompressible fluids
32 
33  Uses the Boussinesq approximation:
34  \f[
35  rho_{k} = 1 - beta(T - T_{ref})
36  \f]
37 
38  where:
39  \f$ rho_{k} \f$ = the effective (driving) density
40  beta = thermal expansion coefficient [1/K]
41  T = temperature [K]
42  \f$ T_{ref} \f$ = reference temperature [K]
43 
44  Valid when:
45  \f[
46  \frac{beta(T - T_{ref})}{rho_{ref}} << 1
47  \f]
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include "fvCFD.H"
54 #include "radiationModel.H"
55 #include "fvOptions.H"
56 #include "simpleControl.H"
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 int main(int argc, char *argv[])
62 {
63  #include "setRootCase.H"
64  #include "createTime.H"
65  #include "createMesh.H"
66 
68 
69  #include "createFields.H"
70 // #include "createIncompressibleRadiationModel.H"
71 
73  (
75  );
76 
78  (
79  "rhoCpRef",
81  1.0
82  );
83 
84  if (radiation->radiation())
85  {
86  IOdictionary transportProperties
87  (
88  IOobject
89  (
90  "transportProperties",
91  runTime.constant(),
92  runTime,
93  IOobject::MUST_READ,
94  IOobject::NO_WRITE,
95  false // Do not register
96  )
97  );
98 
99  dimensionedScalar rhoRef
100  (
101  "rhoRef",
102  dimDensity,
103  transportProperties
104  );
105 
106  dimensionedScalar CpRef
107  (
108  "CpRef",
110  transportProperties
111  );
112 
113  rhoCpRef = rhoRef*CpRef;
114  }
115 
116  //Info << "form the radiation model, get the Qr is" << Qr << endl;
117 
118  #include "createMRF.H"
119  #include "createFvOptions.H"
120  #include "initContinuityErrs.H"
121 
122  turbulence->validate();
123 
124  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126  Info<< "\nStarting time loop\n" << endl;
127 
128  while (simple.loop())
129  {
130  Info<< "Time = " << runTime.timeName() << nl << endl;
131 
132  // Pressure-velocity SIMPLE corrector
133 
134  #include "UEqn.H"
135  #include "TEqn.H"
136  #include "pEqn.H"
137 
138 
140  turbulence->correct();
141 
142  runTime.write();
143 
144  alphaEff.write();
145  alphaTimeAverage.write();
146 
147  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
148  << " ClockTime = " << runTime.elapsedClockTime() << " s"
149  << nl << endl;
150  }
151 //#include "TcoeffEqn.H"
152 
153  Info<< "End\n" << endl;
154 
155  return 0;
156 }
157 
158 
159 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
fvOptions.H
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::simpleControl
SIMPLE control class to supply convergence information/checks for the SIMPLE loop.
Definition: simpleControl.H:46
pEqn.H
Foam::dimEnergy
const dimensionSet dimEnergy
singlePhaseTransportModel.H
Foam::dimDensity
const dimensionSet dimDensity
turbulentTransportModel.H
UEqn.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
simple
Simple relative velocity model.
createFvOptions.H
TEqn.H
rhoCpRef
dimensionedScalar rhoCpRef("rhoCpRef", dimDensity *dimEnergy/dimMass/dimTemperature, 1.0)
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::dimSpecificHeatCapacity
const dimensionSet dimSpecificHeatCapacity(dimGasConstant)
Definition: dimensionSets.H:71
simpleControl.H
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
alphaTimeAverage
volScalarField alphaTimeAverage("alpha", turbulence->nu()/Pr)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
alphaEff
const volScalarField & alphaEff
Definition: setAlphaEff.H:93
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
setRootCase.H
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
T
const volScalarField & T
Definition: createFields.H:25
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
createMesh.H
createTime.H
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::singlePhaseTransportModel::correct
virtual void correct()
Correct the laminar viscosity.
Definition: singlePhaseTransportModel.C:76
fvCFD.H
main
int main(int argc, char *argv[])
Definition: buoyantBoussinesqSimpleFoam.C:58
fixedFluxPressureFvPatchScalarField.H
createFields.H
radiationModel.H