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-2020 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 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  alphaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
58  theta0_(gcpsf.theta0_),
59  uTheta_(gcpsf.uTheta_),
60  thetaA_(gcpsf.thetaA_),
61  thetaR_(gcpsf.thetaR_)
62 {}
63 
64 
67 (
68  const fvPatch& p,
70  const dictionary& dict
71 )
72 :
74  theta0_(dict.lookup<scalar>("theta0")),
75  uTheta_(dict.lookup<scalar>("uTheta")),
76  thetaA_
77  (
78  dict.found("thetaA")
79  ? dict.lookup<scalar>("thetaA")
80  : dict.lookup<scalar>("thetaRec")
81  ),
82  thetaR_
83  (
84  dict.found("thetaR")
85  ? dict.lookup<scalar>("thetaR")
86  : dict.lookup<scalar>("thetaAdv")
87  )
88 {
89  evaluate();
90 }
91 
92 
95 (
97 )
98 :
100  theta0_(gcpsf.theta0_),
101  uTheta_(gcpsf.uTheta_),
102  thetaA_(gcpsf.thetaA_),
103  thetaR_(gcpsf.thetaR_)
104 {}
105 
106 
109 (
112 )
113 :
115  theta0_(gcpsf.theta0_),
116  uTheta_(gcpsf.uTheta_),
117  thetaA_(gcpsf.thetaA_),
118  thetaR_(gcpsf.thetaR_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
126 (
127  const fvPatchVectorField& Up,
128  const fvsPatchVectorField& nHat
129 ) const
130 {
131  if (uTheta_ < small)
132  {
133  return tmp<scalarField>(new scalarField(size(), theta0_));
134  }
135 
136  const vectorField nf(patch().nf());
137 
138  // Calculated the component of the velocity parallel to the wall
139  vectorField Uwall(Up.patchInternalField() - Up);
140  Uwall -= (nf & Uwall)*nf;
141 
142  // Find the direction of the interface parallel to the wall
143  vectorField nWall(nHat - (nf & nHat)*nf);
144 
145  // Normalise nWall
146  nWall /= (mag(nWall) + small);
147 
148  // Calculate Uwall resolved normal to the interface parallel to
149  // the interface
150  scalarField uwall(nWall & Uwall);
151 
152  return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
153 }
154 
155 
157 {
159  writeEntry(os, "theta0", theta0_);
160  writeEntry(os, "uTheta", uTheta_);
161  writeEntry(os, "thetaA", thetaA_);
162  writeEntry(os, "thetaR", thetaR_);
163  writeEntry(os, "value", *this);
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 namespace Foam
170 {
172  (
175  );
176 }
177 
178 
179 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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:194
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 > &)
Abstract base class for alphaContactAngle boundary conditions.
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
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812