cubicEqnI.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29
31{}
32
33
35:
36 VectorSpace<cubicEqn, scalar, 4>(Foam::zero{})
37{}
38
39
41(
42 const scalar a,
43 const scalar b,
44 const scalar c,
45 const scalar d
46)
47{
48 this->v_[A] = a;
49 this->v_[B] = b;
50 this->v_[C] = c;
51 this->v_[D] = d;
52}
53
54
55// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56
57inline Foam::scalar Foam::cubicEqn::a() const
58{
59 return this->v_[A];
60}
61
62
63inline Foam::scalar Foam::cubicEqn::b() const
64{
65 return this->v_[B];
66}
67
68
69inline Foam::scalar Foam::cubicEqn::c() const
70{
71 return this->v_[C];
72}
73
74
75inline Foam::scalar Foam::cubicEqn::d() const
76{
77 return this->v_[D];
78}
79
80
81inline Foam::scalar& Foam::cubicEqn::a()
82{
83 return this->v_[A];
84}
85
86
87inline Foam::scalar& Foam::cubicEqn::b()
88{
89 return this->v_[B];
90}
91
92
93inline Foam::scalar& Foam::cubicEqn::c()
94{
95 return this->v_[C];
96}
97
98
99inline Foam::scalar& Foam::cubicEqn::d()
100{
101 return this->v_[D];
102}
103
104
105inline Foam::scalar Foam::cubicEqn::value(const scalar x) const
106{
107 return x*(x*(x*a() + b()) + c()) + d();
108}
109
110
111inline Foam::scalar Foam::cubicEqn::derivative(const scalar x) const
112{
113 return x*(x*3*a() + 2*b()) + c();
114}
115
116
117inline Foam::scalar Foam::cubicEqn::error(const scalar x) const
118{
119 return
120 SMALL*magSqr(x)*(mag(x*a()) + mag(b()))
121 + SMALL*mag(x)*(mag(x*(x*a() + b())) + mag(c()))
122 + SMALL*(mag(x*(x*(x*a() + b()) + c())) + mag(d()));
123}
124
125
126// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Graphite solid properties.
Definition: C.H:53
Templated vector space.
Definition: VectorSpace.H:79
Container to encapsulate various operations for cubic equation of the forms with real coefficients:
Definition: cubicEqn.H:115
cubicEqn()
Construct null.
Definition: cubicEqnI.H:30
scalar d() const
Definition: cubicEqnI.H:75
scalar value(const scalar x) const
Evaluate the cubic equation at x.
Definition: cubicEqnI.H:105
scalar c() const
Definition: cubicEqnI.H:69
scalar derivative(const scalar x) const
Evaluate the derivative of the cubic equation at x.
Definition: cubicEqnI.H:111
scalar a() const
Definition: cubicEqnI.H:57
scalar b() const
Definition: cubicEqnI.H:63
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:77
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
Namespace for OpenFOAM.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
const dimensionedScalar & D
volScalarField & b
Definition: createFields.H:27