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-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 "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  // Cache calls to patch coupled flags
59  boolList isCoupledPolyPatch(pbm.size(), false);
60  boolList isCoupledFvPatch(fvbm.size(), false);
61  forAll(isCoupledFvPatch, patchi)
62  {
63  isCoupledPolyPatch[patchi] = pbm[patchi].coupled();
64  isCoupledFvPatch[patchi] = fvbm[patchi].coupled();
65  }
66 
67  pf = Zero;
68 
69  // Interpolate from the cells
70  forAll(pointWeights_, pointi)
71  {
72  forAll(pointWeights_[pointi], pointCelli)
73  {
74  const label celli = pointCells[pointi][pointCelli];
75 
76  pf[pointi] += pointWeights_[pointi][pointCelli]*vf[celli];
77  }
78  }
79 
80  // Get the boundary neighbour field
81  const typename VolField<Type>::Boundary vfBnf
82  (
85  );
86 
87  // Interpolate from the boundary faces
88  forAll(boundaryPointWeights_, bPointi)
89  {
90  const label pointi = boundary_.meshPoints()[bPointi];
91 
92  const labelList& pFaces = boundary_.pointFaces()[bPointi];
93 
94  forAll(boundaryPointWeights_[bPointi], bPointFacei)
95  {
96  // FV indices
97  const labelUList patches =
98  mesh().polyBFacePatches()[pFaces[bPointFacei]];
99  const labelUList patchFaces =
100  mesh().polyBFacePatchFaces()[pFaces[bPointFacei]];
101 
102  forAll(boundaryPointWeights_[bPointi][bPointFacei], i)
103  {
104  // If FV coupled only, add the neighbouring cell value
105  if
106  (
107  !isCoupledPolyPatch[patches[i]]
108  && isCoupledFvPatch[patches[i]]
109  )
110  {
111  pf[pointi] +=
112  boundaryPointNbrWeights_[bPointi][bPointFacei][i]
113  *vfBnf[patches[i]][patchFaces[i]];
114  }
115 
116  // If not coupled, add a weight to the boundary value
117  if
118  (
119  !isCoupledPolyPatch[patches[i]]
120  && !isCoupledFvPatch[patches[i]]
121  )
122  {
123  pf[pointi] +=
124  boundaryPointWeights_[bPointi][bPointFacei][i]
125  *vf.boundaryField()[patches[i]][patchFaces[i]];
126  }
127  }
128  }
129  }
130 
131  // Synchronise
133 }
134 
135 
136 template<class Type>
138 (
139  const VolField<Type>& vf,
140  PointField<Type>& pf
141 ) const
142 {
143  interpolateUnconstrained(vf, pf);
144 
145  // Apply constraints
147 }
148 
149 
150 template<class Type>
153 (
154  const VolField<Type>& vf,
155  const word& name,
156  const bool cache
157 ) const
158 {
159  const pointMesh& pm = pointMesh::New(vf.mesh());
160  const objectRegistry& db = pm.thisDb();
161 
162  if (!cache || vf.mesh().changing())
163  {
164  // Delete any old occurrences to avoid double registration
165  if (db.objectRegistry::template foundObject<PointField<Type>>(name))
166  {
167  PointField<Type>& pf =
168  db.objectRegistry::template lookupObjectRef<PointField<Type>>
169  (
170  name
171  );
172 
173  if (pf.ownedByRegistry())
174  {
175  solution::cachePrintMessage("Deleting", name, vf);
176  pf.release();
177  delete &pf;
178  }
179  }
180 
181  tmp<PointField<Type>> tpf
182  (
184  (
185  name,
186  pm,
187  vf.dimensions()
188  )
189  );
190 
191  interpolate(vf, tpf.ref());
192 
193  return tpf;
194  }
195  else
196  {
197  if (!db.objectRegistry::template foundObject<PointField<Type>>(name))
198  {
199  solution::cachePrintMessage("Calculating and caching", name, vf);
200  tmp<PointField<Type>> tpf = interpolate(vf, name, false);
201  PointField<Type>* pfPtr = tpf.ptr();
202  regIOobject::store(pfPtr);
203  return *pfPtr;
204  }
205  else
206  {
207  PointField<Type>& pf =
208  db.objectRegistry::template lookupObjectRef<PointField<Type>>
209  (
210  name
211  );
212 
213  if (pf.upToDate(vf))
214  {
215  solution::cachePrintMessage("Reusing", name, vf);
216  return pf;
217  }
218  else
219  {
220  solution::cachePrintMessage("Deleting", name, vf);
221  pf.release();
222  delete &pf;
223 
224  solution::cachePrintMessage("Recalculating", name, vf);
225  tmp<PointField<Type>> tpf = interpolate(vf, name, false);
226 
227  solution::cachePrintMessage("Storing", name, vf);
228  PointField<Type>* pfPtr = tpf.ptr();
229  regIOobject::store(pfPtr);
230 
231  return *pfPtr;
232  }
233  }
234  }
235 }
236 
237 
238 template<class Type>
241 (
242  const VolField<Type>& vf
243 ) const
244 {
245  return interpolate(vf, "volPointInterpolate(" + vf.name() + ')', false);
246 }
247 
248 
249 template<class Type>
252 (
253  const tmp<VolField<Type>>& tvf
254 ) const
255 {
256  tmp<PointField<Type>> tpf =
257  interpolate(tvf());
258 
259  tvf.clear();
260 
261  return tpf;
262 }
263 
264 
265 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static pointConstraints & New(const word &name, const pointMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
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.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, pointPatchField, pointMesh > > New(const word &name, const Internal &, const PtrList< pointPatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
const labelListList & pointFaces() const
Return point-face addressing.
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: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
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:404
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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.