Table of Contents


Node:Top, Next:, Previous:(dir), Up:(dir)

GiNaC

This is a tutorial that documents GiNaC 0.8.0, an open framework for symbolic computation within the C++ programming language.


Node:Introduction, Next:, Previous:Top, Up:Top

Introduction

The motivation behind GiNaC derives from the observation that most present day computer algebra systems (CAS) are linguistically and semantically impoverished. Although they are quite powerful tools for learning math and solving particular problems they lack modern linguistical structures that allow for the creation of large-scale projects. GiNaC is an attempt to overcome this situation by extending a well established and standardized computer language (C++) by some fundamental symbolic capabilities, thus allowing for integrated systems that embed symbolic manipulations together with more established areas of computer science (like computation-intense numeric applications, graphical interfaces, etc.) under one roof.

The particular problem that led to the writing of the GiNaC framework is still a very active field of research, namely the calculation of higher order corrections to elementary particle interactions. There, theoretical physicists are interested in matching present day theories against experiments taking place at particle accelerators. The computations involved are so complex they call for a combined symbolical and numerical approach. This turned out to be quite difficult to accomplish with the present day CAS we have worked with so far and so we tried to fill the gap by writing GiNaC. But of course its applications are in no way restricted to theoretical physics.

This tutorial is intended for the novice user who is new to GiNaC but already has some background in C++ programming. However, since a hand-made documentation like this one is difficult to keep in sync with the development, the actual documentation is inside the sources in the form of comments. That documentation may be parsed by one of the many Javadoc-like documentation systems. If you fail at generating it you may access it from the GiNaC home page. It is an invaluable resource not only for the advanced user who wishes to extend the system (or chase bugs) but for everybody who wants to comprehend the inner workings of GiNaC. This little tutorial on the other hand only covers the basic things that are unlikely to change in the near future.

License

The GiNaC framework for symbolic computation within the C++ programming language is Copyright © 1999-2001 Johannes Gutenberg University Mainz, Germany.

This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program; see the file COPYING. If not, write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.


Node:A Tour of GiNaC, Next:, Previous:Introduction, Up:Top

A Tour of GiNaC

This quick tour of GiNaC wants to arise your interest in the subsequent chapters by showing off a bit. Please excuse us if it leaves many open questions.


Node:How to use it from within C++, Next:, Previous:A Tour of GiNaC, Up:A Tour of GiNaC

How to use it from within C++

The GiNaC open framework for symbolic computation within the C++ programming language does not try to define a language of its own as conventional CAS do. Instead, it extends the capabilities of C++ by symbolic manipulations. Here is how to generate and print a simple (and rather pointless) bivariate polynomial with some large coefficients:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y");
    ex poly;

    for (int i=0; i<3; ++i)
        poly += factorial(i+16)*pow(x,i)*pow(y,2-i);

    cout << poly << endl;
    return 0;
}

Assuming the file is called hello.cc, on our system we can compile and run it like this:

$ c++ hello.cc -o hello -lcln -lginac
$ ./hello
355687428096000*x*y+20922789888000*y^2+6402373705728000*x^2

(See Package Tools, for tools that help you when creating a software package that uses GiNaC.)

Next, there is a more meaningful C++ program that calls a function which generates Hermite polynomials in a specified free variable.

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

ex HermitePoly(const symbol & x, int n)
{
    ex HKer=exp(-pow(x, 2));
    // uses the identity H_n(x) == (-1)^n exp(x^2) (d/dx)^n exp(-x^2)
    return normal(pow(-1, n) * diff(HKer, x, n) / HKer);
}

int main()
{
    symbol z("z");

    for (int i=0; i<6; ++i)
        cout << "H_" << i << "(z) == " << HermitePoly(z,i) << endl;

    return 0;
}

When run, this will type out

H_0(z) == 1
H_1(z) == 2*z
H_2(z) == 4*z^2-2
H_3(z) == -12*z+8*z^3
H_4(z) == -48*z^2+16*z^4+12
H_5(z) == 120*z-160*z^3+32*z^5

This method of generating the coefficients is of course far from optimal for production purposes.

In order to show some more examples of what GiNaC can do we will now use the ginsh, a simple GiNaC interactive shell that provides a convenient window into GiNaC's capabilities.


Node:What it can do for you, Next:, Previous:How to use it from within C++, Up:A Tour of GiNaC

What it can do for you

After invoking ginsh one can test and experiment with GiNaC's features much like in other Computer Algebra Systems except that it does not provide programming constructs like loops or conditionals. For a concise description of the ginsh syntax we refer to its accompanied man page. Suffice to say that assignments and comparisons in ginsh are written as they are in C, i.e. = assigns and == compares.

It can manipulate arbitrary precision integers in a very fast way. Rational numbers are automatically converted to fractions of coprime integers:

> x=3^150;
369988485035126972924700782451696644186473100389722973815184405301748249
> y=3^149;
123329495011708990974900260817232214728824366796574324605061468433916083
> x/y;
3
> y/x;
1/3

Exact numbers are always retained as exact numbers and only evaluated as floating point numbers if requested. For instance, with numeric radicals is dealt pretty much as with symbols. Products of sums of them can be expanded:

> expand((1+a^(1/5)-a^(2/5))^3);
1+3*a+3*a^(1/5)-5*a^(3/5)-a^(6/5)
> expand((1+3^(1/5)-3^(2/5))^3);
10-5*3^(3/5)
> evalf((1+3^(1/5)-3^(2/5))^3);
0.33408977534118624228

The function evalf that was used above converts any number in GiNaC's expressions into floating point numbers. This can be done to arbitrary predefined accuracy:

> evalf(1/7);
0.14285714285714285714
> Digits=150;
150
> evalf(1/7);
0.1428571428571428571428571428571428571428571428571428571428571428571428
5714285714285714285714285714285714285

Exact numbers other than rationals that can be manipulated in GiNaC include predefined constants like Archimedes' Pi. They can both be used in symbolic manipulations (as an exact number) as well as in numeric expressions (as an inexact number):

> a=Pi^2+x;
x+Pi^2
> evalf(a);
9.869604401089358619+x
> x=2;
2
> evalf(a);
11.869604401089358619

Built-in functions evaluate immediately to exact numbers if this is possible. Conversions that can be safely performed are done immediately; conversions that are not generally valid are not done:

> cos(42*Pi);
1
> cos(acos(x));
x
> acos(cos(x));
acos(cos(x))

(Note that converting the last input to x would allow one to conclude that 42*Pi is equal to 0.)

Linear equation systems can be solved along with basic linear algebra manipulations over symbolic expressions. In C++ GiNaC offers a matrix class for this purpose but we can see what it can do using ginsh's notation of double brackets to type them in:

> lsolve(a+x*y==z,x);
y^(-1)*(z-a);
> lsolve([3*x+5*y == 7, -2*x+10*y == -5], [x, y]);
[x==19/8,y==-1/40]
> M = [[ [[1, 3]], [[-3, 2]] ]];
[[ [[1,3]], [[-3,2]] ]]
> determinant(M);
11
> charpoly(M,lambda);
lambda^2-3*lambda+11

Multivariate polynomials and rational functions may be expanded, collected and normalized (i.e. converted to a ratio of two coprime polynomials):

> a = x^4 + 2*x^2*y^2 + 4*x^3*y + 12*x*y^3 - 3*y^4;
-3*y^4+x^4+12*x*y^3+2*x^2*y^2+4*x^3*y
> b = x^2 + 4*x*y - y^2;
-y^2+x^2+4*x*y
> expand(a*b);
3*y^6+x^6-24*x*y^5+43*x^2*y^4+16*x^3*y^3+17*x^4*y^2+8*x^5*y
> collect(a*b,x);
3*y^6+48*x*y^4+2*x^2*y^2+x^4*(-y^2+x^2+4*x*y)+4*x^3*y*(-y^2+x^2+4*x*y)
> normal(a/b);
3*y^2+x^2

You can differentiate functions and expand them as Taylor or Laurent series in a very natural syntax (the second argument of series is a relation defining the evaluation point, the third specifies the order):

