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