volPointInterpolate.C
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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2020 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 "volPointInterpolation.H"
30 #include "volFields.H"
31 #include "pointFields.H"
32 #include "emptyFvPatch.H"
33 #include "coupledPointPatchField.H"
34 #include "pointConstraints.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 template<class Type>
39 void Foam::volPointInterpolation::pushUntransformedData
40 (
41  List<Type>& pointData
42 ) const
43 {
44  // Transfer onto coupled patch
45  const globalMeshData& gmd = mesh().globalData();
46  const indirectPrimitivePatch& cpp = gmd.coupledPatch();
47  const labelList& meshPoints = cpp.meshPoints();
48 
49  const mapDistribute& slavesMap = gmd.globalCoPointSlavesMap();
50  const labelListList& slaves = gmd.globalCoPointSlaves();
51 
52  List<Type> elems(slavesMap.constructSize());
53  forAll(meshPoints, i)
54  {
55  elems[i] = pointData[meshPoints[i]];
56  }
57 
58  // Combine master data with slave data
59  forAll(slaves, i)
60  {
61  const labelList& slavePoints = slaves[i];
62 
63  // Copy master data to slave slots
64  forAll(slavePoints, j)
65  {
66  elems[slavePoints[j]] = elems[i];
67  }
68  }
69 
70  // Push slave-slot data back to slaves
71  slavesMap.reverseDistribute(elems.size(), elems, false);
72 
73  // Extract back onto mesh
74  forAll(meshPoints, i)
75  {
76  pointData[meshPoints[i]] = elems[i];
77  }
78 }
79 
80 
81 template<class Type>
82 void Foam::volPointInterpolation::addSeparated
83 (
84  GeometricField<Type, pointPatchField, pointMesh>& pf
85 ) const
86 {
87  if (debug)
88  {
89  Pout<< "volPointInterpolation::addSeparated" << endl;
90  }
91 
92  auto& pfi = pf.internalFieldRef();
93  auto& pfbf = pf.boundaryFieldRef();
94 
95  const label startOfRequests = UPstream::nRequests();
96 
97  forAll(pfbf, patchi)
98  {
99  if (pfbf[patchi].coupled())
100  {
101  refCast<coupledPointPatchField<Type>>
102  (pfbf[patchi]).initSwapAddSeparated
103  (
105  pfi
106  );
107  }
108  }
109 
110  // Wait for outstanding requests
111  // (commsType == UPstream::commsTypes::nonBlocking)
112  UPstream::waitRequests(startOfRequests);
113 
114  forAll(pfbf, patchi)
115  {
116  if (pfbf[patchi].coupled())
117  {
118  refCast<coupledPointPatchField<Type>>
119  (pfbf[patchi]).swapAddSeparated
120  (
122  pfi
123  );
124  }
125  }
126 }
127 
128 
129 template<class Type>
131 (
132  const GeometricField<Type, fvPatchField, volMesh>& vf,
133  GeometricField<Type, pointPatchField, pointMesh>& pf
134 ) const
135 {
136  if (debug)
137  {
138  Pout<< "volPointInterpolation::interpolateInternalField("
139  << "const GeometricField<Type, fvPatchField, volMesh>&, "
140  << "GeometricField<Type, pointPatchField, pointMesh>&) : "
141  << "interpolating field " << vf.name()
142  << " from cells to points " << pf.name() << endl;
143  }
144 
145  const labelListList& pointCells = vf.mesh().pointCells();
146 
147  // Multiply volField by weighting factor matrix to create pointField
148  forAll(pointCells, pointi)
149  {
150  if (!isPatchPoint_[pointi])
151  {
152  const scalarList& pw = pointWeights_[pointi];
153  const labelList& ppc = pointCells[pointi];
154 
155  pf[pointi] = Zero;
156 
157  forAll(ppc, pointCelli)
158  {
159  pf[pointi] += pw[pointCelli]*vf[ppc[pointCelli]];
160  }
161  }
162  }
163 }
164 
165 
166 template<class Type>
168 (
169  const DimensionedField<Type, volMesh>& vf,
170  DimensionedField<Type, pointMesh>& pf
171 ) const
172 {
173  if (debug)
174  {
175  Pout<< "volPointInterpolation::interpolateDimensionedInternalField("
176  << "const DimensionedField<Type, volMesh>&, "
177  << "DimensionedField<Type, pointMesh>&) : "
178  << "interpolating field " << vf.name() << " from cells to points "
179  << pf.name() << endl;
180  }
181 
182  const fvMesh& mesh = vf.mesh();
183 
184  const labelListList& pointCells = mesh.pointCells();
185  const pointField& points = mesh.points();
186  const vectorField& cellCentres = mesh.cellCentres();
187 
188  // Re-do weights and interpolation since normal interpolation
189  // pointWeights_ are for non-boundary points only. Not efficient but
190  // then saves on space.
191 
192  // Multiply volField by weighting factor matrix to create pointField
193  scalarField sumW(points.size(), Zero);
194  forAll(pointCells, pointi)
195  {
196  const labelList& ppc = pointCells[pointi];
197 
198  pf[pointi] = Type(Zero);
199 
200  forAll(ppc, pointCelli)
201  {
202  label celli = ppc[pointCelli];
203  scalar pw = 1.0/mag(points[pointi] - cellCentres[celli]);
204 
205  pf[pointi] += pw*vf[celli];
206  sumW[pointi] += pw;
207  }
208  }
209 
210  // Sum collocated contributions
211  pointConstraints::syncUntransformedData(mesh, sumW, plusEqOp<scalar>());
212  pointConstraints::syncUntransformedData(mesh, pf, plusEqOp<Type>());
213 
214  // Normalise
215  forAll(pf, pointi)
216  {
217  scalar s = sumW[pointi];
218  if (s > ROOTVSMALL)
219  {
220  pf[pointi] /= s;
221  }
222  }
223 }
224 
225 
226 template<class Type>
227 Foam::tmp<Foam::Field<Type>> Foam::volPointInterpolation::flatBoundaryField
228 (
229  const GeometricField<Type, fvPatchField, volMesh>& vf
230 ) const
231 {
232  const fvMesh& mesh = vf.mesh();
233  const fvBoundaryMesh& bm = mesh.boundary();
234 
235  tmp<Field<Type>> tboundaryVals
236  (
237  new Field<Type>(mesh.nBoundaryFaces())
238  );
239  Field<Type>& boundaryVals = tboundaryVals.ref();
240 
241  forAll(vf.boundaryField(), patchi)
242  {
243  label bFacei = bm[patchi].patch().start() - mesh.nInternalFaces();
244 
245  if
246  (
247  !isA<emptyFvPatch>(bm[patchi])
248  && !vf.boundaryField()[patchi].coupled()
249  )
250  {
251  SubList<Type>
252  (
253  boundaryVals,
254  vf.boundaryField()[patchi].size(),
255  bFacei
256  ) = vf.boundaryField()[patchi];
257  }
258  else
259  {
260  const polyPatch& pp = bm[patchi].patch();
261 
262  forAll(pp, i)
263  {
264  boundaryVals[bFacei++] = Zero;
265  }
266  }
267  }
269  return tboundaryVals;
270 }
271 
272 
273 template<class Type>
275 (
278 ) const
279 {
280  const primitivePatch& boundary = boundaryPtr_();
281 
282  Field<Type>& pfi = pf.primitiveFieldRef();
283 
284  // Get face data in flat list
285  tmp<Field<Type>> tboundaryVals(flatBoundaryField(vf));
286  const Field<Type>& boundaryVals = tboundaryVals();
287 
288 
289  // Do points on 'normal' patches from the surrounding patch faces
290  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
291 
292  const labelList& mp = boundary.meshPoints();
293 
294  forAll(mp, i)
295  {
296  label pointi = mp[i];
297 
298  if (isPatchPoint_[pointi])
299  {
300  const labelList& pFaces = boundary.pointFaces()[i];
301  const scalarList& pWeights = boundaryPointWeights_[i];
302 
303  Type& val = pfi[pointi];
304 
305  val = Zero;
306  forAll(pFaces, j)
307  {
308  if (boundaryIsPatchFace_[pFaces[j]])
309  {
310  val += pWeights[j]*boundaryVals[pFaces[j]];
311  }
312  }
313  }
314  }
315 
316  // Sum collocated contributions
317  pointConstraints::syncUntransformedData(mesh(), pfi, plusEqOp<Type>());
318 
319  // And add separated contributions
320  addSeparated(pf);
321 
322  // Optionally normalise
323  if (normalisationPtr_)
324  {
325  const scalarField& normalisation = normalisationPtr_();
326  forAll(mp, i)
327  {
328  pfi[mp[i]] *= normalisation[i];
329  }
330  }
331 
332 
333  // Push master data to slaves. It is possible (not sure how often) for
334  // a coupled point to have its master on a different patch so
335  // to make sure just push master data to slaves.
336  pushUntransformedData(pfi);
337 }
338 
339 
340 template<class Type>
342 (
345  const bool overrideFixedValue
346 ) const
347 {
348  interpolateBoundaryField(vf, pf);
349 
350  // Apply constraints
351  const pointConstraints& pcs = pointConstraints::New(pf.mesh());
353  pcs.constrain(pf, overrideFixedValue);
354 }
355 
356 
357 template<class Type>
359 (
362 ) const
363 {
364  if (debug)
365  {
366  Pout<< "volPointInterpolation::interpolate("
367  << "const GeometricField<Type, fvPatchField, volMesh>&, "
368  << "GeometricField<Type, pointPatchField, pointMesh>&) : "
369  << "interpolating field " << vf.name() << " from cells to points "
370  << pf.name() << endl;
371  }
372 
373  interpolateInternalField(vf, pf);
374 
375  // Interpolate to the patches preserving fixed value BCs
376  interpolateBoundaryField(vf, pf, false);
377 }
378 
379 
380 template<class Type>
383 (
385  const wordList& patchFieldTypes
386 ) const
387 {
388  const pointMesh& pm = pointMesh::New(vf.mesh());
389 
390  // Construct tmp<pointField>
392  (
393  IOobject
394  (
395  "volPointInterpolate(" + vf.name() + ')',
396  vf.instance(),
397  pm.thisDb()
398  ),
399  pm,
400  vf.dimensions(),
401  patchFieldTypes
402  );
403 
404  interpolateInternalField(vf, tpf.ref());
405 
406  // Interpolate to the patches overriding fixed value BCs
407  interpolateBoundaryField(vf, tpf.ref(), true);
408 
409  return tpf;
410 }
411 
412 
413 template<class Type>
416 (
418  const wordList& patchFieldTypes
419 ) const
420 {
421  // Construct tmp<pointField>
423  interpolate(tvf(), patchFieldTypes);
424  tvf.clear();
425  return tpf;
426 }
427 
428 
429 template<class Type>
432 (
434  const word& name,
435  const bool cache
436 ) const
437 {
439 
440  const pointMesh& pm = pointMesh::New(vf.mesh());
441  const objectRegistry& db = pm.thisDb();
442 
443  PointFieldType* pfPtr =
444  db.objectRegistry::template getObjectPtr<PointFieldType>(name);
445 
446  if (!cache || vf.mesh().changing())
447  {
448  // Delete any old occurrences to avoid double registration
449  if (pfPtr && pfPtr->ownedByRegistry())
450  {
451  solution::cachePrintMessage("Deleting", name, vf);
452  delete pfPtr;
453  }
454 
456  (
457  IOobject
458  (
459  name,
460  vf.instance(),
461  pm.thisDb()
462  ),
463  pm,
464  vf.dimensions()
465  );
466 
467  interpolate(vf, tpf.ref());
468 
469  return tpf;
470  }
471 
472 
473  if (!pfPtr)
474  {
475  solution::cachePrintMessage("Calculating and caching", name, vf);
476 
477  pfPtr = interpolate(vf, name, false).ptr();
478  regIOobject::store(pfPtr);
479  }
480  else
481  {
482  PointFieldType& pf = *pfPtr;
483 
484  if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
485  {
486  solution::cachePrintMessage("Reusing", name, vf);
487  }
488  else
489  {
490  solution::cachePrintMessage("Updating", name, vf);
491  interpolate(vf, pf);
492  }
493  }
494 
495  return *pfPtr;
496 }
497 
498 
499 template<class Type>
502 (
504 ) const
505 {
506  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
507 }
508 
509 
510 template<class Type>
513 (
515 ) const
516 {
517  // Construct tmp<pointField>
519  interpolate(tvf());
520  tvf.clear();
521  return tpf;
522 }
523 
524 
525 template<class Type>
528 (
530  const word& name,
531  const bool cache
532 ) const
533 {
534  typedef DimensionedField<Type, pointMesh> PointFieldType;
535 
536  const pointMesh& pm = pointMesh::New(vf.mesh());
537  const objectRegistry& db = pm.thisDb();
538 
539 
540  PointFieldType* pfPtr =
541  db.objectRegistry::template getObjectPtr<PointFieldType>(name);
542 
543  if (!cache || vf.mesh().changing())
544  {
545  // Delete any old occurrences to avoid double registration
546  if (pfPtr && pfPtr->ownedByRegistry())
547  {
548  solution::cachePrintMessage("Deleting", name, vf);
549  delete pfPtr;
550  }
551 
553  (
554  IOobject
555  (
556  name,
557  vf.instance(),
558  pm.thisDb()
559  ),
560  pm,
561  vf.dimensions()
562  );
563 
564  interpolateDimensionedInternalField(vf, tpf.ref());
565 
566  return tpf;
567  }
568 
569 
570  if (!pfPtr)
571  {
572  solution::cachePrintMessage("Calculating and caching", name, vf);
573  pfPtr = interpolate(vf, name, false).ptr();
574 
575  regIOobject::store(pfPtr);
576  }
577  else
578  {
579  PointFieldType& pf = *pfPtr;
580 
581  if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
582  {
583  solution::cachePrintMessage("Reusing", name, vf);
584  }
585  else
586  {
587  solution::cachePrintMessage("Updating", name, vf);
588  interpolateDimensionedInternalField(vf, pf);
589  }
590  }
591 
592  return *pfPtr;
593 }
594 
595 
596 template<class Type>
599 (
601 ) const
602 {
603  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
604 }
605 
606 
607 template<class Type>
610 (
612 ) const
613 {
614  // Construct tmp<pointField>
616  tvf.clear();
617  return tpf;
618 }
619 
620 
621 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
faceListList boundary
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
static const pointConstraints & New(const pointMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1538
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:45
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Application of (multi-)patch point constraints.
static void cachePrintMessage(const char *message, const word &name, const FieldType &fld)
Helper for printing cache message.
A list of faces which address into the list of points.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
Generic templated field type.
Definition: Field.H:62
const fvMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:157
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:154
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
void interpolateInternalField(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate internal field from volField to pointField.
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
const labelListList & pointCells() const
const Mesh & mesh() const noexcept
Return mesh.
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:192
int debug
Static debugging option.
void interpolateDimensionedInternalField(const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const
Interpolate dimensioned internal field from cells to points.
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:266
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
"nonBlocking" : (MPI_Isend, MPI_Irecv)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Registry of regIOobjects.
bool coupled
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))
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const dimensionedScalar mp
Proton mass.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127