particleI.H
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-2018, 2020 OpenFOAM Foundation
9  Copyright (C) 2011-2022 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 "polyMesh.H"
30 #include "Time.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 inline Foam::barycentricTensor Foam::particle::stationaryTetTransform() const
35 {
36  const tetPointRef tet = currentTetIndices().tet(mesh_);
37 
38  return barycentricTensor(tet.a(), tet.b(), tet.c(), tet.d());
39 }
40 
41 
42 inline void Foam::particle::movingTetGeometry
43 (
44  const scalar fraction,
45  Pair<vector>& centre,
46  Pair<vector>& base,
47  Pair<vector>& vertex1,
48  Pair<vector>& vertex2
49 ) const
50 {
51  const triFace triIs(currentTetIndices().faceTriIs(mesh_));
52 
53  const pointField& ptsOld = mesh_.oldPoints();
54  const pointField& ptsNew = mesh_.points();
55 
56  // !!! <-- We would be better off using mesh_.cellCentres() here. However,
57  // we need to put a mesh_.oldCellCentres() method in for this to work. The
58  // values obtained from the mesh and those obtained from the cell do not
59  // necessarily match. See mantis #1993.
60  //const vector ccOld = mesh_.cells()[celli_].centre(ptsOld, mesh_.faces());
61  //const vector ccNew = mesh_.cells()[celli_].centre(ptsNew, mesh_.faces());
62 
63  const vector ccOld = mesh_.oldCellCentres()[celli_];
64  const vector ccNew = mesh_.cellCentres()[celli_];
65 
66  // Old and new points and cell centres are not sub-cycled. If we are sub-
67  // cycling, then we have to account for the timestep change here by
68  // modifying the fractions that we take of the old and new geometry.
69  const Pair<scalar> s = stepFractionSpan();
70  const scalar f0 = s[0] + stepFraction_*s[1], f1 = fraction*s[1];
71 
72  centre[0] = ccOld + f0*(ccNew - ccOld);
73  base[0] = ptsOld[triIs[0]] + f0*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
74  vertex1[0] = ptsOld[triIs[1]] + f0*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
75  vertex2[0] = ptsOld[triIs[2]] + f0*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
76 
77  centre[1] = f1*(ccNew - ccOld);
78  base[1] = f1*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
79  vertex1[1] = f1*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
80  vertex2[1] = f1*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
81 }
82 
83 
84 inline Foam::Pair<Foam::barycentricTensor> Foam::particle::movingTetTransform
85 (
86  const scalar fraction
87 ) const
88 {
89  Pair<vector> centre, base, vertex1, vertex2;
90  movingTetGeometry(fraction, centre, base, vertex1, vertex2);
91 
92  return
93  Pair<barycentricTensor>
94  (
95  barycentricTensor(centre[0], base[0], vertex1[0], vertex2[0]),
96  barycentricTensor(centre[1], base[1], vertex1[1], vertex2[1])
97  );
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
103 inline Foam::label Foam::particle::getNewParticleID() const
104 {
105  label id = particleCount_++;
106 
107  if (id == labelMax)
108  {
110  << "Particle counter has overflowed. This might cause problems"
111  << " when reconstructing particle tracks." << endl;
112  }
113  return id;
114 }
115 
117 inline const Foam::polyMesh& Foam::particle::mesh() const noexcept
118 {
119  return mesh_;
120 }
121 
124 {
125  return coordinates_;
126 }
127 
129 inline Foam::label Foam::particle::cell() const noexcept
130 {
131  return celli_;
132 }
133 
135 inline Foam::label& Foam::particle::cell() noexcept
136 {
137  return celli_;
138 }
139 
141 inline Foam::label Foam::particle::tetFace() const noexcept
142 {
143  return tetFacei_;
144 }
145 
147 inline Foam::label& Foam::particle::tetFace() noexcept
148 {
149  return tetFacei_;
150 }
151 
153 inline Foam::label Foam::particle::tetPt() const noexcept
154 {
155  return tetPti_;
156 }
157 
159 inline Foam::label& Foam::particle::tetPt() noexcept
160 {
161  return tetPti_;
162 }
163 
165 inline Foam::label Foam::particle::face() const noexcept
166 {
167  return facei_;
168 }
169 
171 inline Foam::label& Foam::particle::face() noexcept
172 {
173  return facei_;
174 }
175 
177 inline Foam::scalar Foam::particle::stepFraction() const noexcept
178 {
179  return stepFraction_;
180 }
181 
183 inline Foam::scalar& Foam::particle::stepFraction() noexcept
184 {
185  return stepFraction_;
186 }
187 
189 inline Foam::label Foam::particle::origProc() const noexcept
190 {
191  return origProc_;
192 }
193 
195 inline Foam::label& Foam::particle::origProc() noexcept
196 {
197  return origProc_;
198 }
199 
201 inline Foam::label Foam::particle::origId() const noexcept
202 {
203  return origId_;
204 }
205 
207 inline Foam::label& Foam::particle::origId() noexcept
208 {
209  return origId_;
210 }
211 
212 
214 {
215  if (mesh_.time().subCycling())
216  {
217  const TimeState& tsNew = mesh_.time();
218  const TimeState& tsOld = mesh_.time().prevTimeState();
219 
220  const scalar tFrac =
221  (
222  (tsNew.value() - tsNew.deltaTValue())
223  - (tsOld.value() - tsOld.deltaTValue())
224  )/tsOld.deltaTValue();
225 
226  const scalar dtFrac = tsNew.deltaTValue()/tsOld.deltaTValue();
227 
228  return Pair<scalar>(tFrac, dtFrac);
229  }
230 
231  return Pair<scalar>(0, 1);
232 }
233 
234 
235 inline Foam::scalar Foam::particle::currentTimeFraction() const
236 {
237  const Pair<scalar> s = stepFractionSpan();
238 
239  return s[0] + stepFraction_*s[1];
240 }
241 
244 {
245  return tetIndices(celli_, tetFacei_, tetPti_);
246 }
247 
248 
250 {
251  if (mesh_.moving() && stepFraction_ != 1)
252  {
253  return movingTetTransform(0)[0];
254  }
255 
256  return stationaryTetTransform();
257 }
258 
261 {
262  return currentTetIndices().faceTri(mesh_).unitNormal();
263 }
264 
266 inline bool Foam::particle::onFace() const noexcept
267 {
268  return facei_ >= 0;
269 }
270 
272 inline bool Foam::particle::onInternalFace() const noexcept
273 {
274  return onFace() && mesh_.isInternalFace(facei_);
275 }
276 
278 inline bool Foam::particle::onBoundaryFace() const noexcept
279 {
280  return onFace() && !mesh_.isInternalFace(facei_);
281 }
282 
284 inline Foam::label Foam::particle::patch() const
285 {
286  return onFace() ? mesh_.boundaryMesh().whichPatch(facei_) : -1;
287 }
288 
291 {
292  return currentTetTransform() & coordinates_;
293 }
294 
295 
296 inline void Foam::particle::reset()
297 {
298  stepFraction_ = 0;
299  behind_ = 0;
300  nBehind_ = 0;
301 }
302 
303 
305 {
306  if (!onBoundaryFace())
307  {
309  << "Patch data was requested for a particle that isn't on a patch"
310  << exit(FatalError);
311  }
312 
313  if ((mesh_.moving() && stepFraction_ != 1))
314  {
315  Pair<vector> centre, base, vertex1, vertex2;
316  movingTetGeometry(1, centre, base, vertex1, vertex2);
317 
318  n = triPointRef::unitNormal(base[0], vertex1[0], vertex2[0]);
319 
320  // Interpolate the motion of the three face vertices to the current
321  // coordinates
322  U =
323  coordinates_.b()*base[1]
324  + coordinates_.c()*vertex1[1]
325  + coordinates_.d()*vertex2[1];
326 
327  // The movingTetGeometry method gives the motion as a displacement
328  // across the time-step, so we divide by the time-step to get velocity
329  U /= mesh_.time().deltaTValue();
330  }
331  else
332  {
333  n = currentTetIndices().faceTri(mesh_).unitNormal();
334 
335  U = Zero;
336  }
337 }
338 
339 
340 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
label tetFace() const noexcept
Return current tet face particle is in.
Definition: particleI.H:134
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:49
label tetPt() const noexcept
Return current tet face particle is in.
Definition: particleI.H:146
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:47
const barycentric & coordinates() const noexcept
Return current particle coordinates.
Definition: particleI.H:116
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
label origId() const noexcept
Return the particle ID on the originating processor.
Definition: particleI.H:194
bool onFace() const noexcept
Is the particle on a face?
Definition: particleI.H:259
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Pair< scalar > stepFractionSpan() const
Return the step fraction change within the overall time-step.
Definition: particleI.H:206
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
Definition: tetIndicesI.H:129
vector normal() const
The (unit) normal of the tri on tetFacei_ for the current tet.
Definition: particleI.H:253
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: instant.H:46
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particleI.H:297
vector unitNormal() const
Return unit normal.
Definition: triangleI.H:178
void reset()
Reset particle data.
Definition: particleI.H:289
bool onInternalFace() const noexcept
Is the particle on an internal face?
Definition: particleI.H:265
face triFace(3)
tetIndices currentTetIndices() const noexcept
Return indices of the current tet that the particle occupies.
Definition: particleI.H:236
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:242
Vector< scalar > vector
Definition: vector.H:57
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition: particleI.H:170
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:96
BarycentricTensor< scalar > barycentricTensor
A scalar version of the templated BarycentricTensor.
label cell() const noexcept
Return current cell particle is in.
Definition: particleI.H:122
label origProc() const noexcept
Return the originating processor ID.
Definition: particleI.H:182
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:110
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
constexpr label labelMax
Definition: label.H:55
label n
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label face() const noexcept
Return current face particle is on otherwise -1.
Definition: particleI.H:158
bool onBoundaryFace() const noexcept
Is the particle on a boundary face?
Definition: particleI.H:271
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))
vector position() const
Return current particle position.
Definition: particleI.H:283
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:277
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.
Definition: particleI.H:228