> diff(tan(x),x);
tan(x)^2+1
> series(sin(x),x==0,4);
x-1/6*x^3+Order(x^4)
> series(1/tan(x),x==0,4);
x^(-1)-1/3*x+Order(x^2)
> series(tgamma(x),x==0,3);
x^(-1)-Euler+(1/12*Pi^2+1/2*Euler^2)*x+
(-1/3*zeta(3)-1/12*Pi^2*Euler-1/6*Euler^3)*x^2+Order(x^3)
> evalf(");
x^(-1)-0.5772156649015328606+(0.9890559953279725555)*x
-(0.90747907608088628905)*x^2+Order(x^3)
> series(tgamma(2*sin(x)-2),x==Pi/2,6);
-(x-1/2*Pi)^(-2)+(-1/12*Pi^2-1/2*Euler^2-1/240)*(x-1/2*Pi)^2
-Euler-1/12+Order((x-1/2*Pi)^3)

Here we have made use of the ginsh-command " to pop the previously evaluated element from ginsh's internal stack.

If you ever wanted to convert units in C or C++ and found this is cumbersome, here is the solution. Symbolic types can always be used as tags for different types of objects. Converting from wrong units to the metric system is now easy:

> in=.0254*m;
0.0254*m
> lb=.45359237*kg;
0.45359237*kg
> 200*lb/in^2;
140613.91592783185568*kg*m^(-2)


Node:Installation, Next:, Previous:What it can do for you, Up:Top

Installation

GiNaC's installation follows the spirit of most GNU software. It is easily installed on your system by three steps: configuration, build, installation.


Node:Prerequisites, Next:, Previous:Installation, Up:Installation

Prerequisites

In order to install GiNaC on your system, some prerequisites need to be met. First of all, you need to have a C++-compiler adhering to the ANSI-standard ISO/IEC 14882:1998(E). We used GCC for development so if you have a different compiler you are on your own. For the configuration to succeed you need a Posix compliant shell installed in /bin/sh, GNU bash is fine. Perl is needed by the built process as well, since some of the source files are automatically generated by Perl scripts. Last but not least, Bruno Haible's library CLN is extensively used and needs to be installed on your system. Please get it either from ftp://ftp.santafe.edu/pub/gnu/, from GiNaC's FTP site or from Bruno Haible's FTP site (it is covered by GPL) and install it prior to trying to install GiNaC. The configure script checks if it can find it and if it cannot it will refuse to continue.


Node:Configuration, Next:, Previous:Prerequisites, Up:Installation

Configuration

To configure GiNaC means to prepare the source distribution for building. It is done via a shell script called configure that is shipped with the sources and was originally generated by GNU Autoconf. Since a configure script generated by GNU Autoconf never prompts, all customization must be done either via command line parameters or environment variables. It accepts a list of parameters, the complete set of which can be listed by calling it with the --help option. The most important ones will be shortly described in what follows:

In addition, you may specify some environment variables. CXX holds the path and the name of the C++ compiler in case you want to override the default in your path. (The configure script searches your path for c++, g++, gcc, CC, cxx and cc++ in that order.) It may be very useful to define some compiler flags with the CXXFLAGS environment variable, like optimization, debugging information and warning levels. If omitted, it defaults to -g -O2.

The whole process is illustrated in the following two examples. (Substitute setenv VARIABLE value for export VARIABLE=value if the Berkeley C shell is your login shell.)

Here is a simple configuration for a site-wide GiNaC library assuming everything is in default paths:

$ export CXXFLAGS="-Wall -O2"
$ ./configure

And here is a configuration for a private static GiNaC library with several components sitting in custom places (site-wide GCC and private CLN). The compiler is pursuaded to be picky and full assertions and debugging information are switched on:

$ export CXX=/usr/local/gnu/bin/c++
$ export CPPFLAGS="$(CPPFLAGS) -I$(HOME)/include"
$ export CXXFLAGS="$(CXXFLAGS) -DDO_GINAC_ASSERT -ggdb -Wall -ansi -pedantic"
$ export LDFLAGS="$(LDFLAGS) -L$(HOME)/lib"
$ ./configure --disable-shared --prefix=$(HOME)


Node:Building GiNaC, Next:, Previous:Configuration, Up:Installation

Building GiNaC

After proper configuration you should just build the whole library by typing

$ make
at the command prompt and go for a cup of coffee. The exact time it takes to compile GiNaC depends not only on the speed of your machines but also on other parameters, for instance what value for CXXFLAGS you entered. Optimization may be very time-consuming.

Just to make sure GiNaC works properly you may run a collection of regression tests by typing

$ make check

This will compile some sample programs, run them and check the output for correctness. The regression tests fall in three categories. First, the so called exams are performed, simple tests where some predefined input is evaluated (like a pupils' exam). Second, the checks test the coherence of results among each other with possible random input. Third, some timings are performed, which benchmark some predefined problems with different sizes and display the CPU time used in seconds. Each individual test should return a message passed. This is mostly intended to be a QA-check if something was broken during development, not a sanity check of your system. Some of the tests in sections checks and timings may require insane amounts of memory and CPU time. Feel free to kill them if your machine catches fire. Another quite important intent is to allow people to fiddle around with optimization.

Generally, the top-level Makefile runs recursively to the subdirectories. It is therfore safe to go into any subdirectory (doc/, ginsh/, ...) and simply type make target there in case something went wrong.


Node:Installing GiNaC, Next:, Previous:Building GiNaC, Up:Installation

Installing GiNaC

To install GiNaC on your system, simply type

$ make install

As described in the section about configuration the files will be installed in the following directories (the directories will be created if they don't already exist):

For the sake of completeness we will list some other useful make targets: make clean deletes all files generated by make, i.e. all the object files. In addition make distclean removes all files generated by the configuration and make maintainer-clean goes one step further and deletes files that may require special tools to rebuild (like the libtool for instance). Finally make uninstall removes the installed library, header files and documentation1.


Node:Basic Concepts, Next:, Previous:Installing GiNaC, Up:Top

Basic Concepts

This chapter will describe the different fundamental objects that can be handled by GiNaC. But before doing so, it is worthwhile introducing you to the more commonly used class of expressions, representing a flexible meta-class for storing all mathematical objects.


Node:Expressions, Next:, Previous:Basic Concepts, Up:Basic Concepts

Expressions

The most common class of objects a user deals with is the expression ex, representing a mathematical object like a variable, number, function, sum, product, etc... Expressions may be put together to form new expressions, passed as arguments to functions, and so on. Here is a little collection of valid expressions:

ex MyEx1 = 5;                       // simple number
ex MyEx2 = x + 2*y;                 // polynomial in x and y
ex MyEx3 = (x + 1)/(x - 1);         // rational expression
ex MyEx4 = sin(x + 2*y) + 3*z + 41; // containing a function
ex MyEx5 = MyEx4 + 1;               // similar to above

Expressions are handles to other more fundamental objects, that often contain other expressions thus creating a tree of expressions (See Internal Structures, for particular examples). Most methods on ex therefore run top-down through such an expression tree. For example, the method has() scans recursively for occurrences of something inside an expression. Thus, if you have declared MyEx4 as in the example above MyEx4.has(y) will find y inside the argument of sin and hence return true.

The next sections will outline the general picture of GiNaC's class hierarchy and describe the classes of objects that are handled by ex.


Node:The Class Hierarchy, Next:, Previous:Expressions, Up:Basic Concepts

The Class Hierarchy

GiNaC's class hierarchy consists of several classes representing mathematical objects, all of which (except for ex and some helpers) are internally derived from one abstract base class called basic. You do not have to deal with objects of class basic, instead you'll be dealing with symbols, numbers, containers of expressions and so on.

To get an idea about what kinds of symbolic composits may be built we have a look at the most important classes in the class hierarchy and some of the relations among the classes: classhierarchy.png

The abstract classes shown here (the ones without drop-shadow) are of no interest for the user. They are used internally in order to avoid code duplication if two or more classes derived from them share certain features. An example is expairseq, a container for a sequence of pairs each consisting of one expression and a number (numeric). What is visible to the user are the derived classes add and mul, representing sums and products. See Internal Structures, where these two classes are described in more detail. The following table shortly summarizes what kinds of mathematical objects are stored in the different classes:

symbol Algebraic symbols a, x, y...
constant Constants like Pi
numeric All kinds of numbers, 42, 7/3*I, 3.14159...
add Sums like x+y or a-(2*b)+3
mul Products like x*y or 2*a^2*(x+y+z)/b
power Exponentials such as x^2, a^b, sqrt(2) ...
pseries Power Series, e.g. x-1/6*x^3+1/120*x^5+O(x^7)
function A symbolic function like sin(2*x)
lst Lists of expressions [x, 2*y, 3+z]
matrix nxm matrices of expressions
relational A relation like the identity x==y
indexed Indexed object like A_ij
tensor Special tensor like the delta and metric tensors
idx Index of an indexed object
varidx Index with variance


Node:Symbols, Next:, Previous:The Class Hierarchy, Up:Basic Concepts

Symbols

Symbols are for symbolic manipulation what atoms are for chemistry. You can declare objects of class symbol as any other object simply by saying symbol x,y;. There is, however, a catch in here having to do with the fact that C++ is a compiled language. The information about the symbol's name is thrown away by the compiler but at a later stage you may want to print expressions holding your symbols. In order to avoid confusion GiNaC's symbols are able to know their own name. This is accomplished by declaring its name for output at construction time in the fashion symbol x("x");. If you declare a symbol using the default constructor (i.e. without string argument) the system will deal out a unique name. That name may not be suitable for printing but for internal routines when no output is desired it is often enough. We'll come across examples of such symbols later in this tutorial.

This implies that the strings passed to symbols at construction time may not be used for comparing two of them. It is perfectly legitimate to write symbol x("x"),y("x"); but it is likely to lead into trouble. Here, x and y are different symbols and statements like x-y will not be simplified to zero although the output x-x looks funny. Such output may also occur when there are two different symbols in two scopes, for instance when you call a function that declares a symbol with a name already existent in a symbol in the calling function. Again, comparing them (using operator== for instance) will always reveal their difference. Watch out, please.

Although symbols can be assigned expressions for internal reasons, you should not do it (and we are not going to tell you how it is done). If you want to replace a symbol with something else in an expression, you can use the expression's .subs() method (See Substituting Symbols, for more information).


Node:Numbers, Next:, Previous:Symbols, Up:Basic Concepts

Numbers

For storing numerical things, GiNaC uses Bruno Haible's library CLN. The classes therein serve as foundation classes for GiNaC. CLN stands for Class Library for Numbers or alternatively for Common Lisp Numbers. In order to find out more about CLN's internals the reader is refered to the documentation of that library. see , for more information. Suffice to say that it is by itself build on top of another library, the GNU Multiple Precision library GMP, which is an extremely fast library for arbitrary long integers and rationals as well as arbitrary precision floating point numbers. It is very commonly used by several popular cryptographic applications. CLN extends GMP by several useful things: First, it introduces the complex number field over either reals (i.e. floating point numbers with arbitrary precision) or rationals. Second, it automatically converts rationals to integers if the denominator is unity and complex numbers to real numbers if the imaginary part vanishes and also correctly treats algebraic functions. Third it provides good implementations of state-of-the-art algorithms for all trigonometric and hyperbolic functions as well as for calculation of some useful constants.

The user can construct an object of class numeric in several ways. The following example shows the four most important constructors. It uses construction from C-integer, construction of fractions from two integers, construction from C-float and construction from a string:

#include <ginac/ginac.h>
using namespace GiNaC;

int main()
{
    numeric two(2);                       // exact integer 2
    numeric r(2,3);                       // exact fraction 2/3
    numeric e(2.71828);                   // floating point number
    numeric p("3.1415926535897932385");   // floating point number
    // Trott's constant in scientific notation:
    numeric trott("1.0841015122311136151E-2");

    std::cout << two*p << std::endl;  // floating point 6.283...
}

Note that all those constructors are explicit which means you are not allowed to write numeric two=2;. This is because the basic objects to be handled by GiNaC are the expressions ex and we want to keep things simple and wish objects like pow(x,2) to be handled the same way as pow(x,a), which means that we need to allow a general ex as base and exponent. Therefore there is an implicit constructor from C-integers directly to expressions handling numerics at work in most of our examples. This design really becomes convenient when one declares own functions having more than one parameter but it forbids using implicit constructors because that would lead to compile-time ambiguities.

It may be tempting to construct numbers writing numeric r(3/2). This would, however, call C's built-in operator / for integers first and result in a numeric holding a plain integer 1. Never use the operator / on integers unless you know exactly what you are doing! Use the constructor from two integers instead, as shown in the example above. Writing numeric(1)/2 may look funny but works also.

We have seen now the distinction between exact numbers and floating point numbers. Clearly, the user should never have to worry about dynamically created exact numbers, since their `exactness' always determines how they ought to be handled, i.e. how `long' they are. The situation is different for floating point numbers. Their accuracy is controlled by one global variable, called Digits. (For those readers who know about Maple: it behaves very much like Maple's Digits). All objects of class numeric that are constructed from then on will be stored with a precision matching that number of decimal digits:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

void foo()
{
    numeric three(3.0), one(1.0);
    numeric x = one/three;

    cout << "in " << Digits << " digits:" << endl;
    cout << x << endl;
    cout << Pi.evalf() << endl;
}

int main()
{
    foo();
    Digits = 60;
    foo();
    return 0;
}

The above example prints the following output to screen:

in 17 digits:
0.333333333333333333
3.14159265358979324
in 60 digits:
0.333333333333333333333333333333333333333333333333333333333333333333
3.14159265358979323846264338327950288419716939937510582097494459231

It should be clear that objects of class numeric should be used for constructing numbers or for doing arithmetic with them. The objects one deals with most of the time are the polymorphic expressions ex.

Tests on numbers

Once you have declared some numbers, assigned them to expressions and done some arithmetic with them it is frequently desired to retrieve some kind of information from them like asking whether that number is integer, rational, real or complex. For those cases GiNaC provides several useful methods. (Internally, they fall back to invocations of certain CLN functions.)

As an example, let's construct some rational number, multiply it with some multiple of its denominator and test what comes out:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

// some very important constants:
const numeric twentyone(21);
const numeric ten(10);
const numeric five(5);

int main()
{
    numeric answer = twentyone;

    answer /= five;
    cout << answer.is_integer() << endl;  // false, it's 21/5
    answer *= ten;
    cout << answer.is_integer() << endl;  // true, it's 42 now!
}

Note that the variable answer is constructed here as an integer by numeric's copy constructor but in an intermediate step it holds a rational number represented as integer numerator and integer denominator. When multiplied by 10, the denominator becomes unity and the result is automatically converted to a pure integer again. Internally, the underlying CLN is responsible for this behaviour and we refer the reader to CLN's documentation. Suffice to say that the same behaviour applies to complex numbers as well as return values of certain functions. Complex numbers are automatically converted to real numbers if the imaginary part becomes zero. The full set of tests that can be applied is listed in the following table.

Method Returns true if the object is...
.is_zero() ...equal to zero
.is_positive() ...not complex and greater than 0
.is_integer() ...a (non-complex) integer
.is_pos_integer() ...an integer and greater than 0
.is_nonneg_integer() ...an integer and greater equal 0
.is_even() ...an even integer
.is_odd() ...an odd integer
.is_prime() ...a prime integer (probabilistic primality test)
.is_rational() ...an exact rational number (integers are rational, too)
.is_real() ...a real integer, rational or float (i.e. is not complex)
.is_cinteger() ...a (complex) integer (such as 2-3*I)
.is_crational() ...an exact (complex) rational number (such as 2/3+7/2*I)


Node:Constants, Next:, Previous:Numbers, Up:Basic Concepts

Constants

Constants behave pretty much like symbols except that they return some specific number when the method .evalf() is called.

The predefined known constants are:

Name Common Name Numerical Value (to 35 digits)
Pi Archimedes' constant 3.14159265358979323846264338327950288
Catalan Catalan's constant 0.91596559417721901505460351493238411
Euler Euler's (or Euler-Mascheroni) constant 0.57721566490153286060651209008240243


Node:Fundamental containers, Next:, Previous:Constants, Up:Basic Concepts

Fundamental containers: the power, add and mul classes

Simple polynomial expressions are written down in GiNaC pretty much like in other CAS or like expressions involving numerical variables in C. The necessary operators +, -, * and / have been overloaded to achieve this goal. When you run the following code snippet, the constructor for an object of type mul is automatically called to hold the product of a and b and then the constructor for an object of type add is called to hold the sum of that mul object and the number one:

    ...
    symbol a("a"), b("b");
    ex MyTerm = 1+a*b;
    ...

For exponentiation, you have already seen the somewhat clumsy (though C-ish) statement pow(x,2); to represent x squared. This direct construction is necessary since we cannot safely overload the constructor ^ in C++ to construct a power object. If we did, it would have several counterintuitive and undesired effects:

All effects are contrary to mathematical notation and differ from the way most other CAS handle exponentiation, therefore overloading ^ is ruled out for GiNaC's C++ part. The situation is different in ginsh, there the exponentiation-^ exists. (Also note that the other frequently used exponentiation operator ** does not exist at all in C++).

To be somewhat more precise, objects of the three classes described here, are all containers for other expressions. An object of class power is best viewed as a container with two slots, one for the basis, one for the exponent. All valid GiNaC expressions can be inserted. However, basic transformations like simplifying pow(pow(x,2),3) to x^6 automatically are only performed when this is mathematically possible. If we replace the outer exponent three in the example by some symbols a, the simplification is not safe and will not be performed, since a might be 1/2 and x negative.

Objects of type add and mul are containers with an arbitrary number of slots for expressions to be inserted. Again, simple and safe simplifications are carried out like transforming 3*x+4-x to 2*x+4.

The general rule is that when you construct such objects, GiNaC automatically creates them in canonical form, which might differ from the form you typed in your program. This allows for rapid comparison of expressions, since after all a-a is simply zero. Note, that the canonical form is not necessarily lexicographical ordering or in any way easily guessable. It is only guaranteed that constructing the same expression twice, either implicitly or explicitly, results in the same canonical form.


Node:Lists, Next:, Previous:Fundamental containers, Up:Basic Concepts

Lists of expressions

The GiNaC class lst serves for holding a list of arbitrary expressions. These are sometimes used to supply a variable number of arguments of the same type to GiNaC methods such as subs() and to_rational(), so you should have a basic understanding about them.

Lists of up to 15 expressions can be directly constructed from single expressions:

{
    symbol x("x"), y("y");
    lst l(x, 2, y, x+y);
    // now, l is a list holding the expressions 'x', '2', 'y', and 'x+y'
    // ...

Use the nops() method to determine the size (number of expressions) of a list and the op() method to access individual elements:

    // ...
    cout << l.nops() << endl;                   // prints '4'
    cout << l.op(2) << " " << l.op(0) << endl;  // prints 'y x'
    // ...

Finally you can append or prepend an expression to a list with the append() and prepend() methods:

    // ...
    l.append(4*x);   // l is now [x, 2, y, x+y, 4*x]
    l.prepend(0);    // l is now [0, x, 2, y, x+y, 4*x]
}


Node:Mathematical functions, Next:, Previous:Lists, Up:Basic Concepts

Mathematical functions

There are quite a number of useful functions hard-wired into GiNaC. For instance, all trigonometric and hyperbolic functions are implemented (See Built-in Functions, for a complete list).

These functions are all objects of class function. They accept one or more expressions as arguments and return one expression. If the arguments are not numerical, the evaluation of the function may be halted, as it does in the next example, showing how a function returns itself twice and finally an expression that may be really useful:

    ...
    symbol x("x"), y("y");
    ex foo = x+y/2;
    cout << tgamma(foo) << endl;
     // -> tgamma(x+(1/2)*y)
    ex bar = foo.subs(y==1);
    cout << tgamma(bar) << endl;
     // -> tgamma(x+1/2)
    ex foobar = bar.subs(x==7);
    cout << tgamma(foobar) << endl;
     // -> (135135/128)*Pi^(1/2)
    ...

Besides evaluation most of these functions allow differentiation, series expansion and so on. Read the next chapter in order to learn more about this.


Node:Relations, Next:, Previous:Mathematical functions, Up:Basic Concepts

Relations

Sometimes, a relation holding between two expressions must be stored somehow. The class relational is a convenient container for such purposes. A relation is by definition a container for two ex and a relation between them that signals equality, inequality and so on. They are created by simply using the C++ operators ==, !=, <, <=, > and >= between two expressions.

See Mathematical functions, for examples where various applications of the .subs() method show how objects of class relational are used as arguments. There they provide an intuitive syntax for substitutions. They are also used as arguments to the ex::series method, where the left hand side of the relation specifies the variable to expand in and the right hand side the expansion point. They can also be used for creating systems of equations that are to be solved for unknown variables. But the most common usage of objects of this class is rather inconspicuous in statements of the form if (expand(pow(a+b,2))==a*a+2*a*b+b*b) {...}. Here, an implicit conversion from relational to bool takes place. Note, however, that == here does not perform any simplifications, hence expand() must be called explicitly.


Node:Indexed objects, Next:, Previous:Relations, Up:Basic Concepts

Indexed objects

GiNaC allows you to handle expressions containing general indexed objects in arbitrary spaces. It is also able to canonicalize and simplify such expressions and perform symbolic dummy index summations. There are a number of predefined indexed objects provided, like delta and metric tensors.

There are few restrictions placed on indexed objects and their indices and it is easy to construct nonsense expressions, but our intention is to provide a general framework that allows you to implement algorithms with indexed quantities, getting in the way as little as possible.

Indexed quantities and their indices

Indexed expressions in GiNaC are constructed of two special types of objects, index objects and indexed objects.

Note: when printing expressions, covariant indices and indices without variance are denoted .i while contravariant indices are denoted ~i. In the following, we are going to use that notation in the text so instead of A^i_jk we will write A~i.j.k. Index dimensions are not visible in the output.

A simple example shall illustrate the concepts:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
    symbol i_sym("i"), j_sym("j");
    idx i(i_sym, 3), j(j_sym, 3);

    symbol A("A");
    cout << indexed(A, i, j) << endl;
     // -> A.i.j
    ...

The idx constructor takes two arguments, the index value and the index dimension. First we define two index objects, i and j, both with the numeric dimension 3. The value of the index i is the symbol i_sym (which prints as i) and the value of the index j is the symbol j_sym (which prints as j). Next we construct an expression containing one indexed object, A.i.j. It has the symbol A as its base expression and the two indices i and j.

Note the difference between the indices i and j which are of class idx, and the index values which are the sybols i_sym and j_sym. The indices of indexed objects cannot directly be symbols or numbers but must be index objects. For example, the following is not correct and will raise an exception:

symbol i("i"), j("j");
e = indexed(A, i, j); // ERROR: indices must be of type idx

You can have multiple indexed objects in an expression, index values can be numeric, and index dimensions symbolic:

    ...
    symbol B("B"), dim("dim");
    cout << 4 * indexed(A, i)
          + indexed(B, idx(j_sym, 4), idx(2, 3), idx(i_sym, dim)) << endl;
     // -> B.j.2.i+4*A.i
    ...

B has a 4-dimensional symbolic index k, a 3-dimensional numeric index of value 2, and a symbolic index i with the symbolic dimension dim. Note that GiNaC doesn't automatically notify you that the free indices of A and B in the sum don't match (you have to call simplify_indexed() for that, see below).

In fact, base expressions, index values and index dimensions can be arbitrary expressions:

    ...
    cout << indexed(A+B, idx(2*i_sym+1, dim/2)) << endl;
     // -> (B+A).(1+2*i)
    ...

It's also possible to construct nonsense like Pi.sin(x). You will not get an error message from this but you will probably not be able to do anything useful with it.

The methods

ex idx::get_value(void);
ex idx::get_dimension(void);

return the value and dimension of an idx object. If you have an index in an expression, such as returned by calling .op() on an indexed object, you can get a reference to the idx object with the function ex_to_idx() on the expression.

There are also the methods

bool idx::is_numeric(void);
bool idx::is_symbolic(void);
bool idx::is_dim_numeric(void);
bool idx::is_dim_symbolic(void);

for checking whether the value and dimension are numeric or symbolic (non-numeric). Using the info() method of an index (see Information About Expressions) returns information about the index value.

If you need co- and contravariant indices, use the varidx class:

    ...
    symbol mu_sym("mu"), nu_sym("nu");
    varidx mu(mu_sym, 4), nu(nu_sym, 4); // default is contravariant ~mu, ~nu
    varidx mu_co(mu_sym, 4, true);       // covariant index .mu

    cout << indexed(A, mu, nu) << endl;
     // -> A~mu~nu
    cout << indexed(A, mu_co, nu) << endl;
     // -> A.mu~nu
    cout << indexed(A, mu.toggle_variance(), nu) << endl;
     // -> A.mu~nu
    ...

A varidx is an idx with an additional flag that marks it as co- or contravariant. The default is a contravariant (upper) index, but this can be overridden by supplying a third argument to the varidx constructor. The two methods

bool varidx::is_covariant(void);
bool varidx::is_contravariant(void);

allow you to check the variance of a varidx object (use ex_to_varidx() to get the object reference from an expression). There's also the very useful method

ex varidx::toggle_variance(void);

which makes a new index with the same value and dimension but the opposite variance. By using it you only have to define the index once.

Substituting indices

Sometimes you will want to substitute one symbolic index with another symbolic or numeric index, for example when calculating one specific element of a tensor expression. This is done with the .subs() method, as it is done for symbols (see Substituting Symbols).

You have two possibilities here. You can either substitute the whole index by another index or expression:

    ...
    ex e = indexed(A, mu_co);
    cout << e << " becomes " << e.subs(mu_co == nu) << endl;
     // -> A.mu becomes A~nu
    cout << e << " becomes " << e.subs(mu_co == varidx(0, 4)) << endl;
     // -> A.mu becomes A~0
    cout << e << " becomes " << e.subs(mu_co == 0) << endl;
     // -> A.mu becomes A.0
    ...

The third example shows that trying to replace an index with something that is not an index will substitute the index value instead.

Alternatively, you can substitute the symbol of a symbolic index by another expression:

    ...
    ex e = indexed(A, mu_co);
    cout << e << " becomes " << e.subs(mu_sym == nu_sym) << endl;
     // -> A.mu becomes A.nu
    cout << e << " becomes " << e.subs(mu_sym == 0) << endl;
     // -> A.mu becomes A.0
    ...

As you see, with the second method only the value of the index will get substituted. Its other properties, including its dimension, remain unchanged. If you want to change the dimension of an index you have to substitute the whole index by another one with the new dimension.

Finally, substituting the base expression of an indexed object works as expected:

    ...
    ex e = indexed(A, mu_co);
    cout << e << " becomes " << e.subs(A == A+B) << endl;
     // -> A.mu becomes (B+A).mu
    ...

Symmetries

Indexed objects can be declared as being totally symmetric or antisymmetric with respect to their indices. In this case, GiNaC will automatically bring the indices into a canonical order which allows for some immediate simplifications:

    ...
    cout << indexed(A, indexed::symmetric, i, j)
          + indexed(A, indexed::symmetric, j, i) << endl;
     // -> 2*A.j.i
    cout << indexed(B, indexed::antisymmetric, i, j)
          + indexed(B, indexed::antisymmetric, j, j) << endl;
     // -> -B.j.i
    cout << indexed(B, indexed::antisymmetric, i, j)
          + indexed(B, indexed::antisymmetric, j, i) << endl;
     // -> 0
    ...

Dummy indices

GiNaC treats certain symbolic index pairs as dummy indices meaning that a summation over the index range is implied. Symbolic indices which are not dummy indices are called free indices. Numeric indices are neither dummy nor free indices.

To be recognized as a dummy index pair, the two indices must be of the same class and dimension and their value must be the same single symbol (an index like 2*n+1 is never a dummy index). If the indices are of class varidx, they must also be of opposite variance.

The method .get_free_indices() returns a vector containing the free indices of an expression. It also checks that the free indices of the terms of a sum are consistent:

{
    symbol A("A"), B("B"), C("C");

    symbol i_sym("i"), j_sym("j"), k_sym("k"), l_sym("l");
    idx i(i_sym, 3), j(j_sym, 3), k(k_sym, 3), l(l_sym, 3);

    ex e = indexed(A, i, j) * indexed(B, j, k) + indexed(C, k, l, i, l);
    cout << exprseq(e.get_free_indices()) << endl;
     // -> (.i,.k)
     // 'j' and 'l' are dummy indices

    symbol mu_sym("mu"), nu_sym("nu"), rho_sym("rho"), sigma_sym("sigma");
    varidx mu(mu_sym, 4), nu(nu_sym, 4), rho(rho_sym, 4), sigma(sigma_sym, 4);

    e = indexed(A, mu, nu) * indexed(B, nu.toggle_variance(), rho)
      + indexed(C, mu, sigma, rho, sigma.toggle_variance());
    cout << exprseq(e.get_free_indices()) << endl;
     // -> (~mu,~rho)
     // 'nu' is a dummy index, but 'sigma' is not

    e = indexed(A, mu, mu);
    cout << exprseq(e.get_free_indices()) << endl;
     // -> (~mu)
     // 'mu' is not a dummy index because it appears twice with the same
     // variance

    e = indexed(A, mu, nu) + 42;
    cout << exprseq(e.get_free_indices()) << endl; // ERROR
     // this will throw an exception:
     // "add::get_free_indices: inconsistent indices in sum"
}

Simplifying indexed expressions

In addition to the few automatic simplifications that GiNaC performs on indexed expressions (such as re-ordering the indices of symmetric tensors and calculating traces and convolutions of matrices and predefined tensors) there is the method

ex ex::simplify_indexed(void);
ex ex::simplify_indexed(const scalar_products & sp);

that performs some more expensive operations:

The last point is done with the help of the scalar_products class which is used to store scalar products with known values (this is not an arithmetic class, you just pass it to simplify_indexed()):

{
    symbol A("A"), B("B"), C("C"), i_sym("i");
    idx i(i_sym, 3);

    scalar_products sp;
    sp.add(A, B, 0); // A and B are orthogonal
    sp.add(A, C, 0); // A and C are orthogonal
    sp.add(A, A, 4); // A^2 = 4 (A has length 2)

    e = indexed(A + B, i) * indexed(A + C, i);
    cout << e << endl;
     // -> (B+A).i*(A+C).i

    cout << e.expand(expand_options::expand_indexed).simplify_indexed(sp)
         << endl;
     // -> 4+C.i*B.i
}

The scalar_products object sp acts as a storage for the scalar products added to it with the .add() method. This method takes three arguments: the two expressions of which the scalar product is taken, and the expression to replace it with. After sp.add(A, B, 0), simplify_indexed() will replace all scalar products of indexed objects that have the symbols A and B as base expressions with the single value 0. The number, type and dimension of the indices doesn't matter; A~mu~nu*B.mu.nu would also be replaced by 0.

The example above also illustrates a feature of the expand() method: if passed the expand_indexed option it will distribute indices over sums, so (A+B).i becomes A.i+B.i.

Predefined tensors

Some frequently used special tensors such as the delta, epsilon and metric tensors are predefined in GiNaC. They have special properties when contracted with other tensor expressions and some of them have constant matrix representations (they will evaluate to a number when numeric indices are specified).

Delta tensor

The delta tensor takes two indices, is symmetric and has the matrix representation diag(1,1,1,...). It is constructed by the function delta_tensor():

{
    symbol A("A"), B("B");

    idx i(symbol("i"), 3), j(symbol("j"), 3),
        k(symbol("k"), 3), l(symbol("l"), 3);

    ex e = indexed(A, i, j) * indexed(B, k, l)
         * delta_tensor(i, k) * delta_tensor(j, l) << endl;
    cout << e.simplify_indexed() << endl;
     // -> B.i.j*A.i.j

    cout << delta_tensor(i, i) << endl;
     // -> 3
}

General metric tensor

The function metric_tensor() creates a general symmetric metric tensor with two indices that can be used to raise/lower tensor indices. The metric tensor is denoted as g in the output and if its indices are of mixed variance it is automatically replaced by a delta tensor:

{
    symbol A("A");

    varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4);

    ex e = metric_tensor(mu, nu) * indexed(A, nu.toggle_variance(), rho);
    cout << e.simplify_indexed() << endl;
     // -> A~mu~rho

    e = delta_tensor(mu, nu.toggle_variance()) * metric_tensor(nu, rho);
    cout << e.simplify_indexed() << endl;
     // -> g~mu~rho

    e = metric_tensor(mu.toggle_variance(), nu.toggle_variance())
      * metric_tensor(nu, rho);
    cout << e.simplify_indexed() << endl;
     // -> delta.mu~rho

    e = metric_tensor(nu.toggle_variance(), rho.toggle_variance())
      * metric_tensor(mu, nu) * (delta_tensor(mu.toggle_variance(), rho)
        + indexed(A, mu.toggle_variance(), rho));
    cout << e.simplify_indexed() << endl;
     // -> 4+A.rho~rho
}

Minkowski metric tensor

The Minkowski metric tensor is a special metric tensor with a constant matrix representation which is either diag(1, -1, -1, ...) (negative signature, the default) or diag(-1, 1, 1, ...) (positive signature). It is created with the function lorentz_g() (although it is output as eta):

{
    varidx mu(symbol("mu"), 4);

    e = delta_tensor(varidx(0, 4), mu.toggle_variance())
      * lorentz_g(mu, varidx(0, 4));       // negative signature
    cout << e.simplify_indexed() << endl;
     // -> 1

    e = delta_tensor(varidx(0, 4), mu.toggle_variance())
      * lorentz_g(mu, varidx(0, 4), true); // positive signature
    cout << e.simplify_indexed() << endl;
     // -> -1
}

Epsilon tensor

The epsilon tensor is totally antisymmetric, its number of indices is equal to the dimension of the index space (the indices must all be of the same numeric dimension), and eps.1.2.3... (resp. eps~0~1~2...) is defined to be 1. Its behaviour with indices that have a variance also depends on the signature of the metric. Epsilon tensors are output as eps.

There are three functions defined to create epsilon tensors in 2, 3 and 4 dimensions:

ex epsilon_tensor(const ex & i1, const ex & i2);
ex epsilon_tensor(const ex & i1, const ex & i2, const ex & i3);
ex lorentz_eps(const ex & i1, const ex & i2, const ex & i3, const ex & i4, bool pos_sig = false);

The first two functions create an epsilon tensor in 2 or 3 Euclidean dimensions, the last function creates an epsilon tensor in a 4-dimensional Minkowski space (the last bool argument specifies whether the metric has negative or positive signature, as in the case of the Minkowski metric tensor).

Linear algebra

The matrix class can be used with indices to do some simple linear algebra (linear combinations and products of vectors and matrices, traces and scalar products):

{
    idx i(symbol("i"), 2), j(symbol("j"), 2);
    symbol x("x"), y("y");

    matrix A(2, 2, lst(1, 2, 3, 4)), X(2, 1, lst(x, y));

    cout << indexed(A, i, i) << endl;
     // -> 5

    ex e = indexed(A, i, j) * indexed(X, j);
    cout << e.simplify_indexed() << endl;
     // -> [[ [[2*y+x]], [[4*y+3*x]] ]].i

    e = indexed(A, i, j) * indexed(X, i) + indexed(X, j) * 2;
    cout << e.simplify_indexed() << endl;
     // -> [[ [[3*y+3*x,6*y+2*x]] ]].j
}

You can of course obtain the same results with the matrix::add(), matrix::mul() and matrix::trace() methods but with indices you don't have to worry about transposing matrices.

Matrix indices always start at 0 and their dimension must match the number of rows/columns of the matrix. Matrices with one row or one column are vectors and can have one or two indices (it doesn't matter whether it's a row or a column vector). Other matrices must have two indices.

You should be careful when using indices with variance on matrices. GiNaC doesn't look at the variance and doesn't know that F~mu~nu and F.mu.nu are different matrices. In this case you should use only one form for F and explicitly multiply it with a matrix representation of the metric tensor.


Node:Methods and Functions, Next:, Previous:Indexed objects, Up:Top

Methods and Functions

In this chapter the most important algorithms provided by GiNaC will be described. Some of them are implemented as functions on expressions, others are implemented as methods provided by expression objects. If they are methods, there exists a wrapper function around it, so you can alternatively call it in a functional way as shown in the simple example:

    ...
    cout << "As method:   " << sin(1).evalf() << endl;
    cout << "As function: " << evalf(sin(1)) << endl;
    ...

The general rule is that wherever methods accept one or more parameters (arg1, arg2, ...) the order of arguments the function wrapper accepts is the same but preceded by the object to act on (object, arg1, arg2, ...). This approach is the most natural one in an OO model but it may lead to confusion for MapleV users because where they would type A:=x+1; subs(x=2,A); GiNaC would require A=x+1; subs(A,x==2); (after proper declaration of A and x). On the other hand, since MapleV returns 3 on A:=x^2+3; coeff(A,x,0); (GiNaC: A=pow(x,2)+3; coeff(A,x,0);) it is clear that MapleV is not trying to be consistent here. Also, users of MuPAD will in most cases feel more comfortable with GiNaC's convention. All function wrappers are implemented as simple inline functions which just call the corresponding method and are only provided for users uncomfortable with OO who are dead set to avoid method invocations. Generally, nested function wrappers are much harder to read than a sequence of methods and should therefore be avoided if possible. On the other hand, not everything in GiNaC is a method on class ex and sometimes calling a function cannot be avoided.


Node:Information About Expressions, Next:, Previous:Methods and Functions, Up:Methods and Functions

Getting information about expressions

Checking expression types

Sometimes it's useful to check whether a given expression is a plain number, a sum, a polynomial with integer coefficients, or of some other specific type. GiNaC provides two functions for this (the first one is actually a macro):

bool is_ex_of_type(const ex & e, TYPENAME t);
bool ex::info(unsigned flag);

When the test made by is_ex_of_type() returns true, it is safe to call one of the functions ex_to_..., where ... is one of the class names (See The Class Hierarchy, for a list of all classes). For example, assuming e is an ex:

{
    ...
    if (is_ex_of_type(e, numeric))
        numeric n = ex_to_numeric(e);
    ...
}

is_ex_of_type() allows you to check whether the top-level object of an expression e is an instance of the GiNaC class t (See The Class Hierarchy, for a list of all classes). This is most useful, e.g., for checking whether an expression is a number, a sum, or a product:

{
    symbol x("x");
    ex e1 = 42;
    ex e2 = 4*x - 3;
    is_ex_of_type(e1, numeric);  // true
    is_ex_of_type(e2, numeric);  // false
    is_ex_of_type(e1, add);      // false
    is_ex_of_type(e2, add);      // true
    is_ex_of_type(e1, mul);      // false
    is_ex_of_type(e2, mul);      // false
}

The info() method is used for checking certain attributes of expressions. The possible values for the flag argument are defined in ginac/flags.h, the most important being explained in the following table:

Flag Returns true if the object is...
numeric ...a number (same as is_ex_of_type(..., numeric))
real ...a real integer, rational or float (i.e. is not complex)
rational ...an exact rational number (integers are rational, too)
integer ...a (non-complex) integer
crational ...an exact (complex) rational number (such as 2/3+7/2*I)
cinteger ...a (complex) integer (such as 2-3*I)
positive ...not complex and greater than 0
negative ...not complex and less than 0
nonnegative ...not complex and greater than or equal to 0
posint ...an integer greater than 0
negint ...an integer less than 0
nonnegint ...an integer greater than or equal to 0
even ...an even integer
odd ...an odd integer
prime ...a prime integer (probabilistic primality test)
relation ...a relation (same as is_ex_of_type(..., relational))
relation_equal ...a == relation
relation_not_equal ...a != relation
relation_less ...a < relation
relation_less_or_equal ...a <= relation
relation_greater ...a > relation
relation_greater_or_equal ...a >= relation
symbol ...a symbol (same as is_ex_of_type(..., symbol))
list ...a list (same as is_ex_of_type(..., lst))
polynomial ...a polynomial (i.e. only consists of sums and products of numbers and symbols with positive integer powers)
integer_polynomial ...a polynomial with (non-complex) integer coefficients
cinteger_polynomial ...a polynomial with (possibly complex) integer coefficients (such as 2-3*I)
rational_polynomial ...a polynomial with (non-complex) rational coefficients
crational_polynomial ...a polynomial with (possibly complex) rational coefficients (such as 2/3+7/2*I)
rational_function ...a rational function (x+y, z/(x+y))
algebraic ...an algebraic object (sqrt(2), sqrt(x)-1)

Accessing subexpressions

GiNaC provides the two methods

unsigned ex::nops();
ex ex::op(unsigned i);

for accessing the subexpressions in the container-like GiNaC classes like add, mul, lst, and function. nops() determines the number of subexpressions (operands) contained, while op() returns the i-th (0..nops()-1) subexpression. In the case of a power object, op(0) will return the basis and op(1) the exponent. For indexed objects, op(0) is the base expression and op(i), i>0 are the indices.

The left-hand and right-hand side expressions of objects of class relational (and only of these) can also be accessed with the methods

ex ex::lhs();
ex ex::rhs();

Finally, the method

bool ex::has(const ex & other);

checks whether an expression contains the given subexpression other. This only works reliably if other is of an atomic class such as a numeric or a symbol. It is, e.g., not possible to verify that a+b+c contains a+c (or a+b) as a subexpression.

Comparing expressions

Expressions can be compared with the usual C++ relational operators like ==, >, and < but if the expressions contain symbols, the result is usually not determinable and the result will be false, except in the case of the != operator. You should also be aware that GiNaC will only do the most trivial test for equality (subtracting both expressions), so something like (pow(x,2)+x)/x==x+1 will return false.

Actually, if you construct an expression like a == b, this will be represented by an object of the relational class (See Relations.) which is not evaluated until (explicitly or implicitely) cast to a bool.

There are also two methods

bool ex::is_equal(const ex & other);
bool ex::is_zero();

for checking whether one expression is equal to another, or equal to zero, respectively.

Warning: You will also find an ex::compare() method in the GiNaC header files. This method is however only to be used internally by GiNaC to establish a canonical sort order for terms, and using it to compare expressions will give very surprising results.


Node:Substituting Symbols, Next:, Previous:Information About Expressions, Up:Methods and Functions

Substituting symbols

Symbols can be replaced with expressions via the .subs() method:

ex ex::subs(const ex & e);
ex ex::subs(const lst & syms, const lst & repls);

In the first form, subs() accepts a relational of the form symbol == expression or a lst of such relationals. E.g.

{
    symbol x("x"), y("y");
    ex e1 = 2*x^2-4*x+3;
    cout << "e1(7) = " << e1.subs(x == 7) << endl;
    ex e2 = x*y + x;
    cout << "e2(-2, 4) = " << e2.subs(lst(x == -2, y == 4)) << endl;
}

will print 73 and -10, respectively.

If you specify multiple substitutions, they are performed in parallel, so e.g. subs(lst(x == y, y == x)) exchanges x and y.

The second form of subs() takes two lists, one for the symbols and one for the expressions to be substituted (both lists must contain the same number of elements). Using this form, you would write subs(lst(x, y), lst(y, x)) to exchange x and y.


Node:Polynomial Arithmetic, Next:, Previous:Substituting Symbols, Up:Methods and Functions

Polynomial arithmetic

Expanding and collecting

A polynomial in one or more variables has many equivalent representations. Some useful ones serve a specific purpose. Consider for example the trivariate polynomial 4*x*y + x*z + 20*y^2 + 21*y*z + 4*z^2 (written down here in output-style). It is equivalent to the factorized polynomial (x + 5*y + 4*z)*(4*y + z). Other representations are the recursive ones where one collects for exponents in one of the three variable. Since the factors are themselves polynomials in the remaining two variables the procedure can be repeated. In our expample, two possibilities would be (4*y + z)*x + 20*y^2 + 21*y*z + 4*z^2 and 20*y^2 + (21*z + 4*x)*y + 4*z^2 + x*z.

To bring an expression into expanded form, its method

ex ex::expand();

may be called. In our example above, this corresponds to 4*x*y + x*z + 20*y^2 + 21*y*z + 4*z^2. Again, since the canonical form in GiNaC is not easily guessable you should be prepared to see different orderings of terms in such sums!

Another useful representation of multivariate polynomials is as a univariate polynomial in one of the variables with the coefficients being polynomials in the remaining variables. The method collect() accomplishes this task:

ex ex::collect(const ex & s);

Note that the original polynomial needs to be in expanded form in order to be able to find the coefficients properly.

Degree and coefficients

The degree and low degree of a polynomial can be obtained using the two methods

int ex::degree(const ex & s);
int ex::ldegree(const ex & s);

which also work reliably on non-expanded input polynomials (they even work on rational functions, returning the asymptotic degree). To extract a coefficient with a certain power from an expanded polynomial you use

ex ex::coeff(const ex & s, int n);

You can also obtain the leading and trailing coefficients with the methods

ex ex::lcoeff(const ex & s);
ex ex::tcoeff(const ex & s);

which are equivalent to coeff(s, degree(s)) and coeff(s, ldegree(s)), respectively.

An application is illustrated in the next example, where a multivariate polynomial is analyzed:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y");
    ex PolyInp = 4*pow(x,3)*y + 5*x*pow(y,2) + 3*y
                 - pow(x+y,2) + 2*pow(y+2,2) - 8;
    ex Poly = PolyInp.expand();

    for (int i=Poly.ldegree(x); i<=Poly.degree(x); ++i) {
        cout << "The x^" << i << "-coefficient is "
             << Poly.coeff(x,i) << endl;
    }
    cout << "As polynomial in y: "
         << Poly.collect(y) << endl;
}

