verticalDamping.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-2018 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 "verticalDamping.H"
27 #include "fvMesh.H"
28 #include "fvMatrix.H"
29 #include "geometricOneField.H"
30 #include "meshTools.H"
31 #include "Function1.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fv
41 {
42  defineTypeNameAndDebug(verticalDamping, 0);
43  addToRunTimeSelectionTable(option, verticalDamping, dictionary);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::fv::verticalDamping::add
51 (
52  const volVectorField& alphaRhoU,
53  fvMatrix<vector>& eqn,
54  const label fieldi
55 )
56 {
58  mesh_.lookupObject<uniformDimensionedVectorField>("g");
59 
60  const dimensionedSymmTensor lgg(lambda_*sqr(g)/magSqr(g));
61 
62  const DimensionedField<scalar, volMesh>& V = mesh_.V();
63 
64  // Calculate the scale
65  scalarField s(mesh_.nCells(), ramp_.valid() ? 0 : 1);
66  forAll(origins_, i)
67  {
68  const vectorField& c = mesh_.cellCentres();
69  const scalarField x((c - origins_[i]) & directions_[i]);
70  s = max(s, ramp_->value(x));
71  }
72 
73  // Check dimensions
74  eqn.dimensions()
75  - V.dimensions()*(lgg.dimensions() & alphaRhoU.dimensions());
76 
77  // Calculate the force and apply it to the equation
78  vectorField force(cells_.size());
79  forAll(cells_, i)
80  {
81  const label c = cells_[i];
82  force[i] = V[c]*s[c]*(lgg.value() & alphaRhoU[c]);
83  }
84  meshTools::constrainDirection(mesh_, mesh_.solutionD(), force);
85  forAll(cells_, i)
86  {
87  const label c = cells_[i];
88  eqn.source()[c] += force[i];
89  }
90 
91  // Write out the force coefficient for debugging
92  if (debug && mesh_.time().writeTime())
93  {
94  volScalarField forceCoeff
95  (
96  IOobject
97  (
98  type() + ":forceCoeff",
99  mesh_.time().timeName(),
100  mesh_
101  ),
102  mesh_,
103  lambda_*mag(g),
105  );
106  forceCoeff.primitiveFieldRef() *= s;
107  forceCoeff.write();
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  const word& name,
117  const word& modelType,
118  const dictionary& dict,
119  const fvMesh& mesh
120 )
121 :
122  cellSetOption(name, modelType, dict, mesh),
123  lambda_("lambda", dimless/dimTime, coeffs_.lookup("lambda")),
124  ramp_(),
125  origins_(),
126  directions_()
127 {
128  read(dict);
129 }
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 (
136  fvMatrix<vector>& eqn,
137  const label fieldi
138 )
139 {
140  add(eqn.psi(), eqn, fieldi);
141 }
142 
143 
145 (
146  const volScalarField& rho,
147  fvMatrix<vector>& eqn,
148  const label fieldi
149 )
150 {
151  add(rho*eqn.psi(), eqn, fieldi);
152 }
153 
154 
156 (
157  const volScalarField& alpha,
158  const volScalarField& rho,
159  fvMatrix<vector>& eqn,
160  const label fieldi
161 )
162 {
163  add(alpha*rho*eqn.psi(), eqn, fieldi);
164 }
165 
166 
168 {
169  if (cellSetOption::read(dict))
170  {
171  lambda_ =
173  (
174  lambda_.name(),
175  lambda_.dimensions(),
176  coeffs_.lookup(lambda_.name())
177  );
178 
179  const bool foundRamp = coeffs_.found("ramp");
180  const bool foundOgn = coeffs_.found("origin");
181  const bool foundDir = coeffs_.found("direction");
182  const bool foundOgns = coeffs_.found("origins");
183  const bool foundDirs = coeffs_.found("directions");
184  const bool foundAll =
185  foundRamp
186  && (
187  (foundOgn && foundDir && !foundOgns && !foundDirs)
188  || (!foundOgn && !foundDir && foundOgns && foundDirs)
189  );
190  const bool foundAny =
191  foundRamp || foundOgn || foundDir || foundOgns || foundDirs;
192 
193  if (!foundAll)
194  {
195  ramp_ = autoPtr<Function1<scalar>>();
196  origins_.clear();
197  directions_.clear();
198  }
199 
200  if (foundAll)
201  {
202  ramp_ = Function1<scalar>::New("ramp", coeffs_);
203  if (foundOgn)
204  {
205  origins_.setSize(1);
206  directions_.setSize(1);
207  coeffs_.lookup("origin") >> origins_.last();
208  coeffs_.lookup("direction") >> directions_.last();
209  }
210  else
211  {
212  coeffs_.lookup("origins") >> origins_;
213  coeffs_.lookup("directions") >> directions_;
214 
215  if
216  (
217  origins_.size() == 0
218  || directions_.size() == 0
219  || origins_.size() != directions_.size()
220  )
221  {
223  << "The same, non-zero number of origins and "
224  << "directions must be provided" << exit(FatalError);
225  }
226  }
227  forAll(directions_, i)
228  {
229  directions_[i] /= mag(directions_[i]);
230  }
231  }
232 
233  if (!foundAll && foundAny)
234  {
236  << "The ramping specification is incomplete. \"ramp\", "
237  << "\"origin\" and \"direction\" (or \"origins\" and "
238  << "\"directions\"), must all be specified in order to ramp "
239  << "the damping. The damping will be applied uniformly across "
240  << "the cell set." << endl << endl;
241  }
242 
243  fieldNames_ = wordList(1, coeffs_.lookupOrDefault<word>("U", "U"));
244 
245  applied_.setSize(1, false);
246 
247  return true;
248  }
249  else
250  {
251  return false;
252  }
253 }
254 
255 
256 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Source term to momentum equation.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
UniformDimensionedField< vector > uniformDimensionedVectorField
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
static const char *const typeName
Definition: Field.H:94
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:284
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
virtual bool read(const dictionary &dict)
Read source dictionary.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
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
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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:72
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
virtual bool read(const dictionary &dict)
Read dictionary.
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
verticalDamping(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
#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.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Cell-set options abtract 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
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:671
static autoPtr< Function1< Type > > New(const word &entryName, const dictionary &dict)
Selector.
Definition: Function1New.C:32
Namespace for OpenFOAM.