DTRMParticle.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) 2017-2019 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 "DTRMParticle.H"
29#include "constants.H"
31
32using namespace Foam::constant;
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
37(
38 const polyMesh& mesh,
39 const vector& position,
40 const vector& targetPosition,
41 const scalar I,
42 const label cellI,
43 const scalar dA,
44 const label transmissiveId
45)
46:
47 particle(mesh, position, cellI),
48 p0_(position),
49 p1_(targetPosition),
50 I0_(I),
51 I_(I),
52 dA_(dA),
53 transmissiveId_(transmissiveId)
54{}
55
56
58(
59 const polyMesh& mesh,
61 const label celli,
62 const label tetFacei,
63 const label tetPti,
64 const vector& position,
65 const vector& targetPosition,
66 const scalar I,
67 const scalar dA,
68 const label transmissiveId
69)
70:
71 particle(mesh, coordinates, celli, tetFacei, tetPti),
72 p0_(position),
73 p1_(targetPosition),
74 I0_(I),
75 I_(I),
76 dA_(dA),
77 transmissiveId_(transmissiveId)
78{}
79
80
82:
83 particle(p),
84 p0_(p.p0_),
85 p1_(p.p1_),
86 I0_(p.I0_),
87 I_(p.I_),
88 dA_(p.dA_),
89 transmissiveId_(p.transmissiveId_)
90{}
91
92// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93
95(
97 trackingData& td,
98 const scalar trackTime
99)
100{
101 td.switchProcessor = false;
102 td.keepParticle = true;
103
104 do
105 {
106 //Cache old data of particle to use for reflected particle
107 const point pos0 = position();
108 const label cell1 = cell();
109
110 scalar f = 1 - stepFraction();
111 const vector s = p1() - p0() - deviationFromMeshCentre();
112 trackToAndHitFace(f*s, f, spc, td);
113
114 const point p1 = position();
115 vector dsv = p1 - pos0;
116 scalar ds = mag(dsv);
117
118 //const label cell1 = cell();
119
120 //NOTE:
121 // Under the new barocentric tracking alghorithm the newly
122 // inserted particles are tracked to the nearest cell centre first,
123 // then, given the direction, to a face. In both occasions the first call
124 // to trackToAndHitFace returns ds = 0. In this case we do an extra
125 // call to trackToAndHitFace to start the tracking.
126 // This is a temporary fix until the tracking can handle it.
127 if (ds == 0)
128 {
129 trackToAndHitFace(f*s, f, spc, td);
130 dsv = p1 - position();
131 ds = mag(dsv);
132 }
133
134 // Boltzman constant
135 const scalar sigma = physicoChemical::sigma.value();
136
137 label reflectedZoneId = td.relfectedCells()[cell1];
138
139 if
140 (
141 (reflectedZoneId > -1)
142 && (
143 (transmissiveId_ == -1)
144 || (transmissiveId_ != reflectedZoneId)
145 )
146 )
147 {
148 scalar rho(0);
149
150 // Create a new reflected particle when the particles is not
151 // transmissive and larger than an absolute I
152 if (I_ > 0.01*I0_ && ds > 0)
153 {
154 vector pDir = dsv/ds;
155
156 cellPointWeight cpw(mesh(), position(), cell1, face());
157 vector nHat = td.nHatInterp().interpolate(cpw);
158
159 nHat /= (mag(nHat) + ROOTSMALL);
160 scalar cosTheta(-pDir & nHat);
161
162 // Only new incoming rays
163 if (cosTheta > SMALL)
164 {
165 vector newDir =
166 td.reflection()
167 [
168 td.relfectedCells()[cell1]
169 ].R(pDir, nHat);
170
171 // reflectivity
172 rho =
173 min
174 (
175 max
176 (
177 td.reflection()
178 [
179 td.relfectedCells()[cell1]
180 ].rho(cosTheta)
181 , 0.0
182 )
183 , 0.98
184 );
185
186 //scalar delaM = cbrt(mesh().cellVolumes()[cell0]);
187 scalar delaM = cbrt(mesh().cellVolumes()[cell1]);
188
189 const point insertP(position() - pDir*0.1*delaM);
190 label cellI = mesh().findCell(insertP);
191
192 if (cellI > -1)
193 {
194 DTRMParticle* pPtr = new DTRMParticle
195 (
196 mesh(),
197 insertP,
198 insertP + newDir*mesh().bounds().mag(),
199 I_*rho,
200 cellI,
201 dA_,
202 -1
203 );
204
205 // Add to cloud
206 spc.addParticle(pPtr);
207 }
208 }
209 }
210
211 // Change transmissiveId of the particle
212 transmissiveId_ = reflectedZoneId;
213
214 scalar a = td.aInterp().interpolate(pos0, cell1);
215 scalar e = td.eInterp().interpolate(pos0, cell1);
216 scalar E = td.EInterp().interpolate(pos0, cell1);
217 scalar T = td.TInterp().interpolate(pos0, cell1);
218
219
220 // Left intensity after reflection
221 const scalar Itran = I_*(1.0 - rho);
222 const scalar I1 =
223 (
224 Itran
225 + ds*(e*sigma*pow4(T)/mathematical::pi + E)
226 ) / (1 + ds*a);
227
228 td.Q(cell1) += (Itran - max(I1, 0.0))*dA_;
229
230 I_ = I1;
231
232 if (I_ <= 0.01*I0_)
233 {
234 stepFraction() = 1.0;
235 break;
236 }
237 }
238 else
239 {
240 scalar a = td.aInterp().interpolate(pos0, cell1);
241 scalar e = td.eInterp().interpolate(pos0, cell1);
242 scalar E = td.EInterp().interpolate(pos0, cell1);
243 scalar T = td.TInterp().interpolate(pos0, cell1);
244
245 const scalar I1 =
246 (
247 I_
248 + ds*(e*sigma*pow4(T)/mathematical::pi + E)
249 ) / (1 + ds*a);
250
251 td.Q(cell1) += (I_ - max(I1, 0.0))*dA_;
252
253 I_ = I1;
254
255 if ((I_ <= 0.01*I0_))
256 {
257 stepFraction() = 1.0;
258 break;
259 }
260 }
261
262 }while (td.keepParticle && !td.switchProcessor && stepFraction() < 1);
263
264 return td.keepParticle;
265}
266
267
269(
271 trackingData& td
272)
273{
274 td.switchProcessor = true;
275}
276
277
279(
281 trackingData& td
282)
283{
284 td.keepParticle = false;
285}
286
287
289(
291 trackingData& td
292)
293{
294 return false;
295}
296
297
298// ************************************************************************* //
Base cloud calls templated on particle type.
Definition: Cloud.H:68
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:104
Discrete Transfer Radiation Model (DTRM) particle.
Definition: DTRMParticle.H:65
const point & p0() const
Return const access to the initial position.
const point & p1() const
Return const access to the target position.
void hitProcessorPatch(Cloud< DTRMParticle > &, trackingData &td)
Overridable function to handle the particle hitting a processorPatch.
bool hitPatch(Cloud< DTRMParticle > &, trackingData &td)
void hitWallPatch(Cloud< DTRMParticle > &, trackingData &td)
Overridable function to handle the particle hitting a wallPatch.
Foam::cellPointWeight.
const Type & value() const
Return const reference to value.
virtual void move()=0
Base particle class.
Definition: particle.H:79
vector deviationFromMeshCentre() const
Get the displacement from the mesh centre. Used to correct the.
Definition: particle.C:1057
vector position() const
Return current particle position.
Definition: particleI.H:314
label face() const noexcept
Return current face particle is on otherwise -1.
Definition: particleI.H:185
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:137
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition: particleI.H:197
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
label cell() const noexcept
Return current cell particle is in.
Definition: particleI.H:149
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1522
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
PtrList< coordinateSystem > coordinates(solidRegions.size())
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Different types of constants.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
labelList f(nPoints)
volScalarField & e
Definition: createFields.H:11