1 \input texinfo @c -*-texinfo-*-
3 @setfilename ginac.info
4 @settitle GiNaC, an open framework for symbolic computation within the C++ programming language
11 @c I hate putting "@noindent" in front of every paragraph.
19 * ginac: (ginac). C++ library for symbolic computation.
23 This is a tutorial that documents GiNaC @value{VERSION}, an open
24 framework for symbolic computation within the C++ programming language.
26 Copyright (C) 1999-2004 Johannes Gutenberg University Mainz, Germany
28 Permission is granted to make and distribute verbatim copies of
29 this manual provided the copyright notice and this permission notice
30 are preserved on all copies.
33 Permission is granted to process this file through TeX and print the
34 results, provided the printed document carries copying permission
35 notice identical to this one except for the removal of this paragraph
38 Permission is granted to copy and distribute modified versions of this
39 manual under the conditions for verbatim copying, provided that the entire
40 resulting derived work is distributed under the terms of a permission
41 notice identical to this one.
45 @c finalout prevents ugly black rectangles on overfull hbox lines
47 @title GiNaC @value{VERSION}
48 @subtitle An open framework for symbolic computation within the C++ programming language
49 @subtitle @value{UPDATED}
50 @author The GiNaC Group:
51 @author Christian Bauer, Alexander Frink, Richard Kreckel, Jens Vollinga
54 @vskip 0pt plus 1filll
55 Copyright @copyright{} 1999-2004 Johannes Gutenberg University Mainz, Germany
57 Permission is granted to make and distribute verbatim copies of
58 this manual provided the copyright notice and this permission notice
59 are preserved on all copies.
61 Permission is granted to copy and distribute modified versions of this
62 manual under the conditions for verbatim copying, provided that the entire
63 resulting derived work is distributed under the terms of a permission
64 notice identical to this one.
73 @node Top, Introduction, (dir), (dir)
74 @c node-name, next, previous, up
77 This is a tutorial that documents GiNaC @value{VERSION}, an open
78 framework for symbolic computation within the C++ programming language.
81 * Introduction:: GiNaC's purpose.
82 * A Tour of GiNaC:: A quick tour of the library.
83 * Installation:: How to install the package.
84 * Basic Concepts:: Description of fundamental classes.
85 * Methods and Functions:: Algorithms for symbolic manipulations.
86 * Extending GiNaC:: How to extend the library.
87 * A Comparison With Other CAS:: Compares GiNaC to traditional CAS.
88 * Internal Structures:: Description of some internal structures.
89 * Package Tools:: Configuring packages to work with GiNaC.
95 @node Introduction, A Tour of GiNaC, Top, Top
96 @c node-name, next, previous, up
98 @cindex history of GiNaC
100 The motivation behind GiNaC derives from the observation that most
101 present day computer algebra systems (CAS) are linguistically and
102 semantically impoverished. Although they are quite powerful tools for
103 learning math and solving particular problems they lack modern
104 linguistic structures that allow for the creation of large-scale
105 projects. GiNaC is an attempt to overcome this situation by extending a
106 well established and standardized computer language (C++) by some
107 fundamental symbolic capabilities, thus allowing for integrated systems
108 that embed symbolic manipulations together with more established areas
109 of computer science (like computation-intense numeric applications,
110 graphical interfaces, etc.) under one roof.
112 The particular problem that led to the writing of the GiNaC framework is
113 still a very active field of research, namely the calculation of higher
114 order corrections to elementary particle interactions. There,
115 theoretical physicists are interested in matching present day theories
116 against experiments taking place at particle accelerators. The
117 computations involved are so complex they call for a combined symbolical
118 and numerical approach. This turned out to be quite difficult to
119 accomplish with the present day CAS we have worked with so far and so we
120 tried to fill the gap by writing GiNaC. But of course its applications
121 are in no way restricted to theoretical physics.
123 This tutorial is intended for the novice user who is new to GiNaC but
124 already has some background in C++ programming. However, since a
125 hand-made documentation like this one is difficult to keep in sync with
126 the development, the actual documentation is inside the sources in the
127 form of comments. That documentation may be parsed by one of the many
128 Javadoc-like documentation systems. If you fail at generating it you
129 may access it from @uref{http://www.ginac.de/reference/, the GiNaC home
130 page}. It is an invaluable resource not only for the advanced user who
131 wishes to extend the system (or chase bugs) but for everybody who wants
132 to comprehend the inner workings of GiNaC. This little tutorial on the
133 other hand only covers the basic things that are unlikely to change in
137 The GiNaC framework for symbolic computation within the C++ programming
138 language is Copyright @copyright{} 1999-2004 Johannes Gutenberg
139 University Mainz, Germany.
141 This program is free software; you can redistribute it and/or
142 modify it under the terms of the GNU General Public License as
143 published by the Free Software Foundation; either version 2 of the
144 License, or (at your option) any later version.
146 This program is distributed in the hope that it will be useful, but
147 WITHOUT ANY WARRANTY; without even the implied warranty of
148 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
149 General Public License for more details.
151 You should have received a copy of the GNU General Public License
152 along with this program; see the file COPYING. If not, write to the
153 Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston,
157 @node A Tour of GiNaC, How to use it from within C++, Introduction, Top
158 @c node-name, next, previous, up
159 @chapter A Tour of GiNaC
161 This quick tour of GiNaC wants to arise your interest in the
162 subsequent chapters by showing off a bit. Please excuse us if it
163 leaves many open questions.
166 * How to use it from within C++:: Two simple examples.
167 * What it can do for you:: A Tour of GiNaC's features.
171 @node How to use it from within C++, What it can do for you, A Tour of GiNaC, A Tour of GiNaC
172 @c node-name, next, previous, up
173 @section How to use it from within C++
175 The GiNaC open framework for symbolic computation within the C++ programming
176 language does not try to define a language of its own as conventional
177 CAS do. Instead, it extends the capabilities of C++ by symbolic
178 manipulations. Here is how to generate and print a simple (and rather
179 pointless) bivariate polynomial with some large coefficients:
183 #include <ginac/ginac.h>
185 using namespace GiNaC;
189 symbol x("x"), y("y");
192 for (int i=0; i<3; ++i)
193 poly += factorial(i+16)*pow(x,i)*pow(y,2-i);
195 cout << poly << endl;
200 Assuming the file is called @file{hello.cc}, on our system we can compile
201 and run it like this:
204 $ c++ hello.cc -o hello -lcln -lginac
206 355687428096000*x*y+20922789888000*y^2+6402373705728000*x^2
209 (@xref{Package Tools}, for tools that help you when creating a software
210 package that uses GiNaC.)
212 @cindex Hermite polynomial
213 Next, there is a more meaningful C++ program that calls a function which
214 generates Hermite polynomials in a specified free variable.
218 #include <ginac/ginac.h>
220 using namespace GiNaC;
222 ex HermitePoly(const symbol & x, int n)
224 ex HKer=exp(-pow(x, 2));
225 // uses the identity H_n(x) == (-1)^n exp(x^2) (d/dx)^n exp(-x^2)
226 return normal(pow(-1, n) * diff(HKer, x, n) / HKer);
233 for (int i=0; i<6; ++i)
234 cout << "H_" << i << "(z) == " << HermitePoly(z,i) << endl;
240 When run, this will type out
246 H_3(z) == -12*z+8*z^3
247 H_4(z) == -48*z^2+16*z^4+12
248 H_5(z) == 120*z-160*z^3+32*z^5
251 This method of generating the coefficients is of course far from optimal
252 for production purposes.
254 In order to show some more examples of what GiNaC can do we will now use
255 the @command{ginsh}, a simple GiNaC interactive shell that provides a
256 convenient window into GiNaC's capabilities.
259 @node What it can do for you, Installation, How to use it from within C++, A Tour of GiNaC
260 @c node-name, next, previous, up
261 @section What it can do for you
263 @cindex @command{ginsh}
264 After invoking @command{ginsh} one can test and experiment with GiNaC's
265 features much like in other Computer Algebra Systems except that it does
266 not provide programming constructs like loops or conditionals. For a
267 concise description of the @command{ginsh} syntax we refer to its
268 accompanied man page. Suffice to say that assignments and comparisons in
269 @command{ginsh} are written as they are in C, i.e. @code{=} assigns and
272 It can manipulate arbitrary precision integers in a very fast way.
273 Rational numbers are automatically converted to fractions of coprime
278 369988485035126972924700782451696644186473100389722973815184405301748249
280 123329495011708990974900260817232214728824366796574324605061468433916083
287 Exact numbers are always retained as exact numbers and only evaluated as
288 floating point numbers if requested. For instance, with numeric
289 radicals is dealt pretty much as with symbols. Products of sums of them
293 > expand((1+a^(1/5)-a^(2/5))^3);
294 1+3*a+3*a^(1/5)-5*a^(3/5)-a^(6/5)
295 > expand((1+3^(1/5)-3^(2/5))^3);
297 > evalf((1+3^(1/5)-3^(2/5))^3);
298 0.33408977534118624228
301 The function @code{evalf} that was used above converts any number in
302 GiNaC's expressions into floating point numbers. This can be done to
303 arbitrary predefined accuracy:
307 0.14285714285714285714
311 0.1428571428571428571428571428571428571428571428571428571428571428571428
312 5714285714285714285714285714285714285
315 Exact numbers other than rationals that can be manipulated in GiNaC
316 include predefined constants like Archimedes' @code{Pi}. They can both
317 be used in symbolic manipulations (as an exact number) as well as in
318 numeric expressions (as an inexact number):
324 9.869604401089358619+x
328 11.869604401089358619
331 Built-in functions evaluate immediately to exact numbers if
332 this is possible. Conversions that can be safely performed are done
333 immediately; conversions that are not generally valid are not done:
344 (Note that converting the last input to @code{x} would allow one to
345 conclude that @code{42*Pi} is equal to @code{0}.)
347 Linear equation systems can be solved along with basic linear
348 algebra manipulations over symbolic expressions. In C++ GiNaC offers
349 a matrix class for this purpose but we can see what it can do using
350 @command{ginsh}'s bracket notation to type them in:
353 > lsolve(a+x*y==z,x);
355 > lsolve(@{3*x+5*y == 7, -2*x+10*y == -5@}, @{x, y@});
357 > M = [ [1, 3], [-3, 2] ];
361 > charpoly(M,lambda);
363 > A = [ [1, 1], [2, -1] ];
366 [[1,1],[2,-1]]+2*[[1,3],[-3,2]]
369 > B = [ [0, 0, a], [b, 1, -b], [-1/a, 0, 0] ];
370 > evalm(B^(2^12345));
371 [[1,0,0],[0,1,0],[0,0,1]]
374 Multivariate polynomials and rational functions may be expanded,
375 collected and normalized (i.e. converted to a ratio of two coprime
379 > a = x^4 + 2*x^2*y^2 + 4*x^3*y + 12*x*y^3 - 3*y^4;
380 12*x*y^3+2*x^2*y^2+4*x^3*y-3*y^4+x^4
381 > b = x^2 + 4*x*y - y^2;
384 8*x^5*y+17*x^4*y^2+43*x^2*y^4-24*x*y^5+16*x^3*y^3+3*y^6+x^6
386 4*x^3*y-y^2-3*y^4+(12*y^3+4*y)*x+x^4+x^2*(1+2*y^2)
388 12*x*y^3-3*y^4+(-1+2*x^2)*y^2+(4*x+4*x^3)*y+x^2+x^4
393 You can differentiate functions and expand them as Taylor or Laurent
394 series in a very natural syntax (the second argument of @code{series} is
395 a relation defining the evaluation point, the third specifies the
398 @cindex Zeta function
402 > series(sin(x),x==0,4);
404 > series(1/tan(x),x==0,4);
405 x^(-1)-1/3*x+Order(x^2)
406 > series(tgamma(x),x==0,3);
407 x^(-1)-Euler+(1/12*Pi^2+1/2*Euler^2)*x+
408 (-1/3*zeta(3)-1/12*Pi^2*Euler-1/6*Euler^3)*x^2+Order(x^3)
410 x^(-1)-0.5772156649015328606+(0.9890559953279725555)*x
411 -(0.90747907608088628905)*x^2+Order(x^3)
412 > series(tgamma(2*sin(x)-2),x==Pi/2,6);
413 -(x-1/2*Pi)^(-2)+(-1/12*Pi^2-1/2*Euler^2-1/240)*(x-1/2*Pi)^2
414 -Euler-1/12+Order((x-1/2*Pi)^3)
417 Here we have made use of the @command{ginsh}-command @code{%} to pop the
418 previously evaluated element from @command{ginsh}'s internal stack.
420 If you ever wanted to convert units in C or C++ and found this is
421 cumbersome, here is the solution. Symbolic types can always be used as
422 tags for different types of objects. Converting from wrong units to the
423 metric system is now easy:
431 140613.91592783185568*kg*m^(-2)
435 @node Installation, Prerequisites, What it can do for you, Top
436 @c node-name, next, previous, up
437 @chapter Installation
440 GiNaC's installation follows the spirit of most GNU software. It is
441 easily installed on your system by three steps: configuration, build,
445 * Prerequisites:: Packages upon which GiNaC depends.
446 * Configuration:: How to configure GiNaC.
447 * Building GiNaC:: How to compile GiNaC.
448 * Installing GiNaC:: How to install GiNaC on your system.
452 @node Prerequisites, Configuration, Installation, Installation
453 @c node-name, next, previous, up
454 @section Prerequisites
456 In order to install GiNaC on your system, some prerequisites need to be
457 met. First of all, you need to have a C++-compiler adhering to the
458 ANSI-standard @cite{ISO/IEC 14882:1998(E)}. We used GCC for development
459 so if you have a different compiler you are on your own. For the
460 configuration to succeed you need a Posix compliant shell installed in
461 @file{/bin/sh}, GNU @command{bash} is fine. Perl is needed by the built
462 process as well, since some of the source files are automatically
463 generated by Perl scripts. Last but not least, Bruno Haible's library
464 CLN is extensively used and needs to be installed on your system.
465 Please get it either from @uref{ftp://ftp.santafe.edu/pub/gnu/}, from
466 @uref{ftp://ftpthep.physik.uni-mainz.de/pub/gnu/, GiNaC's FTP site} or
467 from @uref{ftp://ftp.ilog.fr/pub/Users/haible/gnu/, Bruno Haible's FTP
468 site} (it is covered by GPL) and install it prior to trying to install
469 GiNaC. The configure script checks if it can find it and if it cannot
470 it will refuse to continue.
473 @node Configuration, Building GiNaC, Prerequisites, Installation
474 @c node-name, next, previous, up
475 @section Configuration
476 @cindex configuration
479 To configure GiNaC means to prepare the source distribution for
480 building. It is done via a shell script called @command{configure} that
481 is shipped with the sources and was originally generated by GNU
482 Autoconf. Since a configure script generated by GNU Autoconf never
483 prompts, all customization must be done either via command line
484 parameters or environment variables. It accepts a list of parameters,
485 the complete set of which can be listed by calling it with the
486 @option{--help} option. The most important ones will be shortly
487 described in what follows:
492 @option{--disable-shared}: When given, this option switches off the
493 build of a shared library, i.e. a @file{.so} file. This may be convenient
494 when developing because it considerably speeds up compilation.
497 @option{--prefix=@var{PREFIX}}: The directory where the compiled library
498 and headers are installed. It defaults to @file{/usr/local} which means
499 that the library is installed in the directory @file{/usr/local/lib},
500 the header files in @file{/usr/local/include/ginac} and the documentation
501 (like this one) into @file{/usr/local/share/doc/GiNaC}.
504 @option{--libdir=@var{LIBDIR}}: Use this option in case you want to have
505 the library installed in some other directory than
506 @file{@var{PREFIX}/lib/}.
509 @option{--includedir=@var{INCLUDEDIR}}: Use this option in case you want
510 to have the header files installed in some other directory than
511 @file{@var{PREFIX}/include/ginac/}. For instance, if you specify
512 @option{--includedir=/usr/include} you will end up with the header files
513 sitting in the directory @file{/usr/include/ginac/}. Note that the
514 subdirectory @file{ginac} is enforced by this process in order to
515 keep the header files separated from others. This avoids some
516 clashes and allows for an easier deinstallation of GiNaC. This ought
517 to be considered A Good Thing (tm).
520 @option{--datadir=@var{DATADIR}}: This option may be given in case you
521 want to have the documentation installed in some other directory than
522 @file{@var{PREFIX}/share/doc/GiNaC/}.
526 In addition, you may specify some environment variables. @env{CXX}
527 holds the path and the name of the C++ compiler in case you want to
528 override the default in your path. (The @command{configure} script
529 searches your path for @command{c++}, @command{g++}, @command{gcc},
530 @command{CC}, @command{cxx} and @command{cc++} in that order.) It may
531 be very useful to define some compiler flags with the @env{CXXFLAGS}
532 environment variable, like optimization, debugging information and
533 warning levels. If omitted, it defaults to @option{-g
534 -O2}.@footnote{The @command{configure} script is itself generated from
535 the file @file{configure.ac}. It is only distributed in packaged
536 releases of GiNaC. If you got the naked sources, e.g. from CVS, you
537 must generate @command{configure} along with the various
538 @file{Makefile.in} by using the @command{autogen.sh} script. This will
539 require a fair amount of support from your local toolchain, though.}
541 The whole process is illustrated in the following two
542 examples. (Substitute @command{setenv @var{VARIABLE} @var{value}} for
543 @command{export @var{VARIABLE}=@var{value}} if the Berkeley C shell is
546 Here is a simple configuration for a site-wide GiNaC library assuming
547 everything is in default paths:
550 $ export CXXFLAGS="-Wall -O2"
554 And here is a configuration for a private static GiNaC library with
555 several components sitting in custom places (site-wide GCC and private
556 CLN). The compiler is persuaded to be picky and full assertions and
557 debugging information are switched on:
560 $ export CXX=/usr/local/gnu/bin/c++
561 $ export CPPFLAGS="$(CPPFLAGS) -I$(HOME)/include"
562 $ export CXXFLAGS="$(CXXFLAGS) -DDO_GINAC_ASSERT -ggdb -Wall -pedantic"
563 $ export LDFLAGS="$(LDFLAGS) -L$(HOME)/lib"
564 $ ./configure --disable-shared --prefix=$(HOME)
568 @node Building GiNaC, Installing GiNaC, Configuration, Installation
569 @c node-name, next, previous, up
570 @section Building GiNaC
571 @cindex building GiNaC
573 After proper configuration you should just build the whole
578 at the command prompt and go for a cup of coffee. The exact time it
579 takes to compile GiNaC depends not only on the speed of your machines
580 but also on other parameters, for instance what value for @env{CXXFLAGS}
581 you entered. Optimization may be very time-consuming.
583 Just to make sure GiNaC works properly you may run a collection of
584 regression tests by typing
590 This will compile some sample programs, run them and check the output
591 for correctness. The regression tests fall in three categories. First,
592 the so called @emph{exams} are performed, simple tests where some
593 predefined input is evaluated (like a pupils' exam). Second, the
594 @emph{checks} test the coherence of results among each other with
595 possible random input. Third, some @emph{timings} are performed, which
596 benchmark some predefined problems with different sizes and display the
597 CPU time used in seconds. Each individual test should return a message
598 @samp{passed}. This is mostly intended to be a QA-check if something
599 was broken during development, not a sanity check of your system. Some
600 of the tests in sections @emph{checks} and @emph{timings} may require
601 insane amounts of memory and CPU time. Feel free to kill them if your
602 machine catches fire. Another quite important intent is to allow people
603 to fiddle around with optimization.
605 By default, the only documentation that will be built is this tutorial
606 in @file{.info} format. To build the GiNaC tutorial and reference manual
607 in HTML, DVI, PostScript, or PDF formats, use one of
616 Generally, the top-level Makefile runs recursively to the
617 subdirectories. It is therefore safe to go into any subdirectory
618 (@code{doc/}, @code{ginsh/}, @dots{}) and simply type @code{make}
619 @var{target} there in case something went wrong.
622 @node Installing GiNaC, Basic Concepts, Building GiNaC, Installation
623 @c node-name, next, previous, up
624 @section Installing GiNaC
627 To install GiNaC on your system, simply type
633 As described in the section about configuration the files will be
634 installed in the following directories (the directories will be created
635 if they don't already exist):
640 @file{libginac.a} will go into @file{@var{PREFIX}/lib/} (or
641 @file{@var{LIBDIR}}) which defaults to @file{/usr/local/lib/}.
642 So will @file{libginac.so} unless the configure script was
643 given the option @option{--disable-shared}. The proper symlinks
644 will be established as well.
647 All the header files will be installed into @file{@var{PREFIX}/include/ginac/}
648 (or @file{@var{INCLUDEDIR}/ginac/}, if specified).
651 All documentation (info) will be stuffed into
652 @file{@var{PREFIX}/share/doc/GiNaC/} (or
653 @file{@var{DATADIR}/doc/GiNaC/}, if @var{DATADIR} was specified).
657 For the sake of completeness we will list some other useful make
658 targets: @command{make clean} deletes all files generated by
659 @command{make}, i.e. all the object files. In addition @command{make
660 distclean} removes all files generated by the configuration and
661 @command{make maintainer-clean} goes one step further and deletes files
662 that may require special tools to rebuild (like the @command{libtool}
663 for instance). Finally @command{make uninstall} removes the installed
664 library, header files and documentation@footnote{Uninstallation does not
665 work after you have called @command{make distclean} since the
666 @file{Makefile} is itself generated by the configuration from
667 @file{Makefile.in} and hence deleted by @command{make distclean}. There
668 are two obvious ways out of this dilemma. First, you can run the
669 configuration again with the same @var{PREFIX} thus creating a
670 @file{Makefile} with a working @samp{uninstall} target. Second, you can
671 do it by hand since you now know where all the files went during
675 @node Basic Concepts, Expressions, Installing GiNaC, Top
676 @c node-name, next, previous, up
677 @chapter Basic Concepts
679 This chapter will describe the different fundamental objects that can be
680 handled by GiNaC. But before doing so, it is worthwhile introducing you
681 to the more commonly used class of expressions, representing a flexible
682 meta-class for storing all mathematical objects.
685 * Expressions:: The fundamental GiNaC class.
686 * Automatic evaluation:: Evaluation and canonicalization.
687 * Error handling:: How the library reports errors.
688 * The Class Hierarchy:: Overview of GiNaC's classes.
689 * Symbols:: Symbolic objects.
690 * Numbers:: Numerical objects.
691 * Constants:: Pre-defined constants.
692 * Fundamental containers:: Sums, products and powers.
693 * Lists:: Lists of expressions.
694 * Mathematical functions:: Mathematical functions.
695 * Relations:: Equality, Inequality and all that.
696 * Matrices:: Matrices.
697 * Indexed objects:: Handling indexed quantities.
698 * Non-commutative objects:: Algebras with non-commutative products.
699 * Hash Maps:: A faster alternative to std::map<>.
703 @node Expressions, Automatic evaluation, Basic Concepts, Basic Concepts
704 @c node-name, next, previous, up
706 @cindex expression (class @code{ex})
709 The most common class of objects a user deals with is the expression
710 @code{ex}, representing a mathematical object like a variable, number,
711 function, sum, product, etc@dots{} Expressions may be put together to form
712 new expressions, passed as arguments to functions, and so on. Here is a
713 little collection of valid expressions:
716 ex MyEx1 = 5; // simple number
717 ex MyEx2 = x + 2*y; // polynomial in x and y
718 ex MyEx3 = (x + 1)/(x - 1); // rational expression
719 ex MyEx4 = sin(x + 2*y) + 3*z + 41; // containing a function
720 ex MyEx5 = MyEx4 + 1; // similar to above
723 Expressions are handles to other more fundamental objects, that often
724 contain other expressions thus creating a tree of expressions
725 (@xref{Internal Structures}, for particular examples). Most methods on
726 @code{ex} therefore run top-down through such an expression tree. For
727 example, the method @code{has()} scans recursively for occurrences of
728 something inside an expression. Thus, if you have declared @code{MyEx4}
729 as in the example above @code{MyEx4.has(y)} will find @code{y} inside
730 the argument of @code{sin} and hence return @code{true}.
732 The next sections will outline the general picture of GiNaC's class
733 hierarchy and describe the classes of objects that are handled by
736 @subsection Note: Expressions and STL containers
738 GiNaC expressions (@code{ex} objects) have value semantics (they can be
739 assigned, reassigned and copied like integral types) but the operator
740 @code{<} doesn't provide a well-defined ordering on them. In STL-speak,
741 expressions are @samp{Assignable} but not @samp{LessThanComparable}.
743 This implies that in order to use expressions in sorted containers such as
744 @code{std::map<>} and @code{std::set<>} you have to supply a suitable
745 comparison predicate. GiNaC provides such a predicate, called
746 @code{ex_is_less}. For example, a set of expressions should be defined
747 as @code{std::set<ex, ex_is_less>}.
749 Unsorted containers such as @code{std::vector<>} and @code{std::list<>}
750 don't pose a problem. A @code{std::vector<ex>} works as expected.
752 @xref{Information About Expressions}, for more about comparing and ordering
756 @node Automatic evaluation, Error handling, Expressions, Basic Concepts
757 @c node-name, next, previous, up
758 @section Automatic evaluation and canonicalization of expressions
761 GiNaC performs some automatic transformations on expressions, to simplify
762 them and put them into a canonical form. Some examples:
765 ex MyEx1 = 2*x - 1 + x; // 3*x-1
766 ex MyEx2 = x - x; // 0
767 ex MyEx3 = cos(2*Pi); // 1
768 ex MyEx4 = x*y/x; // y
771 This behavior is usually referred to as @dfn{automatic} or @dfn{anonymous
772 evaluation}. GiNaC only performs transformations that are
776 at most of complexity
784 algebraically correct, possibly except for a set of measure zero (e.g.
785 @math{x/x} is transformed to @math{1} although this is incorrect for @math{x=0})
788 There are two types of automatic transformations in GiNaC that may not
789 behave in an entirely obvious way at first glance:
793 The terms of sums and products (and some other things like the arguments of
794 symmetric functions, the indices of symmetric tensors etc.) are re-ordered
795 into a canonical form that is deterministic, but not lexicographical or in
796 any other way easy to guess (it almost always depends on the number and
797 order of the symbols you define). However, constructing the same expression
798 twice, either implicitly or explicitly, will always result in the same
801 Expressions of the form 'number times sum' are automatically expanded (this
802 has to do with GiNaC's internal representation of sums and products). For
805 ex MyEx5 = 2*(x + y); // 2*x+2*y
806 ex MyEx6 = z*(x + y); // z*(x+y)
810 The general rule is that when you construct expressions, GiNaC automatically
811 creates them in canonical form, which might differ from the form you typed in
812 your program. This may create some awkward looking output (@samp{-y+x} instead
813 of @samp{x-y}) but allows for more efficient operation and usually yields
814 some immediate simplifications.
816 @cindex @code{eval()}
817 Internally, the anonymous evaluator in GiNaC is implemented by the methods
820 ex ex::eval(int level = 0) const;
821 ex basic::eval(int level = 0) const;
824 but unless you are extending GiNaC with your own classes or functions, there
825 should never be any reason to call them explicitly. All GiNaC methods that
826 transform expressions, like @code{subs()} or @code{normal()}, automatically
827 re-evaluate their results.
830 @node Error handling, The Class Hierarchy, Automatic evaluation, Basic Concepts
831 @c node-name, next, previous, up
832 @section Error handling
834 @cindex @code{pole_error} (class)
836 GiNaC reports run-time errors by throwing C++ exceptions. All exceptions
837 generated by GiNaC are subclassed from the standard @code{exception} class
838 defined in the @file{<stdexcept>} header. In addition to the predefined
839 @code{logic_error}, @code{domain_error}, @code{out_of_range},
840 @code{invalid_argument}, @code{runtime_error}, @code{range_error} and
841 @code{overflow_error} types, GiNaC also defines a @code{pole_error}
842 exception that gets thrown when trying to evaluate a mathematical function
845 The @code{pole_error} class has a member function
848 int pole_error::degree() const;
851 that returns the order of the singularity (or 0 when the pole is
852 logarithmic or the order is undefined).
854 When using GiNaC it is useful to arrange for exceptions to be caught in
855 the main program even if you don't want to do any special error handling.
856 Otherwise whenever an error occurs in GiNaC, it will be delegated to the
857 default exception handler of your C++ compiler's run-time system which
858 usually only aborts the program without giving any information what went
861 Here is an example for a @code{main()} function that catches and prints
862 exceptions generated by GiNaC:
867 #include <ginac/ginac.h>
869 using namespace GiNaC;
877 @} catch (exception &p) @{
878 cerr << p.what() << endl;
886 @node The Class Hierarchy, Symbols, Error handling, Basic Concepts
887 @c node-name, next, previous, up
888 @section The Class Hierarchy
890 GiNaC's class hierarchy consists of several classes representing
891 mathematical objects, all of which (except for @code{ex} and some
892 helpers) are internally derived from one abstract base class called
893 @code{basic}. You do not have to deal with objects of class
894 @code{basic}, instead you'll be dealing with symbols, numbers,
895 containers of expressions and so on.
899 To get an idea about what kinds of symbolic composites may be built we
900 have a look at the most important classes in the class hierarchy and
901 some of the relations among the classes:
903 @image{classhierarchy}
905 The abstract classes shown here (the ones without drop-shadow) are of no
906 interest for the user. They are used internally in order to avoid code
907 duplication if two or more classes derived from them share certain
908 features. An example is @code{expairseq}, a container for a sequence of
909 pairs each consisting of one expression and a number (@code{numeric}).
910 What @emph{is} visible to the user are the derived classes @code{add}
911 and @code{mul}, representing sums and products. @xref{Internal
912 Structures}, where these two classes are described in more detail. The
913 following table shortly summarizes what kinds of mathematical objects
914 are stored in the different classes:
917 @multitable @columnfractions .22 .78
918 @item @code{symbol} @tab Algebraic symbols @math{a}, @math{x}, @math{y}@dots{}
919 @item @code{constant} @tab Constants like
926 @item @code{numeric} @tab All kinds of numbers, @math{42}, @math{7/3*I}, @math{3.14159}@dots{}
927 @item @code{add} @tab Sums like @math{x+y} or @math{a-(2*b)+3}
928 @item @code{mul} @tab Products like @math{x*y} or @math{2*a^2*(x+y+z)/b}
929 @item @code{ncmul} @tab Products of non-commutative objects
930 @item @code{power} @tab Exponentials such as @math{x^2}, @math{a^b},
935 @code{sqrt(}@math{2}@code{)}
938 @item @code{pseries} @tab Power Series, e.g. @math{x-1/6*x^3+1/120*x^5+O(x^7)}
939 @item @code{function} @tab A symbolic function like
946 @item @code{lst} @tab Lists of expressions @{@math{x}, @math{2*y}, @math{3+z}@}
947 @item @code{matrix} @tab @math{m}x@math{n} matrices of expressions
948 @item @code{relational} @tab A relation like the identity @math{x}@code{==}@math{y}
949 @item @code{indexed} @tab Indexed object like @math{A_ij}
950 @item @code{tensor} @tab Special tensor like the delta and metric tensors
951 @item @code{idx} @tab Index of an indexed object
952 @item @code{varidx} @tab Index with variance
953 @item @code{spinidx} @tab Index with variance and dot (used in Weyl-van-der-Waerden spinor formalism)
954 @item @code{wildcard} @tab Wildcard for pattern matching
955 @item @code{structure} @tab Template for user-defined classes
960 @node Symbols, Numbers, The Class Hierarchy, Basic Concepts
961 @c node-name, next, previous, up
963 @cindex @code{symbol} (class)
964 @cindex hierarchy of classes
967 Symbolic indeterminates, or @dfn{symbols} for short, are for symbolic
968 manipulation what atoms are for chemistry.
970 A typical symbol definition looks like this:
975 This definition actually contains three very different things:
977 @item a C++ variable named @code{x}
978 @item a @code{symbol} object stored in this C++ variable; this object
979 represents the symbol in a GiNaC expression
980 @item the string @code{"x"} which is the name of the symbol, used (almost)
981 exclusively for printing expressions holding the symbol
984 Symbols have an explicit name, supplied as a string during construction,
985 because in C++, variable names can't be used as values, and the C++ compiler
986 throws them away during compilation.
988 It is possible to omit the symbol name in the definition:
993 In this case, GiNaC will assign the symbol an internal, unique name of the
994 form @code{symbolNNN}. This won't affect the usability of the symbol but
995 the output of your calculations will become more readable if you give your
996 symbols sensible names (for intermediate expressions that are only used
997 internally such anonymous symbols can be quite useful, however).
999 Now, here is one important property of GiNaC that differentiates it from
1000 other computer algebra programs you may have used: GiNaC does @emph{not} use
1001 the names of symbols to tell them apart, but a (hidden) serial number that
1002 is unique for each newly created @code{symbol} object. In you want to use
1003 one and the same symbol in different places in your program, you must only
1004 create one @code{symbol} object and pass that around. If you create another
1005 symbol, even if it has the same name, GiNaC will treat it as a different
1022 // prints "x^6" which looks right, but...
1024 cout << e.degree(x) << endl;
1025 // ...this doesn't work. The symbol "x" here is different from the one
1026 // in f() and in the expression returned by f(). Consequently, it
1031 One possibility to ensure that @code{f()} and @code{main()} use the same
1032 symbol is to pass the symbol as an argument to @code{f()}:
1034 ex f(int n, const ex & x)
1043 // Now, f() uses the same symbol.
1046 cout << e.degree(x) << endl;
1047 // prints "6", as expected
1051 Another possibility would be to define a global symbol @code{x} that is used
1052 by both @code{f()} and @code{main()}. If you are using global symbols and
1053 multiple compilation units you must take special care, however. Suppose
1054 that you have a header file @file{globals.h} in your program that defines
1055 a @code{symbol x("x");}. In this case, every unit that includes
1056 @file{globals.h} would also get its own definition of @code{x} (because
1057 header files are just inlined into the source code by the C++ preprocessor),
1058 and hence you would again end up with multiple equally-named, but different,
1059 symbols. Instead, the @file{globals.h} header should only contain a
1060 @emph{declaration} like @code{extern symbol x;}, with the definition of
1061 @code{x} moved into a C++ source file such as @file{globals.cpp}.
1063 A different approach to ensuring that symbols used in different parts of
1064 your program are identical is to create them with a @emph{factory} function
1067 const symbol & get_symbol(const string & s)
1069 static map<string, symbol> directory;
1070 map<string, symbol>::iterator i = directory.find(s);
1071 if (i != directory.end())
1074 return directory.insert(make_pair(s, symbol(s))).first->second;
1078 This function returns one newly constructed symbol for each name that is
1079 passed in, and it returns the same symbol when called multiple times with
1080 the same name. Using this symbol factory, we can rewrite our example like
1085 return pow(get_symbol("x"), n);
1092 // Both calls of get_symbol("x") yield the same symbol.
1093 cout << e.degree(get_symbol("x")) << endl;
1098 Instead of creating symbols from strings we could also have
1099 @code{get_symbol()} take, for example, an integer number as its argument.
1100 In this case, we would probably want to give the generated symbols names
1101 that include this number, which can be accomplished with the help of an
1102 @code{ostringstream}.
1104 In general, if you're getting weird results from GiNaC such as an expression
1105 @samp{x-x} that is not simplified to zero, you should check your symbol
1108 As we said, the names of symbols primarily serve for purposes of expression
1109 output. But there are actually two instances where GiNaC uses the names for
1110 identifying symbols: When constructing an expression from a string, and when
1111 recreating an expression from an archive (@pxref{Input/Output}).
1113 In addition to its name, a symbol may contain a special string that is used
1116 symbol x("x", "\\Box");
1119 This creates a symbol that is printed as "@code{x}" in normal output, but
1120 as "@code{\Box}" in LaTeX code (@xref{Input/Output}, for more
1121 information about the different output formats of expressions in GiNaC).
1122 GiNaC automatically creates proper LaTeX code for symbols having names of
1123 greek letters (@samp{alpha}, @samp{mu}, etc.).
1125 @cindex @code{subs()}
1126 Symbols in GiNaC can't be assigned values. If you need to store results of
1127 calculations and give them a name, use C++ variables of type @code{ex}.
1128 If you want to replace a symbol in an expression with something else, you
1129 can invoke the expression's @code{.subs()} method
1130 (@pxref{Substituting Expressions}).
1132 @cindex @code{realsymbol()}
1133 By default, symbols are expected to stand in for complex values, i.e. they live
1134 in the complex domain. As a consequence, operations like complex conjugation,
1135 for example (@pxref{Complex Conjugation}), do @emph{not} evaluate if applied
1136 to such symbols. Likewise @code{log(exp(x))} does not evaluate to @code{x},
1137 because of the unknown imaginary part of @code{x}.
1138 On the other hand, if you are sure that your symbols will hold only real values, you
1139 would like to have such functions evaluated. Therefore GiNaC allows you to specify
1140 the domain of the symbol. Instead of @code{symbol x("x");} you can write
1141 @code{realsymbol x("x");} to tell GiNaC that @code{x} stands in for real values.
1144 @node Numbers, Constants, Symbols, Basic Concepts
1145 @c node-name, next, previous, up
1147 @cindex @code{numeric} (class)
1153 For storing numerical things, GiNaC uses Bruno Haible's library CLN.
1154 The classes therein serve as foundation classes for GiNaC. CLN stands
1155 for Class Library for Numbers or alternatively for Common Lisp Numbers.
1156 In order to find out more about CLN's internals, the reader is referred to
1157 the documentation of that library. @inforef{Introduction, , cln}, for
1158 more information. Suffice to say that it is by itself build on top of
1159 another library, the GNU Multiple Precision library GMP, which is an
1160 extremely fast library for arbitrary long integers and rationals as well
1161 as arbitrary precision floating point numbers. It is very commonly used
1162 by several popular cryptographic applications. CLN extends GMP by
1163 several useful things: First, it introduces the complex number field
1164 over either reals (i.e. floating point numbers with arbitrary precision)
1165 or rationals. Second, it automatically converts rationals to integers
1166 if the denominator is unity and complex numbers to real numbers if the
1167 imaginary part vanishes and also correctly treats algebraic functions.
1168 Third it provides good implementations of state-of-the-art algorithms
1169 for all trigonometric and hyperbolic functions as well as for
1170 calculation of some useful constants.
1172 The user can construct an object of class @code{numeric} in several
1173 ways. The following example shows the four most important constructors.
1174 It uses construction from C-integer, construction of fractions from two
1175 integers, construction from C-float and construction from a string:
1179 #include <ginac/ginac.h>
1180 using namespace GiNaC;
1184 numeric two = 2; // exact integer 2
1185 numeric r(2,3); // exact fraction 2/3
1186 numeric e(2.71828); // floating point number
1187 numeric p = "3.14159265358979323846"; // constructor from string
1188 // Trott's constant in scientific notation:
1189 numeric trott("1.0841015122311136151E-2");
1191 std::cout << two*p << std::endl; // floating point 6.283...
1196 @cindex complex numbers
1197 The imaginary unit in GiNaC is a predefined @code{numeric} object with the
1202 numeric z1 = 2-3*I; // exact complex number 2-3i
1203 numeric z2 = 5.9+1.6*I; // complex floating point number
1207 It may be tempting to construct fractions by writing @code{numeric r(3/2)}.
1208 This would, however, call C's built-in operator @code{/} for integers
1209 first and result in a numeric holding a plain integer 1. @strong{Never
1210 use the operator @code{/} on integers} unless you know exactly what you
1211 are doing! Use the constructor from two integers instead, as shown in
1212 the example above. Writing @code{numeric(1)/2} may look funny but works
1215 @cindex @code{Digits}
1217 We have seen now the distinction between exact numbers and floating
1218 point numbers. Clearly, the user should never have to worry about
1219 dynamically created exact numbers, since their `exactness' always
1220 determines how they ought to be handled, i.e. how `long' they are. The
1221 situation is different for floating point numbers. Their accuracy is
1222 controlled by one @emph{global} variable, called @code{Digits}. (For
1223 those readers who know about Maple: it behaves very much like Maple's
1224 @code{Digits}). All objects of class numeric that are constructed from
1225 then on will be stored with a precision matching that number of decimal
1230 #include <ginac/ginac.h>
1231 using namespace std;
1232 using namespace GiNaC;
1236 numeric three(3.0), one(1.0);
1237 numeric x = one/three;
1239 cout << "in " << Digits << " digits:" << endl;
1241 cout << Pi.evalf() << endl;
1253 The above example prints the following output to screen:
1257 0.33333333333333333334
1258 3.1415926535897932385
1260 0.33333333333333333333333333333333333333333333333333333333333333333334
1261 3.1415926535897932384626433832795028841971693993751058209749445923078
1265 Note that the last number is not necessarily rounded as you would
1266 naively expect it to be rounded in the decimal system. But note also,
1267 that in both cases you got a couple of extra digits. This is because
1268 numbers are internally stored by CLN as chunks of binary digits in order
1269 to match your machine's word size and to not waste precision. Thus, on
1270 architectures with different word size, the above output might even
1271 differ with regard to actually computed digits.
1273 It should be clear that objects of class @code{numeric} should be used
1274 for constructing numbers or for doing arithmetic with them. The objects
1275 one deals with most of the time are the polymorphic expressions @code{ex}.
1277 @subsection Tests on numbers
1279 Once you have declared some numbers, assigned them to expressions and
1280 done some arithmetic with them it is frequently desired to retrieve some
1281 kind of information from them like asking whether that number is
1282 integer, rational, real or complex. For those cases GiNaC provides
1283 several useful methods. (Internally, they fall back to invocations of
1284 certain CLN functions.)
1286 As an example, let's construct some rational number, multiply it with
1287 some multiple of its denominator and test what comes out:
1291 #include <ginac/ginac.h>
1292 using namespace std;
1293 using namespace GiNaC;
1295 // some very important constants:
1296 const numeric twentyone(21);
1297 const numeric ten(10);
1298 const numeric five(5);
1302 numeric answer = twentyone;
1305 cout << answer.is_integer() << endl; // false, it's 21/5
1307 cout << answer.is_integer() << endl; // true, it's 42 now!
1311 Note that the variable @code{answer} is constructed here as an integer
1312 by @code{numeric}'s copy constructor but in an intermediate step it
1313 holds a rational number represented as integer numerator and integer
1314 denominator. When multiplied by 10, the denominator becomes unity and
1315 the result is automatically converted to a pure integer again.
1316 Internally, the underlying CLN is responsible for this behavior and we
1317 refer the reader to CLN's documentation. Suffice to say that
1318 the same behavior applies to complex numbers as well as return values of
1319 certain functions. Complex numbers are automatically converted to real
1320 numbers if the imaginary part becomes zero. The full set of tests that
1321 can be applied is listed in the following table.
1324 @multitable @columnfractions .30 .70
1325 @item @strong{Method} @tab @strong{Returns true if the object is@dots{}}
1326 @item @code{.is_zero()}
1327 @tab @dots{}equal to zero
1328 @item @code{.is_positive()}
1329 @tab @dots{}not complex and greater than 0
1330 @item @code{.is_integer()}
1331 @tab @dots{}a (non-complex) integer
1332 @item @code{.is_pos_integer()}
1333 @tab @dots{}an integer and greater than 0
1334 @item @code{.is_nonneg_integer()}
1335 @tab @dots{}an integer and greater equal 0
1336 @item @code{.is_even()}
1337 @tab @dots{}an even integer
1338 @item @code{.is_odd()}
1339 @tab @dots{}an odd integer
1340 @item @code{.is_prime()}
1341 @tab @dots{}a prime integer (probabilistic primality test)
1342 @item @code{.is_rational()}
1343 @tab @dots{}an exact rational number (integers are rational, too)
1344 @item @code{.is_real()}
1345 @tab @dots{}a real integer, rational or float (i.e. is not complex)
1346 @item @code{.is_cinteger()}
1347 @tab @dots{}a (complex) integer (such as @math{2-3*I})
1348 @item @code{.is_crational()}
1349 @tab @dots{}an exact (complex) rational number (such as @math{2/3+7/2*I})
1353 @subsection Numeric functions
1355 The following functions can be applied to @code{numeric} objects and will be
1356 evaluated immediately:
1359 @multitable @columnfractions .30 .70
1360 @item @strong{Name} @tab @strong{Function}
1361 @item @code{inverse(z)}
1362 @tab returns @math{1/z}
1363 @cindex @code{inverse()} (numeric)
1364 @item @code{pow(a, b)}
1365 @tab exponentiation @math{a^b}
1368 @item @code{real(z)}
1370 @cindex @code{real()}
1371 @item @code{imag(z)}
1373 @cindex @code{imag()}
1374 @item @code{csgn(z)}
1375 @tab complex sign (returns an @code{int})
1376 @item @code{numer(z)}
1377 @tab numerator of rational or complex rational number
1378 @item @code{denom(z)}
1379 @tab denominator of rational or complex rational number
1380 @item @code{sqrt(z)}
1382 @item @code{isqrt(n)}
1383 @tab integer square root
1384 @cindex @code{isqrt()}
1391 @item @code{asin(z)}
1393 @item @code{acos(z)}
1395 @item @code{atan(z)}
1396 @tab inverse tangent
1397 @item @code{atan(y, x)}
1398 @tab inverse tangent with two arguments
1399 @item @code{sinh(z)}
1400 @tab hyperbolic sine
1401 @item @code{cosh(z)}
1402 @tab hyperbolic cosine
1403 @item @code{tanh(z)}
1404 @tab hyperbolic tangent
1405 @item @code{asinh(z)}
1406 @tab inverse hyperbolic sine
1407 @item @code{acosh(z)}
1408 @tab inverse hyperbolic cosine
1409 @item @code{atanh(z)}
1410 @tab inverse hyperbolic tangent
1412 @tab exponential function
1414 @tab natural logarithm
1417 @item @code{zeta(z)}
1418 @tab Riemann's zeta function
1419 @item @code{tgamma(z)}
1421 @item @code{lgamma(z)}
1422 @tab logarithm of gamma function
1424 @tab psi (digamma) function
1425 @item @code{psi(n, z)}
1426 @tab derivatives of psi function (polygamma functions)
1427 @item @code{factorial(n)}
1428 @tab factorial function @math{n!}
1429 @item @code{doublefactorial(n)}
1430 @tab double factorial function @math{n!!}
1431 @cindex @code{doublefactorial()}
1432 @item @code{binomial(n, k)}
1433 @tab binomial coefficients
1434 @item @code{bernoulli(n)}
1435 @tab Bernoulli numbers
1436 @cindex @code{bernoulli()}
1437 @item @code{fibonacci(n)}
1438 @tab Fibonacci numbers
1439 @cindex @code{fibonacci()}
1440 @item @code{mod(a, b)}
1441 @tab modulus in positive representation (in the range @code{[0, abs(b)-1]} with the sign of b, or zero)
1442 @cindex @code{mod()}
1443 @item @code{smod(a, b)}
1444 @tab modulus in symmetric representation (in the range @code{[-iquo(abs(b)-1, 2), iquo(abs(b), 2)]})
1445 @cindex @code{smod()}
1446 @item @code{irem(a, b)}
1447 @tab integer remainder (has the sign of @math{a}, or is zero)
1448 @cindex @code{irem()}
1449 @item @code{irem(a, b, q)}
1450 @tab integer remainder and quotient, @code{irem(a, b, q) == a-q*b}
1451 @item @code{iquo(a, b)}
1452 @tab integer quotient
1453 @cindex @code{iquo()}
1454 @item @code{iquo(a, b, r)}
1455 @tab integer quotient and remainder, @code{r == a-iquo(a, b)*b}
1456 @item @code{gcd(a, b)}
1457 @tab greatest common divisor
1458 @item @code{lcm(a, b)}
1459 @tab least common multiple
1463 Most of these functions are also available as symbolic functions that can be
1464 used in expressions (@pxref{Mathematical functions}) or, like @code{gcd()},
1465 as polynomial algorithms.
1467 @subsection Converting numbers
1469 Sometimes it is desirable to convert a @code{numeric} object back to a
1470 built-in arithmetic type (@code{int}, @code{double}, etc.). The @code{numeric}
1471 class provides a couple of methods for this purpose:
1473 @cindex @code{to_int()}
1474 @cindex @code{to_long()}
1475 @cindex @code{to_double()}
1476 @cindex @code{to_cl_N()}
1478 int numeric::to_int() const;
1479 long numeric::to_long() const;
1480 double numeric::to_double() const;
1481 cln::cl_N numeric::to_cl_N() const;
1484 @code{to_int()} and @code{to_long()} only work when the number they are
1485 applied on is an exact integer. Otherwise the program will halt with a
1486 message like @samp{Not a 32-bit integer}. @code{to_double()} applied on a
1487 rational number will return a floating-point approximation. Both
1488 @code{to_int()/to_long()} and @code{to_double()} discard the imaginary
1489 part of complex numbers.
1492 @node Constants, Fundamental containers, Numbers, Basic Concepts
1493 @c node-name, next, previous, up
1495 @cindex @code{constant} (class)
1498 @cindex @code{Catalan}
1499 @cindex @code{Euler}
1500 @cindex @code{evalf()}
1501 Constants behave pretty much like symbols except that they return some
1502 specific number when the method @code{.evalf()} is called.
1504 The predefined known constants are:
1507 @multitable @columnfractions .14 .30 .56
1508 @item @strong{Name} @tab @strong{Common Name} @tab @strong{Numerical Value (to 35 digits)}
1510 @tab Archimedes' constant
1511 @tab 3.14159265358979323846264338327950288
1512 @item @code{Catalan}
1513 @tab Catalan's constant
1514 @tab 0.91596559417721901505460351493238411
1516 @tab Euler's (or Euler-Mascheroni) constant
1517 @tab 0.57721566490153286060651209008240243
1522 @node Fundamental containers, Lists, Constants, Basic Concepts
1523 @c node-name, next, previous, up
1524 @section Sums, products and powers
1528 @cindex @code{power}
1530 Simple rational expressions are written down in GiNaC pretty much like
1531 in other CAS or like expressions involving numerical variables in C.
1532 The necessary operators @code{+}, @code{-}, @code{*} and @code{/} have
1533 been overloaded to achieve this goal. When you run the following
1534 code snippet, the constructor for an object of type @code{mul} is
1535 automatically called to hold the product of @code{a} and @code{b} and
1536 then the constructor for an object of type @code{add} is called to hold
1537 the sum of that @code{mul} object and the number one:
1541 symbol a("a"), b("b");
1546 @cindex @code{pow()}
1547 For exponentiation, you have already seen the somewhat clumsy (though C-ish)
1548 statement @code{pow(x,2);} to represent @code{x} squared. This direct
1549 construction is necessary since we cannot safely overload the constructor
1550 @code{^} in C++ to construct a @code{power} object. If we did, it would
1551 have several counterintuitive and undesired effects:
1555 Due to C's operator precedence, @code{2*x^2} would be parsed as @code{(2*x)^2}.
1557 Due to the binding of the operator @code{^}, @code{x^a^b} would result in
1558 @code{(x^a)^b}. This would be confusing since most (though not all) other CAS
1559 interpret this as @code{x^(a^b)}.
1561 Also, expressions involving integer exponents are very frequently used,
1562 which makes it even more dangerous to overload @code{^} since it is then
1563 hard to distinguish between the semantics as exponentiation and the one
1564 for exclusive or. (It would be embarrassing to return @code{1} where one
1565 has requested @code{2^3}.)
1568 @cindex @command{ginsh}
1569 All effects are contrary to mathematical notation and differ from the
1570 way most other CAS handle exponentiation, therefore overloading @code{^}
1571 is ruled out for GiNaC's C++ part. The situation is different in
1572 @command{ginsh}, there the exponentiation-@code{^} exists. (Also note
1573 that the other frequently used exponentiation operator @code{**} does
1574 not exist at all in C++).
1576 To be somewhat more precise, objects of the three classes described
1577 here, are all containers for other expressions. An object of class
1578 @code{power} is best viewed as a container with two slots, one for the
1579 basis, one for the exponent. All valid GiNaC expressions can be
1580 inserted. However, basic transformations like simplifying
1581 @code{pow(pow(x,2),3)} to @code{x^6} automatically are only performed
1582 when this is mathematically possible. If we replace the outer exponent
1583 three in the example by some symbols @code{a}, the simplification is not
1584 safe and will not be performed, since @code{a} might be @code{1/2} and
1587 Objects of type @code{add} and @code{mul} are containers with an
1588 arbitrary number of slots for expressions to be inserted. Again, simple
1589 and safe simplifications are carried out like transforming
1590 @code{3*x+4-x} to @code{2*x+4}.
1593 @node Lists, Mathematical functions, Fundamental containers, Basic Concepts
1594 @c node-name, next, previous, up
1595 @section Lists of expressions
1596 @cindex @code{lst} (class)
1598 @cindex @code{nops()}
1600 @cindex @code{append()}
1601 @cindex @code{prepend()}
1602 @cindex @code{remove_first()}
1603 @cindex @code{remove_last()}
1604 @cindex @code{remove_all()}
1606 The GiNaC class @code{lst} serves for holding a @dfn{list} of arbitrary
1607 expressions. They are not as ubiquitous as in many other computer algebra
1608 packages, but are sometimes used to supply a variable number of arguments of
1609 the same type to GiNaC methods such as @code{subs()} and some @code{matrix}
1610 constructors, so you should have a basic understanding of them.
1612 Lists can be constructed by assigning a comma-separated sequence of
1617 symbol x("x"), y("y");
1620 // now, l is a list holding the expressions 'x', '2', 'y', and 'x+y',
1625 There are also constructors that allow direct creation of lists of up to
1626 16 expressions, which is often more convenient but slightly less efficient:
1630 // This produces the same list 'l' as above:
1631 // lst l(x, 2, y, x+y);
1632 // lst l = lst(x, 2, y, x+y);
1636 Use the @code{nops()} method to determine the size (number of expressions) of
1637 a list and the @code{op()} method or the @code{[]} operator to access
1638 individual elements:
1642 cout << l.nops() << endl; // prints '4'
1643 cout << l.op(2) << " " << l[0] << endl; // prints 'y x'
1647 As with the standard @code{list<T>} container, accessing random elements of a
1648 @code{lst} is generally an operation of order @math{O(N)}. Faster read-only
1649 sequential access to the elements of a list is possible with the
1650 iterator types provided by the @code{lst} class:
1653 typedef ... lst::const_iterator;
1654 typedef ... lst::const_reverse_iterator;
1655 lst::const_iterator lst::begin() const;
1656 lst::const_iterator lst::end() const;
1657 lst::const_reverse_iterator lst::rbegin() const;
1658 lst::const_reverse_iterator lst::rend() const;
1661 For example, to print the elements of a list individually you can use:
1666 for (lst::const_iterator i = l.begin(); i != l.end(); ++i)
1671 which is one order faster than
1676 for (size_t i = 0; i < l.nops(); ++i)
1677 cout << l.op(i) << endl;
1681 These iterators also allow you to use some of the algorithms provided by
1682 the C++ standard library:
1686 // print the elements of the list (requires #include <iterator>)
1687 std::copy(l.begin(), l.end(), ostream_iterator<ex>(cout, "\n"));
1689 // sum up the elements of the list (requires #include <numeric>)
1690 ex sum = std::accumulate(l.begin(), l.end(), ex(0));
1691 cout << sum << endl; // prints '2+2*x+2*y'
1695 @code{lst} is one of the few GiNaC classes that allow in-place modifications
1696 (the only other one is @code{matrix}). You can modify single elements:
1700 l[1] = 42; // l is now @{x, 42, y, x+y@}
1701 l.let_op(1) = 7; // l is now @{x, 7, y, x+y@}
1705 You can append or prepend an expression to a list with the @code{append()}
1706 and @code{prepend()} methods:
1710 l.append(4*x); // l is now @{x, 7, y, x+y, 4*x@}
1711 l.prepend(0); // l is now @{0, x, 7, y, x+y, 4*x@}
1715 You can remove the first or last element of a list with @code{remove_first()}
1716 and @code{remove_last()}:
1720 l.remove_first(); // l is now @{x, 7, y, x+y, 4*x@}
1721 l.remove_last(); // l is now @{x, 7, y, x+y@}
1725 You can remove all the elements of a list with @code{remove_all()}:
1729 l.remove_all(); // l is now empty
1733 You can bring the elements of a list into a canonical order with @code{sort()}:
1742 // l1 and l2 are now equal
1746 Finally, you can remove all but the first element of consecutive groups of
1747 elements with @code{unique()}:
1752 l3 = x, 2, 2, 2, y, x+y, y+x;
1753 l3.unique(); // l3 is now @{x, 2, y, x+y@}
1758 @node Mathematical functions, Relations, Lists, Basic Concepts
1759 @c node-name, next, previous, up
1760 @section Mathematical functions
1761 @cindex @code{function} (class)
1762 @cindex trigonometric function
1763 @cindex hyperbolic function
1765 There are quite a number of useful functions hard-wired into GiNaC. For
1766 instance, all trigonometric and hyperbolic functions are implemented
1767 (@xref{Built-in Functions}, for a complete list).
1769 These functions (better called @emph{pseudofunctions}) are all objects
1770 of class @code{function}. They accept one or more expressions as
1771 arguments and return one expression. If the arguments are not
1772 numerical, the evaluation of the function may be halted, as it does in
1773 the next example, showing how a function returns itself twice and
1774 finally an expression that may be really useful:
1776 @cindex Gamma function
1777 @cindex @code{subs()}
1780 symbol x("x"), y("y");
1782 cout << tgamma(foo) << endl;
1783 // -> tgamma(x+(1/2)*y)
1784 ex bar = foo.subs(y==1);
1785 cout << tgamma(bar) << endl;
1787 ex foobar = bar.subs(x==7);
1788 cout << tgamma(foobar) << endl;
1789 // -> (135135/128)*Pi^(1/2)
1793 Besides evaluation most of these functions allow differentiation, series
1794 expansion and so on. Read the next chapter in order to learn more about
1797 It must be noted that these pseudofunctions are created by inline
1798 functions, where the argument list is templated. This means that
1799 whenever you call @code{GiNaC::sin(1)} it is equivalent to
1800 @code{sin(ex(1))} and will therefore not result in a floating point
1801 number. Unless of course the function prototype is explicitly
1802 overridden -- which is the case for arguments of type @code{numeric}
1803 (not wrapped inside an @code{ex}). Hence, in order to obtain a floating
1804 point number of class @code{numeric} you should call
1805 @code{sin(numeric(1))}. This is almost the same as calling
1806 @code{sin(1).evalf()} except that the latter will return a numeric
1807 wrapped inside an @code{ex}.
1810 @node Relations, Matrices, Mathematical functions, Basic Concepts
1811 @c node-name, next, previous, up
1813 @cindex @code{relational} (class)
1815 Sometimes, a relation holding between two expressions must be stored
1816 somehow. The class @code{relational} is a convenient container for such
1817 purposes. A relation is by definition a container for two @code{ex} and
1818 a relation between them that signals equality, inequality and so on.
1819 They are created by simply using the C++ operators @code{==}, @code{!=},
1820 @code{<}, @code{<=}, @code{>} and @code{>=} between two expressions.
1822 @xref{Mathematical functions}, for examples where various applications
1823 of the @code{.subs()} method show how objects of class relational are
1824 used as arguments. There they provide an intuitive syntax for
1825 substitutions. They are also used as arguments to the @code{ex::series}
1826 method, where the left hand side of the relation specifies the variable
1827 to expand in and the right hand side the expansion point. They can also
1828 be used for creating systems of equations that are to be solved for
1829 unknown variables. But the most common usage of objects of this class
1830 is rather inconspicuous in statements of the form @code{if
1831 (expand(pow(a+b,2))==a*a+2*a*b+b*b) @{...@}}. Here, an implicit
1832 conversion from @code{relational} to @code{bool} takes place. Note,
1833 however, that @code{==} here does not perform any simplifications, hence
1834 @code{expand()} must be called explicitly.
1837 @node Matrices, Indexed objects, Relations, Basic Concepts
1838 @c node-name, next, previous, up
1840 @cindex @code{matrix} (class)
1842 A @dfn{matrix} is a two-dimensional array of expressions. The elements of a
1843 matrix with @math{m} rows and @math{n} columns are accessed with two
1844 @code{unsigned} indices, the first one in the range 0@dots{}@math{m-1}, the
1845 second one in the range 0@dots{}@math{n-1}.
1847 There are a couple of ways to construct matrices, with or without preset
1848 elements. The constructor
1851 matrix::matrix(unsigned r, unsigned c);
1854 creates a matrix with @samp{r} rows and @samp{c} columns with all elements
1857 The fastest way to create a matrix with preinitialized elements is to assign
1858 a list of comma-separated expressions to an empty matrix (see below for an
1859 example). But you can also specify the elements as a (flat) list with
1862 matrix::matrix(unsigned r, unsigned c, const lst & l);
1867 @cindex @code{lst_to_matrix()}
1869 ex lst_to_matrix(const lst & l);
1872 constructs a matrix from a list of lists, each list representing a matrix row.
1874 There is also a set of functions for creating some special types of
1877 @cindex @code{diag_matrix()}
1878 @cindex @code{unit_matrix()}
1879 @cindex @code{symbolic_matrix()}
1881 ex diag_matrix(const lst & l);
1882 ex unit_matrix(unsigned x);
1883 ex unit_matrix(unsigned r, unsigned c);
1884 ex symbolic_matrix(unsigned r, unsigned c, const string & base_name);
1885 ex symbolic_matrix(unsigned r, unsigned c, const string & base_name, const string & tex_base_name);
1888 @code{diag_matrix()} constructs a diagonal matrix given the list of diagonal
1889 elements. @code{unit_matrix()} creates an @samp{x} by @samp{x} (or @samp{r}
1890 by @samp{c}) unit matrix. And finally, @code{symbolic_matrix} constructs a
1891 matrix filled with newly generated symbols made of the specified base name
1892 and the position of each element in the matrix.
1894 Matrix elements can be accessed and set using the parenthesis (function call)
1898 const ex & matrix::operator()(unsigned r, unsigned c) const;
1899 ex & matrix::operator()(unsigned r, unsigned c);
1902 It is also possible to access the matrix elements in a linear fashion with
1903 the @code{op()} method. But C++-style subscripting with square brackets
1904 @samp{[]} is not available.
1906 Here are a couple of examples for constructing matrices:
1910 symbol a("a"), b("b");
1924 cout << matrix(2, 2, lst(a, 0, 0, b)) << endl;
1927 cout << lst_to_matrix(lst(lst(a, 0), lst(0, b))) << endl;
1930 cout << diag_matrix(lst(a, b)) << endl;
1933 cout << unit_matrix(3) << endl;
1934 // -> [[1,0,0],[0,1,0],[0,0,1]]
1936 cout << symbolic_matrix(2, 3, "x") << endl;
1937 // -> [[x00,x01,x02],[x10,x11,x12]]
1941 @cindex @code{transpose()}
1942 There are three ways to do arithmetic with matrices. The first (and most
1943 direct one) is to use the methods provided by the @code{matrix} class:
1946 matrix matrix::add(const matrix & other) const;
1947 matrix matrix::sub(const matrix & other) const;
1948 matrix matrix::mul(const matrix & other) const;
1949 matrix matrix::mul_scalar(const ex & other) const;
1950 matrix matrix::pow(const ex & expn) const;
1951 matrix matrix::transpose() const;
1954 All of these methods return the result as a new matrix object. Here is an
1955 example that calculates @math{A*B-2*C} for three matrices @math{A}, @math{B}
1960 matrix A(2, 2), B(2, 2), C(2, 2);
1968 matrix result = A.mul(B).sub(C.mul_scalar(2));
1969 cout << result << endl;
1970 // -> [[-13,-6],[1,2]]
1975 @cindex @code{evalm()}
1976 The second (and probably the most natural) way is to construct an expression
1977 containing matrices with the usual arithmetic operators and @code{pow()}.
1978 For efficiency reasons, expressions with sums, products and powers of
1979 matrices are not automatically evaluated in GiNaC. You have to call the
1983 ex ex::evalm() const;
1986 to obtain the result:
1993 // -> [[1,2],[3,4]]*[[-1,0],[2,1]]-2*[[8,4],[2,1]]
1994 cout << e.evalm() << endl;
1995 // -> [[-13,-6],[1,2]]
2000 The non-commutativity of the product @code{A*B} in this example is
2001 automatically recognized by GiNaC. There is no need to use a special
2002 operator here. @xref{Non-commutative objects}, for more information about
2003 dealing with non-commutative expressions.
2005 Finally, you can work with indexed matrices and call @code{simplify_indexed()}
2006 to perform the arithmetic:
2011 idx i(symbol("i"), 2), j(symbol("j"), 2), k(symbol("k"), 2);
2012 e = indexed(A, i, k) * indexed(B, k, j) - 2 * indexed(C, i, j);
2014 // -> -2*[[8,4],[2,1]].i.j+[[-1,0],[2,1]].k.j*[[1,2],[3,4]].i.k
2015 cout << e.simplify_indexed() << endl;
2016 // -> [[-13,-6],[1,2]].i.j
2020 Using indices is most useful when working with rectangular matrices and
2021 one-dimensional vectors because you don't have to worry about having to
2022 transpose matrices before multiplying them. @xref{Indexed objects}, for
2023 more information about using matrices with indices, and about indices in
2026 The @code{matrix} class provides a couple of additional methods for
2027 computing determinants, traces, characteristic polynomials and ranks:
2029 @cindex @code{determinant()}
2030 @cindex @code{trace()}
2031 @cindex @code{charpoly()}
2032 @cindex @code{rank()}
2034 ex matrix::determinant(unsigned algo=determinant_algo::automatic) const;
2035 ex matrix::trace() const;
2036 ex matrix::charpoly(const ex & lambda) const;
2037 unsigned matrix::rank() const;
2040 The @samp{algo} argument of @code{determinant()} allows to select
2041 between different algorithms for calculating the determinant. The
2042 asymptotic speed (as parametrized by the matrix size) can greatly differ
2043 between those algorithms, depending on the nature of the matrix'
2044 entries. The possible values are defined in the @file{flags.h} header
2045 file. By default, GiNaC uses a heuristic to automatically select an
2046 algorithm that is likely (but not guaranteed) to give the result most
2049 @cindex @code{inverse()} (matrix)
2050 @cindex @code{solve()}
2051 Matrices may also be inverted using the @code{ex matrix::inverse()}
2052 method and linear systems may be solved with:
2055 matrix matrix::solve(const matrix & vars, const matrix & rhs, unsigned algo=solve_algo::automatic) const;
2058 Assuming the matrix object this method is applied on is an @code{m}
2059 times @code{n} matrix, then @code{vars} must be a @code{n} times
2060 @code{p} matrix of symbolic indeterminates and @code{rhs} a @code{m}
2061 times @code{p} matrix. The returned matrix then has dimension @code{n}
2062 times @code{p} and in the case of an underdetermined system will still
2063 contain some of the indeterminates from @code{vars}. If the system is
2064 overdetermined, an exception is thrown.
2067 @node Indexed objects, Non-commutative objects, Matrices, Basic Concepts
2068 @c node-name, next, previous, up
2069 @section Indexed objects
2071 GiNaC allows you to handle expressions containing general indexed objects in
2072 arbitrary spaces. It is also able to canonicalize and simplify such
2073 expressions and perform symbolic dummy index summations. There are a number
2074 of predefined indexed objects provided, like delta and metric tensors.
2076 There are few restrictions placed on indexed objects and their indices and
2077 it is easy to construct nonsense expressions, but our intention is to
2078 provide a general framework that allows you to implement algorithms with
2079 indexed quantities, getting in the way as little as possible.
2081 @cindex @code{idx} (class)
2082 @cindex @code{indexed} (class)
2083 @subsection Indexed quantities and their indices
2085 Indexed expressions in GiNaC are constructed of two special types of objects,
2086 @dfn{index objects} and @dfn{indexed objects}.
2090 @cindex contravariant
2093 @item Index objects are of class @code{idx} or a subclass. Every index has
2094 a @dfn{value} and a @dfn{dimension} (which is the dimension of the space
2095 the index lives in) which can both be arbitrary expressions but are usually
2096 a number or a simple symbol. In addition, indices of class @code{varidx} have
2097 a @dfn{variance} (they can be co- or contravariant), and indices of class
2098 @code{spinidx} have a variance and can be @dfn{dotted} or @dfn{undotted}.
2100 @item Indexed objects are of class @code{indexed} or a subclass. They
2101 contain a @dfn{base expression} (which is the expression being indexed), and
2102 one or more indices.
2106 @strong{Note:} when printing expressions, covariant indices and indices
2107 without variance are denoted @samp{.i} while contravariant indices are
2108 denoted @samp{~i}. Dotted indices have a @samp{*} in front of the index
2109 value. In the following, we are going to use that notation in the text so
2110 instead of @math{A^i_jk} we will write @samp{A~i.j.k}. Index dimensions are
2111 not visible in the output.
2113 A simple example shall illustrate the concepts:
2117 #include <ginac/ginac.h>
2118 using namespace std;
2119 using namespace GiNaC;
2123 symbol i_sym("i"), j_sym("j");
2124 idx i(i_sym, 3), j(j_sym, 3);
2127 cout << indexed(A, i, j) << endl;
2129 cout << index_dimensions << indexed(A, i, j) << endl;
2131 cout << dflt; // reset cout to default output format (dimensions hidden)
2135 The @code{idx} constructor takes two arguments, the index value and the
2136 index dimension. First we define two index objects, @code{i} and @code{j},
2137 both with the numeric dimension 3. The value of the index @code{i} is the
2138 symbol @code{i_sym} (which prints as @samp{i}) and the value of the index
2139 @code{j} is the symbol @code{j_sym} (which prints as @samp{j}). Next we
2140 construct an expression containing one indexed object, @samp{A.i.j}. It has
2141 the symbol @code{A} as its base expression and the two indices @code{i} and
2144 The dimensions of indices are normally not visible in the output, but one
2145 can request them to be printed with the @code{index_dimensions} manipulator,
2148 Note the difference between the indices @code{i} and @code{j} which are of
2149 class @code{idx}, and the index values which are the symbols @code{i_sym}
2150 and @code{j_sym}. The indices of indexed objects cannot directly be symbols
2151 or numbers but must be index objects. For example, the following is not
2152 correct and will raise an exception:
2155 symbol i("i"), j("j");
2156 e = indexed(A, i, j); // ERROR: indices must be of type idx
2159 You can have multiple indexed objects in an expression, index values can
2160 be numeric, and index dimensions symbolic:
2164 symbol B("B"), dim("dim");
2165 cout << 4 * indexed(A, i)
2166 + indexed(B, idx(j_sym, 4), idx(2, 3), idx(i_sym, dim)) << endl;
2171 @code{B} has a 4-dimensional symbolic index @samp{k}, a 3-dimensional numeric
2172 index of value 2, and a symbolic index @samp{i} with the symbolic dimension
2173 @samp{dim}. Note that GiNaC doesn't automatically notify you that the free
2174 indices of @samp{A} and @samp{B} in the sum don't match (you have to call
2175 @code{simplify_indexed()} for that, see below).
2177 In fact, base expressions, index values and index dimensions can be
2178 arbitrary expressions:
2182 cout << indexed(A+B, idx(2*i_sym+1, dim/2)) << endl;
2187 It's also possible to construct nonsense like @samp{Pi.sin(x)}. You will not
2188 get an error message from this but you will probably not be able to do
2189 anything useful with it.
2191 @cindex @code{get_value()}
2192 @cindex @code{get_dimension()}
2196 ex idx::get_value();
2197 ex idx::get_dimension();
2200 return the value and dimension of an @code{idx} object. If you have an index
2201 in an expression, such as returned by calling @code{.op()} on an indexed
2202 object, you can get a reference to the @code{idx} object with the function
2203 @code{ex_to<idx>()} on the expression.
2205 There are also the methods
2208 bool idx::is_numeric();
2209 bool idx::is_symbolic();
2210 bool idx::is_dim_numeric();
2211 bool idx::is_dim_symbolic();
2214 for checking whether the value and dimension are numeric or symbolic
2215 (non-numeric). Using the @code{info()} method of an index (see @ref{Information
2216 About Expressions}) returns information about the index value.
2218 @cindex @code{varidx} (class)
2219 If you need co- and contravariant indices, use the @code{varidx} class:
2223 symbol mu_sym("mu"), nu_sym("nu");
2224 varidx mu(mu_sym, 4), nu(nu_sym, 4); // default is contravariant ~mu, ~nu
2225 varidx mu_co(mu_sym, 4, true); // covariant index .mu
2227 cout << indexed(A, mu, nu) << endl;
2229 cout << indexed(A, mu_co, nu) << endl;
2231 cout << indexed(A, mu.toggle_variance(), nu) << endl;
2236 A @code{varidx} is an @code{idx} with an additional flag that marks it as
2237 co- or contravariant. The default is a contravariant (upper) index, but
2238 this can be overridden by supplying a third argument to the @code{varidx}
2239 constructor. The two methods
2242 bool varidx::is_covariant();
2243 bool varidx::is_contravariant();
2246 allow you to check the variance of a @code{varidx} object (use @code{ex_to<varidx>()}
2247 to get the object reference from an expression). There's also the very useful
2251 ex varidx::toggle_variance();
2254 which makes a new index with the same value and dimension but the opposite
2255 variance. By using it you only have to define the index once.
2257 @cindex @code{spinidx} (class)
2258 The @code{spinidx} class provides dotted and undotted variant indices, as
2259 used in the Weyl-van-der-Waerden spinor formalism:
2263 symbol K("K"), C_sym("C"), D_sym("D");
2264 spinidx C(C_sym, 2), D(D_sym); // default is 2-dimensional,
2265 // contravariant, undotted
2266 spinidx C_co(C_sym, 2, true); // covariant index
2267 spinidx D_dot(D_sym, 2, false, true); // contravariant, dotted
2268 spinidx D_co_dot(D_sym, 2, true, true); // covariant, dotted
2270 cout << indexed(K, C, D) << endl;
2272 cout << indexed(K, C_co, D_dot) << endl;
2274 cout << indexed(K, D_co_dot, D) << endl;
2279 A @code{spinidx} is a @code{varidx} with an additional flag that marks it as
2280 dotted or undotted. The default is undotted but this can be overridden by
2281 supplying a fourth argument to the @code{spinidx} constructor. The two
2285 bool spinidx::is_dotted();
2286 bool spinidx::is_undotted();
2289 allow you to check whether or not a @code{spinidx} object is dotted (use
2290 @code{ex_to<spinidx>()} to get the object reference from an expression).
2291 Finally, the two methods
2294 ex spinidx::toggle_dot();
2295 ex spinidx::toggle_variance_dot();
2298 create a new index with the same value and dimension but opposite dottedness
2299 and the same or opposite variance.
2301 @subsection Substituting indices
2303 @cindex @code{subs()}
2304 Sometimes you will want to substitute one symbolic index with another
2305 symbolic or numeric index, for example when calculating one specific element
2306 of a tensor expression. This is done with the @code{.subs()} method, as it
2307 is done for symbols (see @ref{Substituting Expressions}).
2309 You have two possibilities here. You can either substitute the whole index
2310 by another index or expression:
2314 ex e = indexed(A, mu_co);
2315 cout << e << " becomes " << e.subs(mu_co == nu) << endl;
2316 // -> A.mu becomes A~nu
2317 cout << e << " becomes " << e.subs(mu_co == varidx(0, 4)) << endl;
2318 // -> A.mu becomes A~0
2319 cout << e << " becomes " << e.subs(mu_co == 0) << endl;
2320 // -> A.mu becomes A.0
2324 The third example shows that trying to replace an index with something that
2325 is not an index will substitute the index value instead.
2327 Alternatively, you can substitute the @emph{symbol} of a symbolic index by
2332 ex e = indexed(A, mu_co);
2333 cout << e << " becomes " << e.subs(mu_sym == nu_sym) << endl;
2334 // -> A.mu becomes A.nu
2335 cout << e << " becomes " << e.subs(mu_sym == 0) << endl;
2336 // -> A.mu becomes A.0
2340 As you see, with the second method only the value of the index will get
2341 substituted. Its other properties, including its dimension, remain unchanged.
2342 If you want to change the dimension of an index you have to substitute the
2343 whole index by another one with the new dimension.
2345 Finally, substituting the base expression of an indexed object works as
2350 ex e = indexed(A, mu_co);
2351 cout << e << " becomes " << e.subs(A == A+B) << endl;
2352 // -> A.mu becomes (B+A).mu
2356 @subsection Symmetries
2357 @cindex @code{symmetry} (class)
2358 @cindex @code{sy_none()}
2359 @cindex @code{sy_symm()}
2360 @cindex @code{sy_anti()}
2361 @cindex @code{sy_cycl()}
2363 Indexed objects can have certain symmetry properties with respect to their
2364 indices. Symmetries are specified as a tree of objects of class @code{symmetry}
2365 that is constructed with the helper functions
2368 symmetry sy_none(...);
2369 symmetry sy_symm(...);
2370 symmetry sy_anti(...);
2371 symmetry sy_cycl(...);
2374 @code{sy_none()} stands for no symmetry, @code{sy_symm()} and @code{sy_anti()}
2375 specify fully symmetric or antisymmetric, respectively, and @code{sy_cycl()}
2376 represents a cyclic symmetry. Each of these functions accepts up to four
2377 arguments which can be either symmetry objects themselves or unsigned integer
2378 numbers that represent an index position (counting from 0). A symmetry
2379 specification that consists of only a single @code{sy_symm()}, @code{sy_anti()}
2380 or @code{sy_cycl()} with no arguments specifies the respective symmetry for
2383 Here are some examples of symmetry definitions:
2388 e = indexed(A, i, j);
2389 e = indexed(A, sy_none(), i, j); // equivalent
2390 e = indexed(A, sy_none(0, 1), i, j); // equivalent
2392 // Symmetric in all three indices:
2393 e = indexed(A, sy_symm(), i, j, k);
2394 e = indexed(A, sy_symm(0, 1, 2), i, j, k); // equivalent
2395 e = indexed(A, sy_symm(2, 0, 1), i, j, k); // same symmetry, but yields a
2396 // different canonical order
2398 // Symmetric in the first two indices only:
2399 e = indexed(A, sy_symm(0, 1), i, j, k);
2400 e = indexed(A, sy_none(sy_symm(0, 1), 2), i, j, k); // equivalent
2402 // Antisymmetric in the first and last index only (index ranges need not
2404 e = indexed(A, sy_anti(0, 2), i, j, k);
2405 e = indexed(A, sy_none(sy_anti(0, 2), 1), i, j, k); // equivalent
2407 // An example of a mixed symmetry: antisymmetric in the first two and
2408 // last two indices, symmetric when swapping the first and last index
2409 // pairs (like the Riemann curvature tensor):
2410 e = indexed(A, sy_symm(sy_anti(0, 1), sy_anti(2, 3)), i, j, k, l);
2412 // Cyclic symmetry in all three indices:
2413 e = indexed(A, sy_cycl(), i, j, k);
2414 e = indexed(A, sy_cycl(0, 1, 2), i, j, k); // equivalent
2416 // The following examples are invalid constructions that will throw
2417 // an exception at run time.
2419 // An index may not appear multiple times:
2420 e = indexed(A, sy_symm(0, 0, 1), i, j, k); // ERROR
2421 e = indexed(A, sy_none(sy_symm(0, 1), sy_anti(0, 2)), i, j, k); // ERROR
2423 // Every child of sy_symm(), sy_anti() and sy_cycl() must refer to the
2424 // same number of indices:
2425 e = indexed(A, sy_symm(sy_anti(0, 1), 2), i, j, k); // ERROR
2427 // And of course, you cannot specify indices which are not there:
2428 e = indexed(A, sy_symm(0, 1, 2, 3), i, j, k); // ERROR
2432 If you need to specify more than four indices, you have to use the
2433 @code{.add()} method of the @code{symmetry} class. For example, to specify
2434 full symmetry in the first six indices you would write
2435 @code{sy_symm(0, 1, 2, 3).add(4).add(5)}.
2437 If an indexed object has a symmetry, GiNaC will automatically bring the
2438 indices into a canonical order which allows for some immediate simplifications:
2442 cout << indexed(A, sy_symm(), i, j)
2443 + indexed(A, sy_symm(), j, i) << endl;
2445 cout << indexed(B, sy_anti(), i, j)
2446 + indexed(B, sy_anti(), j, i) << endl;
2448 cout << indexed(B, sy_anti(), i, j, k)
2449 - indexed(B, sy_anti(), j, k, i) << endl;
2454 @cindex @code{get_free_indices()}
2456 @subsection Dummy indices
2458 GiNaC treats certain symbolic index pairs as @dfn{dummy indices} meaning
2459 that a summation over the index range is implied. Symbolic indices which are
2460 not dummy indices are called @dfn{free indices}. Numeric indices are neither
2461 dummy nor free indices.
2463 To be recognized as a dummy index pair, the two indices must be of the same
2464 class and their value must be the same single symbol (an index like
2465 @samp{2*n+1} is never a dummy index). If the indices are of class
2466 @code{varidx} they must also be of opposite variance; if they are of class
2467 @code{spinidx} they must be both dotted or both undotted.
2469 The method @code{.get_free_indices()} returns a vector containing the free
2470 indices of an expression. It also checks that the free indices of the terms
2471 of a sum are consistent:
2475 symbol A("A"), B("B"), C("C");
2477 symbol i_sym("i"), j_sym("j"), k_sym("k"), l_sym("l");
2478 idx i(i_sym, 3), j(j_sym, 3), k(k_sym, 3), l(l_sym, 3);
2480 ex e = indexed(A, i, j) * indexed(B, j, k) + indexed(C, k, l, i, l);
2481 cout << exprseq(e.get_free_indices()) << endl;
2483 // 'j' and 'l' are dummy indices
2485 symbol mu_sym("mu"), nu_sym("nu"), rho_sym("rho"), sigma_sym("sigma");
2486 varidx mu(mu_sym, 4), nu(nu_sym, 4), rho(rho_sym, 4), sigma(sigma_sym, 4);
2488 e = indexed(A, mu, nu) * indexed(B, nu.toggle_variance(), rho)
2489 + indexed(C, mu, sigma, rho, sigma.toggle_variance());
2490 cout << exprseq(e.get_free_indices()) << endl;
2492 // 'nu' is a dummy index, but 'sigma' is not
2494 e = indexed(A, mu, mu);
2495 cout << exprseq(e.get_free_indices()) << endl;
2497 // 'mu' is not a dummy index because it appears twice with the same
2500 e = indexed(A, mu, nu) + 42;
2501 cout << exprseq(e.get_free_indices()) << endl; // ERROR
2502 // this will throw an exception:
2503 // "add::get_free_indices: inconsistent indices in sum"
2507 @cindex @code{simplify_indexed()}
2508 @subsection Simplifying indexed expressions
2510 In addition to the few automatic simplifications that GiNaC performs on
2511 indexed expressions (such as re-ordering the indices of symmetric tensors
2512 and calculating traces and convolutions of matrices and predefined tensors)
2516 ex ex::simplify_indexed();
2517 ex ex::simplify_indexed(const scalar_products & sp);
2520 that performs some more expensive operations:
2523 @item it checks the consistency of free indices in sums in the same way
2524 @code{get_free_indices()} does
2525 @item it tries to give dummy indices that appear in different terms of a sum
2526 the same name to allow simplifications like @math{a_i*b_i-a_j*b_j=0}
2527 @item it (symbolically) calculates all possible dummy index summations/contractions
2528 with the predefined tensors (this will be explained in more detail in the
2530 @item it detects contractions that vanish for symmetry reasons, for example
2531 the contraction of a symmetric and a totally antisymmetric tensor
2532 @item as a special case of dummy index summation, it can replace scalar products
2533 of two tensors with a user-defined value
2536 The last point is done with the help of the @code{scalar_products} class
2537 which is used to store scalar products with known values (this is not an
2538 arithmetic class, you just pass it to @code{simplify_indexed()}):
2542 symbol A("A"), B("B"), C("C"), i_sym("i");
2546 sp.add(A, B, 0); // A and B are orthogonal
2547 sp.add(A, C, 0); // A and C are orthogonal
2548 sp.add(A, A, 4); // A^2 = 4 (A has length 2)
2550 e = indexed(A + B, i) * indexed(A + C, i);
2552 // -> (B+A).i*(A+C).i
2554 cout << e.expand(expand_options::expand_indexed).simplify_indexed(sp)
2560 The @code{scalar_products} object @code{sp} acts as a storage for the
2561 scalar products added to it with the @code{.add()} method. This method
2562 takes three arguments: the two expressions of which the scalar product is
2563 taken, and the expression to replace it with. After @code{sp.add(A, B, 0)},
2564 @code{simplify_indexed()} will replace all scalar products of indexed
2565 objects that have the symbols @code{A} and @code{B} as base expressions
2566 with the single value 0. The number, type and dimension of the indices
2567 don't matter; @samp{A~mu~nu*B.mu.nu} would also be replaced by 0.
2569 @cindex @code{expand()}
2570 The example above also illustrates a feature of the @code{expand()} method:
2571 if passed the @code{expand_indexed} option it will distribute indices
2572 over sums, so @samp{(A+B).i} becomes @samp{A.i+B.i}.
2574 @cindex @code{tensor} (class)
2575 @subsection Predefined tensors
2577 Some frequently used special tensors such as the delta, epsilon and metric
2578 tensors are predefined in GiNaC. They have special properties when
2579 contracted with other tensor expressions and some of them have constant
2580 matrix representations (they will evaluate to a number when numeric
2581 indices are specified).
2583 @cindex @code{delta_tensor()}
2584 @subsubsection Delta tensor
2586 The delta tensor takes two indices, is symmetric and has the matrix
2587 representation @code{diag(1, 1, 1, ...)}. It is constructed by the function
2588 @code{delta_tensor()}:
2592 symbol A("A"), B("B");
2594 idx i(symbol("i"), 3), j(symbol("j"), 3),
2595 k(symbol("k"), 3), l(symbol("l"), 3);
2597 ex e = indexed(A, i, j) * indexed(B, k, l)
2598 * delta_tensor(i, k) * delta_tensor(j, l) << endl;
2599 cout << e.simplify_indexed() << endl;
2602 cout << delta_tensor(i, i) << endl;
2607 @cindex @code{metric_tensor()}
2608 @subsubsection General metric tensor
2610 The function @code{metric_tensor()} creates a general symmetric metric
2611 tensor with two indices that can be used to raise/lower tensor indices. The
2612 metric tensor is denoted as @samp{g} in the output and if its indices are of
2613 mixed variance it is automatically replaced by a delta tensor:
2619 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4);
2621 ex e = metric_tensor(mu, nu) * indexed(A, nu.toggle_variance(), rho);
2622 cout << e.simplify_indexed() << endl;
2625 e = delta_tensor(mu, nu.toggle_variance()) * metric_tensor(nu, rho);
2626 cout << e.simplify_indexed() << endl;
2629 e = metric_tensor(mu.toggle_variance(), nu.toggle_variance())
2630 * metric_tensor(nu, rho);
2631 cout << e.simplify_indexed() << endl;
2634 e = metric_tensor(nu.toggle_variance(), rho.toggle_variance())
2635 * metric_tensor(mu, nu) * (delta_tensor(mu.toggle_variance(), rho)
2636 + indexed(A, mu.toggle_variance(), rho));
2637 cout << e.simplify_indexed() << endl;
2642 @cindex @code{lorentz_g()}
2643 @subsubsection Minkowski metric tensor
2645 The Minkowski metric tensor is a special metric tensor with a constant
2646 matrix representation which is either @code{diag(1, -1, -1, ...)} (negative
2647 signature, the default) or @code{diag(-1, 1, 1, ...)} (positive signature).
2648 It is created with the function @code{lorentz_g()} (although it is output as
2653 varidx mu(symbol("mu"), 4);
2655 e = delta_tensor(varidx(0, 4), mu.toggle_variance())
2656 * lorentz_g(mu, varidx(0, 4)); // negative signature
2657 cout << e.simplify_indexed() << endl;
2660 e = delta_tensor(varidx(0, 4), mu.toggle_variance())
2661 * lorentz_g(mu, varidx(0, 4), true); // positive signature
2662 cout << e.simplify_indexed() << endl;
2667 @cindex @code{spinor_metric()}
2668 @subsubsection Spinor metric tensor
2670 The function @code{spinor_metric()} creates an antisymmetric tensor with
2671 two indices that is used to raise/lower indices of 2-component spinors.
2672 It is output as @samp{eps}:
2678 spinidx A(symbol("A")), B(symbol("B")), C(symbol("C"));
2679 ex A_co = A.toggle_variance(), B_co = B.toggle_variance();
2681 e = spinor_metric(A, B) * indexed(psi, B_co);
2682 cout << e.simplify_indexed() << endl;
2685 e = spinor_metric(A, B) * indexed(psi, A_co);
2686 cout << e.simplify_indexed() << endl;
2689 e = spinor_metric(A_co, B_co) * indexed(psi, B);
2690 cout << e.simplify_indexed() << endl;
2693 e = spinor_metric(A_co, B_co) * indexed(psi, A);
2694 cout << e.simplify_indexed() << endl;
2697 e = spinor_metric(A_co, B_co) * spinor_metric(A, B);
2698 cout << e.simplify_indexed() << endl;
2701 e = spinor_metric(A_co, B_co) * spinor_metric(B, C);
2702 cout << e.simplify_indexed() << endl;
2707 The matrix representation of the spinor metric is @code{[[0, 1], [-1, 0]]}.
2709 @cindex @code{epsilon_tensor()}
2710 @cindex @code{lorentz_eps()}
2711 @subsubsection Epsilon tensor
2713 The epsilon tensor is totally antisymmetric, its number of indices is equal
2714 to the dimension of the index space (the indices must all be of the same
2715 numeric dimension), and @samp{eps.1.2.3...} (resp. @samp{eps~0~1~2...}) is
2716 defined to be 1. Its behavior with indices that have a variance also
2717 depends on the signature of the metric. Epsilon tensors are output as
2720 There are three functions defined to create epsilon tensors in 2, 3 and 4
2724 ex epsilon_tensor(const ex & i1, const ex & i2);
2725 ex epsilon_tensor(const ex & i1, const ex & i2, const ex & i3);
2726 ex lorentz_eps(const ex & i1, const ex & i2, const ex & i3, const ex & i4, bool pos_sig = false);
2729 The first two functions create an epsilon tensor in 2 or 3 Euclidean
2730 dimensions, the last function creates an epsilon tensor in a 4-dimensional
2731 Minkowski space (the last @code{bool} argument specifies whether the metric
2732 has negative or positive signature, as in the case of the Minkowski metric
2737 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4),
2738 sig(symbol("sig"), 4), lam(symbol("lam"), 4), bet(symbol("bet"), 4);
2739 e = lorentz_eps(mu, nu, rho, sig) *
2740 lorentz_eps(mu.toggle_variance(), nu.toggle_variance(), lam, bet);
2741 cout << simplify_indexed(e) << endl;
2742 // -> 2*eta~bet~rho*eta~sig~lam-2*eta~sig~bet*eta~rho~lam
2744 idx i(symbol("i"), 3), j(symbol("j"), 3), k(symbol("k"), 3);
2745 symbol A("A"), B("B");
2746 e = epsilon_tensor(i, j, k) * indexed(A, j) * indexed(B, k);
2747 cout << simplify_indexed(e) << endl;
2748 // -> -B.k*A.j*eps.i.k.j
2749 e = epsilon_tensor(i, j, k) * indexed(A, j) * indexed(A, k);
2750 cout << simplify_indexed(e) << endl;
2755 @subsection Linear algebra
2757 The @code{matrix} class can be used with indices to do some simple linear
2758 algebra (linear combinations and products of vectors and matrices, traces
2759 and scalar products):
2763 idx i(symbol("i"), 2), j(symbol("j"), 2);
2764 symbol x("x"), y("y");
2766 // A is a 2x2 matrix, X is a 2x1 vector
2767 matrix A(2, 2), X(2, 1);
2772 cout << indexed(A, i, i) << endl;
2775 ex e = indexed(A, i, j) * indexed(X, j);
2776 cout << e.simplify_indexed() << endl;
2777 // -> [[2*y+x],[4*y+3*x]].i
2779 e = indexed(A, i, j) * indexed(X, i) + indexed(X, j) * 2;
2780 cout << e.simplify_indexed() << endl;
2781 // -> [[3*y+3*x,6*y+2*x]].j
2785 You can of course obtain the same results with the @code{matrix::add()},
2786 @code{matrix::mul()} and @code{matrix::trace()} methods (@pxref{Matrices})
2787 but with indices you don't have to worry about transposing matrices.
2789 Matrix indices always start at 0 and their dimension must match the number
2790 of rows/columns of the matrix. Matrices with one row or one column are
2791 vectors and can have one or two indices (it doesn't matter whether it's a
2792 row or a column vector). Other matrices must have two indices.
2794 You should be careful when using indices with variance on matrices. GiNaC
2795 doesn't look at the variance and doesn't know that @samp{F~mu~nu} and
2796 @samp{F.mu.nu} are different matrices. In this case you should use only
2797 one form for @samp{F} and explicitly multiply it with a matrix representation
2798 of the metric tensor.
2801 @node Non-commutative objects, Hash Maps, Indexed objects, Basic Concepts
2802 @c node-name, next, previous, up
2803 @section Non-commutative objects
2805 GiNaC is equipped to handle certain non-commutative algebras. Three classes of
2806 non-commutative objects are built-in which are mostly of use in high energy
2810 @item Clifford (Dirac) algebra (class @code{clifford})
2811 @item su(3) Lie algebra (class @code{color})
2812 @item Matrices (unindexed) (class @code{matrix})
2815 The @code{clifford} and @code{color} classes are subclasses of
2816 @code{indexed} because the elements of these algebras usually carry
2817 indices. The @code{matrix} class is described in more detail in
2820 Unlike most computer algebra systems, GiNaC does not primarily provide an
2821 operator (often denoted @samp{&*}) for representing inert products of
2822 arbitrary objects. Rather, non-commutativity in GiNaC is a property of the
2823 classes of objects involved, and non-commutative products are formed with
2824 the usual @samp{*} operator, as are ordinary products. GiNaC is capable of
2825 figuring out by itself which objects commutate and will group the factors
2826 by their class. Consider this example:
2830 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4);
2831 idx a(symbol("a"), 8), b(symbol("b"), 8);
2832 ex e = -dirac_gamma(mu) * (2*color_T(a)) * 8 * color_T(b) * dirac_gamma(nu);
2834 // -> -16*(gamma~mu*gamma~nu)*(T.a*T.b)
2838 As can be seen, GiNaC pulls out the overall commutative factor @samp{-16} and
2839 groups the non-commutative factors (the gammas and the su(3) generators)
2840 together while preserving the order of factors within each class (because
2841 Clifford objects commutate with color objects). The resulting expression is a
2842 @emph{commutative} product with two factors that are themselves non-commutative
2843 products (@samp{gamma~mu*gamma~nu} and @samp{T.a*T.b}). For clarification,
2844 parentheses are placed around the non-commutative products in the output.
2846 @cindex @code{ncmul} (class)
2847 Non-commutative products are internally represented by objects of the class
2848 @code{ncmul}, as opposed to commutative products which are handled by the
2849 @code{mul} class. You will normally not have to worry about this distinction,
2852 The advantage of this approach is that you never have to worry about using
2853 (or forgetting to use) a special operator when constructing non-commutative
2854 expressions. Also, non-commutative products in GiNaC are more intelligent
2855 than in other computer algebra systems; they can, for example, automatically
2856 canonicalize themselves according to rules specified in the implementation
2857 of the non-commutative classes. The drawback is that to work with other than
2858 the built-in algebras you have to implement new classes yourself. Symbols
2859 always commutate and it's not possible to construct non-commutative products
2860 using symbols to represent the algebra elements or generators. User-defined
2861 functions can, however, be specified as being non-commutative.
2863 @cindex @code{return_type()}
2864 @cindex @code{return_type_tinfo()}
2865 Information about the commutativity of an object or expression can be
2866 obtained with the two member functions
2869 unsigned ex::return_type() const;
2870 unsigned ex::return_type_tinfo() const;
2873 The @code{return_type()} function returns one of three values (defined in
2874 the header file @file{flags.h}), corresponding to three categories of
2875 expressions in GiNaC:
2878 @item @code{return_types::commutative}: Commutates with everything. Most GiNaC
2879 classes are of this kind.
2880 @item @code{return_types::noncommutative}: Non-commutative, belonging to a
2881 certain class of non-commutative objects which can be determined with the
2882 @code{return_type_tinfo()} method. Expressions of this category commutate
2883 with everything except @code{noncommutative} expressions of the same
2885 @item @code{return_types::noncommutative_composite}: Non-commutative, composed
2886 of non-commutative objects of different classes. Expressions of this
2887 category don't commutate with any other @code{noncommutative} or
2888 @code{noncommutative_composite} expressions.
2891 The value returned by the @code{return_type_tinfo()} method is valid only
2892 when the return type of the expression is @code{noncommutative}. It is a
2893 value that is unique to the class of the object and usually one of the
2894 constants in @file{tinfos.h}, or derived therefrom.
2896 Here are a couple of examples:
2899 @multitable @columnfractions 0.33 0.33 0.34
2900 @item @strong{Expression} @tab @strong{@code{return_type()}} @tab @strong{@code{return_type_tinfo()}}
2901 @item @code{42} @tab @code{commutative} @tab -
2902 @item @code{2*x-y} @tab @code{commutative} @tab -
2903 @item @code{dirac_ONE()} @tab @code{noncommutative} @tab @code{TINFO_clifford}
2904 @item @code{dirac_gamma(mu)*dirac_gamma(nu)} @tab @code{noncommutative} @tab @code{TINFO_clifford}
2905 @item @code{2*color_T(a)} @tab @code{noncommutative} @tab @code{TINFO_color}
2906 @item @code{dirac_ONE()*color_T(a)} @tab @code{noncommutative_composite} @tab -
2910 Note: the @code{return_type_tinfo()} of Clifford objects is only equal to
2911 @code{TINFO_clifford} for objects with a representation label of zero.
2912 Other representation labels yield a different @code{return_type_tinfo()},
2913 but it's the same for any two objects with the same label. This is also true
2916 A last note: With the exception of matrices, positive integer powers of
2917 non-commutative objects are automatically expanded in GiNaC. For example,
2918 @code{pow(a*b, 2)} becomes @samp{a*b*a*b} if @samp{a} and @samp{b} are
2919 non-commutative expressions).
2922 @cindex @code{clifford} (class)
2923 @subsection Clifford algebra
2925 @cindex @code{dirac_gamma()}
2926 Clifford algebra elements (also called Dirac gamma matrices, although GiNaC
2927 doesn't treat them as matrices) are designated as @samp{gamma~mu} and satisfy
2928 @samp{gamma~mu*gamma~nu + gamma~nu*gamma~mu = 2*eta~mu~nu} where @samp{eta~mu~nu}
2929 is the Minkowski metric tensor. Dirac gammas are constructed by the function
2932 ex dirac_gamma(const ex & mu, unsigned char rl = 0);
2935 which takes two arguments: the index and a @dfn{representation label} in the
2936 range 0 to 255 which is used to distinguish elements of different Clifford
2937 algebras (this is also called a @dfn{spin line index}). Gammas with different
2938 labels commutate with each other. The dimension of the index can be 4 or (in
2939 the framework of dimensional regularization) any symbolic value. Spinor
2940 indices on Dirac gammas are not supported in GiNaC.
2942 @cindex @code{dirac_ONE()}
2943 The unity element of a Clifford algebra is constructed by
2946 ex dirac_ONE(unsigned char rl = 0);
2949 @strong{Note:} You must always use @code{dirac_ONE()} when referring to
2950 multiples of the unity element, even though it's customary to omit it.
2951 E.g. instead of @code{dirac_gamma(mu)*(dirac_slash(q,4)+m)} you have to
2952 write @code{dirac_gamma(mu)*(dirac_slash(q,4)+m*dirac_ONE())}. Otherwise,
2953 GiNaC will complain and/or produce incorrect results.
2955 @cindex @code{dirac_gamma5()}
2956 There is a special element @samp{gamma5} that commutates with all other
2957 gammas, has a unit square, and in 4 dimensions equals
2958 @samp{gamma~0 gamma~1 gamma~2 gamma~3}, provided by
2961 ex dirac_gamma5(unsigned char rl = 0);
2964 @cindex @code{dirac_gammaL()}
2965 @cindex @code{dirac_gammaR()}
2966 The chiral projectors @samp{(1+/-gamma5)/2} are also available as proper
2967 objects, constructed by
2970 ex dirac_gammaL(unsigned char rl = 0);
2971 ex dirac_gammaR(unsigned char rl = 0);
2974 They observe the relations @samp{gammaL^2 = gammaL}, @samp{gammaR^2 = gammaR},
2975 and @samp{gammaL gammaR = gammaR gammaL = 0}.
2977 @cindex @code{dirac_slash()}
2978 Finally, the function
2981 ex dirac_slash(const ex & e, const ex & dim, unsigned char rl = 0);
2984 creates a term that represents a contraction of @samp{e} with the Dirac
2985 Lorentz vector (it behaves like a term of the form @samp{e.mu gamma~mu}
2986 with a unique index whose dimension is given by the @code{dim} argument).
2987 Such slashed expressions are printed with a trailing backslash, e.g. @samp{e\}.
2989 In products of dirac gammas, superfluous unity elements are automatically
2990 removed, squares are replaced by their values, and @samp{gamma5}, @samp{gammaL}
2991 and @samp{gammaR} are moved to the front.
2993 The @code{simplify_indexed()} function performs contractions in gamma strings,
2999 symbol a("a"), b("b"), D("D");
3000 varidx mu(symbol("mu"), D);
3001 ex e = dirac_gamma(mu) * dirac_slash(a, D)
3002 * dirac_gamma(mu.toggle_variance());
3004 // -> gamma~mu*a\*gamma.mu
3005 e = e.simplify_indexed();
3008 cout << e.subs(D == 4) << endl;
3014 @cindex @code{dirac_trace()}
3015 To calculate the trace of an expression containing strings of Dirac gammas
3016 you use one of the functions
3019 ex dirac_trace(const ex & e, const std::set<unsigned char> & rls, const ex & trONE = 4);
3020 ex dirac_trace(const ex & e, const lst & rll, const ex & trONE = 4);
3021 ex dirac_trace(const ex & e, unsigned char rl = 0, const ex & trONE = 4);
3024 These functions take the trace over all gammas in the specified set @code{rls}
3025 or list @code{rll} of representation labels, or the single label @code{rl};
3026 gammas with other labels are left standing. The last argument to
3027 @code{dirac_trace()} is the value to be returned for the trace of the unity
3028 element, which defaults to 4.
3030 The @code{dirac_trace()} function is a linear functional that is equal to the
3031 ordinary matrix trace only in @math{D = 4} dimensions. In particular, the
3032 functional is not cyclic in @math{D != 4} dimensions when acting on
3033 expressions containing @samp{gamma5}, so it's not a proper trace. This
3034 @samp{gamma5} scheme is described in greater detail in
3035 @cite{The Role of gamma5 in Dimensional Regularization}.
3037 The value of the trace itself is also usually different in 4 and in
3038 @math{D != 4} dimensions:
3043 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4);
3044 ex e = dirac_gamma(mu) * dirac_gamma(nu) *
3045 dirac_gamma(mu.toggle_variance()) * dirac_gamma(rho);
3046 cout << dirac_trace(e).simplify_indexed() << endl;
3053 varidx mu(symbol("mu"), D), nu(symbol("nu"), D), rho(symbol("rho"), D);
3054 ex e = dirac_gamma(mu) * dirac_gamma(nu) *
3055 dirac_gamma(mu.toggle_variance()) * dirac_gamma(rho);
3056 cout << dirac_trace(e).simplify_indexed() << endl;
3057 // -> 8*eta~rho~nu-4*eta~rho~nu*D
3061 Here is an example for using @code{dirac_trace()} to compute a value that
3062 appears in the calculation of the one-loop vacuum polarization amplitude in
3067 symbol q("q"), l("l"), m("m"), ldotq("ldotq"), D("D");
3068 varidx mu(symbol("mu"), D), nu(symbol("nu"), D);
3071 sp.add(l, l, pow(l, 2));
3072 sp.add(l, q, ldotq);
3074 ex e = dirac_gamma(mu) *
3075 (dirac_slash(l, D) + dirac_slash(q, D) + m * dirac_ONE()) *
3076 dirac_gamma(mu.toggle_variance()) *
3077 (dirac_slash(l, D) + m * dirac_ONE());
3078 e = dirac_trace(e).simplify_indexed(sp);
3079 e = e.collect(lst(l, ldotq, m));
3081 // -> (8-4*D)*l^2+(8-4*D)*ldotq+4*D*m^2
3085 The @code{canonicalize_clifford()} function reorders all gamma products that
3086 appear in an expression to a canonical (but not necessarily simple) form.
3087 You can use this to compare two expressions or for further simplifications:
3091 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4);
3092 ex e = dirac_gamma(mu) * dirac_gamma(nu) + dirac_gamma(nu) * dirac_gamma(mu);
3094 // -> gamma~mu*gamma~nu+gamma~nu*gamma~mu
3096 e = canonicalize_clifford(e);
3098 // -> 2*ONE*eta~mu~nu
3103 @cindex @code{color} (class)
3104 @subsection Color algebra
3106 @cindex @code{color_T()}
3107 For computations in quantum chromodynamics, GiNaC implements the base elements
3108 and structure constants of the su(3) Lie algebra (color algebra). The base
3109 elements @math{T_a} are constructed by the function
3112 ex color_T(const ex & a, unsigned char rl = 0);
3115 which takes two arguments: the index and a @dfn{representation label} in the
3116 range 0 to 255 which is used to distinguish elements of different color
3117 algebras. Objects with different labels commutate with each other. The
3118 dimension of the index must be exactly 8 and it should be of class @code{idx},
3121 @cindex @code{color_ONE()}
3122 The unity element of a color algebra is constructed by
3125 ex color_ONE(unsigned char rl = 0);
3128 @strong{Note:} You must always use @code{color_ONE()} when referring to
3129 multiples of the unity element, even though it's customary to omit it.
3130 E.g. instead of @code{color_T(a)*(color_T(b)*indexed(X,b)+1)} you have to
3131 write @code{color_T(a)*(color_T(b)*indexed(X,b)+color_ONE())}. Otherwise,
3132 GiNaC may produce incorrect results.
3134 @cindex @code{color_d()}
3135 @cindex @code{color_f()}
3139 ex color_d(const ex & a, const ex & b, const ex & c);
3140 ex color_f(const ex & a, const ex & b, const ex & c);
3143 create the symmetric and antisymmetric structure constants @math{d_abc} and
3144 @math{f_abc} which satisfy @math{@{T_a, T_b@} = 1/3 delta_ab + d_abc T_c}
3145 and @math{[T_a, T_b] = i f_abc T_c}.
3147 @cindex @code{color_h()}
3148 There's an additional function
3151 ex color_h(const ex & a, const ex & b, const ex & c);
3154 which returns the linear combination @samp{color_d(a, b, c)+I*color_f(a, b, c)}.
3156 The function @code{simplify_indexed()} performs some simplifications on
3157 expressions containing color objects:
3162 idx a(symbol("a"), 8), b(symbol("b"), 8), c(symbol("c"), 8),
3163 k(symbol("k"), 8), l(symbol("l"), 8);
3165 e = color_d(a, b, l) * color_f(a, b, k);
3166 cout << e.simplify_indexed() << endl;
3169 e = color_d(a, b, l) * color_d(a, b, k);
3170 cout << e.simplify_indexed() << endl;
3173 e = color_f(l, a, b) * color_f(a, b, k);
3174 cout << e.simplify_indexed() << endl;
3177 e = color_h(a, b, c) * color_h(a, b, c);
3178 cout << e.simplify_indexed() << endl;
3181 e = color_h(a, b, c) * color_T(b) * color_T(c);
3182 cout << e.simplify_indexed() << endl;
3185 e = color_h(a, b, c) * color_T(a) * color_T(b) * color_T(c);
3186 cout << e.simplify_indexed() << endl;
3189 e = color_T(k) * color_T(a) * color_T(b) * color_T(k);
3190 cout << e.simplify_indexed() << endl;
3191 // -> 1/4*delta.b.a*ONE-1/6*T.a*T.b
3195 @cindex @code{color_trace()}
3196 To calculate the trace of an expression containing color objects you use one
3200 ex color_trace(const ex & e, const std::set<unsigned char> & rls);
3201 ex color_trace(const ex & e, const lst & rll);
3202 ex color_trace(const ex & e, unsigned char rl = 0);
3205 These functions take the trace over all color @samp{T} objects in the
3206 specified set @code{rls} or list @code{rll} of representation labels, or the
3207 single label @code{rl}; @samp{T}s with other labels are left standing. For
3212 e = color_trace(4 * color_T(a) * color_T(b) * color_T(c));
3214 // -> -I*f.a.c.b+d.a.c.b
3219 @node Hash Maps, Methods and Functions, Non-commutative objects, Basic Concepts
3220 @c node-name, next, previous, up
3223 @cindex @code{exhashmap} (class)
3225 For your convenience, GiNaC offers the container template @code{exhashmap<T>}
3226 that can be used as a drop-in replacement for the STL
3227 @code{std::map<ex, T, ex_is_less>}, using hash tables to provide faster,
3228 typically constant-time, element look-up than @code{map<>}.
3230 @code{exhashmap<>} supports all @code{map<>} members and operations, with the
3231 following differences:
3235 no @code{lower_bound()} and @code{upper_bound()} methods
3237 no reverse iterators, no @code{rbegin()}/@code{rend()}
3239 no @code{operator<(exhashmap, exhashmap)}
3241 the comparison function object @code{key_compare} is hardcoded to
3244 the constructor @code{exhashmap(size_t n)} allows specifying the minimum
3245 initial hash table size (the actual table size after construction may be
3246 larger than the specified value)
3248 the method @code{size_t bucket_count()} returns the current size of the hash
3251 @code{insert()} and @code{erase()} operations invalidate all iterators
3255 @node Methods and Functions, Information About Expressions, Hash Maps, Top
3256 @c node-name, next, previous, up
3257 @chapter Methods and Functions
3260 In this chapter the most important algorithms provided by GiNaC will be
3261 described. Some of them are implemented as functions on expressions,
3262 others are implemented as methods provided by expression objects. If
3263 they are methods, there exists a wrapper function around it, so you can
3264 alternatively call it in a functional way as shown in the simple
3269 cout << "As method: " << sin(1).evalf() << endl;
3270 cout << "As function: " << evalf(sin(1)) << endl;
3274 @cindex @code{subs()}
3275 The general rule is that wherever methods accept one or more parameters
3276 (@var{arg1}, @var{arg2}, @dots{}) the order of arguments the function
3277 wrapper accepts is the same but preceded by the object to act on
3278 (@var{object}, @var{arg1}, @var{arg2}, @dots{}). This approach is the
3279 most natural one in an OO model but it may lead to confusion for MapleV
3280 users because where they would type @code{A:=x+1; subs(x=2,A);} GiNaC
3281 would require @code{A=x+1; subs(A,x==2);} (after proper declaration of
3282 @code{A} and @code{x}). On the other hand, since MapleV returns 3 on
3283 @code{A:=x^2+3; coeff(A,x,0);} (GiNaC: @code{A=pow(x,2)+3;
3284 coeff(A,x,0);}) it is clear that MapleV is not trying to be consistent
3285 here. Also, users of MuPAD will in most cases feel more comfortable
3286 with GiNaC's convention. All function wrappers are implemented
3287 as simple inline functions which just call the corresponding method and
3288 are only provided for users uncomfortable with OO who are dead set to
3289 avoid method invocations. Generally, nested function wrappers are much
3290 harder to read than a sequence of methods and should therefore be
3291 avoided if possible. On the other hand, not everything in GiNaC is a
3292 method on class @code{ex} and sometimes calling a function cannot be
3296 * Information About Expressions::
3297 * Numerical Evaluation::
3298 * Substituting Expressions::
3299 * Pattern Matching and Advanced Substitutions::
3300 * Applying a Function on Subexpressions::
3301 * Visitors and Tree Traversal::
3302 * Polynomial Arithmetic:: Working with polynomials.
3303 * Rational Expressions:: Working with rational functions.
3304 * Symbolic Differentiation::
3305 * Series Expansion:: Taylor and Laurent expansion.
3307 * Built-in Functions:: List of predefined mathematical functions.
3308 * Multiple polylogarithms::
3309 * Complex Conjugation::
3310 * Built-in Functions:: List of predefined mathematical functions.
3311 * Solving Linear Systems of Equations::
3312 * Input/Output:: Input and output of expressions.
3316 @node Information About Expressions, Numerical Evaluation, Methods and Functions, Methods and Functions
3317 @c node-name, next, previous, up
3318 @section Getting information about expressions
3320 @subsection Checking expression types
3321 @cindex @code{is_a<@dots{}>()}
3322 @cindex @code{is_exactly_a<@dots{}>()}
3323 @cindex @code{ex_to<@dots{}>()}
3324 @cindex Converting @code{ex} to other classes
3325 @cindex @code{info()}
3326 @cindex @code{return_type()}
3327 @cindex @code{return_type_tinfo()}
3329 Sometimes it's useful to check whether a given expression is a plain number,
3330 a sum, a polynomial with integer coefficients, or of some other specific type.
3331 GiNaC provides a couple of functions for this:
3334 bool is_a<T>(const ex & e);
3335 bool is_exactly_a<T>(const ex & e);
3336 bool ex::info(unsigned flag);
3337 unsigned ex::return_type() const;
3338 unsigned ex::return_type_tinfo() const;
3341 When the test made by @code{is_a<T>()} returns true, it is safe to call
3342 one of the functions @code{ex_to<T>()}, where @code{T} is one of the
3343 class names (@xref{The Class Hierarchy}, for a list of all classes). For
3344 example, assuming @code{e} is an @code{ex}:
3349 if (is_a<numeric>(e))
3350 numeric n = ex_to<numeric>(e);
3355 @code{is_a<T>(e)} allows you to check whether the top-level object of
3356 an expression @samp{e} is an instance of the GiNaC class @samp{T}
3357 (@xref{The Class Hierarchy}, for a list of all classes). This is most useful,
3358 e.g., for checking whether an expression is a number, a sum, or a product:
3365 is_a<numeric>(e1); // true
3366 is_a<numeric>(e2); // false
3367 is_a<add>(e1); // false
3368 is_a<add>(e2); // true
3369 is_a<mul>(e1); // false
3370 is_a<mul>(e2); // false
3374 In contrast, @code{is_exactly_a<T>(e)} allows you to check whether the
3375 top-level object of an expression @samp{e} is an instance of the GiNaC
3376 class @samp{T}, not including parent classes.
3378 The @code{info()} method is used for checking certain attributes of
3379 expressions. The possible values for the @code{flag} argument are defined
3380 in @file{ginac/flags.h}, the most important being explained in the following
3384 @multitable @columnfractions .30 .70
3385 @item @strong{Flag} @tab @strong{Returns true if the object is@dots{}}
3386 @item @code{numeric}
3387 @tab @dots{}a number (same as @code{is_a<numeric>(...)})
3389 @tab @dots{}a real integer, rational or float (i.e. is not complex)
3390 @item @code{rational}
3391 @tab @dots{}an exact rational number (integers are rational, too)
3392 @item @code{integer}
3393 @tab @dots{}a (non-complex) integer
3394 @item @code{crational}
3395 @tab @dots{}an exact (complex) rational number (such as @math{2/3+7/2*I})
3396 @item @code{cinteger}
3397 @tab @dots{}a (complex) integer (such as @math{2-3*I})
3398 @item @code{positive}
3399 @tab @dots{}not complex and greater than 0
3400 @item @code{negative}
3401 @tab @dots{}not complex and less than 0
3402 @item @code{nonnegative}
3403 @tab @dots{}not complex and greater than or equal to 0
3405 @tab @dots{}an integer greater than 0
3407 @tab @dots{}an integer less than 0
3408 @item @code{nonnegint}
3409 @tab @dots{}an integer greater than or equal to 0
3411 @tab @dots{}an even integer
3413 @tab @dots{}an odd integer
3415 @tab @dots{}a prime integer (probabilistic primality test)
3416 @item @code{relation}
3417 @tab @dots{}a relation (same as @code{is_a<relational>(...)})
3418 @item @code{relation_equal}
3419 @tab @dots{}a @code{==} relation
3420 @item @code{relation_not_equal}
3421 @tab @dots{}a @code{!=} relation
3422 @item @code{relation_less}
3423 @tab @dots{}a @code{<} relation
3424 @item @code{relation_less_or_equal}
3425 @tab @dots{}a @code{<=} relation
3426 @item @code{relation_greater}
3427 @tab @dots{}a @code{>} relation
3428 @item @code{relation_greater_or_equal}
3429 @tab @dots{}a @code{>=} relation
3431 @tab @dots{}a symbol (same as @code{is_a<symbol>(...)})
3433 @tab @dots{}a list (same as @code{is_a<lst>(...)})
3434 @item @code{polynomial}
3435 @tab @dots{}a polynomial (i.e. only consists of sums and products of numbers and symbols with positive integer powers)
3436 @item @code{integer_polynomial}
3437 @tab @dots{}a polynomial with (non-complex) integer coefficients
3438 @item @code{cinteger_polynomial}
3439 @tab @dots{}a polynomial with (possibly complex) integer coefficients (such as @math{2-3*I})
3440 @item @code{rational_polynomial}
3441 @tab @dots{}a polynomial with (non-complex) rational coefficients
3442 @item @code{crational_polynomial}
3443 @tab @dots{}a polynomial with (possibly complex) rational coefficients (such as @math{2/3+7/2*I})
3444 @item @code{rational_function}
3445 @tab @dots{}a rational function (@math{x+y}, @math{z/(x+y)})
3446 @item @code{algebraic}
3447 @tab @dots{}an algebraic object (@math{sqrt(2)}, @math{sqrt(x)-1})
3451 To determine whether an expression is commutative or non-commutative and if
3452 so, with which other expressions it would commutate, you use the methods
3453 @code{return_type()} and @code{return_type_tinfo()}. @xref{Non-commutative objects},
3454 for an explanation of these.
3457 @subsection Accessing subexpressions
3460 Many GiNaC classes, like @code{add}, @code{mul}, @code{lst}, and
3461 @code{function}, act as containers for subexpressions. For example, the
3462 subexpressions of a sum (an @code{add} object) are the individual terms,
3463 and the subexpressions of a @code{function} are the function's arguments.
3465 @cindex @code{nops()}
3467 GiNaC provides several ways of accessing subexpressions. The first way is to
3472 ex ex::op(size_t i);
3475 @code{nops()} determines the number of subexpressions (operands) contained
3476 in the expression, while @code{op(i)} returns the @code{i}-th
3477 (0..@code{nops()-1}) subexpression. In the case of a @code{power} object,
3478 @code{op(0)} will return the basis and @code{op(1)} the exponent. For
3479 @code{indexed} objects, @code{op(0)} is the base expression and @code{op(i)},
3480 @math{i>0} are the indices.
3483 @cindex @code{const_iterator}
3484 The second way to access subexpressions is via the STL-style random-access
3485 iterator class @code{const_iterator} and the methods
3488 const_iterator ex::begin();
3489 const_iterator ex::end();
3492 @code{begin()} returns an iterator referring to the first subexpression;
3493 @code{end()} returns an iterator which is one-past the last subexpression.
3494 If the expression has no subexpressions, then @code{begin() == end()}. These
3495 iterators can also be used in conjunction with non-modifying STL algorithms.
3497 Here is an example that (non-recursively) prints the subexpressions of a
3498 given expression in three different ways:
3505 for (size_t i = 0; i != e.nops(); ++i)
3506 cout << e.op(i) << endl;
3509 for (const_iterator i = e.begin(); i != e.end(); ++i)
3512 // with iterators and STL copy()
3513 std::copy(e.begin(), e.end(), std::ostream_iterator<ex>(cout, "\n"));
3517 @cindex @code{const_preorder_iterator}
3518 @cindex @code{const_postorder_iterator}
3519 @code{op()}/@code{nops()} and @code{const_iterator} only access an
3520 expression's immediate children. GiNaC provides two additional iterator
3521 classes, @code{const_preorder_iterator} and @code{const_postorder_iterator},
3522 that iterate over all objects in an expression tree, in preorder or postorder,
3523 respectively. They are STL-style forward iterators, and are created with the
3527 const_preorder_iterator ex::preorder_begin();
3528 const_preorder_iterator ex::preorder_end();
3529 const_postorder_iterator ex::postorder_begin();
3530 const_postorder_iterator ex::postorder_end();
3533 The following example illustrates the differences between
3534 @code{const_iterator}, @code{const_preorder_iterator}, and
3535 @code{const_postorder_iterator}:
3539 symbol A("A"), B("B"), C("C");
3540 ex e = lst(lst(A, B), C);
3542 std::copy(e.begin(), e.end(),
3543 std::ostream_iterator<ex>(cout, "\n"));
3547 std::copy(e.preorder_begin(), e.preorder_end(),
3548 std::ostream_iterator<ex>(cout, "\n"));
3555 std::copy(e.postorder_begin(), e.postorder_end(),
3556 std::ostream_iterator<ex>(cout, "\n"));
3565 @cindex @code{relational} (class)
3566 Finally, the left-hand side and right-hand side expressions of objects of
3567 class @code{relational} (and only of these) can also be accessed with the
3576 @subsection Comparing expressions
3577 @cindex @code{is_equal()}
3578 @cindex @code{is_zero()}
3580 Expressions can be compared with the usual C++ relational operators like
3581 @code{==}, @code{>}, and @code{<} but if the expressions contain symbols,
3582 the result is usually not determinable and the result will be @code{false},
3583 except in the case of the @code{!=} operator. You should also be aware that
3584 GiNaC will only do the most trivial test for equality (subtracting both
3585 expressions), so something like @code{(pow(x,2)+x)/x==x+1} will return
3588 Actually, if you construct an expression like @code{a == b}, this will be
3589 represented by an object of the @code{relational} class (@pxref{Relations})
3590 which is not evaluated until (explicitly or implicitly) cast to a @code{bool}.
3592 There are also two methods
3595 bool ex::is_equal(const ex & other);
3599 for checking whether one expression is equal to another, or equal to zero,
3603 @subsection Ordering expressions
3604 @cindex @code{ex_is_less} (class)
3605 @cindex @code{ex_is_equal} (class)
3606 @cindex @code{compare()}
3608 Sometimes it is necessary to establish a mathematically well-defined ordering
3609 on a set of arbitrary expressions, for example to use expressions as keys
3610 in a @code{std::map<>} container, or to bring a vector of expressions into
3611 a canonical order (which is done internally by GiNaC for sums and products).
3613 The operators @code{<}, @code{>} etc. described in the last section cannot
3614 be used for this, as they don't implement an ordering relation in the
3615 mathematical sense. In particular, they are not guaranteed to be
3616 antisymmetric: if @samp{a} and @samp{b} are different expressions, and
3617 @code{a < b} yields @code{false}, then @code{b < a} doesn't necessarily
3620 By default, STL classes and algorithms use the @code{<} and @code{==}
3621 operators to compare objects, which are unsuitable for expressions, but GiNaC
3622 provides two functors that can be supplied as proper binary comparison
3623 predicates to the STL:
3626 class ex_is_less : public std::binary_function<ex, ex, bool> @{
3628 bool operator()(const ex &lh, const ex &rh) const;
3631 class ex_is_equal : public std::binary_function<ex, ex, bool> @{
3633 bool operator()(const ex &lh, const ex &rh) const;
3637 For example, to define a @code{map} that maps expressions to strings you
3641 std::map<ex, std::string, ex_is_less> myMap;
3644 Omitting the @code{ex_is_less} template parameter will introduce spurious
3645 bugs because the map operates improperly.
3647 Other examples for the use of the functors:
3655 std::sort(v.begin(), v.end(), ex_is_less());
3657 // count the number of expressions equal to '1'
3658 unsigned num_ones = std::count_if(v.begin(), v.end(),
3659 std::bind2nd(ex_is_equal(), 1));
3662 The implementation of @code{ex_is_less} uses the member function
3665 int ex::compare(const ex & other) const;
3668 which returns @math{0} if @code{*this} and @code{other} are equal, @math{-1}
3669 if @code{*this} sorts before @code{other}, and @math{1} if @code{*this} sorts
3673 @node Numerical Evaluation, Substituting Expressions, Information About Expressions, Methods and Functions
3674 @c node-name, next, previous, up
3675 @section Numerical Evaluation
3676 @cindex @code{evalf()}
3678 GiNaC keeps algebraic expressions, numbers and constants in their exact form.
3679 To evaluate them using floating-point arithmetic you need to call
3682 ex ex::evalf(int level = 0) const;
3685 @cindex @code{Digits}
3686 The accuracy of the evaluation is controlled by the global object @code{Digits}
3687 which can be assigned an integer value. The default value of @code{Digits}
3688 is 17. @xref{Numbers}, for more information and examples.
3690 To evaluate an expression to a @code{double} floating-point number you can
3691 call @code{evalf()} followed by @code{numeric::to_double()}, like this:
3695 // Approximate sin(x/Pi)
3697 ex e = series(sin(x/Pi), x == 0, 6);
3699 // Evaluate numerically at x=0.1
3700 ex f = evalf(e.subs(x == 0.1));
3702 // ex_to<numeric> is an unsafe cast, so check the type first
3703 if (is_a<numeric>(f)) @{
3704 double d = ex_to<numeric>(f).to_double();
3713 @node Substituting Expressions, Pattern Matching and Advanced Substitutions, Numerical Evaluation, Methods and Functions
3714 @c node-name, next, previous, up
3715 @section Substituting expressions
3716 @cindex @code{subs()}
3718 Algebraic objects inside expressions can be replaced with arbitrary
3719 expressions via the @code{.subs()} method:
3722 ex ex::subs(const ex & e, unsigned options = 0);
3723 ex ex::subs(const exmap & m, unsigned options = 0);
3724 ex ex::subs(const lst & syms, const lst & repls, unsigned options = 0);
3727 In the first form, @code{subs()} accepts a relational of the form
3728 @samp{object == expression} or a @code{lst} of such relationals:
3732 symbol x("x"), y("y");
3734 ex e1 = 2*x^2-4*x+3;
3735 cout << "e1(7) = " << e1.subs(x == 7) << endl;
3739 cout << "e2(-2, 4) = " << e2.subs(lst(x == -2, y == 4)) << endl;