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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
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 
26 Application
27  scalarTransportFoam
28 
29 Group
30  grpBasicSolvers
31 
32 Description
33  Passive scalar transport equation solver.
34 
35  \heading Solver details
36  The equation is given by:
37 
38  \f[
39  \ddt{T} + \div \left(\vec{U} T\right) - \div \left(D_T \grad T \right)
40  = S_{T}
41  \f]
42 
43  Where:
44  \vartable
45  T | Passive scalar
46  D_T | Diffusion coefficient
47  S_T | Source
48  \endvartable
49 
50  \heading Required fields
51  \plaintable
52  T | Passive scalar
53  U | Velocity [m/s]
54  \endplaintable
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #include "fvCFD.H"
59 #include "fvOptions.H"
60 #include "simpleControl.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 int main(int argc, char *argv[])
65 {
66  argList::addNote
67  (
68  "Passive scalar transport equation solver."
69  );
70 
71  #include "addCheckCaseOptions.H"
72  #include "setRootCaseLists.H"
73  #include "createTime.H"
74  #include "createMesh.H"
75 
76  simpleControl simple(mesh);
77 
78  #include "createFields.H"
79 
80  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82  Info<< "\nCalculating scalar transport\n" << endl;
83 
84  #include "CourantNo.H"
85 
86  while (simple.loop())
87  {
88  Info<< "Time = " << runTime.timeName() << nl << endl;
89 
90  while (simple.correctNonOrthogonal())
91  {
93  (
94  fvm::ddt(T)
95  + fvm::div(phi, T)
96  - fvm::laplacian(DT, T)
97  ==
98  fvOptions(T)
99  );
100 
101  TEqn.relax();
102  fvOptions.constrain(TEqn);
103  TEqn.solve();
104  fvOptions.correct(T);
105  }
106 
107  runTime.write();
108  }
109 
110  Info<< "End\n" << endl;
111 
112  return 0;
113 }
114 
115 
116 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
fvOptions.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
simpleControl.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
createMesh.H
Required Variables.
createTime.H
fvCFD.H
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47