Reference documentation for deal.II version 9.1.1
|
Classes | |
class | Expression |
struct | is_sd_number |
struct | is_sd_number< Expression > |
struct | is_sd_number< SymEngine::Expression > |
struct | is_symengine_number |
struct | is_symengine_number< Expression > |
struct | is_symengine_number< SymEngine::Expression > |
Functions | |
static ::ExceptionBase & | ExcSymEngineParserError (std::string arg1) |
Expression | operator & (const Expression &lhs, const Expression &rhs) |
Expression | operator && (const Expression &lhs, const Expression &rhs) |
Power functions | |
Expression | pow (const Expression &base, const Expression &exponent) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | pow (const Expression &base, const NumberType &exponent) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | pow (const NumberType &base, const Expression &exponent) |
Expression | sqrt (const Expression &x) |
Expression | cbrt (const Expression &x) |
Expression | exp (const Expression &exponent) |
Expression | log (const Expression &x) |
Expression | log (const Expression &x, const Expression &base) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | log (const Expression &x, const NumberType &base) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | log (const NumberType &x, const Expression &base) |
Expression | log10 (const Expression &x) |
Trigonometric functions | |
Expression | sin (const Expression &x) |
Expression | cos (const Expression &x) |
Expression | tan (const Expression &x) |
Expression | csc (const Expression &x) |
Expression | sec (const Expression &x) |
Expression | cot (const Expression &x) |
Expression | asin (const Expression &x) |
Expression | acos (const Expression &x) |
Expression | atan (const Expression &x) |
Expression | atan2 (const Expression &y, const Expression &x) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | atan2 (const NumberType &y, const Expression &x) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | atan2 (const Expression &y, const NumberType &x) |
Expression | acsc (const Expression &x) |
Expression | asec (const Expression &x) |
Expression | acot (const Expression &x) |
Hyperbolic trigonometric functions | |
Expression | sinh (const Expression &x) |
Expression | cosh (const Expression &x) |
Expression | tanh (const Expression &x) |
Expression | csch (const Expression &x) |
Expression | sech (const Expression &x) |
Expression | coth (const Expression &x) |
Expression | asinh (const Expression &x) |
Expression | acosh (const Expression &x) |
Expression | atanh (const Expression &x) |
Expression | acsch (const Expression &x) |
Expression | asech (const Expression &x) |
Expression | acoth (const Expression &x) |
Other functions | |
Expression | abs (const Expression &x) |
Expression | fabs (const Expression &x) |
Expression | sign (const Expression &x) |
Expression | copysign (const Expression &value, const Expression &sign) |
Expression | floor (const Expression &x) |
Expression | ceil (const Expression &x) |
Expression | max (const Expression &a, const Expression &b) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | max (const Expression &a, const NumberType &b) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | max (const NumberType &a, const Expression &b) |
Expression | min (const Expression &a, const Expression &b) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | min (const Expression &a, const NumberType &b) |
template<typename NumberType , typename = typename std::enable_if< !std::is_same<NumberType, Expression>::value>::type> | |
Expression | min (const NumberType &a, const Expression &b) |
Expression | erf (const Expression &x) |
Expression | erfc (const Expression &x) |
Bitwise operators | |
std::ostream & | operator<< (std::ostream &stream, const Expression &expression) |
std::istream & | operator>> (std::istream &stream, Expression &expression) |
Comparison operators | |
Expression | operator== (const Expression &lhs, const Expression &rhs) |
Expression | operator!= (const Expression &lhs, const Expression &rhs) |
Expression | operator< (const Expression &lhs, const Expression &rhs) |
Expression | operator> (const Expression &lhs, const Expression &rhs) |
Expression | operator<= (const Expression &lhs, const Expression &rhs) |
Expression | operator>= (const Expression &lhs, const Expression &rhs) |
Logical operators | |
Expression | operator! (const Expression &expression) |
Expression | operator & (const Expression &lhs, const Expression &rhs) |
Expression | operator| (const Expression &lhs, const Expression &rhs) |
Expression | operator^ (const Expression &lhs, const Expression &rhs) |
Expression | operator && (const Expression &lhs, const Expression &rhs) |
Expression | operator|| (const Expression &lhs, const Expression &rhs) |
Mathematical operators | |
Expression | operator+ (Expression lhs, const Expression &rhs) |
Expression | operator- (Expression lhs, const Expression &rhs) |
Expression | operator* (Expression lhs, const Expression &rhs) |
Expression | operator/ (Expression lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator+ (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator+ (const Expression &lhs, const NumberType &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator- (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator- (const Expression &lhs, const NumberType &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator* (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator* (const Expression &lhs, const NumberType &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator/ (const NumberType &lhs, const Expression &rhs) |
template<typename NumberType , typename = typename std::enable_if< std::is_constructible<Expression, NumberType>::value>::type> | |
Expression | operator/ (const Expression &lhs, const NumberType &rhs) |
Symbolic variable creation | |
Expression | make_symbol (const std::string &symbol) |
Expression | make_symbolic_function (const std::string &symbol, const types::symbol_vector &arguments) |
Expression | make_symbolic_function (const std::string &symbol, const types::substitution_map &arguments) |
template<int dim> | |
Tensor< 1, dim, Expression > | make_vector_of_symbols (const std::string &symbol) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | make_tensor_of_symbols (const std::string &symbol) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | make_symmetric_tensor_of_symbols (const std::string &symbol) |
template<int dim> | |
Tensor< 1, dim, Expression > | make_vector_of_symbolic_functions (const std::string &symbol, const types::substitution_map &arguments) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | make_tensor_of_symbolic_functions (const std::string &symbol, const types::substitution_map &arguments) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | make_symmetric_tensor_of_symbolic_functions (const std::string &symbol, const types::substitution_map &arguments) |
Symbolic differentiation | |
Expression | differentiate (const Expression &f, const Expression &x) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | differentiate (const Expression &f, const Tensor< rank, dim, Expression > &T) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | differentiate (const Expression &f, const SymmetricTensor< rank, dim, Expression > &S) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | differentiate (const Tensor< 0, dim, Expression > &f, const Tensor< rank, dim, Expression > &T) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | differentiate (const Tensor< 0, dim, Expression > &f, const SymmetricTensor< rank, dim, Expression > &S) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | differentiate (const Tensor< rank, dim, Expression > &T, const Expression &x) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | differentiate (const SymmetricTensor< rank, dim, Expression > &S, const Expression &x) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | differentiate (const Tensor< rank, dim, Expression > &T, const Tensor< 0, dim, Expression > &x) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | differentiate (const SymmetricTensor< rank, dim, Expression > &S, const Tensor< 0, dim, Expression > &x) |
template<int rank_1, int rank_2, int dim> | |
Tensor< rank_1+rank_2, dim, Expression > | differentiate (const Tensor< rank_1, dim, Expression > &T1, const Tensor< rank_2, dim, Expression > &T2) |
template<int rank_1, int rank_2, int dim> | |
SymmetricTensor< rank_1+rank_2, dim, Expression > | differentiate (const SymmetricTensor< rank_1, dim, Expression > &S1, const SymmetricTensor< rank_2, dim, Expression > &S2) |
template<int rank_1, int rank_2, int dim> | |
Tensor< rank_1+rank_2, dim, Expression > | differentiate (const Tensor< rank_1, dim, Expression > &T, const SymmetricTensor< rank_2, dim, Expression > &S) |
template<int rank_1, int rank_2, int dim> | |
Tensor< rank_1+rank_2, dim, Expression > | differentiate (const SymmetricTensor< rank_1, dim, Expression > &S, const Tensor< rank_2, dim, Expression > &T) |
Symbol map creation and manipulation | |
template<bool ignore_invalid_symbols = false, typename ValueType = double, typename SymbolicType > | |
types::substitution_map | make_symbol_map (const SymbolicType &symbol) |
template<bool ignore_invalid_symbols = false, typename ValueType = double, typename SymbolicType , typename... Args> | |
types::substitution_map | make_symbol_map (const SymbolicType &symbol, const Args &... other_symbols) |
template<bool ignore_invalid_symbols = false, typename ValueType = double> | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const Expression &symbol) |
template<bool ignore_invalid_symbols = false, typename ValueType = double, typename SymbolicType , typename T = typename std::enable_if< !std::is_base_of<Expression, SymbolicType>::value && ::internal::is_explicitly_convertible< SymbolicType, const SymEngine::RCP<const SymEngine::Basic> &>::value>::type> | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const SymbolicType &symbol) |
template<bool ignore_invalid_symbols = false, typename ValueType = double, typename SymbolicType > | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const std::vector< SymbolicType > &symbols) |
template<bool ignore_invalid_symbols = false, typename ValueType = double> | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const types::substitution_map &other_symbols) |
template<bool ignore_invalid_symbols = false, typename ValueType = double, typename SymbolicType , typename... Args> | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const SymbolicType &symbol, const Args &... other_symbols) |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const Expression &symbol, const Expression &value) |
template<typename SymbolicType , typename ValueType , typename T = typename std::enable_if< ::internal::is_explicitly_convertible< SymbolicType, const SymEngine::RCP<const SymEngine::Basic> &>::value && std::is_constructible<SymbolicType, ValueType>::value>::type> | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const SymbolicType &symbol, const ValueType &value) |
template<typename SymbolicType , typename ValueType > | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const std::vector< SymbolicType > &symbols, const std::vector< ValueType > &values) |
template<typename SymbolicType , typename ValueType > | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const std::pair< SymbolicType, ValueType > &symbol_value) |
template<typename SymbolicType , typename ValueType , typename... Args> | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const std::pair< SymbolicType, ValueType > &symbol_value, const Args &... other_symbol_values) |
template<typename SymbolicType , typename ValueType > | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const std::vector< std::pair< SymbolicType, ValueType >> &symbol_values) |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const types::substitution_map &symbol_values) |
template<bool ignore_invalid_symbols = false, typename ValueType = double, int rank, int dim, typename SymbolicType > | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const Tensor< rank, dim, SymbolicType > &symbol_tensor) |
template<bool ignore_invalid_symbols = false, typename ValueType = double, int rank, int dim, typename SymbolicType > | |
void | add_to_symbol_map (types::substitution_map &symbol_map, const SymmetricTensor< rank, dim, SymbolicType > &symbol_tensor) |
template<int rank, int dim, typename SymbolicType , typename ValueType > | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const Tensor< rank, dim, SymbolicType > &symbol_tensor, const Tensor< rank, dim, ValueType > &value_tensor) |
template<int rank, int dim, typename SymbolicType , typename ValueType > | |
void | set_value_in_symbol_map (types::substitution_map &substitution_map, const SymmetricTensor< rank, dim, SymbolicType > &symbol_tensor, const SymmetricTensor< rank, dim, ValueType > &value_tensor) |
Symbol substitution map creation | |
types::substitution_map | make_substitution_map (const Expression &symbol, const Expression &value) |
template<typename ExpressionType , typename ValueType , typename T = typename std::enable_if< ::internal::is_explicitly_convertible< ExpressionType, const SymEngine::RCP<const SymEngine::Basic> &>::value && std::is_constructible<ExpressionType, ValueType>::value>::type> | |
types::substitution_map | make_substitution_map (const ExpressionType &symbol, const ValueType &value) |
template<typename ExpressionType , typename ValueType > | |
types::substitution_map | make_substitution_map (const std::vector< ExpressionType > &symbols, const std::vector< ValueType > &values) |
template<typename ExpressionType , typename ValueType > | |
types::substitution_map | make_substitution_map (const std::pair< ExpressionType, ValueType > &symbol_value) |
template<typename ExpressionType , typename ValueType > | |
types::substitution_map | make_substitution_map (const std::vector< std::pair< ExpressionType, ValueType >> &symbol_values) |
template<typename ExpressionType , typename ValueType , typename... Args> | |
types::substitution_map | make_substitution_map (const std::pair< ExpressionType, ValueType > &symbol_value, const Args &... other_symbol_values) |
template<int rank, int dim, typename ExpressionType , typename ValueType > | |
types::substitution_map | make_substitution_map (const Tensor< rank, dim, ExpressionType > &symbol_tensor, const Tensor< rank, dim, ValueType > &value_tensor) |
template<int rank, int dim, typename ExpressionType , typename ValueType > | |
types::substitution_map | make_substitution_map (const SymmetricTensor< rank, dim, ExpressionType > &symbol_tensor, const SymmetricTensor< rank, dim, ValueType > &value_tensor) |
Symbol substitution map enlargement | |
template<bool ignore_invalid_symbols = false> | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const Expression &symbol, const Expression &value) |
template<bool ignore_invalid_symbols = false, typename ExpressionType , typename ValueType , typename = typename std::enable_if< ::internal::is_explicitly_convertible< ExpressionType, const SymEngine::RCP<const SymEngine::Basic> &>::value && std::is_constructible<ExpressionType, ValueType>::value>::type> | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const ExpressionType &symbol, const ValueType &value) |
template<bool ignore_invalid_symbols = false, typename ExpressionType , typename ValueType , typename = typename std::enable_if< ::internal::is_explicitly_convertible< ExpressionType, const SymEngine::RCP<const SymEngine::Basic> &>::value && std::is_constructible<ExpressionType, ValueType>::value>::type> | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const std::vector< ExpressionType > &symbols, const std::vector< ValueType > &values) |
template<bool ignore_invalid_symbols = false> | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const types::substitution_map &symbol_values) |
template<bool ignore_invalid_symbols = false, typename ExpressionType , typename ValueType > | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const std::pair< ExpressionType, ValueType > &symbol_value) |
template<bool ignore_invalid_symbols = false, typename ExpressionType , typename ValueType > | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const std::vector< std::pair< ExpressionType, ValueType >> &symbol_values) |
template<bool ignore_invalid_symbols = false, typename ExpressionType , typename ValueType , typename... Args> | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const std::pair< ExpressionType, ValueType > &symbol_value, const Args &... other_symbol_values) |
void | merge_substitution_maps (types::substitution_map &substitution_map_out, const types::substitution_map &substitution_map_in) |
template<typename... Args> | |
void | merge_substitution_maps (types::substitution_map &substitution_map_out, const types::substitution_map &substitution_map_in, const Args &... other_substitution_maps_in) |
template<typename... Args> | |
types::substitution_map | merge_substitution_maps (const types::substitution_map &substitution_map_in, const Args &... other_substitution_maps_in) |
template<bool ignore_invalid_symbols = false, int rank, int dim, typename ExpressionType , typename ValueType > | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const Tensor< rank, dim, ExpressionType > &symbol_tensor, const Tensor< rank, dim, ValueType > &value_tensor) |
template<bool ignore_invalid_symbols = false, int rank, int dim, typename ExpressionType , typename ValueType > | |
void | add_to_substitution_map (types::substitution_map &substitution_map, const SymmetricTensor< rank, dim, ExpressionType > &symbol_tensor, const SymmetricTensor< rank, dim, ValueType > &value_tensor) |
Symbol substitution and evaluation | |
types::substitution_map | resolve_explicit_dependencies (const types::substitution_map &substitution_map, const bool force_cyclic_dependency_resolution=false) |
template<typename ExpressionType , typename ValueType > | |
types::substitution_map | resolve_explicit_dependencies (const std::vector< std::pair< ExpressionType, ValueType >> &symbol_values, const bool force_cyclic_dependency_resolution=false) |
Expression | substitute (const Expression &expression, const types::substitution_map &substitution_map) |
template<typename ValueType > | |
Expression | substitute (const Expression &expression, const Expression &symbol, const ValueType &value) |
template<typename ExpressionType , typename... Args> | |
ExpressionType | substitute (const ExpressionType &expression, const Args &... symbol_values) |
template<typename ValueType > | |
ValueType | substitute_and_evaluate (const Expression &expression, const types::substitution_map &substitution_map) |
template<typename ValueType , typename... Args> | |
ValueType | substitute_and_evaluate (const Expression &expression, const Args &... symbol_values) |
template<int rank, int dim> | |
Tensor< rank, dim, Expression > | substitute (const Tensor< rank, dim, Expression > &expression_tensor, const types::substitution_map &substitution_map) |
template<int rank, int dim> | |
SymmetricTensor< rank, dim, Expression > | substitute (const SymmetricTensor< rank, dim, Expression > &expression_tensor, const types::substitution_map &substitution_map) |
template<typename ValueType , int rank, int dim> | |
Tensor< rank, dim, ValueType > | substitute_and_evaluate (const Tensor< rank, dim, Expression > &expression_tensor, const types::substitution_map &substitution_map) |
template<typename ValueType , int rank, int dim> | |
SymmetricTensor< rank, dim, ValueType > | substitute_and_evaluate (const SymmetricTensor< rank, dim, Expression > &expression_tensor, const types::substitution_map &substitution_map) |
Wrappers for symbolic differentiation libraries. Currently there is support for the following libraries:
Expression Differentiation::SD::pow | ( | const Expression & | base, |
const Expression & | exponent | ||
) |
Return a symbolic number that represents a base
value raised to the power of an exponent
.
Mimics the function std::pow(base,exponent)
using the standard math library.
Definition at line 42 of file symengine_math.cc.
Expression Differentiation::SD::pow | ( | const Expression & | base, |
const NumberType & | exponent | ||
) |
Return a symbolic number that represents a base
value raised to the power of an exponent
.
Mimics the function std::pow(base,exponent)
using the standard math library.
This variant is used when the exponent
is not a Expression.
Definition at line 75 of file symengine_math.h.
Expression Differentiation::SD::pow | ( | const NumberType & | base, |
const Expression & | exponent | ||
) |
Return a symbolic number that represents a base
value raised to the power of an exponent
.
Mimics the function std::pow(base,exponent)
using the standard math library.
This variant is used when the base
is not a Expression.
Definition at line 94 of file symengine_math.h.
Expression Differentiation::SD::sqrt | ( | const Expression & | x | ) |
Return a symbolic number that represents the square root of some value x
.
Mimics the function std::sqrt(x)
using the standard math library.
Definition at line 49 of file symengine_math.cc.
Expression Differentiation::SD::cbrt | ( | const Expression & | x | ) |
Return a symbolic number that represents the cubic root of some value x
.
Mimics the function std::cbrt(x)
using the standard math library.
Definition at line 56 of file symengine_math.cc.
Expression Differentiation::SD::exp | ( | const Expression & | exponent | ) |
Return a symbolic number that represents the Euler constant \(e \approx 2.71828\) raised to the given exponent
.
Mimics the function std::exp(exponent)
using the standard math library.
Definition at line 63 of file symengine_math.cc.
Expression Differentiation::SD::log | ( | const Expression & | x | ) |
Return a symbolic number that represents the natural logarithm of a value x
.
Mimics the function std::log(x)
using the standard math library.
Definition at line 70 of file symengine_math.cc.
Expression Differentiation::SD::log | ( | const Expression & | x, |
const Expression & | base | ||
) |
Return a symbolic number that represents the logarithm of a value x
taken with respect to a base
number.
Mimics the function std::log(x,base)
using the standard math library.
Definition at line 77 of file symengine_math.cc.
Expression Differentiation::SD::log | ( | const Expression & | x, |
const NumberType & | base | ||
) |
Return a symbolic number that represents the logarithm of a value x
taken with respect to a base
number.
Mimics the function std::log(x,base)
using the standard math library.
This variant is used when the base
is not a Expression.
Definition at line 160 of file symengine_math.h.
Expression Differentiation::SD::log | ( | const NumberType & | x, |
const Expression & | base | ||
) |
Return a symbolic number that represents the logarithm of a value x
taken with respect to a base
number.
Mimics the function std::log(x,base)
using the standard math library.
This variant is used when the value
is not a Expression.
Definition at line 179 of file symengine_math.h.
Expression Differentiation::SD::log10 | ( | const Expression & | x | ) |
Return a symbolic number that represents the base 10 logarithm of a value x
.
Mimics the function std::log10(x)
using the standard math library.
Definition at line 84 of file symengine_math.cc.
Expression Differentiation::SD::sin | ( | const Expression & | x | ) |
Return a symbolic number that represents the sine function with the given argument x
.
Mimics the function std::sin(x)
using the standard math library.
Definition at line 91 of file symengine_math.cc.
Expression Differentiation::SD::cos | ( | const Expression & | x | ) |
Return a symbolic number that represents the cosine function with the given argument x
.
Mimics the function std::cos(x)
using the standard math library.
Definition at line 98 of file symengine_math.cc.
Expression Differentiation::SD::tan | ( | const Expression & | x | ) |
Return a symbolic number that represents the tangent function with the given argument x
.
Mimics the function std::tan(x)
using the standard math library.
Definition at line 105 of file symengine_math.cc.
Expression Differentiation::SD::csc | ( | const Expression & | x | ) |
Return a symbolic number that represents the cosecant function with the given argument x
.
Mimics the function 1.0/std::sin(x)
using the standard math library.
Definition at line 112 of file symengine_math.cc.
Expression Differentiation::SD::sec | ( | const Expression & | x | ) |
Return a symbolic number that represents the secant function with the given argument x
.
Mimics the function 1.0/std::cos(x)
using the standard math library.
Definition at line 119 of file symengine_math.cc.
Expression Differentiation::SD::cot | ( | const Expression & | x | ) |
Return a symbolic number that represents the cotangent function with the given argument x
.
Mimics the function 1.0/std::tan(x)
using the standard math library.
Definition at line 126 of file symengine_math.cc.
Expression Differentiation::SD::asin | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse sine function with the given argument x
.
Mimics the function std::asin(x)
using the standard math library.
Definition at line 133 of file symengine_math.cc.
Expression Differentiation::SD::acos | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse cosine function with the given argument x
.
Mimics the function std::acos(x)
using the standard math library.
Definition at line 140 of file symengine_math.cc.
Expression Differentiation::SD::atan | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse tangent function with the given argument x
.
Mimics the function std::atan(x)
using the standard math library.
Definition at line 147 of file symengine_math.cc.
Expression Differentiation::SD::atan2 | ( | const Expression & | y, |
const Expression & | x | ||
) |
Return a symbolic number that represents the inverse tangent function with the given arguments x
and y
.
Mimics the function std::atan2(y,x)
using the standard math library.
Definition at line 154 of file symengine_math.cc.
Expression Differentiation::SD::atan2 | ( | const NumberType & | y, |
const Expression & | x | ||
) |
Return a symbolic number that represents the inverse tangent function with the given arguments x
and y
.
Mimics the function std::atan2(y,x)
using the standard math library.
This variant is used when the numerator y
is not a Expression.
Definition at line 319 of file symengine_math.h.
Expression Differentiation::SD::atan2 | ( | const Expression & | y, |
const NumberType & | x | ||
) |
Return a symbolic number that represents the inverse tangent function with the given arguments x
and y
.
Mimics the function std::atan2(y,x)
using the standard math library.
This variant is used when the denominator x
is not a Expression.
Definition at line 339 of file symengine_math.h.
Expression Differentiation::SD::acsc | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse cosecant function with the given argument x
.
Mimics the function 1.0/std::asin(x)
using the standard math library.
Definition at line 161 of file symengine_math.cc.
Expression Differentiation::SD::asec | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse secant function with the given argument x
.
Mimics the function 1.0/std::acos(x)
using the standard math library.
Definition at line 168 of file symengine_math.cc.
Expression Differentiation::SD::acot | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse cotangent function with the given argument x
.
Mimics the function 1.0/std::atan(x)
using the standard math library.
Definition at line 175 of file symengine_math.cc.
Expression Differentiation::SD::sinh | ( | const Expression & | x | ) |
Return a symbolic number that represents the hyperbolic sine function with the given argument x
.
Mimics the function std::sinh(x)
using the standard math library.
Definition at line 182 of file symengine_math.cc.
Expression Differentiation::SD::cosh | ( | const Expression & | x | ) |
Return a symbolic number that represents the hyperbolic cosine function with the given argument x
.
Mimics the function std::cosh(x)
using the standard math library.
Definition at line 189 of file symengine_math.cc.
Expression Differentiation::SD::tanh | ( | const Expression & | x | ) |
Return a symbolic number that represents the hyperbolic tangent function with the given argument x
.
Mimics the function std::tanh(x)
using the standard math library.
Definition at line 196 of file symengine_math.cc.
Expression Differentiation::SD::csch | ( | const Expression & | x | ) |
Return a symbolic number that represents the hyperbolic cosecant function with the given argument x
.
Mimics the function 1.0/std::sinh(x)
using the standard math library.
Definition at line 203 of file symengine_math.cc.
Expression Differentiation::SD::sech | ( | const Expression & | x | ) |
Return a symbolic number that represents the hyperbolic secant function with the given argument x
.
Mimics the function 1.0/std::cosh(x)
using the standard math library.
Definition at line 210 of file symengine_math.cc.
Expression Differentiation::SD::coth | ( | const Expression & | x | ) |
Return a symbolic number that represents the hyperbolic cotangent function with the given argument x
.
Mimics the function 1.0/std::tanh(x)
using the standard math library.
Definition at line 217 of file symengine_math.cc.
Expression Differentiation::SD::asinh | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse hyperbolic sine function with the given argument x
.
Mimics the function std::asinh(x)
using the standard math library.
Definition at line 224 of file symengine_math.cc.
Expression Differentiation::SD::acosh | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse hyperbolic cosine function with the given argument x
.
Mimics the function std::acosh(x)
using the standard math library.
Definition at line 231 of file symengine_math.cc.
Expression Differentiation::SD::atanh | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse hyperbolic tangent function with the given argument x
.
Mimics the function std::atanh(x)
using the standard math library.
Definition at line 238 of file symengine_math.cc.
Expression Differentiation::SD::acsch | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse hyperbolic cosecant function with the given argument x
.
Mimics the function 1.0/std::asin(x)
using the standard math library.
Definition at line 245 of file symengine_math.cc.
Expression Differentiation::SD::asech | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse hyperbolic secant function with the given argument x
.
Mimics the function 1.0/std::acos(x)
using the standard math library.
Definition at line 252 of file symengine_math.cc.
Expression Differentiation::SD::acoth | ( | const Expression & | x | ) |
Return a symbolic number that represents the inverse hyperbolic cotangent function with the given argument x
.
Mimics the function 1.0/std::atan(x)
using the standard math library.
Definition at line 259 of file symengine_math.cc.
Expression Differentiation::SD::abs | ( | const Expression & | x | ) |
Return a symbolic number that represents the absolute value of value x
.
Mimics the function std::abs(x)
using the standard math library.
Definition at line 266 of file symengine_math.cc.
Expression Differentiation::SD::fabs | ( | const Expression & | x | ) |
Return a symbolic number that represents the absolute value of value x
.
Mimics the function std::fabs(x)
using the standard math library.
Definition at line 273 of file symengine_math.cc.
Expression Differentiation::SD::sign | ( | const Expression & | x | ) |
Return a symbolic number that represents the sign of value x
.
Although there is no such function in the standard library, it mimics the function boost::sign(x)
using the boost math library.
Definition at line 280 of file symengine_math.cc.
Expression Differentiation::SD::copysign | ( | const Expression & | value, |
const Expression & | sign | ||
) |
Return a symbolic number that represents the value
of the first argument that takes the sign
of the second argument.
Mimics the function std::copysign(value, sign)
using the standard math library.
Definition at line 287 of file symengine_math.cc.
Expression Differentiation::SD::floor | ( | const Expression & | x | ) |
Return a symbolic number that represents the floor of value x
.
Mimics the function std::floor(x)
using the standard math library.
Definition at line 294 of file symengine_math.cc.
Expression Differentiation::SD::ceil | ( | const Expression & | x | ) |
Return a symbolic number that represents the ceiling of value x
.
Mimics the function std::ceil(x)
using the standard math library.
Definition at line 301 of file symengine_math.cc.
Expression Differentiation::SD::max | ( | const Expression & | a, |
const Expression & | b | ||
) |
Return a symbolic number that represents the maximum of two values a
and b
.
Mimics the function std::max(a,b)
using the standard math library.
Definition at line 308 of file symengine_math.cc.
Expression Differentiation::SD::max | ( | const Expression & | a, |
const NumberType & | b | ||
) |
Return a symbolic number that represents the maximum of two values a
and b
.
Mimics the function std::max(a,b)
using the standard math library.
This variant is used when b
is not a Expression.
Definition at line 608 of file symengine_math.h.
Expression Differentiation::SD::max | ( | const NumberType & | a, |
const Expression & | b | ||
) |
Return a symbolic number that represents the maximum of two values a
and b
.
Mimics the function std::max(a,b)
using the standard math library.
This variant is used when a
is not a Expression.
Definition at line 627 of file symengine_math.h.
Expression Differentiation::SD::min | ( | const Expression & | a, |
const Expression & | b | ||
) |
Return a symbolic number that represents the minimum of two values a
and b
.
Mimics the function std::min(a,b)
using the standard math library.
Definition at line 315 of file symengine_math.cc.
Expression Differentiation::SD::min | ( | const Expression & | a, |
const NumberType & | b | ||
) |
Return a symbolic number that represents the minimum of two values a
and b
.
Mimics the function std::min(a,b)
using the standard math library.
This variant is used when b
is not a Expression.
Definition at line 656 of file symengine_math.h.
Expression Differentiation::SD::min | ( | const NumberType & | a, |
const Expression & | b | ||
) |
Return a symbolic number that represents the minimum of two values a
and b
.
Mimics the function std::min(a,b)
using the standard math library.
This variant is used when a
is not a Expression.
Definition at line 675 of file symengine_math.h.
Expression Differentiation::SD::erf | ( | const Expression & | x | ) |
Return a symbolic number that represents error function with the given argument x
.
Mimics the function std::erf(x)
using the standard math library.
Definition at line 322 of file symengine_math.cc.
Expression Differentiation::SD::erfc | ( | const Expression & | x | ) |
Return a symbolic number that represents complimentary error function with the given argument x
.
Mimics the function std::erfc(x)
using the standard math library.
Definition at line 329 of file symengine_math.cc.
std::ostream & Differentiation::SD::operator<< | ( | std::ostream & | stream, |
const Expression & | expression | ||
) |
Bitwise left shift operator.
This is used to output the expression
to the input stream
.
Definition at line 369 of file symengine_number_types.cc.
std::istream & Differentiation::SD::operator>> | ( | std::istream & | stream, |
Expression & | expression | ||
) |
Bitwise right shift operator.
This is used to read in an expression
from the input stream
.
Definition at line 377 of file symengine_number_types.cc.
Expression Differentiation::SD::operator== | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Equality operator.
Return whether the lhs
is equal to the rhs
.
Definition at line 387 of file symengine_number_types.cc.
Expression Differentiation::SD::operator!= | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Non-equality operator.
Return whether the lhs
is not equal to the rhs
.
Definition at line 394 of file symengine_number_types.cc.
Expression Differentiation::SD::operator< | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Less than operator.
Return whether the lhs
is less than the rhs
.
Definition at line 401 of file symengine_number_types.cc.
Expression Differentiation::SD::operator> | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Greater than operator.
Return whether the lhs
is greater than the rhs
.
Definition at line 408 of file symengine_number_types.cc.
Expression Differentiation::SD::operator<= | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Less than or equals operator.
Return whether the lhs
is less than, or equal to, the rhs
.
Definition at line 415 of file symengine_number_types.cc.
Expression Differentiation::SD::operator>= | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Greater than or equals operator.
Return whether the lhs
is greater than, or equal to, the rhs
.
Definition at line 422 of file symengine_number_types.cc.
Expression Differentiation::SD::operator! | ( | const Expression & | expression | ) |
Logical not operator.
Definition at line 428 of file symengine_number_types.cc.
Expression Differentiation::SD::operator& | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Logical and operator.
Definition at line 440 of file symengine_number_types.cc.
Expression Differentiation::SD::operator| | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Logical inclusive or operator.
Definition at line 457 of file symengine_number_types.cc.
Expression Differentiation::SD::operator^ | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Logical exclusive or (xor) operator.
Definition at line 474 of file symengine_number_types.cc.
Expression Differentiation::SD::operator&& | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
And operator.
Definition at line 491 of file symengine_number_types.cc.
Expression Differentiation::SD::operator|| | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Inclusive or operator.
Definition at line 498 of file symengine_number_types.cc.
Expression Differentiation::SD::operator+ | ( | Expression | lhs, |
const Expression & | rhs | ||
) |
Addition operator.
Return the result of adding the rhs
to the lhs
.
Definition at line 505 of file symengine_number_types.cc.
Expression Differentiation::SD::operator- | ( | Expression | lhs, |
const Expression & | rhs | ||
) |
Subtraction operator.
Return the result of subtracting the rhs
from the lhs
.
Definition at line 513 of file symengine_number_types.cc.
Expression Differentiation::SD::operator* | ( | Expression | lhs, |
const Expression & | rhs | ||
) |
Multiplication operator.
Return the result of multiplying the lhs
by the rhs
.
Definition at line 520 of file symengine_number_types.cc.
Expression Differentiation::SD::operator/ | ( | Expression | lhs, |
const Expression & | rhs | ||
) |
Division operator.
Return the result of dividing the lhs
by the rhs
.
Definition at line 528 of file symengine_number_types.cc.
|
inline |
General addition operator.
Return the result of adding the rhs
to the lhs
. The lhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1020 of file symengine_number_types.h.
|
inline |
General addition operator.
Return the result of adding the rhs
to the lhs
. The rhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1037 of file symengine_number_types.h.
|
inline |
General subtraction operator.
Return the result of subtracting the rhs
from the lhs
. The lhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1054 of file symengine_number_types.h.
|
inline |
General subtraction operator.
Return the result of subtracting the rhs
from the lhs
. The rhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1071 of file symengine_number_types.h.
|
inline |
General multiplication operator.
Return the result of multiplying the lhs
by the rhs
. The lhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1087 of file symengine_number_types.h.
|
inline |
General multiplication operator.
Return the result of multiplying the lhs
by the rhs
. The rhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1103 of file symengine_number_types.h.
|
inline |
General division operator.
Return the result of dividing the lhs
by the rhs
. The lhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1120 of file symengine_number_types.h.
|
inline |
General division operator.
Return the result of dividing the lhs
by the rhs
. The lhs
type may be any supported number type, and the result is promoted to a Expression. The type conversion makes writing scalar expressions using Expression more natural.
Definition at line 1137 of file symengine_number_types.h.
Expression Differentiation::SD::make_symbol | ( | const std::string & | symbol | ) |
Return an Expression representing a scalar symbolic variable with the identifier specified by symbol
.
For example, if the symbol
is the string "x"
then the scalar symbolic variable that is returned represents the scalar \(x\).
[in] | symbol | An identifier (or name) for the returned symbolic variable. |
symbol
.Definition at line 41 of file symengine_scalar_operations.cc.
Expression Differentiation::SD::make_symbolic_function | ( | const std::string & | symbol, |
const types::symbol_vector & | arguments | ||
) |
Return an Expression representing a scalar symbolic function with the identifier specified by symbol
. The function's symbolic dependencies are specified by the input arguments
.
For example, if the symbol
is the string "f"
, and the arguments to the function that is generated are the symbolic variable x
and the symbolic expression y+z
, then the generic symbolic function that is returned represents \(f(x, y+z)\).
[in] | symbol | An identifier (or name) for the returned symbolic function. |
[in] | arguments | A vector of input arguments to the returned symbolic function. |
symbolic
and the number of input arguments equal to the length of arguments
.Definition at line 48 of file symengine_scalar_operations.cc.
Expression Differentiation::SD::make_symbolic_function | ( | const std::string & | symbol, |
const types::substitution_map & | arguments | ||
) |
Return an Expression representing a scalar symbolic function with the identifier specified by symbol
. The function's symbolic dependencies are specified by the keys to the input arguments
map; the values stored in the map are ignored.
For example, if the symbol
is the string "f"
, and the arguments to the function that is generated are the symbolic variable x
and the symbolic expression y+z
, then the generic symbolic function that is returned represents \(f(x, y+z)\).
[in] | symbol | An identifier (or name) for the returned symbolic function. |
[in] | arguments | A map of input arguments to the returned symbolic function. |
symbolic
and the number of input arguments equal to the length of arguments
.Definition at line 56 of file symengine_scalar_operations.cc.
Expression Differentiation::SD::differentiate | ( | const Expression & | f, |
const Expression & | x | ||
) |
Return the symbolic result of computing the partial derivative of the scalar f
with respect to the scalar x
. In most use cases the function or variable f
would be called the dependent variable, and x
the independent variable.
[in] | f | A scalar symbolic function or (dependent) expression. |
[in] | x | A scalar symbolic (independent) variable. |
Definition at line 68 of file symengine_scalar_operations.cc.
types::substitution_map Differentiation::SD::make_symbol_map | ( | const SymbolicType & | symbol | ) |
Return a symbolic map that has a single entry with the key given by the symbol
. It is expected that all entries to be added to the symbolic map are valid symbols or symbolic expressions.
ignore_invalid_symbols | See the add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | Any symbolic type that is understood by the add_to_symbol_map() functions. This includes individual Expression, std::vector<Expression>, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
types::substitution_map Differentiation::SD::make_symbol_map | ( | const SymbolicType & | symbol, |
const Args &... | other_symbols | ||
) |
Return a symbolic map that has the entry keys given by symbol
and all other_symbols
. It is expected that all entries to be added to the symbolic map are valid symbols or symbolic expressions.
With this function it is possible to construct a symbolic map from different types. An example may be as follows:
ignore_invalid_symbols | See the add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | Any symbolic type that is understood by the add_to_symbol_map() functions. This includes individual Expression, std::vector<Expression>, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
Args | A type associated with the parameter pack that contains any number of other SymbolicTypes . All types held by the parameter pack share the same restriction as the SymbolicType documented above. |
void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const Expression & | symbol | ||
) |
A convenience function for adding an empty entry, with the key value given by symbol
, to the symbolic map symbol_map
.
This function is guaranteed to create an ordering that is identical to the typical add_to_substitution_map() call that is used when constructing a map to perform symbol substitution. It exists primarily to create an initial map that can be used in the optimize() call to a BatchOptimizer, specifically if the values that are to be substituted into the map are not known at the time that the symbols used to construct symbolic expressions are defined. This helps one conform to the requirement that the arguments sent into lambda and LLVM JIT compiled functions (created by optimizing symbolic expressions) (i) be the same, and (ii) have a constant ordering.
ignore_invalid_symbols | A template parameter that enforces whether or not the symbol has to be a valid one or not. In the overwhelming majority of cases, the default value of false should be selected, with the result that an exception will be thrown if the input symbolic is, in fact, not a symbolic value or expression. An exceptional case is, for example, when performing symbolic assembly on a finite element level. When extracting the symbolic equivalent of the shape function gradients using FEExtractors, the returned tensor will have some a priori determined zero-valued components. These trivial components are not valid symbols (as they are not symbolic expressions), and we would typically wish to guard against their (erroneous) inclusion. In this scenario, for convenience, one could set ignore_invalid_symbols to true and these zero-valued entries would be skipped over and ignored. |
ValueType
is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const SymbolicType & | symbol | ||
) |
A convenience function for adding an empty entry, with the key value given by symbol
, to the symbolic map symbol_map
.
For more context which this function is used, see the other add_to_symbol_map(types::substitution_map &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
T | An arbitrary type resulting from the application of the SFINAE idiom to selectively specialize this class. The required condition is fulfilled when the SymbolicType can be explicitly converted to a const SymEngine::RCP<const SymEngine::Basic> & . |
void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const std::vector< SymbolicType > & | symbols | ||
) |
A convenience function for adding empty entries, with the key values equal to the entries in symbols
, to the symbolic map symbol_map
. It is expected that all entries in the input vector symbols
be of SymbolicType
, compatible with the other add_to_symbol_map() functions.
For more context which this function is used, see the other add_to_symbol_map(types::substitution_map &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | Any symbolic type that is understood by the add_to_symbol_map() functions. This includes an individual Expression, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const types::substitution_map & | other_symbols | ||
) |
A convenience function for adding empty entries, with the key values equal to the key entries in other_symbols
, to the symbolic map symbol_map
.
For more context which this function is used, see the other add_to_symbol_map(types::substitution_map &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const SymbolicType & | symbol, | ||
const Args &... | other_symbols | ||
) |
A convenience function for adding empty entries, with the key values equal to the entries in symbol
plus other_symbols
, to the symbolic map symbol_map
. It is expected that all entries in symbol
and other_symbols
be of a SymbolicType
, compatible with the other add_to_symbol_map() functions.
For more context which this function is used, see the other add_to_symbol_map(types::substitution_map &, const Expression &) function.
With this function it is possible to add entries from different types to a symbolic map. An example may be as follows:
ignore_invalid_symbols | See the other add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | Any symbolic type that is understood by the add_to_symbol_map() functions. This includes individual Expression, std::vector<Expression>, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
Args | A type associated with the parameter pack that contains any number of other SymbolicTypes . All types held by the parameter pack share the same restriction as the SymbolicType documented above. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const Expression & | symbol, | ||
const Expression & | value | ||
) |
Find the entry for symbol
in the substitution_map
and set its corresponding value
.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
Definition at line 147 of file symengine_scalar_operations.cc.
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const SymbolicType & | symbol, | ||
const ValueType & | value | ||
) |
Find the entry for symbol
in the substitution_map
and set its corresponding value
.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
T | An arbitrary type resulting from the application of the SFINAE idiom to selectively specialize this class. The required condition is fulfilled when the SymbolicType can be explicitly converted to a const SymEngine::RCP<const SymEngine::Basic> & , and it is possible to construct an SymbolicType directly from the ValueType . |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const std::vector< SymbolicType > & | symbols, | ||
const std::vector< ValueType > & | values | ||
) |
Find the entries for symbols
in the substitution_map
and set their corresponding values
.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined SymbolicType can be constructed from. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const std::pair< SymbolicType, ValueType > & | symbol_value | ||
) |
Find the entry for symbols
in the substitution_map
and set their corresponding values
. The modified symbol will have the key given by the first element of symbol_value
and the value given by its second element.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined SymbolicType can be constructed from. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const std::pair< SymbolicType, ValueType > & | symbol_value, | ||
const Args &... | other_symbol_values | ||
) |
Find the entries for symbols
in the substitution_map
and set their corresponding values
, followed by the same operation for the other_symbol_values
.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined SymbolicType can be constructed from. |
Args | A type associated with the parameter pack that contains any number of other pairs of SymbolicTypes and ValueTypes . All types held by the parameter pack share the same restriction as the SymbolicType and ValueType documented above. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const std::vector< std::pair< SymbolicType, ValueType >> & | symbol_values | ||
) |
Find the entries for symbols
in the substitution_map
and set their corresponding values
. The modified symbol will have the key given by the first element of each paired entry in the symbol_values
vector and the value given by its respective second element.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined SymbolicType can be constructed from. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const types::substitution_map & | symbol_values | ||
) |
Find the entries for symbols
in the substitution_map
and set their corresponding values
. The modified symbol will have the key given by the each element the symbol_values
map and the value given by its respective mapped element.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
Definition at line 158 of file symengine_scalar_operations.cc.
types::substitution_map Differentiation::SD::make_substitution_map | ( | const Expression & | symbol, |
const Expression & | value | ||
) |
Return a substitution map that has the entry key given by symbol
and the value given by value
. It is expected that the key entry be valid symbol or symbolic expression.
The values that map to a symbol
would typically be of an arithmetic type. However, in some instances is may be useful to map a symbolic type to another symbolic type (i.e. perform partial substitution). In such a situation the resolve_explicit_dependencies() function may be useful to simplify the final substitution map by resolving all explicit interdependencies between entries in the substitution map.
Definition at line 170 of file symengine_scalar_operations.cc.
types::substitution_map Differentiation::SD::make_substitution_map | ( | const ExpressionType & | symbol, |
const ValueType & | value | ||
) |
Return a substitution map that has the entry key given by symbol
and the value given by value
. It is expected that the key entry be valid symbol or symbolic expression.
The values that map to a symbol
would typically be of a ValueType
(i.e., an arithmetic type). However, in some instances is may be useful to map a symbolic type to another symbolic type (i.e. perform partial substitution). In such a situation the resolve_explicit_dependencies() function may be useful to simplify the final substitution map by resolving all explicit interdependencies between entries in the substitution map.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
T | An arbitrary type resulting from the application of the SFINAE idiom to selectively specialize this class. The required condition is fulfilled when the ExpressionType can be explicitly converted to a const SymEngine::RCP<const SymEngine::Basic> & , and it is possible to construct an ExpressionType directly from the ValueType . |
types::substitution_map Differentiation::SD::make_substitution_map | ( | const std::vector< ExpressionType > & | symbols, |
const std::vector< ValueType > & | values | ||
) |
Return a substitution map that has the entry keys given by symbols
and the values given by values
. It is expected that all key entries be valid symbols or symbolic expressions.
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &, const ValueType &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
types::substitution_map Differentiation::SD::make_substitution_map | ( | const std::pair< ExpressionType, ValueType > & | symbol_value | ) |
Return a substitution map that has the key given by the first entry in symbol_value
, and the value of its second entry. It is expected that the key entry be a valid symbol or symbolic expression.
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &, const ValueType &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
types::substitution_map Differentiation::SD::make_substitution_map | ( | const std::vector< std::pair< ExpressionType, ValueType >> & | symbol_values | ) |
Return a substitution map that has the keys given by the first entry of each element of symbol_values
, and the values given its second entry. It is expected that all key entries be valid symbols or symbolic expressions.
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &, const ValueType &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
types::substitution_map Differentiation::SD::make_substitution_map | ( | const std::pair< ExpressionType, ValueType > & | symbol_value, |
const Args &... | other_symbol_values | ||
) |
Return a substitution map that has the key given by the first entry in symbol_value
, and the value of its second entry, followed by the addition of the other_symbol_values
. It is expected that all key entries be valid symbols or symbolic expressions.
With this function it is possible to construct a symbolic substitution map from different types, so long as there exists a add_to_substitution_map() function with the signature corresponding to the pair types. An example may be as follows:
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &, const ValueType &) function.
ExpressionType | Any symbolic expression type that is understood by the make_substitution_map() functions. This includes individual Expression, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
Args | A type associated with the parameter pack that contains any number of other ExpressionTypes . All types held by the parameter pack share the same restriction as the ExpressionType documented above. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const Expression & | symbol, | ||
const Expression & | value | ||
) |
A convenience function to add an entry to the substitution_map
. The new entry will have the key given by symbol
with its paired value
. Such maps are required to perform substitution of symbolic expressions, with key entries being exchanges with their pair values.
This function is guaranteed to create an ordering that is identical to the typical add_to_symbol_map() call that may be used when initially configuring a BatchOptimizer. This helps one conform to the requirement that the arguments sent into lambda and LLVM JIT compiled functions (created by optimizing symbolic expressions) (i) be the same, and (ii) have a constant ordering.
ignore_invalid_symbols | A template parameter that enforces whether or not the symbol has to be a valid one or not. In the overwhelming majority of cases, the default value of false should be selected, with the result that an exception will be thrown if the input symbolic is, in fact, not a symbolic value or expression. An exceptional case is, for example, when performing symbolic assembly on a finite element level. When extracting the symbolic equivalent of the shape function gradients using FEExtractors, the returned tensor will have some a priori determined zero-valued components. These trivial components are not valid symbols (as they are not symbolic expressions), and we would typically wish to guard against their (erroneous) inclusion. In this scenario, for convenience, one could set ignore_invalid_symbols to true and these zero-valued entries would be skipped over and ignored. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const ExpressionType & | symbol, | ||
const ValueType & | value | ||
) |
A convenience function to add an entry to the substitution_map
. The new entry will have the key given by symbol
with its paired value
.
The ExpressionType
will be used to convert the value
to a compatible SymEngine number type. It is therefore required that the ExpressionType
ValueType
, and thatconst SymEngine::RCP<const SymEngine::Basic> &
.For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const std::vector< ExpressionType > & | symbols, | ||
const std::vector< ValueType > & | values | ||
) |
A convenience function for adding multiple entries to the substitution_map
. The new entries will have the keys given in the symbols
vector, each of which will be paired index-wise with its corresponding element in the values
vector.
The class represented by the ExpressionType
template parameter will be used to convert the p value to a compatible SymEngine number type. It is therefore required that the ExpressionType
ValueType
, and thatconst SymEngine::RCP<const SymEngine::Basic> &
.For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const types::substitution_map & | symbol_values | ||
) |
A convenience function for adding multiple entries to the substitution_map
. The new entries will have the keys given in the symbols
vector, each of which will be paired index-wise with its corresponding element in the values
vector. It is expected that there are no duplicate entries between the two maps.
For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const std::pair< ExpressionType, ValueType > & | symbol_value | ||
) |
A convenience function to add an entry to the substitution_map
. The new entry will have the key given by the first element of symbol_value
and the value given by its second element. It is expected that the key entry be a valid symbol or symbolic expression, and that the paired symbol_value
elements are compatible with the other add_to_substitution_map() functions.
The ExpressionType
and its associated ValueType
need not be scalar types. So, for example, this function could be used to add tensor-valued data to the map in the following way:
For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
ExpressionType | Any symbolic expression type that is understood by the add_to_substitution_map() functions. This includes individual Expression, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const std::vector< std::pair< ExpressionType, ValueType >> & | symbol_values | ||
) |
A convenience function for adding multiple entries to the substitution_map
. The new entries will have the keys given by first entry of each element of symbol_values
, and the values given its second entry. It is expected that the key entry be a valid symbols or symbolic expressions, and that the paired symbol_value
elements are compatible with the other add_to_substitution_map() functions.
The ExpressionType
and its associated ValueType
need not be scalar types. So, for example, this function could be used to add tensor-valued data to the map in the following way:
For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
ExpressionType | Any symbolic expression type that is understood by the add_to_substitution_map() functions. This includes individual Expression, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const std::pair< ExpressionType, ValueType > & | symbol_value, | ||
const Args &... | other_symbol_values | ||
) |
A convenience function for adding multiple entries to the substitution_map
. The new entries will have the keys given by first entry of each element of symbol_values
, and the values given its second entry, along with the addition of the other_symbol_values
. It is expected that the key entry be a valid symbols or symbolic expressions, and that the paired symbol_value
elements are compatible with the other add_to_substitution_map() functions.
With this function it is possible to construct a symbolic substitution map from different types, so long as there exists a add_to_substitution_map() function with the signature corresponding to the pair types. An example may be as follows:
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &, const ValueType &) function.
ExpressionType | Any symbolic expression type that is understood by the add_to_substitution_map() functions. This includes individual Expression, as well as Tensors and SymmetricTensors of Expressions. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
Args | A type associated with the parameter pack that contains any number of other ExpressionTypes . All types held by the parameter pack share the same restriction as the ExpressionType documented above. |
void Differentiation::SD::merge_substitution_maps | ( | types::substitution_map & | substitution_map_out, |
const types::substitution_map & | substitution_map_in | ||
) |
Concatenate two symbolic maps, merging a second map substitution_map_in
in-place into the initial and resultant map substitution_map_out
. The map substitution_map_out
need not initially be empty.
Definition at line 182 of file symengine_scalar_operations.cc.
void Differentiation::SD::merge_substitution_maps | ( | types::substitution_map & | substitution_map_out, |
const types::substitution_map & | substitution_map_in, | ||
const Args &... | other_substitution_maps_in | ||
) |
Concatenate multiple symbolic maps, merging the maps substitution_map_in
and other_substitution_maps_in
, a collection of other maps, in-place into the resultant map substitution_map_out
. The map substitution_map_out
need not initially be empty.
types::substitution_map Differentiation::SD::merge_substitution_maps | ( | const types::substitution_map & | substitution_map_in, |
const Args &... | other_substitution_maps_in | ||
) |
Concatenate multiple symbolic maps, merging the maps substitution_map_in
and other_substitution_maps_in
and returning the result.
types::substitution_map Differentiation::SD::resolve_explicit_dependencies | ( | const types::substitution_map & | substitution_map, |
const bool | force_cyclic_dependency_resolution = false |
||
) |
Return a substitution map that has any explicit interdependencies between the entries of the input substitution_map
resolved.
The force_cyclic_dependency_resolution
flag exists to ensure, if desired, that no cyclic dependencies can exist in the returned map. If a cyclic dependency exists in the input substitution map, substitution_map
, then with this flag set to true
the dependency cycle is broken by a dictionary-ordered substitution. For example, if the substitution map contains two entries map["a"] -> "b"
and map["b"] -> "a"
, then the result of calling this function would be a map with the elements map["a"] -> "a"
and map["b"] -> "a"
.
If one symbol is an explicit function of another, and it is desired that all their values are completely resolved, then it may be necessary to perform substitution a number of times before the result is finalized. This function performs substitution sweeps for a set of symbolic variables until all explicit relationships between the symbols in the map have been resolved. Whether each entry returns a symbolic or real value depends on the nature of the values stored in the substitution map. If the values associated with a key are also symbolic then the returned result may still be symbolic in nature. The terminal result of using the input substitution map, symbol_values
, is then guaranteed to be rendered by a single substitition of the returned dependency-resolved map.
Example: If map["a"] -> 1
and map["b"] -> "a"+ 2
, then then the function \(f(a,b(a)) = a+b\) will be evaluated and the result \(f\vert_{a=1,b=a+2} = 3+a\) is determined upon the completion of the first sweep. A second sweep is therefore necessary to resolve the final symbol, and the returned value is ultimately \(f = [3+a]_{a=1} = 4\). By resolving the explicit relationships between all symbols in the map, we determine that map["a"] -> 1
and map["b"] -> 1 + 2 = 3
and thus, using only one substitution, that \(f = a+b = 1 + 3 = 4\).
Definition at line 206 of file symengine_scalar_operations.cc.
types::substitution_map Differentiation::SD::resolve_explicit_dependencies | ( | const std::vector< std::pair< ExpressionType, ValueType >> & | symbol_values, |
const bool | force_cyclic_dependency_resolution = false |
||
) |
Return a substitution map that has any explicit interdependencies between the entries of the map generated by the paired elements in the symbol_values
vector resolved. The force_cyclic_dependency_resolution
exists to ensure, if desired, that no cyclic dependencies can exist in the returned map.
This function performs substitution sweeps for a set of symbolic variables until all explicit relationships between the symbols in the map have been resolved. Whether each entry returns a symbolic or real value depends on the chosen ValueType and the values represented therein. If the ValueType is another Expression, and they contain symbols then the returned result may still be symbolic in nature.
For an example of what this function does, see the documentation for the other resolve_explicit_dependencies(const types::substitution_map &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
Expression Differentiation::SD::substitute | ( | const Expression & | expression, |
const types::substitution_map & | substitution_map | ||
) |
Perform a single substitution sweep of a set of symbols into the given symbolic expression. The symbols in the expression
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function may be used to give a set of symbolic variables either a numeric interpretation or some symbolic definition.
Examples:
If map["a"] == 1
and map["b"] == "a" + 2
, then the function \(f(a,b(a)) := a+b\) will be evaluated and the result \(f\vert_{a=1,b=a+2} = 3+a\) is returned. This return is because the symbol "a" is substituted throughout the function first, and only then is the symbol "b(a)" substituted, by which time its explicit dependency on "a" cannot be resolved.
map["a"] == "b"+2
and map["b"] == 1
, then the function \(f(a(b),b): = a+b\) will be evaluated and the result \(f\vert_{a=b+2, b} = [b+2+b]_{b=1} = 4\) is returned. This is because the explicitly dependent symbol "a(b)" is substituted first followed by the symbol "b". Definition at line 258 of file symengine_scalar_operations.cc.
Expression Differentiation::SD::substitute | ( | const Expression & | expression, |
const Expression & | symbol, | ||
const ValueType & | value | ||
) |
Perform a substitution of the symbol
into the given expression
. All matches are assigned the corresponding value
. This substitution function may be used to give a set of symbolic variables either a numeric interpretation or some symbolic definition.
For more information regarding the performance of symbolic substitution, see the other substitute(const Expression &, const types::substitution_map &) function.
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type. |
ExpressionType Differentiation::SD::substitute | ( | const ExpressionType & | expression, |
const Args &... | symbol_values | ||
) |
Perform a single substitution sweep of a set of symbols into the given symbolic expression. The symbols in the expression
that correspond to a matching entry key of the symbol_values
vector entry are substituted by the entry's associated value. This substitution function may be used to give a set of symbolic variables either a numeric interpretation or some symbolic definition.
For more information regarding the performance of symbolic substitution, see the other substitute(const Expression &, const types::substitution_map &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, but this may also represent Tensors and SymmetricTensors of Expressions. |
Args | Any symbolic type and value combination that is understood by the make_substitution_map() functions. This includes arguments involving individual Expressions, std::vector<Expression>, as well as Tensors and SymmetricTensors of Expressions. |
ValueType Differentiation::SD::substitute_and_evaluate | ( | const Expression & | expression, |
const types::substitution_map & | substitution_map | ||
) |
Perform a single substitution sweep of a set of symbols into the given symbolic function, and immediately evaluate the result. The symbols in the expression
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function is used to give a set of symbolic variables a numeric interpretation, with the returned result being of the type specified by the ValueType
template argument.
For more information regarding the performance of symbolic substitution, and the outcome of evaluation using a substitution map with cyclic dependencies, see the substitute(const Expression &, const types::substitution_map &) function.
expression
be successfully resolved by the substitution_map
. If only partial substitution is performed, then an error is thrown.ValueType | A type that corresponds to the value that the symbol is to represent. In the context of this particular function, this template parameter is typically arithmetic in nature. |
ValueType Differentiation::SD::substitute_and_evaluate | ( | const Expression & | expression, |
const Args &... | symbol_values | ||
) |
Perform a single substitution sweep of a set of symbols into the given symbolic function, and immediately evaluate the result. The symbols in the expression
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function is used to give a set of symbolic variables a numeric interpretation with the returned result being of the type specified by the ValueType
template argument.
For more information regarding the performance of symbolic substitution, and the outcome of evaluation using a substitution map with cyclic dependencies, see the substitute(const Expression &, const types::substitution_map &) function.
expression
be successfully resolved by the substitution map that is generated with the input collection of symbol_values
. If only partial substitution is performed, then an error is thrown.ValueType | A type that corresponds to the value that the symbol is to represent. In the context of this particular function, this template parameter is typically arithmetic in nature. |
Args | Any symbolic type and value combination that is understood by the make_substitution_map() functions. This includes arguments involving individual Expressions, std::vector<Expression>, as well as Tensors and SymmetricTensors of Expressions. |
Tensor< 1, dim, Expression > Differentiation::SD::make_vector_of_symbols | ( | const std::string & | symbol | ) |
Return a vector of Expressions representing a vectorial symbolic variable with the identifier specified by symbol
.
For example, if the symbol
is the string "v"
then the vectorial symbolic variable that is returned represents the vector \(v\). Each component of \(v\) is prefixed by the given symbol
, and has a suffix that indicates its component index.
dim | The dimension of the returned tensor. |
[in] | symbol | An identifier (or name) for the vector of returned symbolic variables. |
symbol
.Definition at line 237 of file symengine_tensor_operations.cc.
Tensor< rank, dim, Expression > Differentiation::SD::make_tensor_of_symbols | ( | const std::string & | symbol | ) |
Return a tensor of Expressions representing a tensorial symbolic variable with the identifier specified by symbol
.
For example, if the symbol
is the string "T"
then the tensorial symbolic variable that is returned represents the vector \(T\). Each component of \(T\) is prefixed by the given symbol
, and has a suffix that indicates its component indices.
rank | The rank of the returned tensor. |
dim | The dimension of the returned tensor. |
[in] | symbol | An identifier (or name) for the tensor of returned symbolic variables. |
symbol
.Definition at line 254 of file symengine_tensor_operations.cc.
SymmetricTensor< rank, dim, Expression > Differentiation::SD::make_symmetric_tensor_of_symbols | ( | const std::string & | symbol | ) |
Return a symmetric tensor of Expressions representing a tensorial symbolic variable with the identifier specified by symbol
.
For example, if the symbol
is the string "S"
then the tensorial symbolic variable that is returned represents the vector \(S\). Each component of \(S\) is prefixed by the given symbol
, and has a suffix that indicates its component indices.
rank | The rank of the returned tensor. |
dim | The dimension of the returned tensor. |
[in] | symbol | An identifier (or name) for the tensor of returned symbolic variables. |
symbol
.Definition at line 262 of file symengine_tensor_operations.cc.
Tensor< 1, dim, Expression > Differentiation::SD::make_vector_of_symbolic_functions | ( | const std::string & | symbol, |
const types::substitution_map & | arguments | ||
) |
Return a vector of Expression representing a vectorial symbolic function with the identifier specified by symbol
. The functions' symbolic dependencies are specified by the keys to the input arguments
map; the values stored in the map are ignored.
dim | The dimension of the returned tensor. |
[in] | symbol | An identifier (or name) for the vector of returned symbolic functions. |
[in] | arguments | A map of input arguments to the returned symbolic functions. |
symbol
, a suffix that indicates its component index, and the number of input arguments equal to the length of arguments
.Definition at line 245 of file symengine_tensor_operations.cc.
Tensor< rank, dim, Expression > Differentiation::SD::make_tensor_of_symbolic_functions | ( | const std::string & | symbol, |
const types::substitution_map & | arguments | ||
) |
Return a tensor of Expression representing a tensorial symbolic function with the identifier specified by symbol
. The functions' symbolic dependencies are specified by the keys to the input arguments
map; the values stored in the map are ignored.
rank | The rank of the returned tensor. |
dim | The dimension of the returned tensor. |
[in] | symbol | An identifier (or name) for the tensor of returned symbolic functions. |
[in] | arguments | A map of input arguments to the returned symbolic functions. |
symbol
, a suffix that indicates its component indeices, and the number of input arguments equal to the length of arguments
.Definition at line 270 of file symengine_tensor_operations.cc.
SymmetricTensor< rank, dim, Expression > Differentiation::SD::make_symmetric_tensor_of_symbolic_functions | ( | const std::string & | symbol, |
const types::substitution_map & | arguments | ||
) |
Return a symmetric tensor of Expression representing a tensorial symbolic function with the identifier specified by symbol
. The functions' symbolic dependencies are specified by the keys to the input arguments
map; the values stored in the map are ignored.
rank | The rank of the returned tensor. |
dim | The dimension of the returned tensor. |
[in] | symbol | An identifier (or name) for the tensor of returned symbolic functions. |
[in] | arguments | A map of input arguments to the returned symbolic functions. |
symbol
, a suffix that indicates its component indeices, and the number of input arguments equal to the length of arguments
.Definition at line 280 of file symengine_tensor_operations.cc.
Tensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const Expression & | f, |
const Tensor< rank, dim, Expression > & | T | ||
) |
Return the symbolic result of computing the partial derivative of the scalar f
with respect to the tensor T
.
[in] | f | A scalar symbolic function or (dependent) expression. |
[in] | T | A tensor of symbolic (independent) variables. |
SymmetricTensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const Expression & | f, |
const SymmetricTensor< rank, dim, Expression > & | S | ||
) |
Return the symbolic result of computing the partial derivative of the scalar f
with respect to the symmetric tensor S
.
[in] | f | A scalar symbolic function or (dependent) expression. |
[in] | S | A symmetric tensor of symbolic (independent) variables. |
Tensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const Tensor< 0, dim, Expression > & | f, |
const Tensor< rank, dim, Expression > & | T | ||
) |
Return the symbolic result of computing the partial derivative of the rank-0 tensor (or scalar) f
with respect to the tensor T
.
[in] | f | A rank-0 tensor symbolic function or (dependent) expression. |
[in] | T | A tensor of symbolic (independent) variables. |
SymmetricTensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const Tensor< 0, dim, Expression > & | f, |
const SymmetricTensor< rank, dim, Expression > & | S | ||
) |
Return the symbolic result of computing the partial derivative of the rank-0 tensor (or scalar) f
with respect to the symmetric tensor S
.
[in] | f | A rank-0 tensor symbolic function or (dependent) expression. |
[in] | S | A symmetric tensor of symbolic (independent) variables. |
Tensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const Tensor< rank, dim, Expression > & | T, |
const Expression & | x | ||
) |
Return the symbolic result of computing the partial derivative of the tensor T
with respect to the scalar x
.
[in] | T | A tensor of symbolic functions or (dependent) expressions. |
[in] | x | A scalar symbolic (independent) variable. |
SymmetricTensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const SymmetricTensor< rank, dim, Expression > & | S, |
const Expression & | x | ||
) |
Return the symbolic result of computing the partial derivative of the symmetric tensor S
with respect to the scalar x
.
[in] | S | A symmetric tensor of symbolic functions or (dependent) expressions. |
[in] | x | A scalar symbolic (independent) variable. |
Tensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const Tensor< rank, dim, Expression > & | T, |
const Tensor< 0, dim, Expression > & | x | ||
) |
Return the symbolic result of computing the partial derivative of the tensor T
with respect to the rank-0 tensor x
.
[in] | T | A tensor of symbolic functions or (dependent) expressions. |
[in] | x | A rank-0 tensor containing a symbolic (independent) variable. |
SymmetricTensor<rank, dim, Expression> Differentiation::SD::differentiate | ( | const SymmetricTensor< rank, dim, Expression > & | S, |
const Tensor< 0, dim, Expression > & | x | ||
) |
Return the symbolic result of computing the partial derivative of the symmetric tensor S
with respect to the rank-0 tensor x
.
[in] | S | A symmetric tensor of symbolic functions or (dependent) expressions. |
[in] | x | A rank-0 tensor containing a symbolic (independent) variable. |
Tensor<rank_1 + rank_2, dim, Expression> Differentiation::SD::differentiate | ( | const Tensor< rank_1, dim, Expression > & | T1, |
const Tensor< rank_2, dim, Expression > & | T2 | ||
) |
Return the symbolic result of computing the partial derivative of the tensor T1
with respect to the tensor T2
.
[in] | T1 | A tensor of symbolic functions or (dependent) expressions. |
[in] | T2 | A tensor of symbolic (independent) variables. |
SymmetricTensor<rank_1 + rank_2, dim, Expression> Differentiation::SD::differentiate | ( | const SymmetricTensor< rank_1, dim, Expression > & | S1, |
const SymmetricTensor< rank_2, dim, Expression > & | S2 | ||
) |
Return the symbolic result of computing the partial derivative of the symmetric tensor S1
with respect to the symmetric tensor S2
.
[in] | S1 | A symmetric tensor of symbolic functions or (dependent) expressions. |
[in] | S2 | A symmetric tensor of symbolic (independent) variables. |
Tensor<rank_1 + rank_2, dim, Expression> Differentiation::SD::differentiate | ( | const Tensor< rank_1, dim, Expression > & | T, |
const SymmetricTensor< rank_2, dim, Expression > & | S | ||
) |
Return the symbolic result of computing the partial derivative of the tensor T
with respect to the symmetric tensor S
.
[in] | T | A tensor of symbolic functions or (dependent) expressions. |
[in] | S | A symmetric tensor of symbolic (independent) variables. |
Tensor<rank_1 + rank_2, dim, Expression> Differentiation::SD::differentiate | ( | const SymmetricTensor< rank_1, dim, Expression > & | S, |
const Tensor< rank_2, dim, Expression > & | T | ||
) |
Return the symbolic result of computing the partial derivative of the symmetric tensor S
with respect to the tensor T
.
[in] | S | A symmetric tensor of symbolic functions or (dependent) expressions. |
[in] | T | A tensor of symbolic (independent) variables. |
void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const Tensor< rank, dim, SymbolicType > & | symbol_tensor | ||
) |
A convenience function for adding empty entries, with the key values equal to the entries in the symbol_tensor
, to the symbolic map symbol_map
.
For more context which this function is used, see the other add_to_symbol_map(types::substitution_map &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
void Differentiation::SD::add_to_symbol_map | ( | types::substitution_map & | symbol_map, |
const SymmetricTensor< rank, dim, SymbolicType > & | symbol_tensor | ||
) |
A convenience function for adding empty entries, with the key values equal to the entries in the symbol_tensor
, to the symbolic map symbol_map
.
For more context which this function is used, see the other add_to_symbol_map(types::substitution_map &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_symbol_map(types::substitution_map &, const Expression &) function for a detailed discussion on the role of this template argument. |
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. This ValueType is somewhat arbitrary as it is only used to create default-constructed values as entries in the map. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const Tensor< rank, dim, SymbolicType > & | symbol_tensor, | ||
const Tensor< rank, dim, ValueType > & | value_tensor | ||
) |
Find the input symbols
in the substitution_map
and set the entries corresponding to the key values given by symbol_tensor
to the values given by value_tensor
.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::set_value_in_symbol_map | ( | types::substitution_map & | substitution_map, |
const SymmetricTensor< rank, dim, SymbolicType > & | symbol_tensor, | ||
const SymmetricTensor< rank, dim, ValueType > & | value_tensor | ||
) |
Find the input symbols
in the substitution_map
and set the entries corresponding to the key values given by symbol_tensor
to the values given by value_tensor
.
This function may be used to safely transform an existing or null symbolic map (one with uninitialized entries) into one that can be used to conduct symbolic substitution operations (i.e., a substitution map).
SymbolicType | A type that represents a symbolic variable. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined SymbolicType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
types::substitution_map Differentiation::SD::make_substitution_map | ( | const Tensor< rank, dim, ExpressionType > & | symbol_tensor, |
const Tensor< rank, dim, ValueType > & | value_tensor | ||
) |
Return a substitution map that has the entry keys given by the symbol_tensor
and the values given by the value_tensor
. It is expected that all key entries be valid symbols or symbolic expressions.
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &,const ValueType &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
types::substitution_map Differentiation::SD::make_substitution_map | ( | const SymmetricTensor< rank, dim, ExpressionType > & | symbol_tensor, |
const SymmetricTensor< rank, dim, ValueType > & | value_tensor | ||
) |
Return a substitution map that has the entry keys given by the symbol_tensor
and the values given by the value_tensor
. It is expected that all key entries be valid symbols or symbolic expressions.
It is possible to map symbolic types to other symbolic types using this function. For more details on this, see the other make_substitution_map(const Expression &,const ValueType &) function.
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const Tensor< rank, dim, ExpressionType > & | symbol_tensor, | ||
const Tensor< rank, dim, ValueType > & | value_tensor | ||
) |
A convenience function for adding an entry to the substitution_map
. The new entries will have the keys given in the symbol_tensor
with their paired values extracted from the corresponding elements of the value_tensor
.
For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
void Differentiation::SD::add_to_substitution_map | ( | types::substitution_map & | substitution_map, |
const SymmetricTensor< rank, dim, ExpressionType > & | symbol_tensor, | ||
const SymmetricTensor< rank, dim, ValueType > & | value_tensor | ||
) |
A convenience function for adding an entry to the substitution_map
. The new entries will have the keys given in the symbol_tensor
with their paired values extracted from the corresponding elements of the value_tensor
.
For more context which this function is used, see the other add_to_substitution_map(types::substitution_map &,const Expression &, const Expression &) function.
ignore_invalid_symbols | See the other add_to_substitution_map(types::substitution_map &, const Expression &, const Expression &) function for a detailed discussion on the role of this template argument. |
ExpressionType | A type that represents a symbolic expression. The Differentiation::SD::Expression class is often suitable for this purpose, although if the ValueType is not supported by this class then a user-defined ExpressionType should be used. |
ValueType | A type that corresponds to the value that the symbol is to represent. Although it is typically arithmetic in nature, it may also represent another symbolic expression type or be a special type that a user-defined ExpressionType can be constructed from. |
Tensor<rank, dim, Expression> Differentiation::SD::substitute | ( | const Tensor< rank, dim, Expression > & | expression_tensor, |
const types::substitution_map & | substitution_map | ||
) |
Perform a single substitution sweep of a set of symbols into the given tensor of symbolic expressions. The symbols in the expression_tensor
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function may be used to give a set of symbolic variables either a numeric interpretation or some symbolic definition.
For more information regarding the performance of symbolic substitution, and the outcome of evaluation using a substitution map with cyclic dependencies, see the substitute(const Expression &, const types::substitution_map &) function.
SymmetricTensor<rank, dim, Expression> Differentiation::SD::substitute | ( | const SymmetricTensor< rank, dim, Expression > & | expression_tensor, |
const types::substitution_map & | substitution_map | ||
) |
Perform a single substitution sweep of a set of symbols into the given symmetric tensor of symbolic expressions. The symbols in the expression_tensor
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function may be used to give a set of symbolic variables either a numeric interpretation or some symbolic definition.
For more information regarding the performance of symbolic substitution, and the outcome of evaluation using a substitution map with cyclic dependencies, see the substitute(const Expression &, const types::substitution_map &) function.
Tensor<rank, dim, ValueType> Differentiation::SD::substitute_and_evaluate | ( | const Tensor< rank, dim, Expression > & | expression_tensor, |
const types::substitution_map & | substitution_map | ||
) |
Perform a single substitution sweep of a set of symbols into the given tensor of symbolic expressions, and immediately evaluate the tensorial result. The symbols in the expression_tensor
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function is used to give a set of symbolic variables a numeric interpretation with the returned result being of the type specified by the ValueType
template argument.
For more information regarding the performance of symbolic substitution, and the outcome of evaluation using a substitution map with cyclic dependencies, see the substitute(const Expression &, const types::substitution_map &) function.
expression_tensor
be successfully resolved by the substitution_map
. If only partial substitution is performed, then an error is thrown.ValueType | A type that corresponds to the value that the symbol is to represent. In the context of this particular function, this template parameter is typically arithmetic in nature. |
SymmetricTensor<rank, dim, ValueType> Differentiation::SD::substitute_and_evaluate | ( | const SymmetricTensor< rank, dim, Expression > & | expression_tensor, |
const types::substitution_map & | substitution_map | ||
) |
Perform a single substitution sweep of a set of symbols into the given symmetric tensor of symbolic expressions, and immediately evaluate the tensorial result. The symbols in the expression_tensor
that correspond to the entry keys of the substitution_map
are substituted with the map entry's associated value. This substitution function is used to give a set of symbolic variables a numeric interpretation, with the returned result being of the type specified by the ValueType
template argument.
For more information regarding the performance of symbolic substitution, and the outcome of evaluation using a substitution map with cyclic dependencies, see the substitute(const Expression &, const types::substitution_map &) function.
expression_tensor
be successfully resolved by the substitution_map
. If only partial substitution is performed, then an error is thrown.ValueType | A type that corresponds to the value that the symbol is to represent. In the context of this particular function, this template parameter is typically arithmetic in nature. |
Expression Differentiation::SD::operator& | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
Logical and operator.
Definition at line 440 of file symengine_number_types.cc.
Expression Differentiation::SD::operator&& | ( | const Expression & | lhs, |
const Expression & | rhs | ||
) |
And operator.
Definition at line 491 of file symengine_number_types.cc.