When run, it returns an output in the following fashion:

The x^0-coefficient is y^2+11*y
The x^1-coefficient is 5*y^2-2*y
The x^2-coefficient is -1
The x^3-coefficient is 4*y
As polynomial in y: -x^2+(5*x+1)*y^2+(-2*x+4*x^3+11)*y

As always, the exact output may vary between different versions of GiNaC or even from run to run since the internal canonical ordering is not within the user's sphere of influence.

Polynomial division

The two functions

ex quo(const ex & a, const ex & b, const symbol & x);
ex rem(const ex & a, const ex & b, const symbol & x);

compute the quotient and remainder of univariate polynomials in the variable x. The results satisfy a = b*quo(a, b, x) + rem(a, b, x).

The additional function

ex prem(const ex & a, const ex & b, const symbol & x);

computes the pseudo-remainder of a and b which satisfies c*a = b*q + prem(a, b, x), where c = b.lcoeff(x) ^ (a.degree(x) - b.degree(x) + 1).

Exact division of multivariate polynomials is performed by the function

bool divide(const ex & a, const ex & b, ex & q);

If b divides a over the rationals, this function returns true and returns the quotient in the variable q. Otherwise it returns false in which case the value of q is undefined.

Unit, content and primitive part

The methods

ex ex::unit(const symbol & x);
ex ex::content(const symbol & x);
ex ex::primpart(const symbol & x);

