All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
qInterfaceProperties.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of OpenFOAM.
14 
15  OpenFOAM is free software: you can redistribute it and/or modify it
16  under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23  for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27 
28 Class
29  Foam::qInterfaceProperties
30 
31 Group
32  grpQInterfaceProperties
33 
34 Description
35  Contains the interface properties.
36 
37  Properties to aid interFoam:
38  -# Correct the alpha boundary condition for dynamic contact angle.
39  -# Calculate interface curvature.
40 
41 SourceFiles
42  qInterfaceProperties.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef qInterfaceProperties_H
47 #define qInterfaceProperties_H
48 
49 #include "IOdictionary.H"
50 #include "surfaceTensionModel.H"
51 #include "volFields.H"
52 #include "surfaceFields.H"
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 namespace Foam
57 {
58 
59 /*---------------------------------------------------------------------------*\
60  Class qInterfaceProperties Declaration
61 \*---------------------------------------------------------------------------*/
62 
64 {
65  // Private data
66 
67  //- Keep a reference to the transportProperties dictionary
68  const dictionary& transportPropertiesDict_;
69 
70  //- Compression coefficient
71  scalar cAlpha_;
72 
73  //- Surface tension
74  autoPtr<surfaceTensionModel> sigmaPtr_;
75 
76  //- Stabilisation for normalisation of the interface normal
77  const dimensionedScalar deltaN_;
78 
79  const volScalarField& alpha1_;
80  const volVectorField& U_;
81  surfaceScalarField nHatf_;
82  volScalarField K_;
83 
84 
85  // Private Member Functions
86 
87  //- No copy construct
89 
90  //- No copy assignment
91  void operator=(const qInterfaceProperties&) = delete;
92 
93  //- Correction for the boundary condition on the unit normal nHat on
94  // walls to produce the correct contact dynamic angle
95  // calculated from the component of U parallel to the wall
96  void correctContactAngle
97  (
98  surfaceVectorField::Boundary& nHat,
99  const surfaceVectorField::Boundary& gradAlphaf
100  ) const;
101 
102  //- Re-calculate the interface curvature
103  void calculateK();
104 
105 
106 public:
107 
108  // Constructors
109 
110  //- Construct from volume fraction field gamma and IOdictionary
112  (
113  const volScalarField& alpha1,
114  const volVectorField& U,
115  const IOdictionary&
116  );
117 
118 
119  // Member Functions
120 
121  scalar cAlpha() const
122  {
123  return cAlpha_;
124  }
125 
126  const dimensionedScalar& deltaN() const
127  {
128  return deltaN_;
129  }
130 
131  const surfaceScalarField& nHatf() const
132  {
133  return nHatf_;
134  }
136  tmp<volScalarField> sigmaK() const;
137 
138  tmp<surfaceScalarField> surfaceTensionForce() const;
139 
140  //- Indicator of the proximity of the interface
141  // Field values are 1 near and 0 away for the interface.
142  tmp<volScalarField> nearInterface() const;
143 
144  void correct();
146  //- Read transportProperties dictionary
147  bool read();
148 };
149 
150 
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
152 
153 } // End namespace Foam
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 #endif
158 
159 // ************************************************************************* //
const surfaceScalarField & nHatf() const
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
bool read()
Read transportProperties dictionary.
const dimensionedScalar & deltaN() const
tmp< surfaceScalarField > surfaceTensionForce() const
Contains the interface properties.
tmp< volScalarField > sigmaK() const
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.