scalarTransportFoam.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-2018 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 "setRootCaseLists.H"
41  #include "createTime.H"
42  #include "createMesh.H"
43 
44  simpleControl simple(mesh);
45 
46  #include "createFields.H"
47 
48  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50  Info<< "\nCalculating scalar transport\n" << endl;
51 
52  #include "CourantNo.H"
53 
54  while (simple.loop(runTime))
55  {
56  Info<< "Time = " << runTime.timeName() << nl << endl;
57 
58  while (simple.correctNonOrthogonal())
59  {
60  fvScalarMatrix TEqn
61  (
62  fvm::ddt(T)
63  + fvm::div(phi, T)
64  - fvm::laplacian(DT, T)
65  ==
66  fvOptions(T)
67  );
68 
69  TEqn.relax();
70  fvOptions.constrain(TEqn);
71  TEqn.solve();
72  fvOptions.correct(T);
73  }
74 
75  runTime.write();
76  }
77 
78  Info<< "End\n" << endl;
79 
80  return 0;
81 }
82 
83 
84 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
phi
Definition: pEqn.H:104
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:260
const volScalarField & T
messageStream Info
simpleControl simple(mesh)