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 -------------------------------------------------------------------------------
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::fvFieldDecomposer
29 
30 Description
31  Finite Volume volume and surface field decomposer.
32 
33 SourceFiles
34  fvFieldDecomposer.C
35  fvFieldDecomposerCache.C
36  fvFieldDecomposerTemplates.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef Foam_fvFieldDecomposer_H
41 #define Foam_fvFieldDecomposer_H
42 
43 #include "fvMesh.H"
44 #include "fvPatchFieldMapper.H"
45 #include "surfaceFields.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 class IOobjectList;
54 
55 /*---------------------------------------------------------------------------*\
56  Class fvFieldDecomposer Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 {
61 public:
62 
63  // Public Classes
64 
65  //- Patch field decomposer class
67  :
68  public fvPatchFieldMapper
69  {
70  // Private Data
71 
72  labelList directAddressing_;
73 
74  public:
75 
76  // Constructors
77 
78  //- Construct given addressing
80  (
81  const labelUList& addressingSlice,
82  const label addressingOffset
83  );
84 
85 
86  // Member Functions
87 
88  label size() const
89  {
90  return directAddressing_.size();
91  }
92 
93  bool direct() const
94  {
95  return true;
96  }
97 
98  //- Are there unmapped values
99  bool hasUnmapped() const
100  {
101  return false;
102  }
103 
104  const labelUList& directAddressing() const
105  {
106  return directAddressing_;
107  }
108  };
109 
110 
111  //- Processor patch field decomposer class. Maps either owner or
112  // neighbour data (no interpolate anymore - processorFvPatchField
113  // holds neighbour data)
114  class processorVolPatchFieldDecomposer
115  :
116  public fvPatchFieldMapper
117  {
118  // Private Data
119 
120  labelList directAddressing_;
121 
122  public:
123 
124  //- Construct addressing from details
126  (
127  const labelUList& faceOwner,
128  const labelUList& faceNeigbour,
129  const labelUList& addressingSlice
130  );
131 
132  //- Construct given addressing from complete mesh
134  (
135  const polyMesh& mesh,
136  const labelUList& addressingSlice
137  );
138 
139 
140  // Member Functions
141 
142  label size() const
143  {
144  return directAddressing_.size();
145  }
146 
147  bool direct() const
148  {
149  return true;
150  }
151 
152  //- Are there unmapped values
153  bool hasUnmapped() const
154  {
155  return false;
156  }
157 
158  const labelUList& directAddressing() const
159  {
160  return directAddressing_;
161  }
162  };
164 
165  //- Processor patch field decomposer class. Surface field is assumed
166  // to have direction (so manipulates sign when mapping)
168  :
169  public fvPatchFieldMapper
170  {
171  labelListList addressing_;
172  scalarListList weights_;
173 
174  public:
175 
176  //- Construct given addressing
178  (
179  const labelUList& addressingSlice
180  );
181 
182 
183  // Member functions
184 
185  label size() const
186  {
187  return addressing_.size();
188  }
189 
190  bool direct() const
191  {
192  return false;
193  }
194 
195  //- Are there unmapped values
196  bool hasUnmapped() const
197  {
198  return false;
199  }
201  const labelListList& addressing() const
202  {
203  return addressing_;
204  }
206  const scalarListList& weights() const
207  {
208  return weights_;
209  }
210  };
211 
212 
213 private:
214 
215  // Private Data
216 
217  //- Reference to processor mesh
218  const fvMesh& procMesh_;
219 
220  //- Reference to face addressing
221  const labelList& faceAddressing_;
222 
223  //- Reference to cell addressing
224  const labelList& cellAddressing_;
225 
226  //- Reference to boundary addressing
227  const labelList& boundaryAddressing_;
228 
229  //- List of patch field decomposers
230  PtrList<patchFieldDecomposer> patchFieldDecomposerPtrs_;
231 
233  processorVolPatchFieldDecomposerPtrs_;
234 
236  processorSurfacePatchFieldDecomposerPtrs_;
237 
238  PtrList<scalarField> faceSign_;
239 
240 
241  // Private Member Functions
242 
243  //- No copy construct
244  fvFieldDecomposer(const fvFieldDecomposer&) = delete;
245 
246  //- No copy assignment
247  void operator=(const fvFieldDecomposer&) = delete;
248 
249 public:
250 
251  // Public Classes
252  class fieldsCache;
253 
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()
265  (
266  const Foam::zero,
267  const fvMesh& procMesh,
268  const labelList& faceAddressing,
269  const labelList& cellAddressing,
270  const labelList& boundaryAddressing
271  );
272 
273  //- Construct from components using information from the complete mesh
275  (
276  const fvMesh& completeMesh,
277  const fvMesh& procMesh,
278  const labelList& faceAddressing,
279  const labelList& cellAddressing,
280  const labelList& boundaryAddressing
281  );
282 
283  //- Construct from components without the complete mesh
285  (
286  // Information about the complete mesh
287  const List<labelRange>& boundaryRanges,
288  const labelUList& faceOwner,
289  const labelUList& faceNeigbour,
290 
291  // Addressing for processor mesh
292  const fvMesh& procMesh,
293  const labelList& faceAddressing,
294  const labelList& cellAddressing,
295  const labelList& boundaryAddressing
296  );
297 
298 
299  //- Destructor
300  ~fvFieldDecomposer() = default;
301 
302 
303  // Member Functions
304 
305  //- True if no mappers have been allocated
306  bool empty() const;
307 
308  //- Remove all mappers
309  void clear();
310 
311  //- Reset mappers using information from the complete mesh
312  void reset(const fvMesh& completeMesh);
313 
314  //- Reset mapper using information about the complete mesh
315  void reset
316  (
317  const List<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  (
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 PtrList<GeoField>& fields) const;
353 };
354 
355 
356 /*---------------------------------------------------------------------------*\
357  Class fvFieldDecomposer::fieldsCache Declaration
358 \*---------------------------------------------------------------------------*/
359 
361 {
362  // Private Data
363 
364  class privateCache;
365 
366  //- All field and field-field types for lagrangian
367  std::unique_ptr<privateCache> cache_;
368 
369 
370  // Private Member Functions
371 
372  //- No copy construct
373  fieldsCache(const fieldsCache&) = delete;
374 
375  //- No copy assignment
376  void operator=(const fieldsCache&) = delete;
377 
378 
379 public:
380 
381  // Constructors
382 
383  //- Default construct
384  fieldsCache();
385 
386 
387  //- Destructor
388  ~fieldsCache();
389 
390 
391  // Member Functions
392 
393  //- No fields
394  bool empty() const;
395 
396  //- Total number of fields
397  label size() const;
398 
399  //- Clear out
400  void clear();
401 
402 
403  //- Read all fields given mesh and objects
404  void readAllFields
405  (
406  const fvMesh& mesh,
407  const IOobjectList& objects
408  );
409 
410  //- Decompose and write all fields
411  void decomposeAllFields
412  (
413  const fvFieldDecomposer& decomposer,
414  bool report = false
415  ) const;
416 };
418 
419 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
420 
421 } // End namespace Foam
422 
423 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
424 
425 #ifdef NoRepository
427 #endif
428 
429 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
430 
431 #endif
432 
433 // ************************************************************************* //
Foam::surfaceFields.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
rDeltaTY field()
static int verbose_
Output verbosity when writing.
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
bool hasUnmapped() const
Are there unmapped values.
processorSurfacePatchFieldDecomposer(const labelUList &addressingSlice)
Construct given addressing.
label size() const
Total number of fields.
const scalarListList & weights() const
Return the interpolation weights.
void clear()
Remove all mappers.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
void decomposeAllFields(const fvFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
bool direct() const
Is it a direct (non-interpolating) mapper?
const labelUList & directAddressing() const
Return the direct addressing values.
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
bool empty() const
True if no mappers have been allocated.
void reset(const fvMesh &completeMesh)
Reset mappers using information from the complete mesh.
tmp< DimensionedField< Type, volMesh > > decomposeField(const DimensionedField< Type, volMesh > &field) const
Decompose internal field.
void readAllFields(const fvMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
dynamicFvMesh & mesh
~fvFieldDecomposer()=default
Destructor.
bool direct() const
Is it a direct (non-interpolating) mapper?
A FieldMapper for finite-volume patch fields.
Processor patch field decomposer class. Surface field is assumed.
Finite Volume volume and surface field decomposer.
patchFieldDecomposer(const labelUList &addressingSlice, const label addressingOffset)
Construct given addressing.
bool direct() const
Is it a direct (non-interpolating) mapper?
label size() const
The size of the mapper.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
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...
Definition: areaFieldsFwd.H:42
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
const labelUList & directAddressing() const
Return the direct addressing values.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
A class for managing temporary objects.
Definition: HashPtrTable.H:50
processorVolPatchFieldDecomposer(const labelUList &faceOwner, const labelUList &faceNeigbour, const labelUList &addressingSlice)
Construct addressing from details.
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose list of fields.
const labelListList & addressing() const
Return the interpolation addressing.
bool hasUnmapped() const
Are there unmapped values.
Namespace for OpenFOAM.