weighted.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 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::weighted
28
29Group
30 grpFvSurfaceInterpolationSchemes
31
32Description
33 Interpolation scheme class using weights looked-up from the objectRegistry.
34
35SourceFiles
36 weighted.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef weighted_H
41#define weighted_H
42
44#include "volFields.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class weighted Declaration
53\*---------------------------------------------------------------------------*/
54
55template<class Type>
56class weighted
57:
59{
60 // Private member data
61
62 const surfaceScalarField& weights_;
63
64
65 // Private Member Functions
66
67 //- No copy assignment
68 void operator=(const weighted&) = delete;
69
70
71public:
72
73 //- Runtime type information
74 TypeName("weighted");
75
76
77 // Constructors
78
79 //- Construct from weights
81 :
83 weights_(weights)
84 {}
85
86 //- Construct from Istream
87 weighted(const fvMesh& mesh, Istream& is)
88 :
90 weights_
91 (
92 this->mesh().objectRegistry::template
93 lookupObject<const surfaceScalarField>
94 (
95 word(is)
96 )
97 )
98 {}
99
100 //- Construct from faceFlux and Istream
102 (
103 const fvMesh& mesh,
104 const surfaceScalarField&,
105 Istream& is
106 )
107 :
109 weights_
110 (
111 this->mesh().objectRegistry::template
112 lookupObject<const surfaceScalarField>
113 (
114 word(is)
115 )
116 )
117 {}
118
119
120 // Member Functions
121
122 //- Return the interpolation weighting factors
124 (
126 ) const
127 {
128 return weights_;
129 }
130};
131
132
133// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134
135} // End namespace Foam
136
137// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138
139#endif
140
141// ************************************************************************* //
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
Registry of regIOobjects.
Abstract base class for surface interpolation schemes.
const fvMesh & mesh() const
Return mesh reference.
A class for managing temporary objects.
Definition: tmp.H:65
Interpolation scheme class using weights looked-up from the objectRegistry.
Definition: weighted.H:58
weighted(const fvMesh &mesh, const surfaceScalarField &, Istream &is)
Construct from faceFlux and Istream.
Definition: weighted.H:101
TypeName("weighted")
Runtime type information.
weighted(const surfaceScalarField &weights)
Construct from weights.
Definition: weighted.H:79
weighted(const fvMesh &mesh, Istream &is)
Construct from Istream.
Definition: weighted.H:86
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: weighted.H:123
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73