UpwindFitScheme.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 -------------------------------------------------------------------------------
10 License
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 
26 Class
27  Foam::UpwindFitScheme
28 
29 Group
30  grpFvSurfaceInterpolationSchemes
31 
32 Description
33  Upwind biased fit surface interpolation scheme that applies an explicit
34  correction to linear.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef UpwindFitScheme_H
39 #define UpwindFitScheme_H
40 
41 #include "UpwindFitData.H"
42 #include "linear.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class UpwindFitScheme Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class Type, class Polynomial, class Stencil>
54 class UpwindFitScheme
55 :
56  public linear<Type>
57 {
58  // Private Data
59 
60  //- Reference to the surface flux used to choose upwind direction
61  const surfaceScalarField& faceFlux_;
62 
63  //- Factor the fit is allowed to deviate from linear.
64  // This limits the amount of high-order correction and increases
65  // stability on bad meshes
66  const scalar linearLimitFactor_;
67 
68  //- Weights for central stencil
69  const scalar centralWeight_;
70 
71 
72  // Private Member Functions
73 
74  //- No copy construct
75  UpwindFitScheme(const UpwindFitScheme&) = delete;
76 
77  //- No copy assignment
78  void operator=(const UpwindFitScheme&) = delete;
79 
80 
81 public:
82 
83  //- Runtime type information
84  TypeName("UpwindFitScheme");
85 
86 
87  // Constructors
88 
89  //- Construct from mesh and Istream
90  // The name of the flux field is read from the Istream and looked-up
91  // from the mesh objectRegistry
92  UpwindFitScheme(const fvMesh& mesh, Istream& is)
93  :
95  faceFlux_(mesh.lookupObject<surfaceScalarField>(word(is))),
96  linearLimitFactor_(readScalar(is)),
97  centralWeight_(1000)
98  {}
99 
100 
101  //- Construct from mesh, faceFlux and Istream
103  (
104  const fvMesh& mesh,
105  const surfaceScalarField& faceFlux,
106  Istream& is
107  )
108  :
109  linear<Type>(mesh),
110  faceFlux_(faceFlux),
111  linearLimitFactor_(readScalar(is)),
112  centralWeight_(1000)
113  {}
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
126  correction
127  (
129  ) const
130  {
131  const fvMesh& mesh = this->mesh();
134  (
135  mesh,
136  false, //pureUpwind
137  scalar(0.5)
138  );
139 
140  const UpwindFitData<Polynomial>& ufd =
142  (
143  mesh,
144  stencil,
145  true, //calculate as offset to linear
146  linearLimitFactor_,
147  centralWeight_
148  );
149 
150  const List<scalarList>& fo = ufd.owncoeffs();
151  const List<scalarList>& fn = ufd.neicoeffs();
152 
153  return stencil.weightedSum(faceFlux_, vf, fo, fn);
154  }
155 };
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace Foam
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 // Add the patch constructor functions to the hash tables
165 
166 #define makeUpwindFitSurfaceInterpolationTypeScheme\
167 ( \
168  SS, \
169  POLYNOMIAL, \
170  STENCIL, \
171  TYPE \
172 ) \
173  \
174 typedef UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
175  UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
176 defineTemplateTypeNameAndDebugWithName \
177  (UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
178  \
179 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
180 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
181  add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
182  \
183 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
184 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
185  add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
186 
187 #define makeUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
188  \
189 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
190 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
191 makeUpwindFitSurfaceInterpolationTypeScheme \
192 ( \
193  SS, \
194  POLYNOMIAL, \
195  STENCIL, \
196  sphericalTensor \
197 ) \
198 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
199 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 #endif
205 
206 // ************************************************************************* //
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
Central-differencing interpolation scheme class.
Definition: linear.H:51
static const Type & New(const Mesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
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 tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
const List< scalarList > & neicoeffs() const
Return reference to neighbour fit coefficients.
TypeName("UpwindFitScheme")
Runtime type information.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Upwind biased fit surface interpolation scheme that applies an explicit correction to linear...
const List< scalarList > & owncoeffs() const
Return reference to owner fit coefficients.
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:53
Creates upwind stencil by shifting a centred stencil to upwind and downwind faces and optionally remo...
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Namespace for OpenFOAM.