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 -------------------------------------------------------------------------------
10 License
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 
32 using 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,
60  const barycentric& coordinates,
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 
81 Foam::DTRMParticle::DTRMParticle(const DTRMParticle& p)
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 (
96  Cloud<DTRMParticle>& spc,
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  const tetIndices tetIs = this->currentTetIndices();
110 
111  scalar f = 1 - stepFraction();
112  const vector s = p1() - p0() - deviationFromMeshCentre();
113  trackToAndHitFace(f*s, f, spc, td);
114 
115  const point p1 = position();
116  vector dsv = p1 - pos0;
117  scalar ds = mag(dsv);
118 
119  //const label cell1 = cell();
120 
121  //NOTE:
122  // Under the new barocentric tracking alghorithm the newly
123  // inserted particles are tracked to the nearest cell centre first,
124  // then, given the direction, to a face. In both occasions the first call
125  // to trackToAndHitFace returns ds = 0. In this case we do an extra
126  // call to trackToAndHitFace to start the tracking.
127  // This is a temporary fix until the tracking can handle it.
128  if (ds == 0)
129  {
130  trackToAndHitFace(f*s, f, spc, td);
131  dsv = p1 - position();
132  ds = mag(dsv);
133  }
134 
135  // Boltzman constant
136  const scalar sigma = physicoChemical::sigma.value();
137 
138  label reflectedZoneId = td.relfectedCells()[cell1];
139 
140  if
141  (
142  (reflectedZoneId > -1)
143  && (
144  (transmissiveId_ == -1)
145  || (transmissiveId_ != reflectedZoneId)
146  )
147  )
148  {
149  scalar rho(0);
150 
151  // Create a new reflected particle when the particles is not
152  // transmissive and larger than an absolute I
153  if (I_ > 0.01*I0_ && ds > 0)
154  {
155  vector pDir = dsv/ds;
156 
157  cellPointWeight cpw(mesh(), position(), cell1, face());
158  vector nHat = td.nHatInterp().interpolate(cpw);
159 
160  nHat /= (mag(nHat) + ROOTSMALL);
161  scalar cosTheta(-pDir & nHat);
162 
163  // Only new incoming rays
164  if (cosTheta > SMALL)
165  {
166  vector newDir =
167  td.reflection()
168  [
169  td.relfectedCells()[cell1]
170  ].R(pDir, nHat);
171 
172  // reflectivity
173  rho =
174  min
175  (
176  max
177  (
178  td.reflection()
179  [
180  td.relfectedCells()[cell1]
181  ].rho(cosTheta)
182  , 0.0
183  )
184  , 0.98
185  );
186 
187  //scalar delaM = cbrt(mesh().cellVolumes()[cell0]);
188  scalar delaM = cbrt(mesh().cellVolumes()[cell1]);
189 
190  const point insertP(position() - pDir*0.1*delaM);
191  label cellI = mesh().findCell(insertP);
192 
193  if (cellI > -1)
194  {
195  DTRMParticle* pPtr = new DTRMParticle
196  (
197  mesh(),
198  insertP,
199  insertP + newDir*mesh().bounds().mag(),
200  I_*rho,
201  cellI,
202  dA_,
203  -1
204  );
205 
206  // Add to cloud
207  spc.addParticle(pPtr);
208  }
209  }
210  }
211 
212  // Change transmissiveId of the particle
213  transmissiveId_ = reflectedZoneId;
214 
215  scalar a = td.aInterp().interpolate(pos0, cell1);
216  scalar e = td.eInterp().interpolate(pos0, cell1);
217  scalar E = td.EInterp().interpolate(pos0, cell1);
218  scalar T = td.TInterp().interpolate(pos0, cell1);
219 
220 
221  // Left intensity after reflection
222  const scalar Itran = I_*(1.0 - rho);
223  const scalar I1 =
224  (
225  Itran
226  + ds*(e*sigma*pow4(T)/mathematical::pi + E)
227  ) / (1 + ds*a);
228 
229  td.Q(cell1) += (Itran - max(I1, 0.0))*dA_;
230 
231  I_ = I1;
232 
233  if (I_ <= 0.01*I0_)
234  {
235  stepFraction() = 1.0;
236  break;
237  }
238  }
239  else
240  {
241  scalar a = td.aInterp().interpolate(pos0, cell1);
242  scalar e = td.eInterp().interpolate(pos0, cell1);
243  scalar E = td.EInterp().interpolate(pos0, cell1);
244  scalar T = td.TInterp().interpolate(pos0, cell1);
245 
246  const scalar I1 =
247  (
248  I_
249  + ds*(e*sigma*pow4(T)/mathematical::pi + E)
250  ) / (1 + ds*a);
251 
252  td.Q(cell1) += (I_ - max(I1, 0.0))*dA_;
253 
254  I_ = I1;
255 
256  if ((I_ <= 0.01*I0_))
257  {
258  stepFraction() = 1.0;
259  break;
260  }
261  }
262 
263  }while (td.keepParticle && !td.switchProcessor && stepFraction() < 1);
264 
265  return td.keepParticle;
266 }
267 
268 
270 (
271  Cloud<DTRMParticle>&,
272  trackingData& td
273 )
274 {
275  td.switchProcessor = true;
276 }
277 
278 
280 (
281  Cloud<DTRMParticle>&,
282  trackingData& td
283 )
284 {
285  td.keepParticle = false;
286 }
287 
288 
290 (
291  Cloud<DTRMParticle>&,
292  trackingData& td
293 )
294 {
295  return false;
296 }
297 
298 
299 // ************************************************************************* //
p
volScalarField & p
Definition: createFieldRefs.H:8
s
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))
Definition: gmvOutputSpray.H:25
Foam::DTRMParticle::hitWallPatch
void hitWallPatch(Cloud< DTRMParticle > &, trackingData &td)
Foam::constant
Different types of constants.
Definition: atomicConstants.C:31
Foam::DTRMParticle::hitProcessorPatch
void hitProcessorPatch(Cloud< DTRMParticle > &, trackingData &td)
Foam::DTRMParticle::DTRMParticle
DTRMParticle(const polyMesh &mesh, const vector &position, const vector &targetPosition, const scalar I, const label cellI, const scalar dA, const label transmissiveId)
Foam::particle::currentTetIndices
tetIndices currentTetIndices() const
Definition: particleI.H:258
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:181
DTRMParticle.H
Foam::DTRMParticle::move
bool move(Cloud< DTRMParticle > &, trackingData &, const scalar)
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Definition: hashSets.C:26
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:93
Foam::particle::trackToAndHitFace
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Definition: particleTemplates.C:358
Foam::particle::cell
label cell() const
Definition: particleI.H:142
Foam::DTRMParticle::p0
const point & p0() const
Definition: DTRMParticleI.H:103
Foam::DTRMParticle::p1
const point & p1() const
Definition: DTRMParticleI.H:109
rho
rho
Definition: readInitialConditions.H:88
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
Foam::particle::position
vector position() const
Definition: particleI.H:307
Foam::DTRMParticle::hitPatch
bool hitPatch(Cloud< DTRMParticle > &, trackingData &td)
Foam::particle::mesh
const polyMesh & mesh() const
Definition: particleI.H:130
Foam::particle::deviationFromMeshCentre
vector deviationFromMeshCentre() const
Definition: particle.C:1050
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:47
constants.H
Foam::polyMesh::findCell
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Definition: polyMesh.C:1500
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
f
labelList f(nPoints)
physicoChemicalConstants.H
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::barycentric
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:41
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::point
vector point
Point is a vector.
Definition: point.H:37
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:148
Foam::particle::stepFraction
scalar stepFraction() const
Definition: particleI.H:190
Foam::particle::face
label face() const
Definition: particleI.H:178
Foam::I
static const Identity< scalar > I
Definition: Identity.H:89