cellMDLimitedGrad.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-2020 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::cellMDLimitedGrad
26 
27 Description
28  cellMDLimitedGrad gradient scheme applied to a runTime selected base
29  gradient 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 the gradient in each face direction separately.
34 
35 SourceFiles
36  cellMDLimitedGrad.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef cellMDLimitedGrad_H
41 #define cellMDLimitedGrad_H
42 
43 #include "gradScheme.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace fv
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class cellMDLimitedGrad Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 template<class Type>
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 public:
73 
74  //- RunTime type information
75  TypeName("cellMDLimited");
76 
77 
78  // Constructors
79 
80  //- Construct from mesh and schemeData
81  cellMDLimitedGrad(const fvMesh& mesh, Istream& schemeData)
82  :
83  gradScheme<Type>(mesh),
84  basicGradScheme_(fv::gradScheme<Type>::New(mesh, schemeData)),
85  k_(readScalar(schemeData))
86  {
87  if (k_ < 0 || k_ > 1)
88  {
90  (
91  schemeData
92  ) << "coefficient = " << k_
93  << " should be >= 0 and <= 1"
94  << exit(FatalIOError);
95  }
96  }
97 
98  //- Disallow default bitwise copy construction
99  cellMDLimitedGrad(const cellMDLimitedGrad&) = delete;
100 
101 
102  // Member Functions
103 
104  static inline void limitFace
105  (
107  const Type& maxDelta,
108  const Type& minDelta,
109  const vector& dcf
110  );
111 
112  //- Return the gradient of the given field to the gradScheme::grad
113  // for optional caching
114  virtual tmp
115  <
118  > calcGrad
119  (
121  const word& name
122  ) const;
123 
124 
125  // Member Operators
126 
127  //- Disallow default bitwise assignment
128  void operator=(const cellMDLimitedGrad&) = delete;
129 };
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 template<>
136 (
137  vector& g,
138  const scalar& maxDelta,
139  const scalar& minDelta,
140  const vector& dcf
141 )
142 {
143  scalar extrapolate = dcf & g;
144 
145  if (extrapolate > maxDelta)
146  {
147  g = g + dcf*(maxDelta - extrapolate)/magSqr(dcf);
148  }
149  else if (extrapolate < minDelta)
150  {
151  g = g + dcf*(minDelta - extrapolate)/magSqr(dcf);
152  }
153 }
154 
155 
156 template<class Type>
158 (
160  const Type& maxDelta,
161  const Type& minDelta,
162  const vector& dcf
163 )
164 {
165  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
166  {
167  vector gi(g[cmpt], g[cmpt+3], g[cmpt+6]);
169  (
170  gi,
171  maxDelta.component(cmpt),
172  minDelta.component(cmpt),
173  dcf
174  );
175  g[cmpt] = gi.x();
176  g[cmpt+3] = gi.y();
177  g[cmpt+6] = gi.z();
178  }
179 }
180 
181 
182 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
183 
184 template<>
186 (
187  const volScalarField& vsf,
188  const word& name
189 ) const;
190 
191 
192 template<>
194 (
195  const volVectorField& vsf,
196  const word& name
197 ) const;
198 
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 } // End namespace fv
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 } // End namespace Foam
207 
208 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 
210 #endif
211 
212 // ************************************************************************* //
cellMDLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
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
void operator=(const cellMDLimitedGrad &)=delete
Disallow default bitwise assignment.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const Cmpt & z() const
Definition: VectorI.H:87
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
Generic GeometricField class.
const Cmpt & y() const
Definition: VectorI.H:81
Mesh data needed to do the Finite Volume discretisation.
Definition: volMesh.H:53
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 handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static void limitFace(typename outerProduct< vector, Type >::type &g, const Type &maxDelta, const Type &minDelta, const vector &dcf)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:75
const Cmpt & x() const
Definition: VectorI.H:75
cellMDLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
TypeName("cellMDLimited")
RunTime type information.
dimensioned< scalar > magSqr(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:335
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
const Cmpt & component(const direction) const
Definition: VectorSpaceI.H:91
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
const fvMesh & mesh() const
Return mesh reference.
Definition: gradScheme.H:122
Namespace for OpenFOAM.
IOerror FatalIOError