return the unit part, content part, and primitive polynomial of a multivariate polynomial with respect to the variable x (the unit part being the sign of the leading coefficient, the content part being the GCD of the coefficients, and the primitive polynomial being the input polynomial divided by the unit and content parts). The product of unit, content, and primitive part is the original polynomial.

GCD and LCM

The functions for polynomial greatest common divisor and least common multiple have the synopsis

ex gcd(const ex & a, const ex & b);
ex lcm(const ex & a, const ex & b);

The functions gcd() and lcm() accept two expressions a and b as arguments and return a new expression, their greatest common divisor or least common multiple, respectively. If the polynomials a and b are coprime gcd(a,b) returns 1 and lcm(a,b) returns the product of a and b.

#include <ginac/ginac.h>
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y"), z("z");
    ex P_a = 4*x*y + x*z + 20*pow(y, 2) + 21*y*z + 4*pow(z, 2);
    ex P_b = x*y + 3*x*z + 5*pow(y, 2) + 19*y*z + 12*pow(z, 2);

    ex P_gcd = gcd(P_a, P_b);
    // x + 5*y + 4*z
    ex P_lcm = lcm(P_a, P_b);
    // 4*x*y^2 + 13*y*x*z + 20*y^3 + 81*y^2*z + 67*y*z^2 + 3*x*z^2 + 12*z^3
}

