Loading...
Searching...
No Matches
contactAngleForce.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) 2018-2023 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 "contactAngleForce.H"
31#include "fvcGrad.H"
32#include "unitConversion.H"
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace regionModels
40{
41namespace surfaceFilmModels
42{
43
44// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45
47
48// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49
50void contactAngleForce::initialise()
51{
52 const wordRes zeroForcePatches
53 (
54 coeffDict_.getOrDefault<wordRes>("zeroForcePatches", wordRes())
55 );
56
57 if (zeroForcePatches.size())
58 {
60 const scalar dLim = coeffDict_.get<scalar>("zeroForceDistance");
61
62 Info<< " Assigning zero contact force within " << dLim
63 << " of patches:" << endl;
64
65 labelHashSet patchIDs = pbm.patchSet(zeroForcePatches);
66
67 for (const label patchi : patchIDs)
68 {
69 Info<< " " << pbm[patchi].name() << endl;
70 }
71
72 // Temporary implementation until run-time selection covers this case
75 (
77 (
78 "y",
84 ),
86 dimensionedScalar("y", dimLength, GREAT)
87 );
88 dist.correct(y);
89
90 mask_ = pos0(y - dimensionedScalar("dLim", dimLength, dLim));
91 }
92}
93
94
95// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96
97contactAngleForce::contactAngleForce
98(
99 const word& typeName,
100 surfaceFilmRegionModel& film,
101 const dictionary& dict
102)
103:
104 force(typeName, film, dict),
105 Ccf_(coeffDict_.get<scalar>("Ccf")),
106 mask_
107 (
109 (
110 IOobject::scopedName(typeName, "contactForceMask"),
111 filmModel_.time().timeName(),
112 filmModel_.regionMesh(),
113 IOobject::NO_READ,
114 IOobject::NO_WRITE,
115 IOobject::NO_REGISTER
116 ),
117 filmModel_.regionMesh(),
118 dimensionedScalar(word::null, dimless, 1.0)
119 )
121 initialise();
122}
123
124
125// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
128{}
129
130
131// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
132
134{
135 auto tForce = volVectorField::New
136 (
137 IOobject::scopedName(typeName, "contactForce"),
141 );
142 auto& force = tForce.ref().primitiveFieldRef();
143
144 const labelUList& own = filmModel_.regionMesh().owner();
145 const labelUList& nbr = filmModel_.regionMesh().neighbour();
146
147 const scalarField& magSf = filmModel_.magSf();
148
151
152 const tmp<volScalarField> ttheta = theta();
153 const volScalarField& theta = ttheta();
154
155 const volVectorField gradAlpha(fvc::grad(alpha));
156
157 forAll(nbr, facei)
158 {
159 const label cellO = own[facei];
160 const label cellN = nbr[facei];
161
162 label celli = -1;
163 if ((alpha[cellO] > 0.5) && (alpha[cellN] < 0.5))
164 {
165 celli = cellO;
166 }
167 else if ((alpha[cellO] < 0.5) && (alpha[cellN] > 0.5))
168 {
169 celli = cellN;
170 }
171
172 if (celli != -1 && mask_[celli] > 0.5)
173 {
174 const scalar invDx = filmModel_.regionMesh().deltaCoeffs()[facei];
175 const vector n =
176 gradAlpha[celli]/(mag(gradAlpha[celli]) + ROOTVSMALL);
177 const scalar cosTheta = cos(degToRad(theta[celli]));
178 force[celli] += Ccf_*n*sigma[celli]*(1 - cosTheta)/invDx;
179 }
180 }
181
182 forAll(alpha.boundaryField(), patchi)
183 {
184 if (!filmModel_.isCoupledPatch(patchi))
185 {
186 const fvPatchField<scalar>& alphaPf = alpha.boundaryField()[patchi];
187 const fvPatchField<scalar>& maskPf = mask_.boundaryField()[patchi];
188 const fvPatchField<scalar>& sigmaPf = sigma.boundaryField()[patchi];
189 const fvPatchField<scalar>& thetaPf = theta.boundaryField()[patchi];
190 const scalarField& invDx = alphaPf.patch().deltaCoeffs();
191 const labelUList& faceCells = alphaPf.patch().faceCells();
192
193 forAll(alphaPf, facei)
194 {
195 if (maskPf[facei] > 0.5)
196 {
197 label cellO = faceCells[facei];
198
199 if ((alpha[cellO] > 0.5) && (alphaPf[facei] < 0.5))
200 {
201 const vector n =
202 gradAlpha[cellO]
203 /(mag(gradAlpha[cellO]) + ROOTVSMALL);
204 const scalar cosTheta = cos(degToRad(thetaPf[facei]));
205 force[cellO] +=
206 Ccf_*n*sigmaPf[facei]*(1 - cosTheta)/invDx[facei];
207 }
208 }
209 }
210 }
211 }
212
213 force /= magSf;
214
215 if (filmModel_.regionMesh().time().writeTime())
216 {
217 tForce().write();
218 }
219
221
222 tfvm.ref() += tForce;
223
224 return tfvm;
225}
226
227
228// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229
230} // End namespace surfaceFilmModels
231} // End namespace regionModels
232} // End namespace Foam
233
234// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
labelList patchIDs
const polyBoundaryMesh & pbm
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition fvMesh.H:572
const labelUList & neighbour() const
Internal face neighbour.
Definition fvMesh.H:580
const fvPatch & patch() const noexcept
Return the patch.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient except for coupled patches for which the cell-centre to c...
Definition fvPatch.C:183
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
const fvMesh & regionMesh() const
Return the region mesh database.
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
Base-class for film contact angle force models.
virtual tmp< volScalarField > theta() const =0
Return the contact angle field.
virtual tmp< fvVectorMatrix > correct(volVectorField &U)
Correct.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Base class for film (stress-based) force models.
Definition force.H:56
virtual const volScalarField & sigma() const =0
Return the film surface tension [N/m].
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
virtual const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
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
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
Calculate the gradient of the given field.
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar pos0(const dimensionedScalar &ds)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionSet dimForce
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dimensionedScalar cos(const dimensionedScalar &ds)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
volScalarField & alpha
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Unit conversion functions.