Title: | Interface to the 'SymEngine' Library |
---|---|
Description: | Provides an R interface to 'SymEngine' <https://github.com/symengine/>, a standalone 'C++' library for fast symbolic manipulation. The package has functionalities for symbolic computation like calculating exact mathematical expressions, solving systems of linear equations and code generation. |
Authors: | Jialin Ma [cre, aut], Isuru Fernando [aut], Xin Chen [aut] |
Maintainer: | Jialin Ma <[email protected]> |
License: | GPL (>= 2) |
Version: | 0.2.6 |
Built: | 2024-11-05 05:48:25 UTC |
Source: | https://github.com/symengine/symengine.r |
These are S4 methods defined for Basic
, VecBasic
and DenseMatrix
.
## S4 method for signature 'Basic,Basic' e1 == e2 ## S4 method for signature 'Basic,Basic' e1 != e2 ## S4 method for signature 'SymEngineDataType,SymEngineDataType' Arith(e1, e2) ## S4 method for signature 'SymEngineDataType,vector' Arith(e1, e2) ## S4 method for signature 'vector,SymEngineDataType' Arith(e1, e2) ## S4 method for signature 'SymEngineDataType,missing' e1 - e2 ## S4 method for signature 'SymEngineDataType,missing' e1 + e2 ## S4 method for signature 'DenseMatrix,DenseMatrix' x %*% y ## S4 method for signature 'VecBasic,VecBasic' x %*% y ## S4 method for signature 'DenseMatrix,VecBasic' x %*% y ## S4 method for signature 'DenseMatrix,vector' x %*% y ## S4 method for signature 'VecBasic,DenseMatrix' x %*% y ## S4 method for signature 'vector,DenseMatrix' x %*% y ## S4 method for signature 'SymEngineDataType' Math(x) ## S4 method for signature 'SymEngineDataType' sinpi(x) ## S4 method for signature 'SymEngineDataType' cospi(x) ## S4 method for signature 'SymEngineDataType' tanpi(x) ## S4 method for signature 'SymEngineDataType' log(x, base) ## S4 method for signature 'SymEngineDataType' log2(x) ## S4 method for signature 'SymEngineDataType' log10(x) ## S4 method for signature 'SymEngineDataType' log1p(x) ## S4 method for signature 'SymEngineDataType' expm1(x) ## S4 method for signature 'SymEngineDataType' sum(x, ..., na.rm = FALSE) ## S4 method for signature 'SymEngineDataType' prod(x, ..., na.rm = FALSE)
## S4 method for signature 'Basic,Basic' e1 == e2 ## S4 method for signature 'Basic,Basic' e1 != e2 ## S4 method for signature 'SymEngineDataType,SymEngineDataType' Arith(e1, e2) ## S4 method for signature 'SymEngineDataType,vector' Arith(e1, e2) ## S4 method for signature 'vector,SymEngineDataType' Arith(e1, e2) ## S4 method for signature 'SymEngineDataType,missing' e1 - e2 ## S4 method for signature 'SymEngineDataType,missing' e1 + e2 ## S4 method for signature 'DenseMatrix,DenseMatrix' x %*% y ## S4 method for signature 'VecBasic,VecBasic' x %*% y ## S4 method for signature 'DenseMatrix,VecBasic' x %*% y ## S4 method for signature 'DenseMatrix,vector' x %*% y ## S4 method for signature 'VecBasic,DenseMatrix' x %*% y ## S4 method for signature 'vector,DenseMatrix' x %*% y ## S4 method for signature 'SymEngineDataType' Math(x) ## S4 method for signature 'SymEngineDataType' sinpi(x) ## S4 method for signature 'SymEngineDataType' cospi(x) ## S4 method for signature 'SymEngineDataType' tanpi(x) ## S4 method for signature 'SymEngineDataType' log(x, base) ## S4 method for signature 'SymEngineDataType' log2(x) ## S4 method for signature 'SymEngineDataType' log10(x) ## S4 method for signature 'SymEngineDataType' log1p(x) ## S4 method for signature 'SymEngineDataType' expm1(x) ## S4 method for signature 'SymEngineDataType' sum(x, ..., na.rm = FALSE) ## S4 method for signature 'SymEngineDataType' prod(x, ..., na.rm = FALSE)
e1 , e2 , x , y , base , ...
|
Objects. |
na.rm |
Ignored |
==
and !=
will return a logical vector. Other
functions will return a Basic
, VecBasic
or DenseMatrix
.
Miscellaneous S4 methods defined for converting a Basic
or
VecBasic
object to R number/string/language object.
## S4 method for signature 'Basic' as.character(x) ## S4 method for signature 'Basic' as.numeric(x) ## S4 method for signature 'Basic' as.integer(x) ## S4 method for signature 'VecBasic' as.character(x) ## S4 method for signature 'VecBasic' as.numeric(x) ## S4 method for signature 'VecBasic' as.integer(x) as.language(x) ## S4 method for signature 'Basic' as.language(x)
## S4 method for signature 'Basic' as.character(x) ## S4 method for signature 'Basic' as.numeric(x) ## S4 method for signature 'Basic' as.integer(x) ## S4 method for signature 'VecBasic' as.character(x) ## S4 method for signature 'VecBasic' as.numeric(x) ## S4 method for signature 'VecBasic' as.integer(x) as.language(x) ## S4 method for signature 'Basic' as.language(x)
x |
The object to be converted. |
Same as default methods of these generics. as.language()
may return symbol
, integer
, double
or call
.
These are miscellaneous S3/S4 methods defined for DenseMatrix
class.
## S3 method for class 'DenseMatrix' as.matrix(x, ...) ## S4 method for signature 'DenseMatrix' dim(x) ## S4 replacement method for signature 'DenseMatrix' dim(x) <- value ## S4 replacement method for signature 'VecBasic' dim(x) <- value ## S4 replacement method for signature 'Basic' dim(x) <- value ## S4 replacement method for signature 'DenseMatrix' dimnames(x) <- value ## S4 method for signature 'DenseMatrix' dimnames(x) ## S4 method for signature 'DenseMatrix' length(x) ## S4 method for signature 'DenseMatrix,ANY' x[[i, j, ...]] ## S4 replacement method for signature 'DenseMatrix' x[[i, j, ...]] <- value ## S4 method for signature 'DenseMatrix' x[i, j, ..., drop = TRUE] ## S4 replacement method for signature 'DenseMatrix' x[i, j, ...] <- value
## S3 method for class 'DenseMatrix' as.matrix(x, ...) ## S4 method for signature 'DenseMatrix' dim(x) ## S4 replacement method for signature 'DenseMatrix' dim(x) <- value ## S4 replacement method for signature 'VecBasic' dim(x) <- value ## S4 replacement method for signature 'Basic' dim(x) <- value ## S4 replacement method for signature 'DenseMatrix' dimnames(x) <- value ## S4 method for signature 'DenseMatrix' dimnames(x) ## S4 method for signature 'DenseMatrix' length(x) ## S4 method for signature 'DenseMatrix,ANY' x[[i, j, ...]] ## S4 replacement method for signature 'DenseMatrix' x[[i, j, ...]] <- value ## S4 method for signature 'DenseMatrix' x[i, j, ..., drop = TRUE] ## S4 replacement method for signature 'DenseMatrix' x[i, j, ...] <- value
x |
A DenseMatrix object. |
i , j , value , ... , drop
|
Arguments for subsetting, assignment or replacing. |
Same or similar with the generics of these methods.
S3 methods of cbind
and rbind
defined for
DenseMatrix
and VecBasic
.
## S3 method for class 'SymEngineDataType' cbind(..., deparse.level) ## S3 method for class 'SymEngineDataType' rbind(..., deparse.level)
## S3 method for class 'SymEngineDataType' cbind(..., deparse.level) ## S3 method for class 'SymEngineDataType' rbind(..., deparse.level)
... |
DenseMatrix, VecBasic or R objects. |
deparse.level |
Not used. |
DenseMatrix
S4 object.
Generate C/MathML/LaTeX/JavaScript code string from a Basic
or VecBasic
object.
codegen(x, type = c("ccode", "mathml", "latex", "jscode"))
codegen(x, type = c("ccode", "mathml", "latex", "jscode"))
x |
A Basic or a VecBasic object. |
type |
One of "ccode", "mathml", "latex" and "jscode". |
A character vector.
S4 method of D
defined for Basic
. It returns
the derivative of expr
with regards to name
.
name
may be missing if there is only one symbol in
expr
.
## S4 method for signature 'SymEngineDataType' D(expr, name)
## S4 method for signature 'SymEngineDataType' D(expr, name)
expr |
A Basic object. |
name |
A character vector or a Basic object of type Symbol. |
Same type as expr
argument.
expr <- S(~ exp(x)) D(expr) == expr expr <- S(~ x^2 + 2*x + 1) D(expr)
expr <- S(~ exp(x)) D(expr) == expr expr <- S(~ x^2 + 2*x + 1) D(expr)
S4 method of det
defined for DenseMatrix
.
det(x, ...) ## S4 method for signature 'DenseMatrix' det(x, ...)
det(x, ...) ## S4 method for signature 'DenseMatrix' det(x, ...)
x |
A DenseMatrix object. |
... |
Unused. |
A Basic
object.
mat <- Matrix(LETTERS[1:9], 3) det(mat)
mat <- Matrix(LETTERS[1:9], 3) det(mat)
Construct DoubleVisitor
object from Basic
or VecBasic
and use it to numerically evaluate symbolic expressions.
DoubleVisitor( exprs, args, perform_cse = TRUE, llvm_opt_level = if (symengine_have_component("llvm")) 3L else -1L ) visitor_call(visitor, input, do_transpose = FALSE)
DoubleVisitor( exprs, args, perform_cse = TRUE, llvm_opt_level = if (symengine_have_component("llvm")) 3L else -1L ) visitor_call(visitor, input, do_transpose = FALSE)
exprs |
A Basic object or a VecBasic object to be evaluated. |
args |
A VecBasic object indicating order of input arguments. Can be missing. |
perform_cse |
Boolean. |
llvm_opt_level |
Integer. If negative, it will return a |
visitor |
A DoubleVisitor object. |
input |
A numeric matrix. Each row is input value for one argument. |
do_transpose |
Boolean. Matters when |
DoubleVisitor
constructs the visitor and visitor itself is callable.
visitor_call
is the low level function to call the visitor with input.
DoubleVisitor
returns a callable LambdaDoubleVisitor
or
LLVMDoubleVisitor
. visitor_call
returns a numeric vector or matrix.
a <- S("a") b <- S("b") c <- S("c") vec <- c(log(a), log(a)/log(b) + c) func <- DoubleVisitor(vec, args = c(a, b, c)) args(func) ## Use closure func(a = 1:10, b = 10:1, c = 1.43) ## Use visitor_call input <- rbind(a = 1:10, b = 10:1, c = 1.43) visitor_call(func, input, do_transpose = TRUE)
a <- S("a") b <- S("b") c <- S("c") vec <- c(log(a), log(a)/log(b) + c) func <- DoubleVisitor(vec, args = c(a, b, c)) args(func) ## Use closure func(a = 1:10, b = 10:1, c = 1.43) ## Use visitor_call input <- rbind(a = 1:10, b = 10:1, c = 1.43) visitor_call(func, input, do_transpose = TRUE)
This is a wrapper of the odeintr
R package using
symengine objects to specify the ODE system and C code
generation functionality from symengine to generate the
C++ source. The dxdt
function and defined ==
S4 method
allow one to intuitively specify the ODE system with symengine
objects. The ODESystem
will generate C++ source
and compile on the fly with Rcpp. Then predict
can be used to get
results.
dxdt(x) ## S4 method for signature 'DxdtOdeConstructor,ANY' e1 == e2 ODESystem( odesys, ..., method = "rk5_i", atol = 1e-06, rtol = 1e-06, compile = TRUE ) ## S4 method for signature 'ODESystem' predict(object, init, duration, step_size = 1, start = 0)
dxdt(x) ## S4 method for signature 'DxdtOdeConstructor,ANY' e1 == e2 ODESystem( odesys, ..., method = "rk5_i", atol = 1e-06, rtol = 1e-06, compile = TRUE ) ## S4 method for signature 'ODESystem' predict(object, init, duration, step_size = 1, start = 0)
x |
A SymEngine Basic object of type Symbol or a R object
that will be converted to |
e1 |
A DxdtOdeConstructor S4 object which can be returned by 'dxdt'. |
e2 |
A Basic object or an R object that will be converted to 'S(e2)'. |
odesys , ...
|
DxdtOde S4 objects that can be returned with 'dxdt(x) == rhs'. Or 'odesys' can be a list of DxdtOde S4 objects when there is no dot arguments. |
method , atol , rtol
|
Passed to 'odeintr::compile_sys'. |
compile |
Logical, whether to compile the C++ source. Useful if you only want to obtain the code. |
object |
A ODESystem S4 object. |
init |
A numeric vector specifying the initial conditions. It can be named with the variable names or it can be unnamed but in the same of order of equations. |
duration , step_size , start
|
Passed to the function generated by 'odeintr::compile_sys'. |
dxdt
returns a DxdtOdeConstructor S4 object.
S4 method of '==' for "DxdtOdeConstructor" returns a DxdtOde S4 object.
'ODESystem' returns a ODESystem S4 object.
'predict' returns a dataframe.
# A differential equation specified with dxdt and == x <- Symbol("x") eq <- dxdt(x) == 1/exp(x) print(eq) ## Not run: ## Lorenz system use_vars(x, y, z) sigma <- 10 rho <- 28 beta <- 8/3 lorenz_sys <- ODESystem( dxdt(x) == sigma * (y - x), dxdt(y) == (rho - z) * x - y, dxdt(z) == - beta * z + x * y ) res <- predict( lorenz_sys, init = c(x = 1, y = 1, z = 1), duration = 100, step_size = 0.001 ) plot(res[, c(2, 4)], type = 'l', col = "steelblue", main = "Lorenz Attractor") ## End(Not run)
# A differential equation specified with dxdt and == x <- Symbol("x") eq <- dxdt(x) == 1/exp(x) print(eq) ## Not run: ## Lorenz system use_vars(x, y, z) sigma <- 10 rho <- 28 beta <- 8/3 lorenz_sys <- ODESystem( dxdt(x) == sigma * (y - x), dxdt(y) == (rho - z) * x - y, dxdt(z) == - beta * z + x * y ) res <- predict( lorenz_sys, init = c(x = 1, y = 1, z = 1), duration = 100, step_size = 0.001 ) plot(res[, c(2, 4)], type = 'l', col = "steelblue", main = "Lorenz Attractor") ## End(Not run)
This function will evaluate a SymEngine object to its "numerical" form
with given precision. User may further use as.double()
to convert
to R value.
evalf(expr, bits = 53L, complex = FALSE)
evalf(expr, bits = 53L, complex = FALSE)
expr |
A SymEngine object. |
bits |
The precision. |
complex |
Whether or not to be evaluated as a complex number. |
Same type as expr
argument.
expr <- Constant("pi") evalf(expr) as.double(evalf(expr)) == pi
expr <- Constant("pi") evalf(expr) as.double(evalf(expr)) == pi
This function takes a SymEngine object and return its expanded form.
expand(x)
expand(x)
x |
A Basic/VecBasic/DenseMatrix S4 object. |
Same type as input.
expr <- S(~ (x + y) ^ 3) expand(expr)
expr <- S(~ (x + y) ^ 3) expand(expr)
FunctionSymbol
creates a Basic object with type FunctionSymbol
.
Function
returns a generator.
Function(name) FunctionSymbol(name, args)
Function(name) FunctionSymbol(name, args)
name |
Name of the function symbol |
args |
Dependent symbols |
FunctionSymbol
returns a Basic
. Function
returns
a function that will return a Basic
f <- Function("f") a <- Symbol("a") b <- Symbol("b") f(a, b) e <- f(a, f(a + b)) D(e, a) FunctionSymbol("f", c(a,b))
f <- Function("f") a <- Symbol("a") b <- Symbol("b") f(a, b) e <- f(a, f(a + b)) D(e, a) FunctionSymbol("f", c(a,b))
These functions are used to access the underlying properties of a
Basic
object.
get_type(x) get_args(x) get_hash(x) get_str(x) free_symbols(x) function_symbols(x) get_name(x) get_prec(x)
get_type(x) get_args(x) get_hash(x) get_str(x) free_symbols(x) function_symbols(x) get_name(x) get_prec(x)
x |
A Basic object. |
Return the internal type
Return the internal arguments of a Basic object as a VecBasic
Return the hash as a string
Return the string representation of the Basic object
Return free symbols in an expression
Return function symbols in an expression
Return name of a Basic object of type FunctionSymbol
Return precision of a Basic object of type RealMPFR
get_type()
, get_hash()
, get_str()
, get_name()
return a string.
get_args()
, free_symbols()
, function_symbols()
return a VecBasic
S4 object.
get_prec()
returns an integer.
These functions currently use DoubleVisitor
to
convert a Basic/VecBasic object to a DoubleVisitor
which
essentially is a S4 class extending R function.
lambdify(x, args, backend = c("auto", "lambda", "llvm"), perform_cse = TRUE) ## S3 method for class 'BasicOrVecBasic' as.function(x, args, backend = "auto", perform_cse = TRUE, ...)
lambdify(x, args, backend = c("auto", "lambda", "llvm"), perform_cse = TRUE) ## S3 method for class 'BasicOrVecBasic' as.function(x, args, backend = "auto", perform_cse = TRUE, ...)
x |
A Basic object or a VecBasic object. |
args |
A VecBasic object specifying the arguments of the resulted function.
It will be passed to |
backend |
One of "auto", "lambda" and "llvm". If "auto",
|
perform_cse |
Passed to |
... |
Not used |
A DoubleVisitor
S4 object.
These are some special mathematical functions and functions related to number theory.
LCM(a, b) GCD(a, b) nextprime(a) factorial(x) ## S4 method for signature 'SymEngineDataType' factorial(x) choose(n, k) ## S4 method for signature 'SymEngineDataType' choose(n, k) zeta(a) lambertw(a) dirichlet_eta(a) erf(a) erfc(a) ## S4 method for signature 'SymEngineDataType,SymEngineDataType' atan2(y, x) kronecker_delta(x, y) lowergamma(x, a) uppergamma(x, a) ## S4 method for signature 'SymEngineDataType,SymEngineDataType' beta(a, b) ## S4 method for signature 'SymEngineDataType' psigamma(x, deriv = 0L) ## S4 method for signature 'SymEngineDataType' digamma(x) ## S4 method for signature 'SymEngineDataType' trigamma(x)
LCM(a, b) GCD(a, b) nextprime(a) factorial(x) ## S4 method for signature 'SymEngineDataType' factorial(x) choose(n, k) ## S4 method for signature 'SymEngineDataType' choose(n, k) zeta(a) lambertw(a) dirichlet_eta(a) erf(a) erfc(a) ## S4 method for signature 'SymEngineDataType,SymEngineDataType' atan2(y, x) kronecker_delta(x, y) lowergamma(x, a) uppergamma(x, a) ## S4 method for signature 'SymEngineDataType,SymEngineDataType' beta(a, b) ## S4 method for signature 'SymEngineDataType' psigamma(x, deriv = 0L) ## S4 method for signature 'SymEngineDataType' digamma(x) ## S4 method for signature 'SymEngineDataType' trigamma(x)
a , b , x , y , n , k , deriv
|
SymEngine objects ( |
Same type as input.
Miscellaneous S4 methods defined for VecBasic
class.
## S4 method for signature 'VecBasic' length(x) ## S3 method for class 'VecBasic' rep(x, ...) ## S3 method for class 'Basic' rep(x, ...) ## S3 method for class 'VecBasic' unique(x, ...) ## S4 method for signature 'BasicOrVecBasic' c(x, ...) ## S4 method for signature 'VecBasic,numeric' x[[i, j, ...]] ## S4 method for signature 'VecBasic' x[i, j, ..., drop = TRUE] ## S4 replacement method for signature 'VecBasic' x[[i]] <- value ## S4 replacement method for signature 'VecBasic' x[i, j, ...] <- value
## S4 method for signature 'VecBasic' length(x) ## S3 method for class 'VecBasic' rep(x, ...) ## S3 method for class 'Basic' rep(x, ...) ## S3 method for class 'VecBasic' unique(x, ...) ## S4 method for signature 'BasicOrVecBasic' c(x, ...) ## S4 method for signature 'VecBasic,numeric' x[[i, j, ...]] ## S4 method for signature 'VecBasic' x[i, j, ..., drop = TRUE] ## S4 replacement method for signature 'VecBasic' x[[i]] <- value ## S4 replacement method for signature 'VecBasic' x[i, j, ...] <- value
x |
Basic object or Vecbasic object. |
i , j , ... , drop , value
|
Arguments for subsetting or replacing. |
Same or similar to the generics.
This function constructs a symbolic matrix (DenseMatrix
S4 object)
with a similar interface with R's matrix
function.
Matrix(data, nrow = 1L, ncol = 1L, byrow = FALSE)
Matrix(data, nrow = 1L, ncol = 1L, byrow = FALSE)
data |
A R object. |
nrow , ncol
|
Number of rows and columns. |
byrow |
Boolean value. Whether the data should be filled by row or by column. |
DenseMatrix
S4 object.
'S' and 'Basic' converts a R object to a Basic object. 'Symbol', 'Real' and 'Constant' construct a Basic object with type "Symbol", "RealDouble"/"RealMPFR" and "Constant", respectively.
S(x) Basic(x) Symbol(x) Constant(x) Real(x, prec = NULL)
S(x) Basic(x) Symbol(x) Constant(x) Real(x, prec = NULL)
x |
A R object. |
prec |
If supplied, the argument will be parsed as a Basic object of type RealMPFR. |
For double vector, 'S' will check whether it is a whole number – if true, it will be converted to a Integer type. If this behavior is not desired, you can use 'Basic' or 'as(x, "Basic")'.
A Basic
S4 object.
S("(x + y)^2") S(~ (x + y)^2) S(NaN) S(42) Basic(42) as(42, "Basic") pi <- Constant("pi") evalf(pi) if (symengine_have_component("mpfr")) evalf(pi, 300) Real(42) if (symengine_have_component("mpfr")) Real(42, prec = 140)
S("(x + y)^2") S(~ (x + y)^2) S(NaN) S(42) Basic(42) as(42, "Basic") pi <- Constant("pi") evalf(pi) if (symengine_have_component("mpfr")) evalf(pi, 300) Real(42) if (symengine_have_component("mpfr")) Real(42, prec = 140)
Solve system of symbolic equations or solve a polynomial equation. Depending on types of arguments, it supports different modes. See Details and Examples.
solve(a, b, ...) ## S4 method for signature 'DenseMatrix' solve(a, b, ...) ## S4 method for signature 'VecBasic' solve(a, b, ...) ## S4 method for signature 'Basic' solve(a, b, ...)
solve(a, b, ...) ## S4 method for signature 'DenseMatrix' solve(a, b, ...) ## S4 method for signature 'VecBasic' solve(a, b, ...) ## S4 method for signature 'Basic' solve(a, b, ...)
a , b
|
Objects, see details. |
... |
Not used. |
solve
is a generic function dispatched on the class of the first argument.
If a
is a (square) DenseMatrix, it solves the equation
a %*% x = b
for x
. (similar to solve.default()
)
If a
is a DenseMatrix and b
is missing, b
is
taken to be an identity matrix and solve
will return the
inverse of a
. (similar to solve.default()
)
If a
is a VecBasic, it solves the system of linear equations
represented by a
with regards to symbols represented in b
.
If a
is a Basic, it solves the polynomial equation represented by
a with regards to the symbol represented in b
.
A VecBasic
or DenseMatrix
S4 object.
## Inverse of a symbolic matrix mat <- Matrix(c("A", "B", "C", "D"), 2) solve(mat) ## Solve a %*% x == b a <- Matrix(c("a11", "a21", "a12", "a22"), 2) # a is a 2x2 matrix b <- Vector("b1", "b2") # b is a length 2 vector solve(a, b) # Solution of x (2x1 matrix) ## Solve the system of linear equations represented by a with regards to ## symbols in b a <- Vector(~ -2*x + y - 4, # A system of linear equations ~ 3*x + y - 9) b <- Vector(~x, ~y) # Symbols to solve (x and y) solve(a, b) # Solution of x and y
## Inverse of a symbolic matrix mat <- Matrix(c("A", "B", "C", "D"), 2) solve(mat) ## Solve a %*% x == b a <- Matrix(c("a11", "a21", "a12", "a22"), 2) # a is a 2x2 matrix b <- Vector("b1", "b2") # b is a length 2 vector solve(a, b) # Solution of x (2x1 matrix) ## Solve the system of linear equations represented by a with regards to ## symbols in b a <- Vector(~ -2*x + y - 4, # A system of linear equations ~ 3*x + y - 9) b <- Vector(~x, ~y) # Symbols to solve (x and y) solve(a, b) # Solution of x and y
This function will substitute expr
with pairs of
values in the dot arguments. The length of dot arguments must
be a even number.
subs(expr, ...)
subs(expr, ...)
expr |
A |
... |
Pairs of Basic objects or values can be converted to |
Same type as expr
.
symengine
is a R package for symbolic computation.
SymEngine library is a standalone fast symbolic manipulation library written in C++. It allows computation over mathematical expressions in a way which is similar to the traditional manual computations of mathematicians and scientists. The R interface of the library tries to provide a user-friendly way to do symbolic computation in R and can be integrated into other packages to help solve related tasks. The design of the package is somehow similar to the SymPy package in Python. Unlike some other computer algebra systems, it does not invent its own language or domain specific language but uses R language to manipulate the symbolic expressions.
symengine
uses the S4 dispatch system extensively to differentiate between calculation
over normal R objects and symengine objects. For example, the semantics of sin
in
expr <- Symbol("x"); sin(expr)
is different from the sin
used over normal R numbers.
Basic
is simply a S4 class holding a pointer representing a symbolic expression
in symengine. Basic
objects have the same S4 class but can have different
C-level representations which can be accessed via get_type()
.
For example, Basic(~ 1/2)
will have "Rational" type and Basic(1/2)
will have
"RealDouble" type.
A Basic
object will also have a list of associated sub-components
which can be accessed via get_args()
. For example, (expr <- S("x") * 3L * S("a"))
will have type "Mul", and as.list(get_args(expr))
will show the three factors of
the multiplication.
A Basic
object can be constructed via Basic()
, S()
, Symbol()
, Constant()
or
Real()
.
VecBasic and DenseMatrix are S4 classes representing a symbolic vector or matrix.
They can be constructed with Vector()
, V()
, Matrix()
, c()
, rbind()
or cbind()
. For example the following code will construct a 2x3 matrix.
vec <- Vector("a", "b") cbind(vec, vec^2L, c(S("c"), S("d")))
The following functions are expected to work naturally with VecBasic and DenseMatrix classes.
[
, [[
, [<-
and [[<-
for subsetting and assignment.
dim()
, dim<-
, length()
, t()
, det()
, rbind()
, cbind()
, c()
, rep()
%*%
for matrix multiplication
solve(a, b)
: solve a %*% x = b
where a
is a square DenseMatrix and
b
is a VecBasic/DenseMatrix.
solve(a)
: find the inverse of a
where a
is a square DenseMatrix.
solve(a, b)
: solve system of linear equations represented by a
(VecBasic) with
regards to symbols in b
(VecBasic).
Further, the R functions that work on Basic objects (e.g. sin
) are expected work
on VecBasic and DenseMatrix objects as well in a vectorized manner.
The following is a (incomplete) list of functions that are expected to work with
symengine objects. Note that these functions can also be used inside a formula or
R language objects and passed to S or Basic or Vector to construct symengine
objects. For example S(~ sin(x) + 1)
and S(quote(sin(x) + 1))
.
+
, -
, *
, /
, ^
abs
, sqrt
, exp
, expm1
, log
, log10
, log2
, log1p
cos
, cosh
, sin
, sinh
, tan
, tanh
, acos
, acosh
, asin
, asinh
, atan
, atanh
cospi
, sinpi
, tanpi
, gamma
, lgamma
, digamma
, trigamma
lambertw
, zeta
, dirichlet_eta
, erf
, erfc
atan2
, kronecker_delta
, lowergamma
, uppergamma
, psigamma
, beta
Functions to get symengine logo, version and external libraries built with.
symengine_version() symengine_ascii_art() symengine_have_component( which = c("mpfr", "flint", "arb", "mpc", "ecm", "primesieve", "piranha", "boost", "pthread", "llvm") ) symengine_compilation_notes()
symengine_version() symengine_ascii_art() symengine_have_component( which = c("mpfr", "flint", "arb", "mpc", "ecm", "primesieve", "piranha", "boost", "pthread", "llvm") ) symengine_compilation_notes()
which |
A character vector. |
Character vector.
S4 methods of t
defined for Basic
, VecBasic
and DenseMatrix
.
t(x) ## S4 method for signature 'Basic' t(x) ## S4 method for signature 'VecBasic' t(x) ## S4 method for signature 'DenseMatrix' t(x)
t(x) ## S4 method for signature 'Basic' t(x) ## S4 method for signature 'VecBasic' t(x) ## S4 method for signature 'DenseMatrix' t(x)
x |
A SymEngine object. |
A DenseMatrix
S4 object.
This is a convenient way to initialize variables and assign them in the given environment.
use_vars(..., .env = parent.frame(), .quiet = FALSE)
use_vars(..., .env = parent.frame(), .quiet = FALSE)
... |
All the arguments will be quoted and parsed, if a argument is named, the name will be used as the name of variable to assign, otherwise the argument can only be a symbol. |
.env |
Environment to assign. |
.quiet |
Whether to supress the message. |
Invisibly returns a list of assigned variables.
use_vars(x, y, expr = "a + b", p = 3.14) p * x + y expand(expr^2L) rm(x, y, expr, p)
use_vars(x, y, expr = "a + b", p = 3.14) p * x + y expand(expr^2L) rm(x, y, expr, p)
A symbolic vector is represented by VecBasic
S4 class.
Vector
and V
are constructors of VecBasic
.
Vector(x, ...) V(...)
Vector(x, ...) V(...)
x , ...
|
R objects. |
There are some differences between Vector
and V
.
For double values, V
will check whether they are
whole number, and convert them to integer if so.
Vector
will not.
V
does not accept "non-scalar" arguments,
like Vector(c(1,2,3))
.
A VecBasic
.
a <- S("a") b <- S("b") Vector(a, b, a + b, 42L) Vector(list(a, b, 42L)) Vector(1,2,a) V(1,2,a)
a <- S("a") b <- S("b") Vector(a, b, a + b, 42L) Vector(list(a, b, 42L)) Vector(1,2,a) V(1,2,a)