faFieldDecomposerCache.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) 2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "faFieldDecomposer.H"
29 #include "fieldsDistributor.H"
30 #include "areaFields.H"
31 #include "edgeFields.H"
32 #include "IOobjectList.H"
33 #include "PtrListOps.H"
34 
35 // * * * * * * * * * * * * * * * * Declarations * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // All area/edge field types
41 class faFieldDecomposer::fieldsCache::privateCache
42 {
43 public:
44 
45  #undef declareField
46  #define declareField(Type) \
47  PtrList<GeometricField<Type, faPatchField, areaMesh>> Type##AreaFields_; \
48  PtrList<GeometricField<Type, faePatchField, edgeMesh>> Type##EdgeFields_;
49 
50  declareField(scalar);
55  #undef declareField
56 
57  label size() const noexcept
58  {
59  label count = 0;
60 
61  #undef doLocalCode
62  #define doLocalCode(Type) \
63  { \
64  count += Type##AreaFields_.size(); \
65  count += Type##EdgeFields_.size(); \
66  }
67 
68  doLocalCode(scalar);
73  #undef doLocalCode
74 
75  return count;
76  }
77 
78  bool empty() const noexcept { return !size(); }
79 
80  void readAll
81  (
82  const faMesh& mesh,
83  const IOobjectList& objects
84  )
85  {
86  #undef doLocalCode
87  #define doLocalCode(Type) \
88  { \
89  fieldsDistributor::readFields \
90  ( \
91  mesh, \
92  objects, \
93  Type##AreaFields_ \
94  ); \
95  fieldsDistributor::readFields \
96  ( \
97  mesh, \
98  objects, \
99  Type##AreaFields_, \
100  false /* readOldTime = false */ \
101  ); \
102  fieldsDistributor::readFields \
103  ( \
104  mesh, \
105  objects, \
106  Type##EdgeFields_, \
107  false /* readOldTime = false */ \
108  ); \
109  }
110 
111  doLocalCode(scalar);
116 
117  #undef doLocalCode
118  }
119 
120  template<class BoolListType>
121  void readAll
122  (
123  const BoolListType& haveMeshOnProc,
124  const faMeshSubset*& subsetter,
125  const faMesh& mesh,
126  IOobjectList& objects
127  )
128  {
129  #undef doLocalCode
130  #define doLocalCode(Type) \
131  { \
132  fieldsDistributor::readFields \
133  ( \
134  haveMeshOnProc, \
135  subsetter, \
136  mesh, \
137  objects, \
138  Type##AreaFields_ \
139  ); \
140  fieldsDistributor::readFields \
141  ( \
142  haveMeshOnProc, \
143  subsetter, \
144  mesh, \
145  objects, \
146  Type##EdgeFields_ \
147  ); \
148  }
149 
150  doLocalCode(scalar);
155 
156  #undef doLocalCode
157  }
158 
159  template<class GeoField>
160  static void decompose
161  (
162  const faFieldDecomposer& decomposer,
163  const PtrList<GeoField>& fields,
164  bool report
165  )
166  {
167  if (!fields.empty())
168  {
169  if (report)
170  {
171  Info<< " "
173  << "s: "
175  }
176 
177  decomposer.decomposeFields(fields);
178  }
179  }
180 
181  void decomposeAll(const faFieldDecomposer& decomposer, bool report) const
182  {
183  #undef doLocalCode
184  #define doLocalCode(Flavour) \
185  { \
186  decompose(decomposer, scalar##Flavour##Fields_, report); \
187  decompose(decomposer, vector##Flavour##Fields_, report); \
188  decompose(decomposer, sphericalTensor##Flavour##Fields_, report); \
189  decompose(decomposer, symmTensor##Flavour##Fields_, report); \
190  decompose(decomposer, tensor##Flavour##Fields_, report); \
191  }
192 
193  doLocalCode(Area);
194  doLocalCode(Edge);
195 
196  #undef doLocalCode
197  }
198 };
199 
200 } // End namespace Foam
201 
202 
203 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
204 
206 :
207  cache_(new privateCache)
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
212 
213 // Destructor not in header (used incomplete type)
215 {}
216 
217 
218 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 {
222  return (!cache_ || cache_->empty());
223 }
224 
227 {
228  return (cache_ ? cache_->size() : label(0));
229 }
230 
231 
233 {
234  cache_.reset(new privateCache);
235 }
236 
237 
239 (
240  const faMesh& mesh,
241  const IOobjectList& objects
242 )
243 {
244  if (cache_)
245  {
246  cache_->readAll(mesh, objects);
247  }
248 }
249 
250 
252 (
253  const bitSet& haveMeshOnProc,
254  const faMeshSubset* subsetter,
255  const faMesh& mesh,
256  IOobjectList& objects
257 )
258 {
259  if (cache_)
260  {
261  cache_->readAll(haveMeshOnProc, subsetter, mesh, objects);
262  }
263 }
264 
265 
267 (
268  const boolList& haveMeshOnProc,
269  const faMeshSubset* subsetter,
270  const faMesh& mesh,
271  IOobjectList& objects
272 )
273 {
274  if (cache_)
275  {
276  cache_->readAll(haveMeshOnProc, subsetter, mesh, objects);
277  }
278 }
279 
280 
282 (
283  const faFieldDecomposer& decomposer,
284  bool report
285 ) const
286 {
287  if (cache_)
288  {
289  cache_->decomposeAll(decomposer, report);
290  }
291 }
292 
293 
294 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:87
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
#define doLocalCode(Type)
Finite Area area and edge field decomposer.
void readAll(const faMesh &mesh, const IOobjectList &objects)
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
#define declareField(Type)
void readAllFields(const faMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: faMeshSubset.H:61
Functions to operate on Pointer Lists.
void decomposeFields(const PtrList< GeoField > &fields) const
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
dynamicFvMesh & mesh
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static void decompose(const faFieldDecomposer &decomposer, const PtrList< GeoField > &fields, bool report)
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
messageStream Info
Information stream (stdout output on master, null elsewhere)
void decomposeAll(const faFieldDecomposer &decomposer, bool report) const
Tensor of scalars, i.e. Tensor<scalar>.
label size() const
Number of fields.
Namespace for OpenFOAM.
void decomposeAllFields(const faFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225