isotropicDamping.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) 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 "isotropicDamping.H"
27 #include "fvMatrix.H"
28 #include "fvmSup.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fv
36 {
37  defineTypeNameAndDebug(isotropicDamping, 0);
38  addToRunTimeSelectionTable(option, isotropicDamping, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::fv::isotropicDamping::add
46 (
47  const volScalarField::Internal& forceCoeff,
48  fvMatrix<vector>& eqn
49 )
50 {
51  eqn -= fvm::Sp(forceCoeff, eqn.psi());
52  eqn += forceCoeff*value_;
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const word& name,
61  const word& modelType,
62  const dictionary& dict,
63  const fvMesh& mesh
64 )
65 :
66  damping(name, modelType, dict, mesh),
67  value_("value", dimVelocity, coeffs_.lookup("value"))
68 {
69  read(dict);
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 (
77  fvMatrix<vector>& eqn,
78  const label fieldi
79 )
80 {
81  add(this->forceCoeff(), eqn);
82 }
83 
84 
86 (
87  const volScalarField& rho,
88  fvMatrix<vector>& eqn,
89  const label fieldi
90 )
91 {
92  add(rho*forceCoeff(), eqn);
93 }
94 
95 
97 (
98  const volScalarField& alpha,
99  const volScalarField& rho,
100  fvMatrix<vector>& eqn,
101  const label fieldi
102 )
103 {
104  add(alpha()*rho()*this->forceCoeff(), eqn);
105 }
106 
107 
109 {
110  if (damping::read(dict))
111  {
112  value_ =
114  (
115  value_.name(),
116  value_.dimensions(),
117  coeffs_.lookup(value_.name())
118  );
119 
120  return true;
121  }
122  else
123  {
124  return false;
125  }
126 }
127 
128 
129 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: damping.C:116
isotropicDamping(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
virtual bool read(const dictionary &dict)
Read dictionary.
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
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Source term to momentum equation.
Base fvOption for damping functions.
Definition: damping.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.
const dimensionSet dimVelocity