temperatureDependentSurfaceTension.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) 2017-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 
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace surfaceTensionModels
35 {
36  defineTypeNameAndDebug(temperatureDependent, 0);
38  (
39  surfaceTensionModel,
40  temperatureDependent,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict,
52  const fvMesh& mesh
53 )
54 :
55  surfaceTensionModel(mesh),
56  TName_(dict.lookupOrDefault<word>("T", "T")),
57  sigma_(Function1<scalar>::New("sigma", dict))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 {
72  tmp<volScalarField> tsigma
73  (
74  volScalarField::New("sigma", mesh_, dimSigma)
75  );
76  volScalarField& sigma = tsigma.ref();
77 
78  const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
79 
80  sigma.field() = sigma_->value(T.field());
81 
82  volScalarField::Boundary& sigmaBf = sigma.boundaryFieldRef();
83  const volScalarField::Boundary& TBf = T.boundaryField();
84 
85  forAll(sigmaBf, patchi)
86  {
87  sigmaBf[patchi] = sigma_->value(TBf[patchi]);
88  }
89 
90  return tsigma;
91 }
92 
93 
95 (
96  const dictionary& dict
97 )
98 {
99  const dictionary& sigmaDict = surfaceTensionModel::sigmaDict(dict);
100 
101  TName_ = sigmaDict.lookupOrDefault<word>("T", "T");
102  sigma_ = Function1<scalar>::New("sigma", sigmaDict);
103 
104  return true;
105 }
106 
107 
109 (
110  Ostream& os
111 ) const
112 {
114  {
115  os << sigma_() << token::END_STATEMENT << nl;
116  return os.good();
117  }
118  else
119  {
120  return false;
121  }
122 }
123 
124 
125 // ************************************************************************* //
static const dictionary & sigmaDict(const dictionary &dict)
Run-time selectable general function of one variable.
Definition: Function1.H:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
Macros for easy insertion into run-time selection tables.
temperatureDependent(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
virtual bool writeData(Ostream &os) const
Write in dictionary format.
A class for handling words, derived from string.
Definition: word.H:59
virtual bool writeData(Ostream &os) const =0
Write in dictionary format.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
const Field< Type > & field() const
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
virtual bool readDict(const dictionary &dict)
Update surface tension coefficient from given dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
defineTypeNameAndDebug(constant, 0)
virtual tmp< volScalarField > sigma() const
Surface tension coefficient.
addToRunTimeSelectionTable(surfaceTensionModel, constant, dictionary)
A class for managing temporary objects.
Definition: PtrList.H:53
Abstract base-class for surface tension models which return the surface tension coefficient field...
Namespace for OpenFOAM.
static autoPtr< Function1< Type > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32