timeVaryingAlphaContactAngleFvPatchScalarField.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 #include "Time.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
45  t0_(0.0),
46  thetaT0_(0.0),
47  te_(0.0),
48  thetaTe_(0.0)
49 {}
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
62  t0_(gcpsf.t0_),
63  thetaT0_(gcpsf.thetaT0_),
64  te_(gcpsf.te_),
65  thetaTe_(gcpsf.te_)
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
78  t0_(dict.get<scalar>("t0")),
79  thetaT0_(dict.get<scalar>("thetaT0")),
80  te_(dict.get<scalar>("te")),
81  thetaTe_(dict.get<scalar>("thetaTe"))
82 {
83  evaluate();
84 }
85 
86 
89 (
92 )
93 :
95  t0_(gcpsf.t0_),
96  thetaT0_(gcpsf.thetaT0_),
97  te_(gcpsf.te_),
98  thetaTe_(gcpsf.thetaTe_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
106 (
107  const fvPatchVectorField&,
108  const fvsPatchVectorField&
109 ) const
110 {
111  scalar t = patch().boundaryMesh().mesh().time().value();
112  scalar theta0 = thetaT0_;
113 
114  if (t < t0_)
115  {
116  theta0 = thetaT0_;
117  }
118  else if (t > te_)
119  {
120  theta0 = thetaTe_;
121  }
122  else
123  {
124  theta0 = thetaT0_ + (t - t0_)*(thetaTe_ - thetaT0_)/(te_ - t0_);
125  }
126 
127  return tmp<scalarField>::New(size(), theta0);
128 }
129 
130 
132 (
133  Ostream& os
134 ) const
135 {
137  os.writeEntry("t0", t0_);
138  os.writeEntry("thetaT0", thetaT0_);
139  os.writeEntry("te", te_);
140  os.writeEntry("thetaTe", thetaTe_);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146 
147 namespace Foam
148 {
150  (
152  timeVaryingAlphaContactAngleFvPatchScalarField
153  );
154 }
155 
156 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual tmp< scalarField > theta(const fvPatchVectorField &Up, const fvsPatchVectorField &nHat) const
Evaluate and return the time-varying equilibrium contact-angle.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
Macros for easy insertion into run-time selection tables.
A time-varying alphaContactAngle scalar boundary condition (alphaContactAngleTwoPhaseFvPatchScalarFie...
fvPatchField< scalar > fvPatchScalarField
A FieldMapper for finite-volume patch fields.
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.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Abstract base class for two-phase alphaContactAngle boundary conditions.
OBJstream os(runTime.globalPath()/outputName)
timeVaryingAlphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Namespace for OpenFOAM.