NSRDS0.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) 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 \*---------------------------------------------------------------------------*/
25 
26 #include "NSRDS0.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace Function1s
34 {
36 }
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const word& name,
45  const scalar a,
46  const scalar b,
47  const scalar c,
48  const scalar d,
49  const scalar e,
50  const scalar f
51 )
52 :
53  FieldFunction1<scalar, NSRDS0>(name),
54  a_(a),
55  b_(b),
56  c_(c),
57  d_(d),
58  e_(e),
59  f_(f)
60 {}
61 
62 
64 (
65  const word& name,
66  const unitConversions& units,
67  const dictionary& dict
68 )
69 :
70  FieldFunction1<scalar, NSRDS0>(name),
71  a_(dict.lookup<scalar>("a")),
72  b_(dict.lookup<scalar>("b")),
73  c_(dict.lookup<scalar>("c")),
74  d_(dict.lookup<scalar>("d")),
75  e_(dict.lookup<scalar>("e")),
76  f_(dict.lookup<scalar>("f"))
77 {
78  assertNoConvertUnits(typeName, units, dict);
79 }
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
85 (
86  const scalar x1,
87  const scalar x2
88 ) const
89 {
91  return 0;
92 }
93 
94 
96 (
97  const word& name,
98  const scalar a
99 ) const
100 {
101  if (f_ != 0)
102  {
104  << "Integral function of " << typeName << " function requested"
105  << " but the \"f\" coefficient is not zero"
106  << exit(FatalError);
107  }
108 
109  return NSRDS0(name, a, a_, b_/2, c_/3, d_/4, e_/5);
110 }
111 
112 
114 (
115  Ostream& os,
116  const unitConversions& units
117 ) const
118 {
119  writeEntry(os, "a", a_);
120  writeEntry(os, "b", b_);
121  writeEntry(os, "c", c_);
122  writeEntry(os, "d", d_);
123  writeEntry(os, "e", e_);
124  writeEntry(os, "f", f_);
125 }
126 
127 
128 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
NSRDS function number 100.
Definition: NSRDS0.H:72
virtual void write(Ostream &os, const unitConversions &units) const
Write the function coefficients.
Definition: NSRDS0.C:114
NSRDS0(const word &name, const scalar a, const scalar b, const scalar c, const scalar d, const scalar e, const scalar f)
Construct from components.
Definition: NSRDS0.C:43
virtual scalar integral(const scalar x1, const scalar x2) const
Integrate between two scalar values.
Definition: NSRDS0.C:85
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 handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField & b
Definition: createFields.H:25
void writeEntry(Ostream &os, const omega &a)
Definition: omega1.C:97
addScalarFunction1(laminarBL)
const dimensionedScalar e
Elementary charge.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const HashTable< unitConversion > & units()
Get the table of unit conversions.
error FatalError
void assertNoConvertUnits(const word &typeName, const Function1s::unitConversions &units, const dictionary &dict)
Generate an error in an context where unit conversions are not supported.
labelList f(nPoints)
dictionary dict