skewCorrectedSnGrad.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) 2019 Zeljko Tukovic, FSB Zagreb.
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::skewCorrectedSnGrad
29
30Group
31 grpFvSnGradSchemes
32
33Description
34 Surface gradient scheme with skewness and full
35 explicit non-orthogonal corrections.
36
37Usage
38 Minimal example by using \c system/fvSchemes:
39 \verbatim
40 snGradSchemes
41 {
42 snGrad(<term>) skewCorrected;
43 }
44 \endverbatim
45
46SourceFiles
47 skewCorrectedSnGrad.C
48
49\*---------------------------------------------------------------------------*/
50
51#ifndef skewCorrectedSnGrad_H
52#define skewCorrectedSnGrad_H
53
54#include "snGradScheme.H"
55
56// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58namespace Foam
59{
60
61// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62
63namespace fv
64{
65
66/*---------------------------------------------------------------------------*\
67 Class skewCorrectedSnGrad Declaration
68\*---------------------------------------------------------------------------*/
69
70template<class Type>
72:
73 public snGradScheme<Type>
74{
75 // Private Member Functions
76
77 //- No copy assignment
78 void operator=(const skewCorrectedSnGrad&) = delete;
79
80
81public:
82
83 //- Runtime type information
84 TypeName("skewCorrected");
85
86
87 // Constructors
88
89 //- Construct from mesh
91 :
92 snGradScheme<Type>(mesh)
93 {}
94
95 //- Construct from mesh and data stream
97 :
98 snGradScheme<Type>(mesh)
99 {}
100
101
102 //- Destructor
103 virtual ~skewCorrectedSnGrad() = default;
104
105
106 // Member Functions
107
108 //- Return the interpolation weighting factors for the given field
110 (
112 ) const
113 {
114 return this->mesh().nonOrthDeltaCoeffs();
115 }
116
117 //- Return true if this scheme uses an explicit correction
118 virtual bool corrected() const noexcept
119 {
120 return true;
121 }
122
123 //- Return the explicit correction to the skewCorrectedSnGrad
124 //- for the given field using the gradient of the field
127 (
129 ) const;
130
131 //- Return the explicit correction to the skewCorrectedSnGrad
132 //- for the given field using the gradients of the field components
135};
136
137
138// * * * * * * * * Template Member Function Specialisations * * * * * * * * //
139
140template<>
142(
143 const volScalarField& vsf
144) const;
145
146
147template<>
149(
150 const volVectorField& vvf
151) const;
152
153
154// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155
156} // End namespace fv
157
158// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159
160} // End namespace Foam
161
162// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163
164#ifdef NoRepository
165 #include "skewCorrectedSnGrad.C"
166#endif
167
168// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169
170#endif
171
172// ************************************************************************* //
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
Surface gradient scheme with skewness and full explicit non-orthogonal corrections.
virtual ~skewCorrectedSnGrad()=default
Destructor.
skewCorrectedSnGrad(const fvMesh &mesh, Istream &)
Construct from mesh and data stream.
skewCorrectedSnGrad(const fvMesh &mesh)
Construct from mesh.
TypeName("skewCorrected")
Runtime type information.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
virtual bool corrected() const noexcept
Return true if this scheme uses an explicit correction.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fullGradCorrection(const GeometricField< Type, fvPatchField, volMesh > &) const
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
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