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-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 "fvMeshSubset.H"
27 #include "surfaceFields.H"
28 #include "internalFvsPatchField.H"
30 #include "internalFvPatchFields.H"
31 #include "forwardFieldMapper.H"
32 #include "flipOp.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
39 (
40  const VolField<Type>& vf,
41  const fvMesh& sMesh,
42  const labelList& patchMap,
43  const labelList& cellMap,
44  const labelList& faceMap
45 )
46 {
47  // 1. Create the complete field with dummy patch fields
49 
50  forAll(patchFields, patchi)
51  {
52  // Set the first one by hand as it corresponds to the
53  // exposed internal faces. Additional interpolation can be put here
54  // as necessary.
55  if (patchMap[patchi] == -1)
56  {
57  patchFields.set
58  (
59  patchi,
61  (
62  sMesh.boundary()[patchi],
64  )
65  );
66  }
67  else
68  {
69  patchFields.set
70  (
71  patchi,
73  (
75  sMesh.boundary()[patchi],
77  )
78  );
79  }
80  }
81 
82  tmp<VolField<Type>> tresF
83  (
84  new VolField<Type>
85  (
86  IOobject
87  (
88  "subset"+vf.name(),
89  sMesh.time().name(),
90  sMesh,
93  false
94  ),
95  sMesh,
96  vf.dimensions(),
98  patchFields
99  )
100  );
101  VolField<Type>& resF = tresF.ref();
102 
103 
104  // 2. Change the fvPatchFields to the correct type using a mapper
105  // constructor (with reference to the now correct internal field)
106 
107  typename VolField<Type>::
108  Boundary& bf = resF.boundaryFieldRef();
109 
110  forAll(bf, patchi)
111  {
112  if (patchMap[patchi] != -1)
113  {
114  // Construct addressing
115  const fvPatch& subPatch = sMesh.boundary()[patchi];
116  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
117  const label baseStart = basePatch.start();
118  const label baseSize = basePatch.size();
119 
120  labelList directAddressing(subPatch.size());
121 
122  forAll(directAddressing, i)
123  {
124  const label baseFacei = faceMap[subPatch.start() + i];
125 
126  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
127  {
128  directAddressing[i] = baseFacei-baseStart;
129  }
130  else
131  {
132  // Mapped from internal face. Do what? Leave up to
133  // fvPatchField
134  directAddressing[i] = -1;
135  }
136  }
137 
138  bf.set
139  (
140  patchi,
142  (
144  subPatch,
145  resF(),
146  forwardFieldMapper(directAddressing)
147  )
148  );
149  }
150  }
151 
152  return tresF;
153 }
154 
155 
156 template<class Type>
159 (
160  const VolField<Type>& vf
161 ) const
162 {
163  return interpolate
164  (
165  vf,
166  subMesh(),
167  patchMap(),
168  cellMap(),
169  faceMap()
170  );
171 }
172 
173 
174 template<class Type>
177 (
178  const SurfaceField<Type>& sf,
179  const fvMesh& sMesh,
180  const labelList& patchMap,
181  const labelList& cellMap,
182  const labelList& faceMap
183 )
184 {
185  const bool negateIfFlipped = isFlux(sf);
186 
187  // 1. Create the complete field with dummy patch fields
188  PtrList<fvsPatchField<Type>> patchFields(patchMap.size());
189 
190  forAll(patchFields, patchi)
191  {
192  // Set the first one by hand as it corresponds to the
193  // exposed internal faces. Additional interpolation can be put here
194  // as necessary.
195  if (patchMap[patchi] == -1)
196  {
197  patchFields.set
198  (
199  patchi,
201  (
202  sMesh.boundary()[patchi],
204  )
205  );
206  }
207  else
208  {
209  patchFields.set
210  (
211  patchi,
213  (
215  sMesh.boundary()[patchi],
217  )
218  );
219  }
220  }
221 
222  // Create the complete field from the pieces
224  (
226  (
227  IOobject
228  (
229  "subset"+sf.name(),
230  sMesh.time().name(),
231  sMesh,
234  false
235  ),
236  sMesh,
237  sf.dimensions(),
239  (
240  sf.primitiveField(),
242  (
243  faceMap,
244  sMesh.nInternalFaces()
245  )
246  ),
247  patchFields
248  )
249  );
250  SurfaceField<Type>& resF = tresF.ref();
251 
252  // 2. Change the fvsPatchFields to the correct type using a mapper
253  // constructor (with reference to the now correct internal field)
254 
255  typename SurfaceField<Type>::
256  Boundary& bf = resF.boundaryFieldRef();
257 
258  forAll(bf, patchi)
259  {
260  const fvPatch& subPatch = sMesh.boundary()[patchi];
261 
262  labelList directAddressing(subPatch.size(), -1);
263 
264  if (patchMap[patchi] != -1)
265  {
266  // Construct addressing
267  const fvPatch& basePatch = sf.mesh().boundary()[patchMap[patchi]];
268  const label baseStart = basePatch.start();
269  const label baseSize = basePatch.size();
270 
271  forAll(directAddressing, i)
272  {
273  const label baseFacei = faceMap[subPatch.start() + i];
274 
275  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
276  {
277  directAddressing[i] = baseFacei-baseStart;
278  }
279  }
280 
281  bf.set
282  (
283  patchi,
285  (
286  sf.boundaryField()[patchMap[patchi]],
287  subPatch,
288  resF(),
289  forwardFieldMapper(directAddressing)
290  )
291  );
292  }
293 
294  // Map internal face values onto the patch elected to hold
295  // the exposed faces
296  fvsPatchField<Type>& pfld = bf[patchi];
297  const labelUList& fc = bf[patchi].patch().faceCells();
298  const labelList& own = sf.mesh().faceOwner();
299 
300  forAll(pfld, i)
301  {
302  const label baseFacei = faceMap[subPatch.start() + i];
303 
304  if (directAddressing[i] == -1)
305  {
306  if (sf.mesh().isInternalFace(baseFacei))
307  {
308  const Type val = sf.internalField()[baseFacei];
309 
310  if (cellMap[fc[i]] == own[baseFacei] || !negateIfFlipped)
311  {
312  pfld[i] = val;
313  }
314  else
315  {
316  pfld[i] = flipOp()(val);
317  }
318  }
319  else
320  {
321  const label basePatchi =
322  sf.mesh().boundaryMesh().patchIndices()
323  [baseFacei - sf.mesh().nInternalFaces()];
324  const label basePatchFacei =
325  sf.mesh().boundaryMesh()[basePatchi]
326  .whichFace(baseFacei);
327 
328  pfld[i] = sf.boundaryField()[basePatchi][basePatchFacei];
329  }
330  }
331  }
332  }
333 
334  return tresF;
335 }
336 
337 
338 template<class Type>
341 (
342  const SurfaceField<Type>& sf
343 ) const
344 {
345  return interpolate
346  (
347  sf,
348  subMesh(),
349  patchMap(),
350  cellMap(),
351  faceMap()
352  );
353 }
354 
355 
356 template<class Type>
359 (
360  const PointField<Type>& pf,
361  const pointMesh& sMesh,
362  const labelList& patchMap,
363  const labelList& pointMap
364 )
365 {
366  // 1. Create the complete field with dummy patch fields
367  PtrList<pointPatchField<Type>> patchFields(patchMap.size());
368 
369  forAll(patchFields, patchi)
370  {
371  // Set the first one by hand as it corresponds to the
372  // exposed internal faces. Additional interpolation can be put here
373  // as necessary.
374  if (patchMap[patchi] == -1)
375  {
376  patchFields.set
377  (
378  patchi,
380  (
381  sMesh.boundary()[patchi],
383  )
384  );
385  }
386  else
387  {
388  patchFields.set
389  (
390  patchi,
392  (
394  sMesh.boundary()[patchi],
396  )
397  );
398  }
399  }
400 
401  // Create the complete field from the pieces
402  tmp<PointField<Type>> tresF
403  (
404  new PointField<Type>
405  (
406  IOobject
407  (
408  "subset"+pf.name(),
409  sMesh.time().name(),
410  sMesh.thisDb(),
413  false
414  ),
415  sMesh,
416  pf.dimensions(),
417  Field<Type>(pf.primitiveField(), pointMap),
418  patchFields
419  )
420  );
421  PointField<Type>& resF = tresF.ref();
422 
423 
424  // 2. Change the pointPatchFields to the correct type using a mapper
425  // constructor (with reference to the now correct internal field)
426 
427  typename PointField<Type>::
428  Boundary& bf = resF.boundaryFieldRef();
429 
430  forAll(bf, patchi)
431  {
432  // Set the first one by hand as it corresponds to the
433  // exposed internal faces. Additional interpolation can be put here
434  // as necessary.
435  if (patchMap[patchi] != -1)
436  {
437  // Construct addressing
438  const pointPatch& basePatch =
439  pf.mesh().boundary()[patchMap[patchi]];
440 
441  const labelList& meshPoints = basePatch.meshPoints();
442 
443  // Make addressing from mesh to patch point
444  Map<label> meshPointMap(2*meshPoints.size());
445  forAll(meshPoints, localI)
446  {
447  meshPointMap.insert(meshPoints[localI], localI);
448  }
449 
450  // Find which subpatch points originate from which patch point
451  const pointPatch& subPatch = sMesh.boundary()[patchi];
452  const labelList& subMeshPoints = subPatch.meshPoints();
453 
454  // If mapped from outside patch leave handling up to patchField
455  labelList directAddressing(subPatch.size(), -1);
456 
457  forAll(subMeshPoints, localI)
458  {
459  // Get mesh point on original mesh.
460  label meshPointi = pointMap[subMeshPoints[localI]];
461 
462  Map<label>::const_iterator iter = meshPointMap.find(meshPointi);
463 
464  if (iter != meshPointMap.end())
465  {
466  directAddressing[localI] = iter();
467  }
468  }
469 
470  bf.set
471  (
472  patchi,
474  (
475  pf.boundaryField()[patchMap[patchi]],
476  subPatch,
477  resF(),
478  forwardFieldMapper(directAddressing)
479  )
480  );
481  }
482  }
483 
484  return tresF;
485 }
486 
487 
488 template<class Type>
491 (
492  const PointField<Type>& sf
493 ) const
494 {
495  return interpolate
496  (
497  sf,
498  pointMesh::New(subMesh()), // subsetted point mesh
499  patchMap(),
500  pointMap()
501  );
502 }
503 
504 
505 template<class Type>
508 (
510  const fvMesh& sMesh,
511  const labelList& cellMap
512 )
513 {
514  // Create the complete field from the pieces
516  (
518  (
519  IOobject
520  (
521  "subset"+df.name(),
522  sMesh.time().name(),
523  sMesh,
526  ),
527  sMesh,
528  df.dimensions(),
529  Field<Type>(df, cellMap)
530  )
531  );
532 
533  return tresF;
534 }
535 
536 
537 template<class Type>
540 (
542 ) const
543 {
544  return interpolate(df, subMesh(), cellMap());
545 }
546 
547 
548 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
static pointMesh & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
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:83
Generic GeometricField class.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive 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:167
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:62
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
Forward field mapper.
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:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:138
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:132
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
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:53
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:137
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:125
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.