Loading...
Searching...
No Matches
P1.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-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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\*---------------------------------------------------------------------------*/
28
29#include "P1.H"
30#include "fvmLaplacian.H"
31#include "fvmSup.H"
33#include "scatterModel.H"
34#include "constants.H"
36
37using namespace Foam::constant;
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43 namespace radiation
44 {
47 }
48}
49
50
51// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52
53Foam::radiation::P1::P1(const volScalarField& T)
54:
56 G_
57 (
59 (
60 "G",
61 mesh_.time().timeName(),
62 mesh_,
63 IOobject::MUST_READ,
64 IOobject::AUTO_WRITE
65 ),
66 mesh_
67 ),
68 qr_
69 (
71 (
72 "qr",
73 mesh_.time().timeName(),
74 mesh_,
75 IOobject::READ_IF_PRESENT,
76 IOobject::AUTO_WRITE
77 ),
78 mesh_,
80 ),
81 a_
82 (
84 (
85 "a",
86 mesh_.time().timeName(),
87 mesh_,
88 IOobject::NO_READ,
89 IOobject::AUTO_WRITE
90 ),
91 mesh_,
93 ),
94 e_
95 (
97 (
98 "e",
99 mesh_.time().timeName(),
100 mesh_,
101 IOobject::NO_READ,
102 IOobject::NO_WRITE
103 ),
104 mesh_,
106 ),
107 E_
108 (
110 (
111 "E",
112 mesh_.time().timeName(),
113 mesh_,
114 IOobject::NO_READ,
115 IOobject::NO_WRITE
116 ),
117 mesh_,
119 )
120{}
121
122
123Foam::radiation::P1::P1(const dictionary& dict, const volScalarField& T)
124:
126 G_
127 (
129 (
130 "G",
131 mesh_.time().timeName(),
132 mesh_,
133 IOobject::MUST_READ,
134 IOobject::AUTO_WRITE
135 ),
136 mesh_
137 ),
138 qr_
139 (
141 (
142 "qr",
143 mesh_.time().timeName(),
144 mesh_,
145 IOobject::READ_IF_PRESENT,
146 IOobject::AUTO_WRITE
147 ),
148 mesh_,
150 ),
151 a_
152 (
154 (
155 "a",
156 mesh_.time().timeName(),
157 mesh_,
158 IOobject::NO_READ,
159 IOobject::AUTO_WRITE
160 ),
161 mesh_,
163 ),
164 e_
165 (
167 (
168 "e",
169 mesh_.time().timeName(),
170 mesh_,
171 IOobject::NO_READ,
172 IOobject::NO_WRITE
173 ),
174 mesh_,
176 ),
177 E_
178 (
180 (
181 "E",
182 mesh_.time().timeName(),
183 mesh_,
184 IOobject::NO_READ,
185 IOobject::NO_WRITE
186 ),
187 mesh_,
189 )
190{}
191
192
193// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194
196{
198 {
199 // nothing to read
200
201 return true;
202 }
203
204 return false;
205}
206
207
209{
210 a_ = absorptionEmission_->a();
211 e_ = absorptionEmission_->e();
212 E_ = absorptionEmission_->E();
213 const volScalarField sigmaEff(scatter_->sigmaEff());
214
215 const dimensionedScalar a0("a0", a_.dimensions(), ROOTVSMALL);
216
217 // Construct diffusion
218 const auto tgamma = volScalarField::New
219 (
220 "gammaRad",
221 IOobject::REGISTER, // used by boundary conditions
222 1.0/(3.0*a_ + sigmaEff + a0)
223 );
224 const auto& gamma = tgamma();
225
226 // Solve G transport equation
227 solve
228 (
230 - fvm::Sp(a_, G_)
231 ==
232 - 4.0*(e_*physicoChemical::sigma*pow4(T_)) - E_
233 );
234
235 // Calculate radiative heat flux on boundaries.
236 volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
237 const volScalarField::Boundary& GBf = G_.boundaryField();
238 const volScalarField::Boundary& gammaBf = gamma.boundaryField();
239
240 forAll(mesh_.boundaryMesh(), patchi)
241 {
242 if (!GBf[patchi].coupled())
244 qrBf[patchi] = -gammaBf[patchi]*GBf[patchi].snGrad();
245 }
246 }
247}
248
249
251{
253 (
254 "Rp",
263{
264 const volScalarField::Internal& G = G_();
266 const volScalarField::Internal a = absorptionEmission_->aCont()()();
267
268 return a*G - E;
269}
270
271
272Foam::label Foam::radiation::P1::nBands() const
273{
274 return absorptionEmission_->nBands();
275}
276
277
278// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
DimensionedField< scalar, volMesh > Internal
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Works well for combustion applications where optical thickness, tau is large, i.e....
Definition P1.H:61
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant).
Definition P1.C:255
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4).
Definition P1.C:243
label nBands() const
Number of bands.
Definition P1.C:265
bool read()
Read radiation properties dictionary.
Definition P1.C:188
void calculate()
Solve radiation equation(s).
Definition P1.C:201
Top level model for radiation modelling.
const fvMesh & mesh_
Reference to the mesh database.
virtual bool read()=0
Read radiationProperties dictionary.
autoPtr< scatterModel > scatter_
Scatter model.
autoPtr< absorptionEmissionModel > absorptionEmission_
Absorption/emission model.
const volScalarField & T() const noexcept
Return access to the temperature field.
const volScalarField & T_
Reference to the temperature field.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const scalar gamma
Definition EEqn.H:9
bool coupled
Calculate the matrix for the laplacian of the field.
Calculate the finiteVolume matrix for implicit and explicit sources.
word timeName
Definition getTimeIndex.H:3
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for radiation modelling.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar pow4(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define addToRadiationRunTimeSelectionTables(model)
dictionary dict
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299