cellPoint_LagrangianAverage.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) 2025 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 
28 #include "LagrangianFields.H"
29 #include "oneField.H"
30 #include "tetIndices.H"
31 #include "syncTools.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 template<class Type>
36 template<class CellWeight>
39 (
40  const polyMesh& pMesh,
41  const CellWeight& cellWeight
42 )
43 {
44  tmp<Field<scalar>> tresult(new Field<scalar>(pMesh.nCells(), scalar(0)));
45  Field<scalar>& result = tresult.ref();
46 
47  forAll(pMesh.cells(), celli)
48  {
49  const Foam::cell& c = pMesh.cells()[celli];
50 
51  forAll(c, cFacei)
52  {
53  const label facei = c[cFacei];
54  const Foam::face& f = pMesh.faces()[facei];
55 
56  for (label fTrii = 0; fTrii < f.nTriangles(); ++ fTrii)
57  {
58  const tetIndices tetPoints(celli, facei, fTrii + 1);
59  const scalar v = tetPoints.tet(pMesh).mag();
60 
61  result[celli] += v*cellWeight[celli];
62  }
63  }
64  }
65 
66  return tresult;
67 }
68 
69 
70 template<class Type>
71 template<class CellWeight>
74 (
75  const polyMesh& pMesh,
76  const CellWeight& cellWeight
77 )
78 {
79  tmp<Field<scalar>> tresult(new Field<scalar>(pMesh.nPoints(), scalar(0)));
80  Field<scalar>& result = tresult.ref();
81 
82  forAll(pMesh.cells(), celli)
83  {
84  const Foam::cell& c = pMesh.cells()[celli];
85 
86  forAll(c, cFacei)
87  {
88  const label facei = c[cFacei];
89  const Foam::face& f = pMesh.faces()[facei];
90 
91  for (label fTrii = 0; fTrii < f.nTriangles(); ++ fTrii)
92  {
93  const tetIndices tetPoints(celli, facei, fTrii + 1);
94  const scalar v = tetPoints.tet(pMesh).mag();
95  const triFace triPoints = tetPoints.faceTriIs(pMesh);
96 
97  forAll(triPoints, triPointi)
98  {
99  const label pointi = triPoints[triPointi];
100 
101  result[pointi] += v*cellWeight[celli];
102  }
103  }
104  }
105  }
106 
107  syncTools::syncPointList
108  (
109  pMesh,
110  result,
112  scalar(0)
113  );
114 
115  return tresult;
116 }
117 
118 
119 template<class Type>
122 (
123  const polyMesh& pMesh,
124  const Field<scalar>& weightSum
125 )
126 {
127  return weightSum/4;
128 }
129 
130 
131 template<class Type>
134 (
135  const polyMesh& pMesh,
136  const Field<scalar>& weightSum
137 )
138 {
139  return
140  pointVolumeWeightedSum
141  (
142  pMesh,
143  (weightSum/cellVolumeWeightedSum(pMesh, oneField()))()
144  )/4;
145 }
146 
147 
148 template<class Type>
150 {
151  forAll(d.cellAvgCell_, cellAvgi)
152  {
153  d.cellCellAvg_[d.cellAvgCell_[cellAvgi]] = -1;
154  }
155  forAll(d.pointAvgPoint_, pointAvgi)
156  {
157  d.pointPointAvg_[d.pointAvgPoint_[pointAvgi]] = -1;
158  }
159  d.cellAvgCell_.clear();
160  d.pointAvgPoint_.clear();
161  d.cellAvgCount_.clear();
162  d.pointAvgCount_.clear();
163  if (d.cellAvgWeightSumPtr_.valid()) d.cellAvgWeightSumPtr_->clear();
164  if (d.pointAvgWeightSumPtr_.valid()) d.pointAvgWeightSumPtr_->clear();
165  d.cellAvgSum_.clear();
166  d.pointAvgSum_.clear();
167 }
168 
169 
170 template<class Type>
172 (
173  const LagrangianSubSubField<scalar>& weightOrNull,
174  const LagrangianSubSubField<Type>& psiOrWeightPsi,
175  data& d
176 )
177 {
178  const LagrangianSubMesh& subMesh = psiOrWeightPsi.mesh();
179  const LagrangianMesh& mesh = subMesh.mesh();
180 
181  if (notNull(weightOrNull) != d.cellAvgWeightSumPtr_.valid())
182  {
184  << "Inconsistent weight specification for average of field "
185  << psiOrWeightPsi.name() << exit(FatalError);
186  }
187 
188  // Remove from the cell averages
189  forAll(psiOrWeightPsi, subi)
190  {
191  const label i = subMesh.start() + subi;
192 
193  const barycentric& coordinates = mesh.coordinates()[i];
194  const label celli = mesh.celli()[i];
195 
196  const label cellAvgi = d.cellCellAvg_[celli];
197 
198  // Check that this cell can be removed from
199  if (cellAvgi == -1 || d.cellAvgCount_[cellAvgi] == 0)
200  {
202  << "Negative cell count for average of field "
203  << psiOrWeightPsi.name() << exit(FatalError);
204  }
205 
206  // Remove from the cell average
207  d.cellAvgCount_[cellAvgi] --;
208  if (notNull(weightOrNull))
209  {
210  d.cellAvgWeightSumPtr_()[cellAvgi] -=
211  coordinates.a()*weightOrNull[subi];
212  d.cellAvgSum_[cellAvgi] -=
213  coordinates.a()*weightOrNull[subi]*psiOrWeightPsi[subi];
214  }
215  else
216  {
217  d.cellAvgSum_[cellAvgi] -= coordinates.a()*psiOrWeightPsi[subi];
218  }
219  }
220 }
221 
222 
223 template<class Type>
225 (
226  const LagrangianSubSubField<scalar>& weightOrNull,
227  const LagrangianSubSubField<Type>& psiOrWeightPsi,
228  data& d
229 )
230 {
231  const LagrangianSubMesh& subMesh = psiOrWeightPsi.mesh();
232  const LagrangianMesh& mesh = subMesh.mesh();
233 
234  if (notNull(weightOrNull) != d.cellAvgWeightSumPtr_.valid())
235  {
237  << "Inconsistent weight specification for average of field "
238  << psiOrWeightPsi.name() << exit(FatalError);
239  }
240 
241  // Add to the cell averages
242  forAll(psiOrWeightPsi, subi)
243  {
244  const label i = subMesh.start() + subi;
245 
246  const barycentric& coordinates = mesh.coordinates()[i];
247  const label celli = mesh.celli()[i];
248 
249  label cellAvgi = d.cellCellAvg_[celli];
250 
251  // Initialise this cell if it is newly a part of the average
252  if (cellAvgi == -1)
253  {
254  cellAvgi = d.cellAvgCell_.size();
255  d.cellCellAvg_[celli] = cellAvgi;
256  d.cellAvgCell_.append(celli);
257  d.cellAvgCount_.append(label(0));
258  if (notNull(weightOrNull))
259  {
260  d.cellAvgWeightSumPtr_().append(scalar(0));
261  }
262  d.cellAvgSum_.append(pTraits<Type>::zero);
263  }
264 
265  // Add to the cell average
266  d.cellAvgCount_[cellAvgi] ++;
267  if (notNull(weightOrNull))
268  {
269  d.cellAvgWeightSumPtr_()[cellAvgi] +=
270  coordinates.a()*weightOrNull[subi];
271  d.cellAvgSum_[cellAvgi] +=
272  coordinates.a()*weightOrNull[subi]*psiOrWeightPsi[subi];
273  }
274  else
275  {
276  d.cellAvgSum_[cellAvgi] +=
277  coordinates.a()*psiOrWeightPsi[subi];
278  }
279  }
280 }
281 
282 
283 template<class Type>
285 (
286  const LagrangianSubSubField<scalar>& weightOrNull,
287  const LagrangianSubSubField<Type>& psiOrWeightPsi,
288  data& d
289 )
290 {
291  const LagrangianSubMesh& subMesh = psiOrWeightPsi.mesh();
292  const LagrangianMesh& mesh = subMesh.mesh();
293 
294  if (notNull(weightOrNull) != d.cellAvgWeightSumPtr_.valid())
295  {
297  << "Inconsistent weight specification for average of field "
298  << psiOrWeightPsi.name() << exit(FatalError);
299  }
300 
301  // Add to the point averages
302  forAll(psiOrWeightPsi, subi)
303  {
304  const label i = subMesh.start() + subi;
305 
306  const barycentric& coordinates = mesh.coordinates()[i];
307  const label celli = mesh.celli()[i];
308  const label facei = mesh.facei()[i];
309  const label faceTrii = mesh.faceTrii()[i];
310 
311  const triFace triPoints =
312  tetIndices(celli, facei, faceTrii).faceTriIs(mesh.mesh());
313 
314  forAll(triPoints, triPointi)
315  {
316  const label pointi = triPoints[triPointi];
317 
318  label pointAvgi = d.pointPointAvg_[pointi];
319 
320  // Initialise this point if it is newly a part of the average
321  if (pointAvgi == -1)
322  {
323  pointAvgi = d.pointAvgPoint_.size();
324  d.pointPointAvg_[pointi] = pointAvgi;
325  d.pointAvgPoint_.append(pointi);
326  d.pointAvgCount_.append(label(0));
327  if (notNull(weightOrNull))
328  {
329  d.pointAvgWeightSumPtr_().append(scalar(0));
330  }
331  d.pointAvgSum_.append(pTraits<Type>::zero);
332  }
333 
334  // Add to the point average
335  d.pointAvgCount_[pointAvgi] ++;
336  if (notNull(weightOrNull))
337  {
338  d.pointAvgWeightSumPtr_()[pointAvgi] +=
339  coordinates[triPointi + 1]*weightOrNull[subi];
340  d.pointAvgSum_[pointAvgi] +=
341  coordinates[triPointi + 1]
342  *weightOrNull[subi]
343  *psiOrWeightPsi[subi];
344  }
345  else
346  {
347  d.pointAvgSum_[pointAvgi] +=
348  coordinates[triPointi + 1]*psiOrWeightPsi[subi];
349  }
350  }
351  }
352 
353  // Expand to include any points that are being added to remotely
355  (
356  d.pointPointAvg_,
357  d.pointAvgPoint_
358  );
359  d.pointAvgCount_.resize(d.pointAvgPoint_.size(), label(0));
360  if (notNull(weightOrNull))
361  {
362  d.pointAvgWeightSumPtr_().resize(d.pointAvgPoint_.size(), scalar(0));
363  }
364  d.pointAvgSum_.resize(d.pointAvgPoint_.size(), pTraits<Type>::zero);
365 
366  // Sum on coupled points
367  syncTools::syncPointList
368  (
369  mesh.mesh(),
370  d.pointAvgPoint_,
371  d.pointAvgCount_,
372  plusEqOp<label>(),
373  label(0)
374  );
375  if (notNull(weightOrNull))
376  {
377  syncTools::syncPointList
378  (
379  mesh.mesh(),
380  d.pointAvgPoint_,
381  d.pointAvgWeightSumPtr_(),
382  plusEqOp<scalar>(),
383  scalar(0)
384  );
385  }
386  syncTools::syncPointList
387  (
388  mesh.mesh(),
389  d.pointAvgPoint_,
390  d.pointAvgSum_,
391  plusEqOp<Type>(),
392  pTraits<Type>::zero
393  );
394 }
395 
396 
397 template<class Type>
399 (
400  const data& dd,
401  data& d
402 )
403 {
404  // Remove from the point averages
405  forAll(dd.pointAvgPoint_, dPointAvgi)
406  {
407  const label pointi = dd.pointAvgPoint_[dPointAvgi];
408  const label pointAvgi = d.pointPointAvg_[pointi];
409 
410  // Check that this point can be removed from
411  if
412  (
413  pointAvgi == -1
414  || d.pointAvgCount_[pointAvgi] < dd.pointAvgCount_[dPointAvgi]
415  )
416  {
418  << "Negative point count for average "
419  << exit(FatalError);
420  }
421 
422  // Remove from the point average
423  d.pointAvgCount_[pointAvgi] -= dd.pointAvgCount_[dPointAvgi];
424  if (d.pointAvgWeightSumPtr_.valid())
425  {
426  d.pointAvgWeightSumPtr_()[pointAvgi] -=
427  dd.pointAvgWeightSumPtr_()[dPointAvgi];
428  }
429  d.pointAvgSum_[pointAvgi] -= dd.pointAvgSum_[dPointAvgi];
430  }
431 }
432 
433 
434 template<class Type>
436 (
437  const data& dd,
438  data& d
439 )
440 {
441  forAll(dd.pointAvgPoint_, dPointAvgi)
442  {
443  const label pointi = dd.pointAvgPoint_[dPointAvgi];
444 
445  label pointAvgi = d.pointPointAvg_[pointi];
446 
447  // Initialise this point if it is newly a part of the average
448  if (pointAvgi == -1)
449  {
450  pointAvgi = d.pointAvgPoint_.size();
451  d.pointPointAvg_[pointi] = pointAvgi;
452  d.pointAvgPoint_.append(pointi);
453  d.pointAvgCount_.append(label(0));
454  if (d.pointAvgWeightSumPtr_.valid())
455  {
456  d.pointAvgWeightSumPtr_().append(scalar(0));
457  }
458  d.pointAvgSum_.append(pTraits<Type>::zero);
459  }
460 
461  // Add to the point average
462  d.pointAvgCount_[pointAvgi] += dd.pointAvgCount_[dPointAvgi];
463  if (d.pointAvgWeightSumPtr_.valid())
464  {
465  d.pointAvgWeightSumPtr_()[pointAvgi] +=
466  dd.pointAvgWeightSumPtr_()[dPointAvgi];
467  }
468  d.pointAvgSum_[pointAvgi] += dd.pointAvgSum_[dPointAvgi];
469  }
470 }
471 
472 
473 template<class Type>
475 (
476  LagrangianSubField<Type>& result
477 ) const
478 {
479  const LagrangianSubMesh& subMesh = result.mesh();
480  const LagrangianMesh& mesh = subMesh.mesh();
481 
482  forAll(subMesh, subi)
483  {
484  const label i = subMesh.start() + subi;
485 
486  const barycentric& coordinates = mesh.coordinates()[i];
487  const label celli = mesh.celli()[i];
488  const label facei = mesh.facei()[i];
489  const label faceTrii = mesh.faceTrii()[i];
490 
491  const triFace triPoints =
492  tetIndices(celli, facei, faceTrii).faceTriIs(mesh.mesh());
493 
494  const label cellAvgi = data_.cellCellAvg_[celli];
495  const label dCellAvgi = dData_.cellCellAvg_[celli];
496  const bool haveData =
497  cellAvgi != -1 && data_.cellAvgCount_[cellAvgi] != 0;
498  const bool haveDData =
499  dCellAvgi != -1 && dData_.cellAvgCount_[dCellAvgi] != 0;
500 
501  if (!haveData && !haveDData)
502  {
504  << "Interpolated value requested for a cell in which no "
505  << "elements have been averaged"
506  << exit(FatalError);
507  }
508 
509  static const Type& z = pTraits<Type>::zero;
510 
511  Type wr =
512  coordinates.a()
513  *(
514  (haveData ? data_.cellAvgSum_[cellAvgi] : z)
515  + (haveDData ? dData_.cellAvgSum_[dCellAvgi] : z)
516  );
517 
518  scalar w =
519  coordinates.a()
520  *(
521  !cellWeightSumPtr_.valid()
522  ? (haveData ? data_.cellAvgWeightSumPtr_()[cellAvgi] : 0)
523  + (haveDData ? dData_.cellAvgWeightSumPtr_()[dCellAvgi] : 0)
524  : cellWeightSumPtr_()[celli]
525  );
526 
527  forAll(triPoints, triPointi)
528  {
529  const label pointi = triPoints[triPointi];
530 
531  const label pointAvgi = data_.pointPointAvg_[pointi];
532  const label dPointAvgi = dData_.pointPointAvg_[pointi];
533  const bool haveData =
534  pointAvgi != -1 && data_.pointAvgCount_[pointAvgi] != 0;
535  const bool haveDData =
536  dPointAvgi != -1 && dData_.pointAvgCount_[dPointAvgi] != 0;
537 
538  wr +=
539  coordinates[triPointi + 1]
540  *(
541  (haveData ? data_.pointAvgSum_[pointAvgi] : z)
542  + (haveDData ? dData_.pointAvgSum_[dPointAvgi] : z)
543  );
544 
545  w +=
546  coordinates[triPointi + 1]
547  *(
548  !pointWeightSumPtr_.valid()
549  ? (haveData ? data_.pointAvgWeightSumPtr_()[pointAvgi] : 0)
550  + (haveDData ? dData_.pointAvgWeightSumPtr_()[dPointAvgi] : 0)
551  : pointWeightSumPtr_()[pointi]
552  );
553  }
554 
555  result[subi] = wr/w;
556  }
557 }
558 
559 
560 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
561 
562 template<class Type>
564 (
565  const label nCells,
566  const label nPoints,
567  const bool hasWeightSum
568 )
569 :
570  cellCellAvg_(nCells, label(-1)),
571  pointPointAvg_(nPoints, label(-1)),
572  cellAvgCell_(),
573  pointAvgPoint_(),
574  cellAvgCount_(),
575  pointAvgCount_(),
576  cellAvgWeightSumPtr_(hasWeightSum ? new DynamicList<scalar>() : nullptr),
577  pointAvgWeightSumPtr_(hasWeightSum ? new DynamicList<scalar>() : nullptr),
578  cellAvgSum_(),
579  pointAvgSum_()
580 {}
581 
582 
583 template<class Type>
585 (
586  const word& name,
587  const LagrangianMesh& mesh,
588  const dimensionSet& dimensions,
589  const Field<scalar>& weightSum
590 )
591 :
593  cellWeightSumPtr_
594  (
595  !isNull(weightSum)
596  ? cellWeightSum(mesh.mesh(), weightSum).ptr()
597  : nullptr
598  ),
599  pointWeightSumPtr_
600  (
601  !isNull(weightSum)
602  ? pointWeightSum(mesh.mesh(), weightSum).ptr()
603  : nullptr
604  ),
605  data_(mesh.mesh().nCells(), mesh.mesh().nPoints(), isNull(weightSum)),
606  dDataIsValid_(false),
607  dData_(mesh.mesh().nCells(), mesh.mesh().nPoints(), isNull(weightSum))
608 {}
609 
610 
611 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
612 
613 template<class Type>
615 {}
616 
617 
618 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
619 
620 template<class Type>
622 (
623  const LagrangianSubSubField<scalar>& weightOrNull,
624  const LagrangianSubSubField<Type>& psiOrWeightPsi
625 )
626 {
627  if (dDataIsValid_)
628  {
630  << "Cannot remove from an average with a cached difference"
631  << exit(FatalError);
632  }
633 
634  removeFromCells(weightOrNull, psiOrWeightPsi, data_);
635 
636  addToPoints(weightOrNull, psiOrWeightPsi, dData_);
637 
638  removeFromPoints(dData_, data_);
639 
640  clear(dData_);
641 }
642 
643 
644 template<class Type>
646 (
647  const LagrangianSubSubField<scalar>& weightOrNull,
648  const LagrangianSubSubField<Type>& psiOrWeightPsi,
649  const bool cache
650 )
651 {
652  if (dDataIsValid_)
653  {
655  << "Cannot add to an average with a cached difference"
656  << exit(FatalError);
657  }
658 
659  dDataIsValid_ = cache;
660 
661  addToCells(weightOrNull, psiOrWeightPsi, cache ? dData_ : data_);
662 
663  addToPoints(weightOrNull, psiOrWeightPsi, dData_);
664 
665  if (!cache)
666  {
667  addToPoints(dData_, data_);
668 
669  clear(dData_);
670  }
671 }
672 
673 
674 template<class Type>
676 (
677  const LagrangianSubSubField<scalar>& weightOrNull,
678  const LagrangianSubSubField<Type>& psiOrWeightPsi,
679  const bool cache
680 )
681 {
682  if (!dDataIsValid_)
683  {
685  << "Cannot correct an average without a cached difference"
686  << exit(FatalError);
687  }
688 
689  clear(dData_);
690 
691  dDataIsValid_ = cache;
692 
693  addToCells(weightOrNull, psiOrWeightPsi, cache ? dData_ : data_);
694 
695  addToPoints(weightOrNull, psiOrWeightPsi, dData_);
696 
697  if (!cache)
698  {
699  addToPoints(dData_, data_);
700 
701  clear(dData_);
702  }
703 }
704 
705 
706 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
const Cmpt & a() const
Definition: BarycentricI.H:59
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Basic Lagrangian averaging process in which values are averaged in the cells and then interpolated as...
virtual void remove(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi)
Remove weighted values from the average.
virtual void add(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)
Add weighted values to the average.
cellPoint(const word &name, const LagrangianMesh &mesh, const dimensionSet &dimensions, const Field< scalar > &weightSum)
Construct with a name, for a mesh and with given dimensions.
virtual void correct(const LagrangianSubSubField< scalar > &weight, const LagrangianSubSubField< Type > &psi, const bool cache)
Correct weighted values in the average.
Class containing Lagrangian geometry and topology.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
Dimension set for the base types.
Definition: dimensionSet.H:125
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:443
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:53
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
label nCells() const
label nPoints() const
const cellList & cells() const
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
triFace faceTriIs(const polyMesh &mesh) const
Return the indices corresponding to the tri on the face for.
Definition: tetIndicesI.H:67
tetPointRef tet(const polyMesh &mesh, const pointField &meshPoints, const pointField &cellCentres) const
Return the geometry corresponding to this tet and the given.
Definition: tetIndicesI.H:109
scalar mag() const
Return volume.
Definition: tetrahedronI.H:170
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:71
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label nPoints
tUEqn clear()
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar c
Speed of light in a vacuum.
barycentric coordinates(const polyMesh &mesh, const point &position, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the coordinates given the position and tet topology.
Definition: tracking.C:1259
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:45
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
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:64
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
face triFace(3)
labelList f(nPoints)