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 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"
34#include "surfaceFields.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 {
54 if (isA<coupledFacePointPatch>(pp) && !isType<processorPointPatch>(pp))
55 {
56 hasSpecial = true;
57 break;
58 }
59 }
60
61 reduce(hasSpecial, orOp<bool>());
62 return hasSpecial;
63}
64
65
66void Foam::volPointInterpolation::calcBoundaryAddressing()
67{
68 if (debug)
69 {
70 Pout<< "volPointInterpolation::calcBoundaryAddressing() : "
71 << "constructing boundary addressing"
72 << endl;
73 }
74
75 boundaryPtr_.reset
76 (
78 (
79 SubList<face>
80 (
81 mesh().faces(),
82 mesh().nBoundaryFaces(),
83 mesh().nInternalFaces()
84 ),
85 mesh().points()
86 )
87 );
88 const primitivePatch& boundary = boundaryPtr_();
89
90 boundaryIsPatchFace_.setSize(boundary.size());
91 boundaryIsPatchFace_ = false;
92
93 // Store per mesh point whether it is on any 'real' patch. Currently
94 // boolList just so we can use syncUntransformedData (does not take
95 // bitSet. Tbd)
96 boolList isPatchPoint(mesh().nPoints(), false);
97
98 const polyBoundaryMesh& pbm = mesh().boundaryMesh();
99
100 // Get precalculated volField only so we can use coupled() tests for
101 // cyclicAMI
102 const surfaceScalarField& magSf = mesh().magSf();
103
104 forAll(pbm, patchi)
105 {
106 const polyPatch& pp = pbm[patchi];
107
108 if
109 (
110 !isA<emptyPolyPatch>(pp)
111 && !magSf.boundaryField()[patchi].coupled()
112 )
113 {
114 label bFacei = pp.start()-mesh().nInternalFaces();
115
116 forAll(pp, i)
117 {
118 boundaryIsPatchFace_[bFacei] = true;
119
120 const face& f = boundary[bFacei++];
121
122 forAll(f, fp)
123 {
124 isPatchPoint[f[fp]] = true;
125 }
126 }
127 }
128 }
129
130 // Make sure point status is synchronised so even processor that holds
131 // no face of a certain patch still can have boundary points marked.
133 (
134 mesh(),
135 isPatchPoint,
136 orEqOp<bool>()
137 );
138
139 // Convert to bitSet
140 isPatchPoint_.setSize(mesh().nPoints());
141 isPatchPoint_.assign(isPatchPoint);
142
143 if (debug)
144 {
145 label nPatchFace = 0;
146 forAll(boundaryIsPatchFace_, i)
147 {
148 if (boundaryIsPatchFace_[i])
149 {
150 nPatchFace++;
151 }
152 }
153 label nPatchPoint = 0;
154 forAll(isPatchPoint_, i)
155 {
156 if (isPatchPoint_[i])
157 {
158 nPatchPoint++;
159 }
160 }
161 Pout<< "boundary:" << nl
162 << " faces :" << boundary.size() << nl
163 << " of which on proper patch:" << nPatchFace << nl
164 << " points:" << boundary.nPoints() << nl
165 << " of which on proper patch:" << nPatchPoint << endl;
166 }
167}
168
169
170void Foam::volPointInterpolation::makeInternalWeights(scalarField& sumWeights)
171{
172 if (debug)
173 {
174 Pout<< "volPointInterpolation::makeInternalWeights() : "
175 << "constructing weighting factors for internal and non-coupled"
176 << " points." << endl;
177 }
178
179 const pointField& points = mesh().points();
180 const labelListList& pointCells = mesh().pointCells();
181 const vectorField& cellCentres = mesh().cellCentres();
182
183 // Allocate storage for weighting factors
184 pointWeights_.clear();
185 pointWeights_.setSize(points.size());
186
187 // Calculate inverse distances between cell centres and points
188 // and store in weighting factor array
189 forAll(points, pointi)
190 {
191 if (!isPatchPoint_[pointi])
192 {
193 const labelList& pcp = pointCells[pointi];
194
195 scalarList& pw = pointWeights_[pointi];
196 pw.setSize(pcp.size());
197
198 forAll(pcp, pointCelli)
199 {
200 pw[pointCelli] =
201 1.0/mag(points[pointi] - cellCentres[pcp[pointCelli]]);
202
203 sumWeights[pointi] += pw[pointCelli];
204 }
205 }
206 }
207}
208
209
210void Foam::volPointInterpolation::makeBoundaryWeights(scalarField& sumWeights)
211{
212 if (debug)
213 {
214 Pout<< "volPointInterpolation::makeBoundaryWeights() : "
215 << "constructing weighting factors for boundary points." << endl;
216 }
217
218 const pointField& points = mesh().points();
219 const pointField& faceCentres = mesh().faceCentres();
220
221 const primitivePatch& boundary = boundaryPtr_();
222
223 boundaryPointWeights_.clear();
224 boundaryPointWeights_.setSize(boundary.meshPoints().size());
225
226 forAll(boundary.meshPoints(), i)
227 {
228 label pointi = boundary.meshPoints()[i];
229
230 if (isPatchPoint_[pointi])
231 {
232 const labelList& pFaces = boundary.pointFaces()[i];
233
234 scalarList& pw = boundaryPointWeights_[i];
235 pw.setSize(pFaces.size());
236
237 sumWeights[pointi] = 0.0;
238
239 forAll(pFaces, i)
240 {
241 if (boundaryIsPatchFace_[pFaces[i]])
242 {
243 label facei = mesh().nInternalFaces() + pFaces[i];
244
245 pw[i] = 1.0/mag(points[pointi] - faceCentres[facei]);
246 sumWeights[pointi] += pw[i];
247 }
248 else
249 {
250 pw[i] = 0.0;
251 }
252 }
253 }
254 }
255}
256
257
258void Foam::volPointInterpolation::interpolateOne
259(
260 const tmp<scalarField>& tnormalisation,
262) const
263{
264 if (debug)
265 {
266 Pout<< "volPointInterpolation::interpolateOne("
267 << "pointScalarField&) : "
268 << "interpolating oneField"
269 << " from cells to BOUNDARY points "
270 << pf.name() << endl;
271 }
272
273 const primitivePatch& boundary = boundaryPtr_();
274 const labelList& mp = boundary.meshPoints();
275 Field<scalar>& pfi = pf.primitiveFieldRef();
276
277
278 // 1. Interpolate coupled boundary points from cells
279 {
280 forAll(mp, i)
281 {
282 const label pointi = mp[i];
283 if (!isPatchPoint_[pointi])
284 {
285 const scalarList& pw = pointWeights_[pointi];
286
287 scalar& val = pfi[pointi];
288
289 val = Zero;
290 forAll(pw, pointCelli)
291 {
292 val += pw[pointCelli];
293 }
294 }
295 }
296 }
297
298
299 // 2. Interpolate to the patches preserving fixed value BCs
300 {
301 forAll(mp, i)
302 {
303 const label pointi = mp[i];
304
305 if (isPatchPoint_[pointi])
306 {
307 const labelList& pFaces = boundary.pointFaces()[i];
308 const scalarList& pWeights = boundaryPointWeights_[i];
309
310 scalar& val = pfi[pointi];
311
312 val = Zero;
313 forAll(pFaces, j)
314 {
315 if (boundaryIsPatchFace_[pFaces[j]])
316 {
317 val += pWeights[j];
318 }
319 }
320 }
321 }
322
323 // Sum collocated contributions
325 (
326 mesh(),
327 pfi,
328 plusEqOp<scalar>()
329 );
330
331 // And add separated contributions
332 addSeparated(pf);
333
334 // Optionally normalise
335 if (tnormalisation)
336 {
337 const scalarField& normalisation = tnormalisation();
338 forAll(mp, i)
339 {
340 pfi[mp[i]] *= normalisation[i];
341 }
342 }
343 }
344
345 // Apply constraints
346 pointConstraints::New(pf.mesh()).constrain(pf, false);
347}
348
349
350void Foam::volPointInterpolation::makeWeights()
351{
352 if (debug)
353 {
354 Pout<< "volPointInterpolation::makeWeights() : "
355 << "constructing weighting factors"
356 << endl;
357 }
358
359 const pointMesh& pMesh = pointMesh::New(mesh());
360
361 // Update addressing over all boundary faces
362 calcBoundaryAddressing();
363
364
365 // Running sum of weights
366 tmp<pointScalarField> tsumWeights
367 (
369 (
370 IOobject
371 (
372 "volPointSumWeights",
374 mesh()
375 ),
376 pMesh,
378 )
379 );
380 pointScalarField& sumWeights = tsumWeights.ref();
381
382
383 // Create internal weights; add to sumWeights
384 makeInternalWeights(sumWeights);
385
386
387 // Create boundary weights; override sumWeights
388 makeBoundaryWeights(sumWeights);
389
390
391 const primitivePatch& boundary = boundaryPtr_();
392 const labelList& mp = boundary.meshPoints();
393
394
395 // Sum collocated contributions
397 (
398 mesh(),
399 sumWeights,
400 plusEqOp<scalar>()
401 );
402
403
404 // And add separated contributions
405 addSeparated(sumWeights);
406
407
408 // Push master data to slaves. It is possible (not sure how often) for
409 // a coupled point to have its master on a different patch so
410 // to make sure just push master data to slaves. Reuse the syncPointData
411 // structure.
412 pushUntransformedData(sumWeights);
413
414
415 // Normalise internal weights
416 forAll(pointWeights_, pointi)
417 {
418 scalarList& pw = pointWeights_[pointi];
419 // Note:pw only sized for !isPatchPoint
420 forAll(pw, i)
421 {
422 pw[i] /= sumWeights[pointi];
423 }
424 }
425
426 // Normalise boundary weights
427 forAll(mp, i)
428 {
429 const label pointi = mp[i];
430
431 scalarList& pw = boundaryPointWeights_[i];
432 // Note:pw only sized for isPatchPoint
433 forAll(pw, i)
434 {
435 pw[i] /= sumWeights[pointi];
436 }
437 }
438
439
440 // Normalise separated contributions
441 if (hasSeparated_)
442 {
443 if (debug)
444 {
445 Pout<< "volPointInterpolation::makeWeights() : "
446 << "detected separated coupled patches"
447 << " - allocating normalisation" << endl;
448 }
449
450 // Sum up effect of interpolating one (on boundary points only)
451 interpolateOne(tmp<scalarField>(), sumWeights);
452
453 // Store as normalisation factor (on boundary points only)
454 normalisationPtr_ = new scalarField(mp.size());
455 normalisationPtr_.ref() = scalar(1.0);
456 normalisationPtr_.ref() /= scalarField(sumWeights, mp);
457 }
458 else
459 {
460 normalisationPtr_.clear();
461 }
462
463
464 if (debug)
465 {
466 Pout<< "volPointInterpolation::makeWeights() : "
467 << "finished constructing weighting factors"
468 << endl;
469 }
470}
471
472
473// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
474
476:
478 hasSeparated_(hasSeparated(pointMesh::New(vm)))
479{
480 makeWeights();
481}
482
483
484// * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
485
487{}
488
489
490// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
491
493{
494 // Recheck whether has coupled patches
495 hasSeparated_ = hasSeparated(pointMesh::New(mesh()));
496
497 makeWeights();
498}
499
500
502{
503 makeWeights();
504
505 return true;
506}
507
508
510(
511 const volVectorField& vf,
513) const
514{
515 interpolateInternalField(vf, pf);
516
517 // Interpolate to the patches but no constraints
518 interpolateBoundaryField(vf, pf);
519
520 // Apply displacement constraints
522
523 pcs.constrainDisplacement(pf, false);
524}
525
526
527// ************************************************************************* //
const Mesh & mesh() const
Return mesh.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:91
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
void updateMesh()
Update for new mesh topology.
Application of (multi-)patch point constraints.
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:55
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1083
const vectorField & faceCentres() const
const labelListList & pointCells() const
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
Interpolate from cell centres to points (vertices) using inverse distance weighting.
bool movePoints()
Correct weighting factors for moving mesh.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
faceListList boundary
dynamicFvMesh & mesh
const pointField & points
label nPoints
const dimensionedScalar mp
Proton mass.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition: List.H:66
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
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
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< bool > boolList
A List of bools.
Definition: List.H:64
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
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
labelList f(nPoints)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Foam::surfaceFields.