Loading...
Searching...
No Matches
volPointInterpolation.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-2016 OpenFOAM Foundation
9 Copyright (C) 2020-2022 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 "fvMesh.H"
31#include "volFields.H"
32#include "pointFields.H"
33#include "pointConstraints.H"
35#include "processorPointPatch.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
42}
43
44
45// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46
47bool Foam::volPointInterpolation::hasSeparated(const pointMesh& pMesh)
48{
49 const pointBoundaryMesh& pbm = pMesh.boundary();
50
51 bool hasSpecial = false;
52 for (const auto& pp : pbm)
53 {
55 {
56 hasSpecial = true;
57 break;
58 }
59 }
60
61 return returnReduceOr(hasSpecial);
62}
63
64
65void Foam::volPointInterpolation::calcBoundaryAddressing()
66{
67 if (debug)
68 {
69 Pout<< "volPointInterpolation::calcBoundaryAddressing() : "
70 << "constructing boundary addressing"
71 << endl;
72 }
73
75
76 boundaryPtr_.reset(new primitivePatch(pbm.faces(), mesh().points()));
77 const auto& boundary = *boundaryPtr_;
78
79 boundaryIsPatchFace_.resize_nocopy(boundary.size());
80 boundaryIsPatchFace_ = false;
81
82 // Store per mesh point whether it is on any 'real' patch. Currently
83 // boolList just so we can use syncUntransformedData (does not take
84 // bitSet. Tbd)
85 boolList isPatchPoint(mesh().nPoints(), false);
86
87 // Get precalculated volField only so we can use coupled() tests for
88 // cyclicAMI
89 const surfaceScalarField& magSf = mesh().magSf();
90
91 forAll(pbm, patchi)
92 {
93 const polyPatch& pp = pbm[patchi];
94
95 if
96 (
98 && !magSf.boundaryField()[patchi].coupled()
99 )
100 {
101 boundaryIsPatchFace_.set(labelRange(pp.offset(), pp.size()));
102
103 for (const auto& f : pp.faces())
104 {
105 UIndirectList<bool>(isPatchPoint, f) = true;
106 }
107 }
108 }
109
110 // Make sure point status is synchronised so even processor that holds
111 // no face of a certain patch still can have boundary points marked.
113 (
114 mesh(),
115 isPatchPoint,
117 );
118
119 // Convert to bitSet
120 isPatchPoint_.reset();
121 isPatchPoint_.resize(mesh().nPoints());
122 isPatchPoint_.assign(isPatchPoint);
123
124 if (debug)
125 {
126 label nPatchFace = boundaryIsPatchFace_.count();
127 label nPatchPoint = isPatchPoint_.count();
128
129 Pout<< "boundary:" << nl
130 << " faces :" << boundary.size() << nl
131 << " of which on proper patch:" << nPatchFace << nl
132 << " points:" << boundary.nPoints() << nl
133 << " of which on proper patch:" << nPatchPoint << endl;
134 }
135}
136
137
138void Foam::volPointInterpolation::makeInternalWeights(scalarField& sumWeights)
139{
140 if (debug)
141 {
142 Pout<< "volPointInterpolation::makeInternalWeights() : "
143 << "constructing weighting factors for internal and non-coupled"
144 << " points." << endl;
145 }
146
147 const pointField& points = mesh().points();
149 const vectorField& cellCentres = mesh().cellCentres();
150
151 // Allocate storage for weighting factors
152 pointWeights_.clear();
153 pointWeights_.setSize(points.size());
154
155 // Calculate inverse distances between cell centres and points
156 // and store in weighting factor array
157 forAll(points, pointi)
158 {
159 if (!isPatchPoint_[pointi])
160 {
161 const labelList& pcp = pointCells[pointi];
162
163 scalarList& pw = pointWeights_[pointi];
164 pw.setSize(pcp.size());
165
166 forAll(pcp, pointCelli)
167 {
168 pw[pointCelli] =
169 1.0/mag(points[pointi] - cellCentres[pcp[pointCelli]]);
170
171 sumWeights[pointi] += pw[pointCelli];
172 }
173 }
174 }
175}
176
177
178void Foam::volPointInterpolation::makeBoundaryWeights(scalarField& sumWeights)
179{
180 if (debug)
181 {
182 Pout<< "volPointInterpolation::makeBoundaryWeights() : "
183 << "constructing weighting factors for boundary points." << endl;
184 }
185
186 const pointField& points = mesh().points();
187 const pointField& faceCentres = mesh().faceCentres();
188
189 const primitivePatch& boundary = boundaryPtr_();
190
191 boundaryPointWeights_.clear();
192 boundaryPointWeights_.setSize(boundary.meshPoints().size());
193
194 forAll(boundary.meshPoints(), i)
195 {
196 label pointi = boundary.meshPoints()[i];
197
198 if (isPatchPoint_[pointi])
199 {
200 const labelList& pFaces = boundary.pointFaces()[i];
201
202 scalarList& pw = boundaryPointWeights_[i];
203 pw.setSize(pFaces.size());
204
205 sumWeights[pointi] = 0.0;
206
207 forAll(pFaces, i)
208 {
209 if (boundaryIsPatchFace_[pFaces[i]])
210 {
211 label facei = mesh().nInternalFaces() + pFaces[i];
212
213 pw[i] = 1.0/mag(points[pointi] - faceCentres[facei]);
214 sumWeights[pointi] += pw[i];
215 }
216 else
217 {
218 pw[i] = 0.0;
219 }
220 }
221 }
222 }
223}
224
225
226void Foam::volPointInterpolation::interpolateOne
227(
228 const tmp<scalarField>& tnormalisation,
230) const
231{
232 if (debug)
233 {
234 Pout<< "volPointInterpolation::interpolateOne("
235 << "pointScalarField&) : "
236 << "interpolating oneField"
237 << " from cells to BOUNDARY points "
238 << pf.name() << endl;
239 }
240
241 const primitivePatch& boundary = boundaryPtr_();
242 const labelList& mp = boundary.meshPoints();
243 Field<scalar>& pfi = pf.primitiveFieldRef();
244
245
246 // 1. Interpolate coupled boundary points from cells
247 {
248 forAll(mp, i)
249 {
250 const label pointi = mp[i];
251 if (!isPatchPoint_[pointi])
252 {
253 const scalarList& pw = pointWeights_[pointi];
254
255 scalar& val = pfi[pointi];
256
257 val = Zero;
258 forAll(pw, pointCelli)
259 {
260 val += pw[pointCelli];
261 }
262 }
263 }
264 }
265
266
267 // 2. Interpolate to the patches preserving fixed value BCs
268 {
269 forAll(mp, i)
270 {
271 const label pointi = mp[i];
272
273 if (isPatchPoint_[pointi])
274 {
275 const labelList& pFaces = boundary.pointFaces()[i];
276 const scalarList& pWeights = boundaryPointWeights_[i];
277
278 scalar& val = pfi[pointi];
279
280 val = Zero;
281 forAll(pFaces, j)
282 {
283 if (boundaryIsPatchFace_[pFaces[j]])
284 {
285 val += pWeights[j];
286 }
287 }
288 }
289 }
290
291 // Sum collocated contributions
293 (
294 mesh(),
295 pfi,
297 );
298
299 // And add separated contributions
300 addSeparated(pf);
301
302 // Optionally normalise
303 if (tnormalisation)
304 {
305 const scalarField& normalisation = tnormalisation();
306 forAll(mp, i)
307 {
308 pfi[mp[i]] *= normalisation[i];
309 }
310 }
311 }
312
313 // Apply constraints
314 pointConstraints::New(pf.mesh()).constrain(pf, false);
315}
316
317
318void Foam::volPointInterpolation::makeWeights()
319{
320 if (debug)
321 {
322 Pout<< "volPointInterpolation::makeWeights() : "
323 << "constructing weighting factors"
324 << endl;
325 }
326
327 const pointMesh& pMesh = pointMesh::New(mesh());
328
329 // Update addressing over all boundary faces
330 calcBoundaryAddressing();
331
332
333 // Running sum of weights
334 tmp<pointScalarField> tsumWeights
335 (
337 (
339 (
340 "volPointSumWeights",
342 mesh()
343 ),
344 pMesh,
346 )
347 );
348 pointScalarField& sumWeights = tsumWeights.ref();
349
350
351 // Create internal weights; add to sumWeights
352 makeInternalWeights(sumWeights);
353
354
355 // Create boundary weights; override sumWeights
356 makeBoundaryWeights(sumWeights);
357
358
359 const primitivePatch& boundary = boundaryPtr_();
360 const labelList& mp = boundary.meshPoints();
361
362
363 // Sum collocated contributions
365 (
366 mesh(),
367 sumWeights,
369 );
370
371
372 // And add separated contributions
373 addSeparated(sumWeights);
374
375
376 // Push master data to slaves. It is possible (not sure how often) for
377 // a coupled point to have its master on a different patch so
378 // to make sure just push master data to slaves. Reuse the syncPointData
379 // structure.
380 pushUntransformedData(sumWeights);
381
382
383 // Normalise internal weights
384 forAll(pointWeights_, pointi)
385 {
386 scalarList& pw = pointWeights_[pointi];
387 // Note:pw only sized for !isPatchPoint
388 forAll(pw, i)
389 {
390 pw[i] /= sumWeights[pointi];
391 }
392 }
393
394 // Normalise boundary weights
395 forAll(mp, i)
396 {
397 const label pointi = mp[i];
398
399 scalarList& pw = boundaryPointWeights_[i];
400 // Note:pw only sized for isPatchPoint
401 forAll(pw, i)
402 {
403 pw[i] /= sumWeights[pointi];
404 }
405 }
406
407
408 // Normalise separated contributions
409 if (hasSeparated_)
410 {
411 if (debug)
412 {
413 Pout<< "volPointInterpolation::makeWeights() : "
414 << "detected separated coupled patches"
415 << " - allocating normalisation" << endl;
416 }
417
418 // Sum up effect of interpolating one (on boundary points only)
419 interpolateOne(tmp<scalarField>(), sumWeights);
420
421 // Store as normalisation factor (on boundary points only)
422 normalisationPtr_ = new scalarField(mp.size());
423 normalisationPtr_.ref() = scalar(1.0);
424 normalisationPtr_.ref() /= scalarField(sumWeights, mp);
425 }
426 else
427 {
428 normalisationPtr_.clear();
429 }
430
431
432 if (debug)
433 {
434 Pout<< "volPointInterpolation::makeWeights() : "
435 << "finished constructing weighting factors"
437 }
438}
439
440
441// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
442
443Foam::volPointInterpolation::volPointInterpolation(const fvMesh& vm)
444:
445 MeshObject_type(vm),
446 hasSeparated_(hasSeparated(pointMesh::New(vm)))
448 makeWeights();
449}
450
451
452// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
455{}
456
457
458// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
459
461{
462 // Recheck whether has coupled patches
463 hasSeparated_ = hasSeparated(pointMesh::New(mesh()));
464
465 makeWeights();
466}
467
468
471 makeWeights();
472
473 return true;
474}
475
476
478(
479 const volVectorField& vf,
481) const
482{
483 interpolateInternalField(vf, pf);
484
485 // Interpolate to the patches but no constraints
486 interpolateBoundaryField(vf, pf);
487
488 // Apply displacement constraints
490
491 pcs.constrainDisplacement(pf, false);
492}
493
494
495// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
const Mesh & mesh() const noexcept
Return const reference to mesh.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
label size() const noexcept
The number of elements in the list.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
static FOAM_NO_DANGLING_REFERENCE const pointConstraints & New(const pointMesh &mesh, Args &&... args)
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition pointCells.H:55
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.
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const faceList::subList faces() const
Return mesh faces for the entire boundary.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const vectorField & faceCentres() const
const labelListList & pointCells() const
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
A class for managing temporary objects.
Definition tmp.H:75
Interpolate from cell centres to points (vertices) using inverse distance weighting.
bool movePoints()
Correct weighting factors for moving mesh.
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 interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
faceListList boundary
dynamicFvMesh & mesh
const pointField & points
label nPoints
const dimensionedScalar mp
Proton mass.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
GeometricField< vector, pointPatchField, pointMesh > pointVectorField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition List.H:60
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
bool isType(const U &obj)
Check if typeid of the object and Type are identical.
Definition typeInfo.H:112
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
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];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } 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};constexpr 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< DynamicList< face > > pFaces[nBCs]
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.