PureUpwindFitScheme.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-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::PureUpwindFitScheme
28
29Group
30 grpFvSurfaceInterpolationSchemes
31
32Description
33 Upwind biased fit surface interpolation scheme that applies an explicit
34 correction to upwind.
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef PureUpwindFitScheme_H
39#define PureUpwindFitScheme_H
40
41#include "UpwindFitData.H"
42#include "upwind.H"
43#include "Switch.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50/*---------------------------------------------------------------------------*\
51 Class PureUpwindFitScheme Declaration
52\*---------------------------------------------------------------------------*/
53
54template<class Type, class Polynomial, class Stencil>
56:
57 public upwind<Type>
58{
59 // Private Data
60
61 //- Factor the fit is allowed to deviate from linear.
62 // This limits the amount of high-order correction and increases
63 // stability on bad meshes
64 const scalar linearLimitFactor_;
65
66 //- Weights for central stencil
67 const scalar centralWeight_;
68
69
70 // Private Member Functions
71
72 //- No copy construct
74
75 //- No copy assignment
76 void operator=(const PureUpwindFitScheme&) = delete;
77
78
79public:
80
81 //- Runtime type information
82 TypeName("PureUpwindFitScheme");
83
84
85 // Constructors
86
87 //- Construct from mesh and Istream
88 // The name of the flux field is read from the Istream and looked-up
89 // from the mesh objectRegistry
91 :
92 upwind<Type>
93 (
94 mesh,
95 mesh.lookupObject<surfaceScalarField>(word(is))
96 ),
97 linearLimitFactor_(readScalar(is)),
98 centralWeight_(1000)
99 {}
100
101
102 //- Construct from mesh, faceFlux and Istream
104 (
105 const fvMesh& mesh,
106 const surfaceScalarField& faceFlux,
107 Istream& is
108 )
109 :
110 upwind<Type>(mesh, faceFlux),
111 linearLimitFactor_(readScalar(is)),
112 centralWeight_(1000)
113 {}
114
115
116 // Member Functions
117
118 //- Return true if this scheme uses an explicit correction
119 virtual bool corrected() const
120 {
121 return true;
122 }
123
124 //- Return the explicit correction to the face-interpolate
127 (
129 ) const
130 {
131 const fvMesh& mesh = this->mesh();
132
133 // Use the owner/neighbour splitting constructor
135
136 const UpwindFitData<Polynomial>& ufd =
138 (
139 mesh,
140 stencil,
141 false, //offset to upwind
142 linearLimitFactor_,
143 centralWeight_
144 );
145
146 const List<scalarList>& fo = ufd.owncoeffs();
147 const List<scalarList>& fn = ufd.neicoeffs();
148
149 return stencil.weightedSum(this->faceFlux_, vf, fo, fn);
150 }
151};
152
153
154// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155
156} // End namespace Foam
157
158// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159
160// Add the patch constructor functions to the hash tables
162#define makePureUpwindFitSurfaceInterpolationTypeScheme\
163( \
164 SS, \
165 POLYNOMIAL, \
166 STENCIL, \
167 TYPE \
168) \
169 \
170typedef PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
171 PureUpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
172defineTemplateTypeNameAndDebugWithName \
173 (PureUpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
174 \
175surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
176<PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
177 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
178 \
179surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
180<PureUpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
181 add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
183#define makePureUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
184 \
185makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
186makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
187makePureUpwindFitSurfaceInterpolationTypeScheme \
188( \
189 SS, \
190 POLYNOMIAL, \
191 STENCIL, \
192 sphericalTensor \
193) \
194makePureUpwindFitSurfaceInterpolationTypeScheme \
195( \
196 SS, \
197 POLYNOMIAL, \
198 STENCIL, \
199 symmTensor \
200) \
201makePureUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
202
203
204// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205
206#endif
207
208// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
Upwind biased fit surface interpolation scheme that applies an explicit correction to upwind.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
TypeName("PureUpwindFitScheme")
Runtime type information.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
PureUpwindFitScheme(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from mesh, faceFlux and Istream.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Data for the quadratic fit correction interpolation scheme to be used with upwind biased stencil.
Definition: UpwindFitData.H:65
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
const List< scalarList > & neicoeffs() const
Return reference to neighbour fit coefficients.
Creates upwind stencil by shifting a centred stencil to upwind and downwind faces and optionally remo...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > weightedSum(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< scalar > > &ownWeights, const List< List< scalar > > &neiWeights) const
Sum vol field contributions to create face values.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
Definition: tmp.H:65
Upwind differencing scheme class.
Definition: upwind.H:59
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73