volSurfaceMapping.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 Wikki Ltd
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volSurfaceMapping.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 template<class Type>
35 (
37  Field<Type>& result
38 ) const
39 {
40  // The polyPatch/local-face for each of the faceLabels
41  const List<labelPair>& patchFaces = mesh_.whichPatchFaces();
42 
43  // FULLDEBUG: checkSize ?
44  // or simply result.resize(mesh_.nFaces());
45 
46  forAll(patchFaces, i)
47  {
48  const labelPair& patchAndFace = patchFaces[i];
49 
50  if (patchAndFace.first() >= 0)
51  {
52  result[i] = bfld[patchAndFace];
53  }
54  }
55 }
56 
57 
58 template<class Type>
60 (
61  const GeometricBoundaryField<Type, fvPatchField, volMesh>& bfld
62 ) const
63 {
64  auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
65  mapToSurface(bfld, tresult.ref());
66  return tresult;
67 }
68 
69 
70 template<class Type>
72 (
74  Field<Type>& result
75 ) const
76 {
77  mapToSurface(vfld.boundaryField(), result);
78 }
79 
80 
81 template<class Type>
83 (
85  Field<Type>& result
86 ) const
87 {
88  mapToSurface(tvf().boundaryField(), result);
89  tvf.clear();
90 }
91 
92 
93 template<class Type>
95 (
97 ) const
98 {
99  return mapToSurface(vfld.boundaryField());
100 }
101 
102 
103 template<class Type>
105 (
107 ) const
108 {
109  tmp<Field<Type>> tresult(mapToSurface(tvf().boundaryField()));
110  tvf.clear();
111  return tresult;
112 }
113 
114 
115 template<class Type>
117 (
119  Field<Type>& result
120 ) const
121 {
122  const auto& bfld = fld.boundaryField();
123 
124  // The polyPatch/local-face for each of the faceLabels
125  const List<labelPair>& patchFaces = mesh_.whichPatchFaces();
126 
127  // FULLDEBUG: checkSize ?
128  // or simply result.resize(mesh_.nFaces());
129 
130  forAll(patchFaces, i)
131  {
132  const labelPair& patchAndFace = patchFaces[i];
133 
134  if (patchAndFace.first() >= 0)
135  {
136  // Value from boundary
137  result[i] = bfld[patchAndFace];
138  }
139  else if (patchAndFace.second() >= 0)
140  {
141  // Value from internal
142  result[i] = fld[patchAndFace.second()];
143  }
144  }
145 }
146 
147 
148 template<class Type>
150 (
151  const GeometricField<Type, fvsPatchField, surfaceMesh>& fld
152 ) const
153 {
154  auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
155  mapToSurface(fld, tresult.ref());
156  return tresult;
157 }
158 
159 
160 template<class Type>
162 (
164  Field<Type>& result
165 ) const
166 {
167  mapToSurface(tsf(), result);
168  tsf.clear();
169 }
170 
171 
172 template<class Type>
174 (
176 ) const
177 {
178  tmp<Field<Type>> tresult(mapToSurface(tsf()));
179  tsf.clear();
180  return tresult;
181 }
182 
183 
184 template<class Type>
186 (
187  const UPtrList<Field<Type>>& patchFields,
188  Field<Type>& result
189 ) const
190 {
191  // The polyPatch/local-face for each of the faceLabels
192  const List<labelPair>& patchFaces = mesh_.whichPatchFaces();
193 
194  // FULLDEBUG: checkSize ?
195  // or simply result.resize(mesh_.nFaces());
196 
197  forAll(patchFaces, i)
198  {
199  const label patchi = patchFaces[i].first();
200  const label facei = patchFaces[i].second();
201 
202  const auto* pfld = patchFields.get(patchi);
203 
204  if (pfld)
205  {
206  result[i] = (*pfld)[facei];
207  }
208  }
209 }
210 
211 
212 template<class Type>
214 (
215  const PtrMap<Field<Type>>& patchFields,
216  Field<Type>& result
217 ) const
218 {
219  // The polyPatch/local-face for each of the faceLabels
220  const List<labelPair>& patchFaces = mesh_.whichPatchFaces();
221 
222  // FULLDEBUG: checkSize ?
223  // or simply result.resize(mesh_.nFaces());
224 
225  forAll(patchFaces, i)
226  {
227  const label patchi = patchFaces[i].first();
228  const label facei = patchFaces[i].second();
229 
230  const auto* pfld = patchFields.get(patchi);
231 
232  if (pfld)
233  {
234  result[i] = (*pfld)[facei];
235  }
236  }
237 }
238 
239 
240 template<class Type>
242 (
243  const UPtrList<Field<Type>>& patchFields
244 ) const
245 {
246  auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
247  mapToSurface(patchFields, tresult.ref());
248  return tresult;
249 }
250 
251 
252 template<class Type>
254 (
255  const PtrMap<Field<Type>>& patchFields
256 ) const
257 {
258  auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
259  mapToSurface(patchFields, tresult.ref());
260  return tresult;
261 }
262 
263 
264 template<class Type>
266 (
268  Field<Type>& result
269 ) const
270 {
271  PtrList<Field<Type>> patchFields;
272 
273  // All referenced polyPatches (sorted order)
274  const labelList& patches = mesh_.whichPolyPatches();
275 
276  if (!patches.empty())
277  {
278  // maxPolyPatch+1
279  patchFields.resize(patches.last()+1);
280  }
281 
282  // Populate patchInternalField
283  for (const label patchi : patches)
284  {
285  patchFields.set(patchi, bfld[patchi].patchInternalField());
286  }
288  mapToSurface(patchFields, result);
289 }
290 
291 
292 template<class Type>
294 (
296 ) const
297 {
298  auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
299 
300  mapInternalToSurface(bfld, tresult.ref());
302  return tresult;
303 }
304 
305 
306 template<class Type>
308 (
310  Field<Type>& result
311 ) const
312 {
313  mapInternalToSurface(vf.boundaryField(), result);
314 }
315 
316 
317 template<class Type>
319 (
321  Field<Type>& result
322 ) const
323 {
324  mapInternalToSurface(tvf().boundaryField(), result);
325  tvf.clear();
326 }
327 
328 
329 template<class Type>
331 (
333 ) const
334 {
335  return mapInternalToSurface(vfld.boundaryField());
336 }
337 
338 
339 template<class Type>
341 (
343 ) const
344 {
345  tmp<Field<Type>> tresult(mapInternalToSurface(tvf().boundaryField()));
346  tvf.clear();
347  return tresult;
348 }
349 
350 
351 template<class Type>
353 (
356  const label destPatchi
357 ) const
358 {
359  // The polyPatch/local-face for each of the faceLabels
360  const List<labelPair>& patchFaces = mesh_.whichPatchFaces();
361 
362  forAll(patchFaces, i)
363  {
364  const labelPair& patchAndFace = patchFaces[i];
365 
366  if
367  (
368  patchAndFace.first() >= 0
369  && (patchAndFace.first() == destPatchi || destPatchi < 0)
370  )
371  {
372  dest[patchAndFace] = af[i];
373  }
374  }
375 }
376 
377 
378 template<class Type>
380 (
381  const tmp<DimensionedField<Type, areaMesh>>& taf,
382  GeometricBoundaryField<Type, fvPatchField, volMesh>& dest,
383  const label destPatchi
384 ) const
385 {
386  mapToVolume(taf(), dest, destPatchi);
387  taf.clear();
388 }
389 
390 
391 template<class Type>
393 (
396  const label destPatchi
397 ) const
398 {
399  mapToVolume(af, dest.boundaryFieldRef(), destPatchi);
400 }
401 
402 
403 template<class Type>
405 (
408  const label destPatchi
409 ) const
410 {
411  mapToVolume(taf, dest.boundaryFieldRef(), destPatchi);
412 }
413 
414 
415 template<class Type>
417 (
420  const label destPatchi
421 ) const
422 {
423  mapToVolume(taf().internalField(), dest, destPatchi);
424  taf.clear();
425 }
426 
427 
428 template<class Type>
430 (
433  const label destPatchi
434 ) const
435 {
436  mapToVolume(taf().internalField(), dest.boundaryFieldRef(), destPatchi);
437  taf.clear();
438 }
439 
440 
441 template<class Type>
443 (
445  Field<Type>& dest,
446  const label destPatchi
447 ) const
448 {
449  // The polyPatch/local-face for each of the faceLabels
450  const List<labelPair>& patchFaces = mesh_.whichPatchFaces();
451 
452  forAll(patchFaces, i)
453  {
454  const label patchi = patchFaces[i].first();
455  const label facei = patchFaces[i].second();
456 
457  if (patchi >= 0 && patchi == destPatchi)
458  {
459  dest[facei] = af[i];
460  }
461  }
462 }
463 
464 
465 template<class Type>
467 (
468  const tmp<DimensionedField<Type, areaMesh>>& taf,
469  Field<Type>& dest,
470  const label destPatchi
471 ) const
472 {
473  mapToVolumePatch(taf(), dest, destPatchi);
474  taf.clear();
475 }
476 
477 
478 template<class Type>
480 (
482  Field<Type>& dest,
483  const label destPatchi
484 ) const
485 {
486  mapToVolumePatch(taf().internalField(), dest, destPatchi);
487  taf.clear();
488 }
489 
490 
491 // ************************************************************************* //
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
const T & first() const noexcept
Access the first element.
Definition: Pair.H:137
T & last()
Return reference to the last element of the list.
Definition: UPtrList.H:861
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
T & first()
Access first element of the list, position [0].
Definition: UList.H:853
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition: UList.H:782
A HashTable of pointers to objects of type <T> with a label key.
Definition: HashTableFwd.H:42
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Generic templated field type.
Definition: Field.H:62
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
vsm mapToVolume(Cs, Cvf.boundaryFieldRef())
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:95
Generic GeometricBoundaryField class.
Definition: areaFieldsFwd.H:46
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool empty() const noexcept
True if the list is empty (ie, size() is zero)
Definition: UPtrListI.H:99
void mapToVolume(const DimensionedField< Type, areaMesh > &, GeometricBoundaryField< Type, fvPatchField, volMesh > &dest, const label destPatchi=-1) const
Map area field to volume boundary field, optionally restricted to a single destination patch...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
void mapInternalToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map patch internal field to area field.
void mapToVolumePatch(const DimensionedField< Type, areaMesh > &af, Field< Type > &dest, const label destPatchi) const
Map area field to a volume boundary patch.
const polyBoundaryMesh & patches
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const List< labelPair > & whichPatchFaces() const
The polyPatch/local-face for each faceLabels()
Definition: faMeshI.H:145
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const T & second() const noexcept
Access the second element.
Definition: Pair.H:147
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127