faFieldDecomposer.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2021-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 Class
28  Foam::faFieldDecomposer
29 
30 Description
31  Finite Area area and edge field decomposer.
32 
33 Author
34  Zeljko Tukovic, FSB Zagreb
35  Hrvoje Jasak, Wikki Ltd.
36 
37 SourceFiles
38  faFieldDecomposer.C
39  fvFieldDecomposerCache.C
40  faFieldDecomposerTemplates.C
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef Foam_faFieldDecomposer_H
45 #define Foam_faFieldDecomposer_H
46 
47 #include "faMesh.H"
48 #include "faMeshSubset.H"
49 #include "faPatchFieldMapper.H"
50 #include "edgeFields.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 namespace Foam
55 {
56 
57 // Forward Declarations
58 class IOobjectList;
59 
60 /*---------------------------------------------------------------------------*\
61  Class faFieldDecomposer Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 {
66 public:
67 
68  // Public Classes
69 
70  //- Patch field decomposer class
72  :
73  public faPatchFieldMapper
74  {
75  // Private Data
76 
77  label sizeBeforeMapping_;
78  labelList directAddressing_;
79 
80  public:
81 
82  // Constructors
83 
84  //- Construct given addressing
86  (
87  const label sizeBeforeMapping,
88  const labelUList& addressingSlice,
89  const label addressingOffset
90  );
91 
92 
93  // Member functions
94 
95  label size() const
96  {
97  return directAddressing_.size();
98  }
99 
100  virtual label sizeBeforeMapping() const
101  {
102  return sizeBeforeMapping_;
103  }
105  bool direct() const
106  {
107  return true;
108  }
110  virtual bool hasUnmapped() const
111  {
112  return false;
113  }
115  const labelUList& directAddressing() const
116  {
117  return directAddressing_;
118  }
119  };
120 
121 
122  //- Processor patch field decomposer class
123  class processorAreaPatchFieldDecomposer
124  :
125  public faPatchFieldMapper
126  {
127  // Private Data
128 
129  label sizeBeforeMapping_;
130  labelListList addressing_;
131  scalarListList weights_;
132 
133  public:
134 
135  //- Construct addressing from details
137  (
138  const label nTotalFaces,
139  const labelUList& edgeOwner,
140  const labelUList& edgeNeigbour,
141  const labelUList& addressingSlice,
142  const scalarField& edgeWeights = scalarField::null()
143  );
144 
145  //- Construct given addressing from complete mesh
147  (
148  const faMesh& mesh,
149  const labelUList& addressingSlice
150  );
151 
152 
153  // Member Functions
154 
155  label size() const
156  {
157  return addressing_.size();
158  }
159 
160  virtual label sizeBeforeMapping() const
161  {
162  return sizeBeforeMapping_;
163  }
164 
165  bool direct() const
166  {
167  return false;
168  }
169 
170  virtual bool hasUnmapped() const
171  {
172  return false;
173  }
174 
175  const labelListList& addressing() const
176  {
177  return addressing_;
178  }
179 
180  const scalarListList& weights() const
181  {
182  return weights_;
183  }
184  };
186 
187  //- Processor patch field decomposer class
189  :
190  public faPatchFieldMapper
191  {
192  label sizeBeforeMapping_;
193  labelListList addressing_;
194  scalarListList weights_;
196  public:
197 
198  //- Construct given addressing
200  (
201  label sizeBeforeMapping,
202  const labelUList& addressingSlice
203  );
204 
205 
206  // Member Functions
207 
208  label size() const
209  {
210  return addressing_.size();
211  }
212 
213  virtual label sizeBeforeMapping() const
214  {
215  return sizeBeforeMapping_;
216  }
218  bool direct() const
219  {
220  return false;
221  }
223  virtual bool hasUnmapped() const
224  {
225  return false;
226  }
228  const labelListList& addressing() const
229  {
230  return addressing_;
231  }
233  const scalarListList& weights() const
234  {
235  return weights_;
236  }
237  };
238 
239 
240 private:
241 
242  // Private Data
243 
244  //- Reference to processor mesh
245  const faMesh& procMesh_;
246 
247  //- Reference to edge addressing
248  const labelList& edgeAddressing_;
249 
250  //- Reference to face addressing
251  const labelList& faceAddressing_;
252 
253  //- Reference to boundary addressing
254  const labelList& boundaryAddressing_;
255 
256  //- List of patch field decomposers
257  PtrList<patchFieldDecomposer> patchFieldDecomposerPtrs_;
258 
260  processorAreaPatchFieldDecomposerPtrs_;
261 
263  processorEdgePatchFieldDecomposerPtrs_;
264 
265 
266  // Private Member Functions
267 
268  //- No copy construct
269  faFieldDecomposer(const faFieldDecomposer&) = delete;
270 
271  //- No copy assignment
272  void operator=(const faFieldDecomposer&) = delete;
273 
274 
275 public:
276 
277  // Public Classes
278  class fieldsCache;
279 
280 
281  // Constructors
282 
283  //- Construct without mappers, added later with reset()
285  (
286  const Foam::zero,
287  const faMesh& procMesh, // Target mesh
288  const labelList& edgeAddressing,
289  const labelList& faceAddressing,
290  const labelList& boundaryAddressing
291  );
292 
293  //- Construct from components using information from the complete mesh
295  (
296  const faMesh& completeMesh, // Source mesh
297  const faMesh& procMesh, // Target mesh
298  const labelList& edgeAddressing,
299  const labelList& faceAddressing,
300  const labelList& boundaryAddressing
301  );
302 
303  //- Construct from components without the complete mesh
305  (
306  // Information about the complete mesh
307  const label nTotalFaces,
308  const List<labelRange>& boundaryRanges,
309  const labelUList& edgeOwner,
310  const labelUList& edgeNeigbour,
311 
312  // Addressing for processor mesh
313  const faMesh& procMesh, // Target mesh
314  const labelList& edgeAddressing,
315  const labelList& faceAddressing,
316  const labelList& boundaryAddressing
317  );
318 
319 
320  //- Destructor
321  ~faFieldDecomposer() = default;
322 
323 
324  // Member Functions
325 
326  //- True if no mappers have been allocated
327  bool empty() const;
328 
329  //- Remove all mappers
330  void clear();
331 
332  //- Reset mappers using information from the complete mesh
333  void reset(const faMesh& completeMesh);
334 
335  //- Reset mapper using information about the complete mesh
336  void reset
337  (
338  const label nTotalFaces,
339  const List<labelRange>& boundaryRanges,
340  const labelUList& edgeOwner,
341  const labelUList& edgeNeigbour
342  );
343 
344 
345  // Mapping
346 
347  //- Decompose area field
348  template<class Type>
351  (
353  ) const;
354 
355  //- Decompose surface field
356  template<class Type>
359  (
361  ) const;
362 
363  template<class GeoField>
364  void decomposeFields(const PtrList<GeoField>& fields) const;
365 
366 
367  // Reading helpers
368 
369  //- Read the fields and store on the pointer list
370  template
371  <
372  class Type,
373  template<class> class PatchField,
374  class GeoMesh
375  >
376  static void readFields
377  (
378  const typename GeoMesh::Mesh& mesh,
379  const IOobjectList& objects,
381  const bool readOldTime
382  );
383 
384  //- Read fields and store on the pointer list
385  template<class Mesh, class GeoField>
386  static void readFields
387  (
388  const Mesh& mesh,
389  const IOobjectList& objects,
391  );
392 };
393 
394 
395 
396 /*---------------------------------------------------------------------------*\
397  Class faFieldDecomposer::fieldsCache Declaration
398 \*---------------------------------------------------------------------------*/
399 
401 {
402  // Private Data
403 
404  class privateCache;
405 
406  //- All field and field-field types for finiteArea
407  std::unique_ptr<privateCache> cache_;
408 
409 
410  // Private Member Functions
411 
412  //- No copy construct
413  fieldsCache(const fieldsCache&) = delete;
414 
415  //- No copy assignment
416  void operator=(const fieldsCache&) = delete;
417 
418 
419 public:
420 
421  // Constructors
422 
423  //- Default construct
424  fieldsCache();
425 
426 
427  //- Destructor
428  ~fieldsCache();
429 
430 
431  // Member Functions
432 
433  //- No fields
434  bool empty() const;
435 
436  //- Number of fields
437  label size() const;
438 
439  //- Clear out
440  void clear();
441 
442 
443  //- Read all fields given mesh and objects
444  void readAllFields
445  (
446  const faMesh& mesh,
447  const IOobjectList& objects
448  );
449 
450  //- Read all fields given mesh and objects.
451  //- Supports reading/sending fields
452  void readAllFields
453  (
454  const bitSet& haveMeshOnProc,
455  const faMeshSubset* subsetter,
456  const faMesh& mesh,
457  IOobjectList& objects
458  );
459 
460  //- Read all fields given mesh and objects.
461  //- Supports reading/sending fields
462  void readAllFields
463  (
464  const boolList& haveMeshOnProc,
465  const faMeshSubset* subsetter,
466  const faMesh& mesh,
467  IOobjectList& objects
468  );
469 
470  //- Decompose and write all fields
471  void decomposeAllFields
472  (
473  const faFieldDecomposer& decomposer,
474  bool report = false
475  ) const;
476 };
477 
478 
479 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
480 
481 } // End namespace Foam
482 
483 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
484 
485 #ifdef NoRepository
487 #endif
488 
489 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
490 
491 #endif
492 
493 // ************************************************************************* //
bool direct() const
Is it a direct (non-interpolating) mapper?
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
const labelListList & addressing() const
Return the interpolation addressing.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual bool hasUnmapped() const
Any unmapped values?
rDeltaTY field()
Finite Area area and edge field decomposer.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
~faFieldDecomposer()=default
Destructor.
void readAllFields(const faMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: faMeshSubset.H:61
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
const labelListList & addressing() const
Return the interpolation addressing.
processorEdgePatchFieldDecomposer(label sizeBeforeMapping, const labelUList &addressingSlice)
Construct given addressing.
void decomposeFields(const PtrList< GeoField > &fields) const
const scalarListList & weights() const
Return the interpolation weights.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
virtual bool hasUnmapped() const
Any unmapped values?
bool direct() const
Is it a direct (non-interpolating) mapper?
virtual bool hasUnmapped() const
Any unmapped values?
static void readFields(const typename GeoMesh::Mesh &mesh, const IOobjectList &objects, PtrList< GeometricField< Type, PatchField, GeoMesh >> &fields, const bool readOldTime)
Read the fields and store on the pointer list.
void reset(const faMesh &completeMesh)
Reset mappers using information from the complete mesh.
dynamicFvMesh & mesh
patchFieldDecomposer(const label sizeBeforeMapping, const labelUList &addressingSlice, const label addressingOffset)
Construct given addressing.
const scalarListList & weights() const
Return the interpolation weights.
MESH Mesh
Definition: GeoMesh.H:58
void clear()
Remove all mappers.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static const Field< scalar > & null()
Return nullObject reference Field.
Definition: FieldI.H:24
bool empty() const
True if no mappers have been allocated.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
label size() const
The size of the mapper.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:42
A class for managing temporary objects.
Definition: HashPtrTable.H:50
A FieldMapper for finite-area patch fields.
label size() const
Number of fields.
processorAreaPatchFieldDecomposer(const label nTotalFaces, const labelUList &edgeOwner, const labelUList &edgeNeigbour, const labelUList &addressingSlice, const scalarField &edgeWeights=scalarField::null())
Construct addressing from details.
tmp< GeometricField< Type, faPatchField, areaMesh > > decomposeField(const GeometricField< Type, faPatchField, areaMesh > &field) const
Decompose area field.
bool direct() const
Is it a direct (non-interpolating) mapper?
const labelUList & directAddressing() const
Return the direct addressing values.
Namespace for OpenFOAM.
void decomposeAllFields(const faFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.