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 -------------------------------------------------------------------------------
11 License
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 #include "cellCellStencilObject.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fv
40 {
41  defineTypeNameAndDebug(velocityDampingConstraint, 0);
43  (
44  option,
45  velocityDampingConstraint,
46  dictionary
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()
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;
119 
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)
138 {
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 {
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 // ************************************************************************* //
void writeFileHeader(Ostream &os)
Write file header information.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual bool canWriteToFile() const
Flag to allow writing to the file.
Definition: writeFile.C:287
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
Tab [isspace].
Definition: token.H:129
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool read(const dictionary &dict)
Read source dictionary.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Macros for easy insertion into run-time selection tables.
velocityDampingConstraint(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const word name_
Source name.
Definition: fvOption.H:132
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
dimensionedScalar cbrt(const dimensionedScalar &ds)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Constrain vector matrix.
labelList cells_
Set of cells to apply source to.
U
Definition: pEqn.H:72
scalar timeOutputValue() const
Return the current user-time value. (ie, after applying any timeToUserTime() conversion) ...
Definition: TimeStateI.H:24
void addDamping(fvMatrix< vector > &eqn)
Constrain the given velocity fields by a given maximum value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Intermediate abstract class for handling cell-set options for the derived fvOptions.
scalar UMax_
Maximum velocity magnitude.
scalarField & diag()
Definition: lduMatrix.C:197
virtual void writeData(Ostream &os) const
Write data.
Namespace for OpenFOAM.