multivariateSelectionScheme.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
31#include "volFields.H"
32#include "surfaceFields.H"
33#include "upwind.H"
34
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38template<class Type>
40(
41 const fvMesh& mesh,
44 const surfaceScalarField& faceFlux,
45 Istream& schemeData
46)
47:
49 (
50 mesh,
51 fields,
52 faceFlux,
53 schemeData
54 ),
55 schemes_(schemeData),
56 faceFlux_(faceFlux),
57 weights_
58 (
60 (
61 "multivariateWeights",
62 mesh.time().timeName(),
63 mesh
64 ),
65 mesh,
67 )
68{
69 auto iter = this->fields().cbegin();
70
72 (
74 (
75 mesh,
76 faceFlux_,
77 schemes_.lookup(iter()->name())
78 )().limiter(*iter())
79 );
80
81 for (++iter; iter.good(); ++iter)
82 {
83 limiter = min
84 (
85 limiter,
87 (
88 mesh,
89 faceFlux_,
90 schemes_.lookup(iter()->name())
91 )().limiter(*iter())
92 );
93 }
94
95 weights_ =
96 limiter*mesh.surfaceInterpolation::weights()
97 + (scalar(1) - limiter)*upwind<Type>(mesh, faceFlux_).weights();
98}
99
100
101// ************************************************************************* //
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
const_iterator cbegin() const
Iterator to first item in list with const access.
Definition: LList.H:509
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for limited surface interpolation schemes.
Generic multi-variate discretisation scheme class for which any of the NVD, CNVD or NVDV schemes may ...
Abstract base class for multi-variate surface interpolation schemes.
const fieldTable & fields() const
Return fields to be interpolated.
Upwind differencing scheme class.
Definition: upwind.H:59
tmp< surfaceScalarField > weights() const
Return the interpolation weighting factors.
Definition: upwind.H:133
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
const dimensionSet dimless
Dimensionless.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:38
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::surfaceFields.