Fast Numerical Evaluation¶
For many applications such as numerical integration, differential equation approximation, plotting a 3d surface, optimization problems, monte-carlo simulations, etc., one wishes to pass around and evaluate a single algebraic expression many, many times at various floating point values. Doing this via recursive calls over a python representation of the object (even if Maxima or other outside packages are not involved) is extremely inefficient.
Up until now the solution has been to use lambda expressions, but this is neither intuitive, Sage-like, nor efficient (compared to operating on raw C doubles). This module provides a representation of algebraic expression in Reverse Polish Notation, and provides an efficient interpreter on C double values as a callable python object. It does what it can in C, and will call out to Python if necessary.
Essential to the understanding of this class is the distinction
between symbolic expressions and callable symbolic expressions (where
the latter binds argument names to argument positions). The
*vars
parameter passed around encapsulates this information.
See the function fast_float(f, *vars)
to create a fast-callable
version of f.
Note
Sage temporarily has two implementations of this functionality ; one in this file, which will probably be deprecated soon, and one in fast_callable.pyx. The following instructions are for the old implementation; you probably want to be looking at fast_callable.pyx instead.
To provide this interface for a class, implement fast_float_(self, *vars)
. The basic building blocks are
provided by the functions fast_float_constant
(returns a
constant function), fast_float_arg
(selects the n
-th value
when called with \ge_n
arguments), and fast_float_func
which
wraps a callable Python function. These may be combined with the
standard Python arithmetic operators, and support many of the basic
math functions such sqrt
, exp
, and trig functions.
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float
sage: f = fast_float(sqrt(x^7+1), 'x', old=True)
sage: f(1)
1.4142135623730951
sage: f.op_list()
['load 0', 'push 7.0', 'pow', 'push 1.0', 'add', 'call sqrt(1)']
To interpret that last line, we load argument 0 (x
in this case) onto
the stack, push the constant 2.0 onto the stack, call the pow function
(which takes 2 arguments from the stack), push the constant 1.0, add the
top two arguments of the stack, and then call sqrt.
Here we take sin
of the first argument and add it to f
:
sage: from sage.ext.fast_eval import fast_float_arg
sage: g = fast_float_arg(0).sin()
sage: (f+g).op_list()
['load 0', 'push 7.0', 'pow', 'push 1.0', 'add', 'call sqrt(1)', 'load 0', 'call sin(1)', 'add']
AUTHORS:
- Robert Bradshaw (2008-10): Initial version
-
class
sage.ext.fast_eval.
FastDoubleFunc
¶ Bases:
object
This class is for fast evaluation of algebraic expressions over the real numbers (e.g. for plotting). It represents an expression as a stack-based series of operations.
EXAMPLES:
sage: from sage.ext.fast_eval import FastDoubleFunc sage: f = FastDoubleFunc('const', 1.5) # the constant function sage: f() 1.5 sage: g = FastDoubleFunc('arg', 0) # the first argument sage: g(5) 5.0 sage: h = f+g sage: h(17) 18.5 sage: h = h.sin() sage: h(pi/2-1.5) 1.0 sage: h.is_pure_c() True sage: list(h) ['push 1.5', 'load 0', 'add', 'call sin(1)']
We can wrap Python functions too:
sage: h = FastDoubleFunc('callable', lambda x,y: x*x*x - y, g, f) sage: h(10) 998.5 sage: h.is_pure_c() False sage: list(h) ['load 0', 'push 1.5', 'py_call <function <lambda> at 0x...>(2)']
Here’s a more complicated expression:
sage: from sage.ext.fast_eval import fast_float_constant, fast_float_arg sage: a = fast_float_constant(1.5) sage: b = fast_float_constant(3.14) sage: c = fast_float_constant(7) sage: x = fast_float_arg(0) sage: y = fast_float_arg(1) sage: f = a*x^2 + b*x + c - y/sqrt(sin(y)^2+a) sage: f(2,3) 16.846610528508116 sage: f.max_height 4 sage: f.is_pure_c() True sage: list(f) ['push 1.5', 'load 0', 'dup', 'mul', 'mul', 'push 3.14', 'load 0', 'mul', 'add', 'push 7.0', 'add', 'load 1', 'load 1', 'call sin(1)', 'dup', 'mul', 'push 1.5', 'add', 'call sqrt(1)', 'div', 'sub']
AUTHORS:
- Robert Bradshaw
-
abs
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).abs() sage: f(3) 3.0
-
arccos
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).arccos() sage: f(sqrt(3)/2) 0.5235987755982989...
-
arccosh
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).arccosh() sage: f(cosh(5)) 5.0
-
arcsin
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).arcsin() sage: f(0.5) 0.523598775598298...
-
arcsinh
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).arcsinh() sage: f(sinh(5)) 5.0
-
arctan
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).arctan() sage: f(1) 0.785398163397448...
-
arctanh
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).arctanh() sage: abs(f(tanh(0.5)) - 0.5) < 0.0000001 True
-
ceil
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).ceil() sage: f(1.5) 2.0 sage: f(-1.5) -1.0
-
cos
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).cos() sage: f(0) 1.0
-
cosh
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).cosh() sage: f(log(2)) 1.25
-
cot
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).cot() sage: f(pi/4) 1.0...
-
csc
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).csc() sage: f(pi/2) 1.0
-
exp
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).exp() sage: f(1) 2.718281828459045... sage: f(100) 2.6881171418161356e+43
-
floor
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).floor() sage: f(11.5) 11.0 sage: f(-11.5) -12.0
-
is_pure_c
()¶ Returns
True
if this function can be evaluated without any python calls (at any level).EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_constant, fast_float_arg, fast_float_func sage: fast_float_constant(2).is_pure_c() True sage: fast_float_arg(2).sqrt().sin().is_pure_c() True sage: fast_float_func(lambda _: 2).is_pure_c() False
-
log
(base=None)¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).log() sage: f(2) 0.693147180559945... sage: f = fast_float_arg(0).log(2) sage: f(2) 1.0 sage: f = fast_float_arg(0).log(3) sage: f(9) 2.0...
-
max_height
¶
-
nargs
¶
-
nops
¶
-
op_list
()¶ Returns a list of string representations of the operations that make up this expression.
Python and C function calls may be only available by function pointer addresses.
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_constant, fast_float_arg sage: a = fast_float_constant(17) sage: x = fast_float_arg(0) sage: a.op_list() ['push 17.0'] sage: x.op_list() ['load 0'] sage: (a*x).op_list() ['push 17.0', 'load 0', 'mul'] sage: (a+a*x^2).sqrt().op_list() ['push 17.0', 'push 17.0', 'load 0', 'dup', 'mul', 'mul', 'add', 'call sqrt(1)']
-
python_calls
()¶ Returns a list of all python calls used by function.
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_func, fast_float_arg sage: x = fast_float_arg(0) sage: f = fast_float_func(hash, sqrt(x)) sage: f.op_list() ['load 0', 'call sqrt(1)', 'py_call <built-in function hash>(1)'] sage: f.python_calls() [<built-in function hash>]
-
sec
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).sec() sage: f(pi) -1.0
-
sin
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).sin() sage: f(pi/2) 1.0
-
sinh
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).sinh() sage: f(log(2)) 0.75
-
sqrt
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).sqrt() sage: f(4) 2.0
-
tan
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).tan() sage: f(pi/3) 1.73205080756887...
-
tanh
()¶ EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0).tanh() sage: f(0) 0.0
-
sage.ext.fast_eval.
fast_float
(f, old=None, expect_one_var=False, *vars)¶ Tries to create a function that evaluates f quickly using floating-point numbers, if possible. There are two implementations of fast_float in Sage; by default we use the newer, which is slightly faster on most tests.
On failure, returns the input unchanged.
INPUT:
f
– an expressionvars
– the names of the argumentsold
– use the original algorithm for fast_floatexpect_one_var
– don’t give deprecation warning ifvars
is omitted, as long as expression has only one var
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float sage: x,y = var('x,y') sage: f = fast_float(sqrt(x^2+y^2), 'x', 'y') sage: f(3,4) 5.0
Specifying the argument names is essential, as fast_float objects only distinguish between arguments by order.
sage: f = fast_float(x-y, 'x','y') sage: f(1,2) -1.0 sage: f = fast_float(x-y, 'y','x') sage: f(1,2) 1.0
-
sage.ext.fast_eval.
fast_float_arg
(n)¶ Return a fast-to-evaluate argument selector.
INPUT:
n
– the (zero-indexed) argument to select
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_arg sage: f = fast_float_arg(0) sage: f(1,2) 1.0 sage: f = fast_float_arg(1) sage: f(1,2) 2.0
This is all that goes on under the hood:
sage: fast_float_arg(10).op_list() ['load 10']
-
sage.ext.fast_eval.
fast_float_constant
(x)¶ Return a fast-to-evaluate constant function.
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_constant sage: f = fast_float_constant(-2.75) sage: f() -2.75
This is all that goes on under the hood:
sage: fast_float_constant(pi).op_list() ['push 3.1415926535...']
-
sage.ext.fast_eval.
fast_float_func
(f, *args)¶ Returns a wrapper around a python function.
INPUT:
f
– a callable python objectargs
– a list of FastDoubleFunc inputs
EXAMPLES:
sage: from sage.ext.fast_eval import fast_float_func, fast_float_arg sage: f = fast_float_arg(0) sage: g = fast_float_arg(1) sage: h = fast_float_func(lambda x,y: x-y, f, g) sage: h(5, 10) -5.0
This is all that goes on under the hood:
sage: h.op_list() ['load 0', 'load 1', 'py_call <function <lambda> at 0x...>(2)']
-
sage.ext.fast_eval.
is_fast_float
(x)¶