dynamicAlphaContactAngleFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volMesh.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38 )
39 :
41  theta0_(0.0),
42  uTheta_(0.0),
43  thetaA_(0.0),
44  thetaR_(0.0)
45 {}
46 
47 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
57  theta0_(dict.lookup<scalar>("theta0")),
58  uTheta_(dict.lookup<scalar>("uTheta")),
59  thetaA_(dict.lookupBackwardsCompatible<scalar>({"thetaA", "thetaRec"})),
60  thetaR_(dict.lookupBackwardsCompatible<scalar>({"thetaR", "thetaAdv"}))
61 {
62  evaluate();
63 }
64 
65 
68 (
70  const fvPatch& p,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  alphaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
76  theta0_(gcpsf.theta0_),
77  uTheta_(gcpsf.uTheta_),
78  thetaA_(gcpsf.thetaA_),
79  thetaR_(gcpsf.thetaR_)
80 {}
81 
82 
85 (
88 )
89 :
91  theta0_(gcpsf.theta0_),
92  uTheta_(gcpsf.uTheta_),
93  thetaA_(gcpsf.thetaA_),
94  thetaR_(gcpsf.thetaR_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
102 (
103  const fvPatchVectorField& Up,
104  const fvsPatchVectorField& nHat
105 ) const
106 {
107  if (uTheta_ < small)
108  {
109  return tmp<scalarField>(new scalarField(size(), theta0_));
110  }
111 
112  const vectorField nf(patch().nf());
113 
114  // Calculated the component of the velocity parallel to the wall
115  vectorField Uwall(Up.patchInternalField() - Up);
116  Uwall -= (nf & Uwall)*nf;
117 
118  // Find the direction of the interface parallel to the wall
119  vectorField nWall(nHat - (nf & nHat)*nf);
120 
121  // Normalise nWall
122  nWall /= (mag(nWall) + small);
123 
124  // Calculate Uwall resolved normal to the interface parallel to
125  // the interface
126  scalarField uwall(nWall & Uwall);
127 
128  return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
129 }
130 
131 
133 {
135  writeEntry(os, "theta0", theta0_);
136  writeEntry(os, "uTheta", uTheta_);
137  writeEntry(os, "thetaA", thetaA_);
138  writeEntry(os, "thetaR", thetaR_);
139  writeEntry(os, "value", *this);
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 namespace Foam
146 {
148  (
151  );
152 }
153 
154 
155 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
A dynamic alphaContactAngle scalar boundary condition.
Macros for easy insertion into run-time selection tables.
Foam::fvPatchFieldMapper.
dynamicAlphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:174
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
Contact-angle boundary condition for multi-phase interface-capturing simulations. ...
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const
Evaluate and return dynamic contact-angle.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:875
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864