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 Generally, the top-level Makefile runs recursively to the
606 subdirectories. It is therefore safe to go into any subdirectory
607 (@code{doc/}, @code{ginsh/}, @dots{}) and simply type @code{make}
608 @var{target} there in case something went wrong.
611 @node Installing GiNaC, Basic Concepts, Building GiNaC, Installation
612 @c node-name, next, previous, up
613 @section Installing GiNaC
616 To install GiNaC on your system, simply type
622 As described in the section about configuration the files will be
623 installed in the following directories (the directories will be created
624 if they don't already exist):
629 @file{libginac.a} will go into @file{@var{PREFIX}/lib/} (or
630 @file{@var{LIBDIR}}) which defaults to @file{/usr/local/lib/}.
631 So will @file{libginac.so} unless the configure script was
632 given the option @option{--disable-shared}. The proper symlinks
633 will be established as well.
636 All the header files will be installed into @file{@var{PREFIX}/include/ginac/}
637 (or @file{@var{INCLUDEDIR}/ginac/}, if specified).
640 All documentation (HTML and Postscript) will be stuffed into
641 @file{@var{PREFIX}/share/doc/GiNaC/} (or
642 @file{@var{DATADIR}/doc/GiNaC/}, if @var{DATADIR} was specified).
646 For the sake of completeness we will list some other useful make
647 targets: @command{make clean} deletes all files generated by
648 @command{make}, i.e. all the object files. In addition @command{make
649 distclean} removes all files generated by the configuration and
650 @command{make maintainer-clean} goes one step further and deletes files
651 that may require special tools to rebuild (like the @command{libtool}
652 for instance). Finally @command{make uninstall} removes the installed
653 library, header files and documentation@footnote{Uninstallation does not
654 work after you have called @command{make distclean} since the
655 @file{Makefile} is itself generated by the configuration from
656 @file{Makefile.in} and hence deleted by @command{make distclean}. There
657 are two obvious ways out of this dilemma. First, you can run the
658 configuration again with the same @var{PREFIX} thus creating a
659 @file{Makefile} with a working @samp{uninstall} target. Second, you can
660 do it by hand since you now know where all the files went during
664 @node Basic Concepts, Expressions, Installing GiNaC, Top
665 @c node-name, next, previous, up
666 @chapter Basic Concepts
668 This chapter will describe the different fundamental objects that can be
669 handled by GiNaC. But before doing so, it is worthwhile introducing you
670 to the more commonly used class of expressions, representing a flexible
671 meta-class for storing all mathematical objects.
674 * Expressions:: The fundamental GiNaC class.
675 * Automatic evaluation:: Evaluation and canonicalization.
676 * Error handling:: How the library reports errors.
677 * The Class Hierarchy:: Overview of GiNaC's classes.
678 * Symbols:: Symbolic objects.
679 * Numbers:: Numerical objects.
680 * Constants:: Pre-defined constants.
681 * Fundamental containers:: Sums, products and powers.
682 * Lists:: Lists of expressions.
683 * Mathematical functions:: Mathematical functions.
684 * Relations:: Equality, Inequality and all that.
685 * Matrices:: Matrices.
686 * Indexed objects:: Handling indexed quantities.
687 * Non-commutative objects:: Algebras with non-commutative products.
688 * Hash Maps:: A faster alternative to std::map<>.
692 @node Expressions, Automatic evaluation, Basic Concepts, Basic Concepts
693 @c node-name, next, previous, up
695 @cindex expression (class @code{ex})
698 The most common class of objects a user deals with is the expression
699 @code{ex}, representing a mathematical object like a variable, number,
700 function, sum, product, etc@dots{} Expressions may be put together to form
701 new expressions, passed as arguments to functions, and so on. Here is a
702 little collection of valid expressions:
705 ex MyEx1 = 5; // simple number
706 ex MyEx2 = x + 2*y; // polynomial in x and y
707 ex MyEx3 = (x + 1)/(x - 1); // rational expression
708 ex MyEx4 = sin(x + 2*y) + 3*z + 41; // containing a function
709 ex MyEx5 = MyEx4 + 1; // similar to above
712 Expressions are handles to other more fundamental objects, that often
713 contain other expressions thus creating a tree of expressions
714 (@xref{Internal Structures}, for particular examples). Most methods on
715 @code{ex} therefore run top-down through such an expression tree. For
716 example, the method @code{has()} scans recursively for occurrences of
717 something inside an expression. Thus, if you have declared @code{MyEx4}
718 as in the example above @code{MyEx4.has(y)} will find @code{y} inside
719 the argument of @code{sin} and hence return @code{true}.
721 The next sections will outline the general picture of GiNaC's class
722 hierarchy and describe the classes of objects that are handled by
725 @subsection Note: Expressions and STL containers
727 GiNaC expressions (@code{ex} objects) have value semantics (they can be
728 assigned, reassigned and copied like integral types) but the operator
729 @code{<} doesn't provide a well-defined ordering on them. In STL-speak,
730 expressions are @samp{Assignable} but not @samp{LessThanComparable}.
732 This implies that in order to use expressions in sorted containers such as
733 @code{std::map<>} and @code{std::set<>} you have to supply a suitable
734 comparison predicate. GiNaC provides such a predicate, called
735 @code{ex_is_less}. For example, a set of expressions should be defined
736 as @code{std::set<ex, ex_is_less>}.
738 Unsorted containers such as @code{std::vector<>} and @code{std::list<>}
739 don't pose a problem. A @code{std::vector<ex>} works as expected.
741 @xref{Information About Expressions}, for more about comparing and ordering
745 @node Automatic evaluation, Error handling, Expressions, Basic Concepts
746 @c node-name, next, previous, up
747 @section Automatic evaluation and canonicalization of expressions
750 GiNaC performs some automatic transformations on expressions, to simplify
751 them and put them into a canonical form. Some examples:
754 ex MyEx1 = 2*x - 1 + x; // 3*x-1
755 ex MyEx2 = x - x; // 0
756 ex MyEx3 = cos(2*Pi); // 1
757 ex MyEx4 = x*y/x; // y
760 This behavior is usually referred to as @dfn{automatic} or @dfn{anonymous
761 evaluation}. GiNaC only performs transformations that are
765 at most of complexity
773 algebraically correct, possibly except for a set of measure zero (e.g.
774 @math{x/x} is transformed to @math{1} although this is incorrect for @math{x=0})
777 There are two types of automatic transformations in GiNaC that may not
778 behave in an entirely obvious way at first glance:
782 The terms of sums and products (and some other things like the arguments of
783 symmetric functions, the indices of symmetric tensors etc.) are re-ordered
784 into a canonical form that is deterministic, but not lexicographical or in
785 any other way easy to guess (it almost always depends on the number and
786 order of the symbols you define). However, constructing the same expression
787 twice, either implicitly or explicitly, will always result in the same
790 Expressions of the form 'number times sum' are automatically expanded (this
791 has to do with GiNaC's internal representation of sums and products). For
794 ex MyEx5 = 2*(x + y); // 2*x+2*y
795 ex MyEx6 = z*(x + y); // z*(x+y)
799 The general rule is that when you construct expressions, GiNaC automatically
800 creates them in canonical form, which might differ from the form you typed in
801 your program. This may create some awkward looking output (@samp{-y+x} instead
802 of @samp{x-y}) but allows for more efficient operation and usually yields
803 some immediate simplifications.
805 @cindex @code{eval()}
806 Internally, the anonymous evaluator in GiNaC is implemented by the methods
809 ex ex::eval(int level = 0) const;
810 ex basic::eval(int level = 0) const;
813 but unless you are extending GiNaC with your own classes or functions, there
814 should never be any reason to call them explicitly. All GiNaC methods that
815 transform expressions, like @code{subs()} or @code{normal()}, automatically
816 re-evaluate their results.
819 @node Error handling, The Class Hierarchy, Automatic evaluation, Basic Concepts
820 @c node-name, next, previous, up
821 @section Error handling
823 @cindex @code{pole_error} (class)
825 GiNaC reports run-time errors by throwing C++ exceptions. All exceptions
826 generated by GiNaC are subclassed from the standard @code{exception} class
827 defined in the @file{<stdexcept>} header. In addition to the predefined
828 @code{logic_error}, @code{domain_error}, @code{out_of_range},
829 @code{invalid_argument}, @code{runtime_error}, @code{range_error} and
830 @code{overflow_error} types, GiNaC also defines a @code{pole_error}
831 exception that gets thrown when trying to evaluate a mathematical function
834 The @code{pole_error} class has a member function
837 int pole_error::degree() const;
840 that returns the order of the singularity (or 0 when the pole is
841 logarithmic or the order is undefined).
843 When using GiNaC it is useful to arrange for exceptions to be caught in
844 the main program even if you don't want to do any special error handling.
845 Otherwise whenever an error occurs in GiNaC, it will be delegated to the
846 default exception handler of your C++ compiler's run-time system which
847 usually only aborts the program without giving any information what went
850 Here is an example for a @code{main()} function that catches and prints
851 exceptions generated by GiNaC:
856 #include <ginac/ginac.h>
858 using namespace GiNaC;
866 @} catch (exception &p) @{
867 cerr << p.what() << endl;
875 @node The Class Hierarchy, Symbols, Error handling, Basic Concepts
876 @c node-name, next, previous, up
877 @section The Class Hierarchy
879 GiNaC's class hierarchy consists of several classes representing
880 mathematical objects, all of which (except for @code{ex} and some
881 helpers) are internally derived from one abstract base class called
882 @code{basic}. You do not have to deal with objects of class
883 @code{basic}, instead you'll be dealing with symbols, numbers,
884 containers of expressions and so on.
888 To get an idea about what kinds of symbolic composites may be built we
889 have a look at the most important classes in the class hierarchy and
890 some of the relations among the classes:
892 @image{classhierarchy}
894 The abstract classes shown here (the ones without drop-shadow) are of no
895 interest for the user. They are used internally in order to avoid code
896 duplication if two or more classes derived from them share certain
897 features. An example is @code{expairseq}, a container for a sequence of
898 pairs each consisting of one expression and a number (@code{numeric}).
899 What @emph{is} visible to the user are the derived classes @code{add}
900 and @code{mul}, representing sums and products. @xref{Internal
901 Structures}, where these two classes are described in more detail. The
902 following table shortly summarizes what kinds of mathematical objects
903 are stored in the different classes:
906 @multitable @columnfractions .22 .78
907 @item @code{symbol} @tab Algebraic symbols @math{a}, @math{x}, @math{y}@dots{}
908 @item @code{constant} @tab Constants like
915 @item @code{numeric} @tab All kinds of numbers, @math{42}, @math{7/3*I}, @math{3.14159}@dots{}
916 @item @code{add} @tab Sums like @math{x+y} or @math{a-(2*b)+3}
917 @item @code{mul} @tab Products like @math{x*y} or @math{2*a^2*(x+y+z)/b}
918 @item @code{ncmul} @tab Products of non-commutative objects
919 @item @code{power} @tab Exponentials such as @math{x^2}, @math{a^b},
924 @code{sqrt(}@math{2}@code{)}
927 @item @code{pseries} @tab Power Series, e.g. @math{x-1/6*x^3+1/120*x^5+O(x^7)}
928 @item @code{function} @tab A symbolic function like
935 @item @code{lst} @tab Lists of expressions @{@math{x}, @math{2*y}, @math{3+z}@}
936 @item @code{matrix} @tab @math{m}x@math{n} matrices of expressions
937 @item @code{relational} @tab A relation like the identity @math{x}@code{==}@math{y}
938 @item @code{indexed} @tab Indexed object like @math{A_ij}
939 @item @code{tensor} @tab Special tensor like the delta and metric tensors
940 @item @code{idx} @tab Index of an indexed object
941 @item @code{varidx} @tab Index with variance
942 @item @code{spinidx} @tab Index with variance and dot (used in Weyl-van-der-Waerden spinor formalism)
943 @item @code{wildcard} @tab Wildcard for pattern matching
944 @item @code{structure} @tab Template for user-defined classes
949 @node Symbols, Numbers, The Class Hierarchy, Basic Concepts
950 @c node-name, next, previous, up
952 @cindex @code{symbol} (class)
953 @cindex hierarchy of classes
956 Symbolic indeterminates, or @dfn{symbols} for short, are for symbolic
957 manipulation what atoms are for chemistry.
959 A typical symbol definition looks like this:
964 This definition actually contains three very different things:
966 @item a C++ variable named @code{x}
967 @item a @code{symbol} object stored in this C++ variable; this object
968 represents the symbol in a GiNaC expression
969 @item the string @code{"x"} which is the name of the symbol, used (almost)
970 exclusively for printing expressions holding the symbol
973 Symbols have an explicit name, supplied as a string during construction,
974 because in C++, variable names can't be used as values, and the C++ compiler
975 throws them away during compilation.
977 It is possible to omit the symbol name in the definition:
982 In this case, GiNaC will assign the symbol an internal, unique name of the
983 form @code{symbolNNN}. This won't affect the usability of the symbol but
984 the output of your calculations will become more readable if you give your
985 symbols sensible names (for intermediate expressions that are only used
986 internally such anonymous symbols can be quite useful, however).
988 Now, here is one important property of GiNaC that differentiates it from
989 other computer algebra programs you may have used: GiNaC does @emph{not} use
990 the names of symbols to tell them apart, but a (hidden) serial number that
991 is unique for each newly created @code{symbol} object. In you want to use
992 one and the same symbol in different places in your program, you must only
993 create one @code{symbol} object and pass that around. If you create another
994 symbol, even if it has the same name, GiNaC will treat it as a different
1011 // prints "x^6" which looks right, but...
1013 cout << e.degree(x) << endl;
1014 // ...this doesn't work. The symbol "x" here is different from the one
1015 // in f() and in the expression returned by f(). Consequently, it
1020 One possibility to ensure that @code{f()} and @code{main()} use the same
1021 symbol is to pass the symbol as an argument to @code{f()}:
1023 ex f(int n, const ex & x)
1032 // Now, f() uses the same symbol.
1035 cout << e.degree(x) << endl;
1036 // prints "6", as expected
1040 Another possibility would be to define a global symbol @code{x} that is used
1041 by both @code{f()} and @code{main()}. If you are using global symbols and
1042 multiple compilation units you must take special care, however. Suppose
1043 that you have a header file @file{globals.h} in your program that defines
1044 a @code{symbol x("x");}. In this case, every unit that includes
1045 @file{globals.h} would also get its own definition of @code{x} (because
1046 header files are just inlined into the source code by the C++ preprocessor),
1047 and hence you would again end up with multiple equally-named, but different,
1048 symbols. Instead, the @file{globals.h} header should only contain a
1049 @emph{declaration} like @code{extern symbol x;}, with the definition of
1050 @code{x} moved into a C++ source file such as @file{globals.cpp}.
1052 A different approach to ensuring that symbols used in different parts of
1053 your program are identical is to create them with a @emph{factory} function
1056 const symbol & get_symbol(const string & s)
1058 static map<string, symbol> directory;
1059 map<string, symbol>::iterator i = directory.find(s);
1060 if (i != directory.end())
1063 return directory.insert(make_pair(s, symbol(s))).first->second;
1067 This function returns one newly constructed symbol for each name that is
1068 passed in, and it returns the same symbol when called multiple times with
1069 the same name. Using this symbol factory, we can rewrite our example like
1074 return pow(get_symbol("x"), n);
1081 // Both calls of get_symbol("x") yield the same symbol.
1082 cout << e.degree(get_symbol("x")) << endl;
1087 Instead of creating symbols from strings we could also have
1088 @code{get_symbol()} take, for example, an integer number as its argument.
1089 In this case, we would probably want to give the generated symbols names
1090 that include this number, which can be accomplished with the help of an
1091 @code{ostringstream}.
1093 In general, if you're getting weird results from GiNaC such as an expression
1094 @samp{x-x} that is not simplified to zero, you should check your symbol
1097 As we said, the names of symbols primarily serve for purposes of expression
1098 output. But there are actually two instances where GiNaC uses the names for
1099 identifying symbols: When constructing an expression from a string, and when
1100 recreating an expression from an archive (@pxref{Input/Output}).
1102 In addition to its name, a symbol may contain a special string that is used
1105 symbol x("x", "\\Box");
1108 This creates a symbol that is printed as "@code{x}" in normal output, but
1109 as "@code{\Box}" in LaTeX code (@xref{Input/Output}, for more
1110 information about the different output formats of expressions in GiNaC).
1111 GiNaC automatically creates proper LaTeX code for symbols having names of
1112 greek letters (@samp{alpha}, @samp{mu}, etc.).
1114 @cindex @code{subs()}
1115 Symbols in GiNaC can't be assigned values. If you need to store results of
1116 calculations and give them a name, use C++ variables of type @code{ex}.
1117 If you want to replace a symbol in an expression with something else, you
1118 can invoke the expression's @code{.subs()} method
1119 (@pxref{Substituting Expressions}).
1121 @cindex @code{realsymbol()}
1122 By default, symbols are expected to stand in for complex values, i.e. they live
1123 in the complex domain. As a consequence, operations like complex conjugation,
1124 for example (@pxref{Complex Conjugation}), do @emph{not} evaluate if applied
1125 to such symbols. Likewise @code{log(exp(x))} does not evaluate to @code{x},
1126 because of the unknown imaginary part of @code{x}.
1127 On the other hand, if you are sure that your symbols will hold only real values, you
1128 would like to have such functions evaluated. Therefore GiNaC allows you to specify
1129 the domain of the symbol. Instead of @code{symbol x("x");} you can write
1130 @code{realsymbol x("x");} to tell GiNaC that @code{x} stands in for real values.
1133 @node Numbers, Constants, Symbols, Basic Concepts
1134 @c node-name, next, previous, up
1136 @cindex @code{numeric} (class)
1142 For storing numerical things, GiNaC uses Bruno Haible's library CLN.
1143 The classes therein serve as foundation classes for GiNaC. CLN stands
1144 for Class Library for Numbers or alternatively for Common Lisp Numbers.
1145 In order to find out more about CLN's internals, the reader is referred to
1146 the documentation of that library. @inforef{Introduction, , cln}, for
1147 more information. Suffice to say that it is by itself build on top of
1148 another library, the GNU Multiple Precision library GMP, which is an
1149 extremely fast library for arbitrary long integers and rationals as well
1150 as arbitrary precision floating point numbers. It is very commonly used
1151 by several popular cryptographic applications. CLN extends GMP by
1152 several useful things: First, it introduces the complex number field
1153 over either reals (i.e. floating point numbers with arbitrary precision)
1154 or rationals. Second, it automatically converts rationals to integers
1155 if the denominator is unity and complex numbers to real numbers if the
1156 imaginary part vanishes and also correctly treats algebraic functions.
1157 Third it provides good implementations of state-of-the-art algorithms
1158 for all trigonometric and hyperbolic functions as well as for
1159 calculation of some useful constants.
1161 The user can construct an object of class @code{numeric} in several
1162 ways. The following example shows the four most important constructors.
1163 It uses construction from C-integer, construction of fractions from two
1164 integers, construction from C-float and construction from a string:
1168 #include <ginac/ginac.h>
1169 using namespace GiNaC;
1173 numeric two = 2; // exact integer 2
1174 numeric r(2,3); // exact fraction 2/3
1175 numeric e(2.71828); // floating point number
1176 numeric p = "3.14159265358979323846"; // constructor from string
1177 // Trott's constant in scientific notation:
1178 numeric trott("1.0841015122311136151E-2");
1180 std::cout << two*p << std::endl; // floating point 6.283...
1185 @cindex complex numbers
1186 The imaginary unit in GiNaC is a predefined @code{numeric} object with the
1191 numeric z1 = 2-3*I; // exact complex number 2-3i
1192 numeric z2 = 5.9+1.6*I; // complex floating point number
1196 It may be tempting to construct fractions by writing @code{numeric r(3/2)}.
1197 This would, however, call C's built-in operator @code{/} for integers
1198 first and result in a numeric holding a plain integer 1. @strong{Never
1199 use the operator @code{/} on integers} unless you know exactly what you
1200 are doing! Use the constructor from two integers instead, as shown in
1201 the example above. Writing @code{numeric(1)/2} may look funny but works
1204 @cindex @code{Digits}
1206 We have seen now the distinction between exact numbers and floating
1207 point numbers. Clearly, the user should never have to worry about
1208 dynamically created exact numbers, since their `exactness' always
1209 determines how they ought to be handled, i.e. how `long' they are. The
1210 situation is different for floating point numbers. Their accuracy is
1211 controlled by one @emph{global} variable, called @code{Digits}. (For
1212 those readers who know about Maple: it behaves very much like Maple's
1213 @code{Digits}). All objects of class numeric that are constructed from
1214 then on will be stored with a precision matching that number of decimal
1219 #include <ginac/ginac.h>
1220 using namespace std;
1221 using namespace GiNaC;
1225 numeric three(3.0), one(1.0);
1226 numeric x = one/three;
1228 cout << "in " << Digits << " digits:" << endl;
1230 cout << Pi.evalf() << endl;
1242 The above example prints the following output to screen:
1246 0.33333333333333333334
1247 3.1415926535897932385
1249 0.33333333333333333333333333333333333333333333333333333333333333333334
1250 3.1415926535897932384626433832795028841971693993751058209749445923078
1254 Note that the last number is not necessarily rounded as you would
1255 naively expect it to be rounded in the decimal system. But note also,
1256 that in both cases you got a couple of extra digits. This is because
1257 numbers are internally stored by CLN as chunks of binary digits in order
1258 to match your machine's word size and to not waste precision. Thus, on
1259 architectures with different word size, the above output might even
1260 differ with regard to actually computed digits.
1262 It should be clear that objects of class @code{numeric} should be used
1263 for constructing numbers or for doing arithmetic with them. The objects
1264 one deals with most of the time are the polymorphic expressions @code{ex}.
1266 @subsection Tests on numbers
1268 Once you have declared some numbers, assigned them to expressions and
1269 done some arithmetic with them it is frequently desired to retrieve some
1270 kind of information from them like asking whether that number is
1271 integer, rational, real or complex. For those cases GiNaC provides
1272 several useful methods. (Internally, they fall back to invocations of
1273 certain CLN functions.)
1275 As an example, let's construct some rational number, multiply it with
1276 some multiple of its denominator and test what comes out:
1280 #include <ginac/ginac.h>
1281 using namespace std;
1282 using namespace GiNaC;
1284 // some very important constants:
1285 const numeric twentyone(21);
1286 const numeric ten(10);
1287 const numeric five(5);
1291 numeric answer = twentyone;
1294 cout << answer.is_integer() << endl; // false, it's 21/5
1296 cout << answer.is_integer() << endl; // true, it's 42 now!
1300 Note that the variable @code{answer} is constructed here as an integer
1301 by @code{numeric}'s copy constructor but in an intermediate step it
1302 holds a rational number represented as integer numerator and integer
1303 denominator. When multiplied by 10, the denominator becomes unity and
1304 the result is automatically converted to a pure integer again.
1305 Internally, the underlying CLN is responsible for this behavior and we
1306 refer the reader to CLN's documentation. Suffice to say that
1307 the same behavior applies to complex numbers as well as return values of
1308 certain functions. Complex numbers are automatically converted to real
1309 numbers if the imaginary part becomes zero. The full set of tests that
1310 can be applied is listed in the following table.
1313 @multitable @columnfractions .30 .70
1314 @item @strong{Method} @tab @strong{Returns true if the object is@dots{}}
1315 @item @code{.is_zero()}
1316 @tab @dots{}equal to zero
1317 @item @code{.is_positive()}
1318 @tab @dots{}not complex and greater than 0
1319 @item @code{.is_integer()}
1320 @tab @dots{}a (non-complex) integer
1321 @item @code{.is_pos_integer()}
1322 @tab @dots{}an integer and greater than 0
1323 @item @code{.is_nonneg_integer()}
1324 @tab @dots{}an integer and greater equal 0
1325 @item @code{.is_even()}
1326 @tab @dots{}an even integer
1327 @item @code{.is_odd()}
1328 @tab @dots{}an odd integer
1329 @item @code{.is_prime()}
1330 @tab @dots{}a prime integer (probabilistic primality test)
1331 @item @code{.is_rational()}
1332 @tab @dots{}an exact rational number (integers are rational, too)
1333 @item @code{.is_real()}
1334 @tab @dots{}a real integer, rational or float (i.e. is not complex)
1335 @item @code{.is_cinteger()}
1336 @tab @dots{}a (complex) integer (such as @math{2-3*I})
1337 @item @code{.is_crational()}
1338 @tab @dots{}an exact (complex) rational number (such as @math{2/3+7/2*I})
1342 @subsection Numeric functions
1344 The following functions can be applied to @code{numeric} objects and will be
1345 evaluated immediately:
1348 @multitable @columnfractions .30 .70
1349 @item @strong{Name} @tab @strong{Function}
1350 @item @code{inverse(z)}
1351 @tab returns @math{1/z}
1352 @cindex @code{inverse()} (numeric)
1353 @item @code{pow(a, b)}
1354 @tab exponentiation @math{a^b}
1357 @item @code{real(z)}
1359 @cindex @code{real()}
1360 @item @code{imag(z)}
1362 @cindex @code{imag()}
1363 @item @code{csgn(z)}
1364 @tab complex sign (returns an @code{int})
1365 @item @code{numer(z)}
1366 @tab numerator of rational or complex rational number
1367 @item @code{denom(z)}
1368 @tab denominator of rational or complex rational number
1369 @item @code{sqrt(z)}
1371 @item @code{isqrt(n)}
1372 @tab integer square root
1373 @cindex @code{isqrt()}
1380 @item @code{asin(z)}
1382 @item @code{acos(z)}
1384 @item @code{atan(z)}
1385 @tab inverse tangent
1386 @item @code{atan(y, x)}
1387 @tab inverse tangent with two arguments
1388 @item @code{sinh(z)}
1389 @tab hyperbolic sine
1390 @item @code{cosh(z)}
1391 @tab hyperbolic cosine
1392 @item @code{tanh(z)}
1393 @tab hyperbolic tangent
1394 @item @code{asinh(z)}
1395 @tab inverse hyperbolic sine
1396 @item @code{acosh(z)}
1397 @tab inverse hyperbolic cosine
1398 @item @code{atanh(z)}
1399 @tab inverse hyperbolic tangent
1401 @tab exponential function
1403 @tab natural logarithm
1406 @item @code{zeta(z)}
1407 @tab Riemann's zeta function
1408 @item @code{tgamma(z)}
1410 @item @code{lgamma(z)}
1411 @tab logarithm of gamma function
1413 @tab psi (digamma) function
1414 @item @code{psi(n, z)}
1415 @tab derivatives of psi function (polygamma functions)
1416 @item @code{factorial(n)}
1417 @tab factorial function @math{n!}
1418 @item @code{doublefactorial(n)}
1419 @tab double factorial function @math{n!!}
1420 @cindex @code{doublefactorial()}
1421 @item @code{binomial(n, k)}
1422 @tab binomial coefficients
1423 @item @code{bernoulli(n)}
1424 @tab Bernoulli numbers
1425 @cindex @code{bernoulli()}
1426 @item @code{fibonacci(n)}
1427 @tab Fibonacci numbers
1428 @cindex @code{fibonacci()}
1429 @item @code{mod(a, b)}
1430 @tab modulus in positive representation (in the range @code{[0, abs(b)-1]} with the sign of b, or zero)
1431 @cindex @code{mod()}
1432 @item @code{smod(a, b)}
1433 @tab modulus in symmetric representation (in the range @code{[-iquo(abs(b)-1, 2), iquo(abs(b), 2)]})
1434 @cindex @code{smod()}
1435 @item @code{irem(a, b)}
1436 @tab integer remainder (has the sign of @math{a}, or is zero)
1437 @cindex @code{irem()}
1438 @item @code{irem(a, b, q)}
1439 @tab integer remainder and quotient, @code{irem(a, b, q) == a-q*b}
1440 @item @code{iquo(a, b)}
1441 @tab integer quotient
1442 @cindex @code{iquo()}
1443 @item @code{iquo(a, b, r)}
1444 @tab integer quotient and remainder, @code{r == a-iquo(a, b)*b}
1445 @item @code{gcd(a, b)}
1446 @tab greatest common divisor
1447 @item @code{lcm(a, b)}
1448 @tab least common multiple
1452 Most of these functions are also available as symbolic functions that can be
1453 used in expressions (@pxref{Mathematical functions}) or, like @code{gcd()},
1454 as polynomial algorithms.
1456 @subsection Converting numbers
1458 Sometimes it is desirable to convert a @code{numeric} object back to a
1459 built-in arithmetic type (@code{int}, @code{double}, etc.). The @code{numeric}
1460 class provides a couple of methods for this purpose:
1462 @cindex @code{to_int()}
1463 @cindex @code{to_long()}
1464 @cindex @code{to_double()}
1465 @cindex @code{to_cl_N()}
1467 int numeric::to_int() const;
1468 long numeric::to_long() const;
1469 double numeric::to_double() const;
1470 cln::cl_N numeric::to_cl_N() const;
1473 @code{to_int()} and @code{to_long()} only work when the number they are
1474 applied on is an exact integer. Otherwise the program will halt with a
1475 message like @samp{Not a 32-bit integer}. @code{to_double()} applied on a
1476 rational number will return a floating-point approximation. Both
1477 @code{to_int()/to_long()} and @code{to_double()} discard the imaginary
1478 part of complex numbers.
1481 @node Constants, Fundamental containers, Numbers, Basic Concepts
1482 @c node-name, next, previous, up
1484 @cindex @code{constant} (class)
1487 @cindex @code{Catalan}
1488 @cindex @code{Euler}
1489 @cindex @code{evalf()}
1490 Constants behave pretty much like symbols except that they return some
1491 specific number when the method @code{.evalf()} is called.
1493 The predefined known constants are:
1496 @multitable @columnfractions .14 .30 .56
1497 @item @strong{Name} @tab @strong{Common Name} @tab @strong{Numerical Value (to 35 digits)}
1499 @tab Archimedes' constant
1500 @tab 3.14159265358979323846264338327950288
1501 @item @code{Catalan}
1502 @tab Catalan's constant
1503 @tab 0.91596559417721901505460351493238411
1505 @tab Euler's (or Euler-Mascheroni) constant
1506 @tab 0.57721566490153286060651209008240243
1511 @node Fundamental containers, Lists, Constants, Basic Concepts
1512 @c node-name, next, previous, up
1513 @section Sums, products and powers
1517 @cindex @code{power}
1519 Simple rational expressions are written down in GiNaC pretty much like
1520 in other CAS or like expressions involving numerical variables in C.
1521 The necessary operators @code{+}, @code{-}, @code{*} and @code{/} have
1522 been overloaded to achieve this goal. When you run the following
1523 code snippet, the constructor for an object of type @code{mul} is
1524 automatically called to hold the product of @code{a} and @code{b} and
1525 then the constructor for an object of type @code{add} is called to hold
1526 the sum of that @code{mul} object and the number one:
1530 symbol a("a"), b("b");
1535 @cindex @code{pow()}
1536 For exponentiation, you have already seen the somewhat clumsy (though C-ish)
1537 statement @code{pow(x,2);} to represent @code{x} squared. This direct
1538 construction is necessary since we cannot safely overload the constructor
1539 @code{^} in C++ to construct a @code{power} object. If we did, it would
1540 have several counterintuitive and undesired effects:
1544 Due to C's operator precedence, @code{2*x^2} would be parsed as @code{(2*x)^2}.
1546 Due to the binding of the operator @code{^}, @code{x^a^b} would result in
1547 @code{(x^a)^b}. This would be confusing since most (though not all) other CAS
1548 interpret this as @code{x^(a^b)}.
1550 Also, expressions involving integer exponents are very frequently used,
1551 which makes it even more dangerous to overload @code{^} since it is then
1552 hard to distinguish between the semantics as exponentiation and the one
1553 for exclusive or. (It would be embarrassing to return @code{1} where one
1554 has requested @code{2^3}.)
1557 @cindex @command{ginsh}
1558 All effects are contrary to mathematical notation and differ from the
1559 way most other CAS handle exponentiation, therefore overloading @code{^}
1560 is ruled out for GiNaC's C++ part. The situation is different in
1561 @command{ginsh}, there the exponentiation-@code{^} exists. (Also note
1562 that the other frequently used exponentiation operator @code{**} does
1563 not exist at all in C++).
1565 To be somewhat more precise, objects of the three classes described
1566 here, are all containers for other expressions. An object of class
1567 @code{power} is best viewed as a container with two slots, one for the
1568 basis, one for the exponent. All valid GiNaC expressions can be
1569 inserted. However, basic transformations like simplifying
1570 @code{pow(pow(x,2),3)} to @code{x^6} automatically are only performed
1571 when this is mathematically possible. If we replace the outer exponent
1572 three in the example by some symbols @code{a}, the simplification is not
1573 safe and will not be performed, since @code{a} might be @code{1/2} and
1576 Objects of type @code{add} and @code{mul} are containers with an
1577 arbitrary number of slots for expressions to be inserted. Again, simple
1578 and safe simplifications are carried out like transforming
1579 @code{3*x+4-x} to @code{2*x+4}.
1582 @node Lists, Mathematical functions, Fundamental containers, Basic Concepts
1583 @c node-name, next, previous, up
1584 @section Lists of expressions
1585 @cindex @code{lst} (class)
1587 @cindex @code{nops()}
1589 @cindex @code{append()}
1590 @cindex @code{prepend()}
1591 @cindex @code{remove_first()}
1592 @cindex @code{remove_last()}
1593 @cindex @code{remove_all()}
1595 The GiNaC class @code{lst} serves for holding a @dfn{list} of arbitrary
1596 expressions. They are not as ubiquitous as in many other computer algebra
1597 packages, but are sometimes used to supply a variable number of arguments of
1598 the same type to GiNaC methods such as @code{subs()} and some @code{matrix}
1599 constructors, so you should have a basic understanding of them.
1601 Lists can be constructed by assigning a comma-separated sequence of
1606 symbol x("x"), y("y");
1609 // now, l is a list holding the expressions 'x', '2', 'y', and 'x+y',
1614 There are also constructors that allow direct creation of lists of up to
1615 16 expressions, which is often more convenient but slightly less efficient:
1619 // This produces the same list 'l' as above:
1620 // lst l(x, 2, y, x+y);
1621 // lst l = lst(x, 2, y, x+y);
1625 Use the @code{nops()} method to determine the size (number of expressions) of
1626 a list and the @code{op()} method or the @code{[]} operator to access
1627 individual elements:
1631 cout << l.nops() << endl; // prints '4'
1632 cout << l.op(2) << " " << l[0] << endl; // prints 'y x'
1636 As with the standard @code{list<T>} container, accessing random elements of a
1637 @code{lst} is generally an operation of order @math{O(N)}. Faster read-only
1638 sequential access to the elements of a list is possible with the
1639 iterator types provided by the @code{lst} class:
1642 typedef ... lst::const_iterator;
1643 typedef ... lst::const_reverse_iterator;
1644 lst::const_iterator lst::begin() const;
1645 lst::const_iterator lst::end() const;
1646 lst::const_reverse_iterator lst::rbegin() const;
1647 lst::const_reverse_iterator lst::rend() const;
1650 For example, to print the elements of a list individually you can use:
1655 for (lst::const_iterator i = l.begin(); i != l.end(); ++i)
1660 which is one order faster than
1665 for (size_t i = 0; i < l.nops(); ++i)
1666 cout << l.op(i) << endl;
1670 These iterators also allow you to use some of the algorithms provided by
1671 the C++ standard library:
1675 // print the elements of the list (requires #include <iterator>)
1676 std::copy(l.begin(), l.end(), ostream_iterator<ex>(cout, "\n"));
1678 // sum up the elements of the list (requires #include <numeric>)
1679 ex sum = std::accumulate(l.begin(), l.end(), ex(0));
1680 cout << sum << endl; // prints '2+2*x+2*y'
1684 @code{lst} is one of the few GiNaC classes that allow in-place modifications
1685 (the only other one is @code{matrix}). You can modify single elements:
1689 l[1] = 42; // l is now @{x, 42, y, x+y@}
1690 l.let_op(1) = 7; // l is now @{x, 7, y, x+y@}
1694 You can append or prepend an expression to a list with the @code{append()}
1695 and @code{prepend()} methods:
1699 l.append(4*x); // l is now @{x, 7, y, x+y, 4*x@}
1700 l.prepend(0); // l is now @{0, x, 7, y, x+y, 4*x@}
1704 You can remove the first or last element of a list with @code{remove_first()}
1705 and @code{remove_last()}:
1709 l.remove_first(); // l is now @{x, 7, y, x+y, 4*x@}
1710 l.remove_last(); // l is now @{x, 7, y, x+y@}
1714 You can remove all the elements of a list with @code{remove_all()}:
1718 l.remove_all(); // l is now empty
1722 You can bring the elements of a list into a canonical order with @code{sort()}:
1731 // l1 and l2 are now equal
1735 Finally, you can remove all but the first element of consecutive groups of
1736 elements with @code{unique()}:
1741 l3 = x, 2, 2, 2, y, x+y, y+x;
1742 l3.unique(); // l3 is now @{x, 2, y, x+y@}
1747 @node Mathematical functions, Relations, Lists, Basic Concepts
1748 @c node-name, next, previous, up
1749 @section Mathematical functions
1750 @cindex @code{function} (class)
1751 @cindex trigonometric function
1752 @cindex hyperbolic function
1754 There are quite a number of useful functions hard-wired into GiNaC. For
1755 instance, all trigonometric and hyperbolic functions are implemented
1756 (@xref{Built-in Functions}, for a complete list).
1758 These functions (better called @emph{pseudofunctions}) are all objects
1759 of class @code{function}. They accept one or more expressions as
1760 arguments and return one expression. If the arguments are not
1761 numerical, the evaluation of the function may be halted, as it does in
1762 the next example, showing how a function returns itself twice and
1763 finally an expression that may be really useful:
1765 @cindex Gamma function
1766 @cindex @code{subs()}
1769 symbol x("x"), y("y");
1771 cout << tgamma(foo) << endl;
1772 // -> tgamma(x+(1/2)*y)
1773 ex bar = foo.subs(y==1);
1774 cout << tgamma(bar) << endl;
1776 ex foobar = bar.subs(x==7);
1777 cout << tgamma(foobar) << endl;
1778 // -> (135135/128)*Pi^(1/2)
1782 Besides evaluation most of these functions allow differentiation, series
1783 expansion and so on. Read the next chapter in order to learn more about
1786 It must be noted that these pseudofunctions are created by inline
1787 functions, where the argument list is templated. This means that
1788 whenever you call @code{GiNaC::sin(1)} it is equivalent to
1789 @code{sin(ex(1))} and will therefore not result in a floating point
1790 number. Unless of course the function prototype is explicitly
1791 overridden -- which is the case for arguments of type @code{numeric}
1792 (not wrapped inside an @code{ex}). Hence, in order to obtain a floating
1793 point number of class @code{numeric} you should call
1794 @code{sin(numeric(1))}. This is almost the same as calling
1795 @code{sin(1).evalf()} except that the latter will return a numeric
1796 wrapped inside an @code{ex}.
1799 @node Relations, Matrices, Mathematical functions, Basic Concepts
1800 @c node-name, next, previous, up
1802 @cindex @code{relational} (class)
1804 Sometimes, a relation holding between two expressions must be stored
1805 somehow. The class @code{relational} is a convenient container for such
1806 purposes. A relation is by definition a container for two @code{ex} and
1807 a relation between them that signals equality, inequality and so on.
1808 They are created by simply using the C++ operators @code{==}, @code{!=},
1809 @code{<}, @code{<=}, @code{>} and @code{>=} between two expressions.
1811 @xref{Mathematical functions}, for examples where various applications
1812 of the @code{.subs()} method show how objects of class relational are
1813 used as arguments. There they provide an intuitive syntax for
1814 substitutions. They are also used as arguments to the @code{ex::series}
1815 method, where the left hand side of the relation specifies the variable
1816 to expand in and the right hand side the expansion point. They can also
1817 be used for creating systems of equations that are to be solved for
1818 unknown variables. But the most common usage of objects of this class
1819 is rather inconspicuous in statements of the form @code{if
1820 (expand(pow(a+b,2))==a*a+2*a*b+b*b) @{...@}}. Here, an implicit
1821 conversion from @code{relational} to @code{bool} takes place. Note,
1822 however, that @code{==} here does not perform any simplifications, hence
1823 @code{expand()} must be called explicitly.
1826 @node Matrices, Indexed objects, Relations, Basic Concepts
1827 @c node-name, next, previous, up
1829 @cindex @code{matrix} (class)
1831 A @dfn{matrix} is a two-dimensional array of expressions. The elements of a
1832 matrix with @math{m} rows and @math{n} columns are accessed with two
1833 @code{unsigned} indices, the first one in the range 0@dots{}@math{m-1}, the
1834 second one in the range 0@dots{}@math{n-1}.
1836 There are a couple of ways to construct matrices, with or without preset
1837 elements. The constructor
1840 matrix::matrix(unsigned r, unsigned c);
1843 creates a matrix with @samp{r} rows and @samp{c} columns with all elements
1846 The fastest way to create a matrix with preinitialized elements is to assign
1847 a list of comma-separated expressions to an empty matrix (see below for an
1848 example). But you can also specify the elements as a (flat) list with
1851 matrix::matrix(unsigned r, unsigned c, const lst & l);
1856 @cindex @code{lst_to_matrix()}
1858 ex lst_to_matrix(const lst & l);
1861 constructs a matrix from a list of lists, each list representing a matrix row.
1863 There is also a set of functions for creating some special types of
1866 @cindex @code{diag_matrix()}
1867 @cindex @code{unit_matrix()}
1868 @cindex @code{symbolic_matrix()}
1870 ex diag_matrix(const lst & l);
1871 ex unit_matrix(unsigned x);
1872 ex unit_matrix(unsigned r, unsigned c);
1873 ex symbolic_matrix(unsigned r, unsigned c, const string & base_name);
1874 ex symbolic_matrix(unsigned r, unsigned c, const string & base_name, const string & tex_base_name);
1877 @code{diag_matrix()} constructs a diagonal matrix given the list of diagonal
1878 elements. @code{unit_matrix()} creates an @samp{x} by @samp{x} (or @samp{r}
1879 by @samp{c}) unit matrix. And finally, @code{symbolic_matrix} constructs a
1880 matrix filled with newly generated symbols made of the specified base name
1881 and the position of each element in the matrix.
1883 Matrix elements can be accessed and set using the parenthesis (function call)
1887 const ex & matrix::operator()(unsigned r, unsigned c) const;
1888 ex & matrix::operator()(unsigned r, unsigned c);
1891 It is also possible to access the matrix elements in a linear fashion with
1892 the @code{op()} method. But C++-style subscripting with square brackets
1893 @samp{[]} is not available.
1895 Here are a couple of examples for constructing matrices:
1899 symbol a("a"), b("b");
1913 cout << matrix(2, 2, lst(a, 0, 0, b)) << endl;
1916 cout << lst_to_matrix(lst(lst(a, 0), lst(0, b))) << endl;
1919 cout << diag_matrix(lst(a, b)) << endl;
1922 cout << unit_matrix(3) << endl;
1923 // -> [[1,0,0],[0,1,0],[0,0,1]]
1925 cout << symbolic_matrix(2, 3, "x") << endl;
1926 // -> [[x00,x01,x02],[x10,x11,x12]]
1930 @cindex @code{transpose()}
1931 There are three ways to do arithmetic with matrices. The first (and most
1932 direct one) is to use the methods provided by the @code{matrix} class:
1935 matrix matrix::add(const matrix & other) const;
1936 matrix matrix::sub(const matrix & other) const;
1937 matrix matrix::mul(const matrix & other) const;
1938 matrix matrix::mul_scalar(const ex & other) const;
1939 matrix matrix::pow(const ex & expn) const;
1940 matrix matrix::transpose() const;
1943 All of these methods return the result as a new matrix object. Here is an
1944 example that calculates @math{A*B-2*C} for three matrices @math{A}, @math{B}
1949 matrix A(2, 2), B(2, 2), C(2, 2);
1957 matrix result = A.mul(B).sub(C.mul_scalar(2));
1958 cout << result << endl;
1959 // -> [[-13,-6],[1,2]]
1964 @cindex @code{evalm()}
1965 The second (and probably the most natural) way is to construct an expression
1966 containing matrices with the usual arithmetic operators and @code{pow()}.
1967 For efficiency reasons, expressions with sums, products and powers of
1968 matrices are not automatically evaluated in GiNaC. You have to call the
1972 ex ex::evalm() const;
1975 to obtain the result:
1982 // -> [[1,2],[3,4]]*[[-1,0],[2,1]]-2*[[8,4],[2,1]]
1983 cout << e.evalm() << endl;
1984 // -> [[-13,-6],[1,2]]
1989 The non-commutativity of the product @code{A*B} in this example is
1990 automatically recognized by GiNaC. There is no need to use a special
1991 operator here. @xref{Non-commutative objects}, for more information about
1992 dealing with non-commutative expressions.
1994 Finally, you can work with indexed matrices and call @code{simplify_indexed()}
1995 to perform the arithmetic:
2000 idx i(symbol("i"), 2), j(symbol("j"), 2), k(symbol("k"), 2);
2001 e = indexed(A, i, k) * indexed(B, k, j) - 2 * indexed(C, i, j);
2003 // -> -2*[[8,4],[2,1]].i.j+[[-1,0],[2,1]].k.j*[[1,2],[3,4]].i.k
2004 cout << e.simplify_indexed() << endl;
2005 // -> [[-13,-6],[1,2]].i.j
2009 Using indices is most useful when working with rectangular matrices and
2010 one-dimensional vectors because you don't have to worry about having to
2011 transpose matrices before multiplying them. @xref{Indexed objects}, for
2012 more information about using matrices with indices, and about indices in
2015 The @code{matrix} class provides a couple of additional methods for
2016 computing determinants, traces, characteristic polynomials and ranks:
2018 @cindex @code{determinant()}
2019 @cindex @code{trace()}
2020 @cindex @code{charpoly()}
2021 @cindex @code{rank()}
2023 ex matrix::determinant(unsigned algo=determinant_algo::automatic) const;
2024 ex matrix::trace() const;
2025 ex matrix::charpoly(const ex & lambda) const;
2026 unsigned matrix::rank() const;
2029 The @samp{algo} argument of @code{determinant()} allows to select
2030 between different algorithms for calculating the determinant. The
2031 asymptotic speed (as parametrized by the matrix size) can greatly differ
2032 between those algorithms, depending on the nature of the matrix'
2033 entries. The possible values are defined in the @file{flags.h} header
2034 file. By default, GiNaC uses a heuristic to automatically select an
2035 algorithm that is likely (but not guaranteed) to give the result most
2038 @cindex @code{inverse()} (matrix)
2039 @cindex @code{solve()}
2040 Matrices may also be inverted using the @code{ex matrix::inverse()}
2041 method and linear systems may be solved with:
2044 matrix matrix::solve(const matrix & vars, const matrix & rhs, unsigned algo=solve_algo::automatic) const;
2047 Assuming the matrix object this method is applied on is an @code{m}
2048 times @code{n} matrix, then @code{vars} must be a @code{n} times
2049 @code{p} matrix of symbolic indeterminates and @code{rhs} a @code{m}
2050 times @code{p} matrix. The returned matrix then has dimension @code{n}
2051 times @code{p} and in the case of an underdetermined system will still
2052 contain some of the indeterminates from @code{vars}. If the system is
2053 overdetermined, an exception is thrown.
2056 @node Indexed objects, Non-commutative objects, Matrices, Basic Concepts
2057 @c node-name, next, previous, up
2058 @section Indexed objects
2060 GiNaC allows you to handle expressions containing general indexed objects in
2061 arbitrary spaces. It is also able to canonicalize and simplify such
2062 expressions and perform symbolic dummy index summations. There are a number
2063 of predefined indexed objects provided, like delta and metric tensors.
2065 There are few restrictions placed on indexed objects and their indices and
2066 it is easy to construct nonsense expressions, but our intention is to
2067 provide a general framework that allows you to implement algorithms with
2068 indexed quantities, getting in the way as little as possible.
2070 @cindex @code{idx} (class)
2071 @cindex @code{indexed} (class)
2072 @subsection Indexed quantities and their indices
2074 Indexed expressions in GiNaC are constructed of two special types of objects,
2075 @dfn{index objects} and @dfn{indexed objects}.
2079 @cindex contravariant
2082 @item Index objects are of class @code{idx} or a subclass. Every index has
2083 a @dfn{value} and a @dfn{dimension} (which is the dimension of the space
2084 the index lives in) which can both be arbitrary expressions but are usually
2085 a number or a simple symbol. In addition, indices of class @code{varidx} have
2086 a @dfn{variance} (they can be co- or contravariant), and indices of class
2087 @code{spinidx} have a variance and can be @dfn{dotted} or @dfn{undotted}.
2089 @item Indexed objects are of class @code{indexed} or a subclass. They
2090 contain a @dfn{base expression} (which is the expression being indexed), and
2091 one or more indices.
2095 @strong{Note:} when printing expressions, covariant indices and indices
2096 without variance are denoted @samp{.i} while contravariant indices are
2097 denoted @samp{~i}. Dotted indices have a @samp{*} in front of the index
2098 value. In the following, we are going to use that notation in the text so
2099 instead of @math{A^i_jk} we will write @samp{A~i.j.k}. Index dimensions are
2100 not visible in the output.
2102 A simple example shall illustrate the concepts:
2106 #include <ginac/ginac.h>
2107 using namespace std;
2108 using namespace GiNaC;
2112 symbol i_sym("i"), j_sym("j");
2113 idx i(i_sym, 3), j(j_sym, 3);
2116 cout << indexed(A, i, j) << endl;
2118 cout << index_dimensions << indexed(A, i, j) << endl;
2120 cout << dflt; // reset cout to default output format (dimensions hidden)
2124 The @code{idx} constructor takes two arguments, the index value and the
2125 index dimension. First we define two index objects, @code{i} and @code{j},
2126 both with the numeric dimension 3. The value of the index @code{i} is the
2127 symbol @code{i_sym} (which prints as @samp{i}) and the value of the index
2128 @code{j} is the symbol @code{j_sym} (which prints as @samp{j}). Next we
2129 construct an expression containing one indexed object, @samp{A.i.j}. It has
2130 the symbol @code{A} as its base expression and the two indices @code{i} and
2133 The dimensions of indices are normally not visible in the output, but one
2134 can request them to be printed with the @code{index_dimensions} manipulator,
2137 Note the difference between the indices @code{i} and @code{j} which are of
2138 class @code{idx}, and the index values which are the symbols @code{i_sym}
2139 and @code{j_sym}. The indices of indexed objects cannot directly be symbols
2140 or numbers but must be index objects. For example, the following is not
2141 correct and will raise an exception:
2144 symbol i("i"), j("j");
2145 e = indexed(A, i, j); // ERROR: indices must be of type idx
2148 You can have multiple indexed objects in an expression, index values can
2149 be numeric, and index dimensions symbolic:
2153 symbol B("B"), dim("dim");
2154 cout << 4 * indexed(A, i)
2155 + indexed(B, idx(j_sym, 4), idx(2, 3), idx(i_sym, dim)) << endl;
2160 @code{B} has a 4-dimensional symbolic index @samp{k}, a 3-dimensional numeric
2161 index of value 2, and a symbolic index @samp{i} with the symbolic dimension
2162 @samp{dim}. Note that GiNaC doesn't automatically notify you that the free
2163 indices of @samp{A} and @samp{B} in the sum don't match (you have to call
2164 @code{simplify_indexed()} for that, see below).
2166 In fact, base expressions, index values and index dimensions can be
2167 arbitrary expressions:
2171 cout << indexed(A+B, idx(2*i_sym+1, dim/2)) << endl;
2176 It's also possible to construct nonsense like @samp{Pi.sin(x)}. You will not
2177 get an error message from this but you will probably not be able to do
2178 anything useful with it.
2180 @cindex @code{get_value()}
2181 @cindex @code{get_dimension()}
2185 ex idx::get_value();
2186 ex idx::get_dimension();
2189 return the value and dimension of an @code{idx} object. If you have an index
2190 in an expression, such as returned by calling @code{.op()} on an indexed
2191 object, you can get a reference to the @code{idx} object with the function
2192 @code{ex_to<idx>()} on the expression.
2194 There are also the methods
2197 bool idx::is_numeric();
2198 bool idx::is_symbolic();
2199 bool idx::is_dim_numeric();
2200 bool idx::is_dim_symbolic();
2203 for checking whether the value and dimension are numeric or symbolic
2204 (non-numeric). Using the @code{info()} method of an index (see @ref{Information
2205 About Expressions}) returns information about the index value.
2207 @cindex @code{varidx} (class)
2208 If you need co- and contravariant indices, use the @code{varidx} class:
2212 symbol mu_sym("mu"), nu_sym("nu");
2213 varidx mu(mu_sym, 4), nu(nu_sym, 4); // default is contravariant ~mu, ~nu
2214 varidx mu_co(mu_sym, 4, true); // covariant index .mu
2216 cout << indexed(A, mu, nu) << endl;
2218 cout << indexed(A, mu_co, nu) << endl;
2220 cout << indexed(A, mu.toggle_variance(), nu) << endl;
2225 A @code{varidx} is an @code{idx} with an additional flag that marks it as
2226 co- or contravariant. The default is a contravariant (upper) index, but
2227 this can be overridden by supplying a third argument to the @code{varidx}
2228 constructor. The two methods
2231 bool varidx::is_covariant();
2232 bool varidx::is_contravariant();
2235 allow you to check the variance of a @code{varidx} object (use @code{ex_to<varidx>()}
2236 to get the object reference from an expression). There's also the very useful
2240 ex varidx::toggle_variance();
2243 which makes a new index with the same value and dimension but the opposite
2244 variance. By using it you only have to define the index once.
2246 @cindex @code{spinidx} (class)
2247 The @code{spinidx} class provides dotted and undotted variant indices, as
2248 used in the Weyl-van-der-Waerden spinor formalism:
2252 symbol K("K"), C_sym("C"), D_sym("D");
2253 spinidx C(C_sym, 2), D(D_sym); // default is 2-dimensional,
2254 // contravariant, undotted
2255 spinidx C_co(C_sym, 2, true); // covariant index
2256 spinidx D_dot(D_sym, 2, false, true); // contravariant, dotted
2257 spinidx D_co_dot(D_sym, 2, true, true); // covariant, dotted
2259 cout << indexed(K, C, D) << endl;
2261 cout << indexed(K, C_co, D_dot) << endl;
2263 cout << indexed(K, D_co_dot, D) << endl;
2268 A @code{spinidx} is a @code{varidx} with an additional flag that marks it as
2269 dotted or undotted. The default is undotted but this can be overridden by
2270 supplying a fourth argument to the @code{spinidx} constructor. The two
2274 bool spinidx::is_dotted();
2275 bool spinidx::is_undotted();
2278 allow you to check whether or not a @code{spinidx} object is dotted (use
2279 @code{ex_to<spinidx>()} to get the object reference from an expression).
2280 Finally, the two methods
2283 ex spinidx::toggle_dot();
2284 ex spinidx::toggle_variance_dot();
2287 create a new index with the same value and dimension but opposite dottedness
2288 and the same or opposite variance.
2290 @subsection Substituting indices
2292 @cindex @code{subs()}
2293 Sometimes you will want to substitute one symbolic index with another
2294 symbolic or numeric index, for example when calculating one specific element
2295 of a tensor expression. This is done with the @code{.subs()} method, as it
2296 is done for symbols (see @ref{Substituting Expressions}).
2298 You have two possibilities here. You can either substitute the whole index
2299 by another index or expression:
2303 ex e = indexed(A, mu_co);
2304 cout << e << " becomes " << e.subs(mu_co == nu) << endl;
2305 // -> A.mu becomes A~nu
2306 cout << e << " becomes " << e.subs(mu_co == varidx(0, 4)) << endl;
2307 // -> A.mu becomes A~0
2308 cout << e << " becomes " << e.subs(mu_co == 0) << endl;
2309 // -> A.mu becomes A.0
2313 The third example shows that trying to replace an index with something that
2314 is not an index will substitute the index value instead.
2316 Alternatively, you can substitute the @emph{symbol} of a symbolic index by
2321 ex e = indexed(A, mu_co);
2322 cout << e << " becomes " << e.subs(mu_sym == nu_sym) << endl;
2323 // -> A.mu becomes A.nu
2324 cout << e << " becomes " << e.subs(mu_sym == 0) << endl;
2325 // -> A.mu becomes A.0
2329 As you see, with the second method only the value of the index will get
2330 substituted. Its other properties, including its dimension, remain unchanged.
2331 If you want to change the dimension of an index you have to substitute the
2332 whole index by another one with the new dimension.
2334 Finally, substituting the base expression of an indexed object works as
2339 ex e = indexed(A, mu_co);
2340 cout << e << " becomes " << e.subs(A == A+B) << endl;
2341 // -> A.mu becomes (B+A).mu
2345 @subsection Symmetries
2346 @cindex @code{symmetry} (class)
2347 @cindex @code{sy_none()}
2348 @cindex @code{sy_symm()}
2349 @cindex @code{sy_anti()}
2350 @cindex @code{sy_cycl()}
2352 Indexed objects can have certain symmetry properties with respect to their
2353 indices. Symmetries are specified as a tree of objects of class @code{symmetry}
2354 that is constructed with the helper functions
2357 symmetry sy_none(...);
2358 symmetry sy_symm(...);
2359 symmetry sy_anti(...);
2360 symmetry sy_cycl(...);
2363 @code{sy_none()} stands for no symmetry, @code{sy_symm()} and @code{sy_anti()}
2364 specify fully symmetric or antisymmetric, respectively, and @code{sy_cycl()}
2365 represents a cyclic symmetry. Each of these functions accepts up to four
2366 arguments which can be either symmetry objects themselves or unsigned integer
2367 numbers that represent an index position (counting from 0). A symmetry
2368 specification that consists of only a single @code{sy_symm()}, @code{sy_anti()}
2369 or @code{sy_cycl()} with no arguments specifies the respective symmetry for
2372 Here are some examples of symmetry definitions:
2377 e = indexed(A, i, j);
2378 e = indexed(A, sy_none(), i, j); // equivalent
2379 e = indexed(A, sy_none(0, 1), i, j); // equivalent
2381 // Symmetric in all three indices:
2382 e = indexed(A, sy_symm(), i, j, k);
2383 e = indexed(A, sy_symm(0, 1, 2), i, j, k); // equivalent
2384 e = indexed(A, sy_symm(2, 0, 1), i, j, k); // same symmetry, but yields a
2385 // different canonical order
2387 // Symmetric in the first two indices only:
2388 e = indexed(A, sy_symm(0, 1), i, j, k);
2389 e = indexed(A, sy_none(sy_symm(0, 1), 2), i, j, k); // equivalent
2391 // Antisymmetric in the first and last index only (index ranges need not
2393 e = indexed(A, sy_anti(0, 2), i, j, k);
2394 e = indexed(A, sy_none(sy_anti(0, 2), 1), i, j, k); // equivalent
2396 // An example of a mixed symmetry: antisymmetric in the first two and
2397 // last two indices, symmetric when swapping the first and last index
2398 // pairs (like the Riemann curvature tensor):
2399 e = indexed(A, sy_symm(sy_anti(0, 1), sy_anti(2, 3)), i, j, k, l);
2401 // Cyclic symmetry in all three indices:
2402 e = indexed(A, sy_cycl(), i, j, k);
2403 e = indexed(A, sy_cycl(0, 1, 2), i, j, k); // equivalent
2405 // The following examples are invalid constructions that will throw
2406 // an exception at run time.
2408 // An index may not appear multiple times:
2409 e = indexed(A, sy_symm(0, 0, 1), i, j, k); // ERROR
2410 e = indexed(A, sy_none(sy_symm(0, 1), sy_anti(0, 2)), i, j, k); // ERROR
2412 // Every child of sy_symm(), sy_anti() and sy_cycl() must refer to the
2413 // same number of indices:
2414 e = indexed(A, sy_symm(sy_anti(0, 1), 2), i, j, k); // ERROR
2416 // And of course, you cannot specify indices which are not there:
2417 e = indexed(A, sy_symm(0, 1, 2, 3), i, j, k); // ERROR
2421 If you need to specify more than four indices, you have to use the
2422 @code{.add()} method of the @code{symmetry} class. For example, to specify
2423 full symmetry in the first six indices you would write
2424 @code{sy_symm(0, 1, 2, 3).add(4).add(5)}.
2426 If an indexed object has a symmetry, GiNaC will automatically bring the
2427 indices into a canonical order which allows for some immediate simplifications:
2431 cout << indexed(A, sy_symm(), i, j)
2432 + indexed(A, sy_symm(), j, i) << endl;
2434 cout << indexed(B, sy_anti(), i, j)
2435 + indexed(B, sy_anti(), j, i) << endl;
2437 cout << indexed(B, sy_anti(), i, j, k)
2438 - indexed(B, sy_anti(), j, k, i) << endl;
2443 @cindex @code{get_free_indices()}
2445 @subsection Dummy indices
2447 GiNaC treats certain symbolic index pairs as @dfn{dummy indices} meaning
2448 that a summation over the index range is implied. Symbolic indices which are
2449 not dummy indices are called @dfn{free indices}. Numeric indices are neither
2450 dummy nor free indices.
2452 To be recognized as a dummy index pair, the two indices must be of the same
2453 class and their value must be the same single symbol (an index like
2454 @samp{2*n+1} is never a dummy index). If the indices are of class
2455 @code{varidx} they must also be of opposite variance; if they are of class
2456 @code{spinidx} they must be both dotted or both undotted.
2458 The method @code{.get_free_indices()} returns a vector containing the free
2459 indices of an expression. It also checks that the free indices of the terms
2460 of a sum are consistent:
2464 symbol A("A"), B("B"), C("C");
2466 symbol i_sym("i"), j_sym("j"), k_sym("k"), l_sym("l");
2467 idx i(i_sym, 3), j(j_sym, 3), k(k_sym, 3), l(l_sym, 3);
2469 ex e = indexed(A, i, j) * indexed(B, j, k) + indexed(C, k, l, i, l);
2470 cout << exprseq(e.get_free_indices()) << endl;
2472 // 'j' and 'l' are dummy indices
2474 symbol mu_sym("mu"), nu_sym("nu"), rho_sym("rho"), sigma_sym("sigma");
2475 varidx mu(mu_sym, 4), nu(nu_sym, 4), rho(rho_sym, 4), sigma(sigma_sym, 4);
2477 e = indexed(A, mu, nu) * indexed(B, nu.toggle_variance(), rho)
2478 + indexed(C, mu, sigma, rho, sigma.toggle_variance());
2479 cout << exprseq(e.get_free_indices()) << endl;
2481 // 'nu' is a dummy index, but 'sigma' is not
2483 e = indexed(A, mu, mu);
2484 cout << exprseq(e.get_free_indices()) << endl;
2486 // 'mu' is not a dummy index because it appears twice with the same
2489 e = indexed(A, mu, nu) + 42;
2490 cout << exprseq(e.get_free_indices()) << endl; // ERROR
2491 // this will throw an exception:
2492 // "add::get_free_indices: inconsistent indices in sum"
2496 @cindex @code{simplify_indexed()}
2497 @subsection Simplifying indexed expressions
2499 In addition to the few automatic simplifications that GiNaC performs on
2500 indexed expressions (such as re-ordering the indices of symmetric tensors
2501 and calculating traces and convolutions of matrices and predefined tensors)
2505 ex ex::simplify_indexed();
2506 ex ex::simplify_indexed(const scalar_products & sp);
2509 that performs some more expensive operations:
2512 @item it checks the consistency of free indices in sums in the same way
2513 @code{get_free_indices()} does
2514 @item it tries to give dummy indices that appear in different terms of a sum
2515 the same name to allow simplifications like @math{a_i*b_i-a_j*b_j=0}
2516 @item it (symbolically) calculates all possible dummy index summations/contractions
2517 with the predefined tensors (this will be explained in more detail in the
2519 @item it detects contractions that vanish for symmetry reasons, for example
2520 the contraction of a symmetric and a totally antisymmetric tensor
2521 @item as a special case of dummy index summation, it can replace scalar products
2522 of two tensors with a user-defined value
2525 The last point is done with the help of the @code{scalar_products} class
2526 which is used to store scalar products with known values (this is not an
2527 arithmetic class, you just pass it to @code{simplify_indexed()}):
2531 symbol A("A"), B("B"), C("C"), i_sym("i");
2535 sp.add(A, B, 0); // A and B are orthogonal
2536 sp.add(A, C, 0); // A and C are orthogonal
2537 sp.add(A, A, 4); // A^2 = 4 (A has length 2)
2539 e = indexed(A + B, i) * indexed(A + C, i);
2541 // -> (B+A).i*(A+C).i
2543 cout << e.expand(expand_options::expand_indexed).simplify_indexed(sp)
2549 The @code{scalar_products} object @code{sp} acts as a storage for the
2550 scalar products added to it with the @code{.add()} method. This method
2551 takes three arguments: the two expressions of which the scalar product is
2552 taken, and the expression to replace it with. After @code{sp.add(A, B, 0)},
2553 @code{simplify_indexed()} will replace all scalar products of indexed
2554 objects that have the symbols @code{A} and @code{B} as base expressions
2555 with the single value 0. The number, type and dimension of the indices
2556 don't matter; @samp{A~mu~nu*B.mu.nu} would also be replaced by 0.
2558 @cindex @code{expand()}
2559 The example above also illustrates a feature of the @code{expand()} method:
2560 if passed the @code{expand_indexed} option it will distribute indices
2561 over sums, so @samp{(A+B).i} becomes @samp{A.i+B.i}.
2563 @cindex @code{tensor} (class)
2564 @subsection Predefined tensors
2566 Some frequently used special tensors such as the delta, epsilon and metric
2567 tensors are predefined in GiNaC. They have special properties when
2568 contracted with other tensor expressions and some of them have constant
2569 matrix representations (they will evaluate to a number when numeric
2570 indices are specified).
2572 @cindex @code{delta_tensor()}
2573 @subsubsection Delta tensor
2575 The delta tensor takes two indices, is symmetric and has the matrix
2576 representation @code{diag(1, 1, 1, ...)}. It is constructed by the function
2577 @code{delta_tensor()}:
2581 symbol A("A"), B("B");
2583 idx i(symbol("i"), 3), j(symbol("j"), 3),
2584 k(symbol("k"), 3), l(symbol("l"), 3);
2586 ex e = indexed(A, i, j) * indexed(B, k, l)
2587 * delta_tensor(i, k) * delta_tensor(j, l) << endl;
2588 cout << e.simplify_indexed() << endl;
2591 cout << delta_tensor(i, i) << endl;
2596 @cindex @code{metric_tensor()}
2597 @subsubsection General metric tensor
2599 The function @code{metric_tensor()} creates a general symmetric metric
2600 tensor with two indices that can be used to raise/lower tensor indices. The
2601 metric tensor is denoted as @samp{g} in the output and if its indices are of
2602 mixed variance it is automatically replaced by a delta tensor:
2608 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4);
2610 ex e = metric_tensor(mu, nu) * indexed(A, nu.toggle_variance(), rho);
2611 cout << e.simplify_indexed() << endl;
2614 e = delta_tensor(mu, nu.toggle_variance()) * metric_tensor(nu, rho);
2615 cout << e.simplify_indexed() << endl;
2618 e = metric_tensor(mu.toggle_variance(), nu.toggle_variance())
2619 * metric_tensor(nu, rho);
2620 cout << e.simplify_indexed() << endl;
2623 e = metric_tensor(nu.toggle_variance(), rho.toggle_variance())
2624 * metric_tensor(mu, nu) * (delta_tensor(mu.toggle_variance(), rho)
2625 + indexed(A, mu.toggle_variance(), rho));
2626 cout << e.simplify_indexed() << endl;
2631 @cindex @code{lorentz_g()}
2632 @subsubsection Minkowski metric tensor
2634 The Minkowski metric tensor is a special metric tensor with a constant
2635 matrix representation which is either @code{diag(1, -1, -1, ...)} (negative
2636 signature, the default) or @code{diag(-1, 1, 1, ...)} (positive signature).
2637 It is created with the function @code{lorentz_g()} (although it is output as
2642 varidx mu(symbol("mu"), 4);
2644 e = delta_tensor(varidx(0, 4), mu.toggle_variance())
2645 * lorentz_g(mu, varidx(0, 4)); // negative signature
2646 cout << e.simplify_indexed() << endl;
2649 e = delta_tensor(varidx(0, 4), mu.toggle_variance())
2650 * lorentz_g(mu, varidx(0, 4), true); // positive signature
2651 cout << e.simplify_indexed() << endl;
2656 @cindex @code{spinor_metric()}
2657 @subsubsection Spinor metric tensor
2659 The function @code{spinor_metric()} creates an antisymmetric tensor with
2660 two indices that is used to raise/lower indices of 2-component spinors.
2661 It is output as @samp{eps}:
2667 spinidx A(symbol("A")), B(symbol("B")), C(symbol("C"));
2668 ex A_co = A.toggle_variance(), B_co = B.toggle_variance();
2670 e = spinor_metric(A, B) * indexed(psi, B_co);
2671 cout << e.simplify_indexed() << endl;
2674 e = spinor_metric(A, B) * indexed(psi, A_co);
2675 cout << e.simplify_indexed() << endl;
2678 e = spinor_metric(A_co, B_co) * indexed(psi, B);
2679 cout << e.simplify_indexed() << endl;
2682 e = spinor_metric(A_co, B_co) * indexed(psi, A);
2683 cout << e.simplify_indexed() << endl;
2686 e = spinor_metric(A_co, B_co) * spinor_metric(A, B);
2687 cout << e.simplify_indexed() << endl;
2690 e = spinor_metric(A_co, B_co) * spinor_metric(B, C);
2691 cout << e.simplify_indexed() << endl;
2696 The matrix representation of the spinor metric is @code{[[0, 1], [-1, 0]]}.
2698 @cindex @code{epsilon_tensor()}
2699 @cindex @code{lorentz_eps()}
2700 @subsubsection Epsilon tensor
2702 The epsilon tensor is totally antisymmetric, its number of indices is equal
2703 to the dimension of the index space (the indices must all be of the same
2704 numeric dimension), and @samp{eps.1.2.3...} (resp. @samp{eps~0~1~2...}) is
2705 defined to be 1. Its behavior with indices that have a variance also
2706 depends on the signature of the metric. Epsilon tensors are output as
2709 There are three functions defined to create epsilon tensors in 2, 3 and 4
2713 ex epsilon_tensor(const ex & i1, const ex & i2);
2714 ex epsilon_tensor(const ex & i1, const ex & i2, const ex & i3);
2715 ex lorentz_eps(const ex & i1, const ex & i2, const ex & i3, const ex & i4, bool pos_sig = false);
2718 The first two functions create an epsilon tensor in 2 or 3 Euclidean
2719 dimensions, the last function creates an epsilon tensor in a 4-dimensional
2720 Minkowski space (the last @code{bool} argument specifies whether the metric
2721 has negative or positive signature, as in the case of the Minkowski metric
2726 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4),
2727 sig(symbol("sig"), 4), lam(symbol("lam"), 4), bet(symbol("bet"), 4);
2728 e = lorentz_eps(mu, nu, rho, sig) *
2729 lorentz_eps(mu.toggle_variance(), nu.toggle_variance(), lam, bet);
2730 cout << simplify_indexed(e) << endl;
2731 // -> 2*eta~bet~rho*eta~sig~lam-2*eta~sig~bet*eta~rho~lam
2733 idx i(symbol("i"), 3), j(symbol("j"), 3), k(symbol("k"), 3);
2734 symbol A("A"), B("B");
2735 e = epsilon_tensor(i, j, k) * indexed(A, j) * indexed(B, k);
2736 cout << simplify_indexed(e) << endl;
2737 // -> -B.k*A.j*eps.i.k.j
2738 e = epsilon_tensor(i, j, k) * indexed(A, j) * indexed(A, k);
2739 cout << simplify_indexed(e) << endl;
2744 @subsection Linear algebra
2746 The @code{matrix} class can be used with indices to do some simple linear
2747 algebra (linear combinations and products of vectors and matrices, traces
2748 and scalar products):
2752 idx i(symbol("i"), 2), j(symbol("j"), 2);
2753 symbol x("x"), y("y");
2755 // A is a 2x2 matrix, X is a 2x1 vector
2756 matrix A(2, 2), X(2, 1);
2761 cout << indexed(A, i, i) << endl;
2764 ex e = indexed(A, i, j) * indexed(X, j);
2765 cout << e.simplify_indexed() << endl;
2766 // -> [[2*y+x],[4*y+3*x]].i
2768 e = indexed(A, i, j) * indexed(X, i) + indexed(X, j) * 2;
2769 cout << e.simplify_indexed() << endl;
2770 // -> [[3*y+3*x,6*y+2*x]].j
2774 You can of course obtain the same results with the @code{matrix::add()},
2775 @code{matrix::mul()} and @code{matrix::trace()} methods (@pxref{Matrices})
2776 but with indices you don't have to worry about transposing matrices.
2778 Matrix indices always start at 0 and their dimension must match the number
2779 of rows/columns of the matrix. Matrices with one row or one column are
2780 vectors and can have one or two indices (it doesn't matter whether it's a
2781 row or a column vector). Other matrices must have two indices.
2783 You should be careful when using indices with variance on matrices. GiNaC
2784 doesn't look at the variance and doesn't know that @samp{F~mu~nu} and
2785 @samp{F.mu.nu} are different matrices. In this case you should use only
2786 one form for @samp{F} and explicitly multiply it with a matrix representation
2787 of the metric tensor.
2790 @node Non-commutative objects, Hash Maps, Indexed objects, Basic Concepts
2791 @c node-name, next, previous, up
2792 @section Non-commutative objects
2794 GiNaC is equipped to handle certain non-commutative algebras. Three classes of
2795 non-commutative objects are built-in which are mostly of use in high energy
2799 @item Clifford (Dirac) algebra (class @code{clifford})
2800 @item su(3) Lie algebra (class @code{color})
2801 @item Matrices (unindexed) (class @code{matrix})
2804 The @code{clifford} and @code{color} classes are subclasses of
2805 @code{indexed} because the elements of these algebras usually carry
2806 indices. The @code{matrix} class is described in more detail in
2809 Unlike most computer algebra systems, GiNaC does not primarily provide an
2810 operator (often denoted @samp{&*}) for representing inert products of
2811 arbitrary objects. Rather, non-commutativity in GiNaC is a property of the
2812 classes of objects involved, and non-commutative products are formed with
2813 the usual @samp{*} operator, as are ordinary products. GiNaC is capable of
2814 figuring out by itself which objects commute and will group the factors
2815 by their class. Consider this example:
2819 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4);
2820 idx a(symbol("a"), 8), b(symbol("b"), 8);
2821 ex e = -dirac_gamma(mu) * (2*color_T(a)) * 8 * color_T(b) * dirac_gamma(nu);
2823 // -> -16*(gamma~mu*gamma~nu)*(T.a*T.b)
2827 As can be seen, GiNaC pulls out the overall commutative factor @samp{-16} and
2828 groups the non-commutative factors (the gammas and the su(3) generators)
2829 together while preserving the order of factors within each class (because
2830 Clifford objects commute with color objects). The resulting expression is a
2831 @emph{commutative} product with two factors that are themselves non-commutative
2832 products (@samp{gamma~mu*gamma~nu} and @samp{T.a*T.b}). For clarification,
2833 parentheses are placed around the non-commutative products in the output.
2835 @cindex @code{ncmul} (class)
2836 Non-commutative products are internally represented by objects of the class
2837 @code{ncmul}, as opposed to commutative products which are handled by the
2838 @code{mul} class. You will normally not have to worry about this distinction,
2841 The advantage of this approach is that you never have to worry about using
2842 (or forgetting to use) a special operator when constructing non-commutative
2843 expressions. Also, non-commutative products in GiNaC are more intelligent
2844 than in other computer algebra systems; they can, for example, automatically
2845 canonicalize themselves according to rules specified in the implementation
2846 of the non-commutative classes. The drawback is that to work with other than
2847 the built-in algebras you have to implement new classes yourself. Symbols
2848 always commute and it's not possible to construct non-commutative products
2849 using symbols to represent the algebra elements or generators. User-defined
2850 functions can, however, be specified as being non-commutative.
2852 @cindex @code{return_type()}
2853 @cindex @code{return_type_tinfo()}
2854 Information about the commutativity of an object or expression can be
2855 obtained with the two member functions
2858 unsigned ex::return_type() const;
2859 unsigned ex::return_type_tinfo() const;
2862 The @code{return_type()} function returns one of three values (defined in
2863 the header file @file{flags.h}), corresponding to three categories of
2864 expressions in GiNaC:
2867 @item @code{return_types::commutative}: Commutes with everything. Most GiNaC
2868 classes are of this kind.
2869 @item @code{return_types::noncommutative}: Non-commutative, belonging to a
2870 certain class of non-commutative objects which can be determined with the
2871 @code{return_type_tinfo()} method. Expressions of this category commute
2872 with everything except @code{noncommutative} expressions of the same
2874 @item @code{return_types::noncommutative_composite}: Non-commutative, composed
2875 of non-commutative objects of different classes. Expressions of this
2876 category don't commute with any other @code{noncommutative} or
2877 @code{noncommutative_composite} expressions.
2880 The value returned by the @code{return_type_tinfo()} method is valid only
2881 when the return type of the expression is @code{noncommutative}. It is a
2882 value that is unique to the class of the object and usually one of the
2883 constants in @file{tinfos.h}, or derived therefrom.
2885 Here are a couple of examples:
2888 @multitable @columnfractions 0.33 0.33 0.34
2889 @item @strong{Expression} @tab @strong{@code{return_type()}} @tab @strong{@code{return_type_tinfo()}}
2890 @item @code{42} @tab @code{commutative} @tab -
2891 @item @code{2*x-y} @tab @code{commutative} @tab -
2892 @item @code{dirac_ONE()} @tab @code{noncommutative} @tab @code{TINFO_clifford}
2893 @item @code{dirac_gamma(mu)*dirac_gamma(nu)} @tab @code{noncommutative} @tab @code{TINFO_clifford}
2894 @item @code{2*color_T(a)} @tab @code{noncommutative} @tab @code{TINFO_color}
2895 @item @code{dirac_ONE()*color_T(a)} @tab @code{noncommutative_composite} @tab -
2899 Note: the @code{return_type_tinfo()} of Clifford objects is only equal to
2900 @code{TINFO_clifford} for objects with a representation label of zero.
2901 Other representation labels yield a different @code{return_type_tinfo()},
2902 but it's the same for any two objects with the same label. This is also true
2905 A last note: With the exception of matrices, positive integer powers of
2906 non-commutative objects are automatically expanded in GiNaC. For example,
2907 @code{pow(a*b, 2)} becomes @samp{a*b*a*b} if @samp{a} and @samp{b} are
2908 non-commutative expressions).
2911 @cindex @code{clifford} (class)
2912 @subsection Clifford algebra
2914 @cindex @code{dirac_gamma()}
2915 Clifford algebra elements (also called Dirac gamma matrices, although GiNaC
2916 doesn't treat them as matrices) are designated as @samp{gamma~mu} and satisfy
2917 @samp{gamma~mu*gamma~nu + gamma~nu*gamma~mu = 2*eta~mu~nu} where @samp{eta~mu~nu}
2918 is the Minkowski metric tensor. Dirac gammas are constructed by the function
2921 ex dirac_gamma(const ex & mu, unsigned char rl = 0);
2924 which takes two arguments: the index and a @dfn{representation label} in the
2925 range 0 to 255 which is used to distinguish elements of different Clifford
2926 algebras (this is also called a @dfn{spin line index}). Gammas with different
2927 labels commute with each other. The dimension of the index can be 4 or (in
2928 the framework of dimensional regularization) any symbolic value. Spinor
2929 indices on Dirac gammas are not supported in GiNaC.
2931 @cindex @code{dirac_ONE()}
2932 The unity element of a Clifford algebra is constructed by
2935 ex dirac_ONE(unsigned char rl = 0);
2938 @strong{Note:} You must always use @code{dirac_ONE()} when referring to
2939 multiples of the unity element, even though it's customary to omit it.
2940 E.g. instead of @code{dirac_gamma(mu)*(dirac_slash(q,4)+m)} you have to
2941 write @code{dirac_gamma(mu)*(dirac_slash(q,4)+m*dirac_ONE())}. Otherwise,
2942 GiNaC will complain and/or produce incorrect results.
2944 @cindex @code{dirac_gamma5()}
2945 There is a special element @samp{gamma5} that commutes with all other
2946 gammas, has a unit square, and in 4 dimensions equals
2947 @samp{gamma~0 gamma~1 gamma~2 gamma~3}, provided by
2950 ex dirac_gamma5(unsigned char rl = 0);
2953 @cindex @code{dirac_gammaL()}
2954 @cindex @code{dirac_gammaR()}
2955 The chiral projectors @samp{(1+/-gamma5)/2} are also available as proper
2956 objects, constructed by
2959 ex dirac_gammaL(unsigned char rl = 0);
2960 ex dirac_gammaR(unsigned char rl = 0);
2963 They observe the relations @samp{gammaL^2 = gammaL}, @samp{gammaR^2 = gammaR},
2964 and @samp{gammaL gammaR = gammaR gammaL = 0}.
2966 @cindex @code{dirac_slash()}
2967 Finally, the function
2970 ex dirac_slash(const ex & e, const ex & dim, unsigned char rl = 0);
2973 creates a term that represents a contraction of @samp{e} with the Dirac
2974 Lorentz vector (it behaves like a term of the form @samp{e.mu gamma~mu}
2975 with a unique index whose dimension is given by the @code{dim} argument).
2976 Such slashed expressions are printed with a trailing backslash, e.g. @samp{e\}.
2978 In products of dirac gammas, superfluous unity elements are automatically
2979 removed, squares are replaced by their values, and @samp{gamma5}, @samp{gammaL}
2980 and @samp{gammaR} are moved to the front.
2982 The @code{simplify_indexed()} function performs contractions in gamma strings,
2988 symbol a("a"), b("b"), D("D");
2989 varidx mu(symbol("mu"), D);
2990 ex e = dirac_gamma(mu) * dirac_slash(a, D)
2991 * dirac_gamma(mu.toggle_variance());
2993 // -> gamma~mu*a\*gamma.mu
2994 e = e.simplify_indexed();
2997 cout << e.subs(D == 4) << endl;
3003 @cindex @code{dirac_trace()}
3004 To calculate the trace of an expression containing strings of Dirac gammas
3005 you use the function
3008 ex dirac_trace(const ex & e, unsigned char rl = 0, const ex & trONE = 4);
3011 This function takes the trace of all gammas with the specified representation
3012 label; gammas with other labels are left standing. The last argument to
3013 @code{dirac_trace()} is the value to be returned for the trace of the unity
3014 element, which defaults to 4. The @code{dirac_trace()} function is a linear
3015 functional that is equal to the usual trace only in @math{D = 4} dimensions.
3016 In particular, the functional is not cyclic in @math{D != 4} dimensions when
3017 acting on expressions containing @samp{gamma5}, so it's not a proper trace.
3018 This @samp{gamma5} scheme is described in greater detail in
3019 @cite{The Role of gamma5 in Dimensional Regularization}.
3021 The value of the trace itself is also usually different in 4 and in
3022 @math{D != 4} dimensions:
3027 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4), rho(symbol("rho"), 4);
3028 ex e = dirac_gamma(mu) * dirac_gamma(nu) *
3029 dirac_gamma(mu.toggle_variance()) * dirac_gamma(rho);
3030 cout << dirac_trace(e).simplify_indexed() << endl;
3037 varidx mu(symbol("mu"), D), nu(symbol("nu"), D), rho(symbol("rho"), D);
3038 ex e = dirac_gamma(mu) * dirac_gamma(nu) *
3039 dirac_gamma(mu.toggle_variance()) * dirac_gamma(rho);
3040 cout << dirac_trace(e).simplify_indexed() << endl;
3041 // -> 8*eta~rho~nu-4*eta~rho~nu*D
3045 Here is an example for using @code{dirac_trace()} to compute a value that
3046 appears in the calculation of the one-loop vacuum polarization amplitude in
3051 symbol q("q"), l("l"), m("m"), ldotq("ldotq"), D("D");
3052 varidx mu(symbol("mu"), D), nu(symbol("nu"), D);
3055 sp.add(l, l, pow(l, 2));
3056 sp.add(l, q, ldotq);
3058 ex e = dirac_gamma(mu) *
3059 (dirac_slash(l, D) + dirac_slash(q, D) + m * dirac_ONE()) *
3060 dirac_gamma(mu.toggle_variance()) *
3061 (dirac_slash(l, D) + m * dirac_ONE());
3062 e = dirac_trace(e).simplify_indexed(sp);
3063 e = e.collect(lst(l, ldotq, m));
3065 // -> (8-4*D)*l^2+(8-4*D)*ldotq+4*D*m^2
3069 The @code{canonicalize_clifford()} function reorders all gamma products that
3070 appear in an expression to a canonical (but not necessarily simple) form.
3071 You can use this to compare two expressions or for further simplifications:
3075 varidx mu(symbol("mu"), 4), nu(symbol("nu"), 4);
3076 ex e = dirac_gamma(mu) * dirac_gamma(nu) + dirac_gamma(nu) * dirac_gamma(mu);
3078 // -> gamma~mu*gamma~nu+gamma~nu*gamma~mu
3080 e = canonicalize_clifford(e);
3082 // -> 2*ONE*eta~mu~nu
3087 @cindex @code{color} (class)
3088 @subsection Color algebra
3090 @cindex @code{color_T()}
3091 For computations in quantum chromodynamics, GiNaC implements the base elements
3092 and structure constants of the su(3) Lie algebra (color algebra). The base
3093 elements @math{T_a} are constructed by the function
3096 ex color_T(const ex & a, unsigned char rl = 0);
3099 which takes two arguments: the index and a @dfn{representation label} in the
3100 range 0 to 255 which is used to distinguish elements of different color
3101 algebras. Objects with different labels commute with each other. The
3102 dimension of the index must be exactly 8 and it should be of class @code{idx},
3105 @cindex @code{color_ONE()}
3106 The unity element of a color algebra is constructed by
3109 ex color_ONE(unsigned char rl = 0);
3112 @strong{Note:} You must always use @code{color_ONE()} when referring to
3113 multiples of the unity element, even though it's customary to omit it.
3114 E.g. instead of @code{color_T(a)*(color_T(b)*indexed(X,b)+1)} you have to
3115 write @code{color_T(a)*(color_T(b)*indexed(X,b)+color_ONE())}. Otherwise,
3116 GiNaC may produce incorrect results.
3118 @cindex @code{color_d()}
3119 @cindex @code{color_f()}
3123 ex color_d(const ex & a, const ex & b, const ex & c);
3124 ex color_f(const ex & a, const ex & b, const ex & c);
3127 create the symmetric and antisymmetric structure constants @math{d_abc} and
3128 @math{f_abc} which satisfy @math{@{T_a, T_b@} = 1/3 delta_ab + d_abc T_c}
3129 and @math{[T_a, T_b] = i f_abc T_c}.
3131 @cindex @code{color_h()}
3132 There's an additional function
3135 ex color_h(const ex & a, const ex & b, const ex & c);
3138 which returns the linear combination @samp{color_d(a, b, c)+I*color_f(a, b, c)}.
3140 The function @code{simplify_indexed()} performs some simplifications on
3141 expressions containing color objects:
3146 idx a(symbol("a"), 8), b(symbol("b"), 8), c(symbol("c"), 8),
3147 k(symbol("k"), 8), l(symbol("l"), 8);
3149 e = color_d(a, b, l) * color_f(a, b, k);
3150 cout << e.simplify_indexed() << endl;
3153 e = color_d(a, b, l) * color_d(a, b, k);
3154 cout << e.simplify_indexed() << endl;
3157 e = color_f(l, a, b) * color_f(a, b, k);
3158 cout << e.simplify_indexed() << endl;
3161 e = color_h(a, b, c) * color_h(a, b, c);
3162 cout << e.simplify_indexed() << endl;
3165 e = color_h(a, b, c) * color_T(b) * color_T(c);
3166 cout << e.simplify_indexed() << endl;
3169 e = color_h(a, b, c) * color_T(a) * color_T(b) * color_T(c);
3170 cout << e.simplify_indexed() << endl;
3173 e = color_T(k) * color_T(a) * color_T(b) * color_T(k);
3174 cout << e.simplify_indexed() << endl;
3175 // -> 1/4*delta.b.a*ONE-1/6*T.a*T.b
3179 @cindex @code{color_trace()}
3180 To calculate the trace of an expression containing color objects you use the
3184 ex color_trace(const ex & e, unsigned char rl = 0);
3187 This function takes the trace of all color @samp{T} objects with the
3188 specified representation label; @samp{T}s with other labels are left
3189 standing. For example:
3193 e = color_trace(4 * color_T(a) * color_T(b) * color_T(c));
3195 // -> -I*f.a.c.b+d.a.c.b
3200 @node Hash Maps, Methods and Functions, Non-commutative objects, Basic Concepts
3201 @c node-name, next, previous, up
3204 @cindex @code{exhashmap} (class)
3206 For your convenience, GiNaC offers the container template @code{exhashmap<T>}
3207 that can be used as a drop-in replacement for the STL
3208 @code{std::map<ex, T, ex_is_less>}, using hash tables to provide faster,
3209 typically constant-time, element look-up than @code{map<>}.
3211 @code{exhashmap<>} supports all @code{map<>} members and operations, with the
3212 following differences:
3216 no @code{lower_bound()} and @code{upper_bound()} methods
3218 no reverse iterators, no @code{rbegin()}/@code{rend()}
3220 no @code{operator<(exhashmap, exhashmap)}
3222 the comparison function object @code{key_compare} is hardcoded to
3225 the constructor @code{exhashmap(size_t n)} allows specifying the minimum
3226 initial hash table size (the actual table size after construction may be
3227 larger than the specified value)
3229 the method @code{size_t bucket_count()} returns the current size of the hash
3232 @code{insert()} and @code{erase()} operations invalidate all iterators
3236 @node Methods and Functions, Information About Expressions, Hash Maps, Top
3237 @c node-name, next, previous, up
3238 @chapter Methods and Functions
3241 In this chapter the most important algorithms provided by GiNaC will be
3242 described. Some of them are implemented as functions on expressions,
3243 others are implemented as methods provided by expression objects. If
3244 they are methods, there exists a wrapper function around it, so you can
3245 alternatively call it in a functional way as shown in the simple
3250 cout << "As method: " << sin(1).evalf() << endl;
3251 cout << "As function: " << evalf(sin(1)) << endl;
3255 @cindex @code{subs()}
3256 The general rule is that wherever methods accept one or more parameters
3257 (@var{arg1}, @var{arg2}, @dots{}) the order of arguments the function
3258 wrapper accepts is the same but preceded by the object to act on
3259 (@var{object}, @var{arg1}, @var{arg2}, @dots{}). This approach is the
3260 most natural one in an OO model but it may lead to confusion for MapleV
3261 users because where they would type @code{A:=x+1; subs(x=2,A);} GiNaC
3262 would require @code{A=x+1; subs(A,x==2);} (after proper declaration of
3263 @code{A} and @code{x}). On the other hand, since MapleV returns 3 on
3264 @code{A:=x^2+3; coeff(A,x,0);} (GiNaC: @code{A=pow(x,2)+3;
3265 coeff(A,x,0);}) it is clear that MapleV is not trying to be consistent
3266 here. Also, users of MuPAD will in most cases feel more comfortable
3267 with GiNaC's convention. All function wrappers are implemented
3268 as simple inline functions which just call the corresponding method and
3269 are only provided for users uncomfortable with OO who are dead set to
3270 avoid method invocations. Generally, nested function wrappers are much
3271 harder to read than a sequence of methods and should therefore be
3272 avoided if possible. On the other hand, not everything in GiNaC is a
3273 method on class @code{ex} and sometimes calling a function cannot be
3277 * Information About Expressions::
3278 * Numerical Evaluation::
3279 * Substituting Expressions::
3280 * Pattern Matching and Advanced Substitutions::
3281 * Applying a Function on Subexpressions::
3282 * Visitors and Tree Traversal::
3283 * Polynomial Arithmetic:: Working with polynomials.
3284 * Rational Expressions:: Working with rational functions.
3285 * Symbolic Differentiation::
3286 * Series Expansion:: Taylor and Laurent expansion.
3288 * Built-in Functions:: List of predefined mathematical functions.
3289 * Multiple polylogarithms::
3290 * Complex Conjugation::
3291 * Built-in Functions:: List of predefined mathematical functions.
3292 * Solving Linear Systems of Equations::
3293 * Input/Output:: Input and output of expressions.
3297 @node Information About Expressions, Numerical Evaluation, Methods and Functions, Methods and Functions
3298 @c node-name, next, previous, up
3299 @section Getting information about expressions
3301 @subsection Checking expression types
3302 @cindex @code{is_a<@dots{}>()}
3303 @cindex @code{is_exactly_a<@dots{}>()}
3304 @cindex @code{ex_to<@dots{}>()}
3305 @cindex Converting @code{ex} to other classes
3306 @cindex @code{info()}
3307 @cindex @code{return_type()}
3308 @cindex @code{return_type_tinfo()}
3310 Sometimes it's useful to check whether a given expression is a plain number,
3311 a sum, a polynomial with integer coefficients, or of some other specific type.
3312 GiNaC provides a couple of functions for this:
3315 bool is_a<T>(const ex & e);
3316 bool is_exactly_a<T>(const ex & e);
3317 bool ex::info(unsigned flag);
3318 unsigned ex::return_type() const;
3319 unsigned ex::return_type_tinfo() const;
3322 When the test made by @code{is_a<T>()} returns true, it is safe to call
3323 one of the functions @code{ex_to<T>()}, where @code{T} is one of the
3324 class names (@xref{The Class Hierarchy}, for a list of all classes). For
3325 example, assuming @code{e} is an @code{ex}:
3330 if (is_a<numeric>(e))
3331 numeric n = ex_to<numeric>(e);
3336 @code{is_a<T>(e)} allows you to check whether the top-level object of
3337 an expression @samp{e} is an instance of the GiNaC class @samp{T}
3338 (@xref{The Class Hierarchy}, for a list of all classes). This is most useful,
3339 e.g., for checking whether an expression is a number, a sum, or a product:
3346 is_a<numeric>(e1); // true
3347 is_a<numeric>(e2); // false
3348 is_a<add>(e1); // false
3349 is_a<add>(e2); // true
3350 is_a<mul>(e1); // false
3351 is_a<mul>(e2); // false
3355 In contrast, @code{is_exactly_a<T>(e)} allows you to check whether the
3356 top-level object of an expression @samp{e} is an instance of the GiNaC
3357 class @samp{T}, not including parent classes.
3359 The @code{info()} method is used for checking certain attributes of
3360 expressions. The possible values for the @code{flag} argument are defined
3361 in @file{ginac/flags.h}, the most important being explained in the following
3365 @multitable @columnfractions .30 .70
3366 @item @strong{Flag} @tab @strong{Returns true if the object is@dots{}}
3367 @item @code{numeric}
3368 @tab @dots{}a number (same as @code{is_a<numeric>(...)})
3370 @tab @dots{}a real integer, rational or float (i.e. is not complex)
3371 @item @code{rational}
3372 @tab @dots{}an exact rational number (integers are rational, too)
3373 @item @code{integer}
3374 @tab @dots{}a (non-complex) integer
3375 @item @code{crational}
3376 @tab @dots{}an exact (complex) rational number (such as @math{2/3+7/2*I})
3377 @item @code{cinteger}
3378 @tab @dots{}a (complex) integer (such as @math{2-3*I})
3379 @item @code{positive}
3380 @tab @dots{}not complex and greater than 0
3381 @item @code{negative}
3382 @tab @dots{}not complex and less than 0
3383 @item @code{nonnegative}
3384 @tab @dots{}not complex and greater than or equal to 0
3386 @tab @dots{}an integer greater than 0
3388 @tab @dots{}an integer less than 0
3389 @item @code{nonnegint}
3390 @tab @dots{}an integer greater than or equal to 0
3392 @tab @dots{}an even integer
3394 @tab @dots{}an odd integer
3396 @tab @dots{}a prime integer (probabilistic primality test)
3397 @item @code{relation}
3398 @tab @dots{}a relation (same as @code{is_a<relational>(...)})
3399 @item @code{relation_equal}
3400 @tab @dots{}a @code{==} relation
3401 @item @code{relation_not_equal}
3402 @tab @dots{}a @code{!=} relation
3403 @item @code{relation_less}
3404 @tab @dots{}a @code{<} relation
3405 @item @code{relation_less_or_equal}
3406 @tab @dots{}a @code{<=} relation
3407 @item @code{relation_greater}
3408 @tab @dots{}a @code{>} relation
3409 @item @code{relation_greater_or_equal}
3410 @tab @dots{}a @code{>=} relation
3412 @tab @dots{}a symbol (same as @code{is_a<symbol>(...)})
3414 @tab @dots{}a list (same as @code{is_a<lst>(...)})
3415 @item @code{polynomial}
3416 @tab @dots{}a polynomial (i.e. only consists of sums and products of numbers and symbols with positive integer powers)
3417 @item @code{integer_polynomial}
3418 @tab @dots{}a polynomial with (non-complex) integer coefficients
3419 @item @code{cinteger_polynomial}
3420 @tab @dots{}a polynomial with (possibly complex) integer coefficients (such as @math{2-3*I})
3421 @item @code{rational_polynomial}
3422 @tab @dots{}a polynomial with (non-complex) rational coefficients
3423 @item @code{crational_polynomial}
3424 @tab @dots{}a polynomial with (possibly complex) rational coefficients (such as @math{2/3+7/2*I})
3425 @item @code{rational_function}
3426 @tab @dots{}a rational function (@math{x+y}, @math{z/(x+y)})
3427 @item @code{algebraic}
3428 @tab @dots{}an algebraic object (@math{sqrt(2)}, @math{sqrt(x)-1})
3432 To determine whether an expression is commutative or non-commutative and if
3433 so, with which other expressions it would commute, you use the methods
3434 @code{return_type()} and @code{return_type_tinfo()}. @xref{Non-commutative objects},
3435 for an explanation of these.
3438 @subsection Accessing subexpressions
3441 Many GiNaC classes, like @code{add}, @code{mul}, @code{lst}, and
3442 @code{function}, act as containers for subexpressions. For example, the
3443 subexpressions of a sum (an @code{add} object) are the individual terms,
3444 and the subexpressions of a @code{function} are the function's arguments.
3446 @cindex @code{nops()}
3448 GiNaC provides several ways of accessing subexpressions. The first way is to
3453 ex ex::op(size_t i);
3456 @code{nops()} determines the number of subexpressions (operands) contained
3457 in the expression, while @code{op(i)} returns the @code{i}-th
3458 (0..@code{nops()-1}) subexpression. In the case of a @code{power} object,
3459 @code{op(0)} will return the basis and @code{op(1)} the exponent. For
3460 @code{indexed} objects, @code{op(0)} is the base expression and @code{op(i)},
3461 @math{i>0} are the indices.
3464 @cindex @code{const_iterator}
3465 The second way to access subexpressions is via the STL-style random-access
3466 iterator class @code{const_iterator} and the methods
3469 const_iterator ex::begin();
3470 const_iterator ex::end();
3473 @code{begin()} returns an iterator referring to the first subexpression;
3474 @code{end()} returns an iterator which is one-past the last subexpression.
3475 If the expression has no subexpressions, then @code{begin() == end()}. These
3476 iterators can also be used in conjunction with non-modifying STL algorithms.
3478 Here is an example that (non-recursively) prints the subexpressions of a
3479 given expression in three different ways:
3486 for (size_t i = 0; i != e.nops(); ++i)
3487 cout << e.op(i) << endl;
3490 for (const_iterator i = e.begin(); i != e.end(); ++i)
3493 // with iterators and STL copy()
3494 std::copy(e.begin(), e.end(), std::ostream_iterator<ex>(cout, "\n"));
3498 @cindex @code{const_preorder_iterator}
3499 @cindex @code{const_postorder_iterator}
3500 @code{op()}/@code{nops()} and @code{const_iterator} only access an
3501 expression's immediate children. GiNaC provides two additional iterator
3502 classes, @code{const_preorder_iterator} and @code{const_postorder_iterator},
3503 that iterate over all objects in an expression tree, in preorder or postorder,
3504 respectively. They are STL-style forward iterators, and are created with the
3508 const_preorder_iterator ex::preorder_begin();
3509 const_preorder_iterator ex::preorder_end();
3510 const_postorder_iterator ex::postorder_begin();
3511 const_postorder_iterator ex::postorder_end();
3514 The following example illustrates the differences between
3515 @code{const_iterator}, @code{const_preorder_iterator}, and
3516 @code{const_postorder_iterator}:
3520 symbol A("A"), B("B"), C("C");
3521 ex e = lst(lst(A, B), C);
3523 std::copy(e.begin(), e.end(),
3524 std::ostream_iterator<ex>(cout, "\n"));
3528 std::copy(e.preorder_begin(), e.preorder_end(),
3529 std::ostream_iterator<ex>(cout, "\n"));
3536 std::copy(e.postorder_begin(), e.postorder_end(),
3537 std::ostream_iterator<ex>(cout, "\n"));
3546 @cindex @code{relational} (class)
3547 Finally, the left-hand side and right-hand side expressions of objects of
3548 class @code{relational} (and only of these) can also be accessed with the
3557 @subsection Comparing expressions
3558 @cindex @code{is_equal()}
3559 @cindex @code{is_zero()}
3561 Expressions can be compared with the usual C++ relational operators like
3562 @code{==}, @code{>}, and @code{<} but if the expressions contain symbols,
3563 the result is usually not determinable and the result will be @code{false},
3564 except in the case of the @code{!=} operator. You should also be aware that
3565 GiNaC will only do the most trivial test for equality (subtracting both
3566 expressions), so something like @code{(pow(x,2)+x)/x==x+1} will return
3569 Actually, if you construct an expression like @code{a == b}, this will be
3570 represented by an object of the @code{relational} class (@pxref{Relations})
3571 which is not evaluated until (explicitly or implicitly) cast to a @code{bool}.
3573 There are also two methods
3576 bool ex::is_equal(const ex & other);
3580 for checking whether one expression is equal to another, or equal to zero,
3584 @subsection Ordering expressions
3585 @cindex @code{ex_is_less} (class)
3586 @cindex @code{ex_is_equal} (class)
3587 @cindex @code{compare()}
3589 Sometimes it is necessary to establish a mathematically well-defined ordering
3590 on a set of arbitrary expressions, for example to use expressions as keys
3591 in a @code{std::map<>} container, or to bring a vector of expressions into
3592 a canonical order (which is done internally by GiNaC for sums and products).
3594 The operators @code{<}, @code{>} etc. described in the last section cannot
3595 be used for this, as they don't implement an ordering relation in the
3596 mathematical sense. In particular, they are not guaranteed to be
3597 antisymmetric: if @samp{a} and @samp{b} are different expressions, and
3598 @code{a < b} yields @code{false}, then @code{b < a} doesn't necessarily
3601 By default, STL classes and algorithms use the @code{<} and @code{==}
3602 operators to compare objects, which are unsuitable for expressions, but GiNaC
3603 provides two functors that can be supplied as proper binary comparison
3604 predicates to the STL:
3607 class ex_is_less : public std::binary_function<ex, ex, bool> @{
3609 bool operator()(const ex &lh, const ex &rh) const;
3612 class ex_is_equal : public std::binary_function<ex, ex, bool> @{
3614 bool operator()(const ex &lh, const ex &rh) const;
3618 For example, to define a @code{map} that maps expressions to strings you
3622 std::map<ex, std::string, ex_is_less> myMap;
3625 Omitting the @code{ex_is_less} template parameter will introduce spurious
3626 bugs because the map operates improperly.
3628 Other examples for the use of the functors:
3636 std::sort(v.begin(), v.end(), ex_is_less());
3638 // count the number of expressions equal to '1'
3639 unsigned num_ones = std::count_if(v.begin(), v.end(),
3640 std::bind2nd(ex_is_equal(), 1));
3643 The implementation of @code{ex_is_less} uses the member function
3646 int ex::compare(const ex & other) const;
3649 which returns @math{0} if @code{*this} and @code{other} are equal, @math{-1}
3650 if @code{*this} sorts before @code{other}, and @math{1} if @code{*this} sorts
3654 @node Numerical Evaluation, Substituting Expressions, Information About Expressions, Methods and Functions
3655 @c node-name, next, previous, up
3656 @section Numerical Evaluation
3657 @cindex @code{evalf()}
3659 GiNaC keeps algebraic expressions, numbers and constants in their exact form.
3660 To evaluate them using floating-point arithmetic you need to call
3663 ex ex::evalf(int level = 0) const;
3666 @cindex @code{Digits}
3667 The accuracy of the evaluation is controlled by the global object @code{Digits}
3668 which can be assigned an integer value. The default value of @code{Digits}
3669 is 17. @xref{Numbers}, for more information and examples.
3671 To evaluate an expression to a @code{double} floating-point number you can
3672 call @code{evalf()} followed by @code{numeric::to_double()}, like this:
3676 // Approximate sin(x/Pi)
3678 ex e = series(sin(x/Pi), x == 0, 6);
3680 // Evaluate numerically at x=0.1
3681 ex f = evalf(e.subs(x == 0.1));
3683 // ex_to<numeric> is an unsafe cast, so check the type first
3684 if (is_a<numeric>(f)) @{
3685 double d = ex_to<numeric>(f).to_double();
3694 @node Substituting Expressions, Pattern Matching and Advanced Substitutions, Numerical Evaluation, Methods and Functions
3695 @c node-name, next, previous, up
3696 @section Substituting expressions
3697 @cindex @code{subs()}
3699 Algebraic objects inside expressions can be replaced with arbitrary
3700 expressions via the @code{.subs()} method:
3703 ex ex::subs(const ex & e, unsigned options = 0);
3704 ex ex::subs(const exmap & m, unsigned options = 0);
3705 ex ex::subs(const lst & syms, const lst & repls, unsigned options = 0);
3708 In the first form, @code{subs()} accepts a relational of the form
3709 @samp{object == expression} or a @code{lst} of such relationals:
3713 symbol x("x"), y("y");
3715 ex e1 = 2*x^2-4*x+3;
3716 cout << "e1(7) = " << e1.subs(x == 7) << endl;
3720 cout << "e2(-2, 4) = " << e2.subs(lst(x == -2, y == 4)) << endl;
3725 If you specify multiple substitutions, they are performed in parallel, so e.g.
3726 @code{subs(lst(x == y, y == x))} exchanges @samp{x} and @samp{y}.
3728 The second form of @code{subs()} takes an @code{exmap} object which is a
3729 pair associative container that maps expressions to expressions (currently
3730 implemented as a @code{std::map}). This is the most efficient one of the
3731 three @code{subs()} forms and should be used when the number of objects to
3732 be substituted is large or unknown.
3734 Using this form, the second example from above would look like this:
3738 symbol x("x"), y("y");