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-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 
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 {
72  mesh_.lookupObject<heRhoThermopureMixtureliquidProperties>
73  (
75  );
76 
77  tmp<volScalarField> tsigma
78  (
80  (
81  "sigma",
82  mesh_,
83  dimSigma
84  )
85  );
86  volScalarField& sigma = tsigma.ref();
87 
88  const volScalarField& T = thermo.T();
89  const volScalarField& p = thermo.p();
90 
92  const volScalarField::Internal& pi = p;
93  const volScalarField::Internal& Ti = T;
94 
95  forAll(sigmai, celli)
96  {
97  const heRhoThermopureMixtureliquidProperties::transportMixtureType&
98  liquid = thermo.cellTransportMixture(celli);
99 
100  sigmai[celli] = liquid.properties().sigma(pi[celli], Ti[celli]);
101  }
102 
103  volScalarField::Boundary& sigmaBf = sigma.boundaryFieldRef();
104  const volScalarField::Boundary& pBf = p.boundaryField();
105  const volScalarField::Boundary& TBf = T.boundaryField();
106 
107  forAll(sigmaBf, patchi)
108  {
109  scalarField& sigmaPf = sigmaBf[patchi];
110  const scalarField& pPf = pBf[patchi];
111  const scalarField& TPf = TBf[patchi];
112 
113  forAll(sigmaPf, facei)
114  {
115  const heRhoThermopureMixtureliquidProperties::transportMixtureType&
116  liquid = thermo.patchFaceTransportMixture(patchi, facei);
117 
118  sigmaPf[facei] = liquid.properties().sigma(pPf[facei], TPf[facei]);
119  }
120  }
121 
122  return tsigma;
123 }
124 
125 
127 (
128  const dictionary& dict
129 )
130 {
131  return true;
132 }
133 
134 
136 (
137  Ostream& os
138 ) const
139 {
141  {
142  return os.good();
143  }
144  else
145  {
146  return false;
147  }
148 }
149 
150 
151 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
liquidProperties(const dictionary &dict, const fvMesh &mesh)
Construct from dictionary and mesh.
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.
heRhoThermo< rhoThermo::composite, pureMixture< species::thermo< thermophysicalPropertiesSelector< liquidProperties >, sensibleInternalEnergy > >> heRhoThermopureMixtureliquidProperties
Definition: liquidThermo.H:53
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
virtual bool readDict(const dictionary &dict)
Update surface tension coefficient from given dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
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)
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
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.