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) 2016-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2023 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 "limitVelocity.H"
30#include "volFields.H"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace fv
38{
41}
42}
43
44
45// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46
48{
50 writeCommented(os, "Time");
51 writeTabbed(os, "nDampedCells_[count]");
52 writeTabbed(os, "nDampedCells_[%]");
53 writeTabbed(os, "nDampedFaces_[count]");
54 writeTabbed(os, "nDampedFaces_[%]");
55
56 os << endl;
58 writtenHeader_ = true;
59}
60
61
62// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63
65(
66 const word& name,
67 const word& modelType,
68 const dictionary& dict,
69 const fvMesh& mesh
70)
71:
72 fv::cellSetOption(name, modelType, dict, mesh),
73 writeFile(mesh, name, typeName, dict, false),
74 UName_("U"),
75 max_(0)
77 read(dict);
78}
79
80
81// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82
84{
85 if (!(fv::cellSetOption::read(dict) && writeFile::read(dict)))
86 {
87 return false;
88 }
89
90 coeffs_.readEntry("max", max_);
91 coeffs_.readIfPresent("U", UName_);
92
93
94 fieldNames_.resize(1, UName_);
95
97
98
99 if (canResetFile())
100 {
101 resetFile(typeName);
102 }
103
104 if (canWriteHeader())
105 {
106 writeFileHeader(file());
108
109
110 return true;
111}
112
113
115{
116 const scalar maxSqrU = sqr(max_);
117
118 // Count nTotCells ourselves
119 // (maybe only applying on a subset)
120 label nCellsAbove(0);
121 const label nTotCells(returnReduce(cells_.size(), sumOp<label>()));
122
123 vectorField& Uif = U.primitiveFieldRef();
124
125 for (const label celli : cells_)
126 {
127 auto& Uval = Uif[celli];
128
129 const scalar magSqrUi = magSqr(Uval);
130
131 if (magSqrUi > maxSqrU)
132 {
133 Uval *= sqrt(maxSqrU/magSqrUi);
134 ++nCellsAbove;
135 }
136 }
137
138 // Handle boundaries in the case of 'all'
139
140 label nFacesAbove(0);
141 label nTotFaces(0);
142
144 {
145 for (fvPatchVectorField& Up : U.boundaryFieldRef())
146 {
147 if (!Up.fixesValue())
148 {
149 // Do not count patches that fix velocity themselves
150 nTotFaces += Up.size();
151
152 for (auto& Uval : Up)
153 {
154 const scalar magSqrUi = magSqr(Uval);
155
156 if (magSqrUi > maxSqrU)
157 {
158 Uval *= sqrt(maxSqrU/magSqrUi);
159 ++nFacesAbove;
160 }
161 }
162 }
163 }
164 }
165
166 // Percent, max 2 decimal places
167 const auto percent = [](scalar num, label denom) -> scalar
168 {
169 return (denom ? 1e-2*round(1e4*num/denom) : 0);
170 };
171
172
173 reduce(nCellsAbove, sumOp<label>());
174
175 const scalar nCellsAbovePercent = percent(nCellsAbove, nTotCells);
176
177 // Report total numbers and percent
178 Info<< type() << ' ' << name_ << " Limited ";
179
180 Info<< nCellsAbove << " ("
181 << nCellsAbovePercent
182 << "%) of cells";
183
184 reduce(nTotFaces, sumOp<label>());
185 reduce(nFacesAbove, sumOp<label>());
186 scalar nFacesAbovePercent(0);
187 if (nTotFaces)
188 {
189 nFacesAbovePercent = percent(nFacesAbove, nTotFaces);
190
191 Info<< ", " << nFacesAbove << " ("
192 << nFacesAbovePercent
193 << "%) of faces";
194 }
195 Info<< ", with max limit " << max_ << endl;
196
197 if (nCellsAbove || nFacesAbove)
198 {
199 // We've changed internal values so give
200 // boundary conditions opportunity to correct
201 U.correctBoundaryConditions();
202 }
203
204
205 if (canWriteToFile())
206 {
207 file()
208 << mesh_.time().timeOutputValue() << token::TAB
209 << nCellsAbove << token::TAB
210 << nCellsAbovePercent << token::TAB
211 << nFacesAbove << token::TAB
212 << nFacesAbovePercent
213 << endl;
214 }
215}
216
217
218// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
virtual bool canWriteToFile() const
Flag to allow writing to the file.
Definition writeFile.C:292
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
bool writtenHeader_
Flag to identify whether the header has been written.
Definition writeFile.H:157
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void resetFile(const word &name)
Reset internal file pointer to new file with new name.
Definition writeFile.C:165
virtual bool canWriteHeader() const
Flag to allow writing the header.
Definition writeFile.C:304
virtual bool canResetFile() const
Flag to allow resetting the file.
Definition writeFile.C:298
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Intermediate abstract class for handling cell-set options for the derived fvOptions.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
bool useSubMesh() const noexcept
True if sub-selection should be used.
Corrects velocity field (i.e. T) within a specified region by applying a given maximum velocity magni...
word UName_
Name of operand velocity field.
virtual bool read(const dictionary &dict)
Read dictionary.
scalar max_
Maximum velocity magnitude.
limitVelocity(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
void writeFileHeader(Ostream &os)
Write file header information.
virtual void correct(volVectorField &U)
Correct the velocity field.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
@ TAB
Tab [isspace].
Definition token.H:142
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
OBJstream os(runTime.globalPath()/outputName)
Namespace for finite-volume.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
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
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.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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)
fvPatchField< vector > fvPatchVectorField
dictionary dict
volScalarField & e