yade.math module

This python module exposes all C++ math functions for Real and Complex types to python. In fact it sort of duplicates import math, import cmath or import mpmath. Also it facilitates migration of old python scripts to high precision calculations.

This module has following purposes:

  1. To reliably test all C++ math functions of arbitrary precision Real and Complex types against mpmath.
  2. To act as a “migration helper” for python scripts which call python mathematical functions that do not work well with mpmath. As an example see math.linspace below and this merge request
  3. To allow writing python math code in a way that mirrors C++ math code in Yade. As a bonus it will be faster than mpmath because mpmath is a purely python library (which was one of the main difficulties when writing lib/high-precision/ToFromPythonConverter.hpp)
  4. To test Eigen NumTraits
  5. To test CGAL NumTraits

If another C++ math function is needed it should be added to following files:

  1. lib/high-precision/MathFunctions.hpp
  2. py/high-precision/_math.cpp
  3. py/tests/testMath.py
  4. py/tests/testMathHelper.py

If another python math function does not work well with mpmath it should be added below, and original calls to this function should call this function instead, e.g. numpy.linspace(…) is replaced with yade.math.linspace(…).

The RealHP<n> higher precision math functions can be accessed in python by using the .HPn module scope. For example:

import yade.math as mth
mth.HP2.sqrt(2) # produces square root of 2 using RealHP<2> precision
mth.sqrt(2)     # without using HPn module scope it defaults to RealHP<1>
yade.math.Real(arg)

This function is for compatibility of calls like: g = yade.math.toHP1("-9.81"). If yade is compiled with default Real precision set as double, then python won’t accept string arguments as numbers. However when using higher precisions only calls yade.math.toHP1("1.234567890123456789012345678901234567890") do not cut to the first 15 decimal places. The calls such as yade.math.toHP1(1.234567890123456789012345678901234567890) will use default pythondouble conversion and will cut the number to its first 15 digits.

If you are debugging a high precision python script, and have difficulty finding places where such cuts have happened you should use yade.math.toHP1(string) for declaring all python floating point numbers which are physically important in the simulation. This function will throw exception if bad conversion is about to take place.

Also see example high precision check checkGravityRungeKuttaCashKarp54.py.

yade.math.Real1(arg)

This function is for compatibility of calls like: g = yade.math.toHP1("-9.81"). If yade is compiled with default Real precision set as double, then python won’t accept string arguments as numbers. However when using higher precisions only calls yade.math.toHP1("1.234567890123456789012345678901234567890") do not cut to the first 15 decimal places. The calls such as yade.math.toHP1(1.234567890123456789012345678901234567890) will use default pythondouble conversion and will cut the number to its first 15 digits.

If you are debugging a high precision python script, and have difficulty finding places where such cuts have happened you should use yade.math.toHP1(string) for declaring all python floating point numbers which are physically important in the simulation. This function will throw exception if bad conversion is about to take place.

Also see example high precision check checkGravityRungeKuttaCashKarp54.py.

yade.math.degrees(arg)
Returns:arg in radians converted to degrees, using yade.math.Real precision.
yade.math.degreesHP1(arg)[source]
Returns:arg in radians converted to degrees, using yade.math.Real precision.
yade.math.eig(a)[source]

This function calls numpy.linalg.eig(…) or mpmath.eig(…), because numpy.linalg.eig function does not work with mpmath.

yade.math.getRealHPCppDigits10()[source]
Returns:tuple containing amount of decimal digits supported on C++ side by Eigen and CGAL.
yade.math.getRealHPPythonDigits10()[source]
Returns:tuple containing amount of decimal digits supported on python side by yade.minieigenHP.
yade.math.linspace(a, b, num)[source]

This function calls numpy.linspace(…) or mpmath.linspace(…), because numpy.linspace function does not work with mpmath.

yade.math.needsMpmathAtN(N)[source]
Parameters:N – The int N value of RealHP<N> in question. Must be N >= 1.
Returns:True or False with information if using mpmath is necessary to avoid losing precision when working with RealHP<N>.
yade.math.radians(arg)

The default python function import math ; math.radians(arg) only works on 15 digit double precision. If you want to carry on calculations in higher precision it is advisable to use this function yade.math.radiansHP1(arg) instead. It uses full yade Real precision numbers.

NOTE: in the future this function may replace radians(…) function which is called in yade in many scripts, and which in fact is a call to native python math.radians. We only need to find the best backward compatible approach for this. The function yade.math.radiansHP1(arg) will remain as the function which uses native yade Real precision.

yade.math.radiansHP1(arg)[source]

The default python function import math ; math.radians(arg) only works on 15 digit double precision. If you want to carry on calculations in higher precision it is advisable to use this function yade.math.radiansHP1(arg) instead. It uses full yade Real precision numbers.

NOTE: in the future this function may replace radians(…) function which is called in yade in many scripts, and which in fact is a call to native python math.radians. We only need to find the best backward compatible approach for this. The function yade.math.radiansHP1(arg) will remain as the function which uses native yade Real precision.

yade.math.toHP1(arg)[source]

This function is for compatibility of calls like: g = yade.math.toHP1("-9.81"). If yade is compiled with default Real precision set as double, then python won’t accept string arguments as numbers. However when using higher precisions only calls yade.math.toHP1("1.234567890123456789012345678901234567890") do not cut to the first 15 decimal places. The calls such as yade.math.toHP1(1.234567890123456789012345678901234567890) will use default pythondouble conversion and will cut the number to its first 15 digits.

If you are debugging a high precision python script, and have difficulty finding places where such cuts have happened you should use yade.math.toHP1(string) for declaring all python floating point numbers which are physically important in the simulation. This function will throw exception if bad conversion is about to take place.

Also see example high precision check checkGravityRungeKuttaCashKarp54.py.

yade.math.usesHP()[source]
Returns:True if yade is using default Real precision higher than 15 digit (53 bits) double type.
yade._math.CGAL_Is_finite((float)x) → bool

CGAL’s function Is_finite, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:bool indicating if the Real argument is finite.
yade._math.CGAL_Is_valid((float)x) → bool

CGAL’s function Is_valid, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:bool indicating if the Real argument is valid. Checks are performed against NaN and Inf.
yade._math.CGAL_Kth_root((int)arg1, (float)x) → float

CGAL’s function Kth_root, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:Real the k-th root of argument.
yade._math.CGAL_Sgn((float)x) → int

CGAL’s function Sgn, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:sign of the argument, can be -1, 0 or 1. Not very useful in python. In C++ it is useful to obtain a sign of an expression with exact accuracy, CGAL starts using MPFR internally for this when the approximate interval contains zero inside it.
yade._math.CGAL_Sqrt((float)x) → float

CGAL’s function Sqrt, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:Real the square root of argument.
yade._math.CGAL_Square((float)x) → float

CGAL’s function Square, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:Real the argument squared.
yade._math.CGAL_To_interval((float)x) → tuple

CGAL’s function To_interval, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:(double,double) tuple inside which the high-precision Real argument resides.
yade._math.CGAL_simpleTest() → float

Tests a simple CGAL calculation. Distance between plane and point, uses CGAL’s sqrt and pow.

Returns:3.0
yade._math.Catalan([(int)Precision=53]) → float
Returns:Real The catalan constant, exposed to python for testing of eigen numerical traits.
yade._math.Euler([(int)Precision=53]) → float
Returns:Real The Euler–Mascheroni constant, exposed to python for testing of eigen numerical traits.
class yade._math.HP1
AddCost = 1
CGAL_Is_finite((float)x) → bool :

CGAL’s function Is_finite, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:bool indicating if the Real argument is finite.
CGAL_Is_valid((float)x) → bool :

CGAL’s function Is_valid, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:bool indicating if the Real argument is valid. Checks are performed against NaN and Inf.
CGAL_Kth_root((int)arg1, (float)x) → float :

CGAL’s function Kth_root, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:Real the k-th root of argument.
CGAL_Sgn((float)x) → int :

CGAL’s function Sgn, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:sign of the argument, can be -1, 0 or 1. Not very useful in python. In C++ it is useful to obtain a sign of an expression with exact accuracy, CGAL starts using MPFR internally for this when the approximate interval contains zero inside it.
CGAL_Sqrt((float)x) → float :

CGAL’s function Sqrt, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:Real the square root of argument.
CGAL_Square((float)x) → float :

CGAL’s function Square, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:Real the argument squared.
CGAL_To_interval((float)x) → tuple :

CGAL’s function To_interval, as described in CGAL algebraic foundations exposed to python for testing of CGAL numerical traits.

Returns:(double,double) tuple inside which the high-precision Real argument resides.
CGAL_simpleTest() → float :

Tests a simple CGAL calculation. Distance between plane and point, uses CGAL’s sqrt and pow.

Returns:3.0
Catalan([(int)Precision=53]) → float :
Returns:Real The catalan constant, exposed to python for testing of eigen numerical traits.
ComplexAddCost = 2
ComplexMulCost = 6
ComplexReadCost = 2
Euler([(int)Precision=53]) → float :
Returns:Real The Euler–Mascheroni constant, exposed to python for testing of eigen numerical traits.
IsComplex = 0
IsInteger = 0
IsSigned = 1
Log2([(int)Precision=53]) → float :
Returns:Real natural logarithm of 2, exposed to python for testing of eigen numerical traits.
MulCost = 1
Pi([(int)Precision=53]) → float :
Returns:Real The π constant, exposed to python for testing of eigen numerical traits.
ReadCost = 1
RequireInitialization = 0
class Var

The Var class is used to test to/from python converters for arbitrary precision Real

cpl

one Complex variable to test reading from and writing to it.

val

one Real variable for testing.

abs((complex)x) → float :
Returns:the Real absolute value of the Complex argument. Depending on compilation options wraps ::boost::multiprecision::abs(…) or std::abs(…) function.
abs( (float)x) -> float :
return:the Real absolute value of the Real argument. Depending on compilation options wraps ::boost::multiprecision::abs(…) or std::abs(…) function.
acos((complex)x) → complex :
Returns:Complex the arc-cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::acos(…) or std::acos(…) function.
acos( (float)x) -> float :
return:Real the arcus cosine of the argument. Depending on compilation options wraps ::boost::multiprecision::acos(…) or std::acos(…) function.
acosh((complex)x) → complex :
Returns:Complex the arc-hyperbolic cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::acosh(…) or std::acosh(…) function.
acosh( (float)x) -> float :
return:Real the hyperbolic arcus cosine of the argument. Depending on compilation options wraps ::boost::multiprecision::acosh(…) or std::acosh(…) function.
arg((complex)x) → float :
Returns:Real the arg (Phase angle of complex in radians) of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::arg(…) or std::arg(…) function.
asin((complex)x) → complex :
Returns:Complex the arc-sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::asin(…) or std::asin(…) function.
asin( (float)x) -> float :
return:Real the arcus sine of the argument. Depending on compilation options wraps ::boost::multiprecision::asin(…) or std::asin(…) function.
asinh((complex)x) → complex :
Returns:Complex the arc-hyperbolic sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::asinh(…) or std::asinh(…) function.
asinh( (float)x) -> float :
return:Real the hyperbolic arcus sine of the argument. Depending on compilation options wraps ::boost::multiprecision::asinh(…) or std::asinh(…) function.
atan((complex)x) → complex :
Returns:Complex the arc-tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::atan(…) or std::atan(…) function.
atan( (float)x) -> float :
return:Real the arcus tangent of the argument. Depending on compilation options wraps ::boost::multiprecision::atan(…) or std::atan(…) function.
atan2((float)x, (float)y) → float :
Returns:Real the arc tangent of y/x using the signs of the arguments x and y to determine the correct quadrant. Depending on compilation options wraps ::boost::multiprecision::atan2(…) or std::atan2(…) function.
atanh((complex)x) → complex :
Returns:Complex the arc-hyperbolic tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::atanh(…) or std::atanh(…) function.
atanh( (float)x) -> float :
return:Real the hyperbolic arcus tangent of the argument. Depending on compilation options wraps ::boost::multiprecision::atanh(…) or std::atanh(…) function.
cbrt((float)x) → float :
Returns:Real cubic root of the argument. Depending on compilation options wraps ::boost::multiprecision::cbrt(…) or std::cbrt(…) function.
ceil((float)x) → float :
Returns:Real Computes the smallest integer value not less than arg. Depending on compilation options wraps ::boost::multiprecision::ceil(…) or std::ceil(…) function.
conj((complex)x) → complex :
Returns:the complex conjugation a Complex argument. Depending on compilation options wraps ::boost::multiprecision::conj(…) or std::conj(…) function.
cos((complex)x) → complex :
Returns:Complex the cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::cos(…) or std::cos(…) function.
cos( (float)x) -> float :
return:Real the cosine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::cos(…) or std::cos(…) function.
cosh((complex)x) → complex :
Returns:Complex the hyperbolic cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::cosh(…) or std::cosh(…) function.
cosh( (float)x) -> float :
return:Real the hyperbolic cosine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::cosh(…) or std::cosh(…) function.
cylBesselJ((int)k, (float)x) → float :
Returns:Real the Bessel Functions of the First Kind of the order k and the Real argument. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/bessel/bessel_first.html>`__
defprec = 53
dummy_precision() → float :
Returns:similar to the function epsilon, but assumes that last 10% of bits contain the numerical error only. This is sometimes used by Eigen when calling isEqualFuzzy to determine if values differ a lot or if they are vaguely close to each other.
epsilon([(int)Precision=53]) → float :
Returns:Real returns the difference between 1.0 and the next representable value of the Real type. Wraps std::numeric_limits<Real>::epsilon() function.
epsilon( (float)x) -> float :
return:Real returns the difference between 1.0 and the next representable value of the Real type. Wraps std::numeric_limits<Real>::epsilon() function.
erf((float)x) → float :
Returns:Real Computes the error function of argument. Depending on compilation options wraps ::boost::multiprecision::erf(…) or std::erf(…) function.
erfc((float)x) → float :
Returns:Real Computes the complementary error function of argument, that is 1.0-erf(arg). Depending on compilation options wraps ::boost::multiprecision::erfc(…) or std::erfc(…) function.
exp((complex)x) → complex :
Returns:the base e exponential of a Complex argument. Depending on compilation options wraps ::boost::multiprecision::exp(…) or std::exp(…) function.
exp( (float)x) -> float :
return:the base e exponential of a Real argument. Depending on compilation options wraps ::boost::multiprecision::exp(…) or std::exp(…) function.
exp2((float)x) → float :
Returns:the base 2 exponential of a Real argument. Depending on compilation options wraps ::boost::multiprecision::exp2(…) or std::exp2(…) function.
expm1((float)x) → float :
Returns:the base e exponential of a Real argument minus 1.0. Depending on compilation options wraps ::boost::multiprecision::expm1(…) or std::expm1(…) function.
fabs((float)x) → float :
Returns:the Real absolute value of the argument. Depending on compilation options wraps ::boost::multiprecision::abs(…) or std::abs(…) function.
factorial((int)x) → float :
Returns:Real the factorial of the Real argument. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/factorials/sf_factorial.html>`__
floor((float)x) → float :
Returns:Real Computes the largest integer value not greater than arg. Depending on compilation options wraps ::boost::multiprecision::floor(…) or std::floor(…) function.
fma((float)x, (float)y, (float)z) → float :
Returns:Real - computes (x*y) + z as if to infinite precision and rounded only once to fit the result type. Depending on compilation options wraps ::boost::multiprecision::fma(…) or std::fma(…) function.
fmod((float)x, (float)y) → float :
Returns:Real the floating-point remainder of the division operation x/y of the arguments x and y. Depending on compilation options wraps ::boost::multiprecision::fmod(…) or std::fmod(…) function.
frexp((float)x) → tuple :
Returns:tuple of (Real,int), decomposes given floating point Real argument into a normalized fraction and an integral power of two. Depending on compilation options wraps ::boost::multiprecision::frexp(…) or std::frexp(…) function.
fromBits((str)bits[, (int)exp=0[, (int)sign=1]]) → float :
Parameters:
  • bitsstr - a string containing ‘0’, ‘1’ characters.
  • expint - the binary exponent which shifts the bits.
  • signint - the sign, should be -1 or +1, but it is not checked. It multiplies the result when construction from bits is finished.
Returns:

RealHP<N> constructed from string containing ‘0’, ‘1’ bits. This is for debugging purposes, rather slow.

getDecomposedReal((float)x) → dict :
Returns:dict - the dictionary with the debug information how the DecomposedReal class sees this type. This is for debugging purposes, rather slow. Includes result from fpclassify function call, a binary representation and other useful info. See also fromBits.
getDemangledName() → str :
Returns:string - the demangled C++ typnename of RealHP<N>.
getDemangledNameComplex() → str :
Returns:string - the demangled C++ typnename of ComplexHP<N>.
getFloatDistanceULP((float)arg1, (float)arg2) → float :
Returns:an integer value stored in RealHP<N>, the ULP distance calculated by boost::math::float_distance, also see Floating-point Comparison and Prof. Kahan paper about this topic.

Warning

The returned value is the directed distance between two arguments, this means that it can be negative.

getRawBits((float)x) → str :
Returns:string - the raw bits in memory representing this type. Be careful: it only checks the system endianness and either prints bytes in reverse order or not. Does not make any attempts to further interpret the bits of: sign, exponent or significand (on a typical x86 processor they are printed in that order), and different processors might store them differently. It is not useful for types which internally use a pointer because for them this function prints not the floating point number but a pointer. This is for debugging purposes.
hasInfinityNan = True
highest([(int)Precision=53]) → float :
Returns:Real returns the largest finite value of the Real type. Wraps std::numeric_limits<Real>::max() function.
hypot((float)x, (float)y) → float :
Returns:Real the square root of the sum of the squares of x and y, without undue overflow or underflow at intermediate stages of the computation. Depending on compilation options wraps ::boost::multiprecision::hypot(…) or std::hypot(…) function.
ilogb((float)x) → float :
Returns:Real extracts the value of the unbiased exponent from the floating-point argument arg, and returns it as a signed integer value. Depending on compilation options wraps ::boost::multiprecision::ilogb(…) or std::ilogb(…) function.
imag((complex)x) → float :
Returns:the imag part of a Complex argument. Depending on compilation options wraps ::boost::multiprecision::imag(…) or std::imag(…) function.
isApprox((float)a, (float)b, (float)eps) → bool :
Returns:bool, True if a is approximately equal b with provided eps, see also here
isApproxOrLessThan((float)a, (float)b, (float)eps) → bool :
Returns:bool, True if a is approximately less than or equal b with provided eps, see also here
isEqualFuzzy((float)arg1, (float)arg2, (float)arg3) → bool :
Returns:bool, True if the absolute difference between two numbers is smaller than std::numeric_limits<Real>::epsilon()
isMuchSmallerThan((float)a, (float)b, (float)eps) → bool :
Returns:bool, True if a is less than b with provided eps, see also here
isfinite((float)x) → bool :
Returns:bool indicating if the Real argument is Inf. Depending on compilation options wraps ::boost::multiprecision::isfinite(…) or std::isfinite(…) function.
isinf((float)x) → bool :
Returns:bool indicating if the Real argument is Inf. Depending on compilation options wraps ::boost::multiprecision::isinf(…) or std::isinf(…) function.
isnan((float)x) → bool :
Returns:bool indicating if the Real argument is NaN. Depending on compilation options wraps ::boost::multiprecision::isnan(…) or std::isnan(…) function.
laguerre((int)n, (int)m, (float)x) → float :
Returns:Real the Laguerre polynomial of the orders n, m and the Real argument. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/sf_poly/laguerre.html>`__
ldexp((float)x, (int)y) → float :
Returns:Multiplies a floating point value x by the number 2 raised to the exp power. Depending on compilation options wraps ::boost::multiprecision::ldexp(…) or std::ldexp(…) function.
lgamma((float)x) → float :
Returns:Real Computes the natural logarithm of the absolute value of the gamma function of arg. Depending on compilation options wraps ::boost::multiprecision::lgamma(…) or std::lgamma(…) function.
log((complex)x) → complex :
Returns:the Complex natural (base e) logarithm of a complex value z with a branch cut along the negative real axis. Depending on compilation options wraps ::boost::multiprecision::log(…) or std::log(…) function.
log( (float)x) -> float :
return:the Real natural (base e) logarithm of a real value. Depending on compilation options wraps ::boost::multiprecision::log(…) or std::log(…) function.
log10((complex)x) → complex :
Returns:the Complex (base 10) logarithm of a complex value z with a branch cut along the negative real axis. Depending on compilation options wraps ::boost::multiprecision::log10(…) or std::log10(…) function.
log10( (float)x) -> float :
return:the Real decimal (base 10) logarithm of a real value. Depending on compilation options wraps ::boost::multiprecision::log10(…) or std::log10(…) function.
log1p((float)x) → float :
Returns:the Real natural (base e) logarithm of 1+argument. Depending on compilation options wraps ::boost::multiprecision::log1p(…) or std::log1p(…) function.
log2((float)x) → float :
Returns:the Real binary (base 2) logarithm of a real value. Depending on compilation options wraps ::boost::multiprecision::log2(…) or std::log2(…) function.
logb((float)x) → float :
Returns:Extracts the value of the unbiased radix-independent exponent from the floating-point argument arg, and returns it as a floating-point value. Depending on compilation options wraps ::boost::multiprecision::logb(…) or std::logb(…) function.
lowest([(int)Precision=53]) → float :
Returns:Real returns the lowest (negative) finite value of the Real type. Wraps std::numeric_limits<Real>::lowest() function.
max((float)x, (float)y) → float :
Returns:Real larger of the two arguments. Depending on compilation options wraps ::boost::multiprecision::max(…) or std::max(…) function.
max_exp2 = 1024
min((float)x, (float)y) → float :
Returns:Real smaller of the two arguments. Depending on compilation options wraps ::boost::multiprecision::min(…) or std::min(…) function.
modf((float)x) → tuple :
Returns:tuple of (Real,Real), decomposes given floating point Real into integral and fractional parts, each having the same type and sign as x. Depending on compilation options wraps ::boost::multiprecision::modf(…) or std::modf(…) function.
polar((float)x, (float)y) → complex :
Returns:Complex the polar (Complex from polar components) of the Real rho (length), Real theta (angle) arguments in radians. Depending on compilation options wraps ::boost::multiprecision::polar(…) or std::polar(…) function.
pow((complex)x, (complex)pow) → complex :
Returns:the Complex complex arg1 raised to the Complex power arg2. Depending on compilation options wraps ::boost::multiprecision::pow(…) or std::pow(…) function.
pow( (float)x, (float)y) -> float :
return:Real the value of base raised to the power exp. Depending on compilation options wraps ::boost::multiprecision::pow(…) or std::pow(…) function.
proj((complex)x) → complex :
Returns:Complex the proj (projection of the complex number onto the Riemann sphere) of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::proj(…) or std::proj(…) function.
random() → float :
Returns:Real a symmetric random number in interval (-1,1). Used by Eigen.
random( (float)a, (float)b) -> float :
return:Real a random number in interval (a,b). Used by Eigen.
real((complex)x) → float :
Returns:the real part of a Complex argument. Depending on compilation options wraps ::boost::multiprecision::real(…) or std::real(…) function.
remainder((float)x, (float)y) → float :
Returns:Real the IEEE remainder of the floating point division operation x/y. Depending on compilation options wraps ::boost::multiprecision::remainder(…) or std::remainder(…) function.
remquo((float)x, (float)y) → tuple :
Returns:tuple of (Real,long), the floating-point remainder of the division operation x/y as the std::remainder() function does. Additionally, the sign and at least the three of the last bits of x/y are returned, sufficient to determine the octant of the result within a period. Depending on compilation options wraps ::boost::multiprecision::remquo(…) or std::remquo(…) function.
rint((float)x) → float :
Returns:Rounds the floating-point argument arg to an integer value (in floating-point format), using the current rounding mode. Depending on compilation options wraps ::boost::multiprecision::rint(…) or std::rint(…) function.
round((float)x) → float :
Returns:Real the nearest integer value to arg (in floating-point format), rounding halfway cases away from zero, regardless of the current rounding mode.. Depending on compilation options wraps ::boost::multiprecision::round(…) or std::round(…) function.
roundTrip((float)x) → float :
Returns:Real returns the argument x. Can be used to convert type to native RealHP<N> accuracy.
sgn((float)x) → int :
Returns:int the sign of the argument: -1, 0 or 1.
sign((float)x) → int :
Returns:int the sign of the argument: -1, 0 or 1.
sin((complex)x) → complex :
Returns:Complex the sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::sin(…) or std::sin(…) function.
sin( (float)x) -> float :
return:Real the sine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::sin(…) or std::sin(…) function.
sinh((complex)x) → complex :
Returns:Complex the hyperbolic sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::sinh(…) or std::sinh(…) function.
sinh( (float)x) -> float :
return:Real the hyperbolic sine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::sinh(…) or std::sinh(…) function.
smallest_positive() → float :
Returns:Real the smallest number greater than zero. Wraps std::numeric_limits<Real>::min()
sphericalHarmonic((int)l, (int)m, (float)theta, (float)phi) → complex :
Returns:Real the spherical harmonic polynomial of the orders l (unsigned int), m (signed int) and the Real arguments theta and phi. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/sf_poly/sph_harm.html>`__
sqrt((complex)x) → complex :
Returns:the Complex square root of Complex argument. Depending on compilation options wraps ::boost::multiprecision::sqrt(…) or std::sqrt(…) function.
sqrt( (float)x) -> float :
return:Real square root of the argument. Depending on compilation options wraps ::boost::multiprecision::sqrt(…) or std::sqrt(…) function.
squaredNorm((complex)x) → float :
Returns:Real the norm (squared magnitude) of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::norm(…) or std::norm(…) function.
tan((complex)x) → complex :
Returns:Complex the tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::tan(…) or std::tan(…) function.
tan( (float)x) -> float :
return:Real the tangent of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::tan(…) or std::tan(…) function.
tanh((complex)x) → complex :
Returns:Complex the hyperbolic tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::tanh(…) or std::tanh(…) function.
tanh( (float)x) -> float :
return:Real the hyperbolic tangent of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::tanh(…) or std::tanh(…) function.
testArray() → None :

This function tests call to std::vector::data(…) function in order to extract the array.

testCgalNumTraits = True
testConstants() → None :

This function tests lib/high-precision/Constants.hpp, the yade::math::ConstantsHP<N>, former yade::Mathr constants.

tgamma((float)x) → float :
Returns:Real Computes the gamma function of arg. Depending on compilation options wraps ::boost::multiprecision::tgamma(…) or std::tgamma(…) function.
toDouble((float)x) → float :
Returns:float converts Real type to double and returns a native python float.
toHP1((float)x) → float :
Returns:RealHP<1> converted from argument RealHP<1> as a result of static_cast<RealHP<1>>(arg).
toInt((float)x) → int :
Returns:int converts Real type to int and returns a native python int.
toLong((float)x) → int :
Returns:int converts Real type to long int and returns a native python int.
toLongDouble((float)x) → float :
Returns:float converts Real type to long double and returns a native python float.
trunc((float)x) → float :
Returns:Real the nearest integer not greater in magnitude than arg. Depending on compilation options wraps ::boost::multiprecision::trunc(…) or std::trunc(…) function.
yade._math.Log2([(int)Precision=53]) → float
Returns:Real natural logarithm of 2, exposed to python for testing of eigen numerical traits.
yade._math.Pi([(int)Precision=53]) → float
Returns:Real The π constant, exposed to python for testing of eigen numerical traits.
class yade._math.RealHPConfig

RealHPConfig class provides information about RealHP<N> type.

Variables:
  • extraStringDigits10 – this static variable allows to control how many extra digits to use when converting to decimal strings. Assign a different value to it to affect the string conversion done in C++ ↔ python conversions as well as in all other conversions. Be careful, because values smaller than 3 can fail the round trip conversion test.
  • isFloat128Broken – provides runtime information if Yade was compiled with g++ version < 9.2.1 and thus boost::multiprecision::float128 cannot work.
  • isEnabledRealHP – provides runtime information RealHP<N> is available for N higher than 1.
  • workaroundSlowBoostBinFloatboost::multiprecision::cpp_bin_float has some problem that importing it in python is very slow when these functions are exported: erf, erfc, lgamma, tgamma. In such case the python import yade.math can take more than minute. The workaround is to make them unavailable in python for higher N values. See invocation of IfConstexprForSlowFunctions in _math.cpp. This variable contains the highest N in which these functions are available. It equals to highest N when boost::multiprecision::cpp_bin_float is not used.
extraStringDigits10 = 4
getDigits10((int)N) → int :

This is a yade.math.RealHPConfig diagnostic function.

Parameters:Nint - the value of N in RealHP<N>.
Returns:the int representing std::numeric_limits<RealHP<N>>::digits10
getDigits2((int)N) → int :

This is a yade.math.RealHPConfig diagnostic function.

Parameters:Nint - the value of N in RealHP<N>.
Returns:the int representing std::numeric_limits<RealHP<N>>::digits, which corresponds to the number of significand bits used by this type.
getSupportedByEigenCgal() → tuple :
Returns:the tuple containing N from RealHP<N> precisions supported by Eigen and CGAL
getSupportedByMinieigen() → tuple :
Returns:the tuple containing N from RealHP<N> precisions supported by minieigenHP
isEnabledRealHP = False
isFloat128Broken = True
isFloat128Present = False
workaroundSlowBoostBinFloat = 1
class yade._math.Var

The Var class is used to test to/from python converters for arbitrary precision Real

cpl

one Complex variable to test reading from and writing to it.

val

one Real variable for testing.

yade._math.abs((complex)x) → float
return:the Real absolute value of the Complex argument. Depending on compilation options wraps ::boost::multiprecision::abs(…) or std::abs(…) function.
abs( (float)x) → float :
return:the Real absolute value of the Real argument. Depending on compilation options wraps ::boost::multiprecision::abs(…) or std::abs(…) function.
yade._math.acos((complex)x) → complex
return:Complex the arc-cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::acos(…) or std::acos(…) function.
acos( (float)x) → float :
return:Real the arcus cosine of the argument. Depending on compilation options wraps ::boost::multiprecision::acos(…) or std::acos(…) function.
yade._math.acosh((complex)x) → complex
return:Complex the arc-hyperbolic cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::acosh(…) or std::acosh(…) function.
acosh( (float)x) → float :
return:Real the hyperbolic arcus cosine of the argument. Depending on compilation options wraps ::boost::multiprecision::acosh(…) or std::acosh(…) function.
yade._math.arg((complex)x) → float
Returns:Real the arg (Phase angle of complex in radians) of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::arg(…) or std::arg(…) function.
yade._math.asin((complex)x) → complex
return:Complex the arc-sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::asin(…) or std::asin(…) function.
asin( (float)x) → float :
return:Real the arcus sine of the argument. Depending on compilation options wraps ::boost::multiprecision::asin(…) or std::asin(…) function.
yade._math.asinh((complex)x) → complex
return:Complex the arc-hyperbolic sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::asinh(…) or std::asinh(…) function.
asinh( (float)x) → float :
return:Real the hyperbolic arcus sine of the argument. Depending on compilation options wraps ::boost::multiprecision::asinh(…) or std::asinh(…) function.
yade._math.atan((complex)x) → complex
return:Complex the arc-tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::atan(…) or std::atan(…) function.
atan( (float)x) → float :
return:Real the arcus tangent of the argument. Depending on compilation options wraps ::boost::multiprecision::atan(…) or std::atan(…) function.
yade._math.atan2((float)x, (float)y) → float
Returns:Real the arc tangent of y/x using the signs of the arguments x and y to determine the correct quadrant. Depending on compilation options wraps ::boost::multiprecision::atan2(…) or std::atan2(…) function.
yade._math.atanh((complex)x) → complex
return:Complex the arc-hyperbolic tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::atanh(…) or std::atanh(…) function.
atanh( (float)x) → float :
return:Real the hyperbolic arcus tangent of the argument. Depending on compilation options wraps ::boost::multiprecision::atanh(…) or std::atanh(…) function.
yade._math.cbrt((float)x) → float
Returns:Real cubic root of the argument. Depending on compilation options wraps ::boost::multiprecision::cbrt(…) or std::cbrt(…) function.
yade._math.ceil((float)x) → float
Returns:Real Computes the smallest integer value not less than arg. Depending on compilation options wraps ::boost::multiprecision::ceil(…) or std::ceil(…) function.
yade._math.conj((complex)x) → complex
Returns:the complex conjugation a Complex argument. Depending on compilation options wraps ::boost::multiprecision::conj(…) or std::conj(…) function.
yade._math.cos((complex)x) → complex
return:Complex the cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::cos(…) or std::cos(…) function.
cos( (float)x) → float :
return:Real the cosine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::cos(…) or std::cos(…) function.
yade._math.cosh((complex)x) → complex
return:Complex the hyperbolic cosine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::cosh(…) or std::cosh(…) function.
cosh( (float)x) → float :
return:Real the hyperbolic cosine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::cosh(…) or std::cosh(…) function.
yade._math.cylBesselJ((int)k, (float)x) → float
Returns:Real the Bessel Functions of the First Kind of the order k and the Real argument. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/bessel/bessel_first.html>`__
yade._math.dummy_precision() → float
Returns:similar to the function epsilon, but assumes that last 10% of bits contain the numerical error only. This is sometimes used by Eigen when calling isEqualFuzzy to determine if values differ a lot or if they are vaguely close to each other.
yade._math.epsilon([(int)Precision=53]) → float
return:Real returns the difference between 1.0 and the next representable value of the Real type. Wraps std::numeric_limits<Real>::epsilon() function.
epsilon( (float)x) → float :
return:Real returns the difference between 1.0 and the next representable value of the Real type. Wraps std::numeric_limits<Real>::epsilon() function.
yade._math.erf((float)x) → float
Returns:Real Computes the error function of argument. Depending on compilation options wraps ::boost::multiprecision::erf(…) or std::erf(…) function.
yade._math.erfc((float)x) → float
Returns:Real Computes the complementary error function of argument, that is 1.0-erf(arg). Depending on compilation options wraps ::boost::multiprecision::erfc(…) or std::erfc(…) function.
yade._math.exp((complex)x) → complex
return:the base e exponential of a Complex argument. Depending on compilation options wraps ::boost::multiprecision::exp(…) or std::exp(…) function.
exp( (float)x) → float :
return:the base e exponential of a Real argument. Depending on compilation options wraps ::boost::multiprecision::exp(…) or std::exp(…) function.
yade._math.exp2((float)x) → float
Returns:the base 2 exponential of a Real argument. Depending on compilation options wraps ::boost::multiprecision::exp2(…) or std::exp2(…) function.
yade._math.expm1((float)x) → float
Returns:the base e exponential of a Real argument minus 1.0. Depending on compilation options wraps ::boost::multiprecision::expm1(…) or std::expm1(…) function.
yade._math.fabs((float)x) → float
Returns:the Real absolute value of the argument. Depending on compilation options wraps ::boost::multiprecision::abs(…) or std::abs(…) function.
yade._math.factorial((int)x) → float
Returns:Real the factorial of the Real argument. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/factorials/sf_factorial.html>`__
yade._math.floor((float)x) → float
Returns:Real Computes the largest integer value not greater than arg. Depending on compilation options wraps ::boost::multiprecision::floor(…) or std::floor(…) function.
yade._math.fma((float)x, (float)y, (float)z) → float
Returns:Real - computes (x*y) + z as if to infinite precision and rounded only once to fit the result type. Depending on compilation options wraps ::boost::multiprecision::fma(…) or std::fma(…) function.
yade._math.fmod((float)x, (float)y) → float
Returns:Real the floating-point remainder of the division operation x/y of the arguments x and y. Depending on compilation options wraps ::boost::multiprecision::fmod(…) or std::fmod(…) function.
yade._math.frexp((float)x) → tuple
Returns:tuple of (Real,int), decomposes given floating point Real argument into a normalized fraction and an integral power of two. Depending on compilation options wraps ::boost::multiprecision::frexp(…) or std::frexp(…) function.
yade._math.fromBits((str)bits[, (int)exp=0[, (int)sign=1]]) → float
Parameters:
  • bitsstr - a string containing ‘0’, ‘1’ characters.
  • expint - the binary exponent which shifts the bits.
  • signint - the sign, should be -1 or +1, but it is not checked. It multiplies the result when construction from bits is finished.
Returns:

RealHP<N> constructed from string containing ‘0’, ‘1’ bits. This is for debugging purposes, rather slow.

yade._math.getDecomposedReal((float)x) → dict
Returns:dict - the dictionary with the debug information how the DecomposedReal class sees this type. This is for debugging purposes, rather slow. Includes result from fpclassify function call, a binary representation and other useful info. See also fromBits.
yade._math.getDemangledName() → str
Returns:string - the demangled C++ typnename of RealHP<N>.
yade._math.getDemangledNameComplex() → str
Returns:string - the demangled C++ typnename of ComplexHP<N>.
yade._math.getEigenFlags() → dict
Returns:A python dictionary listing flags for all types, see: https://eigen.tuxfamily.org/dox/group__flags.html
yade._math.getEigenStorageOrders() → dict
Returns:A python dictionary listing options for all types, see: https://eigen.tuxfamily.org/dox/group__TopicStorageOrders.html
yade._math.getFloatDistanceULP((float)arg1, (float)arg2) → float
Returns:an integer value stored in RealHP<N>, the ULP distance calculated by boost::math::float_distance, also see Floating-point Comparison and Prof. Kahan paper about this topic.

The returned value is the directed distance between two arguments, this means that it can be negative.

yade._math.getRawBits((float)x) → str
Returns:string - the raw bits in memory representing this type. Be careful: it only checks the system endianness and either prints bytes in reverse order or not. Does not make any attempts to further interpret the bits of: sign, exponent or significand (on a typical x86 processor they are printed in that order), and different processors might store them differently. It is not useful for types which internally use a pointer because for them this function prints not the floating point number but a pointer. This is for debugging purposes.
yade._math.getRealHPErrors((list)testLevelsHP[, (int)testCount=10[, (float)minX=-10.0[, (float)maxX=10.0[, (bool)useRandomArgs=False[, (int)printEveryNth=1000[, (bool)collectArgs=False[, (bool)extraChecks=False]]]]]]]) → dict

Tests mathematical functions against the highest precision in argument testLevelsHP and returns the largest ULP distance found with getFloatDistanceULP. A testCount randomized tries with function arguments in range minX ... maxX are performed on the RealHP<N> types where N is from the list provided in testLevelsHP argument.

Parameters:
  • testLevelsHP – a list of int values consisting of high precision levels N (in RealHP<N>) for which the tests should be done. Must consist at least of two elements so that there is a higher precision type available against which to perform the tests.
  • testCountint - specifies how many randomized tests of each function to perform.
  • minXReal - start of the range in which the random arguments are generated.
  • maxXReal - end of that range.
  • useRandomArgs – If False (default) then minX ... maxX is divided into testCount equidistant points. If True then each call is a random number. This applies only to the first argument of a function, if a function takes more than one argument, then remaining arguments are random - 2D scans are not performed.
  • printEveryNth – will print using LOG_INFO the progress information every Nth step in the testCount loop. To see it e.g. start using yade -f6, also see logger documentation.
  • collectArgs – if True then in returned results will be a longer list of arguments that produce incorrect results.
  • extraChecks – will perform extra checks while executing this funcion. Useful only for debugging of getRealHPErrors.
Returns:

A python dictionary with the largest ULP distance to the correct function value. For each function name there is a dictionary consisting of: how many binary digits (bits) are in the tested RealHP<N> type, the worst arguments for this function, and the ULP distance to the reference value.

The returned ULP error is an absolute value, as opposed to getFloatDistanceULP which is signed.

yade._math.highest([(int)Precision=53]) → float
Returns:Real returns the largest finite value of the Real type. Wraps std::numeric_limits<Real>::max() function.
yade._math.hypot((float)x, (float)y) → float
Returns:Real the square root of the sum of the squares of x and y, without undue overflow or underflow at intermediate stages of the computation. Depending on compilation options wraps ::boost::multiprecision::hypot(…) or std::hypot(…) function.
yade._math.ilogb((float)x) → float
Returns:Real extracts the value of the unbiased exponent from the floating-point argument arg, and returns it as a signed integer value. Depending on compilation options wraps ::boost::multiprecision::ilogb(…) or std::ilogb(…) function.
yade._math.imag((complex)x) → float
Returns:the imag part of a Complex argument. Depending on compilation options wraps ::boost::multiprecision::imag(…) or std::imag(…) function.
yade._math.isApprox((float)a, (float)b, (float)eps) → bool
Returns:bool, True if a is approximately equal b with provided eps, see also here
yade._math.isApproxOrLessThan((float)a, (float)b, (float)eps) → bool
Returns:bool, True if a is approximately less than or equal b with provided eps, see also here
yade._math.isEqualFuzzy((float)arg1, (float)arg2, (float)arg3) → bool
Returns:bool, True if the absolute difference between two numbers is smaller than std::numeric_limits<Real>::epsilon()
yade._math.isMuchSmallerThan((float)a, (float)b, (float)eps) → bool
Returns:bool, True if a is less than b with provided eps, see also here
yade._math.isThisSystemLittleEndian() → bool
Returns:True if this system uses little endian architecture, False otherwise.
yade._math.isfinite((float)x) → bool
Returns:bool indicating if the Real argument is Inf. Depending on compilation options wraps ::boost::multiprecision::isfinite(…) or std::isfinite(…) function.
yade._math.isinf((float)x) → bool
Returns:bool indicating if the Real argument is Inf. Depending on compilation options wraps ::boost::multiprecision::isinf(…) or std::isinf(…) function.
yade._math.isnan((float)x) → bool
Returns:bool indicating if the Real argument is NaN. Depending on compilation options wraps ::boost::multiprecision::isnan(…) or std::isnan(…) function.
yade._math.laguerre((int)n, (int)m, (float)x) → float
Returns:Real the Laguerre polynomial of the orders n, m and the Real argument. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/sf_poly/laguerre.html>`__
yade._math.ldexp((float)x, (int)y) → float
Returns:Multiplies a floating point value x by the number 2 raised to the exp power. Depending on compilation options wraps ::boost::multiprecision::ldexp(…) or std::ldexp(…) function.
yade._math.lgamma((float)x) → float
Returns:Real Computes the natural logarithm of the absolute value of the gamma function of arg. Depending on compilation options wraps ::boost::multiprecision::lgamma(…) or std::lgamma(…) function.
yade._math.log((complex)x) → complex
return:the Complex natural (base e) logarithm of a complex value z with a branch cut along the negative real axis. Depending on compilation options wraps ::boost::multiprecision::log(…) or std::log(…) function.
log( (float)x) → float :
return:the Real natural (base e) logarithm of a real value. Depending on compilation options wraps ::boost::multiprecision::log(…) or std::log(…) function.
yade._math.log10((complex)x) → complex
return:the Complex (base 10) logarithm of a complex value z with a branch cut along the negative real axis. Depending on compilation options wraps ::boost::multiprecision::log10(…) or std::log10(…) function.
log10( (float)x) → float :
return:the Real decimal (base 10) logarithm of a real value. Depending on compilation options wraps ::boost::multiprecision::log10(…) or std::log10(…) function.
yade._math.log1p((float)x) → float
Returns:the Real natural (base e) logarithm of 1+argument. Depending on compilation options wraps ::boost::multiprecision::log1p(…) or std::log1p(…) function.
yade._math.log2((float)x) → float
Returns:the Real binary (base 2) logarithm of a real value. Depending on compilation options wraps ::boost::multiprecision::log2(…) or std::log2(…) function.
yade._math.logb((float)x) → float
Returns:Extracts the value of the unbiased radix-independent exponent from the floating-point argument arg, and returns it as a floating-point value. Depending on compilation options wraps ::boost::multiprecision::logb(…) or std::logb(…) function.
yade._math.lowest([(int)Precision=53]) → float
Returns:Real returns the lowest (negative) finite value of the Real type. Wraps std::numeric_limits<Real>::lowest() function.
yade._math.max((float)x, (float)y) → float
Returns:Real larger of the two arguments. Depending on compilation options wraps ::boost::multiprecision::max(…) or std::max(…) function.
yade._math.min((float)x, (float)y) → float
Returns:Real smaller of the two arguments. Depending on compilation options wraps ::boost::multiprecision::min(…) or std::min(…) function.
yade._math.modf((float)x) → tuple
Returns:tuple of (Real,Real), decomposes given floating point Real into integral and fractional parts, each having the same type and sign as x. Depending on compilation options wraps ::boost::multiprecision::modf(…) or std::modf(…) function.
yade._math.polar((float)x, (float)y) → complex
Returns:Complex the polar (Complex from polar components) of the Real rho (length), Real theta (angle) arguments in radians. Depending on compilation options wraps ::boost::multiprecision::polar(…) or std::polar(…) function.
yade._math.pow((complex)x, (complex)pow) → complex
return:the Complex complex arg1 raised to the Complex power arg2. Depending on compilation options wraps ::boost::multiprecision::pow(…) or std::pow(…) function.
pow( (float)x, (float)y) → float :
return:Real the value of base raised to the power exp. Depending on compilation options wraps ::boost::multiprecision::pow(…) or std::pow(…) function.
yade._math.proj((complex)x) → complex
Returns:Complex the proj (projection of the complex number onto the Riemann sphere) of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::proj(…) or std::proj(…) function.
yade._math.random() → float
return:Real a symmetric random number in interval (-1,1). Used by Eigen.
random( (float)a, (float)b) → float :
return:Real a random number in interval (a,b). Used by Eigen.
yade._math.real((complex)x) → float
Returns:the real part of a Complex argument. Depending on compilation options wraps ::boost::multiprecision::real(…) or std::real(…) function.
yade._math.remainder((float)x, (float)y) → float
Returns:Real the IEEE remainder of the floating point division operation x/y. Depending on compilation options wraps ::boost::multiprecision::remainder(…) or std::remainder(…) function.
yade._math.remquo((float)x, (float)y) → tuple
Returns:tuple of (Real,long), the floating-point remainder of the division operation x/y as the std::remainder() function does. Additionally, the sign and at least the three of the last bits of x/y are returned, sufficient to determine the octant of the result within a period. Depending on compilation options wraps ::boost::multiprecision::remquo(…) or std::remquo(…) function.
yade._math.rint((float)x) → float
Returns:Rounds the floating-point argument arg to an integer value (in floating-point format), using the current rounding mode. Depending on compilation options wraps ::boost::multiprecision::rint(…) or std::rint(…) function.
yade._math.round((float)x) → float
Returns:Real the nearest integer value to arg (in floating-point format), rounding halfway cases away from zero, regardless of the current rounding mode.. Depending on compilation options wraps ::boost::multiprecision::round(…) or std::round(…) function.
yade._math.roundTrip((float)x) → float
Returns:Real returns the argument x. Can be used to convert type to native RealHP<N> accuracy.
yade._math.sgn((float)x) → int
Returns:int the sign of the argument: -1, 0 or 1.
yade._math.sign((float)x) → int
Returns:int the sign of the argument: -1, 0 or 1.
yade._math.sin((complex)x) → complex
return:Complex the sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::sin(…) or std::sin(…) function.
sin( (float)x) → float :
return:Real the sine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::sin(…) or std::sin(…) function.
yade._math.sinh((complex)x) → complex
return:Complex the hyperbolic sine of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::sinh(…) or std::sinh(…) function.
sinh( (float)x) → float :
return:Real the hyperbolic sine of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::sinh(…) or std::sinh(…) function.
yade._math.smallest_positive() → float
Returns:Real the smallest number greater than zero. Wraps std::numeric_limits<Real>::min()
yade._math.sphericalHarmonic((int)l, (int)m, (float)theta, (float)phi) → complex
Returns:Real the spherical harmonic polynomial of the orders l (unsigned int), m (signed int) and the Real arguments theta and phi. See: <https://www.boost.org/doc/libs/1_77_0/libs/math/doc/html/math_toolkit/sf_poly/sph_harm.html>`__
yade._math.sqrt((complex)x) → complex
return:the Complex square root of Complex argument. Depending on compilation options wraps ::boost::multiprecision::sqrt(…) or std::sqrt(…) function.
sqrt( (float)x) → float :
return:Real square root of the argument. Depending on compilation options wraps ::boost::multiprecision::sqrt(…) or std::sqrt(…) function.
yade._math.squaredNorm((complex)x) → float
Returns:Real the norm (squared magnitude) of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::norm(…) or std::norm(…) function.
yade._math.tan((complex)x) → complex
return:Complex the tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::tan(…) or std::tan(…) function.
tan( (float)x) → float :
return:Real the tangent of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::tan(…) or std::tan(…) function.
yade._math.tanh((complex)x) → complex
return:Complex the hyperbolic tangent of the Complex argument in radians. Depending on compilation options wraps ::boost::multiprecision::tanh(…) or std::tanh(…) function.
tanh( (float)x) → float :
return:Real the hyperbolic tangent of the Real argument in radians. Depending on compilation options wraps ::boost::multiprecision::tanh(…) or std::tanh(…) function.
yade._math.testArray() → None

This function tests call to std::vector::data(…) function in order to extract the array.

yade._math.testConstants() → None

This function tests lib/high-precision/Constants.hpp, the yade::math::ConstantsHP<N>, former yade::Mathr constants.

yade._math.testLoopRealHP() → None

This function tests lib/high-precision/Constants.hpp, but the C++ side: all precisions, even those inaccessible from python

yade._math.tgamma((float)x) → float
Returns:Real Computes the gamma function of arg. Depending on compilation options wraps ::boost::multiprecision::tgamma(…) or std::tgamma(…) function.
yade._math.toDouble((float)x) → float
Returns:float converts Real type to double and returns a native python float.
yade._math.toHP1((float)x) → float
Returns:RealHP<1> converted from argument RealHP<1> as a result of static_cast<RealHP<1>>(arg).
yade._math.toInt((float)x) → int
Returns:int converts Real type to int and returns a native python int.
yade._math.toLong((float)x) → int
Returns:int converts Real type to long int and returns a native python int.
yade._math.toLongDouble((float)x) → float
Returns:float converts Real type to long double and returns a native python float.
yade._math.trunc((float)x) → float
Returns:Real the nearest integer not greater in magnitude than arg. Depending on compilation options wraps ::boost::multiprecision::trunc(…) or std::trunc(…) function.