particleTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "IOPosition.H"
27 
28 #include "cyclicPolyPatch.H"
29 #include "cyclicAMIPolyPatch.H"
30 #include "cyclicACMIPolyPatch.H"
32 #include "processorPolyPatch.H"
33 #include "symmetryPlanePolyPatch.H"
34 #include "symmetryPolyPatch.H"
35 #include "wallPolyPatch.H"
36 #include "wedgePolyPatch.H"
37 #include "meshTools.H"
38 
39 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40 
41 template<class TrackCloudType>
42 void Foam::particle::readFields(TrackCloudType& c)
43 {
44  bool valid = c.size();
45 
46  IOobject procIO(c.fieldIOobject("origProcId", IOobject::MUST_READ));
47 
48  bool haveFile = procIO.typeHeaderOk<IOField<label>>(true);
49 
50  IOField<label> origProcId(procIO, valid && haveFile);
51  c.checkFieldIOobject(c, origProcId);
53  (
54  c.fieldIOobject("origId", IOobject::MUST_READ),
55  valid && haveFile
56  );
57  c.checkFieldIOobject(c, origId);
58 
59  label i = 0;
60  forAllIter(typename TrackCloudType, c, iter)
61  {
62  particle& p = iter();
63 
64  p.origProc_ = origProcId[i];
65  p.origId_ = origId[i];
66  i++;
67  }
68 }
69 
70 
71 template<class TrackCloudType>
72 void Foam::particle::writeFields(const TrackCloudType& c)
73 {
74  label np = c.size();
75 
77  ioP.write(np > 0);
78 
80  (
81  c.fieldIOobject("origProcId", IOobject::NO_READ),
82  np
83  );
85  (
86  c.fieldIOobject("origId", IOobject::NO_READ),
87  np
88  );
89 
90  label i = 0;
91  forAllConstIter(typename TrackCloudType, c, iter)
92  {
93  origProc[i] = iter().origProc_;
94  origId[i] = iter().origId_;
95  i++;
96  }
97 
98  origProc.write(np > 0);
99  origId.write(np > 0);
100 }
101 
102 
103 template<class TrackCloudType>
105 (
106  const vector& displacement,
107  const scalar fraction,
108  TrackCloudType& cloud,
109  trackingData& td
110 )
111 {
112  if (debug)
113  {
114  Info << "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
115  }
116 
117  if (onBoundaryFace())
118  {
119  changeToMasterPatch();
120  }
121 
122  hitFaceNoChangeToMasterPatch(displacement, fraction, cloud, td);
123 }
124 
125 
126 template<class TrackCloudType>
128 (
129  const vector& displacement,
130  const scalar fraction,
131  TrackCloudType& cloud,
132  trackingData& td
133 )
134 {
135  if (debug)
136  {
137  Info << "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
138  }
139 
140  typename TrackCloudType::particleType& p =
141  static_cast<typename TrackCloudType::particleType&>(*this);
142  typename TrackCloudType::particleType::trackingData& ttd =
143  static_cast<typename TrackCloudType::particleType::trackingData&>(td);
144 
145  if (!onFace())
146  {
147  return;
148  }
149  else if (onInternalFace())
150  {
151  changeCell();
152  }
153  else if (onBoundaryFace())
154  {
155  if (!p.hitPatch(cloud, ttd))
156  {
157  const polyPatch& patch = mesh_.boundaryMesh()[p.patch()];
158 
159  if (isA<wedgePolyPatch>(patch))
160  {
161  p.hitWedgePatch(cloud, ttd);
162  }
163  else if (isA<symmetryPlanePolyPatch>(patch))
164  {
165  p.hitSymmetryPlanePatch(cloud, ttd);
166  }
167  else if (isA<symmetryPolyPatch>(patch))
168  {
169  p.hitSymmetryPatch(cloud, ttd);
170  }
171  else if (isA<cyclicPolyPatch>(patch))
172  {
173  p.hitCyclicPatch(cloud, ttd);
174  }
175  else if (isA<cyclicACMIPolyPatch>(patch))
176  {
177  p.hitCyclicACMIPatch(displacement, fraction, cloud, ttd);
178  }
179  else if (isA<cyclicAMIPolyPatch>(patch))
180  {
181  p.hitCyclicAMIPatch(displacement, fraction, cloud, ttd);
182  }
183  else if (isA<cyclicRepeatAMIPolyPatch>(patch))
184  {
185  p.hitCyclicRepeatAMIPatch(displacement, fraction, cloud, ttd);
186  }
187  else if (isA<processorPolyPatch>(patch))
188  {
189  p.hitProcessorPatch(cloud, ttd);
190  }
191  else if (isA<wallPolyPatch>(patch))
192  {
193  p.hitWallPatch(cloud, ttd);
194  }
195  else
196  {
197  td.keepParticle = false;
198  }
199  }
200  }
201 }
202 
203 
204 template<class TrackCloudType>
206 (
207  const vector& displacement,
208  const scalar fraction,
209  TrackCloudType& cloud,
210  trackingData& td
211 )
212 {
213  if (debug)
214  {
215  Info << "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
216  }
217 
218  const scalar f = trackToFace(displacement, fraction);
219 
220  hitFace(displacement, fraction, cloud, td);
221 
222  return f;
223 }
224 
225 
226 template<class TrackCloudType>
227 bool Foam::particle::hitPatch(TrackCloudType&, trackingData&)
228 {
229  return false;
230 }
231 
232 
233 template<class TrackCloudType>
234 void Foam::particle::hitWedgePatch(TrackCloudType& cloud, trackingData& td)
235 {
237  << "Hitting a wedge patch should not be possible."
238  << abort(FatalError);
239 
240  hitSymmetryPatch(cloud, td);
241 }
242 
243 
244 template<class TrackCloudType>
246 (
247  TrackCloudType& cloud,
248  trackingData& td
249 )
250 {
251  hitSymmetryPatch(cloud, td);
252 }
253 
254 
255 template<class TrackCloudType>
257 {
258  const vector nf = normal();
259  transformProperties(I - 2.0*nf*nf);
260 }
261 
262 
263 template<class TrackCloudType>
265 {
266  const cyclicPolyPatch& cpp =
267  static_cast<const cyclicPolyPatch&>(mesh_.boundaryMesh()[patch()]);
268  const cyclicPolyPatch& receiveCpp = cpp.neighbPatch();
269  const label receiveFacei = receiveCpp.whichFace(facei_);
270 
271  // Set the topology
272  facei_ = tetFacei_ = cpp.transformGlobalFace(facei_);
273  celli_ = mesh_.faceOwner()[facei_];
274  // See note in correctAfterParallelTransfer for tetPti addressing ...
275  tetPti_ = mesh_.faces()[tetFacei_].size() - 1 - tetPti_;
276 
277  // Reflect to account for the change of triangle orientation in the new cell
278  reflect();
279 
280  // Transform the properties
281  if (!receiveCpp.parallel())
282  {
283  const tensor& T =
284  (
285  receiveCpp.forwardT().size() == 1
286  ? receiveCpp.forwardT()[0]
287  : receiveCpp.forwardT()[receiveFacei]
288  );
290  }
291  else if (receiveCpp.separated())
292  {
293  const vector& s =
294  (
295  (receiveCpp.separation().size() == 1)
296  ? receiveCpp.separation()[0]
297  : receiveCpp.separation()[receiveFacei]
298  );
300  }
301 }
302 
303 
304 template<class TrackCloudType>
306 (
307  const vector& displacement,
308  const scalar fraction,
309  TrackCloudType& cloud,
310  trackingData& td
311 )
312 {
313  const cyclicAMIPolyPatch& cpp =
314  static_cast<const cyclicAMIPolyPatch&>(mesh_.boundaryMesh()[patch()]);
315  const cyclicAMIPolyPatch& receiveCpp = cpp.neighbPatch();
316 
317  if (debug)
318  {
319  Info<< "Particle " << origId() << " crossing AMI from " << cpp.name()
320  << " to " << receiveCpp.name() << endl << endl;
321  }
322 
323  // Get the send patch data
324  vector sendNormal, sendDisplacement;
325  patchData(sendNormal, sendDisplacement);
326 
327  vector pos = position();
328 
329  const labelPair receiveIs =
330  cpp.pointAMIAndFace
331  (
332  cpp.whichFace(facei_),
333  displacement - fraction*sendDisplacement,
334  pos
335  );
336  const label receiveAMIi = receiveIs.first();
337  const label receiveFacei = receiveIs.second();
338 
339  // If the receiving face could not be found then issue a warning and remove
340  // the particle
341  if (receiveFacei < 0)
342  {
343  td.keepParticle = false;
345  << "Particle transfer from " << cyclicAMIPolyPatch::typeName
346  << " patches " << cpp.name() << " to " << receiveCpp.name()
347  << " failed at position " << pos << " and with displacement "
348  << (displacement - fraction*sendDisplacement) << nl
349  << " A receiving face could not be found" << nl
350  << " The particle has been removed" << nl << endl;
351  return;
352  }
353 
354  // Set the topology
355  facei_ = tetFacei_ = receiveFacei + receiveCpp.start();
356 
357  // Locate the particle on the receiving side
358  locate
359  (
360  pos,
361  mesh_.faceOwner()[facei_],
362  false,
363  "Particle crossed between " + cyclicAMIPolyPatch::typeName +
364  " patches " + cpp.name() + " and " + receiveCpp.name() +
365  " to a location outside of the mesh."
366  );
367 
368  // The particle must remain associated with a face for the tracking to
369  // register as incomplete
370  facei_ = tetFacei_;
371 
372  // Transform the properties
373  vector displacementT = displacement;
374 
375  const vectorTensorTransform AMITransform =
376  receiveCpp.owner()
377  ? receiveCpp.AMITransforms()[receiveAMIi]
378  : inv(cpp.AMITransforms()[receiveAMIi]);
379  if (AMITransform.hasR())
380  {
381  transformProperties(AMITransform.R());
382  displacementT = transform(AMITransform.R(), displacementT);
383  }
384  else if (AMITransform.t() != vector::zero)
385  {
386  transformProperties(AMITransform.t());
387  }
388 
389  if (!receiveCpp.parallel())
390  {
391  const tensor& T =
392  (
393  receiveCpp.forwardT().size() == 1
394  ? receiveCpp.forwardT()[0]
395  : receiveCpp.forwardT()[facei_]
396  );
398  displacementT = transform(T, displacementT);
399  }
400  else if (receiveCpp.separated())
401  {
402  const vector& s =
403  (
404  (receiveCpp.separation().size() == 1)
405  ? receiveCpp.separation()[0]
406  : receiveCpp.separation()[facei_]
407  );
409  }
410 
411  // If on a boundary and the displacement points into the receiving face
412  // then issue a warning and remove the particle
413  if (onBoundaryFace())
414  {
415  vector receiveNormal, receiveDisplacement;
416  patchData(receiveNormal, receiveDisplacement);
417 
418  if (((displacementT - fraction*receiveDisplacement)&receiveNormal) > 0)
419  {
420  td.keepParticle = false;
422  << "Particle transfer from " << cyclicAMIPolyPatch::typeName
423  << " patches " << cpp.name() << " to " << receiveCpp.name()
424  << " failed at position " << pos << " and with displacement "
425  << (displacementT - fraction*receiveDisplacement) << nl
426  << " The displacement points into both the source and "
427  << "receiving faces, so the tracking cannot proceed" << nl
428  << " The particle has been removed" << nl << endl;
429  return;
430  }
431  }
432 }
433 
434 
435 template<class TrackCloudType>
437 (
438  const vector& displacement,
439  const scalar fraction,
440  TrackCloudType& cloud,
441  trackingData& td
442 )
443 {
444  typename TrackCloudType::particleType& p =
445  static_cast<typename TrackCloudType::particleType&>(*this);
446 
447  const cyclicACMIPolyPatch& cpp =
448  static_cast<const cyclicACMIPolyPatch&>(mesh_.boundaryMesh()[patch()]);
449 
450  vector patchNormal, patchDisplacement;
451  patchData(patchNormal, patchDisplacement);
452 
453  const label localFacei = cpp.whichFace(facei_);
454 
455  // If the mask is within the patch tolerance at either end, then we can
456  // assume an interaction with the appropriate part of the ACMI pair.
457  const scalar mask = cpp.mask()[localFacei];
458  bool couple = mask >= 1 - cpp.tolerance();
459  bool nonOverlap = mask <= cpp.tolerance();
460 
461  // If the mask is an intermediate value, then we search for a location on
462  // the other side of the AMI. If we can't find a location, then we assume
463  // that we have hit the non-overlap patch.
464  if (!couple && !nonOverlap)
465  {
466  vector pos = position();
467  couple =
468  cpp.pointAMIAndFace
469  (
470  localFacei,
471  displacement - fraction*patchDisplacement,
472  pos
473  ).first() >= 0;
474  nonOverlap = !couple;
475  }
476 
477  if (couple)
478  {
479  p.hitCyclicAMIPatch(displacement, fraction, cloud, td);
480  }
481  else
482  {
483  // Move to the face associated with the non-overlap patch and redo the
484  // face interaction.
485  tetFacei_ = facei_ = cpp.nonOverlapPatch().start() + localFacei;
486  p.hitFaceNoChangeToMasterPatch(displacement, fraction, cloud, td);
487  }
488 }
489 
490 
491 template<class TrackCloudType>
493 (
494  const vector& displacement,
495  const scalar fraction,
496  TrackCloudType& cloud,
497  trackingData& td
498 )
499 {
500  typename TrackCloudType::particleType& p =
501  static_cast<typename TrackCloudType::particleType&>(*this);
502 
503  p.hitCyclicAMIPatch(displacement, fraction, cloud, td);
504 }
505 
506 
507 template<class TrackCloudType>
509 {}
510 
511 
512 template<class TrackCloudType>
514 {}
515 
516 
517 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
void hitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
label origProc() const
Return the originating processor ID.
Definition: particleI.H:173
Cyclic patch for Arbitrarily Coupled Mesh Interface (ACMI)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
virtual bool separated() const
Are the planes separated.
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:266
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Vector-tensor class used to perform translations and rotations in 3D space.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
scalar trackToFace(const vector &displacement, const scalar fraction)
As particle::track, but stops when a face is hit.
Definition: particle.C:644
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
labelPair pointAMIAndFace(const label facei, const vector &n, point &p) const
Return the transform and face indices on neighbour patch which.
scalar trackToAndHitFace(const vector &displacement, 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.
void hitSymmetryPlanePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
virtual void transformProperties(const tensor &T)
Transform the physical properties of the particle.
Definition: particle.C:1046
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a symmetryPatch.
virtual const vectorField & separation() const
If the planes are separated the separation vector.
virtual const tensorField & forwardT() const
Return face transformation tensor.
Base particle class.
Definition: particle.H:84
const polyPatch & nonOverlapPatch() const
Return a const reference to the non-overlapping patch.
const List< vectorTensorTransform > & AMITransforms() const
Return a reference to the AMI transforms.
vector normal() const
Return the normal of the tri on tetFacei_ for the.
Definition: particleI.H:248
void patchData(vector &normal, vector &displacement) const
Get the normal and displacement of the current patch location.
Definition: particleI.H:292
virtual bool parallel() const
Are the cyclic planes parallel.
dimensionedScalar pos(const dimensionedScalar &ds)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
void hitFaceNoChangeToMasterPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
As above, but does not change the master patch. Needed in order for.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
static const Identity< scalar > I
Definition: Identity.H:93
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
Cyclic plane patch.
virtual bool write(const bool write=true) const
Write using setting from DB.
Definition: IOPosition.C:51
void hitCyclicACMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a.
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a cyclicPatch.
Cyclic patch for Arbitrary Mesh Interface (AMI)
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:185
static const char nl
Definition: Ostream.H:265
bool onFace() const
Is the particle on a face?
Definition: particleI.H:254
#define FUNCTION_NAME
labelList f(nPoints)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const Type & second() const
Return second.
Definition: Pair.H:99
void hitWedgePatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wedgePatch.
const scalarField & mask() const
Mask field where 1 = overlap, 0 = no-overlap.
#define WarningInFunction
Report a warning using Foam::Warning.
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:113
virtual const cyclicAMIPolyPatch & neighbPatch() const
Return a reference to the neighbour patch.
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
virtual bool owner() const
Does this side own the patch?
const dimensionedScalar c
Speed of light in a vacuum.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
messageStream Info
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:260
void hitCyclicAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting a cyclicAMIPatch.
label transformGlobalFace(const label facei) const
volScalarField & p
Helper IO class to read and write particle positions.
Definition: Cloud.H:55
const cyclicPolyPatch & neighbPatch() const
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const Type & first() const
Return first.
Definition: Pair.H:87
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
vector position() const
Return current particle position.
Definition: particleI.H:278
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:272
void hitCyclicRepeatAMIPatch(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
static scalar tolerance()
Overlap tolerance.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:380