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.
30 
31 SourceFiles
32  alphaContactAngleFvPatchScalarField.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef alphaContactAngleFvPatchScalarField_H
37 #define alphaContactAngleFvPatchScalarField_H
38 
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class alphaContactAngleFvPatchScalarField Declaration
48 \*---------------------------------------------------------------------------*/
49 
51 :
52  public zeroGradientFvPatchScalarField
53 {
54 public:
55 
57  {
58  // Private Data
59 
60  //- Equilibrium contact angle
61  scalar theta0_;
62 
63  //- Is this contact angle dynamic?
64  bool dynamic_;
65 
66  //- Dynamic contact angle velocity scale
67  scalar uTheta_;
68 
69  //- Limiting advancing contact angle
70  scalar thetaA_;
71 
72  //- Limiting receding contact angle
73  scalar thetaR_;
74 
75 
76  public:
77 
78  // Constructors
79 
80  //- Construct null
82 
83  //- Construct non-dynamic properties from components
84  contactAngleProperties(const scalar theta0);
85 
86  //- Construct dynamic properties from components
88  (
89  const scalar theta0,
90  const scalar uTheta,
91  const scalar thetaA,
92  const scalar thetaR
93  );
94 
95  //- Construct from a dictionary
97 
98 
99  // Member Functions
100 
101  // Access
102 
103  //- Return the equilibrium contact angle theta0
104  inline scalar theta0() const
105  {
106  return theta0_;
107  }
108 
109  //- Return if this contact angle is dynamic
110  inline bool dynamic() const
111  {
112  return dynamic_;
113  }
114 
115  //- Return the dynamic contact angle velocity scale
116  inline scalar uTheta() const
117  {
118  return uTheta_;
119  }
120 
121  //- Return the limiting advancing contact angle
122  inline scalar thetaA() const
123  {
124  return thetaA_;
125  }
126 
127  //- Return the limiting receding contact angle
128  inline scalar thetaR() const
129  {
130  return thetaR_;
131  }
132 
133 
134  //- Return the contact angle properties for the reverse of this
135  // interface (e.g., convert air_water to water_air).
137 
138  //- Write to stream
139  void write(Ostream& os) const;
140 
141 
142  // Member operators
143 
144  //- Check equality between two sets of contact angle properties
145  bool operator==(const contactAngleProperties& thetaProps) const;
146 
147  //- Check inequality between two sets of contact angle properties
148  bool operator!=(const contactAngleProperties& thetaProps) const;
149  };
150 
151 
152 private:
153 
154  // Private Data
155 
156  //- Table of contact angle properties
158 
159 
160 public:
161 
162  //- Runtime type information
163  TypeName("alphaContactAngle");
164 
165 
166  // Constructors
167 
168  //- Construct from patch, internal field and dictionary
170  (
171  const fvPatch&,
173  const dictionary&
174  );
175 
176  //- Construct by mapping given alphaContactAngleFvPatchScalarField
177  // onto a new patch
179  (
181  const fvPatch&,
183  const fvPatchFieldMapper&
184  );
185 
186  //- Construct as copy setting internal field reference
188  (
191  );
192 
193  //- Construct and return a clone setting internal field reference
195  (
197  ) const
198  {
200  (
202  );
203  }
204 
205 
206  // Member Functions
207 
208  //- Return the contact angle properties
210  {
211  return thetaProps_;
212  }
213 
214  //- Write
215  virtual void write(Ostream&) const;
216 };
217 
218 
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
220 
221 } // End namespace Foam
222 
223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
224 
225 #endif
226 
227 // ************************************************************************* //
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.
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:160
Foam::fvPatchFieldMapper.
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