mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
Resolves bug-report http://www.openfoam.org/mantisbt/view.php?id=1938 Because C++ does not support overloading based on the return-type there is a problem defining both const and non-const member functions which are resolved based on the const-ness of the object for which they are called rather than the intent of the programmer declared via the const-ness of the returned type. The issue for the "boundaryField()" member function is that the non-const version increments the event-counter and checks the state of the stored old-time fields in case the returned value is altered whereas the const version has no side-effects and simply returns the reference. If the the non-const function is called within the patch-loop the event-counter may overflow. To resolve this it in necessary to avoid calling the non-const form of "boundaryField()" if the results is not altered and cache the reference outside the patch-loop when mutation of the patch fields is needed. The most straight forward way of resolving this problem is to name the const and non-const forms of the member functions differently e.g. the non-const form could be named: mutableBoundaryField() mutBoundaryField() nonConstBoundaryField() boundaryFieldRef() Given that in C++ a reference is non-const unless specified as const: "T&" vs "const T&" the logical convention would be boundaryFieldRef() boundaryFieldConstRef() and given that the const form which is more commonly used is it could simply be named "boundaryField()" then the logical convention is GeometricBoundaryField& boundaryFieldRef(); inline const GeometricBoundaryField& boundaryField() const; This is also consistent with the new "tmp" class for which non-const access to the stored object is obtained using the ".ref()" member function. This new convention for non-const access to the components of GeometricField will be applied to "dimensionedInternalField()" and "internalField()" in the future, i.e. "dimensionedInternalFieldRef()" and "internalFieldRef()".
229 lines
6.6 KiB
C
229 lines
6.6 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
|
|
\\ / O peration |
|
|
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
|
|
\\/ M anipulation |
|
|
-------------------------------------------------------------------------------
|
|
License
|
|
This file is part of OpenFOAM.
|
|
|
|
OpenFOAM is free software: you can redistribute it and/or modify it
|
|
under the terms of the GNU General Public License as published by
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
(at your option) any later version.
|
|
|
|
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
|
|
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
|
|
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
|
|
for more details.
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "threePhaseInterfaceProperties.H"
|
|
#include "alphaContactAngleFvPatchScalarField.H"
|
|
#include "mathematicalConstants.H"
|
|
#include "surfaceInterpolate.H"
|
|
#include "fvcDiv.H"
|
|
#include "fvcGrad.H"
|
|
#include "fvcSnGrad.H"
|
|
|
|
// * * * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * //
|
|
|
|
const Foam::scalar Foam::threePhaseInterfaceProperties::convertToRad =
|
|
Foam::constant::mathematical::pi/180.0;
|
|
|
|
|
|
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
|
|
|
|
void Foam::threePhaseInterfaceProperties::correctContactAngle
|
|
(
|
|
surfaceVectorField::GeometricBoundaryField& nHatb
|
|
) const
|
|
{
|
|
const volScalarField::GeometricBoundaryField& alpha1 =
|
|
mixture_.alpha1().boundaryField();
|
|
const volScalarField::GeometricBoundaryField& alpha2 =
|
|
mixture_.alpha2().boundaryField();
|
|
const volScalarField::GeometricBoundaryField& alpha3 =
|
|
mixture_.alpha3().boundaryField();
|
|
const volVectorField::GeometricBoundaryField& U =
|
|
mixture_.U().boundaryField();
|
|
|
|
const fvMesh& mesh = mixture_.U().mesh();
|
|
const fvBoundaryMesh& boundary = mesh.boundary();
|
|
|
|
forAll(boundary, patchi)
|
|
{
|
|
if (isA<alphaContactAngleFvPatchScalarField>(alpha1[patchi]))
|
|
{
|
|
const alphaContactAngleFvPatchScalarField& a2cap =
|
|
refCast<const alphaContactAngleFvPatchScalarField>
|
|
(alpha2[patchi]);
|
|
|
|
const alphaContactAngleFvPatchScalarField& a3cap =
|
|
refCast<const alphaContactAngleFvPatchScalarField>
|
|
(alpha3[patchi]);
|
|
|
|
scalarField twoPhaseAlpha2(max(a2cap, scalar(0)));
|
|
scalarField twoPhaseAlpha3(max(a3cap, scalar(0)));
|
|
|
|
scalarField sumTwoPhaseAlpha
|
|
(
|
|
twoPhaseAlpha2 + twoPhaseAlpha3 + SMALL
|
|
);
|
|
|
|
twoPhaseAlpha2 /= sumTwoPhaseAlpha;
|
|
twoPhaseAlpha3 /= sumTwoPhaseAlpha;
|
|
|
|
fvsPatchVectorField& nHatp = nHatb[patchi];
|
|
|
|
scalarField theta
|
|
(
|
|
convertToRad
|
|
* (
|
|
twoPhaseAlpha2*(180 - a2cap.theta(U[patchi], nHatp))
|
|
+ twoPhaseAlpha3*(180 - a3cap.theta(U[patchi], nHatp))
|
|
)
|
|
);
|
|
|
|
vectorField nf(boundary[patchi].nf());
|
|
|
|
// Reset nHatPatch to correspond to the contact angle
|
|
|
|
scalarField a12(nHatp & nf);
|
|
|
|
scalarField b1(cos(theta));
|
|
|
|
scalarField b2(nHatp.size());
|
|
|
|
forAll(b2, facei)
|
|
{
|
|
b2[facei] = cos(acos(a12[facei]) - theta[facei]);
|
|
}
|
|
|
|
scalarField det(1.0 - a12*a12);
|
|
|
|
scalarField a((b1 - a12*b2)/det);
|
|
scalarField b((b2 - a12*b1)/det);
|
|
|
|
nHatp = a*nf + b*nHatp;
|
|
|
|
nHatp /= (mag(nHatp) + deltaN_.value());
|
|
}
|
|
}
|
|
}
|
|
|
|
|
|
void Foam::threePhaseInterfaceProperties::calculateK()
|
|
{
|
|
const volScalarField& alpha1 = mixture_.alpha1();
|
|
|
|
const fvMesh& mesh = alpha1.mesh();
|
|
const surfaceVectorField& Sf = mesh.Sf();
|
|
|
|
// Cell gradient of alpha
|
|
volVectorField gradAlpha(fvc::grad(alpha1));
|
|
|
|
// Interpolated face-gradient of alpha
|
|
surfaceVectorField gradAlphaf(fvc::interpolate(gradAlpha));
|
|
|
|
// Face unit interface normal
|
|
surfaceVectorField nHatfv(gradAlphaf/(mag(gradAlphaf) + deltaN_));
|
|
|
|
correctContactAngle(nHatfv.boundaryFieldRef());
|
|
|
|
// Face unit interface normal flux
|
|
nHatf_ = nHatfv & Sf;
|
|
|
|
// Simple expression for curvature
|
|
K_ = -fvc::div(nHatf_);
|
|
|
|
// Complex expression for curvature.
|
|
// Correction is formally zero but numerically non-zero.
|
|
// volVectorField nHat = gradAlpha/(mag(gradAlpha) + deltaN_);
|
|
// nHat.boundaryField() = nHatfv.boundaryField();
|
|
// K_ = -fvc::div(nHatf_) + (nHat & fvc::grad(nHatfv) & nHat);
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
|
|
|
|
Foam::threePhaseInterfaceProperties::threePhaseInterfaceProperties
|
|
(
|
|
const incompressibleThreePhaseMixture& mixture
|
|
)
|
|
:
|
|
mixture_(mixture),
|
|
cAlpha_
|
|
(
|
|
readScalar
|
|
(
|
|
mixture.U().mesh().solverDict
|
|
(
|
|
mixture_.alpha1().name()
|
|
).lookup("cAlpha")
|
|
)
|
|
),
|
|
sigma12_(mixture.lookup("sigma12")),
|
|
sigma13_(mixture.lookup("sigma13")),
|
|
|
|
deltaN_
|
|
(
|
|
"deltaN",
|
|
1e-8/pow(average(mixture.U().mesh().V()), 1.0/3.0)
|
|
),
|
|
|
|
nHatf_
|
|
(
|
|
IOobject
|
|
(
|
|
"nHatf",
|
|
mixture.alpha1().time().timeName(),
|
|
mixture.alpha1().mesh()
|
|
),
|
|
mixture.alpha1().mesh(),
|
|
dimensionedScalar("nHatf", dimArea, 0.0)
|
|
),
|
|
|
|
K_
|
|
(
|
|
IOobject
|
|
(
|
|
"interfaceProperties:K",
|
|
mixture.alpha1().time().timeName(),
|
|
mixture.alpha1().mesh()
|
|
),
|
|
mixture.alpha1().mesh(),
|
|
dimensionedScalar("K", dimless/dimLength, 0.0)
|
|
)
|
|
{
|
|
calculateK();
|
|
}
|
|
|
|
|
|
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
|
|
|
|
Foam::tmp<Foam::surfaceScalarField>
|
|
Foam::threePhaseInterfaceProperties::surfaceTensionForce() const
|
|
{
|
|
return fvc::interpolate(sigmaK())*fvc::snGrad(mixture_.alpha1());
|
|
}
|
|
|
|
|
|
Foam::tmp<Foam::volScalarField>
|
|
Foam::threePhaseInterfaceProperties::nearInterface() const
|
|
{
|
|
return max
|
|
(
|
|
pos(mixture_.alpha1() - 0.01)*pos(0.99 - mixture_.alpha1()),
|
|
pos(mixture_.alpha2() - 0.01)*pos(0.99 - mixture_.alpha2())
|
|
);
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|