Loading...
Searching...
No Matches
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-2025 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"
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 fvMesh& mesh,
56 const word& defaultAreaName
57)
58:
59 fa::faceSetOption(name, modelType, dict, mesh, defaultAreaName),
60 UName_(coeffs_.getOrDefault<word>("U", "U")),
61 max_(coeffs_.get<scalar>("max"))
62{
63 fieldNames_.resize(1, UName_);
64 applied_.resize(1, false);
65}
66
67
68// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69
71{
73 {
74 coeffs_.readEntry("max", max_);
75
76 return true;
77 }
78
79 return false;
80}
81
82
84{
85 const scalar maxSqrU = sqr(max_);
86
87 // Count nTotFaces ourselves
88 // (maybe only applying on a subset)
89 label nFacesAbove(0);
90 const label nTotFaces(returnReduce(faces_.size(), sumOp<label>()));
91
92 vectorField& Uif = U.primitiveFieldRef();
93
94 for (const label facei : faces_)
95 {
96 auto& Uval = Uif[facei];
97
98 const scalar magSqrUi = magSqr(Uval);
99
100 if (magSqrUi > maxSqrU)
101 {
102 Uval *= sqrt(maxSqrU/max(magSqrUi, SMALL));
103 ++nFacesAbove;
104 }
105 }
106
107 // Handle boundaries in the case of 'all'
108 label nEdgesAbove(0);
109
111 {
112 for (faPatchVectorField& Up : U.boundaryFieldRef())
113 {
114 if (!Up.fixesValue())
115 {
116 for (auto& Uval : Up)
117 {
118 const scalar magSqrUi = magSqr(Uval);
119
120 if (magSqrUi > maxSqrU)
121 {
122 Uval *= sqrt(maxSqrU/max(magSqrUi, SMALL));
123 ++nEdgesAbove;
124 }
125 }
126 }
127 }
128 }
129
130 // Percent, max 2 decimal places
131 const auto percent = [](scalar num, label denom) -> scalar
132 {
133 return (denom ? 1e-2*round(1e4*num/denom) : 0);
134 };
135
136
137 reduce(nFacesAbove, sumOp<label>());
138 reduce(nEdgesAbove, sumOp<label>());
139
140 Info<< type() << ' ' << name_ << " Limited "
141 << nFacesAbove << " ("
142 << percent(nFacesAbove, nTotFaces)
143 << "%) of faces, with max limit " << max_ << endl;
144
145 if (nFacesAbove || nEdgesAbove)
146 {
147 // We've changed internal values so give
148 // boundary conditions opportunity to correct
149 U.correctBoundaryConditions();
150 }
151}
152
153
154// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
labelList faces_
Set of faces to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
bool useSubMesh() const noexcept
True if sub-selection should be used.
faceSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh, const word &defaultAreaName=word())
Construct from components.
Limits the maximum velocity magnitude to the specified max value.
word UName_
Name of operand velocity field.
limitVelocity(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh, const word &defaultAreaName=word())
Construct from components.
virtual void correct(areaVectorField &U)
Correct the velocity field.
virtual bool read(const dictionary &dict)
Read dictionary.
scalar max_
Maximum velocity magnitude.
Base abstract class for handling finite area options (i.e. faOption).
Definition faOption.H:149
const fvMesh & mesh() const noexcept
Return const access to the volume mesh.
Definition faOption.H:385
List< bool > applied_
Applied flag list - corresponds to each fieldNames_ entry.
Definition faOption.H:195
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition faOption.H:190
dictionary coeffs_
Dictionary containing source coefficients.
Definition faOption.H:185
const word & name() const noexcept
The source name.
Definition faOption.H:380
const word name_
Source name.
Definition faOption.H:165
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
Namespace for finite-area.
Definition limitHeight.C:30
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:40
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 POSIX.C:801
faPatchField< vector > faPatchVectorField
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< vector, faPatchField, areaMesh > areaVectorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dictionary dict
volScalarField & e