multivariateUpwind.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) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
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
26Class
27 Foam::multivariateUpwind
28
29Description
30 Upwind-difference form of the multivariate surfaceInterpolationScheme.
31
32SourceFiles
33 multivariateUpwindmake.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef multivariateUpwind_H
38#define multivariateUpwind_H
39
41#include "surfaceFields.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48/*---------------------------------------------------------------------------*\
49 Class multivariateUpwind Declaration
50\*---------------------------------------------------------------------------*/
51
52template<class Type>
54:
56{
57 // Private data
58
59 const surfaceScalarField& faceFlux_;
60
61
62 // Private Member Functions
63
64 //- No copy construct
66
67 //- No copy assignment
68 void operator=(const multivariateUpwind&) = delete;
69
70
71public:
72
73 //- Runtime type information
74 TypeName("upwind");
75
76
77 // Constructors
78
79 //- Construct for field, faceFlux and Istream
81 (
82 const fvMesh& mesh,
85 const surfaceScalarField& faceFlux,
86 Istream& schemeData
87 )
88 :
90 (
91 mesh,
92 fields,
93 faceFlux,
94 schemeData
95 ),
96 faceFlux_(faceFlux)
97 {}
98
99
100 // Member Operators
101
102 //- surfaceInterpolationScheme sub-class returned by operator(field)
103 // which is used as the interpolation scheme for the field
104 class fieldScheme
105 :
106 public multivariateSurfaceInterpolationScheme<Type>::fieldScheme
107 {
108 // Private data
109
110 const surfaceScalarField& faceFlux_;
111
112 public:
113
114 // Constructors
117 (
119 const surfaceScalarField& faceFlux
120 )
121 :
124 faceFlux_(faceFlux)
125 {}
126
127
128 // Member Functions
129
130 //- Return the interpolation weighting factors
132 (
134 ) const
135 {
136 return pos0(faceFlux_);
137 }
138 };
141 (
143 ) const
144 {
146 (
147 new fieldScheme(field, faceFlux_)
148 );
149 }
150};
151
152
153// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
154
155} // End namespace Foam
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159#endif
160
161// ************************************************************************* //
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for multi-variate surface interpolation schemes.
const fieldTable & fields() const
Return fields to be interpolated.
surfaceInterpolationScheme sub-class returned by operator(field)
fieldScheme(const GeometricField< Type, fvPatchField, volMesh > &field, const surfaceScalarField &faceFlux)
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Upwind-difference form of the multivariate surfaceInterpolationScheme.
multivariateUpwind(const fvMesh &mesh, const typename multivariateSurfaceInterpolationScheme< Type >::fieldTable &fields, const surfaceScalarField &faceFlux, Istream &schemeData)
Construct for field, faceFlux and Istream.
TypeName("upwind")
Runtime type information.
A class for managing temporary objects.
Definition: tmp.H:65
rDeltaTY field()
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
Foam::surfaceFields.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73