Commit 054eec50 authored by Henry Weller's avatar Henry Weller
Browse files

surfaceFilmModels::contactAngleForce: Added temperatureDependentContactAngleForce

Created a base-class from contactAngleForce from which the
distributionContactAngleForce (for backward compatibility) and the new
temperatureDependentContactAngleForce are derived:

Description
    Temperature dependent contact angle force

    The contact angle in degrees is specified as a \c Function1 type, to
    enable the use of, e.g.  contant, polynomial, table values.

See also
    Foam::regionModels::surfaceFilmModels::contactAngleForce
    Foam::Function1Types

SourceFiles
    temperatureDependentContactAngleForce.C
parent 5efae4ca
......@@ -115,7 +115,7 @@ Usage
value uniform 0.01;
\endverbatim
SeeAlso
See also
Foam::alphatPhaseChangeJayatillekeWallFunctionFvPatchField
SourceFiles
......
......@@ -63,31 +63,6 @@ void Foam::regionModels::regionModel::constructMeshObjects()
}
void Foam::regionModels::regionModel::constructMeshObjects
(
const dictionary& dict
)
{
// construct region mesh
if (!time_.foundObject<fvMesh>(regionName_))
{
regionMeshPtr_.reset
(
new fvMesh
(
IOobject
(
regionName_,
time_.timeName(),
time_,
IOobject::MUST_READ
)
)
);
}
}
void Foam::regionModels::regionModel::initialise()
{
if (debug)
......@@ -485,7 +460,7 @@ Foam::regionModels::regionModel::regionModel
{
if (active_)
{
constructMeshObjects(dict);
constructMeshObjects();
initialise();
if (readFields)
......
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
......@@ -73,9 +73,6 @@ private:
//- Construct region mesh and fields
void constructMeshObjects();
//- Construct region mesh and dictionary
void constructMeshObjects(const dictionary& dict);
//- Initialise the region
void initialise();
......
......@@ -13,7 +13,9 @@ KINEMATICMODELS=submodels/kinematic
$(KINEMATICMODELS)/force/force/force.C
$(KINEMATICMODELS)/force/force/forceNew.C
$(KINEMATICMODELS)/force/forceList/forceList.C
$(KINEMATICMODELS)/force/contactAngleForce/contactAngleForce.C
$(KINEMATICMODELS)/force/contactAngleForces/contactAngleForce/contactAngleForce.C
$(KINEMATICMODELS)/force/contactAngleForces/distribution/distributionContactAngleForce.C
$(KINEMATICMODELS)/force/contactAngleForces/temperatureDependent/temperatureDependentContactAngleForce.C
$(KINEMATICMODELS)/force/thermocapillaryForce/thermocapillaryForce.C
$(KINEMATICMODELS)/injectionModel/injectionModel/injectionModel.C
......
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
......@@ -42,7 +42,6 @@ namespace surfaceFilmModels
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(contactAngleForce, 0);
addToRunTimeSelectionTable(force, contactAngleForce, dictionary);
// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
......@@ -93,21 +92,13 @@ void contactAngleForce::initialise()
contactAngleForce::contactAngleForce
(
const word& typeName,
surfaceFilmModel& owner,
const dictionary& dict
)
:
force(typeName, owner, dict),
Ccf_(readScalar(coeffDict_.lookup("Ccf"))),
rndGen_(label(0), -1),
distribution_
(
distributionModels::distributionModel::New
(
coeffDict_.subDict("contactAngleDistribution"),
rndGen_
)
),
mask_
(
IOobject
......@@ -163,7 +154,10 @@ tmp<fvVectorMatrix> contactAngleForce::correct(volVectorField& U)
const volScalarField& alpha = owner_.alpha();
const volScalarField& sigma = owner_.sigma();
volVectorField gradAlpha(fvc::grad(alpha));
const tmp<volScalarField> ttheta = theta();
const volScalarField& theta = ttheta();
const volVectorField gradAlpha(fvc::grad(alpha));
forAll(nbr, facei)
{
......@@ -185,14 +179,14 @@ tmp<fvVectorMatrix> contactAngleForce::correct(volVectorField& U)
const scalar invDx = owner_.regionMesh().deltaCoeffs()[facei];
const vector n =
gradAlpha[celli]/(mag(gradAlpha[celli]) + ROOTVSMALL);
scalar theta = cos(degToRad(distribution_->sample()));
force[celli] += Ccf_*n*sigma[celli]*(1.0 - theta)/invDx;
const scalar cosTheta = cos(degToRad(theta[celli]));
force[celli] += Ccf_*n*sigma[celli]*(1 - cosTheta)/invDx;
}
}
forAll(alpha.boundaryField(), patchi)
{
if (!owner().isCoupledPatch(patchi))
if (!owner_.isCoupledPatch(patchi))
{
const fvPatchField<scalar>& alphaf = alpha.boundaryField()[patchi];
const fvPatchField<scalar>& maskf = mask_.boundaryField()[patchi];
......@@ -210,9 +204,9 @@ tmp<fvVectorMatrix> contactAngleForce::correct(volVectorField& U)
const vector n =
gradAlpha[cellO]
/(mag(gradAlpha[cellO]) + ROOTVSMALL);
scalar theta = cos(degToRad(distribution_->sample()));
const scalar cosTheta = cos(degToRad(theta[cellO]));
force[cellO] +=
Ccf_*n*sigma[cellO]*(1.0 - theta)/invDx[facei];
Ccf_*n*sigma[cellO]*(1 - cosTheta)/invDx[facei];
}
}
}
......
......@@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
......@@ -25,7 +25,7 @@ Class
Foam::contactAngleForce
Description
Film contact angle force
Base-class for film contact angle force models.
The effect of the contact angle force can be ignored over a specified
distance from patches.
......@@ -39,8 +39,6 @@ SourceFiles
#define contactAngleForce_H
#include "force.H"
#include "distributionModel.H"
#include "cachedRandom.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
......@@ -59,19 +57,11 @@ class contactAngleForce
:
public force
{
private:
// Private Data
//- Coefficient applied to the contact angle force
scalar Ccf_;
//- Random number generator
cachedRandom rndGen_;
//- Parcel size PDF model
const autoPtr<distributionModels::distributionModel> distribution_;
//- Mask over which force is applied
volScalarField mask_;
......@@ -88,6 +78,12 @@ private:
void operator=(const contactAngleForce&);
protected:
//- Return the contact angle field
virtual tmp<volScalarField> theta() const = 0;
public:
//- Runtime type information
......@@ -99,6 +95,7 @@ public:
//- Construct from surface film model
contactAngleForce
(
const word& typeName,
surfaceFilmModel& owner,
const dictionary& dict
);
......@@ -110,10 +107,8 @@ public:
// Member Functions
// Evolution
//- Correct
virtual tmp<fvVectorMatrix> correct(volVectorField& U);
//- Correct
virtual tmp<fvVectorMatrix> correct(volVectorField& U);
};
......
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 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 "distributionContactAngleForce.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace surfaceFilmModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(distributionContactAngleForce, 0);
addToRunTimeSelectionTable(force, distributionContactAngleForce, dictionary);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
distributionContactAngleForce::distributionContactAngleForce
(
surfaceFilmModel& owner,
const dictionary& dict
)
:
contactAngleForce(typeName, owner, dict),
rndGen_(label(0), -1),
distribution_
(
distributionModels::distributionModel::New
(
coeffDict_.subDict("distribution"),
rndGen_
)
)
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
distributionContactAngleForce::~distributionContactAngleForce()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
tmp<volScalarField> distributionContactAngleForce::theta() const
{
tmp<volScalarField> ttheta
(
new volScalarField
(
IOobject
(
typeName + ":theta",
owner_.time().timeName(),
owner_.regionMesh()
),
owner_.regionMesh(),
dimensionedScalar("0", dimless, 0)
)
);
volScalarField& theta = ttheta.ref();
volScalarField::Internal& thetai = theta.ref();
forAll(thetai, celli)
{
thetai[celli] = distribution_->sample();
}
forAll(theta.boundaryField(), patchi)
{
if (!owner_.isCoupledPatch(patchi))
{
fvPatchField<scalar>& thetaf = theta.boundaryFieldRef()[patchi];
forAll(thetaf, facei)
{
thetaf[facei] = distribution_->sample();
}
}
}
return ttheta;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace surfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 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/>.
Class
Foam::distributionContactAngleForce
Description
PDF distribution based film contact angle force
See also
Foam::regionModels::surfaceFilmModels::contactAngleForce
Foam::distributionModels::distributionModel
SourceFiles
distributionContactAngleForce.C
\*---------------------------------------------------------------------------*/
#ifndef distributionContactAngleForce_H
#define distributionContactAngleForce_H
#include "contactAngleForce.H"
#include "distributionModel.H"
#include "cachedRandom.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace surfaceFilmModels
{
/*---------------------------------------------------------------------------*\
Class distributionContactAngleForce Declaration
\*---------------------------------------------------------------------------*/
class distributionContactAngleForce
:
public contactAngleForce
{
// Private Data
//- Random number generator
cachedRandom rndGen_;
//- Parcel size PDF model
const autoPtr<distributionModels::distributionModel> distribution_;
// Private member functions
//- Disallow default bitwise copy construct
distributionContactAngleForce(const distributionContactAngleForce&);
//- Disallow default bitwise assignment
void operator=(const distributionContactAngleForce&);
protected:
//- Return the contact angle field
virtual tmp<volScalarField> theta() const;
public:
//- Runtime type information
TypeName("distributionContactAngle");
// Constructors
//- Construct from surface film model
distributionContactAngleForce
(
surfaceFilmModel& owner,
const dictionary& dict
);
//- Destructor
virtual ~distributionContactAngleForce();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace surfaceFilmModels
} // End namespace regionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2017 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 "temperatureDependentContactAngleForce.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace regionModels
{
namespace surfaceFilmModels
{
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
defineTypeNameAndDebug(temperatureDependentContactAngleForce, 0);
addToRunTimeSelectionTable
(
force,
temperatureDependentContactAngleForce,
dictionary
);
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
temperatureDependentContactAngleForce::temperatureDependentContactAngleForce
(
surfaceFilmModel& owner,
const dictionary& dict
)
:
contactAngleForce(typeName, owner, dict),
thetaPtr_(Function1<scalar>::New("theta", coeffDict_))
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
temperatureDependentContactAngleForce::~temperatureDependentContactAngleForce()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
tmp<volScalarField> temperatureDependentContactAngleForce::theta() const
{
tmp<volScalarField> ttheta
(
new volScalarField
(
IOobject
(
typeName + ":theta",
owner_.time().timeName(),
owner_.regionMesh()
),
owner_.regionMesh(),
dimensionedScalar("0", dimless, 0)
)
);
volScalarField& theta = ttheta.ref();
const volScalarField& T = owner_.T();
theta.ref().field() = thetaPtr_->value(T());