verticalDamping.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) 2017 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"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace fv
39 {
40  defineTypeNameAndDebug(verticalDamping, 0);
41  addToRunTimeSelectionTable(option, verticalDamping, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::fv::verticalDamping::add
49 (
50  const volVectorField& alphaRhoU,
51  fvMatrix<vector>& eqn,
52  const label fieldi
53 )
54 {
56  mesh_.lookupObject<uniformDimensionedVectorField>("g");
57 
58  const dimensionedSymmTensor lgg(lambda_*sqr(g)/magSqr(g));
59 
60  const DimensionedField<scalar, volMesh>& V = mesh_.V();
61 
62  // Check dimensions
63  eqn.dimensions()
64  - V.dimensions()*(lgg.dimensions() & alphaRhoU.dimensions());
65 
66  // Calculate the force and apply it to the equation
67  vectorField force(cells_.size());
68  forAll(cells_, i)
69  {
70  const label c = cells_[i];
71  force[i] = V[c]*(lgg.value() & alphaRhoU[c]);
72  }
73  meshTools::constrainDirection(mesh_, mesh_.solutionD(), force);
74  forAll(cells_, i)
75  {
76  const label c = cells_[i];
77  eqn.source()[c] += force[i];
78  }
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
85 (
86  const word& name,
87  const word& modelType,
88  const dictionary& dict,
89  const fvMesh& mesh
90 )
91 :
92  cellSetOption(name, modelType, dict, mesh),
93  lambda_("lambda", dimless/dimTime, coeffs_.lookup("lambda"))
94 {
95  read(dict);
96 }
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 (
103  fvMatrix<vector>& eqn,
104  const label fieldi
105 )
106 {
107  add(eqn.psi(), eqn, fieldi);
108 }
109 
110 
112 (
113  const volScalarField& rho,
114  fvMatrix<vector>& eqn,
115  const label fieldi
116 )
117 {
118  add(rho*eqn.psi(), eqn, fieldi);
119 }
120 
121 
123 (
124  const volScalarField& alpha,
125  const volScalarField& rho,
126  fvMatrix<vector>& eqn,
127  const label fieldi
128 )
129 {
130  add(alpha*rho*eqn.psi(), eqn, fieldi);
131 }
132 
133 
135 {
136  if (cellSetOption::read(dict))
137  {
138  lambda_ =
140  (
141  lambda_.name(),
142  lambda_.dimensions(),
143  coeffs_.lookup(lambda_.name())
144  );
145 
146  fieldNames_ = wordList(1, coeffs_.lookupOrDefault<word>("U", "U"));
147 
148  applied_.setSize(1, false);
149 
150  return true;
151  }
152  else
153  {
154  return false;
155  }
156 }
157 
158 
159 // ************************************************************************* //
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
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:284
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
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 > &)
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.
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.
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
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
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:671
Namespace for OpenFOAM.