Loading...
Searching...
No Matches
specularRadiationMixedFvPatchScalarField.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) 2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "radiationModel.H"
29#include "fvDOM.H"
32#include "wedgePolyPatch.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace radiation
40{
41
42// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43
44scalar specularRadiationMixedFvPatchScalarField::azimuthAngle
45(
46 const vector& d
47) const
48{
49 return sign(d.y())*Foam::acos(d.x()/Foam::sqrt(sqr(d.x()) + sqr(d.y())));
50}
51
52
53scalar specularRadiationMixedFvPatchScalarField::polarAngle
54(
55 const vector& d
56) const
57{
58 return Foam::acos(d.z()/mag(d));
59}
60
61
62tmp<scalarField> specularRadiationMixedFvPatchScalarField::interpolateI
63(
64 const fvDOM& dom,
65 const label closestRayi
66) const
67{
68 // Calculate the reflected ray for this ray (KE:p. 2)
69 const vector dAve(normalised(dom.IRay(rayID_).dAve()));
70 const vector dSpe(normalised(dAve - 2*(dAve & n_)*n_));
71
72
73 // Fetch the number of polar and azimuthal segments
74 const label nPolar = dom.nTheta();
75 const label nAzimuth = 4*dom.nPhi();
76
77
78 // Find the neighbouring ray indices of the reflected ray in the east
79 // and west
80 // Go to the east
81 const label polari = std::floor(closestRayi/nAzimuth) + 1;
82 label eastRayi = closestRayi + 1;
83 if (eastRayi == nAzimuth*polari)
84 {
85 eastRayi -= nAzimuth;
86 }
87 // Go to the west
88 label westRayi = closestRayi - 1;
89 if (westRayi < nAzimuth*(polari - 1))
90 {
91 westRayi += nAzimuth;
92 }
93
94 // Find the ray index closest to the reflected ray in the azimuthal
95 // direction
96 label azimuthRayi = -1;
97 bool east = false;
98 const vector dEast(normalised(dom.IRay(eastRayi).dAve()));
99 const vector dWest(normalised(dom.IRay(westRayi).dAve()));
100 if (mag(dSpe - dEast) < mag(dSpe - dWest))
101 {
102 azimuthRayi = eastRayi;
103 east = true;
104 }
105 else
106 {
107 azimuthRayi = westRayi;
108 }
109
110
111 // Find the neighbouring ray indices of the reflected ray in the north
112 // and south
113 // Go to the north
114 label northRayi = closestRayi - nAzimuth;
115 if (northRayi < 0)
116 {
117 // The pole is inside the polar segment - skip
118 northRayi = -1;
119 }
120 // Go to the south
121 label southRayi = closestRayi + nAzimuth;
122 if (southRayi > nPolar*nAzimuth-1)
123 {
124 // The pole is inside the polar segment - skip
125 southRayi = -1;
126 }
127
128
129 // Find the ray index closest to the reflected ray in the polar direction
130 label polarRayi = -1;
131 if (northRayi != -1 && southRayi != -1)
132 {
133 const vector dNorth(normalised(dom.IRay(northRayi).dAve()));
134 const vector dSouth(normalised(dom.IRay(southRayi).dAve()));
135
136 if (mag(dSpe - dNorth) < mag(dSpe - dSouth))
137 {
138 polarRayi = northRayi;
139 }
140 else
141 {
142 polarRayi = southRayi;
143 }
144 }
145 else if (northRayi != -1)
146 {
147 polarRayi = northRayi;
148 }
149 else if (southRayi != -1)
150 {
151 polarRayi = southRayi;
152 }
153
154
155 // Find the ray index neighbouring the azimuthal and polar neighbour rays
156 label cornerRayi = -1;
157 if (polarRayi != -1)
158 {
159 if (east)
160 {
161 cornerRayi = polarRayi + 1;
162
163 if (!(cornerRayi % nAzimuth))
164 {
165 cornerRayi -= nAzimuth;
166 }
167 }
168 else
169 {
170 cornerRayi = polarRayi - 1;
171 if (!(polarRayi % nAzimuth))
172 {
173 cornerRayi += nAzimuth;
174 }
175 }
176 }
177
178
179 // Interpolate the ray intensity of the reflected complementary ray from
180 // the neighbouring rays
181 const label patchi = this->patch().index();
182 auto tIc = tmp<scalarField>::New(this->patch().size(), Zero);
183 auto& Ic = tIc.ref();
184
185
186 if (polarRayi == -1)
187 {
188 // Linear interpolation only in the azimuth direction
189 const vector d1(normalised(dom.IRay(closestRayi).dAve()));
190 const vector d2(normalised(dom.IRay(azimuthRayi).dAve()));
191
192 const scalar phic = azimuthAngle(dSpe);
193 const scalar phi1 = azimuthAngle(d1);
194 const scalar phi2 = azimuthAngle(d2);
195
196 const auto& I1 =
197 dom.IRayLambda(closestRayi, lambdaID_).boundaryField()[patchi];
198 const auto& I2 =
199 dom.IRayLambda(azimuthRayi, lambdaID_).boundaryField()[patchi];
200
201 Ic = lerp(I1, I2, (phic - phi1)/(phi2 - phi1));
202 }
203 else
204 {
205 // Bilinear interpolation in the azimuth and polar directions
206 const vector d1(normalised(dom.IRay(closestRayi).dAve()));
207 const vector d2(normalised(dom.IRay(azimuthRayi).dAve()));
208 const vector d3(normalised(dom.IRay(polarRayi).dAve()));
209 const vector d4(normalised(dom.IRay(cornerRayi).dAve()));
210
211 const scalar phic = azimuthAngle(dSpe);
212 const scalar phi1 = azimuthAngle(d1);
213 const scalar phi2 = azimuthAngle(d2);
214 const scalar phi3 = azimuthAngle(d3);
215 const scalar phi4 = azimuthAngle(d4);
216
217 const scalar thetac = polarAngle(dSpe);
218 const scalar theta1 = polarAngle(d1);
219 const scalar theta3 = polarAngle(d3);
220
221 const auto& I1 =
222 dom.IRayLambda(closestRayi, lambdaID_).boundaryField()[patchi];
223 const auto& I2 =
224 dom.IRayLambda(azimuthRayi, lambdaID_).boundaryField()[patchi];
225 const auto& I3 =
226 dom.IRayLambda(polarRayi, lambdaID_).boundaryField()[patchi];
227 const auto& I4 =
228 dom.IRayLambda(cornerRayi, lambdaID_).boundaryField()[patchi];
229
230 const scalarField Ia(lerp(I1, I2, (phic - phi1)/(phi2 - phi1)));
231 const scalarField Ib(lerp(I3, I4, (phic - phi3)/(phi4 - phi3)));
232
233 Ic = lerp(Ia, Ib, (thetac - theta1)/(theta3 - theta1));
234 }
235
236 return tIc;
237}
238
239
240label specularRadiationMixedFvPatchScalarField::calcComplementaryRayID
241(
242 const fvDOM& dom
243) const
244{
245 vectorList dAve(dom.nRay());
246 forAll(dAve, rayi)
247 {
248 dAve[rayi] = normalised(dom.IRay(rayi).dAve());
249 }
250
251
252 // Check if the ray goes out of the domain, ie. no reflection
253 if ((dAve[rayID_] & n_) > 0)
254 {
255 return -1;
256 }
257
258
259 // Calculate the reflected ray direction for rayID (KE:p. 2)
260 const vector dSpe
261 (
262 normalised(dAve[rayID_] - 2*(dAve[rayID_] & n_)*n_)
263 );
264
265
266 // Find the closest ray to the reflected ray
267 label complementaryRayID = -1;
268 scalar dotProductThisRay = -GREAT;
269 forAll(dAve, i)
270 {
271 const scalar dotProductOtherRay = dAve[i] & dSpe;
272
273 if (dotProductThisRay < dotProductOtherRay)
274 {
275 complementaryRayID = i;
276 dotProductThisRay = dotProductOtherRay;
277 }
278 }
280 return complementaryRayID;
281}
282
283
284// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
285
288(
289 const fvPatch& p,
291)
292:
293 mixedFvPatchScalarField(p, iF),
294 n_(),
295 rayID_(-1),
296 lambdaID_(-1),
297 interpolate_(false)
299 this->refValue() = Zero;
300 this->refGrad() = Zero;
301 this->valueFraction() = Zero;
302}
303
304
307(
309 const fvPatch& p,
311 const fvPatchFieldMapper& mapper
312)
313:
314 mixedFvPatchScalarField(ptf, p, iF, mapper),
315 n_(ptf.n_),
316 rayID_(ptf.rayID_),
317 lambdaID_(ptf.lambdaID_),
318 interpolate_(ptf.interpolate_)
319{}
320
321
324(
325 const fvPatch& p,
327 const dictionary& dict
328)
329:
330 mixedFvPatchScalarField(p, iF),
331 n_(),
332 rayID_(-1),
333 lambdaID_(-1),
334 interpolate_(dict.getOrDefault("interpolate", false))
335{
336 this->refValue() = Zero;
337 this->refGrad() = Zero;
338 this->valueFraction() = Zero;
339
340 if (!this->readValueEntry(dict))
341 {
342 fvPatchScalarField::operator=(this->refValue());
343 }
344
345
346 if (isA<wedgePolyPatch>(p.patch()))
347 {
348 const auto& wp = refCast<const wedgePolyPatch>(p.patch());
349 n_ = wp.n();
350 }
351 else if (isA<symmetryPlanePolyPatch>(p.patch()))
352 {
353 const auto& sp = refCast<const symmetryPlanePolyPatch>(p.patch());
354 n_ = sp.n();
355 }
356 else
357 {
359 << " specularRadiation boundary condition is limited to "
360 << "wedge or symmetry-plane geometries." << nl
361 << exit(FatalError);
362 }
363}
364
365
368(
371)
372:
373 mixedFvPatchScalarField(ptf, iF),
374 n_(ptf.n_),
375 rayID_(ptf.rayID_),
376 lambdaID_(ptf.lambdaID_),
377 interpolate_(ptf.interpolate_)
378{}
379
380
383(
385)
386:
387 mixedFvPatchScalarField(ptf),
388 n_(ptf.n_),
389 rayID_(ptf.rayID_),
390 lambdaID_(ptf.lambdaID_),
391 interpolate_(ptf.interpolate_)
392{}
393
394
395// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
396
398{
399 if (this->updated())
400 {
401 return;
402 }
403
404 const fvDOM& dom = db().lookupObject<fvDOM>("radiationProperties");
405
406 // Get rayID and lambdaID for this ray
407 if (rayID_ == -1 && lambdaID_ == -1)
408 {
409 dom.setRayIdLambdaId(internalField().name(), rayID_, lambdaID_);
410 }
411
412
413 // Find the ray index closest to the reflected ray
414 const label compID = calcComplementaryRayID(dom);
415
416
417 if (compID == -1)
418 {
419 // Apply zero-gradient condition for rays outgoing from the domain
420 this->valueFraction() = 0;
421 }
422 else
423 {
424 // Apply fixed condition for rays incoming to the domain
425 this->valueFraction() = 1;
426
427 if (!interpolate_)
428 {
429 // Fetch the existing ray closest to the reflected ray (KE:Eq. 4)
430 this->refValue() =
431 dom.IRayLambda
432 (
433 compID,
434 lambdaID_
435 ).internalField();
436 }
437 else
438 {
439 // Interpolate the ray intensity from neighbouring rays (KE:p. 2)
440 this->refValue() = interpolateI(dom, compID);
442 }
443
444 mixedFvPatchScalarField::updateCoeffs();
445}
446
447
449{
450 mixedFvPatchScalarField::write(os);
451 os.writeEntryIfDifferent<bool>("interpolate", false, interpolate_);
452 this->writeValueEntry(os);
453}
454
455
456// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
457
459(
462);
463
464
465} // End namespace radiation
466} // End namespace Foam
467
468// ************************************************************************* //
surfaceScalarField phic(mixture.cAlpha() *mag(alphaPhic/mesh.magSf()))
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi2
surfaceScalarField & phi1
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Internal & internalField() const noexcept
Return a const-reference to the dimensioned internal field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void operator=(const UList< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
Definition fvDOM.H:117
const volScalarField & IRayLambda(const label rayI, const label lambdaI) const
Ray intensity for rayI and lambda bandwidth.
Definition fvDOM.H:32
void setRayIdLambdaId(const word &name, label &rayId, label &lambdaId) const
Set the rayId and lambdaId from by decomposing an intensity.
Definition fvDOM.C:738
This boundary condition provides a specular radiation condition for axisymmetric and symmetry-plane f...
specularRadiationMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
auto & name
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for radiation modelling.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
dimensionedScalar sign(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
List< vector > vectorList
List of vector.
Definition vectorList.H:32
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< Type > lerp(const dimensioned< Type > &a, const dimensioned< Type > &b, const scalar t)
dimensionedScalar acos(const dimensionedScalar &ds)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299