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-2021 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 "fvModels.H"
34 #include "fvConstraints.H"
35 #include "simpleControl.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 int main(int argc, char *argv[])
40 {
41  #include "setRootCaseLists.H"
42  #include "createTime.H"
43  #include "createMesh.H"
44 
45  simpleControl simple(mesh);
46 
47  #include "createFields.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< "\nCalculating scalar transport\n" << endl;
52 
53  #include "CourantNo.H"
54 
55  while (simple.loop(runTime))
56  {
57  Info<< "Time = " << runTime.timeName() << nl << endl;
58 
59  fvModels.correct();
60 
61  while (simple.correctNonOrthogonal())
62  {
63  fvScalarMatrix TEqn
64  (
65  fvm::ddt(T)
66  + fvm::div(phi, T)
67  - fvm::laplacian(DT, T)
68  ==
69  fvModels.source(T)
70  );
71 
72  TEqn.relax();
74  TEqn.solve();
76  }
77 
78  runTime.write();
79  }
80 
81  Info<< "End\n" << endl;
82 
83  return 0;
84 }
85 
86 
87 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
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
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
Foam::fvConstraints & fvConstraints
phi
Definition: correctPhi.H:3
static const char nl
Definition: Ostream.H:260
const volScalarField & T
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
messageStream Info
simpleControl simple(mesh)