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-2022 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 {
37  defineTypeNameAndDebug(verticalDamping, 0);
38  addToRunTimeSelectionTable(fvModel, verticalDamping, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void Foam::fv::verticalDamping::add
46 (
47  const volVectorField& alphaRhoU,
48  fvMatrix<vector>& eqn
49 ) const
50 {
52  mesh().lookupObject<uniformDimensionedVectorField>("g");
53 
54  const dimensionedSymmTensor gg(sqr(g)/magSqr(g));
55 
56  eqn -= forceCoeff()*(gg & alphaRhoU());
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 
63 (
64  const word& name,
65  const word& modelType,
66  const dictionary& dict,
67  const fvMesh& mesh
68 )
69 :
70  damping(name, modelType, dict, mesh)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 (
78  fvMatrix<vector>& eqn,
79  const word& fieldName
80 ) const
81 {
82  add(eqn.psi(), eqn);
83 }
84 
85 
87 (
88  const volScalarField& rho,
89  fvMatrix<vector>& eqn,
90  const word& fieldName
91 ) const
92 {
93  add(rho*eqn.psi(), eqn);
94 }
95 
96 
98 (
99  const volScalarField& alpha,
100  const volScalarField& rho,
101  fvMatrix<vector>& eqn,
102  const word& fieldName
103 ) const
104 {
105  add(alpha*rho*eqn.psi(), eqn);
106 }
107 
108 
110 {
111  return true;
112 }
113 
114 
116 {}
117 
118 
120 {}
121 
122 
124 {}
125 
126 
127 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
UniformDimensionedField< vector > uniformDimensionedVectorField
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void addSup(fvMatrix< vector > &eqn, const word &fieldName) const
Source term to momentum equation.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
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.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
verticalDamping(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Base fvModel for damping functions.
Definition: damping.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
virtual bool movePoints()
Update for mesh motion.
Namespace for OpenFOAM.