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-2024 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 {
38  (
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_
58  (
59  Function1<scalar>::New
60  (
61  "sigma",
64  dict
65  )
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
80 {
81  tmp<volScalarField> tsigma
82  (
83  volScalarField::New("sigma", mesh_, dimSigma)
84  );
85  volScalarField& sigma = tsigma.ref();
86 
87  const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
88 
89  sigma.primitiveFieldRef() = sigma_->value(T.primitiveField());
90 
91  volScalarField::Boundary& sigmaBf = sigma.boundaryFieldRef();
92  const volScalarField::Boundary& TBf = T.boundaryField();
93 
94  forAll(sigmaBf, patchi)
95  {
96  sigmaBf[patchi] = sigma_->value(TBf[patchi]);
97  }
98 
99  return tsigma;
100 }
101 
102 
104 (
105  const dictionary& dict
106 )
107 {
108  const dictionary& sigmaDict = surfaceTensionModel::sigmaDict(dict);
109 
110  TName_ = sigmaDict.lookupOrDefault<word>("T", "T");
111  sigma_ =
113  (
114  "sigma",
117  sigmaDict
118  );
119 
120  return true;
121 }
122 
123 
125 (
126  Ostream& os
127 ) const
128 {
130  {
131  os << sigma_() << token::END_STATEMENT << nl;
132  return os.good();
133  }
134  else
135  {
136  return false;
137  }
138 }
139 
140 
141 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Run-time selectable general function of one variable.
Definition: Function1.H:125
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricBoundaryField class.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
const Type & value() const
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Abstract base-class for surface tension models which return the surface tension coefficient field.
static const dictionary & sigmaDict(const dictionary &dict)
virtual bool writeData(Ostream &os) const =0
Write in dictionary format.
virtual bool writeData(Ostream &os) const
Write in dictionary format.
virtual bool readDict(const dictionary &dict)
Update surface tension coefficient from given dictionary.
virtual tmp< volScalarField > sigma() const
Surface tension coefficient.
temperatureDependent(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
@ END_STATEMENT
Definition: token.H:108
A class for handling words, derived from string.
Definition: word.H:62
label patchi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
addToRunTimeSelectionTable(surfaceTensionModel, liquidProperties, dictionary)
defineTypeNameAndDebug(liquidProperties, 0)
Namespace for OpenFOAM.
const dimensionSet dimLength
const dimensionSet dimTemperature
const dimensionSet dimForce
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict