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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "fv.H"
29 #include "limitedSnGrad.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
32 #include "localMax.H"
33 #include "fvcCellReduce.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 namespace fv
43 {
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
47 template<class Type>
49 {}
50 
51 
52 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
53 
54 template<class Type>
57 (
59 ) const
60 {
62  (
63  correctedScheme_().correction(vf)
64  );
65 
67  (
68  min
69  (
70  limitCoeff_
71  *mag(snGradScheme<Type>::snGrad(vf, deltaCoeffs(vf), "SndGrad"))
72  /(
73  (1 - limitCoeff_)*mag(corr)
74  + dimensionedScalar("small", corr.dimensions(), SMALL)
75  ),
76  dimensionedScalar("one", dimless, 1.0)
77  )
78  );
79 
80  if (fv::debug)
81  {
83  << "limiter min: " << min(limiter.primitiveField())
84  << " max: "<< max(limiter.primitiveField())
85  << " avg: " << average(limiter.primitiveField()) << endl;
86 
87 
88  if (fv::debug & 2)
89  {
90  static scalar oldTime = -1;
91  static label subIter = 0;
92  if (vf.mesh().time().value() != oldTime)
93  {
94  oldTime = vf.mesh().time().value();
95  subIter = 0;
96  }
97  else
98  {
99  ++subIter;
100  }
101  word fieldName("limiter_" + Foam::name(subIter));
102 
104  (
105  IOobject
106  (
107  fieldName,
108  vf.mesh().time().timeName(),
109  vf.mesh(),
112  false
113  ),
114  fvc::cellReduce(limiter, minEqOp<scalar>(), scalar(1.0))
115  );
116  Info<< "Writing limiter field to " << volLimiter.objectPath()
117  << endl;
118  volLimiter.write();
119  }
120  }
121 
122  return limiter*corr;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
127 
128 } // End namespace fv
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 } // End namespace Foam
133 
134 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:325
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
limitedSnGrad.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
fv.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
oldTime
Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));if(U.nOldTimes()){ volVectorField *Uold=&U.oldTime();volScalarField *Kold=&K.oldTime();*Kold==0.5 *magSqr(*Uold);while(Uold->nOldTimes()) { Uold=&Uold-> oldTime()
Definition: createK.H:12
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
surfaceFields.H
Foam::surfaceFields.
localMax.H
Foam::correction
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
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:57
Foam::minEqOp
Definition: ops.H:81
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::fvc::cellReduce
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fv::limitedSnGrad::~limitedSnGrad
virtual ~limitedSnGrad()
Destructor.
Definition: limitedSnGrad.C:48
fv
labelList fv(nPoints)
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fv::snGradScheme
Abstract base class for snGrad schemes.
Definition: snGradScheme.H:65
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
fvcCellReduce.H
Construct a volume field from a surface field using a combine operator.
Foam::limiter
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:37
Foam::average
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:328