damping.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2019 OpenFOAM Foundation
6  \\/ M anipulation |
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 
26 #include "damping.H"
27 #include "fvMatrix.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36  defineTypeNameAndDebug(damping, 0);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
42 
44 {
46  (
48  (
49  IOobject
50  (
51  type() + ":forceCoeff",
52  mesh_.time().timeName(),
53  mesh_
54  ),
55  mesh_,
56  dimensionedScalar(lambda_.dimensions(), scale_.valid() ? 0 : 1)
57  )
58  );
59  scalarField& forceCoeff = tforceCoeff.ref();
60 
61  forAll(origins_, i)
62  {
63  const vectorField& c = mesh_.cellCentres();
64  const scalarField x((c - origins_[i]) & directions_[i]);
65  forceCoeff = max(forceCoeff, scale_->value(x));
66  }
67 
68  forceCoeff *= lambda_.value();
69 
70  // Write out the force coefficient for debugging
71  if (debug && mesh_.time().writeTime())
72  {
73  volScalarField vForceCoeff
74  (
75  IOobject
76  (
77  type() + ":forceCoeff",
78  mesh_.time().timeName(),
79  mesh_
80  ),
81  mesh_,
82  lambda_.dimensions(),
84  );
85  vForceCoeff.primitiveFieldRef() = forceCoeff;
86  vForceCoeff.correctBoundaryConditions();
87  vForceCoeff.write();
88  }
89 
90  return tforceCoeff;
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
97 (
98  const word& name,
99  const word& modelType,
100  const dictionary& dict,
101  const fvMesh& mesh
102 )
103 :
104  cellSetOption(name, modelType, dict, mesh),
105  lambda_("lambda", dimless/dimTime, coeffs_.lookup("lambda")),
106  scale_(),
107  origins_(),
108  directions_()
109 {
110  read(dict);
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (cellSetOption::read(dict))
119  {
120  lambda_ =
122  (
123  lambda_.name(),
124  lambda_.dimensions(),
125  coeffs_.lookup(lambda_.name())
126  );
127 
128  const bool foundScale = coeffs_.found("scale");
129  const bool foundOgn = coeffs_.found("origin");
130  const bool foundDir = coeffs_.found("direction");
131  const bool foundOgns = coeffs_.found("origins");
132  const bool foundDirs = coeffs_.found("directions");
133  const bool foundAll =
134  foundScale
135  && (
136  (foundOgn && foundDir && !foundOgns && !foundDirs)
137  || (!foundOgn && !foundDir && foundOgns && foundDirs)
138  );
139  const bool foundAny =
140  foundScale || foundOgn || foundDir || foundOgns || foundDirs;
141 
142  if (!foundAll)
143  {
144  scale_ = autoPtr<Function1<scalar>>();
145  origins_.clear();
146  directions_.clear();
147  }
148 
149  if (foundAll)
150  {
151  scale_ = Function1<scalar>::New("scale", coeffs_);
152  if (foundOgn)
153  {
154  origins_.setSize(1);
155  directions_.setSize(1);
156  coeffs_.lookup("origin") >> origins_.last();
157  coeffs_.lookup("direction") >> directions_.last();
158  }
159  else
160  {
161  coeffs_.lookup("origins") >> origins_;
162  coeffs_.lookup("directions") >> directions_;
163 
164  if
165  (
166  origins_.size() == 0
167  || directions_.size() == 0
168  || origins_.size() != directions_.size()
169  )
170  {
172  << "The same, non-zero number of origins and "
173  << "directions must be provided" << exit(FatalError);
174  }
175  }
176  forAll(directions_, i)
177  {
178  directions_[i] /= mag(directions_[i]);
179  }
180  }
181 
182  if (!foundAll && foundAny)
183  {
185  << "The scaling specification is incomplete. \"scale\", "
186  << "\"origin\" and \"direction\" (or \"origins\" and "
187  << "\"directions\"), must all be specified in order to scale "
188  << "the damping. The damping will be applied uniformly across "
189  << "the cell set." << endl << endl;
190  }
191 
192  fieldNames_ = wordList(1, coeffs_.lookupOrDefault<word>("U", "U"));
193 
194  applied_.setSize(1, false);
195 
196  return true;
197  }
198  else
199  {
200  return false;
201  }
202 }
203 
204 
205 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
damping(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Definition: damping.C:97
virtual bool read(const dictionary &dict)
Read source dictionary.
void clear()
Delete object (if the pointer is valid) and set pointer to.
Definition: autoPtrI.H:126
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: damping.C:116
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
List< word > wordList
A List of words.
Definition: fileName.H:54
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionedScalar c
Speed of light in a vacuum.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
Cell-set options abstract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:69
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
tmp< volScalarField::Internal > forceCoeff() const
Definition: damping.C:43
Namespace for OpenFOAM.