DACE 2.0 API Manual
Differential Algebra Core Engine
|
#include <dace/DA.h>
Public Member Functions | |
DA () | |
Default constructor. More... | |
DA (const DA &da) | |
Copy constructor. More... | |
DA (DA &&da) | |
Move constructor. More... | |
DA (const int i, const double c=1.0) | |
Constructor for DA identities. More... | |
DA (const unsigned int i, const double c=1.0) | |
Constructor for DA identities. More... | |
DA (const double c) | |
Constructor for DA constants. More... | |
~DA () throw () | |
Destructor. More... | |
int | isnan () const |
int | isinf () const |
double | cons () const |
Get constant part of a DA. More... | |
AlgebraicVector< double > | linear () const |
Get linear part of a DA. More... | |
AlgebraicVector< DA > | gradient () const |
Gradient vector with respect to all independent DA variables. More... | |
double | getCoefficient (const std::vector< unsigned int > &jj) const |
Get specific coefficient. More... | |
void | setCoefficient (const std::vector< unsigned int > &jj, const double coeff) |
Set specific coefficient. More... | |
Monomial | getMonomial (const unsigned int npos) const |
Get the Monomial at given position. More... | |
void | getMonomial (const unsigned int npos, Monomial &m) const |
Extract the Monomial at given position. More... | |
std::vector< Monomial > | getMonomials () const |
Get std::vector of all non-zero Monomials. More... | |
DA & | operator= (DA &&da) |
Move assignment from existing DA. More... | |
DA & | operator= (const DA &da) |
Assignment from existing DA. More... | |
DA & | operator= (const double c) |
Assignment from existing double. More... | |
DA & | operator+= (const DA &da) |
Add another DA to this DA. More... | |
DA & | operator+= (const double c) |
Add another double to this DA. More... | |
DA & | operator-= (const DA &da) |
Subtract another DA from this DA. More... | |
DA & | operator-= (const double c) |
Subtract another double from this DA. More... | |
DA & | operator*= (const DA &da) |
Multiply this DA by another DA. More... | |
DA & | operator*= (const double c) |
Multiply this DA by another double. More... | |
DA & | operator/= (const DA &da) |
Divide this DA by another DA. More... | |
DA & | operator/= (const double c) |
Divide this DA by another double. More... | |
DA | operator- () const |
Negate this DA. More... | |
DA | multiplyMonomials (const DA &da) const |
Multiply the DA with the argument da monomial by monomial (i.e. coefficient-wise) More... | |
DA | divide (const unsigned int var, const unsigned int p=1) const |
Divide by an independent variable to some power. More... | |
DA | deriv (const unsigned int i) const |
Derivative with respect to given variable. More... | |
DA | deriv (const std::vector< unsigned int > ind) const |
Derivative with respect to given variables. More... | |
DA | integ (const unsigned int i) const |
Integral with respect to given variable. More... | |
DA | integ (const std::vector< unsigned int > ind) const |
Integral with respect to given variables. More... | |
DA | trim (const unsigned int min, const unsigned int max=DA::getMaxOrder()) const |
Trim the coefficients only include orders between min and max, inclusively. More... | |
DA | trunc () const |
Truncate the constant part to an integer. More... | |
DA | round () const |
Round the constant part to an integer. More... | |
DA | mod (const double p) const |
Modulo of the constant part. More... | |
DA | pow (const int p) const |
Exponentiation to given (integer) power. More... | |
DA | pow (const double p) const |
Exponentiation to given double power. More... | |
DA | root (const int p=2) const |
p-th root More... | |
DA | minv () const |
multiplicative inverse More... | |
DA | sqr () const |
square More... | |
DA | sqrt () const |
square root More... | |
DA | isrt () const |
inverse square root More... | |
DA | cbrt () const |
cubic root More... | |
DA | icrt () const |
inverse cubic root More... | |
DA | hypot (const DA &da) const |
hypotenuse More... | |
DA | exp () const |
exponential More... | |
DA | log () const |
natural logarithm More... | |
DA | logb (const double b=10.0) const |
logarithm with respect to given base More... | |
DA | log10 () const |
10 based logarithm More... | |
DA | log2 () const |
2 based logarithm More... | |
DA | sin () const |
sine More... | |
DA | cos () const |
cosine More... | |
DA | tan () const |
tangent More... | |
DA | asin () const |
arcsine More... | |
DA | acos () const |
arccosine More... | |
DA | atan () const |
arctangent More... | |
DA | atan2 (const DA &da) const |
arctangent in [-pi, pi] More... | |
DA | sinh () const |
hyperbolic sine More... | |
DA | cosh () const |
hyperbolic cosine More... | |
DA | tanh () const |
hyperbolic tangent More... | |
DA | asinh () const |
hyperbolic arcsine More... | |
DA | acosh () const |
hyperbolic arccosine More... | |
DA | atanh () const |
hyperbolic arctangent More... | |
DA | erf () const |
error function More... | |
DA | erfc () const |
complementary error function More... | |
DA | BesselJFunction (const int n) const |
Bessel function J_n. More... | |
DA | BesselYFunction (const int n) const |
Bessel function Y_n. More... | |
DA | BesselIFunction (const int n, const bool scaled=false) const |
Bessel function I_n. More... | |
DA | BesselKFunction (const int n, const bool scaled=false) const |
Bessel function K_n. More... | |
DA | GammaFunction () const |
Gamma function. More... | |
DA | LogGammaFunction () const |
Logartihmic Gamma function. More... | |
DA | PsiFunction (const unsigned int n) const |
Psi function. More... | |
unsigned int | size () const |
Number of non-zero coefficients. More... | |
double | abs () const |
Maximum absolute value of all coefficients. More... | |
double | norm (const unsigned int type=0) const |
Different types of norms over all coefficients. More... | |
std::vector< double > | orderNorm (const unsigned int var=0, const unsigned int type=0) const |
Different types of norms over coefficients of each order separately. More... | |
std::vector< double > | estimNorm (const unsigned int var=0, const unsigned int type=0, const unsigned int nc=DA::getMaxOrder()) const |
Estimate of different types of order sorted norms. More... | |
std::vector< double > | estimNorm (std::vector< double > &err, const unsigned int var=0, const unsigned int type=0, const unsigned int nc=DA::getMaxOrder()) const |
Estimate of different types of order sorted norms with error of the estimate. More... | |
Interval | bound () const |
Estimate range bound over [-1,1] for each independent variable. More... | |
double | convRadius (const double eps, const unsigned int type=1) const |
Estimate the convergence radius of the current DA. More... | |
template<class T > | |
T | eval (const std::vector< T > &args) const |
Evaluation with a vector of arguments (not efficient for repeated evaluation!) More... | |
template<class T > | |
T | eval (const T args[], const unsigned int length) const |
Evaluation with an array of arguments (not efficient for repeated evaluation!) More... | |
template<class T > | |
T | evalScalar (const T &arg) const |
Evaluation with a single arithmetic type T argument (not efficient for repeated evaluation!) More... | |
compiledDA | compile () const |
Compile current DA for efficient repeated evaluation. More... | |
DA | plug (const unsigned int var, const double val=0.0) const |
Partial evaluation to replace given independent DA variable by value val. More... | |
double | evalMonomials (const DA &values) const |
evaluate the DA providing the value of every monomial in da More... | |
DA | replaceVariable (const unsigned int from=0, const unsigned int to=0, const double val=1.0) const |
Replace variable number from by val times variable number to. More... | |
DA | scaleVariable (const unsigned int var=0, const double val=1.0) const |
Scale variable number var by val. More... | |
DA | translateVariable (const unsigned int var=0, const double a=1.0, const double c=0.0) const |
Translate variable number var to a*var + c. More... | |
std::string | toString () const |
Convert to string representation. More... | |
void | write (std::ostream &os) const |
Write binary representation of DA to stream. More... | |
Static Public Member Functions | |
static void | init (const unsigned int ord, const unsigned int nvar) |
DACE initialization. More... | |
static bool | isInitialized () |
Get DACE initialization status. More... | |
static void | version (int &maj, int &min, int &patch) |
DACE core routines version. More... | |
static void | checkVersion () |
Check DACE C++ interface and core routines version for compatibility. More... | |
static double | setEps (const double eps) |
Set truncation epsilon. More... | |
static double | getEps () |
Get truncation epsilon. More... | |
static double | getEpsMac () |
Get the machine epsilon. More... | |
static unsigned int | getMaxOrder () |
Get the maximum order. More... | |
static unsigned int | getMaxVariables () |
Get the maximum number of variables. More... | |
static unsigned int | getMaxMonomials () |
Get the maximum monomials. More... | |
static unsigned int | getTO () |
Get truncation order. More... | |
static unsigned int | setTO (const unsigned int ot=DA::getMaxOrder()) |
Set truncation order. More... | |
static void | pushTO (const unsigned int ot=DA::getMaxOrder()) |
Set truncation order and save previous value. More... | |
static void | popTO () |
Restore truncation order to value before last DA::pushTO() More... | |
static DA | random (const double cm) |
Create new DA filled with random coefficients. More... | |
static DA | identity (const unsigned int var) |
Create the DA identity for independent DA variable var. More... | |
static DA | fromString (const std::string &str) |
Create new DA from string. More... | |
static DA | fromString (const std::vector< std::string > &str) |
Create new DA from vector of strings. More... | |
static DA | read (std::istream &is) |
Create new DA from blob read from binary file. More... | |
static void | memdump () |
Dump memory status to screen (for debugging purposes only) More... | |
Friends | |
class | compiledDA |
class | storedDA |
DA DACE_API | operator+ (const DA &da1, const DA &da2) |
Addition between two DAs. More... | |
DA DACE_API | operator+ (const DA &da, const double c) |
Addition between a DA and a constant. More... | |
DA DACE_API | operator+ (const double c, const DA &da) |
Addition between a constant and a DA. More... | |
DA DACE_API | operator- (const DA &da1, const DA &da2) |
Subtraction between two DAs. More... | |
DA DACE_API | operator- (const DA &da, const double c) |
Subtraction between a DA and a constant. More... | |
DA DACE_API | operator- (const double c, const DA &da) |
Subtraction between a constant and a DA. More... | |
DA DACE_API | operator* (const DA &da1, const DA &da2) |
Multiplication between two DAs. More... | |
DA DACE_API | operator* (const DA &da, const double c) |
Multiplication between a DA and a constant. More... | |
DA DACE_API | operator* (const double c, const DA &da) |
Multiplication between a constant and a DA. More... | |
DA DACE_API | operator/ (const DA &da1, const DA &da2) |
Division between two DAs. More... | |
DA DACE_API | operator/ (const DA &da, const double c) |
Division between a DA and a constant. More... | |
DA DACE_API | operator/ (const double c, const DA &da) |
Division between a constant and a DA. More... | |
DACE_API std::ostream & | operator<< (std::ostream &out, const DA &da) |
Output to C++ stream in text form. More... | |
DACE_API std::istream & | operator>> (std::istream &in, DA &da) |
Input from C++ stream in text form. More... | |
Basic DA class representing a single polynomial.
DACE::DA::DA | ( | ) |
Default constructor.
Create an empty DA object representing the constant zero function.
DACE::DACEException |
DACE::DA::DA | ( | const DA & | da | ) |
Copy constructor.
Create a copy of a DA object.
[in] | da | object to be copied. |
DACE::DACEException |
DACE::DA::DA | ( | DA && | da | ) |
Move constructor.
Move all resources from a DA object to us without copying.
[in] | da | object to be moved. |
DACE::DACEException |
|
explicit |
Constructor for DA identities.
Create a DA object as c times the independent variable number i.
[in] | i | independent variable number (i=0 means the constant part). |
[in] | c | coefficient corresponding to the given independent variable. By default, this value is assumed to be 1.0. |
DACE::DACEException |
|
explicit |
Constructor for DA identities.
Create a DA object as c times the independent variable number i.
[in] | i | independent variable number (i=0 means the constant part). |
[in] | c | coefficient corresponding to the given independent variable. By default, this value is assumed to be 1.0. |
DACE::DACEException |
DACE::DA::DA | ( | const double | c | ) |
Constructor for DA constants.
Create a DA object with the constant part equal to c.
[in] | c | A double value set as constant part of the DA object. |
DACE::DACEException |
DACE::DA::~DA | ( | ) | ||
throw | ( | |||
) |
double DACE::DA::abs | ( | ) | const |
Maximum absolute value of all coefficients.
Compute the max norm of a DA object.
DACE::DACEException |
DA DACE::DA::acos | ( | ) | const |
arccosine
Compute the arccosine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::acosh | ( | ) | const |
hyperbolic arccosine
Compute the hyperbolic arccosine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::asin | ( | ) | const |
arcsine
Compute the arcsine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::asinh | ( | ) | const |
hyperbolic arcsine
Compute the hyperbolic arcsine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::atan | ( | ) | const |
arctangent
Compute the arctangent of a DA object. The result is copied in a new DA object.
DACE::DACEException |
arctangent in [-pi, pi]
Compute the four-quadrant arctangent of Y/X. Y is the current DA object, whereas X is the given da. The result is copied in a new DA object.
[in] | da | DA object |
DACE::DACEException |
DA DACE::DA::atanh | ( | ) | const |
hyperbolic arctangent
Compute the hyperbolic arctangent of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::BesselIFunction | ( | const int | n, |
const bool | scaled = false |
||
) | const |
Bessel function I_n.
Compute the n-th modified Bessel function of first type I_n of a DA object. The result is copied in a new DA object.
[in] | n | order of the Bessel function |
[in] | scaled | if true, the modified Bessel function is scaled by a factor exp(-x), i.e. exp(-x)I_n(x) is returned. |
DACE::DACEException |
DA DACE::DA::BesselJFunction | ( | const int | n | ) | const |
Bessel function J_n.
Compute the n-th Bessel function of first type J_n of a DA object. The result is copied in a new DA object.
[in] | n | order of the Bessel function |
DACE::DACEException |
DA DACE::DA::BesselKFunction | ( | const int | n, |
const bool | scaled = false |
||
) | const |
Bessel function K_n.
Compute the n-th modified Bessel function of second type K_n of a DA object. The result is copied in a new DA object.
[in] | n | order of the Bessel function |
[in] | scaled | if true, the modified Bessel function is scaled by a factor exp(x), i.e. exp(x)K_n(x) is returned. |
DACE::DACEException |
DA DACE::DA::BesselYFunction | ( | const int | n | ) | const |
Bessel function Y_n.
Compute the n-th Bessel function of second type Y_n of a DA object. The result is copied in a new DA object.
[in] | n | order of the Bessel function |
DACE::DACEException |
Interval DACE::DA::bound | ( | ) | const |
Estimate range bound over [-1,1] for each independent variable.
Compute lower and upper bounds of a DA object.
DACE::DACEException |
DA DACE::DA::cbrt | ( | ) | const |
cubic root
Compute the cubic root of a DA object. The result is copied in a new DA object.
DACE::DACEException |
|
static |
Check DACE C++ interface and core routines version for compatibility.
Check the DACE core library version linked to this C++ interface against the interface version and throw an exception if the versions don't match.
DACE::DACEException |
compiledDA DACE::DA::compile | ( | ) | const |
Compile current DA for efficient repeated evaluation.
Compile current DA object and create a compiledDA object.
DACE::DACEException |
double DACE::DA::cons | ( | ) | const |
Get constant part of a DA.
Return the constant part of a DA object.
DACE::DACEException |
double DACE::DA::convRadius | ( | const double | eps, |
const unsigned int | type = 1 |
||
) | const |
Estimate the convergence radius of the current DA.
Estimate the convergence radius of the DA object.
[in] | eps | requested tolerance. |
[in] | type | type of norm (sum norm is used as default) |
DACE::DACEException |
DA DACE::DA::cos | ( | ) | const |
cosine
Compute the cosine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::cosh | ( | ) | const |
hyperbolic cosine
Compute the hyperbolic cosine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::deriv | ( | const unsigned int | i | ) | const |
Derivative with respect to given variable.
Compute the derivative of a DA object with respect to variable i. The result is copied in a new DA object.
[in] | i | variable with respect to which the derivative is calculated. |
DACE::DACEException |
DA DACE::DA::deriv | ( | const std::vector< unsigned int > | ind | ) | const |
Derivative with respect to given variables.
Compute the derivative of a DA object with respect to variables ind. The result is copied in a new DA object.
[in] | ind | vector containing the number of derivatives to take for each independent variable. If ind has fewer entries than there are independent variables, the missing entries are assumed to be zero. If ind has more entries than there are independent variables, extra values are ignored. |
DACE::DACEException |
DA DACE::DA::divide | ( | const unsigned int | var, |
const unsigned int | p = 1 |
||
) | const |
Divide by an independent variable to some power.
Divide by independent variable var raised to power p. The result is copied in a new DA object.
[in] | var | independente variable number to divide by. |
[in] | p | power of the independent variable. |
DACE::DACEException |
DA DACE::DA::erf | ( | ) | const |
error function
Compute the error function of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::erfc | ( | ) | const |
complementary error function
Compute the complementary error function of a DA object. The result is copied in a new DA object.
DACE::DACEException |
std::vector< double > DACE::DA::estimNorm | ( | const unsigned int | var = 0 , |
const unsigned int | type = 0 , |
||
const unsigned int | nc = DA::getMaxOrder() |
||
) | const |
Estimate of different types of order sorted norms.
Estimate different types of order sorted norms for terms of a DA object up to a specified order.
[in] | var | order 0: Terms are sorted by their order (Default) >0: Terms are sorted by the exponent of variable var |
[in] | type | type of norm to be computed. Possible norms are: 0: Max norm (Default) 1: Sum norm >1: Vector norm of given type |
[in] | nc | maximum order to be estimated (Default order = Max order) |
DACE::DACEException |
std::vector< double > DACE::DA::estimNorm | ( | std::vector< double > & | err, |
const unsigned int | var = 0 , |
||
const unsigned int | type = 0 , |
||
const unsigned int | nc = DA::getMaxOrder() |
||
) | const |
Estimate of different types of order sorted norms with error of the estimate.
Estimate different types of order sorted norms for terms of a DA object up to a specified order with error estimates.
[out] | err | returns the amount by which the estimate underestimates the actual ordered norm of the terms in the polynomial up to the minimum of nc or the maximum computation order. |
[in] | var | order 0: Terms are sorted by their order (Default) >0: Terms are sorted by the exponent of variable var |
[in] | type | type of norm to be computed. Possible norms are: 0: Max norm (Default) 1: Sum norm >1: Vector norm of given type |
[in] | nc | maximum order to be estimated (Default order = Max order) |
DACE::DACEException |
T DACE::DA::eval | ( | const std::vector< T > & | args | ) | const |
Evaluation with a vector of arguments (not efficient for repeated evaluation!)
Generic evaluation of the DA with a vector of arithmetic type T arguments.
[in] | args | std::vector<T> of arithmetic type T with which the DA vector is evaluated |
DACE::DACEException |
T DACE::DA::eval | ( | const T | args[], |
const unsigned int | length | ||
) | const |
Evaluation with an array of arguments (not efficient for repeated evaluation!)
Generic evaluation of the DA with an array of arithmetic type T arguments.
[in] | args | array of arithmetic type T with which the DA vector is evaluated. |
[in] | length | number of elements in the array args. |
DACE::DACEException |
double DACE::DA::evalMonomials | ( | const DA & | values | ) | const |
evaluate the DA providing the value of every monomial in da
Evaluates the DA vector using the coefficients in argument values as the values for each monomial. This is equivalent to a monomial-wise dot product of two DA vectors.
[in] | values | DA vector containing the values of each monomial |
DACE::DACEException |
T DACE::DA::evalScalar | ( | const T & | arg | ) | const |
Evaluation with a single arithmetic type T argument (not efficient for repeated evaluation!)
Generic evaluation of the DA with a single arithmetic type T argument.
[in] | arg | single variable of arithmetic type T of the first independent DA variable. |
DA DACE::DA::exp | ( | ) | const |
exponential
Compute the exponential of a DA object. The result is copied in a new DA object.
DACE::DACEException |
|
static |
Create new DA from string.
Convert a string to DA object.
[in] | str | string. |
DACE::DACEException |
|
static |
Create new DA from vector of strings.
Convert a vector of strings to DA object.
[in] | str | vector of strings, each representing one line of the input. |
DACE::DACEException |
DA DACE::DA::GammaFunction | ( | ) | const |
Gamma function.
Compute the Gamma function of a DA object. The result is copied in a new DA object.
DACE::DACEException |
double DACE::DA::getCoefficient | ( | const std::vector< unsigned int > & | jj | ) | const |
Get specific coefficient.
Return a specific coefficient of a DA object.
[in] | jj | vector of the exponents of the coefficient to retrieve. |
DACE::DACEException |
|
static |
Get truncation epsilon.
Return the cutoff value eps currently set for the computations.
DACE::DACEException |
|
static |
Get the machine epsilon.
Return the machine epsilon (pessimistic estimate).
DACE::DACEException |
|
static |
Get the maximum monomials.
Return the maximum number of monomials available with the order and number of variables specified.
DACE::DACEException |
|
static |
Get the maximum order.
Return the maximum order currently set for the computations.
DACE::DACEException |
|
static |
Get the maximum number of variables.
Return the maximum number of variables set for the computations.
DACE::DACEException |
Monomial DACE::DA::getMonomial | ( | const unsigned int | npos | ) | const |
Get the Monomial at given position.
Return the Monomial corresponding to the non-zero coefficient at the given position in the DA object (monomials use one based indexing!).
[in] | npos | position within the DA object. The ordering of the Monomials within a DA object is implementation dependent and does not correspond to the order in which Monomials are listed in the ASCII output routines. |
DACE::DACEException |
void DACE::DA::getMonomial | ( | const unsigned int | npos, |
Monomial & | m | ||
) | const |
Extract the Monomial at given position.
Return the Monomial corresponding to the non-zero coefficient at the given position in the DA object (monomials use one based indexing!).
[in] | npos | position within the DA object. The ordering of the Monomials within a DA object is implementation dependent and does not correspond to the order in which Monomials are listed in the ASCII output routines. |
[out] | m | the monomial object in which to store the corresponding monomial. |
DACE::DACEException |
std::vector< Monomial > DACE::DA::getMonomials | ( | ) | const |
Get std::vector of all non-zero Monomials.
Return a vector of all Monomials in the DA object (differently from getMonomial() where only a single Monomial, corresponding to a specified position in the DA object, is returned).
DACE::DACEException |
|
static |
Get truncation order.
Return the truncation order currently set for the computations. All terms larger than the truncation order are discarded in subsequent operations.
DACE::DACEException |
AlgebraicVector< DA > DACE::DA::gradient | ( | ) | const |
Gradient vector with respect to all independent DA variables.
Compute the gradient of the DA object.
DACE::DACEException |
hypotenuse
Compute the hypotenuse (sqrt(a*a+b*b)) of a DA object and the given DA argument. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::icrt | ( | ) | const |
inverse cubic root
Compute the inverse cubic root of a DA object. The result is copied in a new DA object.
DACE::DACEException |
|
static |
|
static |
DACE initialization.
Initialize the DACE control arrays and set the maximum order and the maximum number of variables.
[in] | ord | order of the Taylor polynomials; |
[in] | nvar | number of variables considered. |
DA DACE::DA::integ | ( | const unsigned int | i | ) | const |
Integral with respect to given variable.
Compute the integral of a DA object with respect to variable i. The result is copied in a new DA object.
[in] | i | variable with respect to which the integral is calculated. |
DACE::DACEException |
DA DACE::DA::integ | ( | const std::vector< unsigned int > | ind | ) | const |
Integral with respect to given variables.
Compute the integral of a DA object with respect to variables i. The result is copied in a new DA object.
[in] | ind | vector containing the number of integrals to take for each independent variable. If ind has fewer entries than there are independent variables, the missing entries are assumed to be zero. If ind has more than there are independent variables, extra values are ignored. |
DACE::DACEException |
int DACE::DA::isinf | ( | ) | const |
Check if a DA object has any INF coefficients.
DACE::DACEException |
|
static |
int DACE::DA::isnan | ( | ) | const |
Check if a DA object has any NAN coefficients.
DACE::DACEException |
DA DACE::DA::isrt | ( | ) | const |
inverse square root
Compute the inverse square root of a DA object. The result is copied in a new DA object.
DACE::DACEException |
AlgebraicVector< double > DACE::DA::linear | ( | ) | const |
Get linear part of a DA.
Return the linear part of a DA object.
DACE::DACEException |
DA DACE::DA::log | ( | ) | const |
natural logarithm
Compute the natural logarithm of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::log10 | ( | ) | const |
10 based logarithm
Compute the 10 based logarithm of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::log2 | ( | ) | const |
2 based logarithm
Compute the 2 based logarithm of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::logb | ( | const double | b = 10.0 | ) | const |
logarithm with respect to given base
Compute the logarithm of a DA object with respect to a given base. The result is copied in a new DA object.
[in] | b | base with respect to which the logarithm is computed (base 10 set as default base). |
DACE::DACEException |
DA DACE::DA::LogGammaFunction | ( | ) | const |
Logartihmic Gamma function.
Compute the Logarithmic Gamma function (i.e. the natural logarithm of Gamma) of a DA object. The result is copied in a new DA object.
DACE::DACEException |
|
static |
Dump memory status to screen (for debugging purposes only)
DA DACE::DA::minv | ( | ) | const |
multiplicative inverse
Compute the multiplicative inverse of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::mod | ( | const double | p | ) | const |
Modulo of the constant part.
Compute the floating-point remainder of c/p (c modulo p), where c is the constant part of the current DA object. The result is copied in a new DA object.
[in] | p | costant with respect to which the modulo function is computed. |
DACE::DACEException |
Multiply the DA with the argument da monomial by monomial (i.e. coefficient-wise)
Multiply the DA vector with another DA vector da monomial by monomial. This is the equivalent of coefficient-wise multiplication (like in DA addition).
[in] | da | DA vector to multiply with coefficient-wise |
DACE::DACEException |
double DACE::DA::norm | ( | const unsigned int | type = 0 | ) | const |
Different types of norms over all coefficients.
Compute different types of norms for a DA object.
[in] | type | type of norm to be computed. Possible norms are: 0: Max norm (Default) 1: Sum norm >1: Vector norm of given type |
DACE::DACEException |
DA & DACE::DA::operator*= | ( | const double | c | ) |
Multiply this DA by another double.
Compute the product between the current DA object and a given constant. The result is directly copied into the current DA object.
[in] | c | constant value to be multiplied. |
DACE::DACEException |
DA & DACE::DA::operator+= | ( | const double | c | ) |
Add another double to this DA.
Compute the sum between the current DA object and a given constant. The result is directly copied into the current DA object.
[in] | c | constant value to be added. |
DACE::DACEException |
DA DACE::DA::operator- | ( | ) | const |
Negate this DA.
Compute the additive inverse of the given DA object. The result is copied in a new DA vector.
DACE::DACEException |
DA & DACE::DA::operator-= | ( | const double | c | ) |
Subtract another double from this DA.
Compute the difference between the current DA object and a given constant. The result is directly copied into the current DA object.
[in] | c | constant value to be subtracted. |
DACE::DACEException |
Compute the division between the current DA object and the given one. The result is directly copied into the current DA object.
DACE::DACEException |
DA & DACE::DA::operator/= | ( | const double | c | ) |
Divide this DA by another double.
Compute the division between the current DA object and a given constant. The result is directly copied into the current DA object.
[in] | c | constant value through which the current DA is divided by. |
DACE::DACEException |
Move assignment from existing DA.
Move all resources from a DA object to us without copying.
[in] | da | object to be moved. |
DACE::DACEException |
DA & DACE::DA::operator= | ( | const double | c | ) |
Assignment from existing double.
Copy a constant polynomial of value c into the DA.
[in] | c | Constant value to be copied. |
DACE::DACEException |
std::vector< double > DACE::DA::orderNorm | ( | const unsigned int | var = 0 , |
const unsigned int | type = 0 |
||
) | const |
Different types of norms over coefficients of each order separately.
Extract different types of order sorted norms from a DA object.
[in] | var | order 0: Terms are sorted by their order (Default) >0: Terms are sorted by the exponent of variable var |
[in] | type | type of norm to be computed. Possible norms are: 0: Max norm (Default) 1: Sum norm >1: Vector norm of given type |
DACE::DACEException |
DA DACE::DA::plug | ( | const unsigned int | var, |
const double | val = 0.0 |
||
) | const |
Partial evaluation to replace given independent DA variable by value val.
Partial evaluation of a DA object. In the DA object, variable var is replaced by the value val. The resulting DA object is returned.
[in] | var | variable number to be replaced |
[in] | val | value by which to replace the variable |
DACE::DACEException |
|
static |
Restore truncation order to value before last DA::pushTO()
Restore the previous truncation order from the truncation order stack. All terms larger than the truncation order are discarded in subsequent operations.
DACE::DACEException |
DA DACE::DA::pow | ( | const int | p | ) | const |
Exponentiation to given (integer) power.
Elevate a DA object to a given integer power. The result is copied in a new DA object.
[in] | p | power to which the DA object is raised. |
DACE::DACEException |
DA DACE::DA::pow | ( | const double | p | ) | const |
Exponentiation to given double power.
Elevate a DA object to a given real power. The constant part must be positive. The result is copied in a new DA object.
[in] | p | power to which the DA object is raised. |
DACE::DACEException |
DA DACE::DA::PsiFunction | ( | const unsigned int | n | ) | const |
Psi function.
Compute the n-th order Psi function, i.e. the (n+1)st derivative of the Logarithmic Gamma function, of a DA object. The result is copied in a new DA object.
DACE::DACEException |
|
static |
Set truncation order and save previous value.
Set a new truncation order, saving the previous one on the truncation order stack. All terms larger than the truncation order are discarded in subsequent operations.
[in] | ot | new truncation order. |
DACE::DACEException |
|
static |
Create new DA filled with random coefficients.
Create a DA object and fill with random entries.
[in] | cm | filling factor (for cm < 0, the DA object is filled with random numbers; for cm > 0, the DA object is filled with weighted decaying numbers). |
DACE::DACEException |
|
static |
Create new DA from blob read from binary file.
Read a binary representation of a DA from is.
DACE::DACEException |
DA DACE::DA::replaceVariable | ( | const unsigned int | from = 0 , |
const unsigned int | to = 0 , |
||
const double | val = 1.0 |
||
) | const |
Replace variable number from by val times variable number to.
Partial evaluation of a DA object. In the DA object, variable from is replaced by the value val times variable to. The resulting DA object is returned.
[in] | from | variable number to be replaced |
[in] | to | variable number to be inserted instead |
[in] | val | value by which to scale the inserted variable |
DACE::DACEException |
DA DACE::DA::root | ( | const int | p = 2 | ) | const |
p-th root
Compute the p-th root of a DA object. The result is copied in a new DA object.
[in] | p | root to be computed. |
DACE::DACEException |
DA DACE::DA::round | ( | ) | const |
Round the constant part to an integer.
Round the constant part of a DA object to an integer. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::scaleVariable | ( | const unsigned int | var = 0 , |
const double | val = 1.0 |
||
) | const |
Scale variable number var by val.
Scaling of an independent variable. In the DA object, variable var is replaced by the value val times var. The resulting DA object is returned.
[in] | var | variable number to be scaled |
[in] | val | value by which to scale the variable |
DACE::DACEException |
void DACE::DA::setCoefficient | ( | const std::vector< unsigned int > & | jj, |
const double | coeff | ||
) |
Set specific coefficient.
Set a specific coefficient into a DA object.
[in] | jj | vector of the exponents of the coefficient to be set. |
[in] | coeff | value to be set as coefficient. |
DACE::DACEException |
|
static |
Set truncation epsilon.
Set the cutoff value eps to a new value and return the previous value.
[in] | eps | new cutoff value. |
DACE::DACEException |
|
static |
Set truncation order.
Set the truncation order ot to a new value and return the previous value. All terms larger than the truncation order are discarded in subsequent operations.
[in] | ot | new truncation order. |
DACE::DACEException |
DA DACE::DA::sin | ( | ) | const |
sine
Compute the sine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::sinh | ( | ) | const |
hyperbolic sine
Compute the hyperbolic sine of a DA object. The result is copied in a new DA object.
DACE::DACEException |
unsigned int DACE::DA::size | ( | ) | const |
Number of non-zero coefficients.
Return the number of non-zero coefficients of a DA object.
DACE::DACEException |
DA DACE::DA::sqr | ( | ) | const |
square
Compute the square of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::sqrt | ( | ) | const |
square root
Compute the square root of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::tan | ( | ) | const |
tangent
Compute the tangent of a DA object. The result is copied in a new DA object.
DACE::DACEException |
DA DACE::DA::tanh | ( | ) | const |
hyperbolic tangent
Compute the hyperbolic tangent of a DA object. The result is copied in a new DA object.
DACE::DACEException |
std::string DACE::DA::toString | ( | ) | const |
Convert to string representation.
Convert DA object to string.
DACE::DACEException |
DA DACE::DA::translateVariable | ( | const unsigned int | var = 0 , |
const double | a = 1.0 , |
||
const double | c = 0.0 |
||
) | const |
Translate variable number var to a*var + c.
Affine translation of an independent variable. In the DA object, variable var is replaced by a*var + c. The resulting DA object is returned.
[in] | var | variable number to be translated |
[in] | a | value by which to scale the variable |
[in] | c | value by which to shift the variable |
DACE::DACEException |
DA DACE::DA::trim | ( | const unsigned int | min, |
const unsigned int | max = DA::getMaxOrder() |
||
) | const |
Trim the coefficients only include orders between min and max, inclusively.
Returns a DA object with all monomials of order less than min and greater than max removed. The result is copied in a new DA object.
[in] | min | The minimum order to keep in the DA object. |
[in] | max | The maximum order to keep in the DA object. |
DACE::DACEException |
DA DACE::DA::trunc | ( | ) | const |
Truncate the constant part to an integer.
Truncate the constant part of a DA object to an integer. The result is copied in a new DA object.
DACE::DACEException |
|
static |
void DACE::DA::write | ( | std::ostream & | os | ) | const |
Write binary representation of DA to stream.
Write a binary representation of the DA as a blob to os.
DACE::DACEException |
|
friend |
Multiplication between a DA and a constant.
Compute the multiplication between a DA object and a given constant. The result is copied in a new DA object.
[in] | da | DA object. |
[in] | c | given constant. |
DACE::DACEException |
Multiplication between a constant and a DA.
Compute the multiplication between a given constant and a DA object. The result is copied in a new DA object.
[in] | c | given constant. |
[in] | da | DA object. |
DACE::DACEException |
Output to C++ stream in text form.
Overload of std::operator<< in iostream.
[in] | out | standard output stream. |
[in] | da | DA object to be printed in the stream |
DACE::DACEException |
Input from C++ stream in text form.
Overload of std::operator>> in iostream. Reads both string and binary DA representations from a file.
[in] | in | standard input stream. |
[in] | da | DA object to be created from the stream |
DACE::DACEException |
|
friend |