Square-free decomposition

GiNaC still lacks proper factorization support. Some form of factorization is, however, easily implemented by noting that factors appearing in a polynomial with power two or more also appear in the derivative and hence can easily be found by computing the GCD of the original polynomial and its derivatives. Any system has an interface for this so called square-free factorization. So we provide one, too:

ex sqrfree(const ex & a, const lst & l = lst());
Here is an example that by the way illustrates how the result may depend on the order of differentiation:
    ...
    symbol x("x"), y("y");
    ex BiVarPol = expand(pow(x-2*y*x,3) * pow(x+y,2) * (x-y));

    cout << sqrfree(BiVarPol, lst(x,y)) << endl;
     // -> (y+x)^2*(-1+6*y+8*y^3-12*y^2)*(y-x)*x^3

    cout << sqrfree(BiVarPol, lst(y,x)) << endl;
     // -> (1-2*y)^3*(y+x)^2*(-y+x)*x^3

    cout << sqrfree(BiVarPol) << endl;
     // -> depending on luck, any of the above
    ...


Node:Rational Expressions, Next:, Previous:Polynomial Arithmetic, Up:Methods and Functions

Rational expressions

The normal method

Some basic form of simplification of expressions is called for frequently. GiNaC provides the method .normal(), which converts a rational function into an equivalent rational function of the form numerator/denominator where numerator and denominator are coprime. If the input expression is already a fraction, it just finds the GCD of numerator and denominator and cancels it, otherwise it performs fraction addition and multiplication.

.normal() can also be used on expressions which are not rational functions as it will replace all non-rational objects (like functions or non-integer powers) by temporary symbols to bring the expression to the domain of rational functions before performing the normalization, and re-substituting these symbols afterwards. This algorithm is also available as a separate method .to_rational(), described below.

This means that both expressions t1 and t2 are indeed simplified in this little program:

#include <ginac/ginac.h>
using namespace GiNaC;

int main()
{
    symbol x("x");
    ex t1 = (pow(x,2) + 2*x + 1)/(x + 1);
    ex t2 = (pow(sin(x),2) + 2*sin(x) + 1)/(sin(x) + 1);
    std::cout << "t1 is " << t1.normal() << std::endl;
    std::cout << "t2 is " << t2.normal() << std::endl;
}

Of course this works for multivariate polynomials too, so the ratio of the sample-polynomials from the section about GCD and LCM above would be normalized to P_a/P_b = (4*y+z)/(y+3*z).

Numerator and denominator

The numerator and denominator of an expression can be obtained with

ex ex::numer();
ex ex::denom();

These functions will first normalize the expression as described above and then return the numerator or denominator, respectively.

Converting to a rational expression

Some of the methods described so far only work on polynomials or rational functions. GiNaC provides a way to extend the domain of these functions to general expressions by using the temporary replacement algorithm described above. You do this by calling

ex ex::to_rational(lst &l);

on the expression to be converted. The supplied lst will be filled with the generated temporary symbols and their replacement expressions in a format that can be used directly for the subs() method. It can also already contain a list of replacements from an earlier application of .to_rational(), so it's possible to use it on multiple expressions and get consistent results.

For example,

{
    symbol x("x");
    ex a = pow(sin(x), 2) - pow(cos(x), 2);
    ex b = sin(x) + cos(x);
    ex q;
    lst l;
    divide(a.to_rational(l), b.to_rational(l), q);
    cout << q.subs(l) << endl;
}

will print sin(x)-cos(x).


Node:Symbolic Differentiation, Next:, Previous:Rational Expressions, Up:Methods and Functions

Symbolic differentiation

GiNaC's objects know how to differentiate themselves. Thus, a polynomial (class add) knows that its derivative is the sum of the derivatives of all the monomials:

#include <ginac/ginac.h>
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y"), z("z");
    ex P = pow(x, 5) + pow(x, 2) + y;

    cout << P.diff(x,2) << endl;  // 20*x^3 + 2
    cout << P.diff(y) << endl;    // 1
    cout << P.diff(z) << endl;    // 0
}

If a second integer parameter n is given, the diff method returns the nth derivative.

If every object and every function is told what its derivative is, all derivatives of composed objects can be calculated using the chain rule and the product rule. Consider, for instance the expression 1/cosh(x). Since the derivative of cosh(x) is sinh(x) and the derivative of pow(x,-1) is -pow(x,-2), GiNaC can readily compute the composition. It turns out that the composition is the generating function for Euler Numbers, i.e. the so called nth Euler number is the coefficient of x^n/n! in the expansion of 1/cosh(x). We may use this identity to code a function that generates Euler numbers in just three lines:

#include <ginac/ginac.h>
using namespace GiNaC;

ex EulerNumber(unsigned n)
{
    symbol x;
    const ex generator = pow(cosh(x),-1);
    return generator.diff(x,n).subs(x==0);
}

int main()
{
    for (unsigned i=0; i<11; i+=2)
        std::cout << EulerNumber(i) << std::endl;
    return 0;
}

When you run it, it produces the sequence 1, -1, 5, -61, 1385, -50521. We increment the loop variable i by two since all odd Euler numbers vanish anyways.


Node:Series Expansion, Next:, Previous:Symbolic Differentiation, Up:Methods and Functions

