uncorrectedSnGrad.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 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::fv::uncorrectedSnGrad
29
30Group
31 grpFvSnGradSchemes
32
33Description
34 Surface gradient scheme with no non-orthogonal correction.
35
36Usage
37 Minimal example by using \c system/fvSchemes:
38 \verbatim
39 snGradSchemes
40 {
41 snGrad(<term>) uncorrected;
42 }
43 \endverbatim
44
45Note
46 - Interpolation weighting factors (i.e. delta coefficients) are based
47 on the \c nonOrthDeltaCoeffs function rather than the \c deltaCoeffs
48 function, which is used by the \c orthogonal scheme.
49
50SourceFiles
51 uncorrectedSnGrad.C
52
53\*---------------------------------------------------------------------------*/
54
55#ifndef uncorrectedSnGrad_H
56#define uncorrectedSnGrad_H
57
58#include "snGradScheme.H"
59
60// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61
62namespace Foam
63{
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67namespace fv
68{
69
70/*---------------------------------------------------------------------------*\
71 Class uncorrectedSnGrad Declaration
72\*---------------------------------------------------------------------------*/
73
74template<class Type>
76:
77 public snGradScheme<Type>
78{
79 // Private Member Functions
80
81 //- No copy assignment
82 void operator=(const uncorrectedSnGrad&) = delete;
83
84
85public:
86
87 //- Runtime type information
88 TypeName("uncorrected");
89
90
91 // Constructors
92
93 //- Construct from mesh
95 :
96 snGradScheme<Type>(mesh)
97 {}
98
99 //- Construct from mesh and data stream
101 :
102 snGradScheme<Type>(mesh)
103 {}
104
105
106 //- Destructor
107 virtual ~uncorrectedSnGrad() = default;
108
109
110 // Member Functions
111
112 //- Return the interpolation weighting factors for the given field
114 (
116 ) const
117 {
118 return this->mesh().nonOrthDeltaCoeffs();
119 }
120
121 //- Return true if this scheme uses an explicit correction
122 virtual bool corrected() const noexcept
123 {
124 return false;
125 }
126
127 //- Return the explicit correction to the uncorrectedSnGrad
128 //- for the given field
131};
132
133
134// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135
136} // End namespace fv
137
138// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140} // End namespace Foam
141
142// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143
144#ifdef NoRepository
145 #include "uncorrectedSnGrad.C"
146#endif
147
148// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149
150#endif
151
152// ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
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
Surface gradient scheme with no non-orthogonal correction.
virtual ~uncorrectedSnGrad()=default
Destructor.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
uncorrectedSnGrad(const fvMesh &mesh)
Construct from mesh.
TypeName("uncorrected")
Runtime type information.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
virtual bool corrected() const noexcept
Return true if this scheme uses an explicit correction.
uncorrectedSnGrad(const fvMesh &mesh, Istream &)
Construct from mesh and data stream.
virtual const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
Definition: tmp.H:65
Namespace for OpenFOAM.
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