timeVaryingAlphaContactAngleFvPatchScalarField.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 #include "Time.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
42  t0_(0.0),
43  thetaT0_(0.0),
44  te_(0.0),
45  thetaTe_(0.0)
46 {}
47 
48 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
58  t0_(dict.lookup<scalar>("t0")),
59  thetaT0_(dict.lookup<scalar>("thetaT0")),
60  te_(dict.lookup<scalar>("te")),
61  thetaTe_(dict.lookup<scalar>("thetaTe"))
62 {
63  evaluate();
64 }
65 
66 
69 (
71  const fvPatch& p,
73  const fvPatchFieldMapper& mapper
74 )
75 :
76  alphaContactAngleFvPatchScalarField(gcpsf, p, iF, mapper),
77  t0_(gcpsf.t0_),
78  thetaT0_(gcpsf.thetaT0_),
79  te_(gcpsf.te_),
80  thetaTe_(gcpsf.te_)
81 {}
82 
83 
86 (
89 )
90 :
92  t0_(gcpsf.t0_),
93  thetaT0_(gcpsf.thetaT0_),
94  te_(gcpsf.te_),
95  thetaTe_(gcpsf.thetaTe_)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
103 (
104  const fvPatchVectorField&,
105  const fvsPatchVectorField&
106 ) const
107 {
108  scalar t = patch().boundaryMesh().mesh().time().value();
109  scalar theta0 = thetaT0_;
110 
111  if (t < t0_)
112  {
113  theta0 = thetaT0_;
114  }
115  else if (t > te_)
116  {
117  theta0 = thetaTe_;
118  }
119  else
120  {
121  theta0 = thetaT0_ + (t - t0_)*(thetaTe_ - thetaT0_)/(te_ - t0_);
122  }
123 
124  return tmp<scalarField>(new scalarField(size(), theta0));
125 }
126 
127 
129 (
130  Ostream& os
131 ) const
132 {
134  writeEntry(os, "t0", t0_);
135  writeEntry(os, "thetaT0", thetaT0_);
136  writeEntry(os, "te", te_);
137  writeEntry(os, "thetaTe", thetaTe_);
138  writeEntry(os, "value", *this);
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
144 namespace Foam
145 {
147  (
150  );
151 }
152 
153 // ************************************************************************* //
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
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:62
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
Macros for easy insertion into run-time selection tables.
A time-varying alphaContactAngle scalar boundary condition (alphaContactAngleFvPatchScalarField) ...
Foam::fvPatchFieldMapper.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
timeVaryingAlphaContactAngleFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
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...
Abstract base class for alphaContactAngle boundary conditions.
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
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:844