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-------------------------------------------------------------------------------
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 "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
37namespace Foam
38{
39
40// All area/edge field types
42{
43public:
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(const faMesh& mesh, const IOobjectList& objects)
81 {
82 #undef doLocalCode
83 #define doLocalCode(Type) \
84 { \
85 fieldsDistributor::readFields \
86 ( \
87 mesh, \
88 objects, \
89 Type##AreaFields_ \
90 ); \
91 fieldsDistributor::readFields \
92 ( \
93 mesh, \
94 objects, \
95 Type##AreaFields_, \
96 false /* readOldTime = false */ \
97 ); \
98 fieldsDistributor::readFields \
99 ( \
100 mesh, \
101 objects, \
102 Type##EdgeFields_, \
103 false /* readOldTime = false */ \
104 ); \
105 }
106
107 doLocalCode(scalar);
112
113 #undef doLocalCode
114 }
115
116 template<class GeoField>
117 static void decompose
118 (
119 const faFieldDecomposer& decomposer,
121 bool report
122 )
123 {
124 if (!fields.empty())
125 {
126 if (report)
127 {
128 Info<< " "
130 << "s: "
132 }
133
134 decomposer.decomposeFields(fields);
135 }
136 }
137
138 void decomposeAll(const faFieldDecomposer& decomposer, bool report) const
139 {
140 #undef doLocalCode
141 #define doLocalCode(Flavour) \
142 { \
143 decompose(decomposer, scalar##Flavour##Fields_, report); \
144 decompose(decomposer, vector##Flavour##Fields_, report); \
145 decompose(decomposer, sphericalTensor##Flavour##Fields_, report); \
146 decompose(decomposer, symmTensor##Flavour##Fields_, report); \
147 decompose(decomposer, tensor##Flavour##Fields_, report); \
148 }
149
150 doLocalCode(Area);
151 doLocalCode(Edge);
152
153 #undef doLocalCode
154 }
155};
156
157} // End namespace Foam
158
159
160// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161
163:
164 cache_(new privateCache)
165{}
166
167
168// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
169
170// Destructor not in header (used incomplete type)
172{}
173
174
175// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176
178{
179 return (!cache_ || cache_->empty());
180}
181
182
184{
185 return (cache_ ? cache_->size() : label(0));
186}
187
188
190{
191 cache_.reset(new privateCache);
192}
193
194
196(
197 const faMesh& mesh,
198 const IOobjectList& objects
199)
200{
201 if (cache_)
202 {
203 cache_->readAll(mesh, objects);
204 }
205}
206
207
209(
210 const faFieldDecomposer& decomposer,
211 bool report
212) const
213{
214 if (cache_)
215 {
216 cache_->decomposeAll(decomposer, report);
217 }
218}
219
220
221// ************************************************************************* //
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 faMesh &mesh, const IOobjectList &objects)
void decomposeAll(const faFieldDecomposer &decomposer, bool report) const
static void decompose(const faFieldDecomposer &decomposer, const PtrList< GeoField > &fields, bool report)
label size() const
Number of fields.
void readAllFields(const faMesh &mesh, const IOobjectList &objects)
Read all fields given mesh and objects.
void decomposeAllFields(const faFieldDecomposer &decomposer, bool report=false) const
Decompose and write all fields.
Finite Area area and edge field decomposer.
void decomposeFields(const PtrList< GeoField > &fields) const
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:100
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
#define doLocalCode(GeoField)