laplacianFoam.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  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 "fvOptions.H"
34 #include "simpleControl.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 int main(int argc, char *argv[])
39 {
40  #include "setRootCaseLists.H"
41 
42  #include "createTime.H"
43  #include "createMesh.H"
44 
45  simpleControl simple(mesh);
46 
47  #include "createFields.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< "\nCalculating temperature distribution\n" << endl;
52 
53  while (simple.loop(runTime))
54  {
55  Info<< "Time = " << runTime.timeName() << nl << endl;
56 
57  while (simple.correctNonOrthogonal())
58  {
60  (
61  fvm::ddt(T) - fvm::laplacian(DT, T)
62  ==
63  fvOptions(T)
64  );
65 
66  fvOptions.constrain(TEqn);
67  TEqn.solve();
68  fvOptions.correct(T);
69  }
70 
71  #include "write.H"
72 
73  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
74  << " ClockTime = " << runTime.elapsedClockTime() << " s"
75  << nl << endl;
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
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
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
static const char nl
Definition: Ostream.H:265
const volScalarField & T
messageStream Info
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
simpleControl simple(mesh)