1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238
|
// Copyright (C) 2005 David Saunders, 2010 LinBox.
// This file is part of LinBox, see COPYING for licence information.
/*
* Coypright (c) LinBox
* ========LICENCE========
* This file is part of the library LinBox.
*
* LinBox is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public
* License as published by the Free Software Foundation; either
* version 2.1 of the License, or (at your option) any later version.
*
* This library is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* Lesser General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with this library; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
* ========LICENCE========
*
*/
/*! \page tuto LinBox tutorial
The person requiring some exact linear algebra computation may exploit LinBox
at any of four general levels. In order from least involved with the details
to most involved, they are:
<ol>
<li>
Access using a LinBox web server.
</li>
The user at this level must attend to preparation of a matrix in a suitable
file format and invoking the webservice. The server itself should provide
adequate documentation for this. Servers are generally available at <a
href="http://www.linalg.org" target="_blank">linalg.org</a>.
<li>
Access using compiled codes. Full programs exist in the examples directory and
can be used through an interface to LinBox in a general purposes system such as
Maple or <a href="http://www.sagemath.org" target="_blank">SAGE</a>.
</li>
The user at this level must see to installation, and then attend to preparation
of her matrix in a suitable file format and to the form of the program or
procedure invocation. A number of programs are available in the examples
directory distributed with LinBox providing for rank, determinant, linear
system solution, etc. The documentation of the examples module should serve
for this level of access.
<!-- what about interfaces maple and gap - can they be documented now? -->
<li>
Use of LinBox as a programmers library for exact linear algebra functions.
</li>
At this level a user must do at least the following:
<ol type="a">
<li>
Choose a field or ring representation and construct a specific field or ring
object \c R.
</li>
<li>
Choose a blackbox or matrix representation suitable to your data and construct
specific matrix \c A over \c R.
</li>
<li>
Call the desired algorithm. The solutions directory is designed to support this by
providing functions with simple problem oriented interfaces (\c rank(), \c
det(), \c solve()), yet allowing some user control of algorithm details. The
programmer may be providing some of the parts such as an application specific
blackbox matrix class.
</li>
</ol>
<li>
Power development.
</li>
Again, this is use of LinBox as a library, but with hands fully on the details.
The programmer at this level apparently needs the best opportunities for high
performance and is willing to use the more complicated internal interfaces to
get it. Direct calls to functions in the algorithms directory and perhaps to
related packages such as \c FFLAS, \c FFPACK, or other components are being
used. The programmer working at this level of intimacy with LinBox is likely
to develop some components that ought to be included in future releases of
LinBox. Such contributions are warmly welcomed. The online documentation
system is intended primarily for the programmer at level 3. Thus documentation
is not generated to the last internal detail. It is supposed that the level 4
(and 3) programmer, when studying internal details, will be best served by
reading the code. It is possible, if you wish, to set <a
href="http://www.doxygen.org/" target="_blank">doxygen</a> parameters so as to
have a documentation listing for every class and function. Good luck
separating wheat from chaff in that case.
</ol>
In this tutorial we will discuss a simple application at level 3, the design
of a program to compute the determinant of a sparse matrix over
\f$\mathbf{GF}(101)\f$. The programmer has 3 major issues: field, matrix,
algorithm.
The basic algorithm choice is between blackbox and elimination algorithms.
If our matrix is \c A, then the algorithm call may look like
\code
det(d, A, Method::Blackbox());
//or
det(d, A, Method::Elimination());
\endcode
To have access to this determinant function just <code>\#include <linbox/solutions/det.h></code>.
The larger and sparser the matrix the more likely that the blackbox algorithm
is fastest. Hybrids are under development to help make this choice, and even
to make it at run time adaptively. The last argument to the solutions
functions is a method object which specifies the algorithm to be used. In each
case there is a default, so this argument is optional. Also, some method
objects may be constructed so as to convey more detail about the algorithm to
use. For example it may help to promise that the matrix is nonsingular or
symmetric or to request a particular algorithm variant. The blackbox algorithm
based fundamentally on Wiedemann's approach is specified by
\c Method::Wiedemann, see \c solutions/methods.h for more details. Specification of
a blocked blackbox algorithm may be advantageous (in the future). Elimination
is likely fastest if the matrix is not extremely sparse or not large or not
subject to rapid fill-in.
Of course, first we must construct the matrix and field. For the field, you
must first choose the class (representation type) of field and then instantiate
your specific field.
\code
#include <givaro/modular.h>
typedef Givaro::Modular<short> Field;
Field F(101);
\endcode
It is a good idea to use a \p typedef, making it easy to change representations
later. The \ref modular field representations are LinBox's own and contain
some useful optimizations. Another possibility is use of <a
href="http://shoup.net/ntl/" target="_blank">NTL</a>'s \c NTL::ZZ_p class. The
LinBox wrapper of that, called \p LinBox::NTL_ZZ_p is found in \c
field/ntl-ZZ_p.h. Or use a <a
href="http://ljk.imag.fr/CASYS/LOGICIELS/givaro/" target="_blank">Givaro</a>
table based representation, \c LinBox::Givaro::GFq in \c field/givaro-gfq.h
...and there are many other options. The program \c tests/test-fields.C will
provide some timing data on the basic functions of each representation. In
%LinBox, a \p Field class and the class representing the field entries are
distinct types. The field object owns the arithmetic operations, whereas the
entry objects may be of a primitive type such as \c short, \c int, \c double.
Each field class \p Fld has the embedded class \c Fld::Element for it's
elements.
\code
Field::Element a, b, c;
F.init(a, 2); F.init(b, 3);
F.mul(c, a, b); // c <- a*b
F.write(cout, c);
\endcode
<!-- put here a diatribe on field and element representations -->
You have seen that the field representations are in subdirectory
\c linbox/field. Similarly the matrix representations are in
\c linbox/blackbox. All of these are are suitable for blackbox
algorithms. Only those deriving from representations in \c linbox/matrix
are suitable for elimination. For a sparse matrix, \c LinBox::TriplesBB is a
good representation choice if a blackbox algorithm is to be used. For a
\f$\{0,1\}-\f$incidence matrix, the class \c LinBox::ZeroOne will be more economical of
space and time. On the other hand, if elimination is to be used, those will
not serve at all, but \c LinBox::SparseMatrix, which allows internal
modification, will do the trick. The ordinary matrix representation is
\c LinBox::BlasMatrix. If your matrix structure allows you to compute
matrix-vector products in a way faster than these classes provide, you can gain
by defining your own matrix class adhering to the \ref blackbox interface. The main
components of the blackbox interface are member functions \c apply() and
\c applyTranspose() functions. \c A.apply(y, x) performs \f$y \gets Ax\f$.
Then there is the question of initializing matrix
entries. Two basic approaches are illustrated here.
\code
TriplesBB A;
A.read("triples-file"); // initialize A from data in file (including shape).
\endcode
Or a program can add entries to a sparse or dense matrix one by one using the \c setEntry
function.
\code
SparseMatrix B(100, 200); // 100 by 200 zero matrix (so far).
for ( ... )
{
Field::Element v;
F.init(v, ...);
B.setEntry(i, j, v); // put nonzero v in i,j position.
}
\endcode
Putting it all together, we may have a program looking like this.
\code
#include <fstream>
#include <givaro/modular.h>
#include <linbox/matrix/dense-matrix.h>
#include <linbox/solutions/det.h>
using namespace LinBox;
int main()
{
typedef Givaro::Modular<double> Field;
Field F(101);
SparseMatrix<Field> A(F);
A.read(std::cin);
Field::Element d;
Method::BlasElimination M;
det(d, A, M);
F.write(std::cout << "the determinant is ", d) << std::endl;
return 0;
}
\endcode
A possible format for the input matrix is that of <a
href="http://hpac.imag.fr">hpac.imag.fr</a>.
This tutorial should continue with a discussion of compilation issues, but it
doesn't. Something useful may be learned from examining the Makefile.am in the examples directory.
*/
// vim:syntax=doxygen
|