Loading...
Searching...
No Matches
fvFieldDecomposer.H
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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2021-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
27Class
28 Foam::fvFieldDecomposer
29
30Description
31 Finite Volume volume and surface field decomposer.
32
33SourceFiles
34 fvFieldDecomposer.cxx
35 fvFieldDecomposer.txx
36 fvFieldDecomposerCache.cxx
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef Foam_fvFieldDecomposer_H
41#define Foam_fvFieldDecomposer_H
42
43#include "fvMesh.H"
44#include "fvPatchFieldMapper.H"
45#include "volFields.H"
46#include "surfaceFields.H"
47
48// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49
50namespace Foam
51{
52
53// Forward Declarations
54class IOobjectList;
56/*---------------------------------------------------------------------------*\
57 Class fvFieldDecomposer Declaration
58\*---------------------------------------------------------------------------*/
59
60class fvFieldDecomposer
61{
62public:
63
64 // Public Classes
65
66 //- Internal caching for field reading
67 class fieldsCache;
68
69 //- Patch field decomposer class
71 :
73 {
74 // Private Data
75
76 labelList directAddressing_;
77
78 public:
79
80 // Constructors
81
82 //- Construct given addressing
84 (
85 const labelUList& addressingSlice,
86 const label addressingOffset
87 );
88
89
90 // Member Functions
91
92 label size() const noexcept
93 {
94 return directAddressing_.size();
95 }
96
97 bool direct() const noexcept
98 {
99 return true;
100 }
101
102 //- Are there unmapped values
103 bool hasUnmapped() const noexcept
104 {
105 return false;
107
109 {
110 return directAddressing_;
112 };
113
114
115 //- Processor patch field decomposer class. Maps either owner or
116 // neighbour data (no interpolate anymore - processorFvPatchField
117 // holds neighbour data)
118 class processorVolPatchFieldDecomposer
119 :
120 public fvPatchFieldMapper
121 {
122 // Private Data
123
124 labelList directAddressing_;
125
126 public:
127
128 //- Construct addressing from details
130 (
131 const labelUList& faceOwner,
132 const labelUList& faceNeigbour,
133 const labelUList& addressingSlice
134 );
135
136 //- Construct given addressing from complete mesh
138 (
139 const polyMesh& mesh,
140 const labelUList& addressingSlice
141 );
142
143
144 // Member Functions
145
146 label size() const noexcept
148 return directAddressing_.size();
149 }
150
151 bool direct() const noexcept
152 {
153 return true;
154 }
155
156 //- Are there unmapped values
157 bool hasUnmapped() const noexcept
158 {
159 return false;
160 }
163 {
164 return directAddressing_;
165 }
166 };
167
168
169 //- Processor patch field decomposer class. Surface field is assumed
170 // to have direction (so manipulates sign when mapping)
172 :
173 public fvPatchFieldMapper
175 labelListList addressing_;
176 scalarListList weights_;
177
178 public:
179
180 //- Construct given addressing
182 (
183 const labelUList& addressingSlice
184 );
185
187 // Member Functions
188
189 label size() const noexcept
190 {
191 return addressing_.size();
192 }
193
194 bool direct() const noexcept
195 {
196 return false;
197 }
199 //- Are there unmapped values
200 bool hasUnmapped() const noexcept
201 {
202 return false;
203 }
204
205 const labelListList& addressing() const noexcept
207 return addressing_;
208 }
209
210 const scalarListList& weights() const noexcept
212 return weights_;
213 }
214 };
215
216
217private:
218
219 // Private Data
220
221 //- Reference to processor mesh
222 const fvMesh& procMesh_;
223
224 //- Reference to face addressing
225 const labelUList& faceAddressing_;
226
227 //- Reference to cell addressing
228 const labelUList& cellAddressing_;
230 //- Reference to boundary addressing
231 const labelUList& boundaryAddressing_;
232
233 //- List of patch field decomposers
234 PtrList<patchFieldDecomposer> patchFieldDecomposerPtrs_;
235
237 processorVolPatchFieldDecomposerPtrs_;
238
240 processorSurfacePatchFieldDecomposerPtrs_;
241
242 PtrList<scalarField> faceSign_;
243
244
245 // Private Member Functions
246
247 //- No copy construct
248 fvFieldDecomposer(const fvFieldDecomposer&) = delete;
249
250 //- No copy assignment
251 void operator=(const fvFieldDecomposer&) = delete;
252
253public:
254
255 // Static Data
256
257 //- Output verbosity when writing
258 static int verbose_;
259
260
261 // Constructors
262
263 //- Construct without mappers, added later with reset()
264 fvFieldDecomposer
265 (
267 const fvMesh& procMesh,
268 const labelUList& faceAddressing,
269 const labelUList& cellAddressing,
270 const labelUList& boundaryAddressing
271 );
272
273 //- Construct from components using information from the complete mesh
274 fvFieldDecomposer
275 (
276 const fvMesh& completeMesh,
277 const fvMesh& procMesh,
278 const labelUList& faceAddressing,
279 const labelUList& cellAddressing,
280 const labelUList& boundaryAddressing
281 );
282
283 //- Construct from components without the complete mesh
284 fvFieldDecomposer
285 (
286 // Information about the complete mesh
287 const UList<labelRange>& boundaryRanges,
288 const labelUList& faceOwner,
289 const labelUList& faceNeigbour,
290
291 // Addressing for processor mesh
292 const fvMesh& procMesh,
293 const labelUList& faceAddressing,
294 const labelUList& cellAddressing,
295 const labelUList& boundaryAddressing
296 );
297
298
299 //- Destructor
300 ~fvFieldDecomposer() = default;
302
303 // Member Functions
304
305 //- True if no mappers have been allocated
306 bool empty() const noexcept;
307
308 //- Remove all mappers
309 void clear();
310
311 //- Reset mappers using information from the complete mesh
312 void reset(const fvMesh& completeMesh);
314 //- Reset mapper using information about the complete mesh
315 void reset
316 (
317 const UList<labelRange>& boundaryRanges,
318 const labelUList& faceOwner,
319 const labelUList& faceNeigbour
320 );
321
322
323 // Mapping
324
325 //- Decompose internal field
326 template<class Type>
329 (
330 const DimensionedField<Type, volMesh>& field
331 ) const;
332
333 //- Decompose volume field
334 template<class Type>
337 (
339 const bool allowUnknownPatchFields = false
340 ) const;
341
342 //- Decompose surface field
343 template<class Type>
346 (
348 ) const;
349
350 //- Decompose list of fields
351 template<class GeoField>
352 void decomposeFields(const UPtrList<GeoField>& fields) const;
353};
354
355
356/*---------------------------------------------------------------------------*\
357 Class fvFieldDecomposer::fieldsCache Declaration
358\*---------------------------------------------------------------------------*/
359
360class fvFieldDecomposer::fieldsCache
362 class privateCache;
363 std::unique_ptr<privateCache> cache_;
364
365 //- No copy construct
366 fieldsCache(const fieldsCache&) = delete;
367
368 //- No copy assignment
369 void operator=(const fieldsCache&) = delete;
370
371
372public:
373
374 // Constructors
375
376 //- Default construct
377 fieldsCache();
378
379
380 //- Destructor
382
383
384 // Member Functions
385
386 //- No fields
387 bool empty() const noexcept;
388
389 //- Total number of fields
390 label size() const noexcept;
392 //- Clear out
393 void clear();
394
395
396 //- Read all fields given mesh and objects
397 void readAllFields
398 (
399 const fvMesh& mesh,
400 const IOobjectList& objects
401 );
403 //- Decompose and write all fields
404 void decomposeAllFields
405 (
406 const fvFieldDecomposer& decomposer,
407 bool report = false
408 ) const;
409};
410
412// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
413
414} // End namespace Foam
415
416// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
417
418#ifdef NoRepository
419 #include "fvFieldDecomposer.txx"
420#endif
421
422// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
423
424#endif
425
426// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
bool empty() const noexcept
No fields.
void readAllFields(const fvMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
void decomposeAllFields(const fvFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
label size() const noexcept
Total number of fields.
bool hasUnmapped() const noexcept
Are there unmapped values.
bool direct() const noexcept
Is it a direct (non-interpolating) mapper?
label size() const noexcept
The size of the mapper.
const labelUList & directAddressing() const noexcept
Return the direct addressing values.
patchFieldDecomposer(const labelUList &addressingSlice, const label addressingOffset)
Construct given addressing.
Processor patch field decomposer class. Surface field is assumed.
bool hasUnmapped() const noexcept
Are there unmapped values.
bool direct() const noexcept
Is it a direct (non-interpolating) mapper?
const scalarListList & weights() const noexcept
Return the interpolation weights.
const labelListList & addressing() const noexcept
Return the interpolation addressing.
label size() const noexcept
The size of the mapper.
processorSurfacePatchFieldDecomposer(const labelUList &addressingSlice)
Construct given addressing.
bool hasUnmapped() const noexcept
Are there unmapped values.
processorVolPatchFieldDecomposer(const labelUList &faceOwner, const labelUList &faceNeigbour, const labelUList &addressingSlice)
Construct addressing from details.
bool direct() const noexcept
Is it a direct (non-interpolating) mapper?
label size() const noexcept
The size of the mapper.
const labelUList & directAddressing() const noexcept
Return the direct addressing values.
processorVolPatchFieldDecomposer(const polyMesh &mesh, const labelUList &addressingSlice)
Construct given addressing from complete mesh.
Finite Volume volume and surface field decomposer.
tmp< DimensionedField< Type, volMesh > > decomposeField(const DimensionedField< Type, volMesh > &field) const
Decompose internal field.
fvFieldDecomposer(const UList< labelRange > &boundaryRanges, const labelUList &faceOwner, const labelUList &faceNeigbour, const fvMesh &procMesh, const labelUList &faceAddressing, const labelUList &cellAddressing, const labelUList &boundaryAddressing)
Construct from components without the complete mesh.
void decomposeFields(const UPtrList< GeoField > &fields) const
Decompose list of fields.
fvFieldDecomposer(const fvMesh &completeMesh, const fvMesh &procMesh, const labelUList &faceAddressing, const labelUList &cellAddressing, const labelUList &boundaryAddressing)
Construct from components using information from the complete mesh.
bool empty() const noexcept
True if no mappers have been allocated.
void reset(const fvMesh &completeMesh)
Reset mappers using information from the complete mesh.
~fvFieldDecomposer()=default
Destructor.
static int verbose_
Output verbosity when writing.
void clear()
Remove all mappers.
fvFieldDecomposer(Foam::zero, const fvMesh &procMesh, const labelUList &faceAddressing, const labelUList &cellAddressing, const labelUList &boundaryAddressing)
Construct without mappers, added later with reset().
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A FieldMapper for finite-volume patch fields.
fvPatchFieldMapper()=default
Default construct.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Mesh data needed to do the Finite Volume discretisation.
Definition surfaceMesh.H:47
A class for managing temporary objects.
Definition tmp.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
rDeltaTY field()
dynamicFvMesh & mesh
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Foam::surfaceFields.