faceCorrectedSnGrad.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::faceCorrectedSnGrad
26 
27 Description
28  Simple central-difference snGrad scheme with non-orthogonal correction.
29 
30 SourceFiles
31  faceCorrectedSnGrad.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef faceCorrectedSnGrad_H
36 #define faceCorrectedSnGrad_H
37 
38 #include "snGradScheme.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace fv
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class faceCorrectedSnGrad Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 template<class Type>
56 :
57  public snGradScheme<Type>
58 {
59  // Private Member Functions
60 
61  //- Disallow default bitwise assignment
62  void operator=(const faceCorrectedSnGrad&) = delete;
63 
64 
65 public:
66 
67  //- Runtime type information
68  TypeName("faceCorrected");
69 
70 
71  // Constructors
72 
73  //- Construct from mesh
75  :
76  snGradScheme<Type>(mesh)
77  {}
78 
79 
80  //- Construct from mesh and data stream
82  :
83  snGradScheme<Type>(mesh)
84  {}
85 
86 
87  //- Destructor
88  virtual ~faceCorrectedSnGrad();
89 
90 
91  // Member Functions
92 
93  //- Return the interpolation weighting factors for the given field
95  (
96  const VolField<Type>&
97  ) const
98  {
99  return this->mesh().nonOrthDeltaCoeffs();
100  }
101 
102  //- Return true if this scheme uses an explicit correction
103  virtual bool corrected() const
104  {
105  return true;
106  }
107 
108  //- Return the explicit correction to the faceCorrectedSnGrad
109  // for the given field using the gradient of the field
112  (
113  const VolField<Type>&
114  ) const;
115 
116  //- Return the explicit correction to the faceCorrectedSnGrad
117  // for the given field using the gradients of the field components
118  virtual tmp<SurfaceField<Type>>
119  correction(const VolField<Type>&) const;
120 };
121 
122 
123 // * * * * * * * * Template Member Function Specialisations * * * * * * * * //
124 
125 template<>
127 (
128  const volScalarField& vsf
129 ) const;
130 
131 
132 template<>
134 (
135  const volVectorField& vvf
136 ) const;
137 
138 
139 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
140 
141 } // End namespace fv
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 } // End namespace Foam
146 
147 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 
149 #ifdef NoRepository
150  #include "faceCorrectedSnGrad.C"
151 #endif
152 
153 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154 
155 #endif
156 
157 // ************************************************************************* //
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:99
Simple central-difference snGrad scheme with non-orthogonal correction.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
tmp< SurfaceField< Type > > fullGradCorrection(const VolField< Type > &) const
Return the explicit correction to the faceCorrectedSnGrad.
virtual tmp< surfaceScalarField > deltaCoeffs(const VolField< Type > &) const
Return the interpolation weighting factors for the given field.
TypeName("faceCorrected")
Runtime type information.
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &) const
Return the explicit correction to the faceCorrectedSnGrad.
virtual ~faceCorrectedSnGrad()
Destructor.
faceCorrectedSnGrad(const fvMesh &mesh)
Construct from mesh.
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:63
const fvMesh & mesh() const
Return mesh reference.
Definition: snGradScheme.H:117
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
labelList fv(nPoints)