volPointInterpolation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "volPointInterpolation.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "pointFields.H"
30 #include "demandDrivenData.H"
31 #include "pointConstraints.H"
32 #include "surfaceFields.H"
33 #include "syncTools.H"
34 #include "cpuTime.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
41 }
42 
43 
44 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
45 
47 :
49  <
50  fvMesh,
53  >(vm),
54  pointWeights_(mesh().points().size()),
55  boundary_
56  (
57  SubList<face>
58  (
59  mesh().faces(),
60  mesh().nFaces() - mesh().nInternalFaces(),
61  mesh().nInternalFaces()
62  ),
63  mesh().points()
64  ),
65  boundaryPointWeights_(boundary_.meshPoints().size()),
66  boundaryPointNbrWeights_(boundary_.meshPoints().size())
67 {
68  if (debug)
69  {
70  Pout<< "volPointInterpolation::volPointInterpolation(const fvMesh&) : "
71  << "constructing weighting factors"
72  << endl;
73  }
74 
75  const pointField& points = mesh().points();
76  const labelListList& pointCells = mesh().pointCells();
77  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
78  const fvBoundaryMesh& fvbm = mesh().boundary();
79 
80  // Cache calls to patch coupled flags
81  boolList isCoupledPolyPatch(pbm.size(), false);
82  boolList isCoupledFvPatch(fvbm.size(), false);
83  forAll(isCoupledFvPatch, patchi)
84  {
85  isCoupledPolyPatch[patchi] = pbm[patchi].coupled();
86  isCoupledFvPatch[patchi] = fvbm[patchi].coupled();
87  }
88 
89  // Determine the factor to which a point has its values set by adjacent
90  // boundary faces, rather than connected cells
91  scalarList pointBoundaryFactor(mesh().nPoints(), Zero);
92  forAll(boundary_.meshPoints(), bPointi)
93  {
94  const label pointi = boundary_.meshPoints()[bPointi];
95 
96  const labelList& pFaces = boundary_.pointFaces()[bPointi];
97 
98  forAll(pFaces, pPointFacei)
99  {
100  // Poly indices
101  const label patchi = pbm.patchIndices()[pFaces[pPointFacei]];
102  const label patchFacei =
103  pbm.patchFaceIndices()[pFaces[pPointFacei]];
104 
105  // FV indices
106  const labelUList patches =
107  mesh().polyBFacePatches()[pFaces[pPointFacei]];
108  const labelUList patchFaces =
109  mesh().polyBFacePatchFaces()[pFaces[pPointFacei]];
110 
111  scalar nonCoupledMagSf = 0;
112  forAll(patches, i)
113  {
114  if
115  (
116  !isCoupledPolyPatch[patches[i]]
117  && !isCoupledFvPatch[patches[i]]
118  )
119  {
120  nonCoupledMagSf += fvbm[patches[i]].magSf()[patchFaces[i]];
121  }
122  }
123 
124  pointBoundaryFactor[pointi] =
125  max
126  (
127  pointBoundaryFactor[pointi],
128  nonCoupledMagSf/pbm[patchi].magFaceAreas()[patchFacei]
129  );
130  }
131  }
133  (
134  mesh(),
135  pointBoundaryFactor,
136  maxEqOp<scalar>(),
137  scalar(0)
138  );
139 
140  // Calculate inverse distances between cell centres and points
141  // and store in the weighting factor array
142  forAll(points, pointi)
143  {
144  if (pointBoundaryFactor[pointi] > 1 - rootSmall) continue;
145 
146  pointWeights_[pointi].setSize(pointCells[pointi].size());
147 
148  const scalar f = pointBoundaryFactor[pointi];
149 
150  forAll(pointCells[pointi], pointCelli)
151  {
152  const label celli = pointCells[pointi][pointCelli];
153 
154  pointWeights_[pointi][pointCelli] =
155  (1 - f)/mag(points[pointi] - mesh().C()[celli]);
156  }
157  }
158 
159  // Get the cell centres on the other side of coupled boundaries
160  typename volVectorField::Boundary CBnf
161  (
162  mesh().boundary(),
165  );
166  forAll(fvbm, patchi)
167  {
168  if
169  (
170  !isCoupledPolyPatch[patchi]
171  && isCoupledFvPatch[patchi]
172  )
173  {
174  CBnf[patchi] = fvbm[patchi].Cn() + fvbm[patchi].delta();
175  }
176  }
177 
178  // Calculate inverse distanced between boundary face centres and points and
179  // store in the boundary weighting factor arrays
180  forAll(boundary_.meshPoints(), bPointi)
181  {
182  const label pointi = boundary_.meshPoints()[bPointi];
183 
184  const labelList& pFaces = boundary_.pointFaces()[bPointi];
185 
186  boundaryPointWeights_[bPointi].resize(pFaces.size());
187  boundaryPointNbrWeights_[bPointi].resize(pFaces.size());
188 
189  forAll(pFaces, bPointFacei)
190  {
191  // Poly indices
192  const label patchi = pbm.patchIndices()[pFaces[bPointFacei]];
193  const label patchFacei =
194  pbm.patchFaceIndices()[pFaces[bPointFacei]];
195 
196  // FV indices
197  const labelUList patches =
198  mesh().polyBFacePatches()[pFaces[bPointFacei]];
199  const labelUList patchFaces =
200  mesh().polyBFacePatchFaces()[pFaces[bPointFacei]];
201 
202  boundaryPointWeights_[bPointi][bPointFacei].resize
203  (
204  patches.size(),
205  0
206  );
207  boundaryPointNbrWeights_[bPointi][bPointFacei].resize
208  (
209  patches.size(),
210  0
211  );
212 
213  forAll(patches, i)
214  {
215  const scalar a =
216  fvbm[patches[i]].magSf()[patchFaces[i]]
217  /pbm[patchi].magFaceAreas()[patchFacei];
218 
219  const scalar f = pointBoundaryFactor[pointi];
220 
221  // If FV coupled only, add a weight to the neighbouring cell.
222  // This is necessary because point synchronisation will not sum
223  // the contributions across this interface as would be the case
224  // with a poly coupled interface.
225  if
226  (
227  !isCoupledPolyPatch[patches[i]]
228  && isCoupledFvPatch[patches[i]]
229  )
230  {
231  const point C = CBnf[patches[i]][patchFaces[i]];
232  boundaryPointNbrWeights_[bPointi][bPointFacei][i] =
233  (1 - f)*a/mag(points[pointi] - C);
234  }
235 
236  // If not coupled, add a weight to the boundary value
237  if
238  (
239  !isCoupledPolyPatch[patches[i]]
240  && !isCoupledFvPatch[patches[i]]
241  )
242  {
243  const point Cf = fvbm[patches[i]].Cf()[patchFaces[i]];
244  boundaryPointWeights_[bPointi][bPointFacei][i] =
245  f*a/mag(points[pointi] - Cf);
246  }
247  }
248  }
249  }
250 
251  // Construct a sum of weights
252  pointScalarField sumWeights
253  (
254  IOobject
255  (
256  "volPointSumWeights",
258  mesh()
259  ),
260  pointMesh::New(mesh()),
262  );
263 
264  // Add the internal weights
265  forAll(pointWeights_, pointi)
266  {
267  forAll(pointWeights_[pointi], i)
268  {
269  sumWeights[pointi] += pointWeights_[pointi][i];
270  }
271  }
272 
273  // Add the boundary weights
274  forAll(boundary_.meshPoints(), bPointi)
275  {
276  const label pointi = boundary_.meshPoints()[bPointi];
277  forAll(boundaryPointWeights_[bPointi], i)
278  {
279  forAll(boundaryPointWeights_[bPointi][i], j)
280  {
281  sumWeights[pointi] += boundaryPointWeights_[bPointi][i][j];
282  sumWeights[pointi] += boundaryPointNbrWeights_[bPointi][i][j];
283  }
284  }
285  }
286 
287  // Synchronise over conformal couplings
289  (
290  mesh(),
291  sumWeights,
293  scalar(0)
294  );
295 
296  // Normalise internal weights
297  forAll(pointWeights_, pointi)
298  {
299  forAll(pointWeights_[pointi], i)
300  {
301  pointWeights_[pointi][i] /= sumWeights[pointi];
302  }
303  }
304 
305  // Normalise boundary weights
306  forAll(boundary_.meshPoints(), bPointi)
307  {
308  const label pointi = boundary_.meshPoints()[bPointi];
309  forAll(boundaryPointWeights_[bPointi], i)
310  {
311  forAll(boundaryPointWeights_[bPointi][i], j)
312  {
313  boundaryPointWeights_[bPointi][i][j] /= sumWeights[pointi];
314  boundaryPointNbrWeights_[bPointi][i][j] /= sumWeights[pointi];
315  }
316  }
317  }
318 }
319 
320 
321 // * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * * //
322 
324 {}
325 
326 
327 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
328 
330 (
331  const volVectorField& vf,
332  pointVectorField& pf
333 ) const
334 {
335  interpolateUnconstrained(vf, pf);
336 
337  // Apply displacement constraints
339 }
340 
341 
342 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Graphite solid properties.
Definition: C.H:51
MeshObject types:
Definition: MeshObjects.H:67
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
void setSize(const label)
Reset size of List.
Definition: List.C:281
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const labelListList & pointFaces() const
Return point-face addressing.
A List obtained as a section of another List.
Definition: SubList.H:56
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Foam::fvBoundaryMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const volVectorField & C() const
Return cell centres.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:911
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:987
void constrainDisplacement(pointVectorField &displacement, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints),.
Foam::polyBoundaryMesh.
const labelList & patchIndices() const
Boundary face patch indices.
const labelList & patchFaceIndices() const
Boundary face patch face indices.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
const labelListList & pointCells() const
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
Interpolate from cell centres to points (vertices) using inverse distance weighting.
volPointInterpolation(const fvMesh &)
Constructor given fvMesh.
void interpolateDisplacement(const volVectorField &, pointVectorField &) const
Interpolate from volField to pointField.
Template functions to aid in the implementation of demand driven data.
label patchi
const pointField & points
label nPoints
const fvPatchList & patches
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
faceListList boundary(nPatches)
labelList f(nPoints)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("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, &mergedCyclicPolyPatch::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:229
Foam::surfaceFields.