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 | Copyright (C) 2011-2015 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 "volPointInterpolation.H"
27 #include "volFields.H"
28 #include "pointFields.H"
29 #include "emptyFvPatch.H"
30 #include "coupledPointPatchField.H"
31 #include "pointConstraints.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 template<class Type>
42 (
44 ) const
45 {
46  // Transfer onto coupled patch
47  const globalMeshData& gmd = mesh().globalData();
48  const indirectPrimitivePatch& cpp = gmd.coupledPatch();
49  const labelList& meshPoints = cpp.meshPoints();
50 
51  const mapDistribute& slavesMap = gmd.globalCoPointSlavesMap();
52  const labelListList& slaves = gmd.globalCoPointSlaves();
53 
54  List<Type> elems(slavesMap.constructSize());
55  forAll(meshPoints, i)
56  {
57  elems[i] = pointData[meshPoints[i]];
58  }
59 
60  // Combine master data with slave data
61  forAll(slaves, i)
62  {
63  const labelList& slavePoints = slaves[i];
64 
65  // Copy master data to slave slots
66  forAll(slavePoints, j)
67  {
68  elems[slavePoints[j]] = elems[i];
69  }
70  }
71 
72  // Push slave-slot data back to slaves
73  slavesMap.reverseDistribute(elems.size(), elems, false);
74 
75  // Extract back onto mesh
76  forAll(meshPoints, i)
77  {
78  pointData[meshPoints[i]] = elems[i];
79  }
80 }
81 
82 
83 template<class Type>
85 (
87 ) const
88 {
89  if (debug)
90  {
91  Pout<< "volPointInterpolation::addSeparated" << endl;
92  }
93 
94  forAll(pf.boundaryField(), patchI)
95  {
96  if (pf.boundaryField()[patchI].coupled())
97  {
98  refCast<coupledPointPatchField<Type> >
99  (pf.boundaryField()[patchI]).initSwapAddSeparated
100  (
102  pf.internalField()
103  );
104  }
105  }
106 
107  // Block for any outstanding requests
109 
110  forAll(pf.boundaryField(), patchI)
111  {
112  if (pf.boundaryField()[patchI].coupled())
113  {
114  refCast<coupledPointPatchField<Type> >
115  (pf.boundaryField()[patchI]).swapAddSeparated
116  (
118  pf.internalField()
119  );
120  }
121  }
122 }
123 
124 
125 template<class Type>
127 (
130 ) const
131 {
132  if (debug)
133  {
134  Pout<< "volPointInterpolation::interpolateInternalField("
135  << "const GeometricField<Type, fvPatchField, volMesh>&, "
136  << "GeometricField<Type, pointPatchField, pointMesh>&) : "
137  << "interpolating field from cells to points"
138  << endl;
139  }
140 
141  const labelListList& pointCells = vf.mesh().pointCells();
142 
143  // Multiply volField by weighting factor matrix to create pointField
144  forAll(pointCells, pointi)
145  {
146  if (!isPatchPoint_[pointi])
147  {
148  const scalarList& pw = pointWeights_[pointi];
149  const labelList& ppc = pointCells[pointi];
150 
151  pf[pointi] = pTraits<Type>::zero;
152 
153  forAll(ppc, pointCelli)
154  {
155  pf[pointi] += pw[pointCelli]*vf[ppc[pointCelli]];
156  }
157  }
158  }
159 }
160 
161 
162 template<class Type>
164 (
167 ) const
168 {
169  if (debug)
170  {
171  Pout<< "volPointInterpolation::interpolateDimensionedInternalField("
172  << "const DimensionedField<Type, volMesh>&, "
173  << "DimensionedField<Type, pointMesh>&) : "
174  << "interpolating field from cells to points"
175  << endl;
176  }
177 
178  const fvMesh& mesh = vf.mesh();
179 
180  const labelListList& pointCells = mesh.pointCells();
181  const pointField& points = mesh.points();
182  const vectorField& cellCentres = mesh.cellCentres();
183 
184  // Re-do weights and interpolation since normal interpolation
185  // pointWeights_ are for non-boundary points only. Not efficient but
186  // then saves on space.
187 
188  // Multiply volField by weighting factor matrix to create pointField
189  scalarField sumW(points.size(), 0.0);
190  forAll(pointCells, pointi)
191  {
192  const labelList& ppc = pointCells[pointi];
193 
194  pf[pointi] = pTraits<Type>::zero;
195 
196  forAll(ppc, pointCelli)
197  {
198  label celli = ppc[pointCelli];
199  scalar pw = 1.0/mag(points[pointi] - cellCentres[celli]);
200 
201  pf[pointi] += pw*vf[celli];
202  sumW[pointi] += pw;
203  }
204  }
205 
206  // Sum collocated contributions
209 
210  // Normalise
211  forAll(pf, pointi)
212  {
213  scalar s = sumW[pointi];
214  if (s > ROOTVSMALL)
215  {
216  pf[pointi] /= s;
217  }
218  }
219 }
220 
221 
222 template<class Type>
224 (
226 ) const
227 {
228  const fvMesh& mesh = vf.mesh();
229  const fvBoundaryMesh& bm = mesh.boundary();
230 
231  tmp<Field<Type> > tboundaryVals
232  (
234  );
235  Field<Type>& boundaryVals = tboundaryVals();
236 
237  forAll(vf.boundaryField(), patchI)
238  {
239  label bFaceI = bm[patchI].patch().start() - mesh.nInternalFaces();
240 
241  if
242  (
243  !isA<emptyFvPatch>(bm[patchI])
244  && !vf.boundaryField()[patchI].coupled()
245  )
246  {
248  (
249  boundaryVals,
250  vf.boundaryField()[patchI].size(),
251  bFaceI
252  ).assign(vf.boundaryField()[patchI]);
253  }
254  else
255  {
256  const polyPatch& pp = bm[patchI].patch();
257 
258  forAll(pp, i)
259  {
260  boundaryVals[bFaceI++] = pTraits<Type>::zero;
261  }
262  }
263  }
264 
265  return tboundaryVals;
266 }
267 
268 
269 template<class Type>
271 (
274 ) const
275 {
276  const primitivePatch& boundary = boundaryPtr_();
277 
278  Field<Type>& pfi = pf.internalField();
279 
280  // Get face data in flat list
281  tmp<Field<Type> > tboundaryVals(flatBoundaryField(vf));
282  const Field<Type>& boundaryVals = tboundaryVals();
283 
284 
285  // Do points on 'normal' patches from the surrounding patch faces
286  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
287 
288  forAll(boundary.meshPoints(), i)
289  {
290  label pointI = boundary.meshPoints()[i];
291 
292  if (isPatchPoint_[pointI])
293  {
294  const labelList& pFaces = boundary.pointFaces()[i];
295  const scalarList& pWeights = boundaryPointWeights_[i];
296 
297  Type& val = pfi[pointI];
298 
299  val = pTraits<Type>::zero;
300  forAll(pFaces, j)
301  {
302  if (boundaryIsPatchFace_[pFaces[j]])
303  {
304  val += pWeights[j]*boundaryVals[pFaces[j]];
305  }
306  }
307  }
308  }
309 
310  // Sum collocated contributions
312 
313  // And add separated contributions
314  addSeparated(pf);
315 
316  // Push master data to slaves. It is possible (not sure how often) for
317  // a coupled point to have its master on a different patch so
318  // to make sure just push master data to slaves.
319  pushUntransformedData(pfi);
320 }
321 
322 
323 template<class Type>
325 (
328  const bool overrideFixedValue
329 ) const
330 {
331  interpolateBoundaryField(vf, pf);
332 
333  // Apply constraints
334  const pointConstraints& pcs = pointConstraints::New(pf.mesh());
335 
336  pcs.constrain(pf, overrideFixedValue);
337 }
338 
339 
340 template<class Type>
342 (
345 ) const
346 {
347  if (debug)
348  {
349  Pout<< "volPointInterpolation::interpolate("
350  << "const GeometricField<Type, fvPatchField, volMesh>&, "
351  << "GeometricField<Type, pointPatchField, pointMesh>&) : "
352  << "interpolating field from cells to points"
353  << endl;
354  }
355 
356  interpolateInternalField(vf, pf);
357 
358  // Interpolate to the patches preserving fixed value BCs
359  interpolateBoundaryField(vf, pf, false);
360 }
361 
362 
363 template<class Type>
366 (
368  const wordList& patchFieldTypes
369 ) const
370 {
371  const pointMesh& pm = pointMesh::New(vf.mesh());
372 
373  // Construct tmp<pointField>
375  (
377  (
378  IOobject
379  (
380  "volPointInterpolate(" + vf.name() + ')',
381  vf.instance(),
382  pm.thisDb()
383  ),
384  pm,
385  vf.dimensions(),
386  patchFieldTypes
387  )
388  );
389 
390  interpolateInternalField(vf, tpf());
391 
392  // Interpolate to the patches overriding fixed value BCs
393  interpolateBoundaryField(vf, tpf(), true);
394 
395  return tpf;
396 }
397 
398 
399 template<class Type>
402 (
404  const wordList& patchFieldTypes
405 ) const
406 {
407  // Construct tmp<pointField>
409  interpolate(tvf(), patchFieldTypes);
410  tvf.clear();
411  return tpf;
412 }
413 
414 
415 template<class Type>
418 (
420  const word& name,
421  const bool cache
422 ) const
423 {
425 
426  const pointMesh& pm = pointMesh::New(vf.mesh());
427  const objectRegistry& db = pm.thisDb();
428 
429  if (!cache || vf.mesh().changing())
430  {
431  // Delete any old occurences to avoid double registration
432  if (db.objectRegistry::template foundObject<PointFieldType>(name))
433  {
434  PointFieldType& pf = const_cast<PointFieldType&>
435  (
436  db.objectRegistry::template lookupObject<PointFieldType>(name)
437  );
438 
439  if (pf.ownedByRegistry())
440  {
441  solution::cachePrintMessage("Deleting", name, vf);
442  pf.release();
443  delete &pf;
444  }
445  }
446 
447 
449  (
451  (
452  IOobject
453  (
454  name,
455  vf.instance(),
456  pm.thisDb()
457  ),
458  pm,
459  vf.dimensions()
460  )
461  );
462 
463  interpolate(vf, tpf());
464 
465  return tpf;
466  }
467  else
468  {
469  if (!db.objectRegistry::template foundObject<PointFieldType>(name))
470  {
471  solution::cachePrintMessage("Calculating and caching", name, vf);
472  tmp<PointFieldType> tpf = interpolate(vf, name, false);
473  PointFieldType* pfPtr = tpf.ptr();
474  regIOobject::store(pfPtr);
475  return *pfPtr;
476  }
477  else
478  {
479  PointFieldType& pf = const_cast<PointFieldType&>
480  (
481  db.objectRegistry::template lookupObject<PointFieldType>(name)
482  );
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  return pf;
494  }
495  }
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  if (!cache || vf.mesh().changing())
540  {
541  // Delete any old occurences to avoid double registration
542  if (db.objectRegistry::template foundObject<PointFieldType>(name))
543  {
544  PointFieldType& pf = const_cast<PointFieldType&>
545  (
546  db.objectRegistry::template lookupObject<PointFieldType>(name)
547  );
548 
549  if (pf.ownedByRegistry())
550  {
551  solution::cachePrintMessage("Deleting", name, vf);
552  pf.release();
553  delete &pf;
554  }
555  }
556 
557 
559  (
561  (
562  IOobject
563  (
564  name,
565  vf.instance(),
566  pm.thisDb()
567  ),
568  pm,
569  vf.dimensions()
570  )
571  );
572 
573  interpolateDimensionedInternalField(vf, tpf());
574 
575  return tpf;
576  }
577  else
578  {
579  if (!db.objectRegistry::template foundObject<PointFieldType>(name))
580  {
581  solution::cachePrintMessage("Calculating and caching", name, vf);
582  tmp<PointFieldType> tpf = interpolate(vf, name, false);
583  PointFieldType* pfPtr = tpf.ptr();
584  regIOobject::store(pfPtr);
585  return *pfPtr;
586  }
587  else
588  {
589  PointFieldType& pf = const_cast<PointFieldType&>
590  (
591  db.objectRegistry::template lookupObject<PointFieldType>(name)
592  );
593 
594  if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
595  {
596  solution::cachePrintMessage("Reusing", name, vf);
597  }
598  else
599  {
600  solution::cachePrintMessage("Updating", name, vf);
601  interpolateDimensionedInternalField(vf, pf);
602  }
603 
604  return pf;
605  }
606  }
607 }
608 
609 
610 template<class Type>
613 (
615 ) const
616 {
617  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
618 }
619 
620 
621 template<class Type>
624 (
626 ) const
627 {
628  // Construct tmp<pointField>
630  tvf.clear();
631  return tpf;
632 }
633 
634 
635 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
636 
637 } // End namespace Foam
638 
639 // ************************************************************************* //
volFields.H
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::globalMeshData::globalCoPointSlavesMap
const mapDistribute & globalCoPointSlavesMap() const
Definition: globalMeshData.C:2363
Foam::volPointInterpolation::flatBoundaryField
tmp< Field< Type > > flatBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Get boundary field in same order as boundary faces. Field is.
Definition: volPointInterpolate.C:224
coupledPointPatchField.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
boundary
faceListList boundary(nPatches)
Foam::UPstream::waitRequests
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:106
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::pointData
Variant of pointEdgePoint with some transported additional data. WIP - should be templated on data li...
Definition: pointData.H:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::solution::cachePrintMessage
static void cachePrintMessage(const char *message, const word &name, const FieldType &vf)
Helper for printing cache message.
Definition: solutionTemplates.C:33
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::MeshObject< pointMesh, UpdateableMeshObject, pointConstraints >::New
static const pointConstraints & New(const pointMesh &mesh)
Definition: MeshObject.C:44
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("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]=cellShape(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
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:52
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::UPstream::nonBlocking
@ nonBlocking
Definition: UPstream.H:68
Foam::interpolate
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::DimensionedField::mesh
const Mesh & mesh() const
Return mesh.
Definition: DimensionedFieldI.H:38
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
pointConstraints.H
Foam::volPointInterpolation::addSeparated
void addSeparated(GeometricField< Type, pointPatchField, pointMesh > &) const
Add separated contributions.
Definition: volPointInterpolate.C:85
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:106
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::plusEqOp
Definition: ops.H:71
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::pointMesh
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
Foam::tmp::ptr
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:146
emptyFvPatch.H
s
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))
Foam::DimensionedField::dimensions
const dimensionSet & dimensions() const
Return dimensions.
Definition: DimensionedFieldI.H:46
Foam::volPointInterpolation::interpolate
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
Definition: volPointInterpolate.C:502
Foam::pointConstraints
Application of (multi-)patch point contraints.
Definition: pointConstraints.H:62
volPointInterpolation.H
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::mapDistributeBase::constructSize
label constructSize() const
Constructed data size.
Definition: mapDistributeBase.H:244
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
Foam::regIOobject::store
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:34
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::pointConstraints::syncUntransformedData
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
Definition: pointConstraintsTemplates.C:39
Foam::List< Type >
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
Foam::globalMeshData::coupledPatch
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Definition: globalMeshData.C:2046
Foam::pointMesh::thisDb
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:118
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::volPointInterpolation::interpolateInternalField
void interpolateInternalField(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate internal field from volField to pointField.
Definition: volPointInterpolate.C:127
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:108
Foam::mapDistribute::reverseDistribute
void reverseDistribute(const label constructSize, List< T > &, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
Definition: mapDistributeTemplates.C:187
Foam::pointConstraints::constrain
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
Definition: pointConstraintsTemplates.C:131
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatchTemplate.C:392
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::volPointInterpolation::interpolateDimensionedInternalField
void interpolateDimensionedInternalField(const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const
Interpolate dimensioned internal field from cells to points.
Definition: volPointInterpolate.C:164
Foam::volPointInterpolation::interpolateBoundaryField
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
Definition: volPointInterpolate.C:271
Foam::volPointInterpolation::pushUntransformedData
void pushUntransformedData(List< Type > &) const
Helper: push master point data to collocated points.
Definition: volPointInterpolate.C:42
pointFields.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::globalMeshData::globalCoPointSlaves
const labelListList & globalCoPointSlaves() const
Definition: globalMeshData.C:2353
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88