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-2016 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 #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  ),
98  sMesh,
99  vf.dimensions(),
100  Field<Type>(vf.primitiveField(), cellMap),
101  patchFields
102  )
103  );
105 
106 
107  // 2. Change the fvPatchFields to the correct type using a mapper
108  // constructor (with reference to the now correct internal field)
109 
111  Boundary& bf = resF.boundaryFieldRef();
112 
113  forAll(bf, patchi)
114  {
115  if (patchMap[patchi] != -1)
116  {
117  // Construct addressing
118  const fvPatch& subPatch = sMesh.boundary()[patchi];
119  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
120  const label baseStart = basePatch.start();
121  const label baseSize = basePatch.size();
122 
123  labelList directAddressing(subPatch.size());
124 
125  forAll(directAddressing, i)
126  {
127  label baseFacei = faceMap[subPatch.start()+i];
128 
129  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
130  {
131  directAddressing[i] = baseFacei-baseStart;
132  }
133  else
134  {
135  // Mapped from internal face. Do what? Leave up to
136  // fvPatchField
137  directAddressing[i] = -1;
138  }
139  }
140 
141  bf.set
142  (
143  patchi,
145  (
146  vf.boundaryField()[patchMap[patchi]],
147  subPatch,
148  resF(),
149  directFvPatchFieldMapper(directAddressing)
150  )
151  );
152  }
153  }
154 
155  return tresF;
156 }
157 
158 
159 template<class Type>
161 (
163 ) const
164 {
165  return interpolate
166  (
167  vf,
168  subMesh(),
169  patchMap(),
170  cellMap(),
171  faceMap()
172  );
173 }
174 
175 
176 template<class Type>
178 (
180  const fvMesh& sMesh,
181  const labelList& patchMap,
182  const labelList& cellMap,
183  const labelList& faceMap,
184  const bool negateIfFlipped
185 )
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"+vf.name(),
230  sMesh.time().timeName(),
231  sMesh,
234  ),
235  sMesh,
236  vf.dimensions(),
238  (
239  vf.primitiveField(),
241  (
242  faceMap,
243  sMesh.nInternalFaces()
244  )
245  ),
246  patchFields
247  )
248  );
250 
251 
252  // 2. Change the fvsPatchFields to the correct type using a mapper
253  // constructor (with reference to the now correct internal field)
254 
256  Boundary& bf = resF.boundaryFieldRef();
257 
258  forAll(bf, patchi)
259  {
260  if (patchMap[patchi] != -1)
261  {
262  // Construct addressing
263  const fvPatch& subPatch = sMesh.boundary()[patchi];
264  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
265  const label baseStart = basePatch.start();
266  const label baseSize = basePatch.size();
267 
268  labelList directAddressing(subPatch.size());
269 
270  forAll(directAddressing, i)
271  {
272  label baseFacei = faceMap[subPatch.start()+i];
273 
274  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
275  {
276  directAddressing[i] = baseFacei-baseStart;
277  }
278  else
279  {
280  // Mapped from internal face. Do what? Leave up to
281  // patchField. This would require also to pass in
282  // original internal field so for now do as postprocessing
283  directAddressing[i] = -1;
284  }
285  }
286 
287  bf.set
288  (
289  patchi,
291  (
292  vf.boundaryField()[patchMap[patchi]],
293  subPatch,
294  resF(),
295  directFvPatchFieldMapper(directAddressing)
296  )
297  );
298 
299 
300  // Postprocess patch field for exposed faces
301 
302  fvsPatchField<Type>& pfld = bf[patchi];
303  const labelUList& fc = bf[patchi].patch().faceCells();
304  const labelList& own = vf.mesh().faceOwner();
305 
306  forAll(pfld, i)
307  {
308  label baseFacei = faceMap[subPatch.start()+i];
309  if (baseFacei < vf.primitiveField().size())
310  {
311  Type val = vf.internalField()[baseFacei];
312 
313  if (cellMap[fc[i]] == own[baseFacei] || !negateIfFlipped)
314  {
315  pfld[i] = val;
316  }
317  else
318  {
319  pfld[i] = flipOp()(val);
320  }
321  }
322  else
323  {
324  // Exposed face from other patch.
325  // Only possible in case of a coupled boundary
326  label patchi = vf.mesh().boundaryMesh().whichPatch
327  (
328  baseFacei
329  );
330  const fvPatch& otherPatch = vf.mesh().boundary()[patchi];
331  label patchFacei = otherPatch.patch().whichFace(baseFacei);
332  pfld[i] = vf.boundaryField()[patchi][patchFacei];
333  }
334  }
335  }
336  }
337 
338  return tresF;
339 }
340 
341 
342 template<class Type>
344 (
346  const bool negateIfFlipped
347 ) const
348 {
349  return interpolate
350  (
351  sf,
352  subMesh(),
353  patchMap(),
354  cellMap(),
355  faceMap(),
356  negateIfFlipped
357  );
358 }
359 
360 
361 template<class Type>
364 (
366  const pointMesh& sMesh,
367  const labelList& patchMap,
368  const labelList& pointMap
369 )
370 {
371  // 1. Create the complete field with dummy patch fields
372  PtrList<pointPatchField<Type>> patchFields(patchMap.size());
373 
374  forAll(patchFields, patchi)
375  {
376  // Set the first one by hand as it corresponds to the
377  // exposed internal faces. Additional interpolation can be put here
378  // as necessary.
379  if (patchMap[patchi] == -1)
380  {
381  patchFields.set
382  (
383  patchi,
385  (
386  sMesh.boundary()[patchi],
388  )
389  );
390  }
391  else
392  {
393  patchFields.set
394  (
395  patchi,
397  (
399  sMesh.boundary()[patchi],
401  )
402  );
403  }
404  }
405 
406  // Create the complete field from the pieces
408  (
410  (
411  IOobject
412  (
413  "subset"+vf.name(),
414  sMesh.time().timeName(),
415  sMesh.thisDb(),
418  ),
419  sMesh,
420  vf.dimensions(),
421  Field<Type>(vf.primitiveField(), pointMap),
422  patchFields
423  )
424  );
426 
427 
428  // 2. Change the pointPatchFields to the correct type using a mapper
429  // constructor (with reference to the now correct internal field)
430 
432  Boundary& bf = resF.boundaryFieldRef();
433 
434  forAll(bf, patchi)
435  {
436  // Set the first one by hand as it corresponds to the
437  // exposed internal faces. Additional interpolation can be put here
438  // as necessary.
439  if (patchMap[patchi] != -1)
440  {
441  // Construct addressing
442  const pointPatch& basePatch =
443  vf.mesh().boundary()[patchMap[patchi]];
444 
445  const labelList& meshPoints = basePatch.meshPoints();
446 
447  // Make addressing from mesh to patch point
448  Map<label> meshPointMap(2*meshPoints.size());
449  forAll(meshPoints, localI)
450  {
451  meshPointMap.insert(meshPoints[localI], localI);
452  }
453 
454  // Find which subpatch points originate from which patch point
455  const pointPatch& subPatch = sMesh.boundary()[patchi];
456  const labelList& subMeshPoints = subPatch.meshPoints();
457 
458  // If mapped from outside patch leave handling up to patchField
459  labelList directAddressing(subPatch.size(), -1);
460 
461  forAll(subMeshPoints, localI)
462  {
463  // Get mesh point on original mesh.
464  label meshPointi = pointMap[subMeshPoints[localI]];
465 
466  Map<label>::const_iterator iter = meshPointMap.find(meshPointi);
467 
468  if (iter != meshPointMap.end())
469  {
470  directAddressing[localI] = iter();
471  }
472  }
473 
474  bf.set
475  (
476  patchi,
478  (
479  vf.boundaryField()[patchMap[patchi]],
480  subPatch,
481  resF(),
482  directPointPatchFieldMapper(directAddressing)
483  )
484  );
485  }
486  }
487 
488  return tresF;
489 }
490 
491 
492 template<class Type>
494 (
496 ) const
497 {
498  return interpolate
499  (
500  sf,
501  pointMesh::New(subMesh()), // subsetted point mesh
502  patchMap(),
503  pointMap()
504  );
505 }
506 
507 
508 template<class Type>
510 (
512  const fvMesh& sMesh,
513  const labelList& cellMap
514 )
515 {
516  // Create the complete field from the pieces
518  (
520  (
521  IOobject
522  (
523  "subset"+df.name(),
524  sMesh.time().timeName(),
525  sMesh,
528  ),
529  sMesh,
530  df.dimensions(),
531  Field<Type>(df, cellMap)
532  )
533  );
534 
535  return tresF;
536 }
537 
538 
539 template<class Type>
541 (
543 ) const
544 {
545  return interpolate(df, subMesh(), cellMap());
546 }
547 
548 
549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
550 
551 } // End namespace Foam
552 
553 // ************************************************************************* //
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:428
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:291
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:110
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:163
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:644
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
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:57
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
This boundary condition is not designed to be evaluated; it is assmued that the value is assigned via...
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:118
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:61
virtual label size() const
Return size.
Definition: fvPatch.H:161
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
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
const Time & time() const
Return time.
Definition: IOobject.C:337
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
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:106
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:377
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545