overLaplacianDyMFoam.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-2015 OpenFOAM Foundation
9  Copyright (C) 2016-2022 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  overLaplacianDyMFoam
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 #include "dynamicFvMesh.H"
61 #include "oversetPatchPhiErr.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 int main(int argc, char *argv[])
66 {
67  argList::addNote
68  (
69  "Overset Laplace equation solver for a scalar quantity."
70  );
71 
72  #include "setRootCaseLists.H"
73  #include "createTime.H"
75 
76  simpleControl simple(mesh);
77 
78  #include "createFields.H"
79  #include "createFvOptions.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  mesh.update();
90 
91  while (simple.correctNonOrthogonal())
92  {
94  (
95  fvm::ddt(T) - fvm::laplacian(DT, T)
96  ==
97  fvOptions(T)
98  );
99 
100  fvOptions.constrain(TEqn);
101  TEqn.solve();
102  fvOptions.correct(T);
103 
105  {
107  }
108  }
109 
110  #include "write.H"
111 
112  runTime.printExecutionTime(Info);
113  }
114 
115  Info<< "End\n" << endl;
116 
117  return 0;
118 }
119 
120 
121 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
oversetPatchErrOutput
fv::options & fvOptions
dynamicFvMesh & mesh
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
oversetPatchPhiErr
Info<< "Reading field T\"<< endl;volScalarField T(IOobject("T", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading transportProperties\"<< endl;IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));Info<< "Reading diffusivity DT\"<< endl;dimensionedScalar DT("DT", dimViscosity, transportProperties);bool oversetPatchErrOutput=simple.dict().getOrDefault("oversetPatchErrOutput", false);tmp< surfaceScalarField > tdummyPhi
Definition: createFields.H:39
Required Classes.
void oversetPatchPhiErr(const fvScalarMatrix &m, const surfaceScalarField &phi)
Prints out sum of m.flux and phi over the fringe faces.
const volScalarField & T
const dictionary & simple
messageStream Info
Information stream (stdout output on master, null elsewhere)
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))