constantSurfaceTension.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 
26 #include "constantSurfaceTension.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace surfaceTensionModels
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const fvMesh& mesh
48 )
49 :
50  surfaceTensionModel(mesh),
51  sigma_("sigma", dimensionSet(1, 0, -2, 0, 0), dict)
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
65 {
66  return volScalarField::New(sigma_.name(), mesh_, sigma_);
67 }
68 
69 
71 {
72  // Handle sub-dictionary format as a special case
73  if (dict.isDict("sigma"))
74  {
75  dict.subDict("sigma").lookup("sigma") >> sigma_;
76  }
77  else
78  {
79  dict.lookup("sigma") >> sigma_;
80  }
81 
82  return true;
83 }
84 
85 
87 {
89  {
90  os << sigma_ << token::END_STATEMENT << nl;
91  return os.good();
92  }
93  else
94  {
95  return false;
96  }
97 }
98 
99 
100 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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
Dimension set for the base types.
Definition: dimensionSet.H:125
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.
virtual bool writeData(Ostream &os) const =0
Write in dictionary format.
Uniform constant surface tension model.
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.
constant(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
A class for managing temporary objects.
Definition: tmp.H:55
@ END_STATEMENT
Definition: token.H:108
addToRunTimeSelectionTable(surfaceTensionModel, liquidProperties, dictionary)
defineTypeNameAndDebug(liquidProperties, 0)
Namespace for OpenFOAM.
static const char nl
Definition: Ostream.H:266
dictionary dict