limitVelocity.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) 2021-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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 "limitVelocity.H"
29#include "areaFields.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace fa
37{
40 (
41 option,
44 );
45}
46}
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
51(
52 const word& name,
53 const word& modelType,
54 const dictionary& dict,
55 const fvPatch& patch
56)
57:
58 faceSetOption(name, modelType, dict, patch),
59 UName_(coeffs_.getOrDefault<word>("U", "U")),
60 max_(coeffs_.get<scalar>("max"))
61{
63 applied_.setSize(1, false);
64}
65
66
67// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68
70{
72 {
73 coeffs_.readEntry("max", max_);
74
75 return true;
76 }
77
78 return false;
79}
80
81
83{
84 const scalar maxSqrU = sqr(max_);
85
86 // Count nTotFaces ourselves
87 // (maybe only applying on a subset)
88 label nFacesAbove(0);
89 const label nTotFaces(returnReduce(faces_.size(), sumOp<label>()));
90
91 vectorField& Uif = U.primitiveFieldRef();
92
93 for (const label facei : faces_)
94 {
95 auto& Uval = Uif[facei];
96
97 const scalar magSqrUi = magSqr(Uval);
98
99 if (magSqrUi > maxSqrU)
100 {
101 Uval *= sqrt(maxSqrU/max(magSqrUi, SMALL));
102 ++nFacesAbove;
103 }
104 }
105
106 // Handle boundaries in the case of 'all'
107 label nEdgesAbove(0);
108
110 {
111 for (faPatchVectorField& Up : U.boundaryFieldRef())
112 {
113 if (!Up.fixesValue())
114 {
115 for (auto& Uval : Up)
116 {
117 const scalar magSqrUi = magSqr(Uval);
118
119 if (magSqrUi > maxSqrU)
120 {
121 Uval *= sqrt(maxSqrU/max(magSqrUi, SMALL));
122 ++nEdgesAbove;
123 }
124 }
125 }
126 }
127 }
128
129 // Percent, max 2 decimal places
130 const auto percent = [](scalar num, label denom) -> scalar
131 {
132 return (denom ? 1e-2*round(1e4*num/denom) : 0);
133 };
134
135
136 reduce(nFacesAbove, sumOp<label>());
137 reduce(nEdgesAbove, sumOp<label>());
138
139 Info<< type() << ' ' << name_ << " Limited "
140 << nFacesAbove << " ("
141 << percent(nFacesAbove, nTotFaces)
142 << "%) of faces, with max limit " << max_ << endl;
143
144 if (nFacesAbove || nEdgesAbove)
145 {
146 // We've changed internal values so give
147 // boundary conditions opportunity to correct
148 U.correctBoundaryConditions();
149 }
150}
151
152
153// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void setSize(const label n)
Alias for resize()
Definition: List.H:218
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual bool read()
Re-read model coefficients if they have changed.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Intermediate abstract class for handling face-set options for the derived faOptions.
bool useSubMesh() const noexcept
True if sub-selection should be used.
Limits the maximum velocity magnitude to the specified max value.
word UName_
Name of operand velocity field.
Base abstract class for handling finite area options (i.e. faOption).
Definition: faOption.H:134
List< bool > applied_
Applied flag list - corresponds to each fieldNames_ entry.
Definition: faOption.H:167
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: faOption.H:164
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce (copy) and return value.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict
volScalarField & e
Definition: createFields.H:11