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-2024 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 "fvMatrix.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::verticalDamping::readCoeffs(const dictionary& dict)
46 {
48 }
49 
50 
51 void Foam::fv::verticalDamping::add
52 (
53  const volVectorField& alphaRhoU,
54  fvMatrix<vector>& eqn
55 ) const
56 {
59 
60  const dimensionedSymmTensor gg(sqr(g)/magSqr(g));
61 
62  eqn -= forceCoeff()*(gg & alphaRhoU());
63 }
64 
65 
66 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
67 
69 (
70  const word& name,
71  const word& modelType,
72  const fvMesh& mesh,
73  const dictionary& dict
74 )
75 :
76  forcing(name, modelType, mesh, dict),
77  UName_(dict.lookupOrDefault<word>("U", "U"))
78 {
79  readCoeffs(coeffs(dict));
81 }
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
87 {
88  return wordList(1, UName_);
89 }
90 
91 
93 (
94  const volVectorField& U,
95  fvMatrix<vector>& eqn
96 ) const
97 {
98  add(U, eqn);
99 }
100 
101 
103 (
104  const volScalarField& rho,
105  const volVectorField& U,
106  fvMatrix<vector>& eqn
107 ) const
108 {
109  add(rho*U, eqn);
110 }
111 
112 
114 (
115  const volScalarField& alpha,
116  const volScalarField& rho,
117  const volVectorField& U,
118  fvMatrix<vector>& eqn
119 ) const
120 {
121  add(alpha*rho*U, eqn);
122 }
123 
124 
126 {
127  return true;
128 }
129 
130 
132 {}
133 
134 
136 {}
137 
138 
140 {}
141 
142 
144 {
145  if (forcing::read(dict))
146  {
147  readCoeffs(coeffs(dict));
148  return true;
149  }
150  else
151  {
152  return false;
153  }
154 }
155 
156 
157 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:96
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
Base fvModel for forcing functions.
Definition: forcing.H:59
void writeForceFields() const
Optionally write the forcing fields:
Definition: forcing.C:214
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: forcing.C:251
void readLambda(const dictionary &dict)
Read the forcing coefficients.
Definition: forcing.C:45
This fvModel applies an explicit forcing force to components of the vector field in the direction of ...
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.
verticalDamping(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
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 mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
UniformDimensionedField< vector > uniformDimensionedVectorField
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void add(LagrangianPatchField< typename typeOfSum< Type1, Type2 >::type > &f, const LagrangianPatchField< Type1 > &f1, const LagrangianPatchField< Type2 > &f2)
labelList fv(nPoints)
dictionary dict