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-2022 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 
87 public:
88 
89  //- RunTime type information
90  TypeName("cellLimited");
91 
92 
93  // Constructors
94 
95  //- Construct from mesh and schemeData
96  cellLimitedGrad(const fvMesh& mesh, Istream& schemeData)
97  :
98  gradScheme<Type>(mesh),
99  Limiter(schemeData),
100  basicGradScheme_(fv::gradScheme<Type>::New(mesh, schemeData)),
101  k_(readScalar(schemeData))
102  {
103  if (k_ < 0 || k_ > 1)
104  {
106  (
107  schemeData
108  ) << "coefficient = " << k_
109  << " should be >= 0 and <= 1"
110  << exit(FatalIOError);
111  }
112  }
113 
114  //- Disallow default bitwise copy construction
115  cellLimitedGrad(const cellLimitedGrad&) = delete;
116 
117 
118  // Member Functions
119 
120  inline void limitFaceCmpt
121  (
122  scalar& limiter,
123  const scalar maxDelta,
124  const scalar minDelta,
125  const scalar extrapolate
126  ) const;
127 
128  inline void limitFace
129  (
130  Type& limiter,
131  const Type& maxDelta,
132  const Type& minDelta,
133  const Type& extrapolate
134  ) const;
135 
136  //- Return the gradient of the given field to the gradScheme::grad
137  // for optional caching
139  calcGrad
140  (
141  const VolField<Type>& vsf,
142  const word& name
143  ) const;
144 
145 
146  // Member Operators
147 
148  //- Disallow default bitwise assignment
149  void operator=(const cellLimitedGrad&) = delete;
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 
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  {
194  limitFaceCmpt
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 // ************************************************************************* //
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:101
cellLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
virtual tmp< VolField< typename outerProduct< vector, Type >::type > > calcGrad(const VolField< Type > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad.
void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate) const
void limitFaceCmpt(scalar &limiter, const scalar maxDelta, const scalar minDelta, const scalar extrapolate) const
TypeName("cellLimited")
RunTime type information.
void operator=(const cellLimitedGrad &)=delete
Disallow default bitwise assignment.
cellLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
Abstract base class for gradient schemes.
Definition: gradScheme.H:63
const fvMesh & mesh() const
Return mesh reference.
Definition: gradScheme.H:122
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: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
void limiter(surfaceScalarField &lambda, 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)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
label & setComponent(label &l, const direction)
Definition: label.H:86
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
IOerror FatalIOError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
uint8_t direction
Definition: direction.H:45
labelList fv(nPoints)