ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / solvers / basic / laplacianFoam / laplacianFoam.C
blob830a266a4017c97186b9c35a81e00e38fe5cd8cc
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Application
25     laplacianFoam
27 Description
28     Solves a simple Laplace equation, e.g. for thermal diffusion in a solid.
30 \*---------------------------------------------------------------------------*/
32 #include "fvCFD.H"
33 #include "simpleControl.H"
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 int main(int argc, char *argv[])
39     #include "setRootCase.H"
41     #include "createTime.H"
42     #include "createMesh.H"
43     #include "createFields.H"
45     simpleControl simple(mesh);
47     // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49     Info<< "\nCalculating temperature distribution\n" << endl;
51     while (simple.loop())
52     {
53         Info<< "Time = " << runTime.timeName() << nl << endl;
55         for (int nonOrth=0; nonOrth<=simple.nNonOrthCorr(); nonOrth++)
56         {
57             solve
58             (
59                 fvm::ddt(T) - fvm::laplacian(DT, T)
60             );
61         }
63         #include "write.H"
65         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
66             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
67             << nl << endl;
68     }
70     Info<< "End\n" << endl;
72     return 0;
76 // ************************************************************************* //