limitedSnGrad.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 | Copyright (C) 2011-2012 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 "fv.H"
27 #include "limitedSnGrad.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "localMax.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace fv
40 {
41 
42 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
43 
44 template<class Type>
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
50 
51 template<class Type>
54 (
56 ) const
57 {
59  (
60  correctedScheme_().correction(vf)
61  );
62 
64  (
65  min
66  (
67  limitCoeff_
68  *mag(snGradScheme<Type>::snGrad(vf, deltaCoeffs(vf), "SndGrad"))
69  /(
70  (1 - limitCoeff_)*mag(corr)
71  + dimensionedScalar("small", corr.dimensions(), SMALL)
72  ),
73  dimensionedScalar("one", dimless, 1.0)
74  )
75  );
76 
77  if (fv::debug)
78  {
79  Info<< "limitedSnGrad :: limiter min: " << min(limiter.internalField())
80  << " max: "<< max(limiter.internalField())
81  << " avg: " << average(limiter.internalField()) << endl;
82  }
83 
84  return limiter*corr;
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89 
90 } // End namespace fv
91 
92 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93 
94 } // End namespace Foam
95 
96 // ************************************************************************* //
volFields.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
limitedSnGrad.H
fv.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::MULES::limiter
void limiter(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
Definition: MULESTemplates.C:159
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
localMax.H
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fv::limitedSnGrad::correction
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the limitedSnGrad.
Definition: limitedSnGrad.C:54
Foam::Info
messageStream Info
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::fv::limitedSnGrad::~limitedSnGrad
virtual ~limitedSnGrad()
Destructor.
Definition: limitedSnGrad.C:45
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
fv
labelList fv(nPoints)
Foam::fv::snGradScheme
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:60
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:335