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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::fv::limitedSnGrad
26 
27 Description
28  Run-time selected snGrad scheme with limited non-orthogonal correction.
29 
30  The limiter is controlled by a coefficient with a value between 0 and 1
31  which when 0 switches the correction off and the scheme behaves as
32  uncorrectedSnGrad, when set to 1 the full correction of the selected scheme
33  is used and when set to 0.5 the limiter is calculated such that the
34  non-orthogonal contribution does not exceed the orthogonal part.
35 
36  Format:
37  limited <corrected scheme> <coefficient>;
38 
39  or
40 
41  limited <coefficient>; // Backward compatibility
42 
43 SourceFiles
44  limitedSnGrad.C
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #ifndef limitedSnGrad_H
49 #define limitedSnGrad_H
50 
51 #include "correctedSnGrad.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace fv
61 {
62 
63 /*---------------------------------------------------------------------------*\
64  Class limitedSnGrad Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 template<class Type>
68 class limitedSnGrad
69 :
70  public snGradScheme<Type>
71 {
72  // Private data
73 
74  tmp<snGradScheme<Type>> correctedScheme_;
75 
76  scalar limitCoeff_;
77 
78 
79  // Private Member Functions
80 
81  //- Disallow default bitwise assignment
82  void operator=(const limitedSnGrad&);
83 
84  //- Lookup function for the corrected to support backward compatibility
85  // of dictionary specification
86  tmp<snGradScheme<Type>> lookupCorrectedScheme(Istream& schemeData)
87  {
88  token nextToken(schemeData);
89 
90  if (nextToken.isNumber())
91  {
92  limitCoeff_ = nextToken.number();
94  (
95  new correctedSnGrad<Type>(this->mesh())
96  );
97  }
98  else
99  {
100  schemeData.putBack(nextToken);
101  tmp<snGradScheme<Type>> tcorrectedScheme
102  (
103  fv::snGradScheme<Type>::New(this->mesh(), schemeData)
104  );
105 
106  schemeData >> limitCoeff_;
107 
108  return tcorrectedScheme;
109  }
110  }
111 
112 
113 public:
114 
115  //- Runtime type information
116  TypeName("limited");
117 
118 
119  // Constructors
120 
121  //- Construct from mesh
122  limitedSnGrad(const fvMesh& mesh)
123  :
124  snGradScheme<Type>(mesh),
125  correctedScheme_(new correctedSnGrad<Type>(this->mesh())),
126  limitCoeff_(1)
127  {}
128 
129 
130  //- Construct from mesh and data stream
131  limitedSnGrad(const fvMesh& mesh, Istream& schemeData)
132  :
133  snGradScheme<Type>(mesh),
134  correctedScheme_(lookupCorrectedScheme(schemeData))
135  {
136  if (limitCoeff_ < 0 || limitCoeff_ > 1)
137  {
139  (
140  schemeData
141  ) << "limitCoeff is specified as " << limitCoeff_
142  << " but should be >= 0 && <= 1"
143  << exit(FatalIOError);
144  }
145  }
146 
147 
148  //- Destructor
149  virtual ~limitedSnGrad();
150 
151 
152  // Member Functions
153 
154  //- Return the interpolation weighting factors for the given field
156  (
158  ) const
159  {
160  return this->mesh().nonOrthDeltaCoeffs();
161  }
162 
163  //- Return true if this scheme uses an explicit correction
164  virtual bool corrected() const
165  {
166  return true;
167  }
168 
169  //- Return the explicit correction to the limitedSnGrad
170  // for the given field
173 };
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace fv
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #ifdef NoRepository
187  #include "limitedSnGrad.C"
188 #endif
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #endif
193 
194 // ************************************************************************* //
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
const fvMesh & mesh() const
Return mesh reference.
Definition: snGradScheme.H:123
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
A token holds items read from Istream.
Definition: token.H:69
void putBack(const token &)
Put back token.
Definition: Istream.C:30
Generic GeometricField class.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the limitedSnGrad.
Definition: limitedSnGrad.C:54
TypeName("limited")
Runtime type information.
limitedSnGrad(const fvMesh &mesh)
Construct from mesh.
Run-time selected snGrad scheme with limited non-orthogonal correction.
Definition: limitedSnGrad.H:67
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
labelList fv(nPoints)
Simple central-difference snGrad scheme with non-orthogonal correction.
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:60
scalar number() const
Definition: tokenI.H:337
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:54
virtual ~limitedSnGrad()
Destructor.
Definition: limitedSnGrad.C:45
bool isNumber() const
Definition: tokenI.H:332
Namespace for OpenFOAM.
IOerror FatalIOError