deferredCorrection.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Class
27  Foam::deferredCorrection
28 
29 Description
30  Deferred correction interpolation scheme wrapper around a run-time
31  selectable base scheme.
32 
33  The characteristics of the base scheme are recovered by applying an
34  explicit correction to the upwind scheme weights.
35 
36 Usage
37  Example of the \c deferredCorrection scheme applied to the \c linear
38  scheme:
39  \verbatim
40  divSchemes
41  {
42  .
43  .
44  div(phi,U) Gauss deferredCorrection linear;
45  .
46  .
47  }
48  \endverbatim
49 
50 SourceFiles
51  deferredCorrection.C
52 
53 SeeAlso
54  Foam::upwind
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #ifndef deferredCorrection_H
59 #define deferredCorrection_H
60 
61 #include "upwind.H"
62 
63 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65 namespace Foam
66 {
67 
68 /*---------------------------------------------------------------------------*\
69  Class deferredCorrection Declaration
70 \*---------------------------------------------------------------------------*/
71 
72 template<class Type>
74 :
75  public upwind<Type>
76 {
77  // Private data
78 
79  //- Base scheme
81 
82 
83  // Private Member Functions
84 
85  //- No copy construct
86  deferredCorrection(const deferredCorrection&) = delete;
87 
88  //- No copy assignment
89  void operator=(const deferredCorrection&) = delete;
90 
91 
92 public:
93 
94  //- Runtime type information
95  TypeName("deferredCorrection");
96 
97 
98  // Constructors
99 
100  //- Construct from Istream.
101  // The name of the flux field is read from the Istream and looked-up
102  // from the mesh objectRegistry
104  (
105  const fvMesh& mesh,
106  Istream& is
107  )
108  :
109  upwind<Type>(mesh, is),
110  tbaseScheme_
111  (
113  )
114  {}
115 
116  //- Construct from faceFlux and Istream
118  (
119  const fvMesh& mesh,
120  const surfaceScalarField& faceFlux,
121  Istream& is
122  )
123  :
124  upwind<Type>(mesh, faceFlux, is),
125  tbaseScheme_
126  (
127  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
128  )
129  {}
130 
131 
132  // Member Functions
133 
134  using upwind<Type>::weights;
135 
136  //- Return true if this scheme uses an explicit correction
137  virtual bool corrected() const
138  {
139  return true;
140  }
141 
142  //- Return the explicit correction to the face-interpolate
144  correction
145  (
147  ) const;
148 };
149 
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace Foam
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 #endif
158 
159 // ************************************************************************* //
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::upwind
Upwind differencing scheme class.
Definition: upwind.H:56
Foam::limitedSurfaceInterpolationScheme::New
static tmp< limitedSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
Definition: limitedSurfaceInterpolationScheme.C:39
Foam::deferredCorrection::corrected
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: deferredCorrection.H:136
Foam::deferredCorrection::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the face-interpolate.
Definition: deferredCorrection.C:35
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:61
upwind.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::deferredCorrection::TypeName
TypeName("deferredCorrection")
Runtime type information.
Foam::deferredCorrection
Deferred correction interpolation scheme wrapper around a run-time selectable base scheme.
Definition: deferredCorrection.H:72
Foam::surfaceInterpolationScheme
Abstract base class for surface interpolation schemes.
Definition: surfaceInterpolationScheme.H:57
Foam::surfaceInterpolationScheme::mesh
const fvMesh & mesh() const
Return mesh reference.
Definition: surfaceInterpolationScheme.H:144
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >