dolfin/math¶
Documentation for C++ code found in dolfin/math/*.h
Contents
Functions¶
between¶
C++ documentation for between
from dolfin/math/basic.h
:
-
bool
dolfin::
between
(double x, std::pair<double, double> range)¶ Check whether x is between x0 and x1 (inclusive, to within DOLFIN_EPS)
Parameters: - x – (double) Value to check
- range – (std::pair<double, double>) Range to check
Returns: bool
ipow¶
C++ documentation for ipow
from dolfin/math/basic.h
:
-
std::size_t
dolfin::
ipow
(std::size_t a, std::size_t n)¶ Return a to the power n. NOTE: Overflow is not checked!
Parameters: - a – (std::size_t) Value
- n – (std::size_t) Power
Returns: std::size_t
near¶
C++ documentation for near
from dolfin/math/basic.h
:
-
bool
dolfin::
near
(double x, double x0, double eps = DOLFIN_EPS)¶ Check whether x is close to x0 (to within DOLFIN_EPS)
Parameters: - x – (double) First value
- x0 – (double) Second value
- eps – (double) Tolerance
Returns: bool
Variables¶
rand_seeded¶
C++ documentation for rand_seeded
from dolfin/math/basic.cpp
:
-
bool
dolfin::
rand_seeded
¶ Flag to determine whether to reseed
dolfin::rand()
. Normally on first call.
Classes¶
Lagrange¶
C++ documentation for Lagrange
from dolfin/math/Lagrange.h
:
-
class
dolfin::
Lagrange
¶ Lagrange
polynomial (basis) with given degree q determined by n = q + 1 nodal points. Example: q = 1 (n = 2)Lagrange
p(1); p.set(0, 0.0); p.set(1, 1.0); It is the callers responsibility that the points are distinct. This creates aLagrange
polynomial (actually twoLagrange
polynomials): p(0,x) = 1 - x (one at x = 0, zero at x = 1) p(1,x) = x (zero at x = 0, one at x = 1)-
dolfin::Lagrange::
Lagrange
(std::size_t q)¶ Constructor.
Parameters: q –
-
std::vector<double>
dolfin::Lagrange::
constants
¶
-
std::size_t
dolfin::Lagrange::
counter
¶
-
double
dolfin::Lagrange::
ddx
(std::size_t i, double x)¶ Return derivate of polynomial i at given point x
Parameters: - i – (std::size_t)
- x – (double)
-
std::size_t
dolfin::Lagrange::
degree
() const¶ Return degree
Returns: std::size_t
-
double
dolfin::Lagrange::
dqdx
(std::size_t i)¶ Return derivative q (a constant) of polynomial
Parameters: i – (std::size_t)
-
double
dolfin::Lagrange::
eval
(std::size_t i, double x)¶ Return value of polynomial i at given point x
Parameters: - i – (std::size_t)
- x – (double)
-
void
dolfin::Lagrange::
init
()¶
-
double
dolfin::Lagrange::
operator()
(std::size_t i, double x)¶ Return value of polynomial i at given point x
Parameters: - i – (std::size_t)
- x – (double)
-
double
dolfin::Lagrange::
point
(std::size_t i) const¶ Return point
Parameters: i – (std::size_t)
-
std::vector<double>
dolfin::Lagrange::
points
¶
-
void
dolfin::Lagrange::
set
(std::size_t i, double x)¶ Specify point
Parameters: - i – (std::size_t)
- x – (double)
-
std::size_t
dolfin::Lagrange::
size
() const¶ Return number of points
Returns: std::size_t
-
std::string
dolfin::Lagrange::
str
(bool verbose) const¶ Return informal string representation (pretty-print)
Parameters: verbose – (bool) Verbosity of output string
-
Legendre¶
C++ documentation for Legendre
from dolfin/math/Legendre.h
:
-
class
dolfin::
Legendre
¶ Interface for computing
Legendre
polynomials via Boost.-
double
dolfin::Legendre::
d2dx
(std::size_t n, double x)¶ Evaluate second derivative of polynomial of order n at point x
Parameters: - n – (std::size_t) Order
- x – (double)
Point
Returns: double
Legendre
polynomial 2nd derivative value at x
-
double