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-2018 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,
40  fieldTable& fields,
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().timeName(),
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 // ************************************************************************* //
Foam::surfaceFields.
tmp< surfaceScalarField > weights() const
Return the interpolation weighting factors.
Definition: upwind.H:115
virtual tmp< surfaceScalarField > limiter(const GeometricField< Type, fvPatchField, volMesh > &) const =0
Return the interpolation weighting factors.
multivariateSelectionScheme(const fvMesh &mesh, const typename multivariateSurfaceInterpolationScheme< Type >::fieldTable &fields, const surfaceScalarField &faceFlux, Istream &schemeData)
Construct for field, faceFlux and Istream.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
const dimensionSet dimless
fvMesh & mesh
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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 complex &)
Return a string representation of a complex.
Definition: complex.C:47
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
Upwind interpolation scheme class.
Definition: upwind.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Abstract base class for multi-variate surface interpolation schemes.
Abstract base class for limited surface interpolation schemes.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98