All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
P1.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-2024 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 \*---------------------------------------------------------------------------*/
25 
26 #include "P1.H"
27 #include "fvmLaplacian.H"
28 #include "fvmSup.H"
30 #include "scatterModel.H"
31 #include "constants.H"
33 
34 using namespace Foam::constant;
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace radiationModels
41 {
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 :
52  radiationModel(typeName, T),
53  G_
54  (
55  IOobject
56  (
57  "G",
58  mesh_.time().name(),
59  mesh_,
60  IOobject::MUST_READ,
61  IOobject::AUTO_WRITE
62  ),
63  mesh_
64  ),
65  qr_
66  (
67  IOobject
68  (
69  "qr",
70  mesh_.time().name(),
71  mesh_,
72  IOobject::READ_IF_PRESENT,
73  IOobject::AUTO_WRITE
74  ),
75  mesh_,
77  ),
78  a_
79  (
80  IOobject
81  (
82  "a",
83  mesh_.time().name(),
84  mesh_,
85  IOobject::NO_READ,
86  IOobject::AUTO_WRITE
87  ),
88  mesh_,
90  ),
91  e_
92  (
93  IOobject
94  (
95  "e",
96  mesh_.time().name(),
97  mesh_,
98  IOobject::NO_READ,
99  IOobject::NO_WRITE
100  ),
101  mesh_,
103  ),
104  E_
105  (
106  IOobject
107  (
108  "E",
109  mesh_.time().name(),
110  mesh_,
111  IOobject::NO_READ,
112  IOobject::NO_WRITE
113  ),
114  mesh_,
116  )
117 {}
118 
119 
121 :
122  radiationModel(typeName, dict, T),
123  G_
124  (
125  IOobject
126  (
127  "G",
128  mesh_.time().name(),
129  mesh_,
130  IOobject::MUST_READ,
131  IOobject::AUTO_WRITE
132  ),
133  mesh_
134  ),
135  qr_
136  (
137  IOobject
138  (
139  "qr",
140  mesh_.time().name(),
141  mesh_,
142  IOobject::READ_IF_PRESENT,
143  IOobject::AUTO_WRITE
144  ),
145  mesh_,
147  ),
148  a_
149  (
150  IOobject
151  (
152  "a",
153  mesh_.time().name(),
154  mesh_,
155  IOobject::NO_READ,
156  IOobject::AUTO_WRITE
157  ),
158  mesh_,
160  ),
161  e_
162  (
163  IOobject
164  (
165  "e",
166  mesh_.time().name(),
167  mesh_,
168  IOobject::NO_READ,
169  IOobject::NO_WRITE
170  ),
171  mesh_,
173  ),
174  E_
175  (
176  IOobject
177  (
178  "E",
179  mesh_.time().name(),
180  mesh_,
181  IOobject::NO_READ,
182  IOobject::NO_WRITE
183  ),
184  mesh_,
186  )
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
191 
193 {}
194 
195 
196 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
197 
199 {
200  if (radiationModel::read())
201  {
202  return true;
203  }
204  else
205  {
206  return false;
207  }
208 }
209 
210 
212 {
213  a_ = absorptionEmission_->a();
214  e_ = absorptionEmission_->e();
215  E_ = absorptionEmission_->E();
216  const volScalarField sigmaEff(scatter_->sigmaEff());
217 
218  const dimensionedScalar a0 ("a0", a_.dimensions(), rootVSmall);
219 
220  // Construct diffusion
221  const volScalarField gamma
222  (
223  IOobject
224  (
225  "gammaRad",
226  G_.mesh().time().name(),
227  G_.mesh(),
230  ),
231  1.0/(3.0*a_ + sigmaEff + a0)
232  );
233 
234  // Solve G transport equation
235  solve
236  (
237  fvm::laplacian(gamma, G_)
238  - fvm::Sp(a_, G_)
239  ==
240  - 4*e_*physicoChemical::sigma*pow4(T_) - E_
241  );
242 
243  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
244 
245  // Calculate radiative heat flux on boundaries.
246  forAll(mesh_.boundaryMesh(), patchi)
247  {
248  if (!G_.boundaryField()[patchi].coupled())
249  {
250  qrBf[patchi] =
251  -gamma.boundaryField()[patchi]
252  *G_.boundaryField()[patchi].snGrad();
253  }
254  }
255 }
256 
257 
259 {
260  return volScalarField::New
261  (
262  "Rp",
263  4.0*absorptionEmission_->eCont()*physicoChemical::sigma
264  );
265 }
266 
267 
270 {
271  const volScalarField::Internal& G = G_();
272  const volScalarField::Internal E = absorptionEmission_->ECont()()();
273  const volScalarField::Internal a = absorptionEmission_->aCont()()();
274 
275  return a*G - E;
276 }
277 
278 
279 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Top level model for radiation modelling.
virtual bool read()=0
Read radiationProperties dictionary.
Works well for combustion applications where optical thickness, tau is large, i.e....
Definition: P1.H:60
P1(const volScalarField &T)
Construct from components.
Definition: P1.C:50
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: P1.C:269
virtual ~P1()
Destructor.
Definition: P1.C:192
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: P1.C:258
bool read()
Read radiation properties dictionary.
Definition: P1.C:198
void calculate()
Solve radiation equation(s)
Definition: P1.C:211
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the matrix for the laplacian of the field.
Calculate the matrix for implicit and explicit sources.
label patchi
const dimensionedScalar a0
Bohr radius: default SI units: [m].
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar G
Newtonian constant of gravitation.
Collection of constants.
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTime
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict