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-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 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"
69 
70 namespace Foam
71 {
72 
73 /*---------------------------------------------------------------------------*\
74  Class wallDampingModel Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 class wallDampingModel
78 :
79  public wallDependentModel
80 {
81 protected:
82 
83  // Protected data
84 
85  //- Interface
86  const dispersedPhaseInterface interface_;
87 
88  //- Diameter coefficient
89  const dimensionedScalar Cd_;
90 
91  //- Distance from wall below which the field is damped to zero
93 
94  //- Set the value to zero in wall-adjacent cells
95  const Switch zeroInNearWallCells_;
96 
97 
98  // Protected member functions
99 
100  //- Return the force limiter field
101  virtual tmp<volScalarField> limiter() const = 0;
102 
103 
104 public:
105 
106  //- Runtime type information
107  TypeName("wallDampingModel");
108 
109 
110  // Declare runtime construction
111 
113  (
116  dictionary,
117  (
118  const dictionary& dict,
119  const phaseInterface& interface
120  ),
121  (dict, interface)
122  );
123 
124 
125  // Constructors
126 
127  //- Construct from a dictionary and an interface
129  (
130  const dictionary& dict,
131  const phaseInterface& interface
132  );
133 
134 
135  //- Destructor
136  virtual ~wallDampingModel();
137 
138 
139  // Selectors
140 
142  (
143  const dictionary& dict,
144  const phaseInterface& interface
145  );
146 
147 
148  // Member Functions
149 
150  //- Return damped coefficient
151  virtual tmp<volScalarField> damping() const;
152 
153  //- Return damped face coefficient
154  virtual tmp<surfaceScalarField> dampingf() const;
155 };
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 } // End namespace Foam
161 
162 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
163 
164 #endif
165 
166 // ************************************************************************* //
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
Wall damping models can be used to filter interfacial models near the walls. This is particularly use...
declareRunTimeSelectionTable(autoPtr, wallDampingModel, dictionary,(const dictionary &dict, const phaseInterface &interface),(dict, interface))
const dimensionedScalar zeroWallDist_
Distance from wall below which the field is damped to zero.
const dispersedPhaseInterface interface_
Interface.
const dimensionedScalar Cd_
Diameter coefficient.
virtual tmp< surfaceScalarField > dampingf() const
Return damped face coefficient.
virtual ~wallDampingModel()
Destructor.
const Switch zeroInNearWallCells_
Set the value to zero in wall-adjacent cells.
static autoPtr< wallDampingModel > New(const dictionary &dict, const phaseInterface &interface)
TypeName("wallDampingModel")
Runtime type information.
virtual tmp< volScalarField > limiter() const =0
Return the force limiter field.
virtual tmp< volScalarField > damping() const
Return damped coefficient.
wallDampingModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Macros to ease declaration of run-time selection tables.
dictionary dict