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 -------------------------------------------------------------------------------
11 License
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_);
156  writeEntry("value", os);
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 namespace Foam
163 {
165  (
168  );
169 }
170 
171 
172 // ************************************************************************* //
Foam::fvPatchField< vector >
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::dynamicAlphaContactAngleFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: dynamicAlphaContactAngleFvPatchScalarField.C:149
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dynamicAlphaContactAngleFvPatchScalarField::theta
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const
Evaluate and return dynamic contact-angle.
Definition: dynamicAlphaContactAngleFvPatchScalarField.C:119
Foam::alphaContactAngleTwoPhaseFvPatchScalarField::write
virtual void write(Ostream &os) const
Write.
Definition: alphaContactAngleTwoPhaseFvPatchScalarField.C:158
volMesh.H
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
fvPatchFieldMapper.H
Foam::alphaContactAngleTwoPhaseFvPatchScalarField
Abstract base class for two-phase alphaContactAngle boundary conditions.
Definition: alphaContactAngleTwoPhaseFvPatchScalarField.H:78
Foam::Field< vector >
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:272
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:233
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
dynamicAlphaContactAngleFvPatchScalarField.H
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::dynamicAlphaContactAngleFvPatchScalarField
A dynamic alphaContactAngle scalar boundary condition (alphaContactAngleTwoPhaseFvPatchScalarField)
Definition: dynamicAlphaContactAngleFvPatchScalarField.H:52
Foam::dynamicAlphaContactAngleFvPatchScalarField::dynamicAlphaContactAngleFvPatchScalarField
dynamicAlphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: dynamicAlphaContactAngleFvPatchScalarField.C:38
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::stringOps::evaluate
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Definition: stringOpsEvaluate.C:37