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 | Copyright (C) 2011-2013 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_(readScalar(dict.lookup("theta0"))),
75  uTheta_(readScalar(dict.lookup("uTheta"))),
76  thetaA_(readScalar(dict.lookup("thetaA"))),
77  thetaR_(readScalar(dict.lookup("thetaR")))
78 {
79  evaluate();
80 }
81 
82 
85 (
87 )
88 :
90  theta0_(gcpsf.theta0_),
91  uTheta_(gcpsf.uTheta_),
92  thetaA_(gcpsf.thetaA_),
93  thetaR_(gcpsf.thetaR_)
94 {}
95 
96 
99 (
102 )
103 :
105  theta0_(gcpsf.theta0_),
106  uTheta_(gcpsf.uTheta_),
107  thetaA_(gcpsf.thetaA_),
108  thetaR_(gcpsf.thetaR_)
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
116 (
117  const fvPatchVectorField& Up,
118  const fvsPatchVectorField& nHat
119 ) const
120 {
121  if (uTheta_ < SMALL)
122  {
123  return tmp<scalarField>(new scalarField(size(), theta0_));
124  }
125 
126  const vectorField nf(patch().nf());
127 
128  // Calculated the component of the velocity parallel to the wall
129  vectorField Uwall(Up.patchInternalField() - Up);
130  Uwall -= (nf & Uwall)*nf;
131 
132  // Find the direction of the interface parallel to the wall
133  vectorField nWall(nHat - (nf & nHat)*nf);
134 
135  // Normalise nWall
136  nWall /= (mag(nWall) + SMALL);
137 
138  // Calculate Uwall resolved normal to the interface parallel to
139  // the interface
140  scalarField uwall(nWall & Uwall);
141 
142  return theta0_ + (thetaA_ - thetaR_)*tanh(uwall/uTheta_);
143 }
144 
145 
147 {
149  os.writeKeyword("theta0") << theta0_ << token::END_STATEMENT << nl;
150  os.writeKeyword("uTheta") << uTheta_ << token::END_STATEMENT << nl;
151  os.writeKeyword("thetaA") << thetaA_ << token::END_STATEMENT << nl;
152  os.writeKeyword("thetaR") << thetaR_ << token::END_STATEMENT << nl;
153  writeEntry("value", os);
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 namespace Foam
160 {
162  (
165  );
166 }
167 
168 
169 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
dimensionedScalar tanh(const dimensionedScalar &ds)
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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 (alphaContactAngleFvPatchScalarField) ...
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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:53
volScalarField scalarField(fieldObject, mesh)
static const char nl
Definition: Ostream.H:262
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:224
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
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
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:576