CentredFitScheme.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::CentredFitScheme
26 
27 Description
28  Centred fit surface interpolation scheme which applies an explicit
29  correction to linear.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #ifndef CentredFitScheme_H
34 #define CentredFitScheme_H
35 
36 #include "CentredFitData.H"
37 #include "linear.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 /*---------------------------------------------------------------------------*\
45  Class CentredFitScheme Declaration
46 \*---------------------------------------------------------------------------*/
47 
48 template<class Type, class Polynomial, class Stencil>
49 class CentredFitScheme
50 :
51  public linear<Type>
52 {
53  // Private Data
54 
55  //- Factor the fit is allowed to deviate from linear.
56  // This limits the amount of high-order correction and increases
57  // stability on bad meshes
58  const scalar linearLimitFactor_;
59 
60  //- Weights for central stencil
61  const scalar centralWeight_;
62 
63 
64 public:
65 
66  //- Runtime type information
67  TypeName("CentredFitScheme");
68 
69 
70  // Constructors
71 
72  //- Construct from mesh and Istream
73  CentredFitScheme(const fvMesh& mesh, Istream& is)
74  :
75  linear<Type>(mesh),
76  linearLimitFactor_(readScalar(is)),
77  centralWeight_(1000)
78  {}
79 
80 
81  //- Construct from mesh, faceFlux and Istream
83  (
84  const fvMesh& mesh,
85  const surfaceScalarField& faceFlux,
86  Istream& is
87  )
88  :
89  linear<Type>(mesh),
90  linearLimitFactor_(readScalar(is)),
91  centralWeight_(1000)
92  {}
93 
94  //- Disallow default bitwise copy construction
95  CentredFitScheme(const CentredFitScheme&) = delete;
96 
97 
98  // Member Functions
99 
100  //- Return true if this scheme uses an explicit correction
101  virtual bool corrected() const
102  {
103  return true;
104  }
105 
106  //- Return the explicit correction to the face-interpolate
109  (
111  ) const
112  {
113  const fvMesh& mesh = this->mesh();
114 
116  (
117  mesh
118  );
119 
120  const CentredFitData<Polynomial>& cfd =
122  (
123  mesh,
124  stencil,
125  linearLimitFactor_,
126  centralWeight_
127  );
128 
129  const List<scalarList>& f = cfd.coeffs();
130 
131  return stencil.weightedSum(vf, f);
132  }
133 
134 
135  // Member Operators
136 
137  //- Disallow default bitwise assignment
138  void operator=(const CentredFitScheme&) = delete;
139 };
140 
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
144 } // End namespace Foam
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 // Add the patch constructor functions to the hash tables
150 #define makeCentredFitSurfaceInterpolationTypeScheme\
151 ( \
152  SS, \
153  POLYNOMIAL, \
154  STENCIL, \
155  TYPE \
156 ) \
157  \
158 typedef CentredFitScheme<TYPE, POLYNOMIAL, STENCIL> \
159  CentredFitScheme##TYPE##POLYNOMIAL##STENCIL##_; \
160 defineTemplateTypeNameAndDebugWithName \
161  (CentredFitScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
162  \
163 surfaceInterpolationScheme<TYPE>::addMeshConstructorToTable \
164 <CentredFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
165  add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
166  \
167 surfaceInterpolationScheme<TYPE>::addMeshFluxConstructorToTable \
168 <CentredFitScheme<TYPE, POLYNOMIAL, STENCIL>> \
169  add##SS##STENCIL##TYPE##MeshFluxConstructorToTable_;
171 #define makeCentredFitSurfaceInterpolationScheme(SS, POLYNOMIAL, STENCIL) \
172  \
173 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
174 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
175 makeCentredFitSurfaceInterpolationTypeScheme \
176 ( \
177  SS, \
178  POLYNOMIAL, \
179  STENCIL, \
180  sphericalTensor \
181 ) \
182 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor)\
183 makeCentredFitSurfaceInterpolationTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
184 
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
188 #endif
189 
190 // ************************************************************************* //
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Centred interpolation 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 List< scalarList > & coeffs() const
Return reference to fit coefficients.
const fvMesh & mesh() const
Return mesh reference.
TypeName("CentredFitScheme")
Runtime type information.
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.
Data for the quadratic fit correction interpolation scheme.
void operator=(const CentredFitScheme &)=delete
Disallow default bitwise assignment.
Centred fit surface interpolation scheme which applies an explicit correction to linear.
CentredFitScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
static CentredFitData< Polynomial > & New(fvMesh &mesh)
Definition: MeshObject.C:54
labelList f(nPoints)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.