Loading...
Searching...
No Matches
velocityDampingConstraint.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) 2015 OpenFOAM Foundation
9 Copyright (C) 2015-2022 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
31#include "fvMatrix.H"
32#include "volFields.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace fv
40{
43 (
44 option,
47 );
48}
49}
50
51
52// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53
55{
56 // Note: we want to add
57 // deltaU/deltaT
58 // where deltaT is a local time scale:
59 // U/(cbrt of volume)
60 // Since directly manipulating the diagonal we multiply by volume.
61
62 const scalarField& vol = mesh_.V();
63 const volVectorField& U = eqn.psi();
64 scalarField& diag = eqn.diag();
65
66 // Count nTotCells ourselves
67 // (maybe only applying on a subset)
68 label nDamped(0);
69 const label nTotCells(returnReduce(cells_.size(), sumOp<label>()));
70
71 for (const label celli : cells_)
72 {
73 const scalar magU = mag(U[celli]);
74 if (magU > UMax_)
75 {
76 const scalar scale = sqr(Foam::cbrt(vol[celli]));
77
78 diag[celli] += C_*scale*(magU-UMax_);
79
80 ++nDamped;
81 }
82 }
83
84 reduce(nDamped, sumOp<label>());
85
86 // Percent, max 2 decimal places
87 const auto percent = [](scalar num, label denom) -> scalar
88 {
89 return (denom ? 1e-2*round(1e4*num/denom) : 0);
90 };
91
92 const scalar nDampedPercent = percent(nDamped, nTotCells);
93
94 Info<< type() << ' ' << name_ << " damped "
95 << nDamped << " ("
96 << nDampedPercent
97 << "%) of cells, with max limit " << UMax_ << endl;
98
99
100 if (canWriteToFile())
101 {
102 file()
103 << mesh_.time().timeOutputValue() << token::TAB
104 << nDamped << token::TAB
105 << nDampedPercent
106 << endl;
107 }
108}
109
110
112{
113 writeHeaderValue(os, "UMax", Foam::name(UMax_));
114 writeCommented(os, "Time");
115 writeTabbed(os, "nDamped_[count]");
116 writeTabbed(os, "nDamped_[%]");
117
118 os << endl;
120 writtenHeader_ = true;
121}
122
123
124// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125
127(
128 const word& name,
129 const word& modelType,
130 const dictionary& dict,
131 const fvMesh& mesh
132)
133:
134 fv::cellSetOption(name, modelType, dict, mesh),
135 writeFile(mesh, name, typeName, dict, false),
136 UMax_(GREAT), // overwritten later
137 C_(1)
139 read(dict);
140}
141
142
143// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
144
146(
147 fvMatrix<vector>& eqn,
148 const label fieldi
149)
150{
151 addDamping(eqn);
152}
153
156{
157 dict_.writeEntry(name_, os);
158}
159
160
162{
163 if (!(fv::cellSetOption::read(dict) && writeFile::read(dict)))
164 {
165 return false;
166 }
167
168 coeffs_.readEntry("UMax", UMax_);
169 coeffs_.readIfPresent("C", C_);
170
171 if (!coeffs_.readIfPresent("UNames", fieldNames_))
172 {
173 fieldNames_.resize(1);
174 fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
175 }
176
178
179
180 if (canResetFile())
181 {
182 resetFile(typeName);
183 }
184
185 if (canWriteHeader())
186 {
187 writeFileHeader(file());
188 }
189
190
191 return true;
192}
193
194
195// ************************************************************************* //
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
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
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.
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 dict_
Top level source dictionary.
Definition fvOption.H:147
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
Constrain given velocity fields to dampen velocity fluctuations exceeding a given value within a spec...
scalar UMax_
Maximum velocity magnitude.
virtual void writeData(Ostream &os) const
Write data.
velocityDampingConstraint(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
void addDamping(fvMatrix< vector > &eqn)
Constrain the given velocity fields by a given maximum value.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Constrain vector matrix.
void writeFileHeader(Ostream &os)
Write file header information.
const scalarField & diag() const
Definition lduMatrix.C:195
@ 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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
dictionary dict
volScalarField & e