particleTemplates.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,2020 OpenFOAM Foundation
9 Copyright (C) 2016-2020,2022 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 "IOPosition.H"
30
31#include "cyclicPolyPatch.H"
32#include "cyclicAMIPolyPatch.H"
33#include "cyclicACMIPolyPatch.H"
34#include "processorPolyPatch.H"
36#include "symmetryPolyPatch.H"
37#include "wallPolyPatch.H"
38#include "wedgePolyPatch.H"
39#include "meshTools.H"
40
41// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42
43template<class Type>
45(
46 Ostream& os,
47 const word& name,
48 const word& delim
49)
50{
52 {
53 os << name;
54 }
55 else
56 {
57 os << '(';
58 for (int i = 0; i < pTraits<Type>::nComponents; ++i)
59 {
60 if (i) os << delim;
61
62 os << name << Foam::name(i);
63 }
64 os << ')';
65 }
66}
67
68
69template<class Type>
71(
72 Ostream& os,
73 const word& name,
74 const Type& value,
75 const bool nameOnly,
76 const word& delim,
77 const wordRes& filters
78)
79{
80 if (!filters.empty() && !filters.match(name))
81 {
82 return;
83 }
84
85 os << delim;
86 if (nameOnly)
87 {
88 writePropertyName<Type>(os, name, delim);
89 }
90 else
91 {
92 os << value;
93 }
94}
95
96
97template<class Type>
99(
100 Ostream& os,
101 const word& name,
102 const Field<Type>& values,
103 const bool nameOnly,
104 const word& delim,
105 const wordRes& filters
106)
107{
108 if (!filters.empty() && !filters.match(name))
109 {
110 return;
111 }
112
113 if (nameOnly)
114 {
115 os << delim;
116 os << "N(";
117 if (values.size())
118 {
119 forAll(values, i)
120 {
121 if (i) os << delim;
122 const word tag(name + Foam::name(i));
123 writePropertyName<Type>(os, tag, delim);
124 }
125 }
126 else
127 {
128 os << name;
129 }
130 os << ')';
131 }
132 else
133 {
134 os << delim << values;
135 }
136}
137
138
139template<class TrackCloudType>
140void Foam::particle::readFields(TrackCloudType& c)
141{
142 const bool valid = c.size();
143
144 IOobject procIO(c.fieldIOobject("origProcId", IOobject::MUST_READ));
145
146 const bool haveFile = procIO.typeHeaderOk<IOField<label>>(true);
147
148 IOField<label> origProcId(procIO, valid && haveFile);
149 c.checkFieldIOobject(c, origProcId);
150
151 IOField<label> origId
152 (
153 c.fieldIOobject("origId", IOobject::MUST_READ),
154 valid && haveFile
155 );
156 c.checkFieldIOobject(c, origId);
157
158 label i = 0;
159 for (particle& p : c)
160 {
161 p.origProc_ = origProcId[i];
162 p.origId_ = origId[i];
163
164 ++i;
165 }
166}
167
168
169template<class TrackCloudType>
170void Foam::particle::writeFields(const TrackCloudType& c)
171{
172 const label np = c.size();
173 const bool valid = np;
174
175 if (writeLagrangianCoordinates)
176 {
178 ioP.write(valid);
179 }
180 else if (!writeLagrangianPositions)
181 {
183 << "Must select coordinates and/or positions" << nl
184 << exit(FatalError);
185 }
186
187 // Optionally write positions file in v1706 format and earlier
188 if (writeLagrangianPositions)
189 {
191 (
192 c,
194 );
195 ioP.write(valid);
196 }
197
198 IOField<label> origProc
199 (
200 c.fieldIOobject("origProcId", IOobject::NO_READ),
201 np
202 );
203 IOField<label> origId
204 (
205 c.fieldIOobject("origId", IOobject::NO_READ),
206 np
207 );
208
209 label i = 0;
210 for (const particle& p : c)
211 {
212 origProc[i] = p.origProc_;
213 origId[i] = p.origId_;
214
215 ++i;
216 }
217
218 origProc.write(valid);
219 origId.write(valid);
220}
221
222
223template<class CloudType>
225{
226 typedef typename CloudType::parcelType parcelType;
227
228 const auto* positionPtr = cloud::findIOPosition(obr);
229
230 const label np = c.size();
231 const label newNp = (positionPtr ? positionPtr->size() : 0);
232
233 // Remove excess parcels
234 for (label i = newNp; i < np; ++i)
235 {
236 parcelType* p = c.last();
237
238 c.deleteParticle(*p);
239 }
240
241 if (newNp)
242 {
243 const auto& position = *positionPtr;
244
245 const auto& origProcId = cloud::lookupIOField<label>("origProc", obr);
246 const auto& origId = cloud::lookupIOField<label>("origId", obr);
247
248 // Create new parcels
249 for (label i = np; i < newNp; ++i)
250 {
251 c.addParticle(new parcelType(c.pMesh(), position[i], -1));
252 }
253
254 label i = 0;
255 for (particle& p : c)
256 {
257 p.origProc_ = origProcId[i];
258 p.origId_ = origId[i];
259
260 if (i < np)
261 {
262 // Use relocate for old particles, not new ones
263 p.relocate(position[i]);
264 }
265
266 ++i;
267 }
268 }
269}
270
271
272template<class CloudType>
274{
275 const label np = c.size();
276
277 auto& origProc = cloud::createIOField<label>("origProc", np, obr);
278 auto& origId = cloud::createIOField<label>("origId", np, obr);
279 auto& position = cloud::createIOField<point>("position", np, obr);
280
281 label i = 0;
282 for (const particle& p : c)
283 {
284 origProc[i] = p.origProc_;
285 origId[i] = p.origId_;
286 position[i] = p.position();
287
288 ++i;
289 }
290}
291
292
293template<class TrackCloudType>
295(
296 const vector& displacement,
297 TrackCloudType& cloud,
298 trackingData& td
299)
300{
301 typename TrackCloudType::particleType& p =
302 static_cast<typename TrackCloudType::particleType&>(*this);
303 typename TrackCloudType::particleType::trackingData& ttd =
304 static_cast<typename TrackCloudType::particleType::trackingData&>(td);
305
306 if (!p.hitPatch(cloud, ttd))
307 {
308 const polyPatch& patch = mesh_.boundaryMesh()[p.patch()];
309
310 if (isA<wedgePolyPatch>(patch))
311 {
312 p.hitWedgePatch(cloud, ttd);
313 }
314 else if (isA<symmetryPlanePolyPatch>(patch))
315 {
316 p.hitSymmetryPlanePatch(cloud, ttd);
317 }
318 else if (isA<symmetryPolyPatch>(patch))
319 {
320 p.hitSymmetryPatch(cloud, ttd);
321 }
322 else if (isA<cyclicPolyPatch>(patch))
323 {
324 p.hitCyclicPatch(cloud, ttd);
325 }
326 else if (isA<cyclicACMIPolyPatch>(patch))
327 {
328 p.hitCyclicACMIPatch(cloud, ttd, displacement);
329 }
330 else if (isA<cyclicAMIPolyPatch>(patch))
331 {
332 p.hitCyclicAMIPatch(cloud, ttd, displacement);
333 }
334 else if (isA<processorPolyPatch>(patch))
335 {
336 p.hitProcessorPatch(cloud, ttd);
337 }
338 else if (isA<wallPolyPatch>(patch))
339 {
340 p.hitWallPatch(cloud, ttd);
341 }
342 else
343 {
344 td.keepParticle = false;
345 }
346 }
347}
348
349
350template<class TrackCloudType>
352(
353 const vector& displacement,
354 TrackCloudType& cloud,
355 trackingData& td
356)
357{
358 if (!onFace())
359 {
360 return;
361 }
362 else if (onInternalFace())
363 {
364 changeCell();
365 }
366 else if (onBoundaryFace())
367 {
368 // If on duplicate baffle (e.g. cyclicACMI) change to lowest numbered
369 // patch ('masterpatch')
370 changeToMasterPatch();
371
372 hitBoundaryFace(displacement, cloud, td);
373 }
374}
375
376
377template<class TrackCloudType>
379(
380 const vector& direction,
381 const scalar fraction,
382 TrackCloudType& cloud,
383 trackingData& td
384)
385{
386 trackToFace(direction, fraction);
387
388 hitFace(direction, cloud, td);
389}
390
391
392template<class TrackCloudType>
394{
395 return false;
396}
397
398
399template<class TrackCloudType>
401{
403 << "Hitting a wedge patch should not be possible."
404 << abort(FatalError);
405
406 hitSymmetryPatch(cloud, td);
407}
408
409
410template<class TrackCloudType>
412(
413 TrackCloudType& cloud,
414 trackingData& td
415)
416{
417 hitSymmetryPatch(cloud, td);
418}
419
420
421template<class TrackCloudType>
423{
424 const vector nf = normal();
425
426 transformProperties(I - 2.0*nf*nf);
427}
428
429
430template<class TrackCloudType>
432{
433 const cyclicPolyPatch& cpp =
434 static_cast<const cyclicPolyPatch&>(mesh_.boundaryMesh()[patch()]);
435 const cyclicPolyPatch& receiveCpp = cpp.neighbPatch();
436 const label receiveFacei = receiveCpp.whichFace(facei_);
437
438 // Set the topology
439 facei_ = tetFacei_ = cpp.transformGlobalFace(facei_);
440 celli_ = mesh_.faceOwner()[facei_];
441 // See note in correctAfterParallelTransfer for tetPti addressing ...
442 tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;
443
444 // Reflect to account for the change of triangle orientation in the new cell
445 reflect();
446
447 // Transform the properties
448 if (!receiveCpp.parallel())
449 {
450 const tensor& T =
451 (
452 receiveCpp.forwardT().size() == 1
453 ? receiveCpp.forwardT()[0]
454 : receiveCpp.forwardT()[receiveFacei]
455 );
456 transformProperties(T);
457 }
458 else if (receiveCpp.separated())
459 {
460 const vector& s =
461 (
462 (receiveCpp.separation().size() == 1)
463 ? receiveCpp.separation()[0]
464 : receiveCpp.separation()[receiveFacei]
465 );
466 transformProperties(-s);
467 }
468}
469
470
471template<class TrackCloudType>
473(
474 TrackCloudType&,
475 trackingData& td,
476 const vector& displacement
477)
478{
479 vector pos = position();
480
481 const cyclicAMIPolyPatch& cpp =
482 static_cast<const cyclicAMIPolyPatch&>(mesh_.boundaryMesh()[patch()]);
483 const cyclicAMIPolyPatch& receiveCpp = cpp.neighbPatch();
484 const label sendFacei = cpp.whichFace(facei_);
485 const label receiveFacei = cpp.pointFace(sendFacei, displacement, pos);
486
487 if (receiveFacei < 0)
488 {
489 // If the patch face of the particle is not known assume that the
490 // particle is lost and mark it to be deleted.
491 td.keepParticle = false;
493 << "Particle lost across " << cyclicAMIPolyPatch::typeName
494 << " patches " << cpp.name() << " and " << receiveCpp.name()
495 << " at position " << pos << endl;
496 }
497
498 // Set the topology
499 facei_ = tetFacei_ = receiveFacei + receiveCpp.start();
500
501 // Locate the particle on the receiving side
502 vector displacementT = displacement;
503 cpp.reverseTransformDirection(displacementT, sendFacei);
504
505 // NOTE: The ray used to find the hit location accross the AMI might not
506 // be consistent in the displacement direction. Therefore a particle can
507 // be looping accross AMI patches indefinitely. Advancing the particle
508 // trajectory inside the cell is a possible solution.
509 const vector dispDir = cpp.fraction()*displacementT;
510 stepFraction_ += cpp.fraction();
511 locate
512 (
513 pos + dispDir,
514 &displacementT,
515 mesh_.faceOwner()[facei_],
516 false,
517 "Particle crossed between " + cyclicAMIPolyPatch::typeName +
518 " patches " + cpp.name() + " and " + receiveCpp.name() +
519 " to a location outside of the mesh."
520 );
521
522 // The particle must remain associated with a face for the tracking to
523 // register as incomplete
524 facei_ = tetFacei_;
525
526 // Transform the properties
527 if (!receiveCpp.parallel())
528 {
529 const tensor& T =
530 (
531 receiveCpp.forwardT().size() == 1
532 ? receiveCpp.forwardT()[0]
533 : receiveCpp.forwardT()[receiveFacei]
534 );
535 transformProperties(T);
536 }
537 else if (receiveCpp.separated())
538 {
539 const vector& s =
540 (
541 (receiveCpp.separation().size() == 1)
542 ? receiveCpp.separation()[0]
543 : receiveCpp.separation()[receiveFacei]
544 );
545 transformProperties(-s);
546 }
547
548 //if (onBoundaryFace())
549 {
550// vector receiveNormal, receiveDisplacement;
551// patchData(receiveNormal, receiveDisplacement);
552//
553// if (((displacementT - fraction*receiveDisplacement)&receiveNormal) > 0)
554// {
555// td.keepParticle = false;
556// WarningInFunction
557// << "Particle transfer from " << cyclicAMIPolyPatch::typeName
558// << " patches " << cpp.name() << " to " << receiveCpp.name()
559// << " failed at position " << pos << " and with displacement "
560// << (displacementT - fraction*receiveDisplacement) << nl
561// << " The displacement points into both the source and "
562// << "receiving faces, so the tracking cannot proceed" << nl
563// << " The particle has been removed" << nl << endl;
564// return;
565// }
566 }
567}
568
569
570template<class TrackCloudType>
572(
573 TrackCloudType& cloud,
574 trackingData& td,
575 const vector& displacement
576)
577{
578 const cyclicACMIPolyPatch& cpp =
579 static_cast<const cyclicACMIPolyPatch&>(mesh_.boundaryMesh()[patch()]);
580
581 const label localFacei = cpp.whichFace(facei_);
582
583 // If the mask is within the patch tolerance at either end, then we can
584 // assume an interaction with the appropriate part of the ACMI pair.
585 const scalar mask = cpp.mask()[localFacei];
586 bool couple = mask >= 1 - cpp.tolerance();
587 bool nonOverlap = mask <= cpp.tolerance();
588
589 // If the mask is an intermediate value, then we search for a location on
590 // the other side of the AMI. If we can't find a location, then we assume
591 // that we have hit the non-overlap patch.
592 if (!couple && !nonOverlap)
593 {
594 vector pos = position();
595 couple = cpp.pointFace(localFacei, displacement, pos) >= 0;
596 nonOverlap = !couple;
597 }
598
599 if (couple)
600 {
601 hitCyclicAMIPatch(cloud, td, displacement);
602 }
603 else
604 {
605 // Move to the face associated with the non-overlap patch and redo the
606 // face interaction.
607 tetFacei_ = facei_ = cpp.nonOverlapPatch().start() + localFacei;
608 // Bypass all the changeToMasterPatch logic so as not to change
609 // the patch back to cyclicACMI
610 hitBoundaryFace(displacement, cloud, td);
611 }
612}
613
614
615template<class TrackCloudType>
617{}
618
619
620template<class TrackCloudType>
622{}
623
624
625// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Generic templated field type.
Definition: Field.H:82
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
Helper IO class to read and write particle coordinates (positions).
Definition: IOPosition.H:56
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: IOPosition.C:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
static const IOField< point > * findIOPosition(const objectRegistry &obr)
Locate the "position" IOField within object registry.
Definition: cloud.H:153
virtual bool separated() const
Are the planes separated.
virtual bool parallel() const
Are the cyclic planes parallel.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
const scalarField & mask() const
Mask field where 1 = overlap(coupled), 0 = no-overlap.
static scalar tolerance()
Overlap tolerance.
Cyclic patch for Arbitrary Mesh Interface (AMI)
label pointFace(const label facei, const vector &n, point &p) const
scalar fraction() const
Particle fraction increase between AMI pathces.
virtual void reverseTransformDirection(vector &d, const label facei) const
Transform a patch-based direction from this side to.
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
Cyclic plane patch.
label transformGlobalFace(const label facei) const
const cyclicPolyPatch & neighbPatch() const
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
Allows specification of different writing frequency of objects registered to the database.
Definition: writeObjects.H:142
Registry of regIOobjects.
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
Base particle class.
Definition: particle.H:79
static void writeProperty(Ostream &os, const word &name, const Type &value, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a.
void hitFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
void hitBoundaryFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Dispatch function for boundary face interaction. Calls one of.
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitCyclicAMIPatch(TrackCloudType &, trackingData &, const vector &)
Overridable function to handle the particle hitting a cyclicAMIPatch.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
static void readObjects(CloudType &c, const objectRegistry &obr)
Read particle fields as objects from the obr registry.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
static void writePropertyName(Ostream &os, const word &name, const word &delim)
Write the name representation to stream.
const word & name() const noexcept
The patch name.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:451
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
virtual bool write(const bool valid=true) const
Write using setting from DB.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
bool match(const std::string &text, bool literal=false) const
Smart match as literal or regex, stopping on the first match.
Definition: wordResI.H:91
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
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))
#define WarningInFunction
Report a warning using Foam::Warning.
dimensionedScalar pos(const dimensionedScalar &ds)
static const Identity< scalar > I
Definition: Identity.H:94
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.