Series expansion

Expressions know how to expand themselves as a Taylor series or (more generally) a Laurent series. As in most conventional Computer Algebra Systems, no distinction is made between those two. There is a class of its own for storing such series (class pseries) and a built-in function (called Order) for storing the order term of the series. As a consequence, if you want to work with series, i.e. multiply two series, you need to call the method ex::series again to convert it to a series object with the usual structure (expansion plus order term). A sample application from special relativity could read:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
    symbol v("v"), c("c");

    ex gamma = 1/sqrt(1 - pow(v/c,2));
    ex mass_nonrel = gamma.series(v==0, 10);

    cout << "the relativistic mass increase with v is " << endl
         << mass_nonrel << endl;

    cout << "the inverse square of this series is " << endl
         << pow(mass_nonrel,-2).series(v==0, 10) << endl;
}

Only calling the series method makes the last output simplify to 1-v^2/c^2+O(v^10), without that call we would just have a long series raised to the power -2.

As another instructive application, let us calculate the numerical value of Archimedes' constant (for which there already exists the built-in constant Pi) using Méchain's amazing formula Pi==16*atan(1/5)-4*atan(1/239). We may expand the arcus tangent around 0 and insert the fractions 1/5 and 1/239. But, as we have seen, a series in GiNaC carries an order term with it and the question arises what the system is supposed to do when the fractions are plugged into that order term. The solution is to use the function series_to_poly() to simply strip the order term off:

#include <ginac/ginac.h>
using namespace GiNaC;

ex mechain_pi(int degr)
{
    symbol x;
    ex pi_expansion = series_to_poly(atan(x).series(x,degr));
    ex pi_approx = 16*pi_expansion.subs(x==numeric(1,5))
                   -4*pi_expansion.subs(x==numeric(1,239));
    return pi_approx;
}

int main()
{
    using std::cout;  // just for fun, another way of...
    using std::endl;  // ...dealing with this namespace std.
    ex pi_frac;
    for (int i=2; i<12; i+=2) {
        pi_frac = mechain_pi(i);
        cout << i << ":\t" << pi_frac << endl
             << "\t" << pi_frac.evalf() << endl;
    }
    return 0;
}

Note how we just called .series(x,degr) instead of .series(x==0,degr). This is a simple shortcut for ex's method series(): if the first argument is a symbol the expression is expanded in that symbol around point 0. When you run this program, it will type out:

2:      3804/1195
        3.1832635983263598326
4:      5359397032/1706489875
        3.1405970293260603143
6:      38279241713339684/12184551018734375
        3.141621029325034425
8:      76528487109180192540976/24359780855939418203125
        3.141591772182177295
10:     327853873402258685803048818236/104359128170408663038552734375
        3.1415926824043995174


Node:Built-in Functions, Next:, Previous:Series Expansion, Up:Methods and Functions

Predefined mathematical functions

GiNaC contains the following predefined mathematical functions:

Name Function
abs(x) absolute value
csgn(x) complex sign
sqrt(x) square root (not a GiNaC function proper but equivalent to pow(x, numeric(1, 2))
sin(x) sine
cos(x) cosine
tan(x) tangent
asin(x) inverse sine
acos(x) inverse cosine
atan(x) inverse tangent
atan2(y, x) inverse tangent with two arguments
sinh(x) hyperbolic sine
cosh(x) hyperbolic cosine
tanh(x) hyperbolic tangent
asinh(x) inverse hyperbolic sine
acosh(x) inverse hyperbolic cosine
atanh(x) inverse hyperbolic tangent
exp(x) exponential function
log(x) natural logarithm
Li2(x) Dilogarithm
zeta(x) Riemann's zeta function
zeta(n, x) derivatives of Riemann's zeta function
tgamma(x) Gamma function
lgamma(x) logarithm of Gamma function
beta(x, y) Beta function (tgamma(x)*tgamma(y)/tgamma(x+y))
psi(x) psi (digamma) function
psi(n, x) derivatives of psi function (polygamma functions)
factorial(n) factorial function
binomial(n, m) binomial coefficients
Order(x) order term function in truncated power series
Derivative(x, l) inert partial differentiation operator (used internally)

For functions that have a branch cut in the complex plane GiNaC follows the conventions for C++ as defined in the ANSI standard as far as possible. In particular: the natural logarithm (log) and the square root (sqrt) both have their branch cuts running along the negative real axis where the points on the axis itself belong to the upper part (i.e. continuous with quadrant II). The inverse trigonometric and hyperbolic functions are not defined for complex arguments by the C++ standard, however. In GiNaC we follow the conventions used by CLN, which in turn follow the carefully designed definitions in the Common Lisp standard. It should be noted that this convention is identical to the one used by the C99 standard and by most serious CAS. It is to be expected that future revisions of the C++ standard incorporate these functions in the complex domain in a manner compatible with C99.


Node:Input/Output, Next:, Previous:Built-in Functions, Up:Methods and Functions

Input and output of expressions

Expression output

The easiest way to print an expression is to write it to a stream:

{
    symbol x("x");
    ex e = 4.5+pow(x,2)*3/2;
    cout << e << endl;    // prints '4.5+3/2*x^2'
    // ...

The output format is identical to the ginsh input syntax and to that used by most computer algebra systems, but not directly pastable into a GiNaC C++ program (note that in the above example, pow(x,2) is printed as x^2).

To print an expression in a way that can be directly used in a C or C++ program, you use the method

void ex::printcsrc(ostream & os, unsigned type, const char *name);

This outputs a line in the form of a variable definition <type> <name> = <expression>. The possible types are defined in ginac/flags.h (csrc_types) and mostly affect the way in which floating point numbers are written:

    // ...
    e.printcsrc(cout, csrc_types::ctype_float, "f");
    e.printcsrc(cout, csrc_types::ctype_double, "d");
    e.printcsrc(cout, csrc_types::ctype_cl_N, "n");
    // ...

The above example will produce (note the x^2 being converted to x*x):

float f = (3.000000e+00/2.000000e+00)*(x*x)+4.500000e+00;
double d = (3.000000e+00/2.000000e+00)*(x*x)+4.500000e+00;
cl_N n = (cl_F("3.0")/cl_F("2.0"))*(x*x)+cl_F("4.5");

Finally, there are the two methods printraw() and printtree() intended for GiNaC developers, that provide a dump of the internal structure of an expression for debugging purposes:

    // ...
    e.printraw(cout); cout << endl << endl;
    e.printtree(cout);
}

produces

ex(+((power(ex(symbol(name=x,serial=1,hash=150875740,flags=11)),ex(numeric(2)),hash=2,flags=3),numeric(3/2)),,hash=0,flags=3))

type=Q25GiNaC3add, hash=0 (0x0), flags=3, nops=2
    power: hash=2 (0x2), flags=3
        x (symbol): serial=1, hash=150875740 (0x8fe2e5c), flags=11
        2 (numeric): hash=2147483714 (0x80000042), flags=11
    3/2 (numeric): hash=2147483745 (0x80000061), flags=11
    -----
    overall_coeff
    4.5L0 (numeric): hash=2147483723 (0x8000004b), flags=11
    =====

The printtree() method is also available in ginsh as the print() function.

Expression input

GiNaC provides no way to directly read an expression from a stream because you will usually want the user to be able to enter something like 2*x+sin(y) and have the x and y correspond to the symbols x and y you defined in your program and there is no way to specify the desired symbols to the >> stream input operator.

Instead, GiNaC lets you construct an expression from a string, specifying the list of symbols to be used:

{
    symbol x("x"), y("y");
    ex e("2*x+sin(y)", lst(x, y));
}

The input syntax is the same as that used by ginsh and the stream output operator <<. The symbols in the string are matched by name to the symbols in the list and if GiNaC encounters a symbol not specified in the list it will throw an exception.

With this constructor, it's also easy to implement interactive GiNaC programs:

#include <iostream>
#include <string>
#include <stdexcept>
#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
     symbol x("x");
     string s;

     cout << "Enter an expression containing 'x': ";
     getline(cin, s);

     try {
         ex e(s, lst(x));
         cout << "The derivative of " << e << " with respect to x is ";
         cout << e.diff(x) << ".\n";
     } catch (exception &p) {
         cerr << p.what() << endl;
     }
}

Archiving

GiNaC allows creating archives of expressions which can be stored to or retrieved from files. To create an archive, you declare an object of class archive and archive expressions in it, giving each expression a unique name:

#include <fstream>
using namespace std;
#include <ginac/ginac.h>
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y"), z("z");

    ex foo = sin(x + 2*y) + 3*z + 41;
    ex bar = foo + 1;

    archive a;
    a.archive_ex(foo, "foo");
    a.archive_ex(bar, "the second one");
    // ...

The archive can then be written to a file:

    // ...
    ofstream out("foobar.gar");
    out << a;
    out.close();
    // ...

The file foobar.gar contains all information that is needed to reconstruct the expressions foo and bar.

The tool viewgar that comes with GiNaC can be used to view the contents of GiNaC archive files:

$ viewgar foobar.gar
foo = 41+sin(x+2*y)+3*z
the second one = 42+sin(x+2*y)+3*z

The point of writing archive files is of course that they can later be read in again:

    // ...
    archive a2;
    ifstream in("foobar.gar");
    in >> a2;
    // ...

And the stored expressions can be retrieved by their name:

    // ...
    lst syms(x, y);

    ex ex1 = a2.unarchive_ex(syms, "foo");
    ex ex2 = a2.unarchive_ex(syms, "the second one");

    cout << ex1 << endl;              // prints "41+sin(x+2*y)+3*z"
    cout << ex2 << endl;              // prints "42+sin(x+2*y)+3*z"
    cout << ex1.subs(x == 2) << endl; // prints "41+sin(2+2*y)+3*z"
}

Note that you have to supply a list of the symbols which are to be inserted in the expressions. Symbols in archives are stored by their name only and if you don't specify which symbols you have, unarchiving the expression will create new symbols with that name. E.g. if you hadn't included x in the syms list above, the ex1.subs(x == 2) statement would have had no effect because the x in ex1 would have been a different symbol than the x which was defined at the beginning of the program, altough both would appear as x when printed.


Node:Extending GiNaC, Next:, Previous:Input/Output, Up:Top

Extending GiNaC

By reading so far you should have gotten a fairly good understanding of GiNaC's design-patterns. From here on you should start reading the sources. All we can do now is issue some recommendations how to tackle GiNaC's many loose ends in order to fulfill everybody's dreams. If you develop some useful extension please don't hesitate to contact the GiNaC authors--they will happily incorporate them into future versions.


Node:What does not belong into GiNaC, Next:, Previous:Extending GiNaC, Up:Extending GiNaC

What doesn't belong into GiNaC

First of all, GiNaC's name must be read literally. It is designed to be a library for use within C++. The tiny ginsh accompanying GiNaC makes this even more clear: it doesn't even attempt to provide a language. There are no loops or conditional expressions in ginsh, it is merely a window into the library for the programmer to test stuff (or to show off). Still, the design of a complete CAS with a language of its own, graphical capabilites and all this on top of GiNaC is possible and is without doubt a nice project for the future.

There are many built-in functions in GiNaC that do not know how to evaluate themselves numerically to a precision declared at runtime (using Digits). Some may be evaluated at certain points, but not generally. This ought to be fixed. However, doing numerical computations with GiNaC's quite abstract classes is doomed to be inefficient. For this purpose, the underlying foundation classes provided by CLN are much better suited.


Node:Symbolic functions, Next:, Previous:What does not belong into GiNaC, Up:Extending GiNaC

Symbolic functions

The easiest and most instructive way to start with is probably to implement your own function. GiNaC's functions are objects of class function. The preprocessor is then used to convert the function names to objects with a corresponding serial number that is used internally to identify them. You usually need not worry about this number. New functions may be inserted into the system via a kind of `registry'. It is your responsibility to care for some functions that are called when the user invokes certain methods. These are usual C++-functions accepting a number of ex as arguments and returning one ex. As an example, if we have a look at a simplified implementation of the cosine trigonometric function, we first need a function that is called when one wishes to eval it. It could look something like this:

static ex cos_eval_method(const ex & x)
{
    // if (!x%(2*Pi)) return 1
    // if (!x%Pi) return -1
    // if (!x%Pi/2) return 0
    // care for other cases...
    return cos(x).hold();
}

