Loading...
Searching...
No Matches
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"
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
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 = volScalarField::New
76 (
77 "sigma",
79 mesh_,
80 dimSigma
81 );
82 auto& sigma = tsigma.ref();
83
84 const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
85
86 sigma.field() = sigma_->value(T.field());
87
88 volScalarField::Boundary& sigmaBf = sigma.boundaryFieldRef();
89 const volScalarField::Boundary& TBf = T.boundaryField();
90
91 forAll(sigmaBf, patchi)
92 {
93 sigmaBf[patchi] = sigma_->value(TBf[patchi]);
94 }
95
96 return tsigma;
97}
98
99
101(
102 const dictionary& dict
103)
104{
106
107 TName_ = sigmaDict.getOrDefault<word>("T", "T");
108 sigma_ = Function1<scalar>::New("sigma", sigmaDict, &mesh_);
109
110 return true;
111}
112
113
115(
116 Ostream& os
117) const
118{
120 {
121 os << sigma_();
122 os.endEntry();
123 return os.good();
124 }
125
126 return false;
127}
128
129
130// ************************************************************************* //
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:92
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
surfaceTensionModel(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
static autoPtr< surfaceTensionModel > New(const dictionary &dict, const phasePair &pair)
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:75
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299