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 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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"
35 #include "symmetryPlanePolyPatch.H"
36 #include "symmetryPolyPatch.H"
37 #include "wallPolyPatch.H"
38 #include "wedgePolyPatch.H"
39 #include "meshTools.H"
40 
41 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
42 
43 template<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 
69 template<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 
97 template<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 
139 template<class TrackCloudType>
140 void 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 
169 template<class TrackCloudType>
170 void Foam::particle::writeFields(const TrackCloudType& c)
171 {
172  const label np = c.size();
173  const bool valid = np;
174 
175  if (writeLagrangianCoordinates)
176  {
177  IOPosition<TrackCloudType> ioP(c);
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  {
190  IOPosition<TrackCloudType> ioP
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 
223 template<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 
272 template<class CloudType>
273 void Foam::particle::writeObjects(const CloudType& c, objectRegistry& obr)
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 
293 template<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 (!onFace())
307  {
308  return;
309  }
310  else if (onInternalFace())
311  {
312  changeCell();
313  }
314  else if (onBoundaryFace())
315  {
316  changeToMasterPatch();
317 
318  if (!p.hitPatch(cloud, ttd))
319  {
320  const polyPatch& patch = mesh_.boundaryMesh()[p.patch()];
321 
322  if (isA<wedgePolyPatch>(patch))
323  {
324  p.hitWedgePatch(cloud, ttd);
325  }
326  else if (isA<symmetryPlanePolyPatch>(patch))
327  {
328  p.hitSymmetryPlanePatch(cloud, ttd);
329  }
330  else if (isA<symmetryPolyPatch>(patch))
331  {
332  p.hitSymmetryPatch(cloud, ttd);
333  }
334  else if (isA<cyclicPolyPatch>(patch))
335  {
336  p.hitCyclicPatch(cloud, ttd);
337  }
338  else if (isA<cyclicACMIPolyPatch>(patch))
339  {
340  p.hitCyclicACMIPatch(cloud, ttd, displacement);
341  }
342  else if (isA<cyclicAMIPolyPatch>(patch))
343  {
344  p.hitCyclicAMIPatch(cloud, ttd, displacement);
345  }
346  else if (isA<processorPolyPatch>(patch))
347  {
348  p.hitProcessorPatch(cloud, ttd);
349  }
350  else if (isA<wallPolyPatch>(patch))
351  {
352  p.hitWallPatch(cloud, ttd);
353  }
354  else
355  {
356  td.keepParticle = false;
357  }
358  }
359  }
360 }
361 
362 
363 template<class TrackCloudType>
365 (
366  const vector& direction,
367  const scalar fraction,
368  TrackCloudType& cloud,
369  trackingData& td
370 )
371 {
372  trackToFace(direction, fraction);
373 
374  hitFace(direction, cloud, td);
375 }
376 
377 
378 template<class TrackCloudType>
379 bool Foam::particle::hitPatch(TrackCloudType&, trackingData&)
380 {
381  return false;
382 }
383 
384 
385 template<class TrackCloudType>
386 void Foam::particle::hitWedgePatch(TrackCloudType& cloud, trackingData& td)
387 {
389  << "Hitting a wedge patch should not be possible."
390  << abort(FatalError);
391 
392  hitSymmetryPatch(cloud, td);
393 }
394 
395 
396 template<class TrackCloudType>
398 (
399  TrackCloudType& cloud,
400  trackingData& td
401 )
402 {
403  hitSymmetryPatch(cloud, td);
404 }
405 
406 
407 template<class TrackCloudType>
408 void Foam::particle::hitSymmetryPatch(TrackCloudType&, trackingData&)
409 {
410  const vector nf = normal();
411 
412  transformProperties(I - 2.0*nf*nf);
413 }
414 
415 
416 template<class TrackCloudType>
417 void Foam::particle::hitCyclicPatch(TrackCloudType&, trackingData&)
418 {
419  const cyclicPolyPatch& cpp =
420  static_cast<const cyclicPolyPatch&>(mesh_.boundaryMesh()[patch()]);
421  const cyclicPolyPatch& receiveCpp = cpp.neighbPatch();
422  const label receiveFacei = receiveCpp.whichFace(facei_);
423 
424  // Set the topology
425  facei_ = tetFacei_ = cpp.transformGlobalFace(facei_);
426  celli_ = mesh_.faceOwner()[facei_];
427  // See note in correctAfterParallelTransfer for tetPti addressing ...
428  tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;
429 
430  // Reflect to account for the change of triangle orientation in the new cell
431  reflect();
432 
433  // Transform the properties
434  if (!receiveCpp.parallel())
435  {
436  const tensor& T =
437  (
438  receiveCpp.forwardT().size() == 1
439  ? receiveCpp.forwardT()[0]
440  : receiveCpp.forwardT()[receiveFacei]
441  );
442  transformProperties(T);
443  }
444  else if (receiveCpp.separated())
445  {
446  const vector& s =
447  (
448  (receiveCpp.separation().size() == 1)
449  ? receiveCpp.separation()[0]
450  : receiveCpp.separation()[receiveFacei]
451  );
452  transformProperties(-s);
453  }
454 }
455 
456 
457 template<class TrackCloudType>
459 (
460  TrackCloudType&,
461  trackingData& td,
462  const vector& displacement
463 )
464 {
465  vector pos = position();
466 
467  const cyclicAMIPolyPatch& cpp =
468  static_cast<const cyclicAMIPolyPatch&>(mesh_.boundaryMesh()[patch()]);
469  const cyclicAMIPolyPatch& receiveCpp = cpp.neighbPatch();
470  const label sendFacei = cpp.whichFace(facei_);
471  const label receiveFacei = cpp.pointFace(sendFacei, displacement, pos);
472 
473  if (receiveFacei < 0)
474  {
475  // If the patch face of the particle is not known assume that the
476  // particle is lost and mark it to be deleted.
477  td.keepParticle = false;
479  << "Particle lost across " << cyclicAMIPolyPatch::typeName
480  << " patches " << cpp.name() << " and " << receiveCpp.name()
481  << " at position " << pos << endl;
482  }
483 
484  // Set the topology
485  facei_ = tetFacei_ = receiveFacei + receiveCpp.start();
486 
487  // Locate the particle on the receiving side
488  vector displacementT = displacement;
489  cpp.reverseTransformDirection(displacementT, sendFacei);
490 
491  // NOTE: The ray used to find the hit location accross the AMI might not
492  // be consistent in the displacement direction. Therefore a particle can
493  // be looping accross AMI patches indefinitely. Advancing the particle
494  // trajectory inside the cell is a possible solution.
495  const vector dispDir = cpp.fraction()*displacementT;
496  stepFraction_ += cpp.fraction();
497  locate
498  (
499  pos + dispDir,
500  &displacementT,
501  mesh_.faceOwner()[facei_],
502  false,
503  "Particle crossed between " + cyclicAMIPolyPatch::typeName +
504  " patches " + cpp.name() + " and " + receiveCpp.name() +
505  " to a location outside of the mesh."
506  );
507 
508  // The particle must remain associated with a face for the tracking to
509  // register as incomplete
510  facei_ = tetFacei_;
511 
512  // Transform the properties
513  if (!receiveCpp.parallel())
514  {
515  const tensor& T =
516  (
517  receiveCpp.forwardT().size() == 1
518  ? receiveCpp.forwardT()[0]
519  : receiveCpp.forwardT()[receiveFacei]
520  );
521  transformProperties(T);
522  }
523  else if (receiveCpp.separated())
524  {
525  const vector& s =
526  (
527  (receiveCpp.separation().size() == 1)
528  ? receiveCpp.separation()[0]
529  : receiveCpp.separation()[receiveFacei]
530  );
531  transformProperties(-s);
532  }
533 
534  //if (onBoundaryFace())
535  {
536 // vector receiveNormal, receiveDisplacement;
537 // patchData(receiveNormal, receiveDisplacement);
538 //
539 // if (((displacementT - fraction*receiveDisplacement)&receiveNormal) > 0)
540 // {
541 // td.keepParticle = false;
542 // WarningInFunction
543 // << "Particle transfer from " << cyclicAMIPolyPatch::typeName
544 // << " patches " << cpp.name() << " to " << receiveCpp.name()
545 // << " failed at position " << pos << " and with displacement "
546 // << (displacementT - fraction*receiveDisplacement) << nl
547 // << " The displacement points into both the source and "
548 // << "receiving faces, so the tracking cannot proceed" << nl
549 // << " The particle has been removed" << nl << endl;
550 // return;
551 // }
552  }
553 }
554 
555 
556 template<class TrackCloudType>
558 (
559  TrackCloudType& cloud,
560  trackingData& td,
561  const vector& displacement
562 )
563 {
564  const cyclicACMIPolyPatch& cpp =
565  static_cast<const cyclicACMIPolyPatch&>(mesh_.boundaryMesh()[patch()]);
566 
567  const label localFacei = cpp.whichFace(facei_);
568 
569  // If the mask is within the patch tolerance at either end, then we can
570  // assume an interaction with the appropriate part of the ACMI pair.
571  const scalar mask = cpp.mask()[localFacei];
572  bool couple = mask >= 1 - cpp.tolerance();
573  bool nonOverlap = mask <= cpp.tolerance();
574 
575  // If the mask is an intermediate value, then we search for a location on
576  // the other side of the AMI. If we can't find a location, then we assume
577  // that we have hit the non-overlap patch.
578  if (!couple && !nonOverlap)
579  {
580  vector pos = position();
581  couple = cpp.pointFace(localFacei, displacement, pos) >= 0;
582  nonOverlap = !couple;
583  }
584 
585  if (couple)
586  {
587  hitCyclicAMIPatch(cloud, td, displacement);
588  }
589  else
590  {
591  // Move to the face associated with the non-overlap patch and redo the
592  // face interaction.
593  tetFacei_ = facei_ = cpp.nonOverlapPatch().start() + localFacei;
594  hitFace(displacement, cloud, td);
595  }
596 }
597 
598 
599 template<class TrackCloudType>
600 void Foam::particle::hitProcessorPatch(TrackCloudType&, trackingData&)
601 {}
602 
603 
604 template<class TrackCloudType>
605 void Foam::particle::hitWallPatch(TrackCloudType&, trackingData&)
606 {}
607 
608 
609 // ************************************************************************* //
Foam::particle::hitProcessorPatch
void hitProcessorPatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:593
Foam::particle::hitCyclicAMIPatch
void hitCyclicAMIPatch(TrackCloudType &, trackingData &, const vector &)
Definition: particleTemplates.C:452
Foam::Tensor
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:268
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
meshTools.H
IOPosition.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::coupledPolyPatch::separation
virtual const vectorField & separation() const
Definition: coupledPolyPatch.H:285
cyclicACMIPolyPatch.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::particle::hitCyclicPatch
void hitCyclicPatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:410
Foam::particle::origId
label origId() const
Definition: particleI.H:214
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:57
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
cyclicPolyPatch.H
Foam::cyclicPolyPatch::transformGlobalFace
label transformGlobalFace(const label facei) const
Definition: cyclicPolyPatch.H:403
Foam::cyclicPolyPatch
Cyclic plane patch.
Definition: cyclicPolyPatch.H:62
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
Definition: HashOps.H:164
Foam::IOobject::typeHeaderOk
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Definition: IOobjectTemplates.C:32
Foam::cyclicPolyPatch::neighbPatch
const cyclicPolyPatch & neighbPatch() const
Definition: cyclicPolyPatch.H:378
wedgePolyPatch.H
Foam::coupledPolyPatch::forwardT
virtual const tensorField & forwardT() const
Definition: coupledPolyPatch.H:297
Foam::particle::hitWedgePatch
void hitWedgePatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:379
Foam::cyclicACMIPolyPatch::tolerance
static scalar tolerance()
Definition: cyclicACMIPolyPatchI.H:58
wallPolyPatch.H
Foam::particle::readObjects
static void readObjects(CloudType &c, const objectRegistry &obr)
Definition: particleTemplates.C:217
Foam::particle::readFields
static void readFields(TrackCloudType &c)
Definition: particleTemplates.C:133
symmetryPolyPatch.H
Foam::dimensioned::name
const word & name() const
Definition: dimensionedType.C:399
Foam::particle::writeFields
static void writeFields(const TrackCloudType &c)
Definition: particleTemplates.C:163
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::wordRes::match
bool match(const std::string &text, bool literal=false) const
Definition: wordResI.H:84
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:56
Foam::particle::writeObjects
static void writeObjects(const CloudType &c, objectRegistry &obr)
Definition: particleTemplates.C:266
Foam::coupledPolyPatch::parallel
virtual bool parallel() const
Definition: coupledPolyPatch.H:291
Foam::cyclicACMIPolyPatch::nonOverlapPatch
const polyPatch & nonOverlapPatch() const
Definition: cyclicACMIPolyPatchI.H:29
Foam::particle::hitPatch
bool hitPatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:372
Foam::coupledPolyPatch::separated
virtual bool separated() const
Definition: coupledPolyPatch.H:279
Foam::particle::trackToAndHitFace
void trackToAndHitFace(const vector &direction, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Definition: particleTemplates.C:358
Foam::Field< Type >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:64
Foam::cyclicAMIPolyPatch::pointFace
label pointFace(const label facei, const vector &n, point &p) const
Definition: cyclicAMIPolyPatch.C:1172
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
Foam::particle::hitCyclicACMIPatch
void hitCyclicACMIPatch(TrackCloudType &, trackingData &, const vector &)
Definition: particleTemplates.C:551
Foam::particle::trackingData::keepParticle
bool keepParticle
Definition: particle.H:99
Foam::CloudType
DSMCCloud< dsmcParcel > CloudType
Definition: makeDSMCParcelBinaryCollisionModels.C:31
Foam::IOPosition
Helper IO class to read and write particle coordinates (positions).
Definition: Cloud.H:52
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Foam::IOPosition::write
virtual bool write(const bool valid=true) const
Definition: IOPosition.C:50
Foam::FatalError
error FatalError
processorPolyPatch.H
Foam::cyclicAMIPolyPatch::neighbPatch
virtual const cyclicAMIPolyPatch & neighbPatch() const
Definition: cyclicAMIPolyPatch.C:963
os
OBJstream os(runTime.globalPath()/outputName)
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:47
Foam::polyPatch::start
label start() const
Definition: polyPatch.H:357
cyclicAMIPolyPatch.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
symmetryPlanePolyPatch.H
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
Foam::particle::writeProperty
static void writeProperty(Ostream &os, const word &name, const Type &value, const bool nameOnly, const word &delim, const wordRes &filters=wordRes::null())
Definition: particleTemplates.C:64
FatalErrorInFunction
#define FatalErrorInFunction
Definition: error.H:465
Foam::cyclicAMIPolyPatch::fraction
scalar fraction() const
Definition: cyclicAMIPolyPatchI.H:55
Foam::polyPatch::whichFace
label whichFace(const label l) const
Definition: polyPatch.H:444
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::foamVersion::patch
const std::string patch
Foam::Vector< scalar >
Foam::pTraits
A traits class, which is primarily used for primitives.
Definition: pTraits.H:50
Foam::cyclicAMIPolyPatch::reverseTransformDirection
virtual void reverseTransformDirection(vector &d, const label facei) const
Definition: cyclicAMIPolyPatch.C:1113
Foam::particle::hitSymmetryPatch
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:401
Foam::particle
Base particle class.
Definition: particle.H:70
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
Foam::direction
uint8_t direction
Definition: direction.H:46
Foam::particle::hitWallPatch
void hitWallPatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:598
Foam::particle::hitSymmetryPlanePatch
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Definition: particleTemplates.C:391
Foam::cloud::findIOPosition
static const IOField< point > * findIOPosition(const objectRegistry &obr)
Definition: cloud.H:149
Foam::constant::universal::c
const dimensionedScalar c
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::DSMCCloud::parcelType
ParcelType parcelType
Definition: DSMCCloud.H:216
Foam::particle::writePropertyName
static void writePropertyName(Ostream &os, const word &name, const word &delim)
Definition: particleTemplates.C:38
Foam::cloud::geometryType::POSITIONS
@ POSITIONS
positions
Foam::patchIdentifier::name
const word & name() const noexcept
Definition: patchIdentifier.H:131
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:52
Foam::cyclicACMIPolyPatch
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI).
Definition: cyclicACMIPolyPatch.H:74
Foam::particle::hitFace
void hitFace(const vector &direction, TrackCloudType &cloud, trackingData &td)
Definition: particleTemplates.C:288
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:184
Foam::cyclicACMIPolyPatch::mask
const scalarField & mask() const
Definition: cyclicACMIPolyPatchI.H:47
Foam::particle::trackingData
Definition: particle.H:89
WarningInFunction
#define WarningInFunction
Definition: messageStream.H:353
Foam::tensor
Tensor< scalar > tensor
Tensor of scalars, i.e. Tensor<scalar>.
Definition: symmTensor.H:57
Foam::I
static const Identity< scalar > I
Definition: Identity.H:89
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:170
Foam::cyclicAMIPolyPatch
Cyclic patch for Arbitrary Mesh Interface (AMI)
Definition: cyclicAMIPolyPatch.H:64