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-2019 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 "emptyFvsPatchField.H"
28 #include "emptyPointPatchField.H"
29 #include "emptyFvPatchFields.H"
31 #include "directPointPatchFieldMapper.H"
32 #include "flipOp.H"
33 
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
40 
41 template<class Type>
42 tmp<GeometricField<Type, fvPatchField, volMesh>> fvMeshSubset::interpolate
43 (
45  const fvMesh& sMesh,
46  const labelList& patchMap,
47  const labelList& cellMap,
48  const labelList& faceMap
49 )
50 {
51  // 1. Create the complete field with dummy patch fields
52  PtrList<fvPatchField<Type>> patchFields(patchMap.size());
53 
54  forAll(patchFields, patchi)
55  {
56  // Set the first one by hand as it corresponds to the
57  // exposed internal faces. Additional interpolation can be put here
58  // as necessary.
59  if (patchMap[patchi] == -1)
60  {
61  patchFields.set
62  (
63  patchi,
65  (
66  sMesh.boundary()[patchi],
68  )
69  );
70  }
71  else
72  {
73  patchFields.set
74  (
75  patchi,
77  (
79  sMesh.boundary()[patchi],
81  )
82  );
83  }
84  }
85 
87  (
89  (
90  IOobject
91  (
92  "subset"+vf.name(),
93  sMesh.time().timeName(),
94  sMesh,
97  false
98  ),
99  sMesh,
100  vf.dimensions(),
101  Field<Type>(vf.primitiveField(), cellMap),
102  patchFields
103  )
104  );
106 
107 
108  // 2. Change the fvPatchFields to the correct type using a mapper
109  // constructor (with reference to the now correct internal field)
110 
112  Boundary& bf = resF.boundaryFieldRef();
113 
114  forAll(bf, patchi)
115  {
116  if (patchMap[patchi] != -1)
117  {
118  // Construct addressing
119  const fvPatch& subPatch = sMesh.boundary()[patchi];
120  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
121  const label baseStart = basePatch.start();
122  const label baseSize = basePatch.size();
123 
124  labelList directAddressing(subPatch.size());
125 
126  forAll(directAddressing, i)
127  {
128  label baseFacei = faceMap[subPatch.start()+i];
129 
130  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
131  {
132  directAddressing[i] = baseFacei-baseStart;
133  }
134  else
135  {
136  // Mapped from internal face. Do what? Leave up to
137  // fvPatchField
138  directAddressing[i] = -1;
139  }
140  }
141 
142  bf.set
143  (
144  patchi,
146  (
147  vf.boundaryField()[patchMap[patchi]],
148  subPatch,
149  resF(),
150  directFvPatchFieldMapper(directAddressing)
151  )
152  );
153  }
154  }
155 
156  return tresF;
157 }
158 
159 
160 template<class Type>
162 (
164 ) const
165 {
166  return interpolate
167  (
168  vf,
169  subMesh(),
170  patchMap(),
171  cellMap(),
172  faceMap()
173  );
174 }
175 
176 
177 template<class Type>
179 (
181  const fvMesh& sMesh,
182  const labelList& patchMap,
183  const labelList& cellMap,
184  const labelList& faceMap,
185  const bool negateIfFlipped
186 )
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"+vf.name(),
231  sMesh.time().timeName(),
232  sMesh,
235  false
236  ),
237  sMesh,
238  vf.dimensions(),
240  (
241  vf.primitiveField(),
243  (
244  faceMap,
245  sMesh.nInternalFaces()
246  )
247  ),
248  patchFields
249  )
250  );
252 
253 
254  // 2. Change the fvsPatchFields to the correct type using a mapper
255  // constructor (with reference to the now correct internal field)
256 
258  Boundary& bf = resF.boundaryFieldRef();
259 
260  forAll(bf, patchi)
261  {
262  if (patchMap[patchi] != -1)
263  {
264  // Construct addressing
265  const fvPatch& subPatch = sMesh.boundary()[patchi];
266  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
267  const label baseStart = basePatch.start();
268  const label baseSize = basePatch.size();
269 
270  labelList directAddressing(subPatch.size());
271 
272  forAll(directAddressing, i)
273  {
274  label baseFacei = faceMap[subPatch.start()+i];
275 
276  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
277  {
278  directAddressing[i] = baseFacei-baseStart;
279  }
280  else
281  {
282  // Mapped from internal face. Do what? Leave up to
283  // patchField. This would require also to pass in
284  // original internal field so for now do as postprocessing
285  directAddressing[i] = -1;
286  }
287  }
288 
289  bf.set
290  (
291  patchi,
293  (
294  vf.boundaryField()[patchMap[patchi]],
295  subPatch,
296  resF(),
297  directFvPatchFieldMapper(directAddressing)
298  )
299  );
300 
301 
302  // Postprocess patch field for exposed faces
303 
304  fvsPatchField<Type>& pfld = bf[patchi];
305  const labelUList& fc = bf[patchi].patch().faceCells();
306  const labelList& own = vf.mesh().faceOwner();
307 
308  forAll(pfld, i)
309  {
310  label baseFacei = faceMap[subPatch.start()+i];
311  if (baseFacei < vf.primitiveField().size())
312  {
313  Type val = vf.internalField()[baseFacei];
314 
315  if (cellMap[fc[i]] == own[baseFacei] || !negateIfFlipped)
316  {
317  pfld[i] = val;
318  }
319  else
320  {
321  pfld[i] = flipOp()(val);
322  }
323  }
324  else
325  {
326  // Exposed face from other patch.
327  // Only possible in case of a coupled boundary
328  label patchi = vf.mesh().boundaryMesh().whichPatch
329  (
330  baseFacei
331  );
332  const fvPatch& otherPatch = vf.mesh().boundary()[patchi];
333  label patchFacei = otherPatch.patch().whichFace(baseFacei);
334  pfld[i] = vf.boundaryField()[patchi][patchFacei];
335  }
336  }
337  }
338  }
339 
340  return tresF;
341 }
342 
343 
344 template<class Type>
346 (
348  const bool negateIfFlipped
349 ) const
350 {
351  return interpolate
352  (
353  sf,
354  subMesh(),
355  patchMap(),
356  cellMap(),
357  faceMap(),
358  negateIfFlipped
359  );
360 }
361 
362 
363 template<class Type>
366 (
368  const pointMesh& sMesh,
369  const labelList& patchMap,
370  const labelList& pointMap
371 )
372 {
373  // 1. Create the complete field with dummy patch fields
374  PtrList<pointPatchField<Type>> patchFields(patchMap.size());
375 
376  forAll(patchFields, patchi)
377  {
378  // Set the first one by hand as it corresponds to the
379  // exposed internal faces. Additional interpolation can be put here
380  // as necessary.
381  if (patchMap[patchi] == -1)
382  {
383  patchFields.set
384  (
385  patchi,
387  (
388  sMesh.boundary()[patchi],
390  )
391  );
392  }
393  else
394  {
395  patchFields.set
396  (
397  patchi,
399  (
401  sMesh.boundary()[patchi],
403  )
404  );
405  }
406  }
407 
408  // Create the complete field from the pieces
410  (
412  (
413  IOobject
414  (
415  "subset"+vf.name(),
416  sMesh.time().timeName(),
417  sMesh.thisDb(),
420  false
421  ),
422  sMesh,
423  vf.dimensions(),
424  Field<Type>(vf.primitiveField(), pointMap),
425  patchFields
426  )
427  );
429 
430 
431  // 2. Change the pointPatchFields to the correct type using a mapper
432  // constructor (with reference to the now correct internal field)
433 
435  Boundary& bf = resF.boundaryFieldRef();
436 
437  forAll(bf, patchi)
438  {
439  // Set the first one by hand as it corresponds to the
440  // exposed internal faces. Additional interpolation can be put here
441  // as necessary.
442  if (patchMap[patchi] != -1)
443  {
444  // Construct addressing
445  const pointPatch& basePatch =
446  vf.mesh().boundary()[patchMap[patchi]];
447 
448  const labelList& meshPoints = basePatch.meshPoints();
449 
450  // Make addressing from mesh to patch point
451  Map<label> meshPointMap(2*meshPoints.size());
452  forAll(meshPoints, localI)
453  {
454  meshPointMap.insert(meshPoints[localI], localI);
455  }
456 
457  // Find which subpatch points originate from which patch point
458  const pointPatch& subPatch = sMesh.boundary()[patchi];
459  const labelList& subMeshPoints = subPatch.meshPoints();
460 
461  // If mapped from outside patch leave handling up to patchField
462  labelList directAddressing(subPatch.size(), -1);
463 
464  forAll(subMeshPoints, localI)
465  {
466  // Get mesh point on original mesh.
467  label meshPointi = pointMap[subMeshPoints[localI]];
468 
469  Map<label>::const_iterator iter = meshPointMap.find(meshPointi);
470 
471  if (iter != meshPointMap.end())
472  {
473  directAddressing[localI] = iter();
474  }
475  }
476 
477  bf.set
478  (
479  patchi,
481  (
482  vf.boundaryField()[patchMap[patchi]],
483  subPatch,
484  resF(),
485  directPointPatchFieldMapper(directAddressing)
486  )
487  );
488  }
489  }
490 
491  return tresF;
492 }
493 
494 
495 template<class Type>
497 (
499 ) const
500 {
501  return interpolate
502  (
503  sf,
504  pointMesh::New(subMesh()), // subsetted point mesh
505  patchMap(),
506  pointMap()
507  );
508 }
509 
510 
511 template<class Type>
513 (
515  const fvMesh& sMesh,
516  const labelList& cellMap
517 )
518 {
519  // Create the complete field from the pieces
521  (
523  (
524  IOobject
525  (
526  "subset"+df.name(),
527  sMesh.time().timeName(),
528  sMesh,
531  ),
532  sMesh,
533  df.dimensions(),
534  Field<Type>(df, cellMap)
535  )
536  );
537 
538  return tresF;
539 }
540 
541 
542 template<class Type>
544 (
546 ) const
547 {
548  return interpolate(df, subMesh(), cellMap());
549 }
550 
551 
552 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
553 
554 } // End namespace Foam
555 
556 // ************************************************************************* //
Class containing functor to negate primitives. Dummy for all other types.
Definition: flipOp.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
const word & name() const
Return name.
Definition: IOobject.H:303
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:149
label nInternalFaces() const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
An empty boundary condition for pointField.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
Generic GeometricField class.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
A calculated boundary condition for pointField.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
static const pointMesh & New(const polyMesh &mesh)
Definition: MeshObject.C:44
const Internal & internalField() const
Return a const-reference to the dimensioned internal field.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
Abstract base class for point-mesh patch fields.
const dimensionSet & dimensions() const
Return dimensions.
A List obtained as a section of another List.
Definition: SubList.H:53
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
Pre-declare SubField and related Field type.
Definition: Field.H:56
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:137
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:115
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual label size() const
Return size.
Definition: fvPatch.H:155
const Mesh & mesh() const
Return mesh.
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
Foam::emptyFvsPatchField.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
const Time & time() const
Return time.
Definition: IOobject.C:328
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This boundary condition provides an &#39;empty&#39; condition for reduced dimensions cases, i.e. 1- and 2-D geometries. Apply this condition to patches whose normal is aligned to geometric directions that do not constitue solution directions.
Foam::calculatedFvsPatchField.
A class for managing temporary objects.
Definition: PtrList.H:53
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:103
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
virtual const labelList & meshPoints() const =0
Return mesh points.
virtual label size() const =0
Return size.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:383
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540