phaseStabilisedSnGrad.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 "phaseStabilisedSnGrad.H"
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "localMin.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace fv
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
43 template<class Type>
45 (
46  const fvMesh& mesh,
47  Istream& schemeData
48 )
49 :
51  correctedScheme_
52  (
53  fv::snGradScheme<Type>::New(this->mesh(), schemeData)
54  )
55 {}
56 
57 
58 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59 
60 template<class Type>
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 template<class Type>
69 (
71 ) const
72 {
73  return correctedScheme_->deltaCoeffs(vf);
74 }
75 
76 
77 template<class Type>
80 (
82 ) const
83 {
85  (
86  correctedScheme_().correction(vf)
87  );
88 
89  const volScalarField& alpha
90  (
91  vf.db().template lookupObject<volScalarField>
92  (
93  IOobject::groupName("alpha", vf.group())
94  )
95  );
96 
98  (
100  );
101 
102  if (fv::debug)
103  {
105  << "limiter min: " << min(limiter.primitiveField())
106  << " max: "<< max(limiter.primitiveField())
107  << " avg: " << average(limiter.primitiveField()) << endl;
108  }
109 
110  return limiter*corr;
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
115 
116 } // End namespace fv
117 
118 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
119 
120 } // End namespace Foam
121 
122 // ************************************************************************* //
Foam::surfaceFields.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Generic GeometricField class.
fvMesh & mesh
dimensionedScalar pos(const dimensionedScalar &ds)
labelList fv(nPoints)
static word groupName(Name name, const word &group)
phaseStabilisedSnGrad(const fvMesh &mesh, Istream &schemeData)
Construct from mesh and data stream.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the phaseStabilisedSnGrad.
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)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
const Mesh & mesh() const
Return mesh.
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:60
virtual ~phaseStabilisedSnGrad()
Destructor.
virtual tmp< surfaceScalarField > deltaCoeffs(const GeometricField< Type, fvPatchField, volMesh > &vf) const
Return the interpolation weighting factors for the given field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
A class for managing temporary objects.
Definition: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
Local minimum interpolation scheme in which the face value is set to the minimum of the two neighbour...
Definition: localMin.H:52
Namespace for OpenFOAM.
#define InfoInFunction
Report an information message using Foam::Info.