limitedSnGrad.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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::fv::limitedSnGrad
29 
30 Group
31  grpFvSnGradSchemes
32 
33 Description
34  Surface gradient scheme with limited explicit non-orthogonal correction.
35 
36  The limiter is controlled by a coefficient with a value between 0 and 1
37  which when 0 switches the correction off and the scheme behaves as
38  \c uncorrected \c snGrad, when set to 1 the full correction of the
39  selected scheme is used and the scheme behaves as \c corrected \c snGrad,
40  and when set to 0.5 the limiter is calculated such that the non-orthogonal
41  component does not exceed the orthogonal component.
42 
43 Usage
44  Minimal example by using \c system/fvSchemes:
45  \verbatim
46  snGradSchemes
47  {
48  snGrad(<term>) limited <corrected scheme> <coefficient>;
49 
50  // Backward compatibility
51  snGrad(<term>) limited <coefficient>;
52  }
53  \endverbatim
54 
55 SourceFiles
56  limitedSnGrad.C
57 
58 \*---------------------------------------------------------------------------*/
59 
60 #ifndef limitedSnGrad_H
61 #define limitedSnGrad_H
62 
63 #include "correctedSnGrad.H"
64 
65 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67 namespace Foam
68 {
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace fv
73 {
74 
75 /*---------------------------------------------------------------------------*\
76  Class limitedSnGrad Declaration
77 \*---------------------------------------------------------------------------*/
78 
79 template<class Type>
80 class limitedSnGrad
81 :
82  public snGradScheme<Type>
83 {
84  // Private Data
85 
86  //- Type of correction scheme
87  tmp<snGradScheme<Type>> correctedScheme_;
88 
89  //- Limiter coefficient
90  scalar limitCoeff_;
91 
92 
93  // Private Member Functions
94 
95  //- No copy assignment
96  void operator=(const limitedSnGrad&) = delete;
97 
98  //- Lookup function for the corrected to support backward compatibility
99  //- of dictionary specification
100  tmp<snGradScheme<Type>> lookupCorrectedScheme(Istream& schemeData)
101  {
102  token nextToken(schemeData);
103 
104  if (nextToken.isNumber())
105  {
106  limitCoeff_ = nextToken.number();
107  return tmp<snGradScheme<Type>>
108  (
109  new correctedSnGrad<Type>(this->mesh())
110  );
111  }
112  else
113  {
114  schemeData.putBack(nextToken);
115  tmp<snGradScheme<Type>> tcorrectedScheme
116  (
117  fv::snGradScheme<Type>::New(this->mesh(), schemeData)
118  );
119 
120  schemeData >> limitCoeff_;
121 
122  return tcorrectedScheme;
123  }
124  }
125 
126 
127 public:
128 
129  //- Runtime type information
130  TypeName("limited");
131 
132 
133  // Constructors
134 
135  //- Construct from mesh
136  limitedSnGrad(const fvMesh& mesh)
137  :
138  snGradScheme<Type>(mesh),
139  correctedScheme_(new correctedSnGrad<Type>(this->mesh())),
140  limitCoeff_(1)
141  {}
142 
143  //- Construct from mesh and data stream
144  limitedSnGrad(const fvMesh& mesh, Istream& schemeData)
145  :
146  snGradScheme<Type>(mesh),
147  correctedScheme_(lookupCorrectedScheme(schemeData))
148  {
149  if (limitCoeff_ < 0 || limitCoeff_ > 1)
150  {
151  FatalIOErrorInFunction(schemeData)
152  << "limitCoeff is specified as " << limitCoeff_
153  << " but should be >= 0 && <= 1"
154  << exit(FatalIOError);
155  }
156  }
157 
158 
159  //- Destructor
160  virtual ~limitedSnGrad() = default;
161 
162 
163  // Member Functions
164 
165  //- Return the interpolation weighting factors for the given field
167  (
169  ) const
170  {
171  return this->mesh().nonOrthDeltaCoeffs();
172  }
173 
174  //- Return true if this scheme uses an explicit correction
175  virtual bool corrected() const noexcept
176  {
177  return true;
178  }
179 
180  //- Return the explicit correction to the limitedSnGrad
181  //- for the given field
184 };
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace fv
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace Foam
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 #ifdef NoRepository
198  #include "limitedSnGrad.C"
199 #endif
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 #endif
204 
205 // ************************************************************************* //
Foam::surfaceInterpolation::nonOrthDeltaCoeffs
virtual const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
Definition: surfaceInterpolation.C:125
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::fv::limitedSnGrad::corrected
virtual bool corrected() const noexcept
Return true if this scheme uses an explicit correction.
Definition: limitedSnGrad.H:174
Foam::FatalIOError
IOerror FatalIOError
Foam::token
A token holds an item read from Istream.
Definition: token.H:68
Foam::fv::limitedSnGrad::TypeName
TypeName("limited")
Runtime type information.
Foam::fv::limitedSnGrad::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Definition: limitedSnGrad.C:51
Foam::token::isNumber
bool isNumber() const noexcept
Token is LABEL, FLOAT or DOUBLE.
Definition: tokenI.H:587
Foam::fv::limitedSnGrad::~limitedSnGrad
virtual ~limitedSnGrad()=default
Destructor.
limitedSnGrad.C
Foam::token::number
scalar number() const
Return label, float or double value.
Definition: tokenI.H:593
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
Foam::fv::limitedSnGrad::limitedSnGrad
limitedSnGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and data stream.
Definition: limitedSnGrad.H:143
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
correctedSnGrad.H
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
fv
labelList fv(nPoints)
Foam::fv::snGradScheme::mesh
const fvMesh & mesh() const
Return const reference to mesh.
Definition: snGradScheme.H:139
Foam::fv::limitedSnGrad
Surface gradient scheme with limited explicit non-orthogonal correction.
Definition: limitedSnGrad.H:79
Foam::Istream::putBack
void putBack(const token &tok)
Put back a token. Only a single put back is permitted.
Definition: Istream.C:70
Foam::fv::correctedSnGrad
Surface gradient scheme with full explicit non-orthogonal correction.
Definition: correctedSnGrad.H:69
Foam::fv::snGradScheme
Abstract base class for runtime selected snGrad surface normal gradient schemes.
Definition: snGradScheme.H:76
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::fv::limitedSnGrad::limitedSnGrad
limitedSnGrad(const fvMesh &mesh)
Construct from mesh.
Definition: limitedSnGrad.H:135
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::fv::limitedSnGrad::deltaCoeffs
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
Definition: limitedSnGrad.H:166