The last line returns cos(x) if we don't know what else to do and stops a potential recursive evaluation by saying .hold(), which sets a flag to the expression signaling that it has been evaluated. We should also implement a method for numerical evaluation and since we are lazy we sweep the problem under the rug by calling someone else's function that does so, in this case the one in class numeric:

static ex cos_evalf(const ex & x)
{
    return cos(ex_to_numeric(x));
}

Differentiation will surely turn up and so we need to tell cos what the first derivative is (higher derivatives (.diff(x,3) for instance are then handled automatically by basic::diff and ex::diff):

static ex cos_deriv(const ex & x, unsigned diff_param)
{
    return -sin(x);
}

The second parameter is obligatory but uninteresting at this point. It specifies which parameter to differentiate in a partial derivative in case the function has more than one parameter and its main application is for correct handling of the chain rule. For Taylor expansion, it is enough to know how to differentiate. But if the function you want to implement does have a pole somewhere in the complex plane, you need to write another method for Laurent expansion around that point.

Now that all the ingredients for cos have been set up, we need to tell the system about it. This is done by a macro and we are not going to descibe how it expands, please consult your preprocessor if you are curious:

REGISTER_FUNCTION(cos, eval_func(cos_eval).
                       evalf_func(cos_evalf).
                       derivative_func(cos_deriv));

The first argument is the function's name used for calling it and for output. The second binds the corresponding methods as options to this object. Options are separated by a dot and can be given in an arbitrary order. GiNaC functions understand several more options which are always specified as .option(params), for example a method for series expansion .series_func(cos_series). Again, if no series expansion method is given, GiNaC defaults to simple Taylor expansion, which is correct if there are no poles involved as is the case for the cos function. The way GiNaC handles poles in case there are any is best understood by studying one of the examples, like the Gamma (tgamma) function for instance. (In essence the function first checks if there is a pole at the evaluation point and falls back to Taylor expansion if there isn't. Then, the pole is regularized by some suitable transformation.) Also, the new function needs to be declared somewhere. This may also be done by a convenient preprocessor macro:

DECLARE_FUNCTION_1P(cos)

