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  value_ =
49  (
50  value_.name(),
51  value_.dimensions(),
52  coeffs().lookup(value_.name())
53  );
54 }
55 
56 
57 void Foam::fv::isotropicDamping::add
58 (
59  const volScalarField::Internal& forceCoeff,
60  fvMatrix<vector>& eqn
61 ) const
62 {
63  eqn -= fvm::Sp(forceCoeff, eqn.psi());
64  eqn += forceCoeff*value_;
65 }
66 
67 
68 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69 
71 (
72  const word& name,
73  const word& modelType,
74  const fvMesh& mesh,
75  const dictionary& dict
76 )
77 :
78  forcing(name, modelType, mesh, dict),
79  UName_(coeffs().lookupOrDefault<word>("U", "U")),
80  value_("value", dimVelocity, vector::uniform(NaN))
81 {
82  readCoeffs();
83 }
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
90  return wordList(1, UName_);
91 }
92 
93 
95 (
96  fvMatrix<vector>& eqn,
97  const word& fieldName
98 ) const
99 {
100  add(this->forceCoeff(), eqn);
101 }
102 
103 
105 (
106  const volScalarField& rho,
107  fvMatrix<vector>& eqn,
108  const word& fieldName
109 ) const
110 {
111  add(rho*forceCoeff(), eqn);
112 }
113 
114 
116 (
117  const volScalarField& alpha,
118  const volScalarField& rho,
119  fvMatrix<vector>& eqn,
120  const word& fieldName
121 ) const
122 {
123  add(alpha()*rho()*this->forceCoeff(), eqn);
124 }
125 
126 
128 {
129  return true;
130 }
131 
132 
134 {}
135 
136 
138 {}
139 
140 
142 {}
143 
144 
146 {
147  if (forcing::read(dict))
148  {
149  readCoeffs();
150  return true;
151  }
152  else
153  {
154  return false;
155  }
156 }
157 
158 
159 // ************************************************************************* //
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:160
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:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
Base fvModel for forcing functions.
Definition: forcing.H:59
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: forcing.C:223
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 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.
virtual void addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
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.
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 > &)
static Type NaN()
Return a primitive with all components set to NaN.
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
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
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict