laplacianFoam.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-2015 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  laplacianFoam
26 
27 Group
28  grpBasicSolvers
29 
30 Description
31  Laplace equation solver for a scalar quantity.
32 
33  \heading Solver details
34  The solver is applicable to, e.g. for thermal diffusion in a solid. The
35  equation is given by:
36 
37  \f[
38  \ddt{T} = \div \left( D_T \grad T \right)
39  \f]
40 
41  Where:
42  \vartable
43  T | Scalar field which is solved for, e.g. temperature
44  D_T | Diffusion coefficient
45  \endvartable
46 
47  \heading Required fields
48  \plaintable
49  T | Scalar field which is solved for, e.g. temperature
50  \endplaintable
51 
52 \*---------------------------------------------------------------------------*/
53 
54 #include "fvCFD.H"
55 #include "simpleControl.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 int main(int argc, char *argv[])
60 {
61  #include "setRootCase.H"
62 
63  #include "createTime.H"
64  #include "createMesh.H"
65 
66  simpleControl simple(mesh);
67 
68  #include "createFields.H"
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nCalculating temperature distribution\n" << endl;
73 
74  while (simple.loop())
75  {
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  while (simple.correctNonOrthogonal())
79  {
80  solve
81  (
82  fvm::ddt(T) - fvm::laplacian(DT, T)
83  );
84  }
85 
86  #include "write.H"
87 
88  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
89  << " ClockTime = " << runTime.elapsedClockTime() << " s"
90  << nl << endl;
91  }
92 
93  Info<< "End\n" << endl;
94 
95  return 0;
96 }
97 
98 
99 // ************************************************************************* //
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
simple
Simple relative velocity model.
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
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
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