temperatureDependentSurfaceTension.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) 2017 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30#include "volFields.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace surfaceTensionModels
38{
41 (
45 );
46}
47}
48
49
50// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51
53(
54 const dictionary& dict,
55 const fvMesh& mesh
56)
57:
59 TName_(dict.getOrDefault<word>("T", "T")),
60 sigma_(Function1<scalar>::New("sigma", dict, &mesh))
61{}
62
63
64// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
65
67{}
68
69
70// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
71
74{
75 auto tsigma = tmp<volScalarField>::New
76 (
78 (
79 "sigma",
80 mesh_.time().timeName(),
81 mesh_,
84 false
85 ),
86 mesh_,
87 dimSigma
88 );
89 auto& sigma = tsigma.ref();
90
91 const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
92
93 sigma.field() = sigma_->value(T.field());
94
95 volScalarField::Boundary& sigmaBf = sigma.boundaryFieldRef();
96 const volScalarField::Boundary& TBf = T.boundaryField();
97
98 forAll(sigmaBf, patchi)
99 {
100 sigmaBf[patchi] = sigma_->value(TBf[patchi]);
101 }
102
103 return tsigma;
104}
105
106
108(
109 const dictionary& dict
110)
111{
113
114 TName_ = sigmaDict.getOrDefault<word>("T", "T");
115 sigma_ = Function1<scalar>::New("sigma", sigmaDict, &mesh_);
116
117 return true;
118}
119
120
122(
123 Ostream& os
124) const
125{
126 if (surfaceTensionModel::writeData(os))
127 {
128 os << sigma_() << token::END_STATEMENT << nl;
129 return os.good();
130 }
131
132 return false;
133}
134
135
136// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
const Type & value() const
Return const reference to value.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
static const dictionary & sigmaDict(const dictionary &dict)
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.
A class for managing temporary objects.
Definition: tmp.H:65
@ END_STATEMENT
End entry [isseparator].
Definition: token.H:154
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333