Loading...
Searching...
No Matches
correctedSnGrad.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::correctedSnGrad
29
30Group
31 grpFvSnGradSchemes
32
33Description
34 Surface gradient scheme with full explicit non-orthogonal correction.
35
36Usage
37 Minimal example by using \c system/fvSchemes:
38 \verbatim
39 snGradSchemes
40 {
41 snGrad(<term>) corrected;
42 }
43 \endverbatim
44
45SourceFiles
46 correctedSnGrad.C
47
48\*---------------------------------------------------------------------------*/
49
50#ifndef correctedSnGrad_H
51#define correctedSnGrad_H
52
53#include "snGradScheme.H"
54
55// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57namespace Foam
58{
59
60// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61
62namespace fv
63{
64
65/*---------------------------------------------------------------------------*\
66 Class correctedSnGrad Declaration
67\*---------------------------------------------------------------------------*/
68
69template<class Type>
71:
72 public snGradScheme<Type>
73{
74 // Private Member Functions
75
76 //- No copy assignment
77 void operator=(const correctedSnGrad&) = delete;
78
79
80public:
81
82 //- Runtime type information
83 TypeName("corrected");
84
85
86 // Constructors
87
88 //- Construct from mesh
90 :
91 snGradScheme<Type>(mesh)
92 {}
93
94 //- Construct from mesh and data stream
96 :
97 snGradScheme<Type>(mesh)
98 {}
99
100
101 //- Destructor
102 virtual ~correctedSnGrad() = default;
103
104
105 // Member Functions
106
107 //- Return the interpolation weighting factors for the given field
109 (
111 ) const
112 {
113 return this->mesh().nonOrthDeltaCoeffs();
114 }
116 //- Return true if this scheme uses an explicit correction
117 virtual bool corrected() const noexcept
118 {
119 return true;
120 }
121
122 //- Return the explicit correction to the correctedSnGrad
123 //- for the given field using the gradient of the field
128 ) const;
129
130 //- Return the explicit correction to the correctedSnGrad
131 //- for the given field using the gradients of the field components
134};
135
136
137// * * * * * * * * Template Member Function Specialisations * * * * * * * * //
138
139template<>
141(
142 const volScalarField& vsf
143) const;
144
145
146template<>
148(
149 const volVectorField& vvf
150) const;
151
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155} // End namespace fv
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159} // End namespace Foam
161// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
162
163#ifdef NoRepository
164 #include "correctedSnGrad.C"
165#endif
166
167// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168
169#endif
170
171// ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
TypeName("corrected")
Runtime type information.
correctedSnGrad(const fvMesh &mesh, Istream &)
Construct from mesh and data stream.
correctedSnGrad(const fvMesh &mesh)
Construct from mesh.
virtual ~correctedSnGrad()=default
Destructor.
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
Return the explicit correction to the correctedSnGrad for the given field using the gradients of the ...
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > fullGradCorrection(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the correctedSnGrad for the given field using the gradient of the f...
Abstract base class for runtime selected snGrad surface normal gradient schemes.
const fvMesh & mesh() const
Return const reference to mesh.
virtual const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
Definition tmp.H:75
dynamicFvMesh & mesh
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const direction noexcept
Definition scalarImpl.H:265
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68