alphaContactAngleFvPatchScalarField.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) 2011-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::alphaContactAngleFvPatchScalarField
26 
27 Description
28  Contact-angle boundary condition for multi-phase interface-capturing
29  simulations. Sets of coefficients are given for the contact angle with each
30  other phase. These coefficients can specify either a constant or a dynamic
31  contact angle.
32 
33 Usage
34  \table
35  Property | Description | Required | Default value
36  theta0 | Equilibrium contact angle | yes |
37  uTheta | Velocity scale | no | none
38  thetaA | Limiting advancing contact angle | if uTheta | none
39  thetaR | Limiting receding contact angle | if uTheta | none
40  \endtable
41 
42  Example of the boundary condition specification:
43  \verbatim
44  <patchName>
45  {
46  type alphaContactAngle;
47 
48  contactAngleProperties
49  {
50  // Constant contact angle with air
51  air
52  {
53  theta0 90;
54  }
55 
56  // Dynamic contact angle with water
57  oil
58  {
59  theta0 70;
60  uTheta 1;
61  thetaA 100;
62  thetaR 50;
63  }
64  }
65 
66  value uniform 0;
67  }
68  \endverbatim
69 
70 SourceFiles
71  alphaContactAngleFvPatchScalarField.C
72 
73 \*---------------------------------------------------------------------------*/
74 
75 #ifndef alphaContactAngleFvPatchScalarField_H
76 #define alphaContactAngleFvPatchScalarField_H
77 
79 
80 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
81 
82 namespace Foam
83 {
84 
85 /*---------------------------------------------------------------------------*\
86  Class alphaContactAngleFvPatchScalarField Declaration
87 \*---------------------------------------------------------------------------*/
88 
89 class alphaContactAngleFvPatchScalarField
90 :
91  public zeroGradientFvPatchScalarField
92 {
93 public:
94 
95  class contactAngleProperties
96  {
97  // Private Data
98 
99  //- Equilibrium contact angle
100  scalar theta0_;
101 
102  //- Is this contact angle dynamic?
103  bool dynamic_;
104 
105  //- Dynamic contact angle velocity scale
106  scalar uTheta_;
107 
108  //- Limiting advancing contact angle
109  scalar thetaA_;
110 
111  //- Limiting receding contact angle
112  scalar thetaR_;
113 
114 
115  public:
116 
117  // Constructors
118 
119  //- Construct null
121 
122  //- Construct non-dynamic properties from components
123  contactAngleProperties(const scalar theta0);
124 
125  //- Construct dynamic properties from components
127  (
128  const scalar theta0,
129  const scalar uTheta,
130  const scalar thetaA,
131  const scalar thetaR
132  );
133 
134  //- Construct from a dictionary
136 
137 
138  // Member Functions
139 
140  // Access
141 
142  //- Return the equilibrium contact angle theta0
143  inline scalar theta0() const
144  {
145  return theta0_;
146  }
147 
148  //- Return if this contact angle is dynamic
149  inline bool dynamic() const
150  {
151  return dynamic_;
152  }
153 
154  //- Return the dynamic contact angle velocity scale
155  inline scalar uTheta() const
156  {
157  return uTheta_;
158  }
159 
160  //- Return the limiting advancing contact angle
161  inline scalar thetaA() const
162  {
163  return thetaA_;
164  }
165 
166  //- Return the limiting receding contact angle
167  inline scalar thetaR() const
168  {
169  return thetaR_;
170  }
171 
172 
173  //- Return the contact angle properties for the reverse of this
174  // interface (e.g., convert air_water to water_air).
176 
177  //- Write to stream
178  void write(Ostream& os) const;
179 
180 
181  // Member operators
182 
183  //- Check equality between two sets of contact angle properties
184  bool operator==(const contactAngleProperties& thetaProps) const;
185 
186  //- Check inequality between two sets of contact angle properties
187  bool operator!=(const contactAngleProperties& thetaProps) const;
188  };
189 
190 
191 private:
192 
193  // Private Data
194 
195  //- Table of contact angle properties
197 
198 
199 public:
200 
201  //- Runtime type information
202  TypeName("alphaContactAngle");
203 
204 
205  // Constructors
206 
207  //- Construct from patch, internal field and dictionary
209  (
210  const fvPatch&,
212  const dictionary&
213  );
214 
215  //- Construct by mapping given alphaContactAngleFvPatchScalarField
216  // onto a new patch
218  (
220  const fvPatch&,
222  const fieldMapper&
223  );
224 
225  //- Construct as copy setting internal field reference
227  (
230  );
231 
232  //- Construct and return a clone setting internal field reference
234  (
236  ) const
237  {
239  (
241  );
242  }
243 
244 
245  // Member Functions
246 
247  //- Return the contact angle properties
249  {
250  return thetaProps_;
251  }
252 
253  //- Write
254  virtual void write(Ostream&) const;
255 };
256 
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 } // End namespace Foam
261 
262 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263 
264 #endif
265 
266 // ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An STL-conforming hash table.
Definition: HashTable.H:127
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
contactAngleProperties reversed() const
Return the contact angle properties for the reverse of this.
bool operator==(const contactAngleProperties &thetaProps) const
Check equality between two sets of contact angle properties.
bool operator!=(const contactAngleProperties &thetaProps) const
Check inequality between two sets of contact angle properties.
scalar theta0() const
Return the equilibrium contact angle theta0.
scalar thetaR() const
Return the limiting receding contact angle.
scalar thetaA() const
Return the limiting advancing contact angle.
scalar uTheta() const
Return the dynamic contact angle velocity scale.
Contact-angle boundary condition for multi-phase interface-capturing simulations. Sets of coefficient...
const HashTable< contactAngleProperties > & thetaProps() const
Return the contact angle properties.
alphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
TypeName("alphaContactAngle")
Runtime type information.
virtual tmp< fvPatchScalarField > clone(const DimensionedField< scalar, volMesh > &iF) const
Construct and return a clone setting internal field reference.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
dictionary dict