cellLimitedGrad.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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::cellLimitedGrad
26 
27 Description
28  cellLimitedGrad gradient scheme applied to a runTime selected base gradient
29  scheme.
30 
31  The scalar limiter based on limiting the extrapolated face values
32  between the maximum and minimum cell and cell neighbour values and is
33  applied to all components of the gradient.
34 
35 SourceFiles
36  cellLimitedGrad.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef cellLimitedGrad_H
41 #define cellLimitedGrad_H
42 
43 #include "gradScheme.H"
44 #include "Field.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 namespace fv
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class cellLimitedGrad Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 template<class Type, class Limiter>
59 class cellLimitedGrad
60 :
61  public fv::gradScheme<Type>,
62  public Limiter
63 {
64  // Private Data
65 
66  tmp<fv::gradScheme<Type>> basicGradScheme_;
67 
68  //- Limiter coefficient
69  const scalar k_;
70 
71 
72  // Private Member Functions
73 
74  void limitGradient
75  (
76  const Field<scalar>& limiter,
77  Field<vector>& gIf
78  ) const;
79 
80  void limitGradient
81  (
82  const Field<vector>& limiter,
83  Field<tensor>& gIf
84  ) const;
85 
86  //- Disallow default bitwise copy construct
88 
89  //- Disallow default bitwise assignment
90  void operator=(const cellLimitedGrad&);
91 
92 
93 public:
94 
95  //- RunTime type information
96  TypeName("cellLimited");
97 
98 
99  // Constructors
100 
101  //- Construct from mesh and schemeData
102  cellLimitedGrad(const fvMesh& mesh, Istream& schemeData)
103  :
104  gradScheme<Type>(mesh),
105  Limiter(schemeData),
106  basicGradScheme_(fv::gradScheme<Type>::New(mesh, schemeData)),
107  k_(readScalar(schemeData))
108  {
109  if (k_ < 0 || k_ > 1)
110  {
112  (
113  schemeData
114  ) << "coefficient = " << k_
115  << " should be >= 0 and <= 1"
116  << exit(FatalIOError);
117  }
118  }
119 
120 
121  // Member Functions
122 
123  inline void limitFaceCmpt
124  (
125  scalar& limiter,
126  const scalar maxDelta,
127  const scalar minDelta,
128  const scalar extrapolate
129  ) const;
130 
131  inline void limitFace
132  (
133  Type& limiter,
134  const Type& maxDelta,
135  const Type& minDelta,
136  const Type& extrapolate
137  ) const;
138 
139  //- Return the gradient of the given field to the gradScheme::grad
140  // for optional caching
141  virtual tmp
142  <
145  > calcGrad
146  (
148  const word& name
149  ) const;
150 };
151 
152 
153 // * * * * * * * * * * * * Inline Member Function * * * * * * * * * * * * * //
154 
155 template<class Type, class Limiter>
157 (
158  scalar& limiter,
159  const scalar maxDelta,
160  const scalar minDelta,
161  const scalar extrapolate
162 ) const
163 {
164  scalar r = 1;
165 
166  if (extrapolate > small)
167  {
168  r = maxDelta/extrapolate;
169  }
170  else if (extrapolate < -small)
171  {
172  r = minDelta/extrapolate;
173  }
174  else
175  {
176  return;
177  }
178 
179  limiter = min(limiter, Limiter::limiter(r));
180 }
181 
182 
183 template<class Type, class Limiter>
185 (
186  Type& limiter,
187  const Type& maxDelta,
188  const Type& minDelta,
189  const Type& extrapolate
190 ) const
191 {
192  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
193  {
195  (
196  setComponent(limiter, cmpt),
197  component(maxDelta, cmpt),
198  component(minDelta, cmpt),
199  component(extrapolate, cmpt)
200  );
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 } // End namespace fv
208 
209 } // End namespace Foam
210 
211 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 
213 #ifdef NoRepository
214  #include "cellLimitedGrad.C"
215 #endif
216 
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 
219 #endif
220 
221 // ************************************************************************* //
TypeName("cellLimited")
RunTime type information.
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
uint8_t direction
Definition: direction.H:45
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
Generic GeometricField class.
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Abstract base class for gradient schemes.
Definition: gradScheme.H:60
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
void limitFaceCmpt(scalar &limiter, const scalar maxDelta, const scalar minDelta, const scalar extrapolate) const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
static tmp< gradScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return a pointer to a new gradScheme created on freestore.
Definition: gradScheme.C:34
A class for managing temporary objects.
Definition: PtrList.H:53
const fvMesh & mesh() const
Return mesh reference.
Definition: gradScheme.H:122
label & setComponent(label &l, const direction)
Definition: label.H:79
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.
void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate) const
IOerror FatalIOError