fvFieldDecomposerCache.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-------------------------------------------------------------------------------
10License
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 "fvFieldDecomposer.H"
29#include "fieldsDistributor.H"
30#include "volFields.H"
31#include "surfaceFields.H"
32#include "IOobjectList.H"
33#include "PtrListOps.H"
34
35// * * * * * * * * * * * * * * * * Declarations * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39
40// All volume internal types
42{
43public:
44
45 #undef declareField
46 #define declareField(Type) \
47 PtrList<DimensionedField<Type, volMesh>> Type##DimFields_; \
48 PtrList<GeometricField<Type, fvPatchField, volMesh>> Type##VolFields_; \
49 PtrList<GeometricField<Type, fvsPatchField, surfaceMesh>> Type##SurfFields_;
50
51 declareField(scalar);
56 #undef declareField
57
58 label size() const noexcept
59 {
60 label count = 0;
61
62 #undef doLocalCode
63 #define doLocalCode(Type) \
64 { \
65 count += Type##DimFields_.size(); \
66 count += Type##VolFields_.size(); \
67 count += Type##SurfFields_.size(); \
68 }
69
70 doLocalCode(scalar);
75 #undef doLocalCode
76
77 return count;
78 }
79
80 bool empty() const noexcept { return !size(); }
81
82 void readAll(const fvMesh& mesh, const IOobjectList& objects)
83 {
84 #undef doLocalCode
85 #define doLocalCode(Type) \
86 { \
87 fieldsDistributor::readFields \
88 ( \
89 mesh, \
90 objects, \
91 Type##DimFields_ \
92 ); \
93 fieldsDistributor::readFields \
94 ( \
95 mesh, \
96 objects, \
97 Type##VolFields_, \
98 false /* readOldTime = false */ \
99 ); \
100 fieldsDistributor::readFields \
101 ( \
102 mesh, \
103 objects, \
104 Type##SurfFields_, \
105 false /* readOldTime = false */ \
106 ); \
107 }
108
109 doLocalCode(scalar);
114
115 #undef doLocalCode
116 }
117
118 template<class GeoField>
119 static void decompose
120 (
121 const fvFieldDecomposer& decomposer,
123 bool report
124 )
125 {
126 if (!fields.empty())
127 {
128 if (report)
129 {
130 Info<< " "
132 << "s: "
134 }
135
136 decomposer.decomposeFields(fields);
137 }
138 }
139
141 (
142 const fvFieldDecomposer& decomposer,
143 bool report
144 ) const
145 {
146 #undef doLocalCode
147 #define doLocalCode(Flavour) \
148 { \
149 decompose(decomposer, scalar##Flavour##Fields_, report); \
150 decompose(decomposer, vector##Flavour##Fields_, report); \
151 decompose(decomposer, sphericalTensor##Flavour##Fields_, report); \
152 decompose(decomposer, symmTensor##Flavour##Fields_, report); \
153 decompose(decomposer, tensor##Flavour##Fields_, report); \
154 }
155
156 doLocalCode(Vol);
157 doLocalCode(Surf);
158 doLocalCode(Dim);
159
160 #undef doLocalCode
161 }
162};
163
164} // End namespace Foam
165
166
167// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168
170:
171 cache_(new privateCache)
172{}
173
174
175// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
176
177// Destructor not in header (used incomplete type)
179{}
180
181
182// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183
185{
186 return (!cache_ || cache_->empty());
187}
188
189
191{
192 return (cache_ ? cache_->size() : label(0));
193}
194
195
197{
198 cache_.reset(new privateCache);
199}
200
201
203(
204 const fvMesh& mesh,
205 const IOobjectList& objects
206)
207{
208 if (cache_)
209 {
210 cache_->readAll(mesh, objects);
211 }
212}
213
214
216(
217 const fvFieldDecomposer& decomposer,
218 bool report
219) const
220{
221 if (cache_)
222 {
223 cache_->decomposeAll(decomposer, report);
224 }
225}
226
227
228// ************************************************************************* //
Functions to operate on Pointer Lists.
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
void readAll(const fvMesh &mesh, const IOobjectList &objects)
void decomposeAll(const fvFieldDecomposer &decomposer, bool report) const
static void decompose(const fvFieldDecomposer &decomposer, const PtrList< GeoField > &fields, bool report)
label size() const
Total number of 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.
Finite Volume volume and surface field decomposer.
void decomposeFields(const PtrList< GeoField > &fields) const
Decompose list of fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A traits class, which is primarily used for primitives.
Definition: pTraits.H:59
dynamicFvMesh & mesh
#define declareField(Type)
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Namespace for OpenFOAM.
messageStream Info
Information stream (stdout output on master, null elsewhere)
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
const direction noexcept
Definition: Scalar.H:223
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::surfaceFields.
#define doLocalCode(GeoField)