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