scalarTransportFoam.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-2016 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  scalarTransportFoam
26 
27 Description
28  Solves the steady or transient transport equation for a passive scalar.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "fvOptions.H"
34 #include "simpleControl.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 int main(int argc, char *argv[])
39 {
40  #include "setRootCase.H"
41  #include "createTime.H"
42  #include "createMesh.H"
43 
44  simpleControl simple(mesh);
45 
46  #include "createFields.H"
47  #include "createFvOptions.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< "\nCalculating scalar transport\n" << endl;
52 
53  #include "CourantNo.H"
54 
55  while (simple.loop())
56  {
57  Info<< "Time = " << runTime.timeName() << nl << endl;
58 
59  while (simple.correctNonOrthogonal())
60  {
62  (
63  fvm::ddt(T)
64  + fvm::div(phi, T)
65  - fvm::laplacian(DT, T)
66  ==
67  fvOptions(T)
68  );
69 
70  TEqn.relax();
71  fvOptions.constrain(TEqn);
72  TEqn.solve();
73  fvOptions.correct(T);
74  }
75 
76  runTime.write();
77  }
78 
79  Info<< "End\n" << endl;
80 
81  return 0;
82 }
83 
84 
85 // ************************************************************************* //
surfaceScalarField & phi
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
fv::options & fvOptions
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:262
const volScalarField & T
const dictionary & simple
messageStream Info
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))