deferred.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) 2021 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::deferred
26 
27 Description
28  Deferred correction interpolation scheme derived from upwind
29  which returns upwind weighting factors and an explicit correction obtained
30  from the specified scheme.
31 
32  This ensures that the transport matrix generated is at least diagonally
33  equal and avoids the need for relaxation of the matrix (which can affect
34  conservation) for transient running.
35 
36 SourceFiles
37  deferred.C
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #ifndef deferred_H
42 #define deferred_H
43 
44 #include "upwind.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class deferred Declaration
53 \*---------------------------------------------------------------------------*/
54 
55 template<class Type>
56 class deferred
57 :
58  public upwind<Type>
59 {
60  // Private member data
61 
62  //- Scheme
64 
65 
66 public:
67 
68  //- Runtime type information
69  TypeName("deferred");
70 
71 
72  // Constructors
73 
74  //- Construct from Istream.
75  // The name of the flux field is read from the Istream and looked-up
76  // from the mesh objectRegistry
78  (
79  const fvMesh& mesh,
80  Istream& schemeData
81  )
82  :
83  upwind<Type>(mesh, schemeData),
84  tScheme_
85  (
86  surfaceInterpolationScheme<Type>::New(mesh, schemeData)
87  )
88  {}
89 
90  //- Construct from faceFlux and Istream
92  (
93  const fvMesh& mesh,
94  const surfaceScalarField& faceFlux,
95  Istream& schemeData
96  )
97  :
98  upwind<Type>(mesh, faceFlux, schemeData),
99  tScheme_
100  (
102  (
103  mesh,
104  faceFlux,
105  schemeData
106  )
107  )
108  {}
109 
110  //- Disallow default bitwise copy construction
111  deferred(const deferred&) = delete;
112 
113 
114  // Member Functions
115 
116  //- Return true if this scheme uses an explicit correction
117  virtual bool corrected() const
118  {
119  return true;
120  }
121 
122  //- Return the explicit correction to the face-interpolate
125  (
127  ) const
128  {
130  (
131  "deferred::correction(" + vf.name() + ')',
132  tScheme_->interpolate(vf)
134  );
135  }
136 
137 
138  // Member Operators
139 
140  //- Disallow default bitwise assignment
141  void operator=(const deferred&) = delete;
142 };
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 } // End namespace Foam
148 
149 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
150 
151 #endif
152 
153 // ************************************************************************* //
const word & name() const
Return name.
Definition: IOobject.H:303
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
static tmp< limitedSurfaceInterpolationScheme< Type > > New(const fvMesh &mesh, Istream &schemeData)
Return new tmp interpolation scheme.
virtual bool corrected() const
Return true if this scheme uses an explicit correction.
Definition: deferred.H:116
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the explicit correction to the face-interpolate.
Definition: deferred.H:124
deferred(const fvMesh &mesh, Istream &schemeData)
Construct from Istream.
Definition: deferred.H:77
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const tmp< surfaceScalarField > &, const tmp< surfaceScalarField > &)
Return the face-interpolate of the given cell field.
void operator=(const deferred &)=delete
Disallow default bitwise assignment.
Upwind differencing scheme class.
Definition: upwind.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class for managing temporary objects.
Definition: PtrList.H:53
Deferred correction interpolation scheme derived from upwind which returns upwind weighting factors a...
Definition: deferred.H:55
Abstract base class for surface interpolation schemes.
TypeName("deferred")
Runtime type information.
Namespace for OpenFOAM.