APIdiffCoef.H
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) 2011-2024 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 Class
25  Foam::Function2s::APIdiffCoefFunc
26 
27 Description
28  API function for vapour mass diffusivity
29 
30  Source:
31  \verbatim
32  API (American Petroleum Institute)
33  Technical Data Book
34  \endverbatim
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef APIdiffCoef_H
39 #define APIdiffCoef_H
40 
41 #include "Function2.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 namespace Function2s
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class APIdiffCoef Declaration
52 \*---------------------------------------------------------------------------*/
53 
54 class APIdiffCoef
55 :
56  public FieldFunction2<scalar, APIdiffCoef>
57 {
58  // Private Data
59 
60  // API vapour mass diffusivity function coefficients
61  scalar a_, b_, wf_, wa_;
62 
63  // Helper variables
64  scalar alpha_, beta_;
65 
66 
67 public:
68 
69  //- Runtime type information
70  TypeName("APIdiffCoef");
71 
72 
73  // Constructors
74 
75  //- Construct from components
77  (
78  const word& name,
79  const scalar a,
80  const scalar b,
81  const scalar wf,
82  const scalar wa
83  );
84 
85  //- Construct from name and dictionary
87  (
88  const word& name,
89  const unitConversions& units,
90  const dictionary& dict
91  );
92 
93  //- Construct and return a clone
94  virtual tmp<Function2<scalar>> clone() const
95  {
96  return tmp<Function2<scalar>>(new APIdiffCoef(*this));
97  }
98 
99 
100  // Member Functions
101 
102  //- Inherit base class value methods to prevent clang warnings
104 
105  //- API vapour mass diffusivity function using properties from
106  // construction
107  virtual scalar value(scalar p, scalar T) const
108  {
109  return 3.6059e-3*(pow(1.8*T, 1.75))*alpha_/(p*beta_);
110  }
111 
112  //- API vapour mass diffusivity function using properties from
113  // construction - with specified binary pair
114  scalar value(scalar p, scalar T, scalar Wa) const
115  {
116  const scalar alphaBinary = sqrt(1/wf_ + 1/Wa);
117  return 3.6059e-3*(pow(1.8*T, 1.75))*alphaBinary/(p*beta_);
118  }
119 
120  //- Write the function coefficients
121  virtual void write(Ostream& os, const unitConversions& units) const;
122 };
123 
124 
125 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
126 
127 } // End namespace Function2s
128 } // End namespace Foam
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 #endif
133 
134 // ************************************************************************* //
const word & name() const
Return the name of the entry.
Definition: Function2.C:78
virtual scalar value(scalar p, scalar T) const
API vapour mass diffusivity function using properties from.
Definition: APIdiffCoef.H:106
TypeName("APIdiffCoef")
Runtime type information.
virtual void write(Ostream &os, const unitConversions &units) const
Write the function coefficients.
Definition: APIdiffCoef.C:83
APIdiffCoef(const word &name, const scalar a, const scalar b, const scalar wf, const scalar wa)
Construct from components.
Definition: APIdiffCoef.C:43
virtual tmp< Function2< scalar > > clone() const
Construct and return a clone.
Definition: APIdiffCoef.H:93
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
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
volScalarField & b
Definition: createFields.H:25
Namespace for OpenFOAM.
const HashTable< unitConversion > & units()
Get the table of unit conversions.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & p