streamLineParticle.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 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 "streamLineParticle.H"
31#include "vectorFieldIOField.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
36(
37 const trackingData& td,
38 const point& position,
39 const label celli,
40 const label facei
41)
42{
43 if (celli == -1)
44 {
46 << "Cell:" << celli << abort(FatalError);
47 }
48
49 sampledScalars_.setSize(td.vsInterp_.size());
50 forAll(td.vsInterp_, scalari)
51 {
52 sampledScalars_[scalari].append
53 (
54 td.vsInterp_[scalari].interpolate
55 (
57 celli,
58 facei
59 )
60 );
61 }
62
63 sampledVectors_.setSize(td.vvInterp_.size());
64 forAll(td.vvInterp_, vectori)
65 {
66 sampledVectors_[vectori].append
67 (
68 td.vvInterp_[vectori].interpolate
69 (
71 celli,
72 facei
73 )
74 );
75 }
76
77 const DynamicList<vector>& U = sampledVectors_[td.UIndex_];
78
79 return U.last();
80}
81
82
83// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84
86(
87 const polyMesh& mesh,
88 const vector& position,
89 const label celli,
90 const bool trackForward,
91 const label lifeTime
92)
93:
94 particle(mesh, position, celli),
95 trackForward_(trackForward),
96 lifeTime_(lifeTime)
97{}
98
99
101(
102 const polyMesh& mesh,
103 Istream& is,
104 bool readFields,
105 bool newFormat
106)
107:
108 particle(mesh, is, readFields, newFormat)
109{
110 if (readFields)
111 {
112 List<scalarList> sampledScalars;
113 List<vectorList> sampledVectors;
114
115 is >> trackForward_ >> lifeTime_
116 >> sampledPositions_ >> sampledScalars
117 >> sampledVectors;
118
119 sampledScalars_.setSize(sampledScalars.size());
120 forAll(sampledScalars, i)
121 {
122 sampledScalars_[i].transfer(sampledScalars[i]);
123 }
124 sampledVectors_.setSize(sampledVectors.size());
125 forAll(sampledVectors, i)
126 {
127 sampledVectors_[i].transfer(sampledVectors[i]);
128 }
129 }
130
132}
133
134
136(
137 const streamLineParticle& p
138)
139:
140 particle(p),
141 trackForward_(p.trackForward_),
142 lifeTime_(p.lifeTime_),
143 sampledPositions_(p.sampledPositions_),
144 sampledScalars_(p.sampledScalars_),
145 sampledVectors_(p.sampledVectors_)
146{}
147
148
149// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150
152(
154 trackingData& td,
155 const scalar
156)
157{
158 td.switchProcessor = false;
159 td.keepParticle = true;
160
161 const scalar maxDt = mesh().bounds().mag();
162
163 while (td.keepParticle && !td.switchProcessor && lifeTime_ > 0)
164 {
165 scalar dt = maxDt;
166
167 // Cross cell in steps:
168 // - at subiter 0 calculate dt to cross cell in nSubCycle steps
169 // - at the last subiter do all of the remaining track
170 for (label subIter = 0; subIter < max(1, td.nSubCycle_); subIter++)
171 {
172 --lifeTime_;
173
174 // Store current position and sampled velocity.
175 sampledPositions_.append(position());
176 vector U = interpolateFields(td, position(), cell(), face());
177
178 if (!trackForward_)
179 {
180 U = -U;
181 }
182
183 scalar magU = mag(U);
184
185 if (magU < SMALL)
186 {
187 // Stagnant particle. Might as well stop
188 lifeTime_ = 0;
189 break;
190 }
191
192 U /= magU;
193
194 if (td.trackLength_ < GREAT)
195 {
196 // No sub-cycling. Track a set length on each step.
197 dt = td.trackLength_;
198 }
199 else if (subIter == 0)
200 {
201 // Sub-cycling. Cross the cell in nSubCycle steps.
202 particle copy(*this);
203 copy.trackToFace(maxDt*U, 1);
204 dt *= (copy.stepFraction() - stepFraction())/td.nSubCycle_;
205 }
206 else if (subIter == td.nSubCycle_ - 1)
207 {
208 // Sub-cycling. Track the whole cell on the last step.
209 dt = maxDt;
210 }
211
212 trackToAndHitFace(dt*U, 0, cloud, td);
213
214 if
215 (
216 onFace()
217 || !td.keepParticle
218 || td.switchProcessor
219 || lifeTime_ == 0
220 )
221 {
222 break;
223 }
224 }
225 }
226
227 if (!td.keepParticle || lifeTime_ == 0)
228 {
229 if (lifeTime_ == 0)
230 {
231 // Failure exit. Particle stagnated or it's life ran out.
232 if (debug)
233 {
234 Pout<< "streamLineParticle: Removing stagnant particle:"
235 << position() << " sampled positions:"
236 << sampledPositions_.size() << endl;
237 }
238 td.keepParticle = false;
239 }
240 else
241 {
242 // Normal exit. Store last position and fields
243 sampledPositions_.append(position());
244 interpolateFields(td, position(), cell(), face());
245
246 if (debug)
247 {
248 Pout<< "streamLineParticle: Removing particle:" << position()
249 << " sampled positions:" << sampledPositions_.size()
250 << endl;
251 }
252 }
253
254 // Transfer particle data into trackingData.
255 td.allPositions_.append(vectorList());
256 vectorList& top = td.allPositions_.last();
257 top.transfer(sampledPositions_);
258
259 forAll(sampledScalars_, i)
260 {
261 td.allScalars_[i].append(scalarList());
262 scalarList& top = td.allScalars_[i].last();
263 top.transfer(sampledScalars_[i]);
264 }
265 forAll(sampledVectors_, i)
266 {
267 td.allVectors_[i].append(vectorList());
268 vectorList& top = td.allVectors_[i].last();
269 top.transfer(sampledVectors_[i]);
270 }
271 }
272
273 return td.keepParticle;
274}
275
276
278{
279 // Disable generic patch interaction
280 return false;
281}
282
283
285(
287 trackingData& td
288)
289{
290 // Remove particle
291 td.keepParticle = false;
292}
293
294
296(
298 trackingData& td
299)
300{
301 // Remove particle
302 td.keepParticle = false;
303}
304
305
307(
309 trackingData& td
310)
311{
312 // Remove particle
313 td.keepParticle = false;
314}
315
316
318(
320 trackingData& td
321)
322{
323 // Remove particle
324 td.keepParticle = false;
325}
326
327
329(
331 trackingData& td,
332 const vector&
333)
334{
335 // Remove particle
336 td.keepParticle = false;
337}
338
339
341(
343 trackingData& td,
344 const vector&
345)
346{
347 // Remove particle
348 td.keepParticle = false;
349}
350
351
353(
355 trackingData& td
356)
357{
358 // Switch particle
359 td.switchProcessor = true;
360}
361
362
364(
366 trackingData& td
367)
368{
369 // Remove particle
370 td.keepParticle = false;
371}
372
373
375{
376 const bool valid = c.size();
377
379
380 IOField<label> lifeTime
381 (
382 c.fieldIOobject("lifeTime", IOobject::MUST_READ),
383 valid
384 );
385 c.checkFieldIOobject(c, lifeTime);
386
387 vectorFieldIOField sampledPositions
388 (
389 c.fieldIOobject("sampledPositions", IOobject::MUST_READ),
390 valid
391 );
392 c.checkFieldIOobject(c, sampledPositions);
393
394 label i = 0;
395 for (streamLineParticle& p : c)
396 {
397 p.lifeTime_ = lifeTime[i];
398 p.sampledPositions_.transfer(sampledPositions[i]);
399 ++i;
400 }
401}
402
403
405{
407
408 const label np = c.size();
409
410 IOField<label> lifeTime
411 (
412 c.fieldIOobject("lifeTime", IOobject::NO_READ),
413 np
414 );
415 vectorFieldIOField sampledPositions
416 (
417 c.fieldIOobject("sampledPositions", IOobject::NO_READ),
418 np
419 );
420
421 label i = 0;
422 for (const streamLineParticle& p : c)
423 {
424 lifeTime[i] = p.lifeTime_;
425 sampledPositions[i] = p.sampledPositions_;
426 ++i;
427 }
428
429 lifeTime.write(np > 0);
430 sampledPositions.write(np > 0);
431}
432
433
434// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
435
437{
438 os << static_cast<const particle&>(p)
439 << token::SPACE << p.trackForward_
440 << token::SPACE << p.lifeTime_
441 << token::SPACE << p.sampledPositions_
442 << token::SPACE << p.sampledScalars_
443 << token::SPACE << p.sampledVectors_;
444
446 return os;
447}
448
449
450// ************************************************************************* //
Base cloud calls templated on particle type.
Definition: Cloud.H:68
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void transfer(List< T > &list)
Definition: List.C:447
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
virtual void move()=0
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
void interpolateFields()
Linearly interpolate volume fields to generate surface fields.
bool switchProcessor
Flag to switch processor.
Definition: particle.H:102
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
Base particle class.
Definition: particle.H:79
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but also stops on internal faces.
Definition: particle.C:658
vector position() const
Return current particle position.
Definition: particleI.H:314
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition: particleI.H:197
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
virtual bool write(const bool valid=true) const
Write using setting from DB.
A Cloud of streamLine particles.
DynamicList< vectorList > & allPositions_
List< DynamicList< vectorList > > & allVectors_
List< DynamicList< scalarList > > & allScalars_
Particle class that samples fields as it passes through. Used in streamline calculation.
void hitCyclicACMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
void hitSymmetryPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitSymmetryPlanePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitCyclicPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
void hitCyclicAMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamLineParticle > &)
Read.
void hitWedgePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
void hitProcessorPatch(streamLineParticleCloud &, trackingData &)
bool hitPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitWallPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
@ SPACE
Space [isspace].
Definition: token.H:125
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:64
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333