fvMeshSubsetInterpolate.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 "fvMeshSubset.H"
27 #include "surfaceFields.H"
28 #include "internalFvsPatchField.H"
30 #include "internalFvPatchFields.H"
32 #include "directPointPatchFieldMapper.H"
33 #include "flipOp.H"
34 
35 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
36 
37 template<class Type>
40 (
41  const VolField<Type>& vf,
42  const fvMesh& sMesh,
43  const labelList& patchMap,
44  const labelList& cellMap,
45  const labelList& faceMap
46 )
47 {
48  // 1. Create the complete field with dummy patch fields
50 
51  forAll(patchFields, patchi)
52  {
53  // Set the first one by hand as it corresponds to the
54  // exposed internal faces. Additional interpolation can be put here
55  // as necessary.
56  if (patchMap[patchi] == -1)
57  {
58  patchFields.set
59  (
60  patchi,
62  (
63  sMesh.boundary()[patchi],
65  )
66  );
67  }
68  else
69  {
70  patchFields.set
71  (
72  patchi,
74  (
76  sMesh.boundary()[patchi],
78  )
79  );
80  }
81  }
82 
83  tmp<VolField<Type>> tresF
84  (
85  new VolField<Type>
86  (
87  IOobject
88  (
89  "subset"+vf.name(),
90  sMesh.time().name(),
91  sMesh,
94  false
95  ),
96  sMesh,
97  vf.dimensions(),
99  patchFields
100  )
101  );
102  VolField<Type>& resF = tresF.ref();
103 
104 
105  // 2. Change the fvPatchFields to the correct type using a mapper
106  // constructor (with reference to the now correct internal field)
107 
108  typename VolField<Type>::
109  Boundary& bf = resF.boundaryFieldRef();
110 
111  forAll(bf, patchi)
112  {
113  if (patchMap[patchi] != -1)
114  {
115  // Construct addressing
116  const fvPatch& subPatch = sMesh.boundary()[patchi];
117  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
118  const label baseStart = basePatch.start();
119  const label baseSize = basePatch.size();
120 
121  labelList directAddressing(subPatch.size());
122 
123  forAll(directAddressing, i)
124  {
125  const label baseFacei = faceMap[subPatch.start() + i];
126 
127  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
128  {
129  directAddressing[i] = baseFacei-baseStart;
130  }
131  else
132  {
133  // Mapped from internal face. Do what? Leave up to
134  // fvPatchField
135  directAddressing[i] = -1;
136  }
137  }
138 
139  bf.set
140  (
141  patchi,
143  (
145  subPatch,
146  resF(),
147  directFvPatchFieldMapper(directAddressing)
148  )
149  );
150  }
151  }
152 
153  return tresF;
154 }
155 
156 
157 template<class Type>
160 (
161  const VolField<Type>& vf
162 ) const
163 {
164  return interpolate
165  (
166  vf,
167  subMesh(),
168  patchMap(),
169  cellMap(),
170  faceMap()
171  );
172 }
173 
174 
175 template<class Type>
178 (
179  const SurfaceField<Type>& sf,
180  const fvMesh& sMesh,
181  const labelList& patchMap,
182  const labelList& cellMap,
183  const labelList& faceMap
184 )
185 {
186  const bool negateIfFlipped = isFlux(sf);
187 
188  // 1. Create the complete field with dummy patch fields
189  PtrList<fvsPatchField<Type>> patchFields(patchMap.size());
190 
191  forAll(patchFields, patchi)
192  {
193  // Set the first one by hand as it corresponds to the
194  // exposed internal faces. Additional interpolation can be put here
195  // as necessary.
196  if (patchMap[patchi] == -1)
197  {
198  patchFields.set
199  (
200  patchi,
202  (
203  sMesh.boundary()[patchi],
205  )
206  );
207  }
208  else
209  {
210  patchFields.set
211  (
212  patchi,
214  (
216  sMesh.boundary()[patchi],
218  )
219  );
220  }
221  }
222 
223  // Create the complete field from the pieces
225  (
227  (
228  IOobject
229  (
230  "subset"+sf.name(),
231  sMesh.time().name(),
232  sMesh,
235  false
236  ),
237  sMesh,
238  sf.dimensions(),
240  (
241  sf.primitiveField(),
243  (
244  faceMap,
245  sMesh.nInternalFaces()
246  )
247  ),
248  patchFields
249  )
250  );
251  SurfaceField<Type>& resF = tresF.ref();
252 
253  // 2. Change the fvsPatchFields to the correct type using a mapper
254  // constructor (with reference to the now correct internal field)
255 
256  typename SurfaceField<Type>::
257  Boundary& bf = resF.boundaryFieldRef();
258 
259  forAll(bf, patchi)
260  {
261  const fvPatch& subPatch = sMesh.boundary()[patchi];
262 
263  labelList directAddressing(subPatch.size(), -1);
264 
265  if (patchMap[patchi] != -1)
266  {
267  // Construct addressing
268  const fvPatch& basePatch = sf.mesh().boundary()[patchMap[patchi]];
269  const label baseStart = basePatch.start();
270  const label baseSize = basePatch.size();
271 
272  forAll(directAddressing, i)
273  {
274  const label baseFacei = faceMap[subPatch.start() + i];
275 
276  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
277  {
278  directAddressing[i] = baseFacei-baseStart;
279  }
280  }
281 
282  bf.set
283  (
284  patchi,
286  (
287  sf.boundaryField()[patchMap[patchi]],
288  subPatch,
289  resF(),
290  directFvPatchFieldMapper(directAddressing)
291  )
292  );
293  }
294 
295  // Map internal face values onto the patch elected to hold
296  // the exposed faces
297  fvsPatchField<Type>& pfld = bf[patchi];
298  const labelUList& fc = bf[patchi].patch().faceCells();
299  const labelList& own = sf.mesh().faceOwner();
300 
301  forAll(pfld, i)
302  {
303  const label baseFacei = faceMap[subPatch.start() + i];
304 
305  if (directAddressing[i] == -1)
306  {
307  if (sf.mesh().isInternalFace(baseFacei))
308  {
309  const Type val = sf.internalField()[baseFacei];
310 
311  if (cellMap[fc[i]] == own[baseFacei] || !negateIfFlipped)
312  {
313  pfld[i] = val;
314  }
315  else
316  {
317  pfld[i] = flipOp()(val);
318  }
319  }
320  else
321  {
322  const label basePatchi =
323  sf.mesh().boundaryMesh().patchID()
324  [baseFacei - sf.mesh().nInternalFaces()];
325  const label basePatchFacei =
326  sf.mesh().boundaryMesh()[basePatchi]
327  .whichFace(baseFacei);
328 
329  pfld[i] = sf.boundaryField()[basePatchi][basePatchFacei];
330  }
331  }
332  }
333  }
334 
335  return tresF;
336 }
337 
338 
339 template<class Type>
342 (
343  const SurfaceField<Type>& sf
344 ) const
345 {
346  return interpolate
347  (
348  sf,
349  subMesh(),
350  patchMap(),
351  cellMap(),
352  faceMap()
353  );
354 }
355 
356 
357 template<class Type>
360 (
361  const PointField<Type>& pf,
362  const pointMesh& sMesh,
363  const labelList& patchMap,
364  const labelList& pointMap
365 )
366 {
367  // 1. Create the complete field with dummy patch fields
368  PtrList<pointPatchField<Type>> patchFields(patchMap.size());
369 
370  forAll(patchFields, patchi)
371  {
372  // Set the first one by hand as it corresponds to the
373  // exposed internal faces. Additional interpolation can be put here
374  // as necessary.
375  if (patchMap[patchi] == -1)
376  {
377  patchFields.set
378  (
379  patchi,
381  (
382  sMesh.boundary()[patchi],
384  )
385  );
386  }
387  else
388  {
389  patchFields.set
390  (
391  patchi,
393  (
395  sMesh.boundary()[patchi],
397  )
398  );
399  }
400  }
401 
402  // Create the complete field from the pieces
403  tmp<PointField<Type>> tresF
404  (
405  new PointField<Type>
406  (
407  IOobject
408  (
409  "subset"+pf.name(),
410  sMesh.time().name(),
411  sMesh.thisDb(),
414  false
415  ),
416  sMesh,
417  pf.dimensions(),
418  Field<Type>(pf.primitiveField(), pointMap),
419  patchFields
420  )
421  );
422  PointField<Type>& resF = tresF.ref();
423 
424 
425  // 2. Change the pointPatchFields to the correct type using a mapper
426  // constructor (with reference to the now correct internal field)
427 
428  typename PointField<Type>::
429  Boundary& bf = resF.boundaryFieldRef();
430 
431  forAll(bf, patchi)
432  {
433  // Set the first one by hand as it corresponds to the
434  // exposed internal faces. Additional interpolation can be put here
435  // as necessary.
436  if (patchMap[patchi] != -1)
437  {
438  // Construct addressing
439  const pointPatch& basePatch =
440  pf.mesh().boundary()[patchMap[patchi]];
441 
442  const labelList& meshPoints = basePatch.meshPoints();
443 
444  // Make addressing from mesh to patch point
445  Map<label> meshPointMap(2*meshPoints.size());
446  forAll(meshPoints, localI)
447  {
448  meshPointMap.insert(meshPoints[localI], localI);
449  }
450 
451  // Find which subpatch points originate from which patch point
452  const pointPatch& subPatch = sMesh.boundary()[patchi];
453  const labelList& subMeshPoints = subPatch.meshPoints();
454 
455  // If mapped from outside patch leave handling up to patchField
456  labelList directAddressing(subPatch.size(), -1);
457 
458  forAll(subMeshPoints, localI)
459  {
460  // Get mesh point on original mesh.
461  label meshPointi = pointMap[subMeshPoints[localI]];
462 
463  Map<label>::const_iterator iter = meshPointMap.find(meshPointi);
464 
465  if (iter != meshPointMap.end())
466  {
467  directAddressing[localI] = iter();
468  }
469  }
470 
471  bf.set
472  (
473  patchi,
475  (
476  pf.boundaryField()[patchMap[patchi]],
477  subPatch,
478  resF(),
479  directPointPatchFieldMapper(directAddressing)
480  )
481  );
482  }
483  }
484 
485  return tresF;
486 }
487 
488 
489 template<class Type>
492 (
493  const PointField<Type>& sf
494 ) const
495 {
496  return interpolate
497  (
498  sf,
499  pointMesh::New(subMesh()), // subsetted point mesh
500  patchMap(),
501  pointMap()
502  );
503 }
504 
505 
506 template<class Type>
509 (
511  const fvMesh& sMesh,
512  const labelList& cellMap
513 )
514 {
515  // Create the complete field from the pieces
517  (
519  (
520  IOobject
521  (
522  "subset"+df.name(),
523  sMesh.time().name(),
524  sMesh,
527  ),
528  sMesh,
529  df.dimensions(),
530  Field<Type>(df, cellMap)
531  )
532  );
533 
534  return tresF;
535 }
536 
537 
538 template<class Type>
541 (
543 ) const
544 {
545  return interpolate(df, subMesh(), cellMap());
546 }
547 
548 
549 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:318
const word & name() const
Return name.
Definition: IOobject.H:310
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
A List obtained as a section of another List.
Definition: SubList.H:56
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
Foam::calculatedFvsPatchField.
A calculated boundary condition for pointField.
const word & name() const
Return const reference to name.
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:51
const labelList & faceMap() const
Return face map.
const labelList & cellMap() const
Return cell map.
const labelList & patchMap() const
Return patch map.
static tmp< VolField< Type > > interpolate(const VolField< Type > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:157
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:151
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
Constraint fvPatchField to hold values for internal face exposed by sub-setting.
Constraint fvsPatchField to hold values for internal face exposed by sub-setting.
Constraint pointPatchField to hold values for internal face exposed by sub-setting.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:116
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:104
Abstract base class for point-mesh patch fields.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:57
virtual const labelList & meshPoints() const =0
Return mesh points.
virtual label size() const =0
Return size.
label nInternalFaces() const
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
label patchi
volScalarField sf(fieldObject, mesh)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
Foam::surfaceFields.