UpwindFitScheme.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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::UpwindFitScheme
26 
27 Description
28  Upwind biased fit surface interpolation scheme that applies an explicit
29  correction to linear.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef UpwindFitScheme_H
34 #define UpwindFitScheme_H
35 
36 #include "UpwindFitData.H"
37 #include "linear.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class UpwindFitScheme Declaration
46 \*---------------------------------------------------------------------------*/
47 
48 template<class Type, class Polynomial, class Stencil>
49 class UpwindFitScheme
50 :
51  public linear<Type>
52 {
53  // Private Data
54 
55  //- Reference to the surface flux used to choose upwind direction
56  const surfaceScalarField& faceFlux_;
57 
58  //- Factor the fit is allowed to deviate from linear.
59  // This limits the amount of high-order correction and increases
60  // stability on bad meshes
61  const scalar linearLimitFactor_;
62 
63  //- Weights for central stencil
64  const scalar centralWeight_;
65 
66 
67 public:
68 
69  //- Runtime type information
70  TypeName("UpwindFitScheme");
71 
72 
73  // Constructors
74 
75  //- Construct from mesh and Istream
76  // The name of the flux field is read from the Istream and looked-up
77  // from the mesh objectRegistry
78  UpwindFitScheme(const fvMesh& mesh, Istream& is)
79  :
80  linear<Type>(mesh),
81  faceFlux_(mesh.lookupObject<surfaceScalarField>(word(is))),
82  linearLimitFactor_(readScalar(is)),
83  centralWeight_(1000)
84  {}
85 
86 
87  //- Construct from mesh, faceFlux and Istream
89  (
90  const fvMesh& mesh,
91  const surfaceScalarField& faceFlux,
92  Istream& is
93  )
94  :
95  linear<Type>(mesh),
96  faceFlux_(faceFlux),
97  linearLimitFactor_(readScalar(is)),
98  centralWeight_(1000)
99  {}
100 
101  //- Disallow default bitwise copy construction
102  UpwindFitScheme(const UpwindFitScheme&) = delete;
103 
104 
105  // Member Functions
106 
107  //- Return true if this scheme uses an explicit correction
108  virtual bool corrected() const
109  {
110  return true;
111  }
112 
113  //- Return the explicit correction to the face-interpolate
116  (
118  ) const
119  {
120  const fvMesh& mesh = this->mesh();
121 
123  (
124  mesh,
125  false, // pureUpwind
126  scalar(0.5)
127  );
128 
129  const UpwindFitData<Polynomial>& ufd =
131  (
132  mesh,
133  stencil,
134  true, // calculate as offset to linear
135  linearLimitFactor_,
136  centralWeight_
137  );
138 
139  const List<scalarList>& fo = ufd.owncoeffs();
140  const List<scalarList>& fn = ufd.neicoeffs();
141 
142  return stencil.weightedSum(faceFlux_, vf, fo, fn);
143  }
144 
145 
146  // Member Operators
147 
148  //- Disallow default bitwise assignment
149  void operator=(const UpwindFitScheme&) = delete;
150 };
151 
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 } // End namespace Foam
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 // Add the patch constructor functions to the hash tables
161 #define makeUpwindFitSurfaceInterpolationTypeScheme\
162 ( \
163  SS, \
164  POLYNOMIAL, \
165  STENCIL, \
166  TYPE \
167 ) \
168  \
169 typedef UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL> \
170  UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
171 defineTemplateTypeNameAndDebugWithName \
172  (UpwindFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
173  \
174 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
175 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
176  add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
177  \
178 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
179 <UpwindFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
180  add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
182 #define makeUpwindFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
183  \
184 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
185 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
186 makeUpwindFitSurfaceInterpolationTypeScheme \
187 ( \
188  SS, \
189  POLYNOMIAL, \
190  STENCIL, \
191  sphericalTensor \
192 ) \
193 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
194 makeUpwindFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 #endif
200 
201 // ************************************************************************* //
Central-differencing interpolation scheme class.
Definition: linear.H:50
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.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
Upwind biased fit surface interpolation scheme that applies an explicit correction to linear...
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 successful.
Definition: doubleScalar.H:75
void operator=(const UpwindFitScheme &)=delete
Disallow default bitwise assignment.
UpwindFitScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
static UpwindFitData< Polynomial > & New(fvMesh &mesh)
Definition: MeshObject.C:54
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...
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: PtrList.H:53
Namespace for OpenFOAM.