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 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "polyMesh.H"
29 #include "Time.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::particle::stationaryTetGeometry
34 (
35  vector& centre,
36  vector& base,
37  vector& vertex1,
38  vector& vertex2
39 ) const
40 {
41  const triFace triIs(currentTetIndices().faceTriIs(mesh_));
42  const vectorField& ccs = mesh_.cellCentres();
43  const pointField& pts = mesh_.points();
44 
45  centre = ccs[celli_];
46  base = pts[triIs[0]];
47  vertex1 = pts[triIs[1]];
48  vertex2 = pts[triIs[2]];
49 }
50 
51 
52 inline Foam::barycentricTensor Foam::particle::stationaryTetTransform() const
53 {
54  vector centre, base, vertex1, vertex2;
55  stationaryTetGeometry(centre, base, vertex1, vertex2);
56 
57  return barycentricTensor(centre, base, vertex1, vertex2);
58 }
59 
60 
61 inline void Foam::particle::movingTetGeometry
62 (
63  const scalar fraction,
64  Pair<vector>& centre,
65  Pair<vector>& base,
66  Pair<vector>& vertex1,
67  Pair<vector>& vertex2
68 ) const
69 {
70  const triFace triIs(currentTetIndices().faceTriIs(mesh_));
71  const pointField& ptsOld = mesh_.oldPoints();
72  const pointField& ptsNew = mesh_.points();
73 
74  // !!! <-- We would be better off using mesh_.cellCentres() here. However,
75  // we need to put a mesh_.oldCellCentres() method in for this to work. The
76  // values obtained from the mesh and those obtained from the cell do not
77  // necessarily match. See mantis #1993.
78  const vector ccOld = mesh_.cells()[celli_].centre(ptsOld, mesh_.faces());
79  const vector ccNew = mesh_.cells()[celli_].centre(ptsNew, mesh_.faces());
80 
81  // Old and new points and cell centres are not sub-cycled. If we are sub-
82  // cycling, then we have to account for the timestep change here by
83  // modifying the fractions that we take of the old and new geometry.
84  const Pair<scalar> s = stepFractionSpan();
85  const scalar f0 = s[0] + stepFraction_*s[1], f1 = fraction*s[1];
86 
87  centre[0] = ccOld + f0*(ccNew - ccOld);
88  base[0] = ptsOld[triIs[0]] + f0*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
89  vertex1[0] = ptsOld[triIs[1]] + f0*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
90  vertex2[0] = ptsOld[triIs[2]] + f0*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
91 
92  centre[1] = f1*(ccNew - ccOld);
93  base[1] = f1*(ptsNew[triIs[0]] - ptsOld[triIs[0]]);
94  vertex1[1] = f1*(ptsNew[triIs[1]] - ptsOld[triIs[1]]);
95  vertex2[1] = f1*(ptsNew[triIs[2]] - ptsOld[triIs[2]]);
96 }
97 
98 
99 inline Foam::Pair<Foam::barycentricTensor> Foam::particle::movingTetTransform
100 (
101  const scalar fraction
102 ) const
103 {
104  Pair<vector> centre, base, vertex1, vertex2;
105  movingTetGeometry(fraction, centre, base, vertex1, vertex2);
106 
107  return
108  Pair<barycentricTensor>
109  (
110  barycentricTensor(centre[0], base[0], vertex1[0], vertex2[0]),
111  barycentricTensor(centre[1], base[1], vertex1[1], vertex2[1])
112  );
113 }
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  label id = particleCount_++;
121 
122  if (id == labelMax)
123  {
125  << "Particle counter has overflowed. This might cause problems"
126  << " when reconstructing particle tracks." << endl;
127  }
128  return id;
129 }
130 
131 
133 {
134  return mesh_;
135 }
136 
137 
139 {
140  return coordinates_;
141 }
142 
143 
145 {
146  return celli_;
147 }
148 
149 
151 {
152  return celli_;
153 }
154 
155 
157 {
158  return tetFacei_;
159 }
160 
161 
163 {
164  return tetFacei_;
165 }
166 
167 
169 {
170  return tetPti_;
171 }
172 
173 
175 {
176  return tetPti_;
177 }
178 
179 
181 {
182  return facei_;
183 }
184 
185 
187 {
188  return facei_;
189 }
190 
191 
192 inline Foam::scalar Foam::particle::stepFraction() const
193 {
194  return stepFraction_;
195 }
196 
197 
198 inline Foam::scalar& Foam::particle::stepFraction()
199 {
200  return stepFraction_;
201 }
202 
203 
205 {
206  return origProc_;
207 }
208 
209 
211 {
212  return origProc_;
213 }
214 
215 
217 {
218  return origId_;
219 }
220 
221 
223 {
224  return origId_;
225 }
226 
227 
229 {
230  if (mesh_.time().subCycling())
231  {
232  const TimeState& tsNew = mesh_.time();
233  const TimeState& tsOld = mesh_.time().prevTimeState();
234 
235  const scalar tFrac =
236  (
237  (tsNew.value() - tsNew.deltaTValue())
238  - (tsOld.value() - tsOld.deltaTValue())
239  )/tsOld.deltaTValue();
240 
241  const scalar dtFrac = tsNew.deltaTValue()/tsOld.deltaTValue();
242 
243  return Pair<scalar>(tFrac, dtFrac);
244  }
245  else
246  {
247  return Pair<scalar>(0, 1);
248  }
249 }
250 
251 
252 inline Foam::scalar Foam::particle::currentTimeFraction() const
253 {
254  const Pair<scalar> s = stepFractionSpan();
255 
256  return s[0] + stepFraction_*s[1];
257 }
258 
259 
261 {
262  return tetIndices(celli_, tetFacei_, tetPti_);
263 }
264 
265 
267 {
268  if (mesh_.moving())
269  {
270  return movingTetTransform(0)[0];
271  }
272  else
273  {
274  return stationaryTetTransform();
275  }
276 }
277 
278 
280 {
281  return currentTetIndices().faceTri(mesh_).unitNormal();
282 }
283 
284 
285 inline bool Foam::particle::onFace() const
286 {
287  return facei_ >= 0;
288 }
289 
290 
292 {
293  return onFace() && mesh_.isInternalFace(facei_);
294 }
295 
296 
298 {
299  return onFace() && !mesh_.isInternalFace(facei_);
300 }
301 
302 
304 {
305  return onFace() ? mesh_.boundaryMesh().whichPatch(facei_) : -1;
306 }
307 
308 
310 {
311  return currentTetTransform() & coordinates_;
312 }
313 
314 
316 {
317  if (!onBoundaryFace())
318  {
320  << "Patch data was requested for a particle that isn't on a patch"
321  << exit(FatalError);
322  }
323 
324  if (mesh_.moving())
325  {
326  Pair<vector> centre, base, vertex1, vertex2;
327  movingTetGeometry(1, centre, base, vertex1, vertex2);
328 
329  n = triPointRef(base[0], vertex1[0], vertex2[0]).unitNormal();
330 
331  // Interpolate the motion of the three face vertices to the current
332  // coordinates
333  U =
334  coordinates_.b()*base[1]
335  + coordinates_.c()*vertex1[1]
336  + coordinates_.d()*vertex2[1];
337 
338  // The movingTetGeometry method gives the motion as a displacement
339  // across the time-step, so we divide by the time-step to get velocity
340  U /= mesh_.time().deltaTValue();
341  }
342  else
343  {
344  vector centre, base, vertex1, vertex2;
345  stationaryTetGeometry(centre, base, vertex1, vertex2);
346 
347  n = triPointRef(base, vertex1, vertex2).unitNormal();
348 
349  U = Zero;
350  }
351 }
352 
353 
354 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1038
Foam::labelMax
constexpr label labelMax
Definition: label.H:65
Foam::barycentricTensor
BarycentricTensor< scalar > barycentricTensor
A scalar version of the templated BarycentricTensor.
Definition: barycentricTensor.H:48
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
Foam::particle::onFace
bool onFace() const
Is the particle on a face?
Definition: particleI.H:285
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::particle::origId
label origId() const
Return the particle ID on the originating processor.
Definition: particleI.H:216
Foam::particle::stepFraction
scalar stepFraction() const
Return the fraction of time-step completed.
Definition: particleI.H:192
Foam::particle::tetFace
label tetFace() const
Return current tet face particle is in.
Definition: particleI.H:156
Foam::particle::getNewParticleID
label getNewParticleID() const
Get unique particle creation id.
Definition: particleI.H:118
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::particle::currentTetIndices
tetIndices currentTetIndices() const
Return the indices of the current tet that the.
Definition: particleI.H:260
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
polyMesh.H
Foam::particle::onInternalFace
bool onInternalFace() const
Is the particle on an internal face?
Definition: particleI.H:291
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::particle::coordinates
const barycentric & coordinates() const
Return current particle coordinates.
Definition: particleI.H:138
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::TimeState
The time value with time-stepping information, user-defined remapping, etc.
Definition: TimeState.H:51
Foam::label
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:62
Foam::Vector::centre
const Vector< Cmpt > & centre(const Foam::List< Vector< Cmpt >> &) const
Return *this (used for point which is a typedef to Vector<scalar>.
Definition: VectorI.H:119
Foam::particle::position
vector position() const
Return current particle position.
Definition: particleI.H:309
Foam::Barycentric< scalar >
Foam::particle::patch
label patch() const
Return the index of patch that the particle is on.
Definition: particleI.H:303
Foam::particle::mesh
const polyMesh & mesh() const
Return the mesh database.
Definition: particleI.H:132
Foam::particle::currentTimeFraction
scalar currentTimeFraction() const
Return the current fraction within the timestep. This differs.
Definition: particleI.H:252
Foam::particle::origProc
label origProc() const
Return the originating processor ID.
Definition: particleI.H:204
Foam::FatalError
error FatalError
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::particle::onBoundaryFace
bool onBoundaryFace() const
Is the particle on a boundary face?
Definition: particleI.H:297
Time.H
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:83
U
U
Definition: pEqn.H:72
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:176
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Pair
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: Pair.H:54
Foam::Vector< scalar >
Foam::TimeState::deltaTValue
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:42
Foam::particle::stepFractionSpan
Pair< scalar > stepFractionSpan() const
Return the step fraction change within the overall time-step.
Definition: particleI.H:228
Foam::BarycentricTensor
Templated 4x3 tensor derived from VectorSpace. Has 12 components. Can represent a barycentric transfo...
Definition: BarycentricTensor.H:57
Foam::particle::tetPt
label tetPt() const
Return current tet face particle is in.
Definition: particleI.H:168
Foam::particle::normal
vector normal() const
The (unit) normal of the tri on tetFacei_ for the current tet.
Definition: particleI.H:279
Foam::particle::cell
label cell() const
Return current cell particle is in.
Definition: particleI.H:144
Foam::particle::face
label face() const
Return current face particle is on otherwise -1.
Definition: particleI.H:180
triFace
face triFace(3)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::triPointRef
triangle< point, const point & > triPointRef
Definition: triangle.H:78
Foam::particle::patchData
void patchData(vector &n, vector &U) const
Get the normal and velocity of the current patch location.
Definition: particleI.H:315
Foam::particle::currentTetTransform
barycentricTensor currentTetTransform() const
Return the current tet transformation tensor.
Definition: particleI.H:266