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-2022 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 phaseInterface& interface
47 )
48 :
49  saturationModel(dict, interface),
50  function_(Function1<scalar>::New("function", dict))
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
64 (
65  const volScalarField& T
66 ) const
67 {
69  return volScalarField::null();
70 }
71 
72 
75 (
76  const volScalarField& T
77 ) const
78 {
80  return volScalarField::null();
81 }
82 
83 
86 (
87  const volScalarField& T
88 ) const
89 {
91  return volScalarField::null();
92 }
93 
94 
97 (
98  const volScalarField& p
99 ) const
100 {
101  tmp<volScalarField> tTsat
102  (
104  (
105  "Tsat",
106  p.mesh(),
108  )
109  );
110 
111  volScalarField& Tsat = tTsat.ref();
112 
113  forAll(Tsat, celli)
114  {
115  Tsat[celli] = function_->value(p[celli]);
116  }
117 
118  volScalarField::Boundary& TsatBf = Tsat.boundaryFieldRef();
119 
120  forAll(Tsat.boundaryField(), patchi)
121  {
122  scalarField& Tsatp = TsatBf[patchi];
123  const scalarField& pp = p.boundaryField()[patchi];
124 
125  forAll(Tsatp, facei)
126  {
127  Tsatp[facei] = function_->value(pp[facei]);
128 
129  }
130  }
131 
132  return tTsat;
133 }
134 
135 
136 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volScalarField > pSatPrime(const volScalarField &T) const
Saturation pressure derivative w.r.t. temperature.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of the boundary field.
function1(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
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:58
virtual ~function1()
Destructor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > pSat(const volScalarField &T) const
Saturation pressure.
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:353
virtual tmp< volScalarField > lnPSat(const volScalarField &T) const
Natural log of the saturation pressure.
const dimensionSet dimTemperature
Namespace for OpenFOAM.