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