multivariateSelectionScheme.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "upwind.H"
31 
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  const fvMesh& mesh,
41  const surfaceScalarField& faceFlux,
42  Istream& schemeData
43 )
44 :
46  (
47  mesh,
48  fields,
49  faceFlux,
50  schemeData
51  ),
52  schemes_(schemeData),
53  faceFlux_(faceFlux),
54  weights_
55  (
56  IOobject
57  (
58  "multivariateWeights",
59  mesh.time().name(),
60  mesh
61  ),
62  mesh,
63  dimless
64  )
65 {
67  fieldTable::const_iterator iter = this->fields().begin();
68 
70  (
72  (
73  mesh,
74  faceFlux_,
75  schemes_.lookup(iter()->name())
76  )().limiter(*iter())
77  );
78 
79  for (++iter; iter != this->fields().end(); ++iter)
80  {
81  limiter = min
82  (
83  limiter,
85  (
86  mesh,
87  faceFlux_,
88  schemes_.lookup(iter()->name())
89  )().limiter(*iter())
90  );
91  }
92 
93  weights_ =
94  limiter*mesh.surfaceInterpolation::weights()
95  + (scalar(1) - limiter)*upwind<Type>(mesh, faceFlux_).weights();
96 }
97 
98 
99 // ************************************************************************* //
Generic GeometricField class.
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract base class for limited surface interpolation schemes.
multivariateSelectionScheme(const fvMesh &mesh, const typename multivariateSurfaceInterpolationScheme< Type >::fieldTable &fields, const surfaceScalarField &faceFlux, Istream &schemeData)
Construct for field, faceFlux and Istream.
Abstract base class for multi-variate surface interpolation schemes.
const fieldTable & fields() const
Return fields to be interpolated.
Upwind interpolation scheme class.
Definition: upwind.H:54
tmp< surfaceScalarField > weights() const
Return the interpolation weighting factors.
Definition: upwind.H:115
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
void limiter(surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
Foam::surfaceFields.