Loading...
Searching...
No Matches
stabilityBlendingFactor.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) 2018-2020 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::functionObjects::stabilityBlendingFactor
28
29Group
30 grpFieldFunctionObjects
31
32Description
33 Computes the \c stabilityBlendingFactor to be used by the
34 local blended convection scheme. The output is a surface field weight
35 between 0-1.
36
37 The weight of a blended scheme, i.e. \c w, is given by a function of
38 the blending factor, \c f:
39
40 \f[
41 w = f_{scheme_1} + (1 - f_{scheme_2})
42 \f]
43
44 The factor is calculated based on six criteria:
45 \verbatim
46 1. mesh non-orthogonality field
47 2. magnitude of cell centres gradient
48 3. convergence rate of residuals
49 4. faceWeight
50 5. skewness
51 6. Courant number
52 \endverbatim
53
54 The user can enable them individually.
55
56 For option 1, the following relation is used, where \f$\phi_1\f$ is
57 the non-orthogonality:
58 \f[
59 fNon =
60 min
61 (
62 max
63 (
64 0.0,
65 (\phi_1 - max(\phi_1))
66 /(min(\phi_1) - max(\phi_1))
67 ),
68 1.0
69 )
70 \f]
71
72 For option 2, the following relation is used, where \f$\phi_2\f$ is
73 the magnitude of cell centres gradient (Note that \f$\phi_2 = 3\f$
74 for orthogonal meshes):
75
76 \f[
77 fMagGradCc =
78 min
79 (
80 max
81 (
82 0.0,
83 (\phi_2 - max(\phi_2))
84 / (min(\phi_2) - max(\phi_2))
85 ),
86 1.0
87 )
88 \f]
89
90 For option 3, a PID control is used in order to control residual
91 unbounded fluctuations for individual cells.
92
93 \f[
94 factor =
95 P*residual
96 + I*residualIntegral
97 + D*residualDifferential
98 \f]
99
100 where \c P, \c I and \c D are user inputs.
101
102 The following relation is used:
103 \f[
104 fRes = (factor - meanRes)/(maxRes*meanRes);
105 \f]
106
107 where
108 \vartable
109 meanRes | Average(residual)
110 maxRes | User input
111 \endvartable
112
113 Note that \f$f_{Res}\f$ will blend more towards one as
114 the cell residual is larger then the domain mean residuals.
115
116
117 For option 4, the following relation is used, where \f$\phi_4\f$ is
118 the face weight (Note that \f$\phi_4 = 0.5\f$ for orthogonal meshes):
119
120 \f[
121 ffaceWeight = min
122 (
123 max
124 (
125 0.0,
126 (min(\phi_4) - \phi_4)
127 / (min(\phi_4) - max(\phi_4))
128 ),
129 1.0
130 )
131 \f]
132
133
134 For option 5, the following relation is used, where \f$\phi_5\f$ is
135 the cell skewness:
136
137 \f[
138 fskewness =
139 min
140 (
141 max
142 (
143 0.0,
144 (\phi_5 - max(\phi_5))
145 / (min(\phi_5) - max(\phi_5))
146 ),
147 1.0
148 )
149 \f]
150
151
152 For option 6, the following relation is used:
153
154 \f[
155 fCoWeight = clamp((Co - Co1)/(Co2 - Co1), zero_one{});
156 \f]
157
158 where
159 \vartable
160 Co1 | Courant number below which scheme2 is used
161 Co2 | Courant number above which scheme1 is used
162 \endvartable
163
164 The final factor is determined by:
165
166 \f[
167 f = max(fNon, fMagGradCc, fRes, ffaceWeight, fskewness, fCoWeight)
168 \f]
169
170 An indicator (volume) field, named \c blendedIndicator
171 is generated if the log flag is on:
172 - 1 represent scheme1 as active,
173 - 0 represent scheme2 as active.
174
175 Additional reporting is written to the standard output, providing
176 statistics as to the number of cells used by each scheme.
177
178 Operands:
179 \table
180 Operand | Type | Location
181 input | - | -
182 output file | dat | postProcessing/<FO>/<time>/file
183 output field | volScalarField | <time>/outputField
184 \endtable
185
186Usage
187 Minimal example by using \c system/controlDict.functions:
188 \verbatim
189 stabilityBlendingFactorFO
190 {
191 // Mandatory entries
192 type stabilityBlendingFactor;
193 libs (fieldFunctionObjects);
194 field <word>; // U;
195 result <word>; // UBlendingFactor;
196
197 // Optional entries
198 tolerance <scalar>;
199
200 // Conditional entries
201 // Any of the options can be chosen in combinations
202
203 // Option-1
204 switchNonOrtho true;
205 nonOrthogonality nonOrthoAngle;
206 maxNonOrthogonality 20;
207 minNonOrthogonality 60;
208
209 // Option-2
210 switchGradCc true;
211 maxGradCc 3;
212 minGradCc 4;
213
214 // Option-3
215 switchResiduals true;
216 maxResidual 10;
217 residual initialResidual:p;
218 P 1.5;
219 I 0;
220 D 0.5;
221
222 // Option-4
223 switchFaceWeight true;
224 maxFaceWeight 0.3;
225 minFaceWeight 0.2;
226
227 // Option-5
228 switchSkewness true;
229 maxSkewness 2;
230 minSkewness 3;
231
232 // Option-6
233 switchCo true;
234 U U;
235 Co1 1;
236 Co2 10;
237
238 // Inherited entries
239 ...
240 }
241 \endverbatim
242
243 Example of function object specification to calculate the \c residuals used
244 by \c stabilityBlendingFactor. The following writes 'initialResidual:p'
245 field
246 \verbatim
247 residuals
248 {
249 type residuals;
250 libs (utilityFunctionObjects);
251 writeFields true;
252 writeControl writeTime;
253 fields (p);
254 }
255 \endverbatim
256
257 where the entries mean:
258 \table
259 Property | Description | Type | Reqd | Deflt
260 type | Type name: stabilityBlendingFactor | word | yes | -
261 libs | Library name: fieldFunctionObjects | word | yes | -
262 field | Name of operand field | word | yes | -
263 result | Name of surface field to be used in the localBlended scheme <!--
264 --> | word | yes
265 switchNonOrtho | Select non-orthogonal method | bool | no | false
266 nonOrthogonality | Name of the non-orthogonal field <!--
267 --> | word | no | nonOrthoAngle
268 maxNonOrthogonality| Maximum non-orthogonal for scheme2 | scalar | no | 20
269 minNonOrthogonality| Minimum non-orthogonal for scheme1 | scalar | no | 60
270 switchGradCc | Select cell centre gradient method | bool | no | false
271 maxGradCc| Maximum gradient for scheme2 | scalar | no | 2
272 minGradCc| Minimum gradient for scheme1 | scalar | no | 4
273 switchResiduals | Select residual evolution method | bool | no | false
274 residual | Name of the residual field | word | no | initialResidual:p
275 maxResidual| Maximum residual-mean ratio for scheme1 | scalar | no | 10
276 P | Proportional factor for PID | scalar | no | 3
277 I | Integral factor for PID | scalar | no | 0
278 D | Differential factor for PID | scalar | no | 0.25
279 switchFaceWeight | Select face weight method | bool | no | false
280 faceWeight | Name of the faceWeight field | word | no | faceWeight
281 maxFaceWeight | Maximum face weight for scheme1 | scalar | no | 0.2
282 minFaceWeight | Minimum face weight for scheme2 | scalar | no | 0.3
283 switchSkewness | Select skewness method | bool | no | false
284 skewness | Name of the skewness field | word | no | skewness
285 maxSkewness | Maximum skewness for scheme2 | scalar | no | 2
286 minSkewness | Minimum skewness for scheme1 | scalar | no | 3
287 switchCo | Select Co blended method | bool | no | false
288 U | Name of the flux field for Co blended | word | no | U
289 Co1 | Courant number below which scheme2 is used | scalar | no | 1
290 Co2 | Courant number above which scheme1 is used | scalar | no | 10
291 tolerance | Tolerance for number of blended cells | scalar | no | 0.001
292 \endtable
293
294 The \c result entry is the field which is read by the \c localBlended scheme
295 specified in \c fvSchemes. This name is determined by the \c localBlended
296 class.
297
298 The inherited entries are elaborated in:
299 - \link fieldExpression.H \endlink
300 - \link writeFile.H \endlink
301
302SourceFiles
303 stabilityBlendingFactor.C
304
305\*---------------------------------------------------------------------------*/
306
307#ifndef Foam_functionObjects_stabilityBlendingFactor_H
308#define Foam_functionObjects_stabilityBlendingFactor_H
309
310#include "fieldExpression.H"
311#include "writeFile.H"
312#include "volFields.H"
313
314// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315
316namespace Foam
317{
318namespace functionObjects
319{
320
321/*---------------------------------------------------------------------------*\
322 Class stabilityBlendingFactor Declaration
323\*---------------------------------------------------------------------------*/
324
326:
327 public fieldExpression,
328 public writeFile
329{
330 // Private Member Data
331
332 // Switches
333
334 //- Switch for non-orthogonality
335 Switch nonOrthogonality_;
336
337 //- Switch for grad of cell centres
338 Switch gradCc_;
339
340 //- Switch for residuals
341 Switch residuals_;
342
343 //- Switch for face weight
344 Switch faceWeight_;
345
346 //- Switch for skewness
347 Switch skewness_;
348
349 //- Switch for Co
350 Switch Co_;
351
352
353 // Lower and upper limits
354
355 //- Maximum non-orthogonality for fully scheme 2
356 scalar maxNonOrthogonality_;
357
358 //- Minimum non-orthogonality for fully scheme 1
359 scalar minNonOrthogonality_;
360
361 //- Maximum gradcc for fully scheme 2
362 scalar maxGradCc_;
363
364 //- Minimum gradcc for fully scheme 1
365 scalar minGradCc_;
366
367 //- Maximum ratio to average residual for scheme 2
368 scalar maxResidual_;
369
370 //- Minimum face weight for fully scheme 2
371 scalar minFaceWeight_;
372
373 //- Maximum face weight for fully scheme 1
374 scalar maxFaceWeight_;
375
376 //- Maximum skewness for fully scheme 2
377 scalar maxSkewness_;
378
379 //- Minimum skewness for fully scheme 1
380 scalar minSkewness_;
381
382 //- Maximum Co for fully scheme 2
383 scalar Co1_;
384
385 //- Minimum Co for fully scheme 1
386 scalar Co2_;
387
388
389 // File names
390
391 //- Name of the non-orthogonality field
392 word nonOrthogonalityName_;
393
394 //- Name of the face weight field
395 word faceWeightName_;
396
397 //- Name of the skewnes field
398 word skewnessName_;
399
400 //- Name of the residual field
401 word residualName_;
402
403 //- Name of the U used for Co based blended
404 word UName_;
405
406
407 //- Tolerance used when calculating the number of blended cells
408 scalar tolerance_;
409
410
411 //- Error fields
412 scalarField error_;
413 scalarField errorIntegral_;
414 scalarField oldError_;
415 scalarField oldErrorIntegral_;
416
417 //- Proportional gain
418 scalar P_;
419
420 //- Integral gain
421 scalar I_;
422
423 //- Derivative gain
424 scalar D_;
425
426
427 // Private Member Functions
428
429 //- Init fields
430 bool init(bool first);
431
432 //- Return access to the indicator field
433 volScalarField& indicator();
434
435 //- Calculate statistics
436 void calcStats(label&, label&, label&) const ;
437
438 //- Calculate the blending factor field and return true if successful
439 virtual bool calc();
440
441
442protected:
443
444 // Protected Member Functions
445
446 //- Write the file header
447 virtual void writeFileHeader(Ostream& os) const;
448
449
450public:
451
452 //- Runtime type information
453 TypeName("stabilityBlendingFactor");
454
455
456 // Constructors
457
458 //- Construct from name, Time and dictionary
460 (
461 const word& name,
462 const Time& runTime,
463 const dictionary& dict
464 );
465
466 //- No copy construct
468
469 //- No copy assignment
470 void operator=(const stabilityBlendingFactor&) = delete;
471
472
473 //- Destructor
474 virtual ~stabilityBlendingFactor() = default;
475
476
477 // Member Functions
478
479 //- Read the function-object dictionary
480 virtual bool read(const dictionary& dict);
481
482 //- Write the function-object results
483 virtual bool write();
484};
485
486
487// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
488
489} // End namespace functionObjects
490} // End namespace Foam
491
492// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
493
494#endif
495
496// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const word & name() const noexcept
Return the name of this functionObject.
Intermediate class for handling field expression function objects (e.g. blendingFactor etc....
fieldExpression(const word &name, const Time &runTime, const dictionary &dict, const word &fieldName=word::null, const word &resultName=word::null)
Construct from name, Time and dictionary.
virtual bool calc()=0
Calculate the components of the field and return true if successful.
Computes the stabilityBlendingFactor to be used by the local blended convection scheme....
TypeName("stabilityBlendingFactor")
Runtime type information.
virtual ~stabilityBlendingFactor()=default
Destructor.
stabilityBlendingFactor(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
stabilityBlendingFactor(const stabilityBlendingFactor &)=delete
No copy construct.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual void writeFileHeader(Ostream &os) const
Write the file header.
virtual bool write()
Write the function-object results.
void operator=(const stabilityBlendingFactor &)=delete
No copy assignment.
Base class for writing single files from the function objects.
Definition writeFile.H:113
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
A class for handling words, derived from Foam::string.
Definition word.H:66
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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