faceMDLimitedGrad.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::faceMDLimitedGrad
26 
27 Description
28  faceMDLimitedGrad gradient scheme applied to a runTime selected
29  base gradient scheme.
30 
31  The scalar limiter based on limiting the extrapolated face values
32  between the face-neighbour cell values and is applied to the gradient
33  in each face direction separately.
34 
35 SourceFiles
36  faceMDLimitedGrad.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef faceMDLimitedGrad_H
41 #define faceMDLimitedGrad_H
42 
43 #include "gradScheme.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace fv
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class faceMDLimitedGrad 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  // 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 construction
84  faceMDLimitedGrad(const faceMDLimitedGrad&) = delete;
85 
86  //- Disallow default bitwise assignment
87  void operator=(const faceMDLimitedGrad&) = delete;
88 
89 
90 public:
91 
92  //- RunTime type information
93  TypeName("faceMDLimited");
94 
95 
96  // Constructors
97 
98  //- Construct from mesh and schemeData
99  faceMDLimitedGrad(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
122  calcGrad
123  (
124  const VolField<Type>& vsf,
125  const word& name
126  ) const;
127 };
128 
129 
130 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
131 
132 template<>
134 (
135  const volScalarField& vsf,
136  const word& name
137 ) const;
138 
139 
140 template<>
142 (
143  const volVectorField& vsf,
144  const word& name
145 ) const;
146 
147 
148 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149 
150 } // End namespace fv
151 
152 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 
154 } // End namespace Foam
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 #endif
159 
160 // ************************************************************************* //
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
faceMDLimitedGrad gradient scheme applied to a runTime selected base gradient scheme.
TypeName("faceMDLimited")
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.
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
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
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)