function1.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 
26 #include "function1.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace saturationModels
34 {
35  defineTypeNameAndDebug(function1, 0);
36  addToRunTimeSelectionTable(saturationModel, function1, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const dictionary& dict,
46  const objectRegistry& db
47 )
48 :
49  saturationModel(db),
50  function_
51  (
52  Function1<scalar>::New("function", dict)
53  )
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
67 (
68  const volScalarField& T
69 ) const
70 {
72  return volScalarField::null();
73 }
74 
75 
78 (
79  const volScalarField& T
80 ) const
81 {
83  return volScalarField::null();
84 }
85 
86 
89 (
90  const volScalarField& T
91 ) const
92 {
94  return volScalarField::null();
95 }
96 
97 
100 (
101  const volScalarField& p
102 ) const
103 {
104  tmp<volScalarField> tTsat
105  (
106  new volScalarField
107  (
108  IOobject
109  (
110  "Tsat",
111  p.mesh().time().timeName(),
112  p.mesh(),
115  ),
116  p.mesh(),
118  )
119  );
120 
121  volScalarField& Tsat = tTsat.ref();
122 
123  forAll(Tsat, celli)
124  {
125  Tsat[celli] = function_->value(p[celli]);
126  }
127 
128  volScalarField::Boundary& TsatBf = Tsat.boundaryFieldRef();
129 
130  forAll(Tsat.boundaryField(), patchi)
131  {
132  scalarField& Tsatp = TsatBf[patchi];
133  const scalarField& pp = p.boundaryField()[patchi];
134 
135  forAll(Tsatp, facei)
136  {
137  Tsatp[facei] = function_->value(pp[facei]);
138 
139  }
140  }
141 
142  return tTsat;
143 }
144 
145 
146 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivetive w.r.t. temperature.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:209
virtual tmp< volScalarField > Tsat(const volScalarField &p) const
Saturation temperature.
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual ~function1()
Destructor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
function1(const dictionary &dict, const objectRegistry &db)
Construct from a dictionary.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Internal & ref()
Return a reference to the dimensioned internal field.
static const GeometricField< scalar, fvPatchField, volMesh > & null()
Return a null geometric field.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
Namespace for OpenFOAM.