linearEqnI.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 -------------------------------------------------------------------------------
10 License
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<linearEqn, scalar, 2>(Foam::zero())
37 {}
38 
39 
40 inline Foam::linearEqn::linearEqn(const scalar a, const scalar b)
41 {
42  this->v_[A] = a;
43  this->v_[B] = b;
44 }
45 
46 
47 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48 
49 inline Foam::scalar Foam::linearEqn::a() const
50 {
51  return this->v_[A];
52 }
53 
54 
55 inline Foam::scalar Foam::linearEqn::b() const
56 {
57  return this->v_[B];
58 }
59 
60 
61 inline Foam::scalar& Foam::linearEqn::a()
62 {
63  return this->v_[A];
64 }
65 
66 
67 inline Foam::scalar& Foam::linearEqn::b()
68 {
69  return this->v_[B];
70 }
71 
72 
73 inline Foam::scalar Foam::linearEqn::value(const scalar x) const
74 {
75  return x*a() + b();
76 }
77 
78 
79 inline Foam::scalar Foam::linearEqn::derivative(const scalar x) const
80 {
81  return a();
82 }
83 
84 
85 inline Foam::scalar Foam::linearEqn::error(const scalar x) const
86 {
87  return SMALL*(mag(x*a()) + mag(b()));
88 }
89 
90 
92 {
93  /*
94 
95  This function solves a linear equation of the following form:
96 
97  a*x + b = 0
98  x + B = 0
99 
100  */
101 
102  const scalar a = this->a();
103  const scalar b = this->b();
104 
105  if (a == 0)
106  {
107  return Roots<1>(roots::nan, 0);
108  }
109 
110  if (mag(b/VGREAT) >= mag(a))
111  {
112  return Roots<1>(sign(a) == sign(b) ? roots::negInf : roots::posInf, 0);
113  }
114 
115  return Roots<1>(roots::real, - b/a);
116 }
117 
118 
119 // ************************************************************************* //
Foam::linearEqn::a
scalar a() const
Definition: linearEqnI.H:49
Foam::roots::posInf
Definition: Roots.H:58
Foam::roots::nan
Definition: Roots.H:60
Foam::linearEqn::derivative
scalar derivative(const scalar x) const
Evaluate the derivative at x.
Definition: linearEqnI.H:79
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
Foam::linearEqn::b
scalar b() const
Definition: linearEqnI.H:55
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::sign
dimensionedScalar sign(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:166
Foam::linearEqn
Linear equation of the form a*x + b = 0.
Definition: linearEqn.H:50
Foam::VectorSpace
Templated vector space.
Definition: VectorSpace.H:56
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::linearEqn::value
scalar value(const scalar x) const
Evaluate at x.
Definition: linearEqnI.H:73
Foam::linearEqn::linearEqn
linearEqn()
Construct null.
Definition: linearEqnI.H:30
Foam::roots::negInf
Definition: Roots.H:59
Foam::roots::real
Definition: Roots.H:56
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::linearEqn::roots
Roots< 1 > roots() const
Get the roots.
Definition: linearEqnI.H:91
Foam::linearEqn::error
scalar error(const scalar x) const
Estimate the error of evaluation at x.
Definition: linearEqnI.H:85
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::Roots
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition: Roots.H:70
Foam::zero
A class representing the concept of 0 (zero), which can be used to avoid manipulating objects that ar...
Definition: zero.H:61