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  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "fv.H"
30 #include "limitedSnGrad.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "localMax.H"
34 #include "fvcCellReduce.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace fv
44 {
45 
46 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
47 
48 template<class Type>
51 (
53 ) const
54 {
56  (
57  correctedScheme_().correction(vf)
58  );
59 
61  (
62  min
63  (
64  limitCoeff_
65  *mag(snGradScheme<Type>::snGrad(vf, deltaCoeffs(vf), "SndGrad"))
66  /(
67  (1 - limitCoeff_)*mag(corr)
68  + dimensionedScalar("small", corr.dimensions(), SMALL)
69  ),
70  dimensionedScalar("one", dimless, 1.0)
71  )
72  );
73 
74  if (fv::debug)
75  {
77  << "limiter min: " << min(limiter.primitiveField())
78  << " max: "<< max(limiter.primitiveField())
79  << " avg: " << average(limiter.primitiveField()) << endl;
80 
81 
82  if (fv::debug & 2)
83  {
84  static scalar oldTime = -1;
85  static label subIter = 0;
86  if (vf.mesh().time().value() != oldTime)
87  {
88  oldTime = vf.mesh().time().value();
89  subIter = 0;
90  }
91  else
92  {
93  ++subIter;
94  }
95  word fieldName("limiter_" + Foam::name(subIter));
96 
98  (
99  IOobject
100  (
101  fieldName,
102  vf.mesh().time().timeName(),
103  vf.mesh().thisDb(),
107  ),
109  );
110  Info<< "Writing limiter field to " << volLimiter.objectPath()
111  << endl;
112  volLimiter.write();
113  }
114  }
115 
116  return limiter*corr;
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace fv
123 
124 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
125 
126 } // End namespace Foam
127 
128 // ************************************************************************* //
Foam::surfaceFields.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1)
Info<< "Creating field kinetic energy K\"<< 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
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Construct a volume field from a surface field using a combine operator.
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:31
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:284
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the explicit correction to the limitedSnGrad for the given field.
Definition: limitedSnGrad.C:44
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const Mesh & mesh() const noexcept
Return mesh.
int debug
Static debugging option.
Abstract base class for runtime selected snGrad surface normal gradient schemes.
Definition: snGradScheme.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Do not request registration (bool: false)
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
#define InfoInFunction
Report an information message using Foam::Info.