wallDampingModel.H
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) 2015-2020 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 Class
25  Foam::wallDampingModel
26 
27 Description
28  Wall damping models can be used to filter interfacial models near the walls.
29  This is particularly useful for the lift force because of its dependence on
30  the velocity gradient.
31 
32  All damping functions accept the following parameters:
33 
34  - Cd: A coefficient for filtering the distance from the wall based on the
35  dispersed phase diameter. This can be useful to correct gradient
36  sampling error when the dispersed phase diameter is significantly
37  larger than near wall mesh resolution.
38  - zeroWallDist: A constant offset from the wall for the zero point of
39  the damping function. Below this distance, the damping will reduce the
40  value to zero.
41  - zeroInNearWallCells: A switch which sets the value to zero in near wall
42  cells regardless of the other parameters. This is recommended to be set
43  if a lift force is applied together with turbulent wall functions.
44 
45 Usage
46  \table
47  Property | Description | Required | Default value
48  Cd | Diameter coefficient | yes | none
49  zeroWallDist | Offset from wall | no | 0
50  zeroInNearWallCells | Zero near wall cells | no | no
51  \endtable
52 
53 SourceFiles
54  wallDampingModel.C
55  wallDampingModelNew.C
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #ifndef wallDampingModel_H
60 #define wallDampingModel_H
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 #include "wallDependentModel.H"
65 #include "volFields.H"
66 #include "dictionary.H"
67 #include "runTimeSelectionTables.H"
68 
69 namespace Foam
70 {
71 
72 class phasePair;
73 
74 /*---------------------------------------------------------------------------*\
75  Class wallDampingModel Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class wallDampingModel
79 :
80  public wallDependentModel
81 {
82 protected:
83 
84  // Protected data
85 
86  //- Phase pair
87  const phasePair& pair_;
88 
89  //- Diameter coefficient
90  const dimensionedScalar Cd_;
91 
92  //- Distance from wall below which the field is damped to zero
94 
95  //- Set the value to zero in wall-adjacent cells
96  const Switch zeroInNearWallCells_;
97 
98 
99  // Protected member functions
100 
101  //- Return the force limiter field
102  virtual tmp<volScalarField> limiter() const = 0;
103 
104 
105 public:
107  //- Runtime type information
108  TypeName("wallDampingModel");
110 
111  // Declare runtime construction
114  (
117  dictionary,
118  (
119  const dictionary& dict,
120  const phasePair& pair
121  ),
122  (dict, pair)
123  );
124 
125 
126  // Static Data Members
127 
128  //- Coefficient dimensions
129  static const dimensionSet dimF;
130 
131 
132  // Constructors
133 
134  //- Construct from components
136  (
137  const dictionary& dict,
138  const phasePair& pair
139  );
140 
141 
142  //- Destructor
143  virtual ~wallDampingModel();
144 
145 
146  // Selectors
147 
149  (
150  const dictionary& dict,
151  const phasePair& pair
152  );
153 
154 
155  // Member Functions
156 
157  //- Return damped coefficient
158  virtual tmp<volScalarField> damping() const;
159 
160  //- Return damped face coefficient
161  virtual tmp<surfaceScalarField> dampingf() const;
162 };
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 } // End namespace Foam
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 #endif
172 
173 // ************************************************************************* //
declareRunTimeSelectionTable(autoPtr, wallDampingModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair))
Wall damping models can be used to filter interfacial models near the walls. This is particularly use...
const dimensionedScalar Cd_
Diameter coefficient.
dictionary dict
virtual tmp< volScalarField > limiter() const =0
Return the force limiter field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual tmp< surfaceScalarField > dampingf() const
Return damped face coefficient.
const dimensionedScalar zeroWallDist_
Distance from wall below which the field is damped to zero.
static const dimensionSet dimF
Coefficient dimensions.
const Switch zeroInNearWallCells_
Set the value to zero in wall-adjacent cells.
Dimension set for the base types.
Definition: dimensionSet.H:120
virtual ~wallDampingModel()
Destructor.
TypeName("wallDampingModel")
Runtime type information.
static autoPtr< wallDampingModel > New(const dictionary &dict, const phasePair &pair)
wallDampingModel(const dictionary &dict, const phasePair &pair)
Construct from components.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual tmp< volScalarField > damping() const
Return damped coefficient.
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: PtrList.H:53
const phasePair & pair_
Phase pair.
Namespace for OpenFOAM.