Loading...
Searching...
No Matches
interfaceProperties.C
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-2017 OpenFOAM Foundation
9 Copyright (C) 2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "interfaceProperties.H"
32#include "surfaceInterpolate.H"
33#include "fvcDiv.H"
34#include "fvcGrad.H"
35#include "fvcSnGrad.H"
36#include "fvcAverage.H"
37#include "unitConversion.H"
38
39// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
40
41// Correction for the boundary condition on the unit normal nHat on
42// walls to produce the correct contact angle.
43
44// The dynamic contact angle is calculated from the component of the
45// velocity on the direction of the interface, parallel to the wall.
46
47void Foam::interfaceProperties::correctContactAngle
48(
49 surfaceVectorField::Boundary& nHatb,
50 const surfaceVectorField::Boundary& gradAlphaf
51) const
52{
53 const fvMesh& mesh = alpha1_.mesh();
54 const volScalarField::Boundary& abf = alpha1_.boundaryField();
55
56 const fvBoundaryMesh& boundary = mesh.boundary();
57
58 forAll(boundary, patchi)
59 {
61 {
62 alphaContactAngleTwoPhaseFvPatchScalarField& acap =
63 const_cast<alphaContactAngleTwoPhaseFvPatchScalarField&>
64 (
66 (
67 abf[patchi]
68 )
69 );
70
71 fvsPatchVectorField& nHatp = nHatb[patchi];
72 const scalarField theta
73 (
74 degToRad() * acap.theta(U_.boundaryField()[patchi], nHatp)
75 );
76
77 const vectorField nf
78 (
79 boundary[patchi].nf()
80 );
81
82 // Reset nHatp to correspond to the contact angle
83
84 const scalarField a12(nHatp & nf);
85 const scalarField b1(cos(theta));
86
87 scalarField b2(nHatp.size());
88 forAll(b2, facei)
89 {
90 b2[facei] = cos(acos(a12[facei]) - theta[facei]);
91 }
92
93 const scalarField det(1.0 - a12*a12);
94
95 scalarField a((b1 - a12*b2)/det);
96 scalarField b((b2 - a12*b1)/det);
97
98 nHatp = a*nf + b*nHatp;
99 nHatp /= (mag(nHatp) + deltaN_.value());
100
101 acap.gradient() = (nf & nHatp)*mag(gradAlphaf[patchi]);
102 acap.evaluate();
103 }
104 }
105}
106
107
108void Foam::interfaceProperties::calculateK()
109{
110 const fvMesh& mesh = alpha1_.mesh();
111 const surfaceVectorField& Sf = mesh.Sf();
112
113 // Cell gradient of alpha
114 tmp<volVectorField> tgradAlpha;
115 if (nAlphaSmoothCurvature_ < 1)
116 {
117 tgradAlpha = fvc::grad(alpha1_, "nHat");
118 }
119 else
120 {
121 // Smooth interface curvature to reduce spurious currents
122 auto talpha1L = tmp<volScalarField>::New(alpha1_);
123 auto& alpha1L = talpha1L.ref();
124
125 for (int i = 0; i < nAlphaSmoothCurvature_; ++i)
126 {
127 alpha1L = fvc::average(fvc::interpolate(alpha1L));
128 }
129
130 tgradAlpha = fvc::grad(talpha1L, "nHat");
131 }
132
133 // Interpolated face-gradient of alpha
134 surfaceVectorField gradAlphaf(fvc::interpolate(tgradAlpha));
135
136 //gradAlphaf -=
137 // (mesh.Sf()/mesh.magSf())
138 // *(fvc::snGrad(alpha1_) - (mesh.Sf() & gradAlphaf)/mesh.magSf());
139
140 // Face unit interface normal
141 surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
142 // surfaceVectorField nHatfv
143 // (
144 // (gradAlphaf + deltaN_*vector(0, 0, 1)
145 // *sign(gradAlphaf.component(vector::Z)))/(mag(gradAlphaf) + deltaN_)
146 // );
147 correctContactAngle(nHatfv.boundaryFieldRef(), gradAlphaf.boundaryField());
148
149 // Face unit interface normal flux
150 nHatf_ = nHatfv & Sf;
151
152 // Simple expression for curvature
153 K_ = -fvc::div(nHatf_);
154
155 // Complex expression for curvature.
156 // Correction is formally zero but numerically non-zero.
157 /*
158 volVectorField nHat(gradAlpha/(mag(gradAlpha) + deltaN_));
159 forAll(nHat.boundaryField(), patchi)
160 {
161 nHat.boundaryFieldRef()[patchi] = nHatfv.boundaryField()[patchi];
162 }
163
164 K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
165 */
166}
167
168
169// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170
171Foam::interfaceProperties::interfaceProperties
172(
173 const volScalarField& alpha1,
174 const volVectorField& U,
175 const IOdictionary& dict
176)
177:
178 transportPropertiesDict_(dict),
179 nAlphaSmoothCurvature_
180 (
181 alpha1.mesh().solverDict(alpha1.name()).
182 getOrDefault<int>("nAlphaSmoothCurvature", 0)
183 ),
184 cAlpha_
185 (
186 alpha1.mesh().solverDict(alpha1.name()).get<scalar>("cAlpha")
187 ),
188
189 sigmaPtr_(surfaceTensionModel::New(dict, alpha1.mesh())),
190
191 deltaN_
192 (
193 "deltaN",
194 1e-8/cbrt(average(alpha1.mesh().V()))
195 ),
196
197 alpha1_(alpha1),
198 U_(U),
199
200 nHatf_
201 (
203 (
204 "nHatf",
205 alpha1_.time().timeName(),
206 alpha1_.mesh()
207 ),
208 alpha1_.mesh(),
210 ),
211
212 K_
213 (
215 (
216 "interfaceProperties:K",
217 alpha1_.time().timeName(),
218 alpha1_.mesh()
219 ),
220 alpha1_.mesh(),
222 )
223{
224 calculateK();
225}
226
227
228// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
229
232{
233 return sigmaPtr_->sigma()*K_;
234}
235
236
242
243
246{
247 return pos0(alpha1_ - 0.01)*pos0(0.99 - alpha1_);
248}
249
252{
253 calculateK();
254}
255
256
258{
259 alpha1_.mesh().solverDict(alpha1_.name()).readEntry("cAlpha", cAlpha_);
260 sigmaPtr_->readDict(transportPropertiesDict_);
261
262 return true;
263}
264
265
266// ************************************************************************* //
const volScalarField & alpha1
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const surfaceVectorField & Sf() const
Return cell face area vectors.
tmp< surfaceScalarField > surfaceTensionForce() const
tmp< volScalarField > nearInterface() const
Indicator of the proximity of the interface.
tmp< volScalarField > sigmaK() const
bool read()
Read transportProperties dictionary.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
U
Definition pEqn.H:72
faceListList boundary
dynamicFvMesh & mesh
Area-weighted average a surfaceField creating a volField.
Calculate the divergence of the given field.
Calculate the gradient of the given field.
Calculate the snGrad of the given volField.
word timeName
Definition getTimeIndex.H:3
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > average(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Area-weighted average a surfaceField creating a volField.
Definition fvcAverage.C:41
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvcSnGrad.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar det(const dimensionedSphericalTensor &dt)
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< Type > average(const DimensionedField< Type, GeoMesh > &f1, const label comm)
fvsPatchField< vector > fvsPatchVectorField
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
dictionary dict
volScalarField & b
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.