CentredFitSnGradScheme.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::CentredFitSnGradScheme
29
30Group
31 grpFvSnGradSchemes
32
33Description
34 Centred fit snGrad scheme which applies an explicit correction to snGrad
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef CentredFitSnGradScheme_H
39#define CentredFitSnGradScheme_H
40
42#include "snGradScheme.H"
43
44// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46namespace Foam
47{
48
49class fvMesh;
50
51// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53namespace fv
54{
55/*---------------------------------------------------------------------------*\
56 Class CentredFitSnGradSnGradScheme Declaration
57\*---------------------------------------------------------------------------*/
58
59template<class Type, class Polynomial, class Stencil>
61:
62 public snGradScheme<Type>
63{
64 // Private Data
65
66 //- Factor the fit is allowed to deviate from linear.
67 // This limits the amount of high-order correction and increases
68 // stability on bad meshes
69 const scalar linearLimitFactor_;
70
71 //- Weights for central stencil
72 const scalar centralWeight_;
73
74
75 // Private Member Functions
76
77 //- No copy construct
79
80 //- No copy assignment
81 void operator=(const CentredFitSnGradScheme&) = delete;
82
83
84public:
85
86 //- Runtime type information
87 TypeName("CentredFitSnGradScheme");
88
89
90 // Constructors
91
92 //- Construct from mesh and Istream
94 :
95 snGradScheme<Type>(mesh),
96 linearLimitFactor_(readScalar(is)),
97 centralWeight_(1000)
98 {}
99
100
101 // Member Functions
102
103 //- Return the interpolation weighting factors for the given field
105 (
107 ) const
108 {
109 return this->mesh().nonOrthDeltaCoeffs();
110 }
111
112 //- Return true if this scheme uses an explicit correction
113 virtual bool corrected() const noexcept
114 {
115 return true;
116 }
117
118 //- Return the explicit correction to the face-interpolate
121 (
123 ) const
124 {
125 const fvMesh& mesh = this->mesh();
126
128 (
129 mesh
130 );
131
134 (
135 mesh,
136 stencil,
137 linearLimitFactor_,
138 centralWeight_
139 );
140
142 (
143 stencil.weightedSum(vf, cfd.coeffs())
144 );
145
146 sft.ref().dimensions() /= dimLength;
147
148 return sft;
149 }
150};
151
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155} // End namespace fv
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159} // End namespace Foam
160
161// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162
163// Add the patch constructor functions to the hash tables
165#define makeCentredFitSnGradTypeScheme(SS, POLYNOMIAL, STENCIL, TYPE) \
166 typedef Foam::fv::CentredFitSnGradScheme \
167 <Foam::TYPE, Foam::POLYNOMIAL, Foam::STENCIL> \
168 CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_; \
169 \
170 defineTemplateTypeNameAndDebugWithName \
171 (CentredFitSnGradScheme##TYPE##POLYNOMIAL##STENCIL##_, #SS, 0); \
172 \
173 namespace Foam \
174 { \
175 namespace fv \
176 { \
177 snGradScheme<TYPE>::addMeshConstructorToTable \
178 <CentredFitSnGradScheme<TYPE, POLYNOMIAL, STENCIL>> \
179 add##SS##STENCIL##TYPE##MeshConstructorToTable_; \
180 } \
181 }
183#define makeCentredFitSnGradScheme(SS, POLYNOMIAL, STENCIL) \
184 \
185 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,scalar) \
186 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,vector) \
187 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,sphericalTensor) \
188 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,symmTensor) \
189 makeCentredFitSnGradTypeScheme(SS,POLYNOMIAL,STENCIL,tensor)
190
191
192// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193
194#endif
195
196// ************************************************************************* //
Data for centred fit snGrad schemes.
const List< scalarList > & coeffs() const
Return reference to fit coefficients.
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
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.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
TypeName("CentredFitSnGradScheme")
Runtime type information.
virtual bool corrected() const noexcept
Return true if this scheme uses an explicit correction.
CentredFitSnGradScheme(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
Abstract base class for runtime selected snGrad surface normal gradient schemes.
Definition: snGradScheme.H:79
const fvMesh & mesh() const
Return const reference to mesh.
Definition: snGradScheme.H:139
virtual const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
Namespace for OpenFOAM.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
const direction noexcept
Definition: Scalar.H:223
labelList fv(nPoints)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73