Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
src
finiteVolume
interpolation
surfaceInterpolation
limitedSchemes
blended
blended.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
License
9
This file is part of OpenFOAM.
10
11
OpenFOAM is free software: you can redistribute it and/or modify it
12
under the terms of the GNU General Public License as published by
13
the Free Software Foundation, either version 3 of the License, or
14
(at your option) any later version.
15
16
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19
for more details.
20
21
You should have received a copy of the GNU General Public License
22
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23
24
Class
25
Foam::blended
26
27
Description
28
linear/upwind blended differencing scheme.
29
30
SourceFiles
31
blended.C
32
33
\*---------------------------------------------------------------------------*/
34
35
#ifndef blended_H
36
#define blended_H
37
38
#include <
finiteVolume/limitedSurfaceInterpolationScheme.H
>
39
40
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42
namespace
Foam
43
{
44
45
/*---------------------------------------------------------------------------*\
46
Class blended Declaration
47
\*---------------------------------------------------------------------------*/
48
49
template
<
class
Type>
50
class
blended
51
:
52
public
limitedSurfaceInterpolationScheme
<Type>
53
{
54
// Private data
55
56
const
scalar blendingFactor_;
57
58
59
// Private Member Functions
60
61
//- Disallow default bitwise copy construct
62
blended
(
const
blended
&);
63
64
//- Disallow default bitwise assignment
65
void
operator=(
const
blended
&);
66
67
68
public
:
69
70
//- Runtime type information
71
TypeName
(
"blended"
);
72
73
74
// Constructors
75
76
//- Construct from mesh, faceFlux and blendingFactor
77
blended
78
(
79
const
fvMesh
&
mesh
,
80
const
surfaceScalarField
& faceFlux,
81
const
scalar blendingFactor
82
)
83
:
84
limitedSurfaceInterpolationScheme<Type>
(
mesh
, faceFlux),
85
blendingFactor_(blendingFactor)
86
{}
87
88
//- Construct from mesh and Istream.
89
// The name of the flux field is read from the Istream and looked-up
90
// from the mesh objectRegistry
91
blended
92
(
93
const
fvMesh
& mesh,
94
Istream
& is
95
)
96
:
97
limitedSurfaceInterpolationScheme<Type>
(
mesh
, is),
98
blendingFactor_(
readScalar
(is))
99
{}
100
101
//- Construct from mesh, faceFlux and Istream
102
blended
103
(
104
const
fvMesh
& mesh,
105
const
surfaceScalarField
& faceFlux,
106
Istream
& is
107
)
108
:
109
limitedSurfaceInterpolationScheme<Type>
(
mesh
, faceFlux),
110
blendingFactor_(
readScalar
(is))
111
{}
112
113
114
// Member Functions
115
116
//- Return the interpolation limiter
117
virtual
tmp<surfaceScalarField>
limiter
118
(
119
const
GeometricField<Type, fvPatchField, volMesh>
&
120
)
const
121
{
122
return
tmp<surfaceScalarField>
123
(
124
new
surfaceScalarField
125
(
126
IOobject
127
(
128
"blendedLimiter"
,
129
this->
mesh
().time().
timeName
(),
130
this->
mesh
()
131
),
132
this->
mesh
(),
133
dimensionedScalar
134
(
135
"blendedLimiter"
,
136
dimless
,
137
1 - blendingFactor_
138
)
139
)
140
);
141
}
142
143
//- Return the interpolation weighting factors
144
virtual
tmp<surfaceScalarField>
weights
145
(
146
const
GeometricField<Type, fvPatchField, volMesh>
& vf
147
)
const
148
{
149
return
150
blendingFactor_*this->
mesh
().surfaceInterpolation::weights()
151
+ (1 - blendingFactor_)*
pos
(this->
faceFlux_
);
152
}
153
};
154
155
156
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158
}
// End namespace Foam
159
160
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162
#endif
163
164
// ************************ vim: set sw=4 sts=4 et: ************************ //