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-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::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
115  calcGrad
116  (
117  const VolField<Type>& vsf,
118  const word& name
119  ) const;
120 
121 
122  // Member Operators
123 
124  //- Disallow default bitwise assignment
125  void operator=(const cellMDLimitedGrad&) = delete;
126 };
127 
128 
129 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
130 
131 template<>
133 (
134  vector& g,
135  const scalar& maxDelta,
136  const scalar& minDelta,
137  const vector& dcf
138 )
139 {
140  scalar extrapolate = dcf & g;
141 
142  if (extrapolate > maxDelta)
143  {
144  g = g + dcf*(maxDelta - extrapolate)/magSqr(dcf);
145  }
146  else if (extrapolate < minDelta)
147  {
148  g = g + dcf*(minDelta - extrapolate)/magSqr(dcf);
149  }
150 }
151 
152 
153 template<class Type>
155 (
157  const Type& maxDelta,
158  const Type& minDelta,
159  const vector& dcf
160 )
161 {
162  for (direction cmpt=0; cmpt<Type::nComponents; cmpt++)
163  {
164  vector gi(g[cmpt], g[cmpt+3], g[cmpt+6]);
166  (
167  gi,
168  maxDelta.component(cmpt),
169  minDelta.component(cmpt),
170  dcf
171  );
172  g[cmpt] = gi.x();
173  g[cmpt+3] = gi.y();
174  g[cmpt+6] = gi.z();
175  }
176 }
177 
178 
179 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
180 
181 template<>
183 (
184  const volScalarField& vsf,
185  const word& name
186 ) const;
187 
188 
189 template<>
191 (
192  const volVectorField& vsf,
193  const word& name
194 ) const;
195 
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 } // End namespace fv
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 } // End namespace Foam
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 #endif
208 
209 // ************************************************************************* //
Generic GeometricField class.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
const Cmpt & component(const direction) const
Definition: VectorSpaceI.H:91
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
cellMDLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
static void limitFace(typename outerProduct< vector, Type >::type &g, const Type &maxDelta, const Type &minDelta, const vector &dcf)
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.
cellMDLimitedGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and schemeData.
void operator=(const cellMDLimitedGrad &)=delete
Disallow default bitwise assignment.
TypeName("cellMDLimited")
RunTime type information.
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
typeOfRank< typename pTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank) >::type type
Definition: products.H:90
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
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
dimensioned< scalar > magSqr(const dimensioned< Type > &)
uint8_t direction
Definition: direction.H:45
labelList fv(nPoints)