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 | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvMatrix.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37  defineTypeNameAndDebug(velocityDampingConstraint, 0);
39  (
40  option,
41  velocityDampingConstraint,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
51 {
52  // Note: we want to add
53  // deltaU/deltaT
54  // where deltaT is a local time scale:
55  // U/(cbrt of volume)
56  // Since directly manipulating the diagonal we multiply by volume.
57 
58  const scalarField& vol = mesh_.V();
59  const volVectorField& U = eqn.psi();
60  scalarField& diag = eqn.diag();
61 
62  label nDamped = 0;
63 
64  forAll(U, cellI)
65  {
66  scalar magU = mag(U[cellI]);
67  if (magU > UMax_)
68  {
69  scalar scale = sqr(Foam::cbrt(vol[cellI]));
70 
71  diag[cellI] += scale*(magU-UMax_);
72 
73  nDamped++;
74  }
75  }
76 
77  reduce(nDamped, sumOp<label>());
78 
79  Info<< type() << " " << name_ << " damped "
80  << nDamped << " ("
81  << 100*scalar(nDamped)/mesh_.globalData().nTotalCells()
82  << "%) of cells" << endl;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const word& name,
91  const word& modelType,
92  const dictionary& dict,
93  const fvMesh& mesh
94 )
95 :
96  cellSetOption(name, modelType, dict, mesh)
97 {
98  read(dict);
99 }
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 (
106  fvMatrix<vector>& eqn,
107  const label fieldI
108 )
109 {
110  addDamping(eqn);
111 }
112 
113 
115 {
116  os << indent << name_ << endl;
117  dict_.write(os);
118 }
119 
120 
122 {
124  {
125  UMax_ = readScalar(coeffs_.lookup("UMax"));
126 
127  if (coeffs_.found("UNames"))
128  {
129  coeffs_.lookup("UNames") >> fieldNames_;
130  }
131  else if (coeffs_.found("UName"))
132  {
133  word UName(coeffs_.lookup("UName"));
134  fieldNames_ = wordList(1, UName);
135  }
136  else
137  {
138  fieldNames_ = wordList(1, "U");
139  }
140 
141  applied_.setSize(fieldNames_.size(), false);
142 
143  return true;
144  }
145  else
146  {
147  return false;
148  }
149 }
150 
151 
152 // ************************************************************************* //
volFields.H
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fv::cellSetOption
Cell-set options abtract base class. Provides a base set of controls, e.g.
Definition: cellSetOption.H:71
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::diag
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
Definition: pointPatchFieldFunctions.H:262
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::globalMeshData::nTotalCells
label nTotalCells() const
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:394
Foam::fv::velocityDampingConstraint::velocityDampingConstraint
velocityDampingConstraint(const velocityDampingConstraint &)
Disallow default bitwise copy construct.
fvMatrix.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::fv::option::name_
const word name_
Source name.
Definition: fvOption.H:72
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
velocityDampingConstraint.H
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:78
Foam::fv::velocityDampingConstraint::addDamping
void addDamping(fvMatrix< vector > &eqn)
Definition: velocityDampingConstraint.C:50
U
U
Definition: pEqn.H:46
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:54
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fv::velocityDampingConstraint::writeData
virtual void writeData(Ostream &) const
Write data.
Definition: velocityDampingConstraint.C:114
Foam::fv::velocityDampingConstraint::UMax_
scalar UMax_
Maximum velocity.
Definition: velocityDampingConstraint.H:80
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::Info
messageStream Info
Foam::fv::cellSetOption::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: cellSetOptionIO.C:30
Foam::fv::velocityDampingConstraint::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: velocityDampingConstraint.C:121
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
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:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sumOp
Definition: ops.H:162
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::fv::velocityDampingConstraint::constrain
virtual void constrain(fvMatrix< vector > &eqn, const label fieldI)
Constrain vector matrix.
Definition: velocityDampingConstraint.C:105
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1143
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47