uncorrectedSnGrad.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-2018 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::uncorrectedSnGrad
26 
27 Description
28  Simple central-difference snGrad scheme using the non-orthogonal mesh
29  delta-coefficients but without non-orthogonal correction.
30 
31 SourceFiles
32  uncorrectedSnGrad.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef uncorrectedSnGrad_H
37 #define uncorrectedSnGrad_H
38 
39 #include "snGradScheme.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace fv
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class uncorrectedSnGrad Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class Type>
57 :
58  public snGradScheme<Type>
59 {
60  // Private Member Functions
61 
62  //- Disallow default bitwise assignment
63  void operator=(const uncorrectedSnGrad&);
64 
65 
66 public:
67 
68  //- Runtime type information
69  TypeName("uncorrected");
70 
71 
72  // Constructors
73 
74  //- Construct from mesh
76  :
77  snGradScheme<Type>(mesh)
78  {}
79 
80 
81  //- Construct from mesh and data stream
83  :
84  snGradScheme<Type>(mesh)
85  {}
86 
87 
88  //- Destructor
89  virtual ~uncorrectedSnGrad();
90 
91 
92  // Member Functions
93 
94  //- Return the interpolation weighting factors for the given field
96  (
98  ) const
99  {
100  return this->mesh().nonOrthDeltaCoeffs();
101  }
102 
103  //- Return true if this scheme uses an explicit correction
104  virtual bool corrected() const
105  {
106  return false;
107  }
108 
109  //- Return the explicit correction to the uncorrectedSnGrad
110  // for the given field
113 };
114 
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117 
118 } // End namespace fv
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace Foam
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 #ifdef NoRepository
127  #include "uncorrectedSnGrad.C"
128 #endif
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 #endif
133 
134 // ************************************************************************* //
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Simple central-difference snGrad scheme without non-orthogonal correction.
Generic GeometricField class.
const surfaceScalarField & nonOrthDeltaCoeffs() const
Return reference to non-orthogonal cell-centre difference.
virtual ~uncorrectedSnGrad()
Destructor.
labelList fv(nPoints)
Simple central-difference snGrad scheme using the non-orthogonal mesh delta-coefficients but without ...
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:60
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the uncorrectedSnGrad.
const fvMesh & mesh() const
Return mesh reference.
Definition: snGradScheme.H:123
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
uncorrectedSnGrad(const fvMesh &mesh)
Construct from mesh.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors for the given field.
TypeName("uncorrected")
Runtime type information.
Namespace for OpenFOAM.