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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::blended
29 
30 Group
31  grpFvLimitedSurfaceInterpolationSchemes
32 
33 Description
34  linear/upwind blended differencing scheme.
35 
36 SourceFiles
37  blended.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef blended_H
42 #define blended_H
43 
45 #include "blendedSchemeBase.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class blended Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 template<class Type>
57 class blended
58 :
60  public blendedSchemeBase<Type>
61 {
62  // Private data
63 
64  const scalar blendingFactor_;
65 
66 
67  // Private Member Functions
68 
69  //- No copy construct
70  blended(const blended&) = delete;
71 
72  //- No copy assignment
73  void operator=(const blended&) = delete;
74 
75 
76 public:
77 
78  //- Runtime type information
79  TypeName("blended");
80 
81 
82  // Constructors
83 
84  //- Construct from mesh, faceFlux and blendingFactor
86  (
87  const fvMesh& mesh,
88  const surfaceScalarField& faceFlux,
89  const scalar blendingFactor
90  )
91  :
92  limitedSurfaceInterpolationScheme<Type>(mesh, faceFlux),
93  blendingFactor_(blendingFactor)
94  {}
95 
96  //- Construct from mesh and Istream.
97  // The name of the flux field is read from the Istream and looked-up
98  // from the mesh objectRegistry
100  (
101  const fvMesh& mesh,
102  Istream& is
103  )
104  :
106  blendingFactor_(readScalar(is))
107  {}
108 
109  //- Construct from mesh, faceFlux and Istream
111  (
112  const fvMesh& mesh,
113  const surfaceScalarField& faceFlux,
114  Istream& is
115  )
116  :
117  limitedSurfaceInterpolationScheme<Type>(mesh, faceFlux),
118  blendingFactor_(readScalar(is))
119  {}
120 
121 
122  //- Destructor
123  virtual ~blended() = default;
124 
125 
126  // Member Functions
127 
128  //- Return the face-based blending factor
130  (
132  ) const
133  {
135  (
137  (
138  IOobject
139  (
140  vf.name() + "BlendingFactor",
141  this->mesh().time().timeName(),
142  this->mesh(),
145  false
146  ),
147  this->mesh(),
149  (
150  "blendingFactor",
151  dimless,
152  blendingFactor_
153  )
154  )
155  );
156  }
157 
158  //- Return the interpolation limiter
160  (
162  ) const
163  {
165  (
167  (
168  IOobject
169  (
170  "blendedLimiter",
171  this->mesh().time().timeName(),
172  this->mesh(),
175  false
176  ),
177  this->mesh(),
179  (
180  "blendedLimiter",
181  dimless,
182  1 - blendingFactor_
183  )
184  )
185  );
186  }
187 
188  //- Return the interpolation weighting factors
190  (
192  ) const
193  {
194  return
195  blendingFactor_*this->mesh().surfaceInterpolation::weights()
196  + (1 - blendingFactor_)*pos0(this->faceFlux_);
197  }
198 };
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 } // End namespace Foam
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 #endif
208 
209 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
blendedSchemeBase.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::blended::weights
virtual tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors.
Definition: blended.H:189
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::blended::TypeName
TypeName("blended")
Runtime type information.
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::blended::blendingFactor
virtual tmp< surfaceScalarField > blendingFactor(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the face-based blending factor.
Definition: blended.H:129
Foam::blended::~blended
virtual ~blended()=default
Destructor.
Foam::blended
linear/upwind blended differencing scheme.
Definition: blended.H:56
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::blended::limiter
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation limiter.
Definition: blended.H:159
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::limitedSurfaceInterpolationScheme::faceFlux_
const surfaceScalarField & faceFlux_
Definition: limitedSurfaceInterpolationScheme.H:74
Foam::limitedSurfaceInterpolationScheme
Abstract base class for limited surface interpolation schemes.
Definition: limitedSurfaceInterpolationScheme.H:54
limitedSurfaceInterpolationScheme.H
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::surfaceInterpolationScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: surfaceInterpolationScheme.H:144
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::blendedSchemeBase
Base class for blended schemes to provide access to the blending factor surface field.
Definition: blendedSchemeBase.H:55
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189