faceLimitedGrad.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::faceLimitedGrad
26 
27 Description
28  faceLimitedGrad 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 face-neighbour cell values and is applied to all components
33  of the gradient.
34 
35 SourceFiles
36  faceLimitedGrad.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef faceLimitedGrad_H
41 #define faceLimitedGrad_H
42 
43 #include "gradScheme.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace fv
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class faceLimitedGrad Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class Type>
60 class faceLimitedGrad
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  inline void limitFace
75  (
76  scalar& limiter,
77  const scalar maxDelta,
78  const scalar minDelta,
79  const scalar extrapolate
80  ) const;
81 
82 
83  //- Disallow default bitwise copy construct
85 
86  //- Disallow default bitwise assignment
87  void operator=(const faceLimitedGrad&);
88 
89 
90 public:
91 
92  //- RunTime type information
93  TypeName("faceLimited");
94 
95 
96  // Constructors
97 
98  //- Construct from mesh and schemeData
99  faceLimitedGrad(const fvMesh& mesh, Istream& schemeData)
100  :
101  gradScheme<Type>(mesh),
102  basicGradScheme_(fv::gradScheme<Type>::New(mesh, schemeData)),
103  k_(readScalar(schemeData))
104  {
105  if (k_ < 0 || k_ > 1)
106  {
108  (
109  schemeData
110  ) << "coefficient = " << k_
111  << " should be >= 0 and <= 1"
112  << exit(FatalIOError);
113  }
114  }
115 
116 
117  // Member Functions
118 
119  //- Return the gradient of the given field to the gradScheme::grad
120  // for optional caching
121  virtual tmp
122  <
126  (
128  const word& name
129  ) const
130  {
131  return grad(vsf);
132  }
133 };
134 
135 
136 // * * * * * * * * * * * * Inline Member Function * * * * * * * * * * * * * //
137 
138 template<class Type>
140 (
141  scalar& limiter,
142  const scalar maxDelta,
143  const scalar minDelta,
144  const scalar extrapolate
145 ) const
146 {
147  if (extrapolate > maxDelta + vSmall)
148  {
149  limiter = min(limiter, maxDelta/extrapolate);
150  }
151  else if (extrapolate < minDelta - vSmall)
152  {
153  limiter = min(limiter, minDelta/extrapolate);
154  }
155 }
156 
157 
158 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
159 
160 template<>
162 (
163  const volScalarField& vsf,
164  const word& name
165 ) const;
166 
167 
168 template<>
170 (
171  const volVectorField& vsf,
172  const word& name
173 ) const;
174 
175 
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177 
178 } // End namespace fv
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #endif
187 
188 // ************************************************************************* //
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
TypeName("faceLimited")
RunTime type information.
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
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.
A class for managing temporary objects.
Definition: PtrList.H:53
faceLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
const fvMesh & mesh() const
Return mesh reference.
Definition: gradScheme.H:122
Namespace for OpenFOAM.
IOerror FatalIOError
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvPatchField, volMesh > &, const word &name) const
Calculate and return the grad of the given field.