weightedFlux.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) 2019 Norbert Weber, HZDR
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::weightedFlux
28
29Description
30 Weighted flux interpolation scheme class.
31
32 This scheme is used to compute fluxes with variable diffusivity or
33 conductivity, as e.g.
34 - a thermal flux: lambda*grad(T)
35 - a mass flux: D*grad(u)
36 - an electric current: -sigma*grad(potential)
37
38 When using the Gauss theorem to compute a gradient, cell centred values
39 need to be interpolated to the faces. Using this scheme, temperature (T)
40 is weighted by thermal conductivity when being interpolated. Similarly,
41 velocity is weighted by diffusivity (D) and the electric potential by
42 the electric conductivity (sigma). Lambda, D or sigma are read from the
43 object registry - the names need to be specified in fvSchemes as e.g.
44
45 \verbatim
46 gradSchemes
47 {
48 grad(T) Gauss weightedFlux lambda;
49 grad(u) Gauss weightedFlux D;
50 grad(potential) Gauss weightedFlux sigma;
51 }
52 \endverbatim
53
54 For more details, see equation 16 and 17 in
55 \verbatim
56 Weber, N., Beckstein, P., Galindo, V., Starace, M. & Weier, T. (2018).
57 Electro-vortex flow simulation using coupled meshes.
58 Computers and Fluids 168, 101-109.
59 doi:10.1016/j.compfluid.2018.03.047
60 https://arxiv.org/pdf/1707.06546.pdf
61 \endverbatim
62
63Note
64 For support, contact Norbert.Weber@hzdr.de
65
66SourceFiles
67 weightedFlux.C
68
69\*---------------------------------------------------------------------------*/
70
71#ifndef weightedFlux_H
72#define weightedFlux_H
73
75#include "volFields.H"
76
77// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78
79namespace Foam
80{
81
82/*---------------------------------------------------------------------------*\
83 Class weightedFlux Declaration
84\*---------------------------------------------------------------------------*/
85
86template<class Type>
87class weightedFlux
88:
90{
91 // Private Data
92
93 //- Const reference to step-wise pre-gradient factor field
94 const volScalarField& sigma_;
95
96 // Demand-driven data
97
98 //- Face to owner cell distance
99 mutable surfaceScalarField* oDelta_;
100
101 //- Face to neighbour cell distance
102 mutable surfaceScalarField* nDelta_;
103
104
105 // Private Member Functions
106
107 //- Compute face-owner and face-neighbour distance
108 void makeDeltas() const;
109
110 //- No copy assignment
111 void operator=(const weightedFlux&) = delete;
112
113
114protected:
115
116 // Protected Member Functions
117
118 // Storage management
119
120 //- Clear all fields
121 void clearOut();
122
123
124public:
125
126 //- Runtime type information
127 TypeName("weightedFlux");
128
129
130 // Constructors
131
132 //- Construct from Istream.
133 // The name of the flux field is read from the Istream and looked-up
134 // from the mesh objectRegistry
136 (
137 const fvMesh& mesh,
138 Istream& is
139 )
140 :
142 sigma_(this->mesh().objectRegistry::template
143 lookupObject<volScalarField>(word(is))),
144 oDelta_(nullptr),
145 nDelta_(nullptr)
146 {}
147
148 //- Construct from faceFlux and Istream
150 (
151 const fvMesh& mesh,
152 const surfaceScalarField& faceFlux,
153 Istream& is
154 )
155 :
157 sigma_(this->mesh().objectRegistry::template
158 lookupObject<volScalarField>(word(is))),
159 oDelta_(nullptr),
160 nDelta_(nullptr)
161 {}
162
163
164 //- Destructor
165 virtual ~weightedFlux();
166
167
168 // Member Functions
169
170 //- Return the interpolation weighting factors
172 (
174 ) const
175 {
176 return this->mesh().surfaceInterpolation::weights();
177 }
178
179 //- Return the distance between face and owner cell
180 const surfaceScalarField& oDelta() const
181 {
182 if (!oDelta_)
183 {
184 makeDeltas();
185 }
186
187 return *oDelta_;
188 }
189
190 //- Return the distance between face and neighbour cell
191 const surfaceScalarField& nDelta() const
192 {
193 if (!nDelta_)
194 {
195 makeDeltas();
196 }
197
198 return *nDelta_;
199 }
200
201 //- Interpolate the cell values to faces
204 (
206 ) const;
207};
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212} // End namespace Foam
213
214// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215
216#endif
217
218// ************************************************************************* //
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
Weighted flux interpolation scheme class.
Definition: weightedFlux.H:89
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Interpolate the cell values to faces.
Definition: weightedFlux.C:154
virtual ~weightedFlux()
Destructor.
Definition: weightedFlux.C:44
const surfaceScalarField & oDelta() const
Return the distance between face and owner cell.
Definition: weightedFlux.H:179
const surfaceScalarField & nDelta() const
Return the distance between face and neighbour cell.
Definition: weightedFlux.H:190
weightedFlux(const fvMesh &mesh, const surfaceScalarField &faceFlux, Istream &is)
Construct from faceFlux and Istream.
Definition: weightedFlux.H:149
TypeName("weightedFlux")
Runtime type information.
tmp< surfaceScalarField > weights(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
Definition: weightedFlux.H:171
void clearOut()
Clear all fields.
Definition: weightedFlux.C:34
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