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-2024 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 "tracking.H"
28 #include "IOPosition.H"
29 
30 #include "cyclicPolyPatch.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 
47  (
48  c.fieldIOobject("origProcId", IOobject::MUST_READ)
49  );
50 
51  bool haveFile = procIO.headerOk();
52 
53  IOField<label> origProcId(procIO, valid && haveFile);
54  c.checkFieldIOobject(c, origProcId);
56  (
57  c.fieldIOobject("origId", IOobject::MUST_READ),
58  valid && haveFile
59  );
60  c.checkFieldIOobject(c, origId);
61 
62  label i = 0;
63  forAllIter(typename TrackCloudType, c, iter)
64  {
65  particle& p = iter();
66 
67  p.origProc_ = origProcId[i];
68  p.origId_ = origId[i];
69  i++;
70  }
71 }
72 
73 
74 template<class TrackCloudType>
75 void Foam::particle::writeFields(const TrackCloudType& c)
76 {
77  label np = c.size();
78 
80  ioP.write(np > 0);
81 
82  IOField<label> origProc
83  (
84  c.fieldIOobject("origProcId", IOobject::NO_READ),
85  np
86  );
87  IOField<label> origId
88  (
89  c.fieldIOobject("origId", IOobject::NO_READ),
90  np
91  );
92 
93  label i = 0;
94  forAllConstIter(typename TrackCloudType, c, iter)
95  {
96  origProc[i] = iter().origProc_;
97  origId[i] = iter().origId_;
98  i++;
99  }
100 
101  origProc.write(np > 0);
102  origId.write(np > 0);
103 }
104 
105 
106 template<class TrackCloudType>
108 (
109  TrackCloudType& cloud,
110  trackingData& td
111 )
112 {
113  if (td.sendFromPatch == patch(td.mesh))
114  {
115  prepareForProcessorTransfer(td);
116  }
117  else
118  {
119  prepareForNonConformalCyclicTransfer
120  (
121  td.mesh,
122  td.sendFromPatch,
123  td.sendToPatchFace,
124  td.sendToPosition
125  );
126  }
127 }
128 
129 
130 template<class TrackCloudType>
132 (
133  TrackCloudType& cloud,
134  trackingData& td
135 )
136 {
137  const polyPatch& pp = td.mesh.boundaryMesh()[td.sendToPatch];
138 
139  if (isA<processorPolyPatch>(pp))
140  {
141  correctAfterProcessorTransfer(td);
142  }
143  else if (isA<nonConformalCyclicPolyPatch>(pp))
144  {
145  correctAfterNonConformalCyclicTransfer
146  (
147  td.mesh,
148  td.sendToPatch,
150  );
151  }
152  else
153  {
155  << "Transfer patch type not recognised"
156  << exit(FatalError);
157  }
158 }
159 
160 
161 template<class TrackCloudType>
163 (
164  const vector& displacement,
165  const scalar fraction,
166  TrackCloudType& cloud,
167  trackingData& td
168 )
169 {
170  typename TrackCloudType::particleType& p =
171  static_cast<typename TrackCloudType::particleType&>(*this);
172  typename TrackCloudType::particleType::trackingData& ttd =
173  static_cast<typename TrackCloudType::particleType::trackingData&>(td);
174 
175  if (!onFace())
176  {
177  return;
178  }
179  else if (onInternalFace(td.mesh))
180  {
182  (
183  td.mesh,
184  coordinates_,
185  celli_,
186  tetFacei_,
187  tetPti_
188  );
189  }
190  else if (onBoundaryFace(td.mesh))
191  {
192  forAll(cloud.patchNonConformalCyclicPatches()[p.patch(td.mesh)], i)
193  {
194  if
195  (
196  p.hitNonConformalCyclicPatch
197  (
198  displacement,
199  fraction,
200  cloud.patchNonConformalCyclicPatches()[p.patch(td.mesh)][i],
201  cloud,
202  ttd
203  )
204  )
205  {
206  return;
207  }
208  }
209 
210  if (!p.hitPatch(cloud, ttd))
211  {
212  const polyPatch& patch = td.mesh.boundaryMesh()[p.patch(td.mesh)];
213 
214  if (isA<wedgePolyPatch>(patch))
215  {
216  p.hitWedgePatch(cloud, ttd);
217  }
218  else if (isA<symmetryPlanePolyPatch>(patch))
219  {
220  p.hitSymmetryPlanePatch(cloud, ttd);
221  }
222  else if (isA<symmetryPolyPatch>(patch))
223  {
224  p.hitSymmetryPatch(cloud, ttd);
225  }
226  else if (isA<cyclicPolyPatch>(patch))
227  {
228  p.hitCyclicPatch(cloud, ttd);
229  }
230  else if (isA<processorPolyPatch>(patch))
231  {
232  p.hitProcessorPatch(cloud, ttd);
233  }
234  else if (isA<wallPolyPatch>(patch))
235  {
236  p.hitWallPatch(cloud, ttd);
237  }
238  else
239  {
240  p.hitBasicPatch(cloud, ttd);
241  }
242  }
243  }
244 }
245 
246 
247 template<class TrackCloudType>
249 (
250  const vector& displacement,
251  const scalar fraction,
252  TrackCloudType& cloud,
253  trackingData& td
254 )
255 {
256  const scalar f = trackToFace(td.mesh, displacement, fraction);
257 
258  hitFace(displacement, fraction, cloud, td);
259 
260  return f;
261 }
262 
263 
264 template<class TrackCloudType>
265 bool Foam::particle::hitPatch(TrackCloudType&, trackingData&)
266 {
267  return false;
268 }
269 
270 
271 template<class TrackCloudType>
273 {
275  << "Hitting a wedge patch should not be possible."
276  << abort(FatalError);
277 
278  hitSymmetryPatch(cloud, td);
279 }
280 
281 
282 template<class TrackCloudType>
284 (
285  TrackCloudType& cloud,
286  trackingData& td
287 )
288 {
289  hitSymmetryPatch(cloud, td);
290 }
291 
292 
293 template<class TrackCloudType>
295 {
296  const vector nf = normal(td.mesh);
297  transformProperties(transformer::rotation(I - 2.0*nf*nf));
298 }
299 
300 
301 template<class TrackCloudType>
303 {
304  const cyclicPolyPatch& cpp =
305  static_cast<const cyclicPolyPatch&>
306  (
307  td.mesh.boundaryMesh()[patch(td.mesh)]
308  );
309  const cyclicPolyPatch& receiveCpp = cpp.nbrPatch();
310 
311  tracking::crossCyclic(cpp, coordinates_, celli_, tetFacei_, tetPti_);
312 
313  // Remain on the face
314  facei_ = tetFacei_;
315 
316  // Transform the properties
317  if (receiveCpp.transform().transformsPosition())
318  {
319  transformProperties(receiveCpp.transform());
320  }
321 }
322 
323 
324 template<class TrackCloudType>
326 (
327  const vector& displacement,
328  const scalar fraction,
329  const label patchi,
330  TrackCloudType& cloud,
331  trackingData& td
332 )
333 {
334  const nonConformalCyclicPolyPatch& nccpp =
335  static_cast<const nonConformalCyclicPolyPatch&>
336  (td.mesh.boundaryMesh()[patchi]);
337 
338  const point sendPos = position(td.mesh);
339 
340  // Get the send patch data
341  vector sendNormal, sendDisplacement;
342  patchData(td.mesh, sendNormal, sendDisplacement);
343 
344  // Project the particle through the non-conformal patch
345  point receivePos;
346  const remote receiveProcFace =
347  nccpp.ray
348  (
349  stepFraction_,
350  nccpp.origPatch().whichFace(facei_),
351  sendPos,
352  displacement - fraction*sendDisplacement,
353  receivePos
354  );
355 
356  // If we didn't hit anything then this particle is assumed to project to
357  // the orig patch, or another different non-conformal patch. Return, so
358  // these can be tried.
359  if (receiveProcFace.proci == -1) return false;
360 
361  // If we are transferring between processes then mark as such and return.
362  // The cloud will handle all processor transfers as a single batch.
363  if (receiveProcFace.proci != Pstream::myProcNo())
364  {
365  td.sendToProc = receiveProcFace.proci;
366  td.sendFromPatch = nccpp.index();
367  td.sendToPatch = nccpp.nbrPatchIndex();
368  td.sendToPatchFace = receiveProcFace.elementi;
369  td.sendToPosition = receivePos;
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  receivePos
381  );
382  correctAfterNonConformalCyclicTransfer
383  (
384  td.mesh,
385  nccpp.nbrPatchIndex(),
387  );
388 
389  return true;
390 }
391 
392 
393 template<class TrackCloudType>
395 {
396  td.sendToProc = cloud.patchNbrProc()[patch(td.mesh)];
397  td.sendFromPatch = patch(td.mesh);
398  td.sendToPatch = cloud.patchNbrProcPatch()[patch(td.mesh)];
399  td.sendToPatchFace =
400  td.mesh.boundaryMesh()[patch(td.mesh)].whichFace(face());
401 }
402 
403 
404 template<class TrackCloudType>
406 {}
407 
408 
409 template<class TrackCloudType>
410 void Foam::particle::hitBasicPatch(TrackCloudType&, trackingData& td)
411 {
412  td.keepParticle = false;
413 }
414 
415 
416 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
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:59
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
Cyclic plane patch.
const cyclicPolyPatch & nbrPatch() const
virtual const transformer & transform() const
Return transformation between the coupled patches.
virtual label nbrPatchIndex() const
Neighbour patchID.
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.
labelList patchNLocateBoundaryHits
Number of boundary hits that occurred during locate executions.
Definition: particle.H:125
label sendToProc
Processor to send the particle to. -1 indicates that this.
Definition: particle.H:109
label sendFromPatch
Patch from which to send the particle.
Definition: particle.H:112
vector sendToPosition
Position to which to send.
Definition: particle.H:121
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
label sendToPatch
Patch to which to send the particle.
Definition: particle.H:115
const polyMesh & mesh
Reference to the mesh.
Definition: particle.H:102
label sendToPatchFace
Patch face to which to send the particle.
Definition: particle.H:118
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:108
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.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
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:70
label proci
Processor index.
Definition: remote.H:67
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:530
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
bool valid(const PtrList< ModelType > &l)
const dimensionedScalar c
Speed of light in a vacuum.
void crossInternalFace(const polyMesh &mesh, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Cross an internal face.
Definition: tracking.C:1681
void crossCyclic(const cyclicPolyPatch &inPatch, barycentric &coordinates, label &celli, label &facei, label &faceTrii)
Cross a cyclic patch.
Definition: tracking.C:1736
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const Identity< scalar > I
Definition: Identity.H:93
error FatalError
labelList f(nPoints)
volScalarField & p