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-2022 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 "particle.H"
27 #include "IOPosition.H"
28 
29 #include "cyclicPolyPatch.H"
31 #include "processorPolyPatch.H"
32 #include "symmetryPlanePolyPatch.H"
33 #include "symmetryPolyPatch.H"
34 #include "wallPolyPatch.H"
35 #include "wedgePolyPatch.H"
36 #include "meshTools.H"
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 template<class TrackCloudType>
41 void Foam::particle::readFields(TrackCloudType& c)
42 {
43  bool valid = c.size();
44 
46  (
47  c.fieldIOobject("origProcId", IOobject::MUST_READ)
48  );
49 
50  bool haveFile = procIO.headerOk();
51 
52  IOField<label> origProcId(procIO, valid && haveFile);
53  c.checkFieldIOobject(c, origProcId);
55  (
56  c.fieldIOobject("origId", IOobject::MUST_READ),
57  valid && haveFile
58  );
59  c.checkFieldIOobject(c, origId);
60 
61  label i = 0;
62  forAllIter(typename TrackCloudType, c, iter)
63  {
64  particle& p = iter();
65 
66  p.origProc_ = origProcId[i];
67  p.origId_ = origId[i];
68  i++;
69  }
70 }
71 
72 
73 template<class TrackCloudType>
74 void Foam::particle::writeFields(const TrackCloudType& c)
75 {
76  label np = c.size();
77 
79  ioP.write(np > 0);
80 
81  IOField<label> origProc
82  (
83  c.fieldIOobject("origProcId", IOobject::NO_READ),
84  np
85  );
86  IOField<label> origId
87  (
88  c.fieldIOobject("origId", IOobject::NO_READ),
89  np
90  );
91 
92  label i = 0;
93  forAllConstIter(typename TrackCloudType, c, iter)
94  {
95  origProc[i] = iter().origProc_;
96  origId[i] = iter().origId_;
97  i++;
98  }
99 
100  origProc.write(np > 0);
101  origId.write(np > 0);
102 }
103 
104 
105 template<class TrackCloudType>
107 (
108  TrackCloudType& cloud,
109  trackingData& td
110 )
111 {
112  if (td.sendFromPatch == patch(td.mesh))
113  {
114  prepareForProcessorTransfer(td);
115  }
116  else
117  {
118  prepareForNonConformalCyclicTransfer
119  (
120  td.mesh,
121  td.sendFromPatch,
122  td.sendToPatchFace
123  );
124  }
125 }
126 
127 
128 template<class TrackCloudType>
130 (
131  TrackCloudType& cloud,
132  trackingData& td
133 )
134 {
135  const polyPatch& pp = td.mesh.boundaryMesh()[td.sendToPatch];
136 
137  if (isA<processorPolyPatch>(pp))
138  {
139  correctAfterProcessorTransfer(td);
140  }
141  else if (isA<nonConformalCyclicPolyPatch>(pp))
142  {
143  correctAfterNonConformalCyclicTransfer(td.mesh, td.sendToPatch);
144  }
145  else
146  {
148  << "Transfer patch type not recognised"
149  << exit(FatalError);
150  }
151 }
152 
153 
154 template<class TrackCloudType>
156 (
157  const vector& displacement,
158  const scalar fraction,
159  TrackCloudType& cloud,
160  trackingData& td
161 )
162 {
163  if (debug)
164  {
165  Info << "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
166  }
167 
168  typename TrackCloudType::particleType& p =
169  static_cast<typename TrackCloudType::particleType&>(*this);
170  typename TrackCloudType::particleType::trackingData& ttd =
171  static_cast<typename TrackCloudType::particleType::trackingData&>(td);
172 
173  if (!onFace())
174  {
175  return;
176  }
177  else if (onInternalFace(td.mesh))
178  {
179  changeCell(td.mesh);
180  }
181  else if (onBoundaryFace(td.mesh))
182  {
183  forAll(cloud.patchNonConformalCyclicPatches()[p.patch(td.mesh)], i)
184  {
185  if
186  (
187  p.hitNonConformalCyclicPatch
188  (
189  displacement,
190  fraction,
191  cloud.patchNonConformalCyclicPatches()[p.patch(td.mesh)][i],
192  cloud,
193  ttd
194  )
195  )
196  {
197  return;
198  }
199  }
200 
201  if (!p.hitPatch(cloud, ttd))
202  {
203  const polyPatch& patch = td.mesh.boundaryMesh()[p.patch(td.mesh)];
204 
205  if (isA<wedgePolyPatch>(patch))
206  {
207  p.hitWedgePatch(cloud, ttd);
208  }
209  else if (isA<symmetryPlanePolyPatch>(patch))
210  {
211  p.hitSymmetryPlanePatch(cloud, ttd);
212  }
213  else if (isA<symmetryPolyPatch>(patch))
214  {
215  p.hitSymmetryPatch(cloud, ttd);
216  }
217  else if (isA<cyclicPolyPatch>(patch))
218  {
219  p.hitCyclicPatch(cloud, ttd);
220  }
221  else if (isA<processorPolyPatch>(patch))
222  {
223  p.hitProcessorPatch(cloud, ttd);
224  }
225  else if (isA<wallPolyPatch>(patch))
226  {
227  p.hitWallPatch(cloud, ttd);
228  }
229  else
230  {
231  p.hitBasicPatch(cloud, ttd);
232  }
233  }
234  }
235 }
236 
237 
238 template<class TrackCloudType>
240 (
241  const vector& displacement,
242  const scalar fraction,
243  TrackCloudType& cloud,
244  trackingData& td
245 )
246 {
247  if (debug)
248  {
249  Info << "Particle " << origId() << nl << FUNCTION_NAME << nl << endl;
250  }
251 
252  const scalar f = trackToFace(td.mesh, displacement, fraction);
253 
254  hitFace(displacement, fraction, cloud, td);
255 
256  return f;
257 }
258 
259 
260 template<class TrackCloudType>
261 bool Foam::particle::hitPatch(TrackCloudType&, trackingData&)
262 {
263  return false;
264 }
265 
266 
267 template<class TrackCloudType>
269 {
271  << "Hitting a wedge patch should not be possible."
272  << abort(FatalError);
273 
274  hitSymmetryPatch(cloud, td);
275 }
276 
277 
278 template<class TrackCloudType>
280 (
281  TrackCloudType& cloud,
282  trackingData& td
283 )
284 {
285  hitSymmetryPatch(cloud, td);
286 }
287 
288 
289 template<class TrackCloudType>
291 {
292  const vector nf = normal(td.mesh);
293  transformProperties(transformer::rotation(I - 2.0*nf*nf));
294 }
295 
296 
297 template<class TrackCloudType>
299 {
300  const cyclicPolyPatch& cpp =
301  static_cast<const cyclicPolyPatch&>
302  (
303  td.mesh.boundaryMesh()[patch(td.mesh)]
304  );
305  const cyclicPolyPatch& receiveCpp = cpp.nbrPatch();
306 
307  // Set the topology
308  facei_ = tetFacei_ = cpp.transformGlobalFace(facei_);
309  celli_ = td.mesh.faceOwner()[facei_];
310 
311  // See note in correctAfterParallelTransfer for tetPti addressing ...
312  tetPti_ = td.mesh.faces()[tetFacei_].size() - 1 - tetPti_;
313 
314  // Reflect to account for the change of triangle orientation in the new cell
315  reflect();
316 
317  // Transform the properties
318  if (receiveCpp.transform().transformsPosition())
319  {
320  transformProperties(receiveCpp.transform());
321  }
322 }
323 
324 
325 template<class TrackCloudType>
327 (
328  const vector& displacement,
329  const scalar fraction,
330  const label patchi,
331  TrackCloudType& cloud,
332  trackingData& td
333 )
334 {
335  const nonConformalCyclicPolyPatch& nccpp =
336  static_cast<const nonConformalCyclicPolyPatch&>
337  (td.mesh.boundaryMesh()[patchi]);
338 
339  const point sendPos = position(td.mesh);
340 
341  // Get the send patch data
342  vector sendNormal, sendDisplacement;
343  patchData(td.mesh, sendNormal, sendDisplacement);
344 
345  // Project the particle through the non-conformal patch
346  point receivePos;
347  const remote receiveProcFace =
348  nccpp.ray
349  (
350  stepFraction_,
351  nccpp.origPatch().whichFace(facei_),
352  sendPos,
353  displacement - fraction*sendDisplacement,
354  receivePos
355  );
356 
357  // If we didn't hit anything then this particle is assumed to project to
358  // the orig patch, or another different non-conformal patch. Return, so
359  // these can be tried.
360  if (receiveProcFace.proci == -1) return false;
361 
362  // If we are transferring between processes then mark as such and return.
363  // The cloud will handle all processor transfers as a single batch.
364  if (receiveProcFace.proci != Pstream::myProcNo())
365  {
366  td.sendFromPatch = nccpp.index();
367  td.sendToProc = receiveProcFace.proci;
368  td.sendToPatch = nccpp.nbrPatchID();
369  td.sendToPatchFace = receiveProcFace.elementi;
370 
371  return true;
372  }
373 
374  // If both sides are on the same process, then do the local transfer
375  prepareForNonConformalCyclicTransfer
376  (
377  td.mesh,
378  nccpp.index(),
379  receiveProcFace.elementi
380  );
381  correctAfterNonConformalCyclicTransfer
382  (
383  td.mesh,
384  nccpp.nbrPatchID()
385  );
386 
387  return true;
388 }
389 
390 
391 template<class TrackCloudType>
393 {
394  td.sendToProc = cloud.patchNbrProc()[patch(td.mesh)];
395  td.sendFromPatch = patch(td.mesh);
396  td.sendToPatch = cloud.patchNbrProcPatch()[patch(td.mesh)];
397  td.sendToPatchFace =
398  td.mesh.boundaryMesh()[patch(td.mesh)].whichFace(face());
399 }
400 
401 
402 template<class TrackCloudType>
404 {}
405 
406 
407 template<class TrackCloudType>
408 void Foam::particle::hitBasicPatch(TrackCloudType&, trackingData& td)
409 {
410  td.keepParticle = false;
411 }
412 
413 
414 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
A primitive field of type <Type> with automated input and output.
Definition: IOField.H:53
Helper IO class to read and write particle positions.
Definition: IOPosition.H:60
virtual bool write(const bool write=true) const
Write using setting from DB.
Definition: IOPosition.C:52
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
Cyclic plane patch.
const cyclicPolyPatch & nbrPatch() const
virtual const transformer & transform() const
Return transformation between the coupled patches.
virtual label nbrPatchID() const
Neighbour patchID.
label transformGlobalFace(const label facei) const
For a given patch face index, return the corresponding index of the.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Non-conformal cyclic poly patch. As nonConformalCoupledPolyPatch, but the neighbouring patch is local...
remote ray(const scalar fraction, const label origFacei, const vector &p, const vector &n, point &nbrP) const
Compute a ray intersection across the coupling.
const polyPatch & origPatch() const
Original patch.
label sendToProc
Processor to send the particle to. -1 indicates that this.
Definition: particle.H:113
label sendFromPatch
Patch from which to send the particle.
Definition: particle.H:116
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:109
label sendToPatch
Patch to which to send the particle.
Definition: particle.H:119
const polyMesh & mesh
Reference to the mesh.
Definition: particle.H:106
label sendToPatchFace
Patch face to which to send the particle.
Definition: particle.H:122
Base particle class.
Definition: particle.H:83
void prepareForParallelTransfer(TrackCloudType &, trackingData &)
Make changes prior to a parallel transfer. Runs either.
void hitBasicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a basic.
static void readFields(TrackCloudType &c)
Read the fields associated with the owner cloud.
void hitCyclicPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
void hitProcessorPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
scalar trackToAndHitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Convenience function. Combines trackToFace and hitFace.
bool hitPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitFace(const vector &displacement, const scalar fraction, TrackCloudType &cloud, trackingData &td)
Hit the current face. If the current face is internal than this.
void hitSymmetryPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a.
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 correctAfterParallelTransfer(TrackCloudType &, trackingData &)
Make changes following a parallel transfer. Runs either.
bool hitNonConformalCyclicPatch(const vector &displacement, const scalar fraction, const label patchi, TrackCloudType &cloud, trackingData &td)
Overridable function to handle the particle hitting an.
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:188
void hitWallPatch(TrackCloudType &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
static void writeFields(const TrackCloudType &c)
Write the fields associated with the owner cloud.
label index() const
Return the index of this patch in the boundaryMesh.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1387
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:360
virtual bool write(const bool write=true) const
Write using setting from DB.
Struct for keeping processor, element (cell, face, point) index.
Definition: remote.H:57
label elementi
Element index.
Definition: remote.H:66
label proci
Processor index.
Definition: remote.H:63
static transformer rotation(const tensor &T)
Construct a pure rotation transformer.
Definition: transformerI.H:43
bool transformsPosition() const
Return true if the transformer transforms a point.
Definition: transformerI.H:146
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
#define FUNCTION_NAME
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
static const Identity< scalar > I
Definition: Identity.H:93
error FatalError
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
volScalarField & p