Loading...
Searching...
No Matches
dynamicAlphaContactAngleFvPatchScalarField.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-2013 OpenFOAM Foundation
9 Copyright (C) 2019 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
31#include "fvPatchFieldMapper.H"
32#include "volMesh.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
38(
39 const fvPatch& p,
41)
42:
44 theta0_(0.0),
45 uTheta_(0.0),
46 thetaA_(0.0),
47 thetaR_(0.0)
48{}
49
50
53(
55 const fvPatch& p,
57 const fvPatchFieldMapper& mapper
58)
59:
61 theta0_(gcpsf.theta0_),
62 uTheta_(gcpsf.uTheta_),
63 thetaA_(gcpsf.thetaA_),
64 thetaR_(gcpsf.thetaR_)
65{}
66
67
70(
71 const fvPatch& p,
73 const dictionary& dict
74)
75:
77 theta0_(dict.get<scalar>("theta0")),
78 uTheta_(dict.get<scalar>("uTheta")),
79 thetaA_(dict.get<scalar>("thetaA")),
80 thetaR_(dict.get<scalar>("thetaR"))
81{
82 evaluate();
83}
84
85
88(
90)
91:
93 theta0_(gcpsf.theta0_),
94 uTheta_(gcpsf.uTheta_),
95 thetaA_(gcpsf.thetaA_),
96 thetaR_(gcpsf.thetaR_)
97{}
98
99
102(
105)
106:
108 theta0_(gcpsf.theta0_),
109 uTheta_(gcpsf.uTheta_),
110 thetaA_(gcpsf.thetaA_),
111 thetaR_(gcpsf.thetaR_)
112{}
113
114
115// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116
119(
120 const fvPatchVectorField& Up,
121 const fvsPatchVectorField& nHat
122) const
123{
124 if (uTheta_ < SMALL)
125 {
126 return tmp<scalarField>::New(size(), theta0_);
127 }
128
129 const vectorField nf(patch().nf());
130
131 // Calculated the component of the velocity parallel to the wall
132 vectorField Uwall(Up.patchInternalField() - Up);
133 Uwall -= (nf & Uwall)*nf;
134
135 // Find the direction of the interface parallel to the wall
136 vectorField nWall(nHat - (nf & nHat)*nf);
137
138 // Normalise nWall
139 nWall /= (mag(nWall) + SMALL);
140
141 // Calculate Uwall resolved normal to the interface parallel to
142 // the interface
143 scalarField uwall(nWall & Uwall);
144
145 return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
146}
147
148
150{
152 os.writeEntry("theta0", theta0_);
153 os.writeEntry("uTheta", uTheta_);
154 os.writeEntry("thetaA", thetaA_);
155 os.writeEntry("thetaR", thetaR_);
158
159
160// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161
162namespace Foam
163{
165 (
167 dynamicAlphaContactAngleFvPatchScalarField
168 );
169}
170
171
172// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Abstract base class for two-phase alphaContactAngle boundary conditions.
alphaContactAngleTwoPhaseFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::buffered)
Evaluate the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A dynamic alphaContactAngle scalar boundary condition (alphaContactAngleTwoPhaseFvPatchScalarField).
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const
Evaluate and return dynamic contact-angle.
dynamicAlphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
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
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
dimensionedScalar tanh(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
fvsPatchField< vector > fvsPatchVectorField
fvPatchField< vector > fvPatchVectorField
fvPatchField< scalar > fvPatchScalarField
dictionary dict