cellLimitedGrad.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::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 minumum 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 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace fv
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class cellLimitedGrad Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class Type>
60 class cellLimitedGrad
61 :
62  public fv::gradScheme<Type>
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  //- Disallow default bitwise copy construct
76 
77  //- Disallow default bitwise assignment
78  void operator=(const cellLimitedGrad&);
79 
80 
81 public:
82 
83  //- RunTime type information
84  TypeName("cellLimited");
85 
86 
87  // Constructors
88 
89  //- Construct from mesh and schemeData
90  cellLimitedGrad(const fvMesh& mesh, Istream& schemeData)
91  :
92  gradScheme<Type>(mesh),
93  basicGradScheme_(fv::gradScheme<Type>::New(mesh, schemeData)),
94  k_(readScalar(schemeData))
95  {
96  if (k_ < 0 || k_ > 1)
97  {
99  (
100  schemeData
101  ) << "coefficient = " << k_
102  << " should be >= 0 and <= 1"
103  << exit(FatalIOError);
104  }
105  }
106 
107 
108  // Member Functions
109 
110  static inline void limitFace
111  (
112  Type& limiter,
113  const Type& maxDelta,
114  const Type& minDelta,
115  const Type& extrapolate
116  );
117 
118  //- Return the gradient of the given field to the gradScheme::grad
119  // for optional caching
120  virtual tmp
121  <
124  > calcGrad
125  (
127  const word& name
128  ) const;
129 };
130 
131 
132 // * * * * * * * * * * * * Inline Member Function * * * * * * * * * * * * * //
133 
134 template<>
136 (
137  scalar& limiter,
138  const scalar& maxDelta,
139  const scalar& minDelta,
140  const scalar& extrapolate
141 )
142 {
143  if (extrapolate > maxDelta + VSMALL)
144  {
145  limiter = min(limiter, maxDelta/extrapolate);
146  }
147  else if (extrapolate < minDelta - VSMALL)
148  {
149  limiter = min(limiter, minDelta/extrapolate);
150  }
151 }
152 
153 
154 template<class Type>
156 (
157  Type& limiter,
158  const Type& maxDelta,
159  const Type& minDelta,
160  const Type& extrapolate
161 )
162 {
163  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
164  {
166  (
167  limiter.component(cmpt),
168  maxDelta.component(cmpt),
169  minDelta.component(cmpt),
170  extrapolate.component(cmpt)
171  );
172  }
173 }
174 
175 
176 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
177 
178 template<>
180 (
181  const volScalarField& vsf,
182  const word& name
183 ) const;
184 
185 
186 template<>
188 (
189  const volVectorField& vsf,
190  const word& name
191 ) const;
192 
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace fv
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 } // End namespace Foam
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 #endif
205 
206 // ************************************************************************* //
TypeName("cellLimited")
RunTime type information.
uint8_t direction
Definition: direction.H:46
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 scalar psiMax, const scalar psiMin)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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.
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:65
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
const fvMesh & mesh() const
Return mesh reference.
Definition: gradScheme.H:122
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 succesful.
Definition: doubleScalar.H:63
static void limitFace(Type &limiter, const Type &maxDelta, const Type &minDelta, const Type &extrapolate)
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
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:54
Namespace for OpenFOAM.
IOerror FatalIOError