Loading...
Searching...
No Matches
zoneBlended.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) 2024 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
26Class
27 Foam::zoneBlended
28
29Group
30 grpFvSurfaceInterpolationSchemes
31
32Description
33 Multi-faceZone based blending differencing scheme.
34
35 Schemes are set in dictonary format according to:
36
37 \verbatim
38 divSchemes
39 {
40 .
41 .
42 div(phi,U) Gauss zoneBlended
43 {
44 default defaultScheme;
45 faceZone1 scheme1;
46 faceZone2 scheme2;
47 ...
48 faceZoneN schemeN;
49 }
50 .
51 .
52 }
53 \endverbatim
54
55 The \c default entry specifies the background scheme; additional schemes
56 can be set per \c faceZone, e.g. \c scheme1 is applied to \c facZone1,
57 \c scheme2 is applied to \c facZone2 etc.
58
59
60Usage
61 Example of the \c zoneBlended scheme to use \c linearUpwind as the
62 background scheme and \c upwind in \c faceZone1:
63
64 \verbatim
65 divSchemes
66 {
67 .
68 .
69 div(phi,U) Gauss zoneBlended
70 {
71 default linearUpwind grad(U);
72 faceZone1 upwind;
73 };
74 .
75 .
76 }
77 \endverbatim
78
79SourceFiles
80 zoneBlended.C
81
82\*---------------------------------------------------------------------------*/
83
84#ifndef Foam_zoneBlended_H
85#define Foam_zoneBlended_H
86
88#include "blendedSchemeBase.H"
89#include "surfaceInterpolate.H"
90#include "UIndirectList.H"
91
92// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93
94namespace Foam
95{
96
97/*---------------------------------------------------------------------------*\
98 Class zoneBlended Declaration
99\*---------------------------------------------------------------------------*/
100
101template<class Type>
102class zoneBlended
103:
104 public surfaceInterpolationScheme<Type>,
105 public blendedSchemeBase<Type>
106{
109
110 // Private data
111
112 //- Face zones
113 wordList zoneNames_;
114
115 //- Schemes
116 // Note: 0 index used to describe default/background scheme
118
119 //- Corrected flag - true if any of the schemes has corrected() set
120 bool corrected_;
121
122
123 // Private Member Functions
124
125 //- Set the lists of face zone names and schemes
126 void setSchemes(const fvMesh& mesh, const dictionary& dict);
127
128 //- Set the lists of face zone names and schemes
129 void setSchemes
130 (
131 const fvMesh& mesh,
132 const surfaceScalarField& faceFlux,
133 const dictionary& dict
134 );
135
136 //- Retrieve a scheme from the list
137 const surfaceInterpolationScheme<Type>& scheme
138 (
139 const label schemei
140 ) const
141 {
142 return schemePtrs_[schemei]();
143 }
144
145 //- Set destination values from source values for a face zone
146 template<class FieldType>
147 void setFaceZoneValues
148 (
149 FieldType& dest,
150 const FieldType& src,
151 const faceZone& fz
152 ) const;
153
154 //- Set destination values to zero for a face zone
155 template<class FieldType>
156 void zeroFaceZoneValues(FieldType& dest, const faceZone& fz) const;
157
158 //- No copy construct
159 zoneBlended(const zoneBlended&) = delete;
160
161 //- No copy assignment
162 void operator=(const zoneBlended&) = delete;
163
164
165public:
166
167 //- Runtime type information
168 TypeName("zoneBlended");
169
170
171 // Constructors
172
173 //- Construct from mesh and Istream.
174 // The name of the flux field is read from the Istream and looked-up
175 // from the mesh objectRegistry
176 zoneBlended(const fvMesh& mesh, Istream& is)
177 :
178 surfaceInterpolationScheme<Type>(mesh),
179 zoneNames_(),
180 schemePtrs_(),
181 corrected_(false)
182 {
183 const dictionary dict(is);
184
185 setSchemes(mesh, dict);
187
188
189 //- Construct from mesh, faceFlux and Istream
191 (
192 const fvMesh& mesh,
193 const surfaceScalarField& faceFlux,
194 Istream& is
195 )
196 :
197 surfaceInterpolationScheme<Type>(mesh),
198 zoneNames_(),
199 schemePtrs_(),
200 corrected_(false)
201 {
202 const dictionary dict(is);
203
204 setSchemes(mesh, faceFlux, dict);
205 }
206
207
208 // Member Functions
209
210 //- Return the face-based blending factor
212 (
213 const VolumeField& vf
214 ) const
215 {
216 const auto& mesh = vf.mesh();
217 auto tbf = surfaceScalarField::New("blendingFactor", mesh, dimless);
218 auto& bf = tbf.ref();
219 auto& bbf = bf.boundaryFieldRef();
220 bf = 0.0;
221
222 const auto& pbm = mesh.boundaryMesh();
223 const auto& zones = mesh.faceZones();
224
225 // Use blending factor to show different zones
226 for (label zonei=1; zonei<zoneNames_.size(); ++zonei)
227 {
228 const word& name = zoneNames_[zonei];
229 const auto& fz = zones[name];
230
231 for (const label facei : fz)
232 {
233 if (mesh.isInternalFace(facei))
234 {
235 bf[facei] = zonei;
237 else
238 {
239 const labelPair pf = pbm.whichPatchFace(facei);
240 auto& pbf = bbf[pf.first()];
241 if (pbf.size())
242 {
243 pbf[pf.second()] = zonei;
244 }
245 }
246 }
247 }
248
249 return tbf;
250 }
251
252
253 //- Return the interpolation weighting factors
254 tmp<surfaceScalarField> weights(const VolumeField& vf) const
255 {
256 const auto& mesh = vf.mesh();
257 auto tweights =
258 surfaceScalarField::New("weights", vf.mesh(), dimless);
259 auto& weights = tweights.ref();
260
261 // Set default scheme weights
262 weights = this->scheme(0).weights(vf);
263
264 // Set face zone weights
265 const auto& zones = mesh.faceZones();
266
267 for (label schemei=1; schemei<schemePtrs_.size(); ++schemei)
268 {
269 const auto& scheme = this->scheme(schemei);
270
271 auto tschemeWeights = scheme.weights(vf);
272 const auto& schemeWeights = tschemeWeights();
273 const auto& fz = zones[zoneNames_[schemei]];
274
275 setFaceZoneValues(weights, schemeWeights, fz);
276 }
277
278 return tweights;
279 }
280
282 //- Return the face-interpolate of the given cell field
283 // with explicit correction
284 tmp<SurfaceField> interpolate(const VolumeField& vf) const
285 {
286 return
288 (
289 vf,
290 weights(vf)
291 );
292 }
293
294
295 //- Return true if this scheme uses an explicit correction
296 virtual bool corrected() const
297 {
298 return corrected_;
299 }
300
301
302 //- Return the explicit correction to the face-interpolate
303 //- for the given field
304 virtual tmp<SurfaceField> correction
305 (
306 const VolumeField& vf
307 ) const
308 {
309 const auto& mesh = vf.mesh();
310 auto tcorr =
311 SurfaceField::New("correction", vf.mesh(), vf.dimensions());
312 auto& corr = tcorr.ref();
313 corr = Zero;
315 // Set default scheme correction
316 const auto& scheme0 = this->scheme(0);
317 if (scheme0.corrected())
318 {
319 corr = scheme0.correction(vf);
320 }
321
322 // Only the default scheme exists - can exit early
323 if (schemePtrs_.size() == 1) return tcorr;
324
325 // Set correction field in faceZones
326 const auto& zones = mesh.faceZones();
327
328 for (label schemei=1; schemei<schemePtrs_.size(); ++schemei)
329 {
330 const auto& scheme = this->scheme(schemei);
331
332 if (scheme.corrected())
333 {
334 auto tschemeCorr = scheme.correction(vf);
335 const auto& schemeCorr = tschemeCorr();
336 const auto& fz = zones[zoneNames_[schemei]];
337
338 setFaceZoneValues(corr, schemeCorr, fz);
339 }
340 else
341 {
342 if (scheme0.corrected())
343 {
344 // Remove correction from base scheme face zone faces
345 const auto& fz = zones[zoneNames_[schemei]];
346 zeroFaceZoneValues(corr, fz);
347 }
348 }
349 }
350
351 return tcorr;
352 }
353};
354
355
356// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
357
358} // End namespace Foam
359
360// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
361
362#ifdef NoRepository
363 #include "zoneBlendedTemplates.C"
364#endif
365
366// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
367
368#endif
369
370// ************************************************************************* //
const polyBoundaryMesh & pbm
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvsPatchField< scalar >::calculatedType())
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
const T & first() const noexcept
Access the first element.
Definition Pair.H:137
const T & second() const noexcept
Access the second element.
Definition Pair.H:147
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
blendedSchemeBase()=default
Constructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A subset of mesh faces organised as a primitive patch.
Definition faceZone.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
Multi-faceZone based blending differencing scheme.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
tmp< surfaceScalarField > weights(const VolumeField &vf) const
Return the interpolation weighting factors.
zoneBlended(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
TypeName("zoneBlended")
Runtime type information.
tmp< SurfaceField > interpolate(const VolumeField &vf) const
Return the face-interpolate of the given cell field.
zoneBlended(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
virtual tmp< SurfaceField > correction(const VolumeField &vf) const
Return the explicit correction to the face-interpolate for the given field.
virtual tmp< surfaceScalarField > blendingFactor(const VolumeField &vf) const
Return the face-based blending factor.
dynamicFvMesh & mesh
Namespace for OpenFOAM.
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
const dimensionSet dimless
Dimensionless.
GeometricField< Type, fvPatchField, volMesh > VolumeField
A volume field for a given type.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68