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 Group
28  grpBasicSolvers
29 
30 Description
31  Passive scalar transport equation solver.
32 
33  \heading Solver details
34  The equation is given by:
35 
36  \f[
37  \ddt{T} + \div \left(\vec{U} T\right) - \div \left(D_T \grad T \right)
38  = S_{T}
39  \f]
40 
41  Where:
42  \vartable
43  T | Passive scalar
44  D_T | Diffusion coefficient
45  S_T | Source
46  \endvartable
47 
48  \heading Required fields
49  \plaintable
50  T | Passive scalar
51  U | Velocity [m/s]
52  \endplaintable
53 
54 \*---------------------------------------------------------------------------*/
55 
56 #include "fvCFD.H"
57 #include "fvOptions.H"
58 #include "simpleControl.H"
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 int main(int argc, char *argv[])
63 {
64  #include "setRootCase.H"
65  #include "createTime.H"
66  #include "createMesh.H"
67 
68  simpleControl simple(mesh);
69 
70  #include "createFields.H"
71  #include "createFvOptions.H"
72 
73  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74 
75  Info<< "\nCalculating scalar transport\n" << endl;
76 
77  #include "CourantNo.H"
78 
79  while (simple.loop())
80  {
81  Info<< "Time = " << runTime.timeName() << nl << endl;
82 
83  while (simple.correctNonOrthogonal())
84  {
86  (
87  fvm::ddt(T)
88  + fvm::div(phi, T)
89  - fvm::laplacian(DT, T)
90  ==
91  fvOptions(T)
92  );
93 
94  fvOptions.constrain(TEqn);
95 
96  TEqn.solve();
97 
98  }
99 
100  runTime.write();
101  }
102 
103  Info<< "End\n" << endl;
104 
105  return 0;
106 }
107 
108 
109 // ************************************************************************* //
fvOptions.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
simple
Simple relative velocity model.
createFvOptions.H
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
simpleControl.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
TEqn
fvScalarMatrix TEqn(fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
setRootCase.H
T
const volScalarField & T
Definition: createFields.H:25
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
createMesh.H
createTime.H
fvCFD.H
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16