liquidPropertiesSurfaceTension.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-2018 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 "liquidThermo.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace surfaceTensionModels
35 {
36  defineTypeNameAndDebug(liquidProperties, 0);
38  (
39  surfaceTensionModel,
40  liquidProperties,
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  phaseName_(dict.lookup("phase"))
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
70 {
73  (
75  );
76 
77  const Foam::liquidProperties& liquid = thermo.mixture().properties();
78 
79  tmp<volScalarField> tsigma
80  (
82  (
83  "sigma",
84  mesh_,
85  dimSigma
86  )
87  );
88  volScalarField& sigma = tsigma.ref();
89 
90  const volScalarField& T = thermo.T();
91  const volScalarField& p = thermo.p();
92 
94  const volScalarField::Internal& pi = p;
95  const volScalarField::Internal& Ti = T;
96 
97  forAll(sigmai, celli)
98  {
99  sigmai[celli] = liquid.sigma(pi[celli], Ti[celli]);
100  }
101 
102  volScalarField::Boundary& sigmaBf = sigma.boundaryFieldRef();
103  const volScalarField::Boundary& pBf = p.boundaryField();
104  const volScalarField::Boundary& TBf = T.boundaryField();
105 
106  forAll(sigmaBf, patchi)
107  {
108  scalarField& sigmaPf = sigmaBf[patchi];
109  const scalarField& pPf = pBf[patchi];
110  const scalarField& TPf = TBf[patchi];
111 
112  forAll(sigmaPf, facei)
113  {
114  sigmaPf[facei] = liquid.sigma(pPf[facei], TPf[facei]);
115  }
116  }
117 
118  return tsigma;
119 }
120 
121 
123 (
124  const dictionary& dict
125 )
126 {
127  return true;
128 }
129 
130 
132 (
133  Ostream& os
134 ) const
135 {
137  {
138  return os.good();
139  }
140  else
141  {
142  return false;
143  }
144 }
145 
146 
147 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:79
liquidProperties(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].
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual bool writeData(Ostream &os) const
Write in dictionary format.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
virtual bool readDict(const dictionary &dict)
Update surface tension coefficient from given dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
heRhoThermo< rhoThermo, pureMixture< species::thermo< thermophysicalPropertiesSelector< liquidProperties >, sensibleInternalEnergy > >> heRhoThermopureMixtureliquidProperties
Definition: liquidThermo.H:53
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
The thermophysical properties of a liquid.
virtual bool writeData(Ostream &os) const =0
Write in dictionary format.
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
defineTypeNameAndDebug(constant, 0)
addToRunTimeSelectionTable(surfaceTensionModel, constant, dictionary)
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > sigma() const
Surface tension coefficient.
Namespace for OpenFOAM.