volPointInterpolationTemplates.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-2022 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 "volFields.H"
28 #include "pointFields.H"
29 #include "surfaceFields.H"
30 #include "emptyFvPatch.H"
31 #include "coupledPointPatchField.H"
32 #include "pointConstraints.H"
33 #include "UCompactListList.H"
34 #include "syncTools.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 template<class Type>
40 (
41  const VolField<Type>& vf,
43 ) const
44 {
45  if (debug)
46  {
47  Pout<< "volPointInterpolation::interpolateUnconstrained("
48  << "const VolField<Type>&, "
49  << "PointField<Type>&) : "
50  << "interpolating field from cells to points"
51  << endl;
52  }
53 
54  const labelListList& pointCells = mesh().pointCells();
55  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
56  const fvBoundaryMesh& fvbm = mesh().boundary();
57 
58  const primitivePatch& boundary = boundaryPtr_();
59 
60  // Cache calls to patch coupled flags
61  boolList isCoupledPolyPatch(pbm.size(), false);
62  boolList isCoupledFvPatch(fvbm.size(), false);
63  forAll(isCoupledFvPatch, patchi)
64  {
65  isCoupledPolyPatch[patchi] = pbm[patchi].coupled();
66  isCoupledFvPatch[patchi] = fvbm[patchi].coupled();
67  }
68 
69  pf = Zero;
70 
71  // Interpolate from the cells
72  forAll(pointWeights_, pointi)
73  {
74  forAll(pointWeights_[pointi], pointCelli)
75  {
76  const label celli = pointCells[pointi][pointCelli];
77 
78  pf[pointi] += pointWeights_[pointi][pointCelli]*vf[celli];
79  }
80  }
81 
82  // Get the boundary neighbour field
83  const typename VolField<Type>::Boundary vfBnf
84  (
87  );
88 
89  // Interpolate from the boundary faces
90  forAll(boundaryPointWeights_, bPointi)
91  {
92  const label pointi = boundary.meshPoints()[bPointi];
93 
94  const labelList& pFaces = boundary.pointFaces()[bPointi];
95 
96  forAll(boundaryPointWeights_[bPointi], bPointFacei)
97  {
98  // FV indices
99  const labelUList patches =
100  mesh().polyBFacePatches()[pFaces[bPointFacei]];
101  const labelUList patchFaces =
102  mesh().polyBFacePatchFaces()[pFaces[bPointFacei]];
103 
104  forAll(boundaryPointWeights_[bPointi][bPointFacei], i)
105  {
106  // If FV coupled only, add the neighbouring cell value
107  if
108  (
109  !isCoupledPolyPatch[patches[i]]
110  && isCoupledFvPatch[patches[i]]
111  )
112  {
113  pf[pointi] +=
114  boundaryPointNbrWeights_[bPointi][bPointFacei][i]
115  *vfBnf[patches[i]][patchFaces[i]];
116  }
117 
118  // If not coupled, add a weight to the boundary value
119  if
120  (
121  !isCoupledPolyPatch[patches[i]]
122  && !isCoupledFvPatch[patches[i]]
123  )
124  {
125  pf[pointi] +=
126  boundaryPointWeights_[bPointi][bPointFacei][i]
127  *vf.boundaryField()[patches[i]][patchFaces[i]];
128  }
129  }
130  }
131  }
132 
133  // Synchronise
135 }
136 
137 
138 template<class Type>
140 (
141  const VolField<Type>& vf,
142  PointField<Type>& pf
143 ) const
144 {
145  interpolateUnconstrained(vf, pf);
146 
147  // Apply constraints
149 }
150 
151 
152 template<class Type>
155 (
156  const VolField<Type>& vf,
157  const word& name,
158  const bool cache
159 ) const
160 {
161  const pointMesh& pm = pointMesh::New(vf.mesh());
162  const objectRegistry& db = pm.thisDb();
163 
164  if (!cache || vf.mesh().changing())
165  {
166  // Delete any old occurrences to avoid double registration
167  if (db.objectRegistry::template foundObject<PointField<Type>>(name))
168  {
169  PointField<Type>& pf =
170  db.objectRegistry::template lookupObjectRef<PointField<Type>>
171  (
172  name
173  );
174 
175  if (pf.ownedByRegistry())
176  {
177  solution::cachePrintMessage("Deleting", name, vf);
178  pf.release();
179  delete &pf;
180  }
181  }
182 
183  tmp<PointField<Type>> tpf
184  (
186  (
187  name,
188  pm,
189  vf.dimensions()
190  )
191  );
192 
193  interpolate(vf, tpf.ref());
194 
195  return tpf;
196  }
197  else
198  {
199  if (!db.objectRegistry::template foundObject<PointField<Type>>(name))
200  {
201  solution::cachePrintMessage("Calculating and caching", name, vf);
202  tmp<PointField<Type>> tpf = interpolate(vf, name, false);
203  PointField<Type>* pfPtr = tpf.ptr();
204  regIOobject::store(pfPtr);
205  return *pfPtr;
206  }
207  else
208  {
209  PointField<Type>& pf =
210  db.objectRegistry::template lookupObjectRef<PointField<Type>>
211  (
212  name
213  );
214 
215  if (pf.upToDate(vf))
216  {
217  solution::cachePrintMessage("Reusing", name, vf);
218  return pf;
219  }
220  else
221  {
222  solution::cachePrintMessage("Deleting", name, vf);
223  pf.release();
224  delete &pf;
225 
226  solution::cachePrintMessage("Recalculating", name, vf);
227  tmp<PointField<Type>> tpf = interpolate(vf, name, false);
228 
229  solution::cachePrintMessage("Storing", name, vf);
230  PointField<Type>* pfPtr = tpf.ptr();
231  regIOobject::store(pfPtr);
232 
233  return *pfPtr;
234  }
235  }
236  }
237 }
238 
239 
240 template<class Type>
243 (
244  const VolField<Type>& vf
245 ) const
246 {
247  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
248 }
249 
250 
251 template<class Type>
254 (
255  const tmp<VolField<Type>>& tvf
256 ) const
257 {
258  tmp<PointField<Type>> tpf =
259  interpolate(tvf());
260 
261  tvf.clear();
262 
263  return tpf;
264 }
265 
266 
267 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
Generic GeometricField class.
static tmp< GeometricField< Type, pointPatchField, pointMesh > > New(const word &name, const Internal &, const PtrList< pointPatchField< Type >> &)
Return a temporary field constructed from name,.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A list of faces which address into the list of points.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
Foam::fvBoundaryMesh.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
const UCompactListList< label > & polyBFacePatches() const
Return poly-bFace-patch addressing.
Definition: fvMesh.C:951
const UCompactListList< label > & polyBFacePatchFaces() const
Return poly-bFace-patch-face addressing.
Definition: fvMesh.C:1027
Traits class for primitives.
Definition: pTraits.H:53
void constrain(PointField< Type > &pf, const bool overrideValue=false) const
Apply boundary conditions (single-patch constraints) and.
Foam::polyBoundaryMesh.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
const labelListList & pointCells() const
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
static void cachePrintMessage(const char *message, const word &name, const FieldType &vf)
Helper for printing cache message.
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
A class for managing temporary objects.
Definition: tmp.H:55
tmp< PointField< Type > > interpolate(const VolField< Type > &) const
Interpolate volField using inverse distance weighting.
void interpolateUnconstrained(const VolField< Type > &, PointField< Type > &) const
Interpolate from volField to pointField.
label patchi
const fvPatchList & patches
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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:251
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
faceListList boundary(nPatches)
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.