Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
33template<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
58template<class Type>
59Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface
60(
62) const
63{
64 auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
65 mapToSurface(bfld, tresult.ref());
66 return tresult;
67}
68
69
70template<class Type>
72(
74 Field<Type>& result
75) const
76{
77 mapToSurface(vfld.boundaryField(), result);
78}
79
80
81template<class Type>
83(
85 Field<Type>& result
86) const
88 mapToSurface(tvf().boundaryField(), result);
89 tvf.clear();
90}
91
92
93template<class Type>
95(
97) const
98{
99 return mapToSurface(vfld.boundaryField());
100}
101
102
103template<class Type>
105(
107) const
108{
109 tmp<Field<Type>> tresult(mapToSurface(tvf().boundaryField()));
110 tvf.clear();
111 return tresult;
112}
113
114
115template<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
148template<class Type>
149Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface
150(
152) const
153{
154 auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
155 mapToSurface(fld, tresult.ref());
156 return tresult;
157}
158
159
160template<class Type>
162(
164 Field<Type>& result
165) const
167 mapToSurface(tsf(), result);
168 tsf.clear();
169}
170
171
172template<class Type>
174(
176) const
177{
179 tsf.clear();
180 return tresult;
181}
182
183
184template<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
212template<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
240template<class Type>
241Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface
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
252template<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
264template<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());
287
288 mapToSurface(patchFields, result);
289}
290
291
292template<class Type>
294(
295 const GeometricBoundaryField<Type, fvPatchField, volMesh>& bfld
296) const
297{
298 auto tresult = tmp<Field<Type>>::New(mesh_.nFaces(), Zero);
299
300 mapInternalToSurface(bfld, tresult.ref());
301
302 return tresult;
303}
304
305
306template<class Type>
308(
310 Field<Type>& result
311) const
312{
314}
315
316
317template<class Type>
319(
321 Field<Type>& result
322) const
324 mapInternalToSurface(tvf().boundaryField(), result);
325 tvf.clear();
326}
327
328
329template<class Type>
331(
333) const
334{
335 return mapInternalToSurface(vfld.boundaryField());
336}
337
338
339template<class Type>
341(
343) const
344{
345 tmp<Field<Type>> tresult(mapInternalToSurface(tvf().boundaryField()));
346 tvf.clear();
347 return tresult;
348}
349
350
351template<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
378template<class Type>
380(
383 const label destPatchi
384) const
386 mapToVolume(taf(), dest, destPatchi);
387 taf.clear();
388}
389
390
391template<class Type>
393(
396 const label destPatchi
397) const
398{
399 mapToVolume(af, dest.boundaryFieldRef(), destPatchi);
400}
401
402
403template<class Type>
405(
408 const label destPatchi
409) const
410{
411 mapToVolume(taf, dest.boundaryFieldRef(), destPatchi);
412}
413
414
415template<class Type>
417(
420 const label destPatchi
421) const
423 mapToVolume(taf().internalField(), dest, destPatchi);
424 taf.clear();
425}
426
427
428template<class Type>
430(
433 const label destPatchi
434) const
436 mapToVolume(taf().internalField(), dest.boundaryFieldRef(), destPatchi);
437 taf.clear();
438}
439
440
441template<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
465template<class Type>
467(
469 Field<Type>& dest,
470 const label destPatchi
471) const
473 mapToVolumePatch(taf(), dest, destPatchi);
474 taf.clear();
475}
476
477
478template<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// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
Generic GeometricBoundaryField class.
Generic GeometricField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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:171
void resize(const label newLen)
Adjust size of PtrList.
Definition PtrList.C:124
A HashTable of pointers to objects of type <T> with a label key.
Definition PtrMap.H:49
T & first()
Access first element of the list, position [0].
Definition UList.H:957
bool get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition UList.H:868
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
A class for managing temporary objects.
Definition tmp.H:75
void mapInternalToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map patch internal field to area field.
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.
void mapToVolumePatch(const DimensionedField< Type, areaMesh > &af, Field< Type > &dest, const label destPatchi) const
Map area field to a volume boundary patch.
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
const polyBoundaryMesh & patches
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
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.
List< label > labelList
A List of labels.
Definition List.H:62
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
vsm mapToVolume(Cs, Cvf.boundaryFieldRef())
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299