phaseStabilisedTemplates.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) 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 
26 #include "phaseStabilised.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class Type>
32 (
33  const surfaceScalarField& faceFlux
34 ) const
35 {
36  return faceFlux.db().lookupObject<volScalarField>
37  (
38  IOobject::groupName("alpha", faceFlux.group())
39  );
40 }
41 
42 
43 template<class Type>
46 {
47  return pos(upwind_.interpolate(alpha_) - 1e-3);
48 }
49 
50 
51 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
52 
53 template<class Type>
55 (
56  const fvMesh& mesh,
57  Istream& is
58 )
59 :
60  surfaceInterpolationScheme<Type>(mesh),
61  faceFlux_
62  (
63  mesh.lookupObject<surfaceScalarField>
64  (
65  word(is)
66  )
67  ),
68  tScheme_
69  (
70  surfaceInterpolationScheme<Type>::New(mesh, faceFlux_, is)
71  ),
72  upwind_(mesh, faceFlux_),
73  alpha_(lookupAlpha(faceFlux_))
74 {}
75 
76 
77 template<class Type>
79 (
80  const fvMesh& mesh,
81  const surfaceScalarField& faceFlux,
82  Istream& is
83 )
84 :
85  surfaceInterpolationScheme<Type>(mesh),
86  faceFlux_(faceFlux),
87  tScheme_
88  (
89  surfaceInterpolationScheme<Type>::New(mesh, faceFlux, is)
90  ),
91  upwind_(mesh, faceFlux_),
92  alpha_(lookupAlpha(faceFlux_))
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
98 template<class Type>
100 (
101  const VolField<Type>& vf
102 ) const
103 {
104  const surfaceScalarField lambdaf(this->lambdaf());
105  return lambdaf*tScheme_().weights(vf) + (1 - lambdaf)*upwind_.weights();
106 }
107 
108 
109 template<class Type>
112 (
113  const VolField<Type>& vf
114 ) const
115 {
116  if (tScheme_().corrected())
117  {
118  return lambdaf()*tScheme_().correction(vf);
119  }
120  else
121  {
122  return tmp<SurfaceField<Type>>
123  (
124  nullptr
125  );
126  }
127 }
128 
129 
130 // ************************************************************************* //
Generic GeometricField class.
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:101
Phase-stabilised interpolation scheme.
phaseStabilised(const fvMesh &mesh, Istream &is)
Construct from mesh and Istream.
tmp< surfaceScalarField > weights(const VolField< Type > &vf) const
Return the interpolation weighting factors.
virtual tmp< SurfaceField< Type > > correction(const VolField< Type > &vf) const
Return the explicit correction to the face-interpolate.
Abstract base class for surface interpolation schemes.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar e
Elementary charge.
dimensionedScalar pos(const dimensionedScalar &ds)
SurfaceField< scalar > surfaceScalarField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)