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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "volPointInterpolation.H"
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "pointFields.H"
32 #include "demandDrivenData.H"
33 #include "pointConstraints.H"
34 #include "surfaceFields.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(volPointInterpolation, 0);
41 }
42 
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::volPointInterpolation::calcBoundaryAddressing()
47 {
48  if (debug)
49  {
50  Pout<< "volPointInterpolation::calcBoundaryAddressing() : "
51  << "constructing boundary addressing"
52  << endl;
53  }
54 
55  boundaryPtr_.reset
56  (
57  new primitivePatch
58  (
59  SubList<face>
60  (
61  mesh().faces(),
62  mesh().nBoundaryFaces(),
63  mesh().nInternalFaces()
64  ),
65  mesh().points()
66  )
67  );
68  const primitivePatch& boundary = boundaryPtr_();
69 
70  boundaryIsPatchFace_.setSize(boundary.size());
71  boundaryIsPatchFace_ = false;
72 
73  isPatchPoint_.setSize(mesh().nPoints());
74  isPatchPoint_ = false;
75 
76  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
77 
78  // Get precalculated volField only so we can use coupled() tests for
79  // cyclicAMI
80  const surfaceScalarField& magSf = mesh().magSf();
81 
82  forAll(pbm, patchi)
83  {
84  const polyPatch& pp = pbm[patchi];
85 
86  if
87  (
88  !isA<emptyPolyPatch>(pp)
89  && !magSf.boundaryField()[patchi].coupled()
90  )
91  {
92  label bFacei = pp.start()-mesh().nInternalFaces();
93 
94  forAll(pp, i)
95  {
96  boundaryIsPatchFace_[bFacei] = true;
97 
98  const face& f = boundary[bFacei++];
99 
100  forAll(f, fp)
101  {
102  isPatchPoint_[f[fp]] = true;
103  }
104  }
105  }
106  }
107 
108  // Make sure point status is synchronised so even processor that holds
109  // no face of a certain patch still can have boundary points marked.
110  if (debug)
111  {
112  boolList oldData(isPatchPoint_);
113 
115  (
116  mesh(),
117  isPatchPoint_,
118  orEqOp<bool>()
119  );
120 
121  forAll(isPatchPoint_, pointi)
122  {
123  if (isPatchPoint_[pointi] != oldData[pointi])
124  {
125  Pout<< "volPointInterpolation::calcBoundaryAddressing():"
126  << " added dangling mesh point:" << pointi
127  << " at:" << mesh().points()[pointi]
128  << endl;
129  }
130  }
131 
132  label nPatchFace = 0;
133  forAll(boundaryIsPatchFace_, i)
134  {
135  if (boundaryIsPatchFace_[i])
136  {
137  nPatchFace++;
138  }
139  }
140  label nPatchPoint = 0;
141  forAll(isPatchPoint_, i)
142  {
143  if (isPatchPoint_[i])
144  {
145  nPatchPoint++;
146  }
147  }
148  Pout<< "boundary:" << nl
149  << " faces :" << boundary.size() << nl
150  << " of which on proper patch:" << nPatchFace << nl
151  << " points:" << boundary.nPoints() << nl
152  << " of which on proper patch:" << nPatchPoint << endl;
153  }
154 }
155 
156 
157 void Foam::volPointInterpolation::makeInternalWeights(scalarField& sumWeights)
158 {
159  if (debug)
160  {
161  Pout<< "volPointInterpolation::makeInternalWeights() : "
162  << "constructing weighting factors for internal and non-coupled"
163  << " points." << endl;
164  }
165 
166  const pointField& points = mesh().points();
167  const labelListList& pointCells = mesh().pointCells();
168  const vectorField& cellCentres = mesh().cellCentres();
169 
170  // Allocate storage for weighting factors
171  pointWeights_.clear();
172  pointWeights_.setSize(points.size());
173 
174  // Calculate inverse distances between cell centres and points
175  // and store in weighting factor array
176  forAll(points, pointi)
177  {
178  if (!isPatchPoint_[pointi])
179  {
180  const labelList& pcp = pointCells[pointi];
181 
182  scalarList& pw = pointWeights_[pointi];
183  pw.setSize(pcp.size());
184 
185  forAll(pcp, pointCelli)
186  {
187  pw[pointCelli] =
188  1.0/mag(points[pointi] - cellCentres[pcp[pointCelli]]);
189 
190  sumWeights[pointi] += pw[pointCelli];
191  }
192  }
193  }
194 }
195 
196 
197 void Foam::volPointInterpolation::makeBoundaryWeights(scalarField& sumWeights)
198 {
199  if (debug)
200  {
201  Pout<< "volPointInterpolation::makeBoundaryWeights() : "
202  << "constructing weighting factors for boundary points." << endl;
203  }
204 
205  const pointField& points = mesh().points();
206  const pointField& faceCentres = mesh().faceCentres();
207 
208  const primitivePatch& boundary = boundaryPtr_();
209 
210  boundaryPointWeights_.clear();
211  boundaryPointWeights_.setSize(boundary.meshPoints().size());
212 
213  forAll(boundary.meshPoints(), i)
214  {
215  label pointi = boundary.meshPoints()[i];
216 
217  if (isPatchPoint_[pointi])
218  {
219  const labelList& pFaces = boundary.pointFaces()[i];
220 
221  scalarList& pw = boundaryPointWeights_[i];
222  pw.setSize(pFaces.size());
223 
224  sumWeights[pointi] = 0.0;
225 
226  forAll(pFaces, i)
227  {
228  if (boundaryIsPatchFace_[pFaces[i]])
229  {
230  label facei = mesh().nInternalFaces() + pFaces[i];
231 
232  pw[i] = 1.0/mag(points[pointi] - faceCentres[facei]);
233  sumWeights[pointi] += pw[i];
234  }
235  else
236  {
237  pw[i] = 0.0;
238  }
239  }
240  }
241  }
242 }
243 
244 
245 void Foam::volPointInterpolation::makeWeights()
246 {
247  if (debug)
248  {
249  Pout<< "volPointInterpolation::makeWeights() : "
250  << "constructing weighting factors"
251  << endl;
252  }
253 
254  // Update addressing over all boundary faces
255  calcBoundaryAddressing();
256 
257 
258  // Running sum of weights
259  pointScalarField sumWeights
260  (
261  IOobject
262  (
263  "volPointSumWeights",
265  mesh()
266  ),
267  pointMesh::New(mesh()),
269  );
270 
271 
272  // Create internal weights; add to sumWeights
273  makeInternalWeights(sumWeights);
274 
275 
276  // Create boundary weights; override sumWeights
277  makeBoundaryWeights(sumWeights);
278 
279 
280  //forAll(boundary.meshPoints(), i)
281  //{
282  // label pointi = boundary.meshPoints()[i];
283  //
284  // if (isPatchPoint_[pointi])
285  // {
286  // Pout<< "Calculated Weight at boundary point:" << i
287  // << " at:" << mesh().points()[pointi]
288  // << " sumWeight:" << sumWeights[pointi]
289  // << " from:" << boundaryPointWeights_[i]
290  // << endl;
291  // }
292  //}
293 
294 
295  // Sum collocated contributions
297  (
298  mesh(),
299  sumWeights,
300  plusEqOp<scalar>()
301  );
302 
303  // And add separated contributions
304  addSeparated(sumWeights);
305 
306  // Push master data to slaves. It is possible (not sure how often) for
307  // a coupled point to have its master on a different patch so
308  // to make sure just push master data to slaves. Reuse the syncPointData
309  // structure.
310  pushUntransformedData(sumWeights);
311 
312 
313  // Normalise internal weights
314  forAll(pointWeights_, pointi)
315  {
316  scalarList& pw = pointWeights_[pointi];
317  // Note:pw only sized for !isPatchPoint
318  forAll(pw, i)
319  {
320  pw[i] /= sumWeights[pointi];
321  }
322  }
323 
324  // Normalise boundary weights
325  const primitivePatch& boundary = boundaryPtr_();
326 
327  forAll(boundary.meshPoints(), i)
328  {
329  label pointi = boundary.meshPoints()[i];
330 
331  scalarList& pw = boundaryPointWeights_[i];
332  // Note:pw only sized for isPatchPoint
333  forAll(pw, i)
334  {
335  pw[i] /= sumWeights[pointi];
336  }
337  }
338 
339 
340  if (debug)
341  {
342  Pout<< "volPointInterpolation::makeWeights() : "
343  << "finished constructing weighting factors"
344  << endl;
345  }
346 }
347 
348 
349 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
350 
351 Foam::volPointInterpolation::volPointInterpolation(const fvMesh& vm)
352 :
354 {
355  makeWeights();
356 }
357 
358 
359 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
360 
362 {}
363 
364 
365 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
366 
368 {
369  makeWeights();
370 }
371 
372 
374 {
375  makeWeights();
376 
377  return true;
378 }
379 
380 
382 (
383  const volVectorField& vf,
384  pointVectorField& pf
385 ) const
386 {
387  interpolateInternalField(vf, pf);
388 
389  // Interpolate to the patches but no constraints
390  interpolateBoundaryField(vf, pf);
391 
392  // Apply displacement constraints
393  const pointConstraints& pcs = pointConstraints::New(pf.mesh());
394 
395  pcs.constrainDisplacement(pf, false);
396 }
397 
398 
399 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
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:1038
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::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.
Definition: zero.H:128
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:167
Foam::MeshObject< polyMesh, UpdateableMeshObject, pointMesh >::New
static const pointMesh & New(const polyMesh &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:72
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
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
An Ostream wrapper for parallel output to std::cout.
Foam::UpdateableMeshObject
Definition: MeshObject.H:241
Foam::volPointInterpolation::updateMesh
void updateMesh(const mapPolyMesh &)
Update mesh topology using the morph engine.
Definition: volPointInterpolation.C:367
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
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:347
Foam::MeshObject< fvMesh, UpdateableMeshObject, volPointInterpolation >::mesh
const fvMesh & mesh() const
Definition: MeshObject.H:122
pointConstraints.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::volPointInterpolation::~volPointInterpolation
~volPointInterpolation()
Destructor.
Definition: volPointInterpolation.C:361
Foam::pointScalarField
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
Definition: pointFields.H:51
Foam::pointConstraints
Application of (multi-)patch point contraints.
Definition: pointConstraints.H:64
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:176
Foam::nl
constexpr char nl
Definition: Ostream.H:372
f
labelList f(nPoints)
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
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:373
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:145
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:382
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
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< face, SubList, const pointField & > primitivePatch
Addressing for a faceList slice.
Definition: primitivePatch.H:47
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