fvMeshSubsetInterpolate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2013 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"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
39 
40 template<class Type>
41 tmp<GeometricField<Type, fvPatchField, volMesh> > fvMeshSubset::interpolate
42 (
44  const fvMesh& sMesh,
45  const labelList& patchMap,
46  const labelList& cellMap,
47  const labelList& faceMap
48 )
49 {
50  // 1. Create the complete field with dummy patch fields
51  PtrList<fvPatchField<Type> > patchFields(patchMap.size());
52 
53  forAll(patchFields, patchI)
54  {
55  // Set the first one by hand as it corresponds to the
56  // exposed internal faces. Additional interpolation can be put here
57  // as necessary.
58  if (patchMap[patchI] == -1)
59  {
60  patchFields.set
61  (
62  patchI,
64  (
65  sMesh.boundary()[patchI],
67  )
68  );
69  }
70  else
71  {
72  patchFields.set
73  (
74  patchI,
76  (
78  sMesh.boundary()[patchI],
80  )
81  );
82  }
83  }
84 
86  (
88  (
89  IOobject
90  (
91  "subset"+vf.name(),
92  sMesh.time().timeName(),
93  sMesh,
96  ),
97  sMesh,
98  vf.dimensions(),
99  Field<Type>(vf.internalField(), cellMap),
100  patchFields
101  )
102  );
104 
105 
106  // 2. Change the fvPatchFields to the correct type using a mapper
107  // constructor (with reference to the now correct internal field)
108 
110  GeometricBoundaryField& bf = resF.boundaryField();
111 
112  forAll(bf, patchI)
113  {
114  if (patchMap[patchI] != -1)
115  {
116  // Construct addressing
117  const fvPatch& subPatch = sMesh.boundary()[patchI];
118  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchI]];
119  const label baseStart = basePatch.start();
120  const label baseSize = basePatch.size();
121 
122  labelList directAddressing(subPatch.size());
123 
124  forAll(directAddressing, i)
125  {
126  label baseFaceI = faceMap[subPatch.start()+i];
127 
128  if (baseFaceI >= baseStart && baseFaceI < baseStart+baseSize)
129  {
130  directAddressing[i] = baseFaceI-baseStart;
131  }
132  else
133  {
134  // Mapped from internal face. Do what? Leave up to
135  // fvPatchField
136  directAddressing[i] = -1;
137  }
138  }
139 
140  bf.set
141  (
142  patchI,
144  (
145  vf.boundaryField()[patchMap[patchI]],
146  subPatch,
147  resF.dimensionedInternalField(),
148  directFvPatchFieldMapper(directAddressing)
149  )
150  );
151  }
152  }
153 
154  return tresF;
155 }
156 
157 
158 template<class Type>
160 (
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>
177 (
179  const fvMesh& sMesh,
180  const labelList& patchMap,
181  const labelList& faceMap
182 )
183 {
184  // 1. Create the complete field with dummy patch fields
185  PtrList<fvsPatchField<Type> > patchFields(patchMap.size());
186 
187  forAll(patchFields, patchI)
188  {
189  // Set the first one by hand as it corresponds to the
190  // exposed internal faces. Additional interpolation can be put here
191  // as necessary.
192  if (patchMap[patchI] == -1)
193  {
194  patchFields.set
195  (
196  patchI,
198  (
199  sMesh.boundary()[patchI],
201  )
202  );
203  }
204  else
205  {
206  patchFields.set
207  (
208  patchI,
210  (
212  sMesh.boundary()[patchI],
214  )
215  );
216  }
217  }
218 
219  // Create the complete field from the pieces
221  (
223  (
224  IOobject
225  (
226  "subset"+vf.name(),
227  sMesh.time().timeName(),
228  sMesh,
231  ),
232  sMesh,
233  vf.dimensions(),
235  (
236  vf.internalField(),
238  (
239  faceMap,
240  sMesh.nInternalFaces()
241  )
242  ),
243  patchFields
244  )
245  );
247 
248 
249  // 2. Change the fvsPatchFields to the correct type using a mapper
250  // constructor (with reference to the now correct internal field)
251 
253  GeometricBoundaryField& bf = resF.boundaryField();
254 
255  forAll(bf, patchI)
256  {
257  if (patchMap[patchI] != -1)
258  {
259  // Construct addressing
260  const fvPatch& subPatch = sMesh.boundary()[patchI];
261  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchI]];
262  const label baseStart = basePatch.start();
263  const label baseSize = basePatch.size();
264 
265  labelList directAddressing(subPatch.size());
266 
267  forAll(directAddressing, i)
268  {
269  label baseFaceI = faceMap[subPatch.start()+i];
270 
271  if (baseFaceI >= baseStart && baseFaceI < baseStart+baseSize)
272  {
273  directAddressing[i] = baseFaceI-baseStart;
274  }
275  else
276  {
277  // Mapped from internal face. Do what? Leave up to
278  // patchField. This would require also to pass in
279  // original internal field so for now do as postprocessing
280  directAddressing[i] = -1;
281  }
282  }
283 
284  bf.set
285  (
286  patchI,
288  (
289  vf.boundaryField()[patchMap[patchI]],
290  subPatch,
292  directFvPatchFieldMapper(directAddressing)
293  )
294  );
295 
296 
297  // Postprocess patch field for exposed faces
298 
299  fvsPatchField<Type>& pfld = bf[patchI];
300 
301  forAll(pfld, i)
302  {
303  label baseFaceI = faceMap[subPatch.start()+i];
304  if (baseFaceI < vf.internalField().size())
305  {
306  // Exposed internal face
307  pfld[i] = vf.internalField()[baseFaceI];
308  }
309  else
310  {
311  // Exposed face from other patch.
312  // Only possible in case of a coupled boundary
313  label patchI = vf.mesh().boundaryMesh().whichPatch
314  (
315  baseFaceI
316  );
317  const fvPatch& otherPatch = vf.mesh().boundary()[patchI];
318  label patchFaceI = otherPatch.patch().whichFace(baseFaceI);
319  pfld[i] = vf.boundaryField()[patchI][patchFaceI];
320  }
321  }
322  }
323  }
324 
325  return tresF;
326 }
327 
328 
329 template<class Type>
331 (
333 ) const
334 {
335  return interpolate
336  (
337  sf,
338  subMesh(),
339  patchMap(),
340  faceMap()
341  );
342 }
343 
344 
345 template<class Type>
348 (
350  const pointMesh& sMesh,
351  const labelList& patchMap,
352  const labelList& pointMap
353 )
354 {
355  // 1. Create the complete field with dummy patch fields
356  PtrList<pointPatchField<Type> > patchFields(patchMap.size());
357 
358  forAll(patchFields, patchI)
359  {
360  // Set the first one by hand as it corresponds to the
361  // exposed internal faces. Additional interpolation can be put here
362  // as necessary.
363  if (patchMap[patchI] == -1)
364  {
365  patchFields.set
366  (
367  patchI,
369  (
370  sMesh.boundary()[patchI],
372  )
373  );
374  }
375  else
376  {
377  patchFields.set
378  (
379  patchI,
381  (
383  sMesh.boundary()[patchI],
385  )
386  );
387  }
388  }
389 
390  // Create the complete field from the pieces
392  (
394  (
395  IOobject
396  (
397  "subset"+vf.name(),
398  sMesh.time().timeName(),
399  sMesh.thisDb(),
402  ),
403  sMesh,
404  vf.dimensions(),
405  Field<Type>(vf.internalField(), pointMap),
406  patchFields
407  )
408  );
410 
411 
412  // 2. Change the pointPatchFields to the correct type using a mapper
413  // constructor (with reference to the now correct internal field)
414 
416  GeometricBoundaryField& bf = resF.boundaryField();
417 
418  forAll(bf, patchI)
419  {
420  // Set the first one by hand as it corresponds to the
421  // exposed internal faces. Additional interpolation can be put here
422  // as necessary.
423  if (patchMap[patchI] != -1)
424  {
425  // Construct addressing
426  const pointPatch& basePatch =
427  vf.mesh().boundary()[patchMap[patchI]];
428 
429  const labelList& meshPoints = basePatch.meshPoints();
430 
431  // Make addressing from mesh to patch point
432  Map<label> meshPointMap(2*meshPoints.size());
433  forAll(meshPoints, localI)
434  {
435  meshPointMap.insert(meshPoints[localI], localI);
436  }
437 
438  // Find which subpatch points originate from which patch point
439  const pointPatch& subPatch = sMesh.boundary()[patchI];
440  const labelList& subMeshPoints = subPatch.meshPoints();
441 
442  // If mapped from outside patch leave handling up to patchField
443  labelList directAddressing(subPatch.size(), -1);
444 
445  forAll(subMeshPoints, localI)
446  {
447  // Get mesh point on original mesh.
448  label meshPointI = pointMap[subMeshPoints[localI]];
449 
450  Map<label>::const_iterator iter = meshPointMap.find(meshPointI);
451 
452  if (iter != meshPointMap.end())
453  {
454  directAddressing[localI] = iter();
455  }
456  }
457 
458  bf.set
459  (
460  patchI,
462  (
463  vf.boundaryField()[patchMap[patchI]],
464  subPatch,
465  resF.dimensionedInternalField(),
466  directPointPatchFieldMapper(directAddressing)
467  )
468  );
469  }
470  }
471 
472  return tresF;
473 }
474 
475 
476 template<class Type>
478 (
480 ) const
481 {
482  return interpolate
483  (
484  sf,
485  pointMesh::New(subMesh()), // subsetted point mesh
486  patchMap(),
487  pointMap()
488  );
489 }
490 
491 
492 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
493 
494 } // End namespace Foam
495 
496 // ************************************************************************* //
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.
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:118
virtual label size() const =0
Return size.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
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.
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const labelList & meshPoints() const =0
Return mesh points.
bool set(const label) const
Is element set.
Definition: PtrListI.H:107
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
static const pointMesh & New(const polyMesh &mesh)
Foam::emptyFvsPatchField.
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
InternalField & internalField()
Return internal field.
DimensionedInternalField & dimensionedInternalField()
Return dimensioned internal field.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
const Mesh & mesh() const
Return mesh.
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Abstract base class for point-mesh patch fields.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
A calculated boundary condition for pointField.
#define forAll(list, i)
Definition: UList.H:421
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:106
const dimensionSet & dimensions() const
Return dimensions.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Pre-declare SubField and related Field type.
Definition: Field.H:57
const word & name() const
Return name.
Definition: IOobject.H:260
An empty boundary condition for pointField.
label nInternalFaces() const
Generic GeometricField class.
const Time & time() const
Return time.
Definition: IOobject.C:251
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
A List obtained as a section of another List.
Definition: SubList.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::calculatedFvsPatchField.
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:139
virtual label size() const
Return size.
Definition: fvPatch.H:161
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
A class for managing temporary objects.
Definition: PtrList.H:118
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:552