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 Description
28  Solves a simple Laplace equation, e.g. for thermal diffusion in a solid.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "simpleControl.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 int main(int argc, char *argv[])
38 {
39  #include "setRootCase.H"
40 
41  #include "createTime.H"
42  #include "createMesh.H"
43 
44  simpleControl simple(mesh);
45 
46  #include "createFields.H"
47 
48  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50  Info<< "\nCalculating temperature distribution\n" << endl;
51 
52  while (simple.loop())
53  {
54  Info<< "Time = " << runTime.timeName() << nl << endl;
55 
56  while (simple.correctNonOrthogonal())
57  {
58  solve
59  (
60  fvm::ddt(T) - fvm::laplacian(DT, T)
61  );
62  }
63 
64  #include "write.H"
65 
66  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
67  << " ClockTime = " << runTime.elapsedClockTime() << " s"
68  << nl << endl;
69  }
70 
71  Info<< "End\n" << endl;
72 
73  return 0;
74 }
75 
76 
77 // ************************************************************************* //
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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
rhoEqn solve()
static const char nl
Definition: Ostream.H:262
const volScalarField & T
const dictionary & simple
messageStream Info