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-2023 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 {
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::fv::isotropicDamping::readCoeffs()
46 {
47  readLambda();
48 
49  value_ =
51  (
52  value_.name(),
53  value_.dimensions(),
54  coeffs().lookup(value_.name())
55  );
56 }
57 
58 
59 void Foam::fv::isotropicDamping::add
60 (
61  const volScalarField::Internal& forceCoeff,
62  fvMatrix<vector>& eqn
63 ) const
64 {
65  eqn -= fvm::Sp(forceCoeff, eqn.psi());
66  eqn += forceCoeff*value_;
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
73 (
74  const word& name,
75  const word& modelType,
76  const fvMesh& mesh,
77  const dictionary& dict
78 )
79 :
80  forcing(name, modelType, mesh, dict),
81  UName_(coeffs().lookupOrDefault<word>("U", "U")),
82  value_("value", dimVelocity, vector::uniform(NaN))
83 {
84  readCoeffs();
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
93  return wordList(1, UName_);
94 }
95 
96 
98 (
99  const volVectorField& U,
100  fvMatrix<vector>& eqn
101 ) const
102 {
103  add(this->forceCoeff(), eqn);
104 }
105 
106 
108 (
109  const volScalarField& rho,
110  const volVectorField& U,
111  fvMatrix<vector>& eqn
112 ) const
113 {
114  add(rho*this->forceCoeff(), eqn);
115 }
116 
117 
119 (
120  const volScalarField& alpha,
121  const volScalarField& rho,
122  const volVectorField& U,
123  fvMatrix<vector>& eqn
124 ) const
125 {
126  add(alpha()*rho()*this->forceCoeff(), eqn);
127 }
128 
129 
131 {
132  return true;
133 }
134 
135 
137 {}
138 
139 
141 {}
142 
143 
145 {}
146 
147 
149 {
150  if (forcing::read(dict))
151  {
152  readCoeffs();
153  return true;
154  }
155  else
156  {
157  return false;
158  }
159 }
160 
161 
162 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dimensionSet & dimensions() const
Return const reference to dimensions.
const word & name() const
Return const reference to name.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
Base fvModel for forcing functions.
Definition: forcing.H:59
void readLambda()
Read the forcing coefficients.
Definition: forcing.C:45
void writeForceFields() const
Optionally write the forcing fields:
Definition: forcing.C:214
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: forcing.C:251
This fvModel applies an implicit forcing force to all components of the vector field to relax the fie...
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
isotropicDamping(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Calculate the matrix for implicit and explicit sources.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
const dimensionSet dimVelocity
labelList fv(nPoints)
dictionary dict