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