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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  laplacianFoam
29 
30 Group
31  grpBasicSolvers
32 
33 Description
34  Laplace equation solver for a scalar quantity.
35 
36  \heading Solver details
37  The solver is applicable to, e.g. for thermal diffusion in a solid. The
38  equation is given by:
39 
40  \f[
41  \ddt{T} = \div \left( D_T \grad T \right)
42  \f]
43 
44  Where:
45  \vartable
46  T | Scalar field which is solved for, e.g. temperature
47  D_T | Diffusion coefficient
48  \endvartable
49 
50  \heading Required fields
51  \plaintable
52  T | Scalar field which is solved for, e.g. temperature
53  \endplaintable
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #include "fvCFD.H"
58 #include "fvOptions.H"
59 #include "simpleControl.H"
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 int main(int argc, char *argv[])
64 {
65  argList::addNote
66  (
67  "Laplace equation solver for a scalar quantity."
68  );
69 
70  #include "postProcess.H"
71 
72  #include "addCheckCaseOptions.H"
73  #include "setRootCaseLists.H"
74  #include "createTime.H"
75  #include "createMesh.H"
76 
77  simpleControl simple(mesh);
78 
79  #include "createFields.H"
80 
81  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
82 
83  Info<< "\nCalculating temperature distribution\n" << endl;
84 
85  while (simple.loop())
86  {
87  Info<< "Time = " << runTime.timeName() << nl << endl;
88 
89  while (simple.correctNonOrthogonal())
90  {
92  (
93  fvm::ddt(T) - fvm::laplacian(DT, T)
94  ==
95  fvOptions(T)
96  );
97 
98  fvOptions.constrain(TEqn);
99  TEqn.solve();
100  fvOptions.correct(T);
101  }
102 
103  #include "write.H"
104 
105  runTime.printExecutionTime(Info);
106  }
107 
108  Info<< "End\n" << endl;
109 
110  return 0;
111 }
112 
113 
114 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
simple
const dictionary & simple
Definition: readFluidMultiRegionSIMPLEControls.H:1
fvOptions.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
TEqn
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
simpleControl.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
createMesh.H
Required Variables.
createTime.H
fvCFD.H
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47