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-------------------------------------------------------------------------------
11License
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
30#include "volFields.H"
31#include "pointFields.H"
32#include "emptyFvPatch.H"
34#include "pointConstraints.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38template<class Type>
39void 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
81template<class Type>
82void 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.ref();
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 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
128template<class Type>
130(
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
165template<class Type>
167(
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
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
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
225template<class Type>
226Foam::tmp<Foam::Field<Type>> Foam::volPointInterpolation::flatBoundaryField
227(
229) const
230{
231 const fvMesh& mesh = vf.mesh();
232 const fvBoundaryMesh& bm = mesh.boundary();
233
234 tmp<Field<Type>> tboundaryVals
235 (
237 );
238 Field<Type>& boundaryVals = tboundaryVals.ref();
239
240 forAll(vf.boundaryField(), patchi)
241 {
242 label bFacei = bm[patchi].patch().start() - mesh.nInternalFaces();
243
244 if
245 (
246 !isA<emptyFvPatch>(bm[patchi])
247 && !vf.boundaryField()[patchi].coupled()
248 )
249 {
251 (
252 boundaryVals,
253 vf.boundaryField()[patchi].size(),
254 bFacei
255 ) = vf.boundaryField()[patchi];
256 }
257 else
258 {
259 const polyPatch& pp = bm[patchi].patch();
260
261 forAll(pp, i)
262 {
263 boundaryVals[bFacei++] = Zero;
264 }
265 }
266 }
267
268 return tboundaryVals;
269}
270
271
272template<class Type>
274(
277) const
278{
279 const primitivePatch& boundary = boundaryPtr_();
280
281 Field<Type>& pfi = pf.primitiveFieldRef();
282
283 // Get face data in flat list
284 tmp<Field<Type>> tboundaryVals(flatBoundaryField(vf));
285 const Field<Type>& boundaryVals = tboundaryVals();
286
287
288 // Do points on 'normal' patches from the surrounding patch faces
289 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
290
291 const labelList& mp = boundary.meshPoints();
292
293 forAll(mp, i)
294 {
295 label pointi = mp[i];
296
297 if (isPatchPoint_[pointi])
298 {
299 const labelList& pFaces = boundary.pointFaces()[i];
300 const scalarList& pWeights = boundaryPointWeights_[i];
301
302 Type& val = pfi[pointi];
303
304 val = Zero;
305 forAll(pFaces, j)
306 {
307 if (boundaryIsPatchFace_[pFaces[j]])
308 {
309 val += pWeights[j]*boundaryVals[pFaces[j]];
310 }
311 }
312 }
313 }
314
315 // Sum collocated contributions
317
318 // And add separated contributions
319 addSeparated(pf);
320
321 // Optionally normalise
322 if (normalisationPtr_)
323 {
324 const scalarField& normalisation = normalisationPtr_();
325 forAll(mp, i)
326 {
327 pfi[mp[i]] *= normalisation[i];
328 }
329 }
330
331
332 // Push master data to slaves. It is possible (not sure how often) for
333 // a coupled point to have its master on a different patch so
334 // to make sure just push master data to slaves.
335 pushUntransformedData(pfi);
336}
337
338
339template<class Type>
341(
344 const bool overrideFixedValue
345) const
346{
347 interpolateBoundaryField(vf, pf);
348
349 // Apply constraints
351
352 pcs.constrain(pf, overrideFixedValue);
353}
354
355
356template<class Type>
358(
361) const
362{
363 if (debug)
364 {
365 Pout<< "volPointInterpolation::interpolate("
366 << "const GeometricField<Type, fvPatchField, volMesh>&, "
367 << "GeometricField<Type, pointPatchField, pointMesh>&) : "
368 << "interpolating field " << vf.name() << " from cells to points "
369 << pf.name() << endl;
370 }
371
372 interpolateInternalField(vf, pf);
373
374 // Interpolate to the patches preserving fixed value BCs
375 interpolateBoundaryField(vf, pf, false);
376}
377
378
379template<class Type>
382(
384 const wordList& patchFieldTypes
385) const
386{
387 const pointMesh& pm = pointMesh::New(vf.mesh());
388
389 // Construct tmp<pointField>
391 (
393 (
394 "volPointInterpolate(" + vf.name() + ')',
395 vf.instance(),
396 pm.thisDb()
397 ),
398 pm,
399 vf.dimensions(),
400 patchFieldTypes
401 );
402
403 interpolateInternalField(vf, tpf.ref());
404
405 // Interpolate to the patches overriding fixed value BCs
406 interpolateBoundaryField(vf, tpf.ref(), true);
407
408 return tpf;
409}
410
411
412template<class Type>
415(
417 const wordList& patchFieldTypes
418) const
419{
420 // Construct tmp<pointField>
422 interpolate(tvf(), patchFieldTypes);
423 tvf.clear();
424 return tpf;
425}
426
427
428template<class Type>
431(
433 const word& name,
434 const bool cache
435) const
436{
438
439 const pointMesh& pm = pointMesh::New(vf.mesh());
440 const objectRegistry& db = pm.thisDb();
441
442 PointFieldType* pfPtr =
443 db.objectRegistry::template getObjectPtr<PointFieldType>(name);
444
445 if (!cache || vf.mesh().changing())
446 {
447 // Delete any old occurrences to avoid double registration
448 if (pfPtr && pfPtr->ownedByRegistry())
449 {
450 solution::cachePrintMessage("Deleting", name, vf);
451 delete pfPtr;
452 }
453
455 (
457 (
458 name,
459 vf.instance(),
460 pm.thisDb()
461 ),
462 pm,
463 vf.dimensions()
464 );
465
466 interpolate(vf, tpf.ref());
467
468 return tpf;
469 }
470
471
472 if (!pfPtr)
473 {
474 solution::cachePrintMessage("Calculating and caching", name, vf);
475
476 pfPtr = interpolate(vf, name, false).ptr();
477 regIOobject::store(pfPtr);
478 }
479 else
480 {
481 PointFieldType& pf = *pfPtr;
482
483 if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
484 {
485 solution::cachePrintMessage("Reusing", name, vf);
486 }
487 else
488 {
489 solution::cachePrintMessage("Updating", name, vf);
490 interpolate(vf, pf);
491 }
492 }
493
494 return *pfPtr;
495}
496
497
498template<class Type>
501(
503) const
504{
505 return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
506}
507
508
509template<class Type>
512(
514) const
515{
516 // Construct tmp<pointField>
518 interpolate(tvf());
519 tvf.clear();
520 return tpf;
521}
522
523
524template<class Type>
527(
529 const word& name,
530 const bool cache
531) const
532{
533 typedef DimensionedField<Type, pointMesh> PointFieldType;
534
535 const pointMesh& pm = pointMesh::New(vf.mesh());
536 const objectRegistry& db = pm.thisDb();
537
538
539 PointFieldType* pfPtr =
540 db.objectRegistry::template getObjectPtr<PointFieldType>(name);
541
542 if (!cache || vf.mesh().changing())
543 {
544 // Delete any old occurrences to avoid double registration
545 if (pfPtr && pfPtr->ownedByRegistry())
546 {
547 solution::cachePrintMessage("Deleting", name, vf);
548 delete pfPtr;
549 }
550
552 (
554 (
555 name,
556 vf.instance(),
557 pm.thisDb()
558 ),
559 pm,
560 vf.dimensions()
561 );
562
563 interpolateDimensionedInternalField(vf, tpf.ref());
564
565 return tpf;
566 }
567
568
569 if (!pfPtr)
570 {
571 solution::cachePrintMessage("Calculating and caching", name, vf);
572 pfPtr = interpolate(vf, name, false).ptr();
573
574 regIOobject::store(pfPtr);
575 }
576 else
577 {
578 PointFieldType& pf = *pfPtr;
579
580 if (pf.upToDate(vf)) //TBD: , vf.mesh().points()))
581 {
582 solution::cachePrintMessage("Reusing", name, vf);
583 }
584 else
585 {
586 solution::cachePrintMessage("Updating", name, vf);
587 interpolateDimensionedInternalField(vf, pf);
588 }
589 }
590
591 return *pfPtr;
592}
593
594
595template<class Type>
598(
600) const
601{
602 return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
603}
604
605
606template<class Type>
609(
611) const
612{
613 // Construct tmp<pointField>
615 tvf.clear();
616 return tpf;
617}
618
619
620// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:70
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:90
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:100
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
Registry of regIOobjects.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:59
Application of (multi-)patch point constraints.
void constrain(GeometricField< Type, pointPatchField, pointMesh > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
static void syncUntransformedData(const polyMesh &mesh, List< Type > &pointData, const CombineOp &cop)
Helper: sync data on collocated points only.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:126
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
const labelListList & pointCells() const
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
bool interpolate() const noexcept
Same as isPointData()
static void cachePrintMessage(const char *message, const word &name, const FieldType &vf)
Helper for printing cache message.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
void interpolateBoundaryField(const GeometricField< Type, fvPatchField, volMesh > &vf, GeometricField< Type, pointPatchField, pointMesh > &pf) const
Interpolate boundary field without applying constraints/boundary.
void interpolateInternalField(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate internal field from volField to pointField.
void interpolateDimensionedInternalField(const DimensionedField< Type, volMesh > &vf, DimensionedField< Type, pointMesh > &pf) const
Interpolate dimensioned internal field from cells to points.
A class for handling words, derived from Foam::string.
Definition: word.H:68
faceListList boundary
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
dynamicFvMesh & mesh
const pointField & points
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))
List< label > labelList
A List of labels.
Definition: List.H:66
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
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
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333