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-2021 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 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace fv
39 {
40  defineTypeNameAndDebug(velocityDampingConstraint, 0);
42  (
43  option,
44  velocityDampingConstraint,
45  dictionary
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  label nDamped = 0;
66 
67  forAll(U, cellI)
68  {
69  const scalar magU = mag(U[cellI]);
70  if (magU > UMax_)
71  {
72  const scalar scale = sqr(Foam::cbrt(vol[cellI]));
73 
74  diag[cellI] += scale*(magU-UMax_);
75 
76  ++nDamped;
77  }
78  }
79 
80  reduce(nDamped, sumOp<label>());
81 
82  Info<< type() << " " << name_ << " damped "
83  << nDamped << " ("
84  << 100*scalar(nDamped)/mesh_.globalData().nTotalCells()
85  << "%) of cells" << endl;
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
92 (
93  const word& name,
94  const word& modelType,
95  const dictionary& dict,
96  const fvMesh& mesh
97 )
98 :
99  fv::cellSetOption(name, modelType, dict, mesh)
100 {
101  read(dict);
102 }
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 (
109  fvMatrix<vector>& eqn,
110  const label fieldi
111 )
112 {
113  addDamping(eqn);
114 }
115 
116 
118 {
119  dict_.writeEntry(name_, os);
120 }
121 
122 
124 {
126  {
127  coeffs_.readEntry("UMax", UMax_);
128 
129  if (!coeffs_.readIfPresent("UNames", fieldNames_))
130  {
131  fieldNames_.resize(1);
132  fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
133  }
134 
136 
137  return true;
138  }
139 
140  return false;
141 }
142 
143 
144 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::cellSetOption
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Definition: cellSetOption.H:163
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:287
Foam::fv::option::resetApplied
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fv::option::name_
const word name_
Source name.
Definition: fvOption.H:133
velocityDampingConstraint.H
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Foam::sumOp
Definition: ops.H:213
Foam::fv::velocityDampingConstraint::addDamping
void addDamping(fvMatrix< vector > &eqn)
Constrain the given velocity fields by a given maximum value.
Definition: velocityDampingConstraint.C:53
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::Field< scalar >
Foam::fv::velocityDampingConstraint::UMax_
scalar UMax_
Maximum velocity magnitude.
Definition: velocityDampingConstraint.H:157
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::fv::velocityDampingConstraint::writeData
virtual void writeData(Ostream &os) const
Write data.
Definition: velocityDampingConstraint.C:117
Foam::fv::cellSetOption::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: cellSetOption.C:255
Foam::fv::velocityDampingConstraint::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: velocityDampingConstraint.C:123
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:399
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
fv
labelList fv(nPoints)
U
U
Definition: pEqn.H:72
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::fv::velocityDampingConstraint::constrain
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Constrain vector matrix.
Definition: velocityDampingConstraint.C:108
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::globalMeshData::nTotalCells
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:371
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fv::velocityDampingConstraint::velocityDampingConstraint
velocityDampingConstraint(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: velocityDampingConstraint.C:92
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179