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"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace fv
39{
42 (
43 option,
46 );
47}
48}
49
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
54{
55 // Note: we want to add
56 // deltaU/deltaT
57 // where deltaT is a local time scale:
58 // U/(cbrt of volume)
59 // Since directly manipulating the diagonal we multiply by volume.
60
61 const scalarField& vol = mesh_.V();
62 const volVectorField& U = eqn.psi();
63 scalarField& diag = eqn.diag();
64
65 // Count nTotCells ourselves
66 // (maybe only applying on a subset)
67 label nDamped(0);
68 const label nTotCells(returnReduce(cells_.size(), sumOp<label>()));
69
70 for (label celli : cells_)
71 {
72 const scalar magU = mag(U[celli]);
73 if (magU > UMax_)
74 {
75 const scalar scale = sqr(Foam::cbrt(vol[celli]));
76
77 diag[celli] += scale*(magU-UMax_);
78
79 ++nDamped;
80 }
81 }
82
83 reduce(nDamped, sumOp<label>());
84
85 // Percent, max 2 decimal places
86 const auto percent = [](scalar num, label denom) -> scalar
87 {
88 return (denom ? 1e-2*round(1e4*num/denom) : 0);
89 };
90
91 Info<< type() << ' ' << name_ << " damped "
92 << nDamped << " ("
93 << percent(nDamped, nTotCells)
94 << "%) of cells, with max limit " << UMax_ << endl;
95}
96
97
98// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99
101(
102 const word& name,
103 const word& modelType,
104 const dictionary& dict,
105 const fvMesh& mesh
106)
107:
108 fv::cellSetOption(name, modelType, dict, mesh)
109{
110 read(dict);
111}
112
113
114// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115
117(
118 fvMatrix<vector>& eqn,
119 const label fieldi
120)
121{
122 addDamping(eqn);
123}
124
125
127{
128 dict_.writeEntry(name_, os);
129}
130
131
133{
135 {
136 coeffs_.readEntry("UMax", UMax_);
137
138 if (!coeffs_.readIfPresent("UNames", fieldNames_))
139 {
140 fieldNames_.resize(1);
141 fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
142 }
143
145
146 return true;
147 }
148
149 return false;
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.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
virtual bool read()
Re-read model coefficients if they have changed.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:412
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
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.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:127
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
const word name_
Source name.
Definition: fvOption.H:133
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.
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.
scalarField & diag()
Definition: lduMatrix.C:192
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
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
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
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar cbrt(const dimensionedScalar &ds)
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.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
labelList fv(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11