relaxedSnGrad.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) 2021 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::fv::relaxedSnGrad
28
29Group
30 grpFvSnGradSchemes
31
32Description
33 Surface gradient scheme with under-/over-relaxed
34 full or limited explicit non-orthogonal correction.
35
36Usage
37 Minimal example by using \c system/fvSchemes:
38 \verbatim
39 snGradSchemes
40 {
41 snGrad(<term>) relaxed;
42 }
43 \endverbatim
44
45 and by using \c system/fvSolution:
46 \verbatim
47 relaxationFactors
48 {
49 fields
50 {
51 snGrad(<term>) <relaxation factor>;
52 }
53 }
54 \endverbatim
55
56
57SourceFiles
58 relaxedSnGrad.C
59
60\*---------------------------------------------------------------------------*/
61
62#ifndef relaxedSnGrad_H
63#define relaxedSnGrad_H
64
65#include "correctedSnGrad.H"
66
67// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69namespace Foam
70{
71
72// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73
74namespace fv
75{
76
77/*---------------------------------------------------------------------------*\
78 Class relaxedSnGrad Declaration
79\*---------------------------------------------------------------------------*/
80
81template<class Type>
82class relaxedSnGrad
83:
84 public snGradScheme<Type>
85{
86 // Private Data
87
88 //- Type of correction scheme
89 tmp<snGradScheme<Type>> correctedScheme_;
90
91
92 // Private Member Functions
93
94 //- No copy assignment
95 void operator=(const relaxedSnGrad&) = delete;
96
97
98public:
99
100 //- Runtime type information
101 TypeName("relaxed");
102
103
104 // Constructors
105
106 //- Construct from mesh
108 :
109 snGradScheme<Type>(mesh),
110 correctedScheme_(new correctedSnGrad<Type>(this->mesh()))
111 {}
112
113 //- Construct from mesh and data stream
114 relaxedSnGrad(const fvMesh& mesh, Istream& schemeData)
115 :
116 snGradScheme<Type>(mesh),
117 correctedScheme_(new correctedSnGrad<Type>(this->mesh()))
118 {}
119
120
121 //- Destructor
122 virtual ~relaxedSnGrad() = default;
123
124
125 // Member Functions
126
127 //- Return the interpolation weighting factors for the given field
129 (
131 ) const
132 {
133 return this->mesh().nonOrthDeltaCoeffs();
134 }
135
136 //- Return true if this scheme uses an explicit correction
137 virtual bool corrected() const noexcept
138 {
139 return true;
140 }
141
142 //- Return the explicit correction to the relaxedSnGrad
143 //- for the given field using the gradients of the field components
146};
147
148// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149
150} // End namespace fv
151
152// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153
154} // End namespace Foam
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158#ifdef NoRepository
159 #include "relaxedSnGrad.C"
160#endif
161
162// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163
164#endif
165
166// ************************************************************************* //
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 full explicit non-orthogonal correction.
Surface gradient scheme with under-/over-relaxed full or limited explicit non-orthogonal correction.
Definition: relaxedSnGrad.H:84
relaxedSnGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and data stream.
relaxedSnGrad(const fvMesh &mesh)
Construct from mesh.
virtual ~relaxedSnGrad()=default
Destructor.
TypeName("relaxed")
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
Definition: relaxedSnGrad.C:38
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