The suffix _1P stands for one parameter. Of course, this implementation of cos is very incomplete and lacks several safety mechanisms. Please, have a look at the real implementation in GiNaC. (By the way: in case you are worrying about all the macros above we can assure you that functions are GiNaC's most macro-intense classes. We have done our best to avoid macros where we can.)


Node:Adding classes, Next:, Previous:Symbolic functions, Up:Extending GiNaC

Adding classes

If you are doing some very specialized things with GiNaC you may find that you have to implement your own algebraic classes to fit your needs. This section will explain how to do this by giving the example of a simple 'string' class. After reading this section you will know how to properly declare a GiNaC class and what the minimum required member functions are that you have to implement. We only cover the implementation of a 'leaf' class here (i.e. one that doesn't contain subexpressions). Creating a container class like, for example, a class representing tensor products is more involved but this section should give you enough information so you can consult the source to GiNaC's predefined classes if you want to implement something more complicated.

GiNaC's run-time type information system

All algebraic classes (that is, all classes that can appear in expressions) in GiNaC are direct or indirect subclasses of the class basic. So a basic * (which is essentially what an ex is) represents a generic pointer to an algebraic class. Occasionally it is necessary to find out what the class of an object pointed to by a basic * really is. Also, for the unarchiving of expressions it must be possible to find the unarchive() function of a class given the class name (as a string). A system that provides this kind of information is called a run-time type information (RTTI) system. The C++ language provides such a thing (see the standard header file <typeinfo>) but for efficiency reasons GiNaC implements its own, simpler RTTI.

The RTTI in GiNaC is based on two mechanisms:

The disadvantage of this proprietary RTTI implementation is that there's a little more to do when implementing new classes (C++'s RTTI works more or less automatic) but don't worry, most of the work is simplified by macros.

A minimalistic example

Now we will start implementing a new class mystring that allows placing character strings in algebraic expressions (this is not very useful, but it's just an example). This class will be a direct subclass of basic. You can use this sample implementation as a starting point for your own classes.

The code snippets given here assume that you have included some header files as follows:

#include <iostream>
#include <string>
#include <stdexcept>
using namespace std;

#include <ginac/ginac.h>
using namespace GiNaC;

The first thing we have to do is to define a tinfo_key for our new class. This can be any arbitrary unsigned number that is not already taken by one of the existing classes but it's better to come up with something that is unlikely to clash with keys that might be added in the future. The numbers in tinfos.h are modeled somewhat after the class hierarchy which is not a requirement but we are going to stick with this scheme:

const unsigned TINFO_mystring = 0x42420001U;

Now we can write down the class declaration. The class stores a C++ string and the user shall be able to construct a mystring object from a C or C++ string:

class mystring : public basic
{
    GINAC_DECLARE_REGISTERED_CLASS(mystring, basic)

public:
    mystring(const string &s);
    mystring(const char *s);

private:
    string str;
};

GIANC_IMPLEMENT_REGISTERED_CLASS(mystring, basic)

The GINAC_DECLARE_REGISTERED_CLASS and GINAC_IMPLEMENT_REGISTERED_CLASS macros are defined in registrar.h. They take the name of the class and its direct superclass as arguments and insert all required declarations for the RTTI system. The GINAC_DECLARE_REGISTERED_CLASS should be the first line after the opening brace of the class definition. The GINAC_IMPLEMENT_REGISTERED_CLASS may appear anywhere else in the source (at global scope, of course, not inside a function).

GINAC_DECLARE_REGISTERED_CLASS contains, among other things the declarations of the default and copy constructor, the destructor, the assignment operator and a couple of other functions that are required. It also defines a type inherited which refers to the superclass so you don't have to modify your code every time you shuffle around the class hierarchy. GINAC_IMPLEMENT_REGISTERED_CLASS implements the copy constructor, the destructor and the assignment operator.

Now there are nine member functions we have to implement to get a working class:

Let's proceed step-by-step. The default constructor looks like this:

mystring::mystring() : inherited(TINFO_mystring)
{
    // dynamically allocate resources here if required
}

The golden rule is that in all constructors you have to set the tinfo_key member to the TINFO_* value of your class. Otherwise it will be set by the constructor of the superclass and all hell will break loose in the RTTI. For your convenience, the basic class provides a constructor that takes a tinfo_key value, which we are using here (remember that in our case inherited = basic). If the superclass didn't have such a constructor, we would have to set the tinfo_key to the right value manually.

In the default constructor you should set all other member variables to reasonable default values (we don't need that here since our str member gets set to an empty string automatically). The constructor(s) are of course also the right place to allocate any dynamic resources you require.

Next, the destroy() function:

void mystring::destroy(bool call_parent)
{
    // free dynamically allocated resources here if required
    if (call_parent)
        inherited::destroy(call_parent);
}

This function is where we free all dynamically allocated resources. We don't have any so we're not doing anything here, but if we had, for example, used a C-style char * to store our string, this would be the place to delete[] the string storage. If call_parent is true, we have to call the destroy() function of the superclass after we're done (to mimic C++'s automatic invocation of superclass destructors where destroy() is called from outside a destructor).

The copy() function just copies over the member variables from another object:

void mystring::copy(const mystring &other)
{
    inherited::copy(other);
    str = other.str;
}

We can simply overwrite the member variables here. There's no need to worry about dynamically allocated storage. The assignment operator (which is automatically defined by GINAC_IMPLEMENT_REGISTERED_CLASS, as you recall) calls destroy() before it calls copy(). You have to explicitly call the copy() function of the superclass here so all the member variables will get copied.

Next are the three functions for archiving. You have to implement them even if you don't plan to use archives, but the minimum required implementation is really simple. First, the archiving function:

void mystring::archive(archive_node &n) const
{
    inherited::archive(n);
    n.add_string("string", str);
}

The only thing that is really required is calling the archive() function of the superclass. Optionally, you can store all information you deem necessary for representing the object into the passed archive_node. We are just storing our string here. For more information on how the archiving works, consult the archive.h header file.

The unarchiving constructor is basically the inverse of the archiving function:

mystring::mystring(const archive_node &n, const lst &sym_lst) : inherited(n, sym_lst)
{
    n.find_string("string", str);
}

If you don't need archiving, just leave this function empty (but you must invoke the unarchiving constructor of the superclass). Note that we don't have to set the tinfo_key here because it is done automatically by the unarchiving constructor of the basic class.

Finally, the unarchiving function:

ex mystring::unarchive(const archive_node &n, const lst &sym_lst)
{
    return (new mystring(n, sym_lst))->setflag(status_flags::dynallocated);
}

You don't have to understand how exactly this works. Just copy these four lines into your code literally (replacing the class name, of course). It calls the unarchiving constructor of the class and unless you are doing something very special (like matching archive_nodes to global objects) you don't need a different implementation. For those who are interested: setting the dynallocated flag puts the object under the control of GiNaC's garbage collection. It will get deleted automatically once it is no longer referenced.

Our compare_same_type() function uses a provided function to compare the string members:

int mystring::compare_same_type(const basic &other) const
{
    const mystring &o = static_cast<const mystring &>(other);
    int cmpval = str.compare(o.str);
    if (cmpval == 0)
        return 0;
    else if (cmpval < 0)
        return -1;
    else
        return 1;
}

Although this function takes a basic &, it will always be a reference to an object of exactly the same class (objects of different classes are not comparable), so the cast is safe. If this function returns 0, the two objects are considered equal (in the sense that A-B=0), so you should compare all relevant member variables.

Now the only thing missing is our two new constructors:

mystring::mystring(const string &s) : inherited(TINFO_mystring), str(s)
{
    // dynamically allocate resources here if required
}

mystring::mystring(const char *s) : inherited(TINFO_mystring), str(s)
{
    // dynamically allocate resources here if required
}

No surprises here. We set the str member from the argument and remember to pass the right tinfo_key to the basic constructor.

That's it! We now have a minimal working GiNaC class that can store strings in algebraic expressions. Let's confirm that the RTTI works:

ex e = mystring("Hello, world!");
cout << is_ex_of_type(e, mystring) << endl;
 // -> 1 (true)

cout << e.bp->class_name() << endl;
 // -> mystring

Obviously it does. Let's see what the expression e looks like:

cout << e << endl;
 // -> [mystring object]

Hm, not exactly what we expect, but of course the mystring class doesn't yet know how to print itself. This is done in the print() member function. Let's say that we wanted to print the string surrounded by double quotes:

class mystring : public basic
{
    ...
public:
    void print(ostream &os, unsigned upper_precedence) const;
    ...
};

void mystring::print(ostream &os, unsigned upper_precedence) const
{
    os << '\"' << str << '\"';
}

The upper_precedence argument is only required for container classes to correctly parenthesize the output. Let's try again to print the expression:

cout << e << endl;
 // -> "Hello, world!"

Much better. The mystring class can be used in arbitrary expressions:

e += mystring("GiNaC rulez");
cout << e << endl;
 // -> "GiNaC rulez"+"Hello, world!"

(note that GiNaC's automatic term reordering is in effect here), or even

e = pow(mystring("One string"), 2*sin(Pi-mystring("Another string")));
cout << e << endl;
 // -> "One string"^(2*sin(-"Another string"+Pi))

Whether this makes sense is debatable but remember that this is only an example. At least it allows you to implement your own symbolic algorithms for your objects.

Note that GiNaC's algebraic rules remain unchanged:

e = mystring("Wow") * mystring("Wow");
cout << e << endl;
 // -> "Wow"^2

e = pow(mystring("First")-mystring("Second"), 2);
cout << e.expand() << endl;
 // -> -2*"First"*"Second"+"First"^2+"Second"^2

There's no way to, for example, make GiNaC's add class perform string concatenation. You would have to implement this yourself.

Automatic evaluation

When dealing with objects that are just a little more complicated than the simple string objects we have implemented, chances are that you will want to have some automatic simplifications or canonicalizations performed on them. This is done in the evaluation member function eval(). Let's say that we wanted all strings automatically converted to lowercase with non-alphabetic characters stripped, and empty strings removed:

class mystring : public basic
{
    ...
public:
    ex eval(int level = 0) const;
    ...
};

ex mystring::eval(int level) const
{
    string new_str;
    for (int i=0; i<str.length(); i++) {
        char c = str[i];
        if (c >= 'A' && c <= 'Z')
            new_str += tolower(c);
        else if (c >= 'a' && c <= 'z')
            new_str += c;
    }

    if (new_str.length() == 0)
        return _ex0();
    else
        return mystring(new_str).hold();
}

The level argument is used to limit the recursion depth of the evaluation. We don't have any subexpressions in the mystring class so we are not concerned with this. If we had, we would call the eval() functions of the subexpressions with level - 1 as the argument if level != 1. The hold() member function sets a flag in the object that prevents further evaluation. Otherwise we might end up in an endless loop. When you want to return the object unmodified, use return this->hold();.

Let's confirm that it works:

ex e = mystring("Hello, world!") + mystring("!?#");
cout << e << endl;
 // -> "helloworld"

e = mystring("Wow!") + mystring("WOW") + mystring(" W ** o ** W");
cout << e << endl;
 // -> 3*"wow"

Other member functions

We have implemented only a small set of member functions to make the class work in the GiNaC framework. For a real algebraic class, there are probably some more functions that you will want to re-implement, such as evalf(), series() or op(). Have a look at basic.h or the header file of the class you want to make a subclass of to see what's there. You can, of course, also add your own new member functions. In this case you will probably want to define a little helper function like

inline const mystring &ex_to_mystring(const ex &e)
{
    return static_cast<const mystring &>(*e.bp);
}

that let's you get at the object inside an expression (after you have verified that the type is correct) so you can call member functions that are specific to the class.

That's it. May the source be with you!


Node:A Comparison With Other CAS, Next:, Previous:Adding classes, Up:Top

A Comparison With Other CAS

This chapter will give you some information on how GiNaC compares to other, traditional Computer Algebra Systems, like Maple, Mathematica or Reduce, where it has advantages and disadvantages over these systems.


Node:Advantages, Next:, Previous:A Comparison With Other CAS, Up:A Comparison With Other CAS

Advantages

GiNaC has several advantages over traditional Computer Algebra Systems, like


Node:Disadvantages, Next:, Previous:Advantages, Up:A Comparison With Other CAS

Disadvantages

Of course it also has some disadvantages:


Node:Why C++?, Next:, Previous:Disadvantages, Up:A Comparison With Other CAS

Why C++?

Why did we choose to implement GiNaC in C++ instead of Java or any other language? C++ is not perfect: type checking is not strict (casting is possible), separation between interface and implementation is not complete, object oriented design is not enforced. The main reason is the often scolded feature of operator overloading in C++. While it may be true that operating on classes with a + operator is rarely meaningful, it is perfectly suited for algebraic expressions. Writing 3x+5y as 3*x+5*y instead of x.times(3).plus(y.times(5)) looks much more natural. Furthermore, the main developers are more familiar with C++ than with any other programming language.


Node:Internal Structures, Next:, Previous:Why C++?, Up:Top

Internal Structures


Node:Expressions are reference counted, Next:, Previous:Internal Structures, Up:Internal Structures

Expressions are reference counted

An expression is extremely light-weight since internally it works like a handle to the actual representation and really holds nothing more than a pointer to some other object. What this means in practice is that whenever you create two ex and set the second equal to the first no copying process is involved. Instead, the copying takes place as soon as you try to change the second. Consider the simple sequence of code:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y"), z("z");
    ex e1, e2;

    e1 = sin(x + 2*y) + 3*z + 41;
    e2 = e1;                // e2 points to same object as e1
    cout << e2 << endl;     // prints sin(x+2*y)+3*z+41
    e2 += 1;                // e2 is copied into a new object
    cout << e2 << endl;     // prints sin(x+2*y)+3*z+42
}

The line e2 = e1; creates a second expression pointing to the object held already by e1. The time involved for this operation is therefore constant, no matter how large e1 was. Actual copying, however, must take place in the line e2 += 1; because e1 and e2 are not handles for the same object any more. This concept is called copy-on-write semantics. It increases performance considerably whenever one object occurs multiple times and represents a simple garbage collection scheme because when an ex runs out of scope its destructor checks whether other expressions handle the object it points to too and deletes the object from memory if that turns out not to be the case. A slightly less trivial example of differentiation using the chain-rule should make clear how powerful this can be:

#include <ginac/ginac.h>
using namespace std;
using namespace GiNaC;

int main()
{
    symbol x("x"), y("y");

    ex e1 = x + 3*y;
    ex e2 = pow(e1, 3);
    ex e3 = diff(sin(e2), x);   // first derivative of sin(e2) by x
    cout << e1 << endl          // prints x+3*y
         << e2 << endl          // prints (x+3*y)^3
         << e3 << endl;         // prints 3*(x+3*y)^2*cos((x+3*y)^3)
}

Here, e1 will actually be referenced three times while e2 will be referenced two times. When the power of an expression is built, that expression needs not be copied. Likewise, since the derivative of a power of an expression can be easily expressed in terms of that expression, no copying of e1 is involved when e3 is constructed. So, when e3 is constructed it will print as 3*(x+3*y)^2*cos((x+3*y)^3) but the argument of cos() only holds a reference to e2 and the factor in front is just 3*e1^2.

As a user of GiNaC, you cannot see this mechanism of copy-on-write semantics. When you insert an expression into a second expression, the result behaves exactly as if the contents of the first expression were inserted. But it may be useful to remember that this is not what happens. Knowing this will enable you to write much more efficient code. If you still have an uncertain feeling with copy-on-write semantics, we recommend you have a look at the C++-FAQ lite by Marshall Cline. Chapter 16 covers this issue and presents an implementation which is pretty close to the one in GiNaC.


Node:Internal representation of products and sums, Next:, Previous:Expressions are reference counted, Up:Internal Structures

Internal representation of products and sums

Although it should be completely transparent for the user of GiNaC a short discussion of this topic helps to understand the sources and also explain performance to a large degree. Consider the unexpanded symbolic expression 2*d^3*(4*a+5*b-3) which could naively be represented by a tree of linear containers for addition and multiplication, one container for exponentiation with base and exponent and some atomic leaves of symbols and numbers in this fashion: repnaive.png

However, doing so results in a rather deeply nested tree which will quickly become inefficient to manipulate. We can improve on this by representing the sum as a sequence of terms, each one being a pair of a purely numeric multiplicative coefficient and its rest. In the same spirit we can store the multiplication as a sequence of terms, each having a numeric exponent and a possibly complicated base, the tree becomes much more flat: reppair.png

The number 3 above the symbol d shows that mul objects are treated similarly where the coefficients are interpreted as exponents now. Addition of sums of terms or multiplication of products with numerical exponents can be coded to be very efficient with such a pair-wise representation. Internally, this handling is performed by most CAS in this way. It typically speeds up manipulations by an order of magnitude. The overall multiplicative factor 2 and the additive term -3 look somewhat out of place in this representation, however, since they are still carrying a trivial exponent and multiplicative factor 1 respectively. Within GiNaC, this is avoided by adding a field that carries an overall numeric coefficient. This results in the realistic picture of internal representation for 2*d^3*(4*a+5*b-3): repreal.png

This also allows for a better handling of numeric radicals, since sqrt(2) can now be carried along calculations. Now it should be clear, why both classes add and mul are derived from the same abstract class: the data representation is the same, only the semantics differs. In the class hierarchy, methods for polynomial expansion and the like are reimplemented for add and mul, but the data structure is inherited from expairseq.


Node:Package Tools, Next:, Previous:Internal representation of products and sums, Up:Top

Package Tools

If you are creating a software package that uses the GiNaC library, setting the correct command line options for the compiler and linker can be difficult. GiNaC includes two tools to make this process easier.


Node:ginac-config, Next:, Previous:Package Tools, Up:Package Tools

ginac-config

ginac-config is a shell script that you can use to determine the compiler and linker command line options required to compile and link a program with the GiNaC library.

ginac-config takes the following flags:

--version
Prints out the version of GiNaC installed.
--cppflags
Prints '-I' flags pointing to the installed header files.
--libs
Prints out the linker flags necessary to link a program against GiNaC.
--prefix[=PREFIX]
If PREFIX is specified, overrides the configured value of $prefix. (And of exec-prefix, unless --exec-prefix is also specified) Otherwise, prints out the configured value of $prefix.
--exec-prefix[=PREFIX]
If PREFIX is specified, overrides the configured value of $exec_prefix. Otherwise, prints out the configured value of $exec_prefix.

Typically, ginac-config will be used within a configure script, as described below. It, however, can also be used directly from the command line using backquotes to compile a simple program. For example:

c++ -o simple `ginac-config --cppflags` simple.cpp `ginac-config --libs`

This command line might expand to (for example):

cc -o simple -I/usr/local/include simple.cpp -L/usr/local/lib \
  -lginac -lcln -lstdc++

Not only is the form using ginac-config easier to type, it will work on any system, no matter how GiNaC was configured.


Node:AM_PATH_GINAC, Next:, Previous:ginac-config, Up:Package Tools

AM_PATH_GINAC

For packages configured using GNU automake, GiNaC also provides a macro to automate the process of checking for GiNaC.

AM_PATH_GINAC([MINIMUM-VERSION, [ACTION-IF-FOUND [, ACTION-IF-NOT-FOUND]]])

This macro:

This macro is in file ginac.m4 which is installed in $datadir/aclocal. Note that if automake was installed with a different --prefix than GiNaC, you will either have to manually move ginac.m4 to automake's $datadir/aclocal, or give aclocal the -I option when running it.


Node:Configure script options, Next:, Previous:AM_PATH_GINAC, Up:AM_PATH_GINAC

Configuring a package that uses AM_PATH_GINAC

Simply make sure that ginac-config is in your path, and run the configure script.

Notes:

Advanced note:


Node:Example package, Next:, Previous:Configure script options, Up:AM_PATH_GINAC

Example of a package using AM_PATH_GINAC

The following shows how to build a simple package using automake and the AM_PATH_GINAC macro. The program used here is simple.cpp:

#include <ginac/ginac.h>

int main(void)
{
    GiNaC::symbol x("x");
    GiNaC::ex a = GiNaC::sin(x);
    std::cout << "Derivative of " << a
              << " is " << a.diff(x) << std::endl;
    return 0;
}

You should first read the introductory portions of the automake Manual, if you are not already familiar with it.

Two files are needed, configure.in, which is used to build the configure script:

dnl Process this file with autoconf to produce a configure script.
AC_INIT(simple.cpp)
AM_INIT_AUTOMAKE(simple.cpp, 1.0.0)

AC_PROG_CXX
AC_PROG_INSTALL
AC_LANG_CPLUSPLUS

AM_PATH_GINAC(0.7.0, [
  LIBS="$LIBS $GINACLIB_LIBS"
  CPPFLAGS="$CPPFLAGS $GINACLIB_CPPFLAGS"
], AC_MSG_ERROR([need to have GiNaC installed]))

AC_OUTPUT(Makefile)

The only command in this which is not standard for automake is the AM_PATH_GINAC macro.

That command does the following: If a GiNaC version greater or equal than 0.7.0 is found, then it adds $GINACLIB_LIBS to $LIBS and $GINACLIB_CPPFLAGS to $CPPFLAGS. Otherwise, it dies with the error message `need to have GiNaC installed'

And the Makefile.am, which will be used to build the Makefile.

## Process this file with automake to produce Makefile.in
bin_PROGRAMS = simple
simple_SOURCES = simple.cpp

This Makefile.am, says that we are building a single executable, from a single sourcefile simple.cpp. Since every program we are building uses GiNaC we simply added the GiNaC options to $LIBS and $CPPFLAGS, but in other circumstances, we might want to specify them on a per-program basis: for instance by adding the lines:

simple_LDADD = $(GINACLIB_LIBS)
INCLUDES = $(GINACLIB_CPPFLAGS)

to the Makefile.am.

To try this example out, create a new directory and add the three files above to it.

Now execute the following commands:

$ automake --add-missing
$ aclocal
$ autoconf

You now have a package that can be built in the normal fashion

$ ./configure
$ make
$ make install


Node:Bibliography, Next:, Previous:Example package, Up:Top

Bibliography


Node:Concept Index, Previous:Bibliography, Up:Top

Concept Index


Footnotes

  1. Uninstallation does not work after you have called make distclean since the Makefile is itself generated by the configuration from Makefile.in and hence deleted by make distclean. There are two obvious ways out of this dilemma. First, you can run the configuration again with the same PREFIX thus creating a Makefile with a working uninstall target. Second, you can do it by hand since you now know where all the files went during installation.

  2. This is because CLN uses PROVIDE/REQUIRE like macros to let the compiler gather all static initializations, which works for GNU C++ only.