Loading...
Searching...
No Matches
multivariateScheme.C
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 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "volFields.H"
30#include "surfaceFields.H"
31#include "upwind.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35template<class Type, class Scheme>
36Foam::multivariateScheme<Type, Scheme>::multivariateScheme
37(
38 const fvMesh& mesh,
39 const typename multivariateSurfaceInterpolationScheme<Type>::
41 const surfaceScalarField& faceFlux,
42 Istream& schemeData
43)
44:
45 multivariateSurfaceInterpolationScheme<Type>
46 (
47 mesh,
48 fields,
49 faceFlux,
50 schemeData
51 ),
52 Scheme::LimiterType(schemeData),
53 faceFlux_(faceFlux),
54 weights_
55 (
57 (
58 "multivariateWeights",
59 mesh.time().timeName(),
60 mesh
61 ),
62 mesh,
64 )
65{
66 auto iter = this->fields().cbegin();
67
69 (
70 Scheme(mesh, faceFlux_, *this).limiter(*iter())
71 );
72
73 for (++iter; iter.good(); ++iter)
74 {
75 limiter = min
76 (
77 limiter,
78 Scheme(mesh, faceFlux_, *this).limiter(*iter())
79 );
80 }
81
82 weights_ =
83 limiter*mesh.surfaceInterpolation::weights()
84 + (scalar(1) - limiter)*upwind<Type>(mesh, faceFlux_).weights();
85}
86
87
88// ************************************************************************* //
const_iterator cbegin() const
const_iterator set to the beginning of the HashTable
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const fieldTable & fields() const
Return fields to be interpolated.
Upwind differencing scheme class.
Definition upwind.H:55
tmp< surfaceScalarField > weights() const
Return the interpolation weighting factors.
Definition upwind.H:145
word timeName
Definition getTimeIndex.H:3
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition faNVDscheme.C:31
Foam::surfaceFields.