CoBlended.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) 2012-2016 OpenFOAM Foundation
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::CoBlended
28
29Group
30 grpFvSurfaceInterpolationSchemes
31
32Description
33 Two-scheme Courant number based blending differencing scheme.
34
35 Similar to localBlended but uses a blending factor computed from the
36 face-based Courant number and the lower and upper Courant number limits
37 supplied:
38 \f[
39 weight = 1 - max(min((Co - Co1)/(Co2 - Co1), 1), 0)
40 \f]
41 where
42 \vartable
43 Co1 | Courant number below which scheme1 is used
44 Co2 | Courant number above which scheme2 is used
45 \endvartable
46
47 The weight applies to the first scheme and 1-weight to the second scheme.
48
49Usage
50 Example of the CoBlended scheme specification using LUST for Courant numbers
51 less than 1 and linearUpwind for Courant numbers greater than 10:
52 \verbatim
53 divSchemes
54 {
55 .
56 .
57 div(phi,U) Gauss CoBlended 1 LUST grad(U) 10 linearUpwind grad(U);
58 .
59 .
60 }
61 \endverbatim
62
63SourceFiles
64 CoBlended.C
65
66\*---------------------------------------------------------------------------*/
67
68#ifndef CoBlended_H
69#define CoBlended_H
70
72#include "blendedSchemeBase.H"
73#include "surfaceInterpolate.H"
74
75// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77namespace Foam
78{
79
80/*---------------------------------------------------------------------------*\
81 Class CoBlended Declaration
82\*---------------------------------------------------------------------------*/
83
84template<class Type>
85class CoBlended
86:
87 public surfaceInterpolationScheme<Type>,
88 public blendedSchemeBase<Type>
89{
90 // Private data
91
92 //- Courant number below which scheme1 is used
93 const scalar Co1_;
94
95 //- Scheme 1
97
98 //- Courant number above which scheme2 is used
99 const scalar Co2_;
100
101 //- Scheme 2
103
104 //- The face-flux used to compute the face Courant number
105 const surfaceScalarField& faceFlux_;
106
107
108 // Private Member Functions
109
110 //- No copy construct
111 CoBlended(const CoBlended&) = delete;
112
113 //- No copy assignment
114 void operator=(const CoBlended&) = delete;
115
116
117public:
118
119 //- Runtime type information
120 TypeName("CoBlended");
121
122
123 // Constructors
124
125 //- Construct from mesh and Istream.
126 // The name of the flux field is read from the Istream and looked-up
127 // from the mesh objectRegistry
129 (
130 const fvMesh& mesh,
131 Istream& is
132 )
133 :
135 Co1_(readScalar(is)),
136 tScheme1_
137 (
139 ),
140 Co2_(readScalar(is)),
141 tScheme2_
142 (
144 ),
145 faceFlux_
146 (
147 mesh.lookupObject<surfaceScalarField>(word(is))
148 )
149 {
150 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
151 {
153 << "coefficients = " << Co1_ << " and " << Co2_
154 << " should be > 0 and Co2 > Co1"
155 << exit(FatalIOError);
156 }
157 }
158
159
160 //- Construct from mesh, faceFlux and Istream
162 (
163 const fvMesh& mesh,
164 const surfaceScalarField& faceFlux,
165 Istream& is
166 )
167 :
169 Co1_(readScalar(is)),
170 tScheme1_
171 (
172 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
173 ),
174 Co2_(readScalar(is)),
175 tScheme2_
176 (
177 surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
178 ),
179 faceFlux_(faceFlux)
180 {
181 if (Co1_ < 0 || Co2_ < 0 || Co1_ >= Co2_)
182 {
184 << "coefficients = " << Co1_ << " and " << Co2_
185 << " should be > 0 and Co2 > Co1"
186 << exit(FatalIOError);
187 }
188 }
189
190
191 // Member Functions
192
193 //- Return the face-based blending factor
195 (
197 ) const
198 {
199 const fvMesh& mesh = this->mesh();
200 tmp<surfaceScalarField> tUflux = faceFlux_;
202 if (faceFlux_.dimensions() == dimDensity*dimVelocity*dimArea)
203 {
204 // Currently assume that the density field
205 // corresponding to the mass-flux is named "rho"
206 const volScalarField& rho =
207 mesh.objectRegistry::template lookupObject<volScalarField>
208 ("rho");
209
210 tUflux = faceFlux_/fvc::interpolate(rho);
211 }
212 else if (faceFlux_.dimensions() != dimVelocity*dimArea)
213 {
215 << "dimensions of faceFlux are not correct"
216 << exit(FatalError);
217 }
218
219 return tmp<surfaceScalarField>
220 (
222 (
223 vf.name() + "BlendingFactor",
224 scalar(1)
225 - max
226 (
227 min
228 (
229 (
231 *mag(tUflux)/mesh.magSf()
232 - Co1_
233 )/(Co2_ - Co1_),
234 scalar(1)
235 ),
236 scalar(0)
237 )
238 )
239 );
240 }
241
242
243 //- Return the interpolation weighting factors
244 tmp<surfaceScalarField>
245 weights
246 (
247 const GeometricField<Type, fvPatchField, volMesh>& vf
248 ) const
249 {
251
252 return
253 bf*tScheme1_().weights(vf)
254 + (scalar(1) - bf)*tScheme2_().weights(vf);
255 }
256
257
258 //- Return the face-interpolate of the given cell field
259 // with explicit correction
262 (
264 ) const
265 {
267
268 return
269 bf*tScheme1_().interpolate(vf)
270 + (scalar(1) - bf)*tScheme2_().interpolate(vf);
271 }
272
273
274 //- Return true if this scheme uses an explicit correction
275 virtual bool corrected() const
276 {
277 return tScheme1_().corrected() || tScheme2_().corrected();
278 }
279
280
281 //- Return the explicit correction to the face-interpolate
282 // for the given field
285 (
287 ) const
288 {
290
291 if (tScheme1_().corrected())
292 {
293 if (tScheme2_().corrected())
294 {
295 return
296 (
297 bf
298 * tScheme1_().correction(vf)
299 + (scalar(1) - bf)
300 * tScheme2_().correction(vf)
301 );
302 }
303 else
304 {
305 return
306 (
307 bf
308 * tScheme1_().correction(vf)
309 );
310 }
311 }
312 else if (tScheme2_().corrected())
313 {
314 return
315 (
316 (scalar(1) - bf)
317 * tScheme2_().correction(vf)
318 );
319 }
320 else
321 {
322 return tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
323 (
324 nullptr
325 );
326 }
327 }
328};
329
330
331// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332
333} // End namespace Foam
334
335// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336
337#endif
338
339// ************************************************************************* //
Two-scheme Courant number based blending differencing scheme.
Definition: CoBlended.H:96
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: CoBlended.H:282
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: CoBlended.H:253
TypeName("CoBlended")
Runtime type information.
CoBlended(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
Definition: CoBlended.H:169
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: CoBlended.H:292
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-interpolate of the given cell field.
Definition: CoBlended.H:269
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Definition: CoBlended.H:202
const dimensionSet & dimensions() const
Return dimensions.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
dimensionedScalar deltaT() const
Return time step.
Definition: TimeStateI.H:55
Base class for blended schemes to provide access to the blending factor surface field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
static tmp< surfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73