Loading...
Searching...
No Matches
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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32:
33 VectorSpace<linearEqn, scalar, 2>(Foam::zero{})
34{}
35
36
37inline Foam::linearEqn::linearEqn(const scalar a, const scalar b)
38{
39 this->v_[A] = a;
40 this->v_[B] = b;
41}
42
43
44// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46inline Foam::scalar Foam::linearEqn::value(const scalar x) const
47{
48 return x*a() + b();
49}
50
52inline Foam::scalar Foam::linearEqn::derivative(const scalar x) const
53{
54 return a();
55}
56
58inline Foam::scalar Foam::linearEqn::error(const scalar x) const
59{
60 return SMALL*(mag(x*a()) + mag(b()));
61}
62
63
65{
66 const scalar a = this->a();
67 const scalar b = this->b();
68
69 if (mag(a) < VSMALL)
70 {
71 return Roots<1>(roots::nan, 0);
72 }
73 else if (mag(b/VGREAT) >= mag(a))
74 {
75 return Roots<1>(sign(a) == sign(b) ? roots::negInf : roots::posInf, 0);
76 }
77
78 return Roots<1>(roots::real, -b/a);
79}
80
81
82// ************************************************************************* //
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition Roots.H:71
linearEqn()=default
Default construct.
scalar value(const scalar x) const
Evaluate the linear equation at x.
Definition linearEqnI.H:39
scalar derivative(const scalar x) const
Evaluate the derivative of the linear equation at x.
Definition linearEqnI.H:45
scalar a() const noexcept
Definition linearEqn.H:90
Roots< 1 > roots() const
Return the real root of the linear equation.
Definition linearEqnI.H:57
scalar b() const noexcept
Definition linearEqn.H:91
scalar error(const scalar x) const
Estimate the error of evaluation of the linear equation at x.
Definition linearEqnI.H:51
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
Namespace for OpenFOAM.
dimensionedScalar sign(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
volScalarField & b