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-------------------------------------------------------------------------------
11License
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
38namespace Foam
39{
40
41// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42
43namespace fv
44{
45
46// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
47
48template<class Type>
49tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
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 (
100 (
101 fieldName,
102 vf.mesh().time().timeName(),
103 vf.mesh(),
106 false
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// ************************************************************************* //
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
fileName objectPath() const
The complete path + object name.
Definition: IOobjectI.H:214
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > correction(const GeometricField< Type, fvPatchField, volMesh > &) const
Definition: limitedSnGrad.C:51
Abstract base class for runtime selected snGrad surface normal gradient schemes.
Definition: snGradScheme.H:79
virtual bool write(const bool valid=true) const
Write using setting from DB.
A class for handling words, derived from Foam::string.
Definition: word.H:68
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()
Construct a volume field from a surface field using a combine operator.
#define InfoInFunction
Report an information message using Foam::Info.
tmp< GeometricField< Type, fvPatchField, volMesh > > cellReduce(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const CombineOp &cop, const Type &nullValue=pTraits< Type >::zero)
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &df)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
tmp< areaScalarField > limiter(const areaScalarField &phi)
Definition: faNVDscheme.C:38
labelList fv(nPoints)
Foam::surfaceFields.