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-2024 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 (non-blocking)
111  UPstream::waitRequests(startOfRequests);
112 
113  forAll(pfbf, patchi)
114  {
115  if (pfbf[patchi].coupled())
116  {
117  refCast<coupledPointPatchField<Type>>
118  (pfbf[patchi]).swapAddSeparated
119  (
121  pfi
122  );
123  }
124  }
125 }
126 
127 
128 template<class Type>
130 (
131  const GeometricField<Type, fvPatchField, volMesh>& vf,
132  GeometricField<Type, pointPatchField, pointMesh>& pf
133 ) const
134 {
135  if (debug)
136  {
137  Pout<< "volPointInterpolation::interpolateInternalField("
138  << "const GeometricField<Type, fvPatchField, volMesh>&, "
139  << "GeometricField<Type, pointPatchField, pointMesh>&) : "
140  << "interpolating field " << vf.name()
141  << " from cells to points " << pf.name() << endl;
142  }
143 
144  const labelListList& pointCells = vf.mesh().pointCells();
145 
146  // Multiply volField by weighting factor matrix to create pointField
147  forAll(pointCells, pointi)
148  {
149  if (!isPatchPoint_[pointi])
150  {
151  const scalarList& pw = pointWeights_[pointi];
152  const labelList& ppc = pointCells[pointi];
153 
154  pf[pointi] = Zero;
155 
156  forAll(ppc, pointCelli)
157  {
158  pf[pointi] += pw[pointCelli]*vf[ppc[pointCelli]];
159  }
160  }
161  }
162 }
163 
164 
165 template<class Type>
167 (
168  const DimensionedField<Type, volMesh>& vf,
169  DimensionedField<Type, pointMesh>& pf
170 ) const
171 {
172  if (debug)
173  {
174  Pout<< "volPointInterpolation::interpolateDimensionedInternalField("
175  << "const DimensionedField<Type, volMesh>&, "
176  << "DimensionedField<Type, pointMesh>&) : "
177  << "interpolating field " << vf.name() << " from cells to points "
178  << pf.name() << endl;
179  }
180 
181  const fvMesh& mesh = vf.mesh();
182 
183  const labelListList& pointCells = mesh.pointCells();
184  const pointField& points = mesh.points();
185  const vectorField& cellCentres = mesh.cellCentres();
186 
187  // Re-do weights and interpolation since normal interpolation
188  // pointWeights_ are for non-boundary points only. Not efficient but
189  // then saves on space.
190 
191  // Multiply volField by weighting factor matrix to create pointField
192  scalarField sumW(points.size(), Zero);
193  forAll(pointCells, pointi)
194  {
195  const labelList& ppc = pointCells[pointi];
196 
197  pf[pointi] = Type(Zero);
198 
199  forAll(ppc, pointCelli)
200  {
201  label celli = ppc[pointCelli];
202  scalar pw = 1.0/mag(points[pointi] - cellCentres[celli]);
203 
204  pf[pointi] += pw*vf[celli];
205  sumW[pointi] += pw;
206  }
207  }
208 
209  // Sum collocated contributions
210  pointConstraints::syncUntransformedData(mesh, sumW, plusEqOp<scalar>());
211  pointConstraints::syncUntransformedData(mesh, pf, plusEqOp<Type>());
212 
213  // Normalise
214  forAll(pf, pointi)
215  {
216  scalar s = sumW[pointi];
217  if (s > ROOTVSMALL)
218  {
219  pf[pointi] /= s;
220  }
221  }
222 }
223 
224 
225 template<class Type>
227 Foam::volPointInterpolation::flatBoundaryField
228 (
229  const GeometricField<Type, fvPatchField, volMesh>& vf
230 ) const
231 {
232  const polyBoundaryMesh& pbm = vf.mesh().boundaryMesh();
233 
234  auto tboundaryVals = tmp<Field<Type>>::New(pbm.nFaces(), Foam::zero{});
235  auto& values = tboundaryVals.ref();
236 
237  forAll(vf.boundaryField(), patchi)
238  {
239  const auto& pp = pbm[patchi];
240  const auto& pfld = vf.boundaryField()[patchi];
241 
242  // Note: restrict transcribing to actual size of the patch field
243  // - handles "empty" patch type etc.
244 
245  SubList<Type> slice(values, pfld.size(), pp.offset());
246 
247  if
248  (
249  !isA<emptyPolyPatch>(pp)
250  && !pfld.coupled()
251  )
252  {
253  slice = pfld;
254  }
255  }
257  return tboundaryVals;
258 }
259 
260 
261 template<class Type>
263 (
266 ) const
267 {
268  const primitivePatch& boundary = boundaryPtr_();
269 
270  Field<Type>& pfi = pf.primitiveFieldRef();
271 
272  // Get face data in flat list
273  tmp<Field<Type>> tboundaryVals(flatBoundaryField(vf));
274  const Field<Type>& boundaryVals = tboundaryVals();
275 
276 
277  // Do points on 'normal' patches from the surrounding patch faces
278  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
279 
280  const labelList& mp = boundary.meshPoints();
281 
282  forAll(mp, i)
283  {
284  label pointi = mp[i];
285 
286  if (isPatchPoint_[pointi])
287  {
288  const labelList& pFaces = boundary.pointFaces()[i];
289  const scalarList& pWeights = boundaryPointWeights_[i];
290 
291  Type& val = pfi[pointi];
292 
293  val = Zero;
294  forAll(pFaces, j)
295  {
296  if (boundaryIsPatchFace_[pFaces[j]])
297  {
298  val += pWeights[j]*boundaryVals[pFaces[j]];
299  }
300  }
301  }
302  }
303 
304  // Sum collocated contributions
305  pointConstraints::syncUntransformedData(mesh(), pfi, plusEqOp<Type>());
306 
307  // And add separated contributions
308  addSeparated(pf);
309 
310  // Optionally normalise
311  if (normalisationPtr_)
312  {
313  const scalarField& normalisation = normalisationPtr_();
314  forAll(mp, i)
315  {
316  pfi[mp[i]] *= normalisation[i];
317  }
318  }
319 
320 
321  // Push master data to slaves. It is possible (not sure how often) for
322  // a coupled point to have its master on a different patch so
323  // to make sure just push master data to slaves.
324  pushUntransformedData(pfi);
325 }
326 
327 
328 template<class Type>
330 (
333  const bool overrideFixedValue
334 ) const
335 {
336  interpolateBoundaryField(vf, pf);
337 
338  // Apply constraints
339  const pointConstraints& pcs = pointConstraints::New(pf.mesh());
341  pcs.constrain(pf, overrideFixedValue);
342 }
343 
344 
345 template<class Type>
347 (
350 ) const
351 {
352  if (debug)
353  {
354  Pout<< "volPointInterpolation::interpolate("
355  << "const GeometricField<Type, fvPatchField, volMesh>&, "
356  << "GeometricField<Type, pointPatchField, pointMesh>&) : "
357  << "interpolating field " << vf.name() << " from cells to points "
358  << pf.name() << endl;
359  }
360 
361  interpolateInternalField(vf, pf);
362 
363  // Interpolate to the patches preserving fixed value BCs
364  interpolateBoundaryField(vf, pf, false);
365 }
366 
367 
368 template<class Type>
371 (
373  const wordList& patchFieldTypes
374 ) const
375 {
376  const pointMesh& pm = pointMesh::New(vf.mesh());
377 
378  // Construct tmp<pointField>
380  (
381  IOobject
382  (
383  "volPointInterpolate(" + vf.name() + ')',
384  vf.instance(),
385  pm.thisDb()
386  ),
387  pm,
388  vf.dimensions(),
389  patchFieldTypes
390  );
391 
392  interpolateInternalField(vf, tpf.ref());
393 
394  // Interpolate to the patches overriding fixed value BCs
395  interpolateBoundaryField(vf, tpf.ref(), true);
396 
397  return tpf;
398 }
399 
400 
401 template<class Type>
404 (
406  const wordList& patchFieldTypes
407 ) const
408 {
409  // Construct tmp<pointField>
411  interpolate(tvf(), patchFieldTypes);
412  tvf.clear();
413  return tpf;
414 }
415 
416 
417 template<class Type>
420 (
422  const word& name,
423  const bool cache
424 ) const
425 {
427 
428  const pointMesh& pm = pointMesh::New(vf.mesh());
429  const objectRegistry& db = pm.thisDb();
430 
431  PointFieldType* pfPtr =
432  db.objectRegistry::template getObjectPtr<PointFieldType>(name);
433 
434  if (!cache || vf.mesh().changing())
435  {
436  // Delete any old occurrences to avoid double registration
437  if (pfPtr && pfPtr->ownedByRegistry())
438  {
439  solution::cachePrintMessage("Deleting", name, vf);
440  delete pfPtr;
441  }
442 
444  (
445  IOobject
446  (
447  name,
448  vf.instance(),
449  pm.thisDb()
450  ),
451  pm,
452  vf.dimensions()
453  );
454 
455  interpolate(vf, tpf.ref());
456 
457  return tpf;
458  }
459 
460 
461  if (!pfPtr)
462  {
463  solution::cachePrintMessage("Calculating and caching", name, vf);
464 
465  pfPtr = interpolate(vf, name, false).ptr();
466  regIOobject::store(pfPtr);
467  }
468  else
469  {
470  PointFieldType& pf = *pfPtr;
471 
472  if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
473  {
474  solution::cachePrintMessage("Reusing", name, vf);
475  }
476  else
477  {
478  solution::cachePrintMessage("Updating", name, vf);
479  interpolate(vf, pf);
480  }
481  }
482 
483  return *pfPtr;
484 }
485 
486 
487 template<class Type>
490 (
492 ) const
493 {
494  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
495 }
496 
497 
498 template<class Type>
501 (
503 ) const
504 {
505  // Construct tmp<pointField>
507  interpolate(tvf());
508  tvf.clear();
509  return tpf;
510 }
511 
512 
513 template<class Type>
516 (
518  const word& name,
519  const bool cache
520 ) const
521 {
522  typedef DimensionedField<Type, pointMesh> PointFieldType;
523 
524  const pointMesh& pm = pointMesh::New(vf.mesh());
525  const objectRegistry& db = pm.thisDb();
526 
527 
528  PointFieldType* pfPtr =
529  db.objectRegistry::template getObjectPtr<PointFieldType>(name);
530 
531  if (!cache || vf.mesh().changing())
532  {
533  // Delete any old occurrences to avoid double registration
534  if (pfPtr && pfPtr->ownedByRegistry())
535  {
536  solution::cachePrintMessage("Deleting", name, vf);
537  delete pfPtr;
538  }
539 
541  (
542  IOobject
543  (
544  name,
545  vf.instance(),
546  pm.thisDb()
547  ),
548  pm,
549  vf.dimensions()
550  );
551 
552  interpolateDimensionedInternalField(vf, tpf.ref());
553 
554  return tpf;
555  }
556 
557 
558  if (!pfPtr)
559  {
560  solution::cachePrintMessage("Calculating and caching", name, vf);
561 
562  pfPtr = interpolate(vf, name, false).ptr();
563  regIOobject::store(pfPtr);
564  }
565  else
566  {
567  PointFieldType& pf = *pfPtr;
568 
569  if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
570  {
571  solution::cachePrintMessage("Reusing", name, vf);
572  }
573  else
574  {
575  solution::cachePrintMessage("Updating", name, vf);
576  interpolateDimensionedInternalField(vf, pf);
577  }
578  }
579 
580  return *pfPtr;
581 }
582 
583 
584 template<class Type>
587 (
589 ) const
590 {
591  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
592 }
593 
594 
595 template<class Type>
598 (
600 ) const
601 {
602  // Construct tmp<pointField>
604  tvf.clear();
605  return tpf;
606 }
607 
608 
609 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
faceListList boundary
const polyBoundaryMesh & pbm
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 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.
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:134
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1561
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.
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
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
const polyMesh & mesh() const noexcept
Return the mesh reference.
Generic templated field type.
Definition: Field.H:62
const fvMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:255
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
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:164
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.
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
label nFaces() const noexcept
The number of boundary faces in the underlying mesh.
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...
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
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))
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
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