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-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::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 public:
84 
85  //- RunTime type information
86  TypeName("faceLimited");
87 
88 
89  // Constructors
90 
91  //- Construct from mesh and schemeData
92  faceLimitedGrad(const fvMesh& mesh, Istream& schemeData)
93  :
94  gradScheme<Type>(mesh),
95  basicGradScheme_(fv::gradScheme<Type>::New(mesh, schemeData)),
96  k_(readScalar(schemeData))
97  {
98  if (k_ < 0 || k_ > 1)
99  {
101  (
102  schemeData
103  ) << "coefficient = " << k_
104  << " should be >= 0 and <= 1"
105  << exit(FatalIOError);
106  }
107  }
108 
109  //- Disallow default bitwise copy construction
110  faceLimitedGrad(const faceLimitedGrad&) = delete;
111 
112 
113  // Member Functions
114 
115  //- Return the gradient of the given field to the gradScheme::grad
116  // for optional caching
118  calcGrad
119  (
120  const VolField<Type>& vsf,
121  const word& name
122  ) const
123  {
124  return grad(vsf);
125  }
126 
127 
128  // Member Operators
129 
130  //- Disallow default bitwise assignment
131  void operator=(const faceLimitedGrad&) = delete;
132 };
133 
134 
135 // * * * * * * * * * * * * Inline Member Function * * * * * * * * * * * * * //
136 
137 template<class Type>
139 (
140  scalar& limiter,
141  const scalar maxDelta,
142  const scalar minDelta,
143  const scalar extrapolate
144 ) const
145 {
146  if (extrapolate > maxDelta + vSmall)
147  {
148  limiter = min(limiter, maxDelta/extrapolate);
149  }
150  else if (extrapolate < minDelta - vSmall)
151  {
152  limiter = min(limiter, minDelta/extrapolate);
153  }
154 }
155 
156 
157 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
158 
159 template<>
161 (
162  const volScalarField& vsf,
163  const word& name
164 ) const;
165 
166 
167 template<>
169 (
170  const volVectorField& vsf,
171  const word& name
172 ) const;
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace fv
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace Foam
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 #endif
186 
187 // ************************************************************************* //
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
faceLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
TypeName("faceLimited")
RunTime type information.
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.
faceLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
void operator=(const faceLimitedGrad &)=delete
Disallow default bitwise assignment.
Abstract base class for gradient schemes.
Definition: gradScheme.H:63
const fvMesh & mesh() const
Return mesh reference.
Definition: gradScheme.H:122
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const VolField< Type > &, const word &name) const
Calculate and return the grad of the given field.
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
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
labelList fv(nPoints)