All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CentredFitSnGradScheme.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) 2013-2019 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::fv::CentredFitSnGradScheme
26 
27 Description
28  Centred fit snGrad scheme which applies an explicit correction to snGrad
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #ifndef CentredFitSnGradScheme_H
33 #define CentredFitSnGradScheme_H
34 
35 #include "CentredFitSnGradData.H"
36 #include "snGradScheme.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 
43 class fvMesh;
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace fv
48 {
49 /*---------------------------------------------------------------------------*\
50  Class CentredFitSnGradSnGradScheme Declaration
51 \*---------------------------------------------------------------------------*/
52 
53 template<class Type, class Polynomial, class Stencil>
55 :
56  public snGradScheme<Type>
57 {
58  // Private Data
59 
60  //- Factor the fit is allowed to deviate from linear.
61  // This limits the amount of high-order correction and increases
62  // stability on bad meshes
63  const scalar linearLimitFactor_;
64 
65  //- Weights for central stencil
66  const scalar centralWeight_;
67 
68 
69  // Private Member Functions
70 
71  //- Disallow default bitwise copy construction
73 
74  //- Disallow default bitwise assignment
75  void operator=(const CentredFitSnGradScheme&) = delete;
76 
77 
78 public:
79 
80  //- Runtime type information
81  TypeName("CentredFitSnGradScheme");
82 
83 
84  // Constructors
85 
86  //- Construct from mesh and Istream
88  :
89  snGradScheme<Type>(mesh),
90  linearLimitFactor_(readScalar(is)),
91  centralWeight_(1000)
92  {}
93 
94 
95  // Member Functions
96 
97  //- Return the interpolation weighting factors for the given field
99  (
101  ) const
102  {
103  return this->mesh().nonOrthDeltaCoeffs();
104  }
105 
106  //- Return true if this scheme uses an explicit correction
107  virtual bool corrected() const
108  {
109  return true;
110  }
111 
112  //- Return the explicit correction to the face-interpolate
115  (
117  ) const
118  {
119  const fvMesh& mesh = this->mesh();
120 
122  (
123  mesh
124  );
125 
128  (
129  mesh,
130  stencil,
131  linearLimitFactor_,
132  centralWeight_
133  );
134 
136  (
137  stencil.weightedSum(vf, cfd.coeffs())
138  );
139 
140  sft.ref().dimensions() /= dimLength;
141 
142  return sft;
143  }
144 };
145 
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 } // End namespace fv
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace Foam
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 // Add the patch constructor functions to the hash tables
159 #define makeCentredFitSnGradTypeScheme(SS, POLYNOMIAL, STENCIL, TYPE) \
160  typedef Foam::fv::CentredFitSnGradScheme \
161  <Foam::TYPE, Foam::POLYNOMIAL, Foam::STENCIL> \
162  CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_; \
163  \
164  defineTemplateTypeNameAndDebugWithName \
165  (CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
166  \
167  namespace Foam \
168  { \
169  namespace fv \
170  { \
171  snGradScheme<TYPE>::addMeshConstructorToTable \
172  <CentredFitSnGradScheme<TYPE, POLYNOMIAL, STENCIL>> \
173  add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
174  } \
175  }
177 #define makeCentredFitSnGradScheme(SS, POLYNOMIAL, STENCIL) \
178  \
179  makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
180  makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
181  makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,sphericalTensor) \
182  makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
183  makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
184 
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #endif
189 
190 // ************************************************************************* //
const List< scalarList > & coeffs() const
Return reference to fit coefficients.
TypeName("CentredFitSnGradScheme")
Runtime type information.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Generic GeometricField class.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > weightedSum(const GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< scalar >> &stencilWeights) const
Sum vol field contributions to create face values.
static const CentredFitSnGradData< Polynomial > & New(const fvMesh &mesh)
Definition: MeshObject.C:44
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
labelList fv(nPoints)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
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.
Centred fit snGrad scheme which applies an explicit correction to snGrad.
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:60
const fvMesh & mesh() const
Return mesh reference.
Definition: snGradScheme.H:117
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Data for centred fit snGrad schemes.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.