1%-------------------------------------------------------------------------------
2% The UserGuide.stex file.  Processed into UserGuide.tex via sed.
3%-------------------------------------------------------------------------------
4
5\documentclass[11pt]{article}
6
7\newcommand{\m}[1]{{\bf{#1}}}       % for matrices and vectors
8\newcommand{\tr}{^{\sf T}}          % transpose
9\newcommand{\he}{^{\sf H}}          % complex conjugate transpose
10\newcommand{\implies}{\rightarrow}
11
12\topmargin 0in
13\textheight 8.5in
14\oddsidemargin 0pt
15\evensidemargin 0pt
16\textwidth 6.5in
17
18\begin{document}
19
20\author{Timothy A. Davis \\
21DrTimothyAldenDavis@gmail.com, http://www.suitesparse.com}
22\title{UMFPACK User Guide}
23\date{VERSION 5.7.9, Oct 20, 2019}
24\maketitle
25
26%-------------------------------------------------------------------------------
27\begin{abstract}
28    UMFPACK is a set of routines for solving unsymmetric sparse linear
29    systems, $\m{Ax}=\m{b}$, using the Unsymmetric MultiFrontal method
30    and direct sparse LU factorization.  It is written in ANSI/ISO C, with a
31    MATLAB interface.  UMFPACK relies on the Level-3 Basic
32    Linear Algebra Subprograms (dense matrix multiply) for its performance.
33    This code works on Windows and many versions of Unix (Sun Solaris,
34    Red Hat Linux, IBM AIX, SGI IRIX, and Compaq Alpha).
35\end{abstract}
36%-------------------------------------------------------------------------------
37
38Technical Report TR-04-003 (revised)
39
40Copyright\copyright 1995-2018 by Timothy A. Davis.
41All Rights Reserved.
42UMFPACK is available under alternate licences; contact T. Davis for details.
43
44{\bf UMFPACK License:} See UMFPACK/Doc/License.txt for the license.
45
46{\bf Availability:}
47    http://www.suitesparse.com
48
49{\bf Acknowledgments:}
50
51    This work was supported by the National Science Foundation, under
52    grants DMS-9504974, DMS-9803599, and CCR-0203270.
53    The upgrade to Version 4.1 and the inclusion of the
54    symmetric and 2-by-2 pivoting strategies
55    were done while the author was on sabbatical at
56    Stanford University and Lawrence Berkeley National Laboratory.
57
58%-------------------------------------------------------------------------------
59\newpage
60%-------------------------------------------------------------------------------
61
62\tableofcontents
63
64%-------------------------------------------------------------------------------
65\newpage
66\section{Overview}
67%-------------------------------------------------------------------------------
68
69UMFPACK\footnote{Pronounced with two syllables: umph-pack}
70is a set of routines for solving systems of linear
71equations, $\m{Ax}=\m{b}$, when $\m{A}$ is sparse and unsymmetric.  It is based
72on the Unsymmetric-pattern MultiFrontal method \cite{DavisDuff97,DavisDuff99}.
73UMFPACK factorizes
74$\m{PAQ}$, $\m{PRAQ}$, or $\m{PR}^{-1}\m{AQ}$ into the product $\m{LU}$,
75where $\m{L}$ and $\m{U}$
76are lower and upper triangular, respectively, $\m{P}$ and $\m{Q}$ are
77permutation matrices, and $\m{R}$ is a diagonal matrix of row scaling factors
78(or $\m{R}=\m{I}$ if row-scaling is not used).  Both $\m{P}$ and $\m{Q}$ are
79chosen to reduce fill-in (new nonzeros in $\m{L}$ and $\m{U}$ that are not
80present in $\m{A}$).  The permutation $\m{P}$ has the dual role of reducing
81fill-in and maintaining numerical accuracy (via relaxed partial pivoting
82and row interchanges).
83
84The sparse matrix $\m{A}$ can be square or rectangular, singular
85or non-singular, and real or complex (or any combination).  Only square
86matrices $\m{A}$ can be used to solve $\m{Ax}=\m{b}$ or related systems.
87Rectangular matrices can only be factorized.
88
89UMFPACK first finds a column pre-ordering that reduces fill-in, without regard
90to numerical values.  It scales and analyzes the matrix, and then automatically
91selects one of two strategies for pre-ordering the rows and columns:
92{\em unsymmetric}
93and {\em symmetric}.  These strategies are described below.
94
95First, all pivots with zero Markowitz cost are eliminated and placed in the
96LU factors.  The remaining submatrix $\m{S}$ is then analyzed.
97The following rules are applied, and the first one that matches defines
98the strategy.
99
100\begin{itemize}
101\item Rule 1: $\m{A}$ rectangular $\implies$ unsymmetric.
102\item Rule 2:
103    If the zero-Markowitz elimination results in a rectangular $\m{S}$,
104    or an $\m{S}$ whose diagonal has not been preserved, the
105    unsymmetric strategy is used.
106\item The symmetry $\sigma_1$ of $\m{S}$ is computed.  It is defined as
107    the number of {\em matched} off-diagonal entries, divided by the
108    total number of off-diagonal entries.  An entry $s_{ij}$ is matched
109    if $s_{ji}$ is also an entry.  They need not be numerically equal.
110    An {\em entry} is a value in $\m{A}$ which is present
111    in the input data structure.  All nonzeros are entries, but some entries
112    may be numerically zero.
113    Let $d$ be the number of nonzero entries on the diagonal of $\m{S}$.
114    Let $\m{S}$ be $\nu$-by-$\nu$.
115    Rule 3: $(\sigma_1 \ge 0.5) \:\wedge\: (d \ge 0.9 \nu) \implies$ symmetric.
116    The matrix has a nearly symmetric nonzero pattern (50\% or more),
117    and a mostly-zero-free diagonal (90\% or more nonzero).
118\item Rule 4:
119    Otherwise, the unsymmetric strategy is used.
120\end{itemize}
121
122Each strategy is described below:
123\begin{itemize}
124\item {\em unsymmetric}:
125The column pre-ordering of $\m{S}$ is computed by a modified version of COLAMD
126\cite{DavisGilbertLarimoreNg00_algo,DavisGilbertLarimoreNg00}.
127The method finds a symmetric permutation $\m{Q}$ of the matrix $\m{S}\tr\m{S}$
128(without forming $\m{S}\tr\m{S}$ explicitly).  This is a good choice for
129$\m{Q}$, since the Cholesky factors of $\m{(SQ)\tr(SQ)}$ are an upper bound (in
130terms of nonzero pattern) of the factor $\m{U}$ for the unsymmetric LU
131factorization ($\m{PSQ}=\m{LU}$) regardless of the choice of $\m{P}$
132\cite{GeorgeNg85,GeorgeNg87,GilbertNg93}.  This modified version of
133COLAMD also computes the column elimination tree and post-orders the
134tree.  It finds the upper bound on the number of nonzeros in L and U.
135It also has a different threshold for determining dense rows and columns.
136During factorization, the column pre-ordering can be modified.
137Columns within a single super-column can be reshuffled, to reduce fill-in.
138Threshold partial pivoting is used with no preference given to the diagonal
139entry.  Within a given pivot column $j$, an entry $a_{ij}$ can be chosen if
140$|a_{ij}| \ge 0.1 \max |a_{*j}|$.  Among those numerically acceptable
141entries, the sparsest row $i$ is chosen as the pivot row.
142
143\item {\em symmetric}:
144The column ordering is computed from AMD
145\cite{AmestoyDavisDuff96,AmestoyDavisDuff03},
146applied to the pattern of $\m{S}+\m{S}\tr$
147followed by a post-ordering of the supernodal elimination
148tree of $\m{S}+\m{S}\tr$.
149No modification of the column pre-ordering is made during numerical
150factorization.  Threshold partial pivoting is used, with a strong
151preference given to the diagonal entry.  The diagonal entry is chosen if
152$a_{jj} \ge 0.001 \max |a_{*j}|$.  Otherwise, a sparse row is selected,
153using the same method used by the unsymmetric strategy.
154
155\end{itemize}
156
157The symmetric strategy, and their automatic selection,
158are new to Version 4.1.  Version 4.0 only used the unsymmetric strategy.
159Versions 4.1 through 5.3 included a {\em 2-by-2} ordering strategy,
160but this option has been disabled in Version 5.4.
161
162Once the strategy is selected,
163the factorization of the matrix $\m{A}$ is broken down into the factorization
164of a sequence of dense rectangular frontal matrices.  The frontal matrices are
165related to each other by a supernodal column elimination tree, in which each
166node in the tree represents one frontal matrix.  This analysis phase also
167determines upper bounds on the memory usage, the floating-point operation count,
168and the number of nonzeros in the LU factors.
169
170UMFPACK factorizes each {\em chain} of frontal matrices in a single working
171array, similar to how the unifrontal method \cite{dusc:96} factorizes the whole
172matrix.  A chain of frontal matrices is a sequence of fronts where the parent
173of front $i$ is $i$+1 in the supernodal column elimination tree.  For the
174nonsingular matrices factorized with the unsymmetric strategy, there are
175exactly the same number of chains as there are leaves in the supernodal
176column elimination tree.  UMFPACK is an
177outer-product based, right-looking method.  At the $k$-th step of Gaussian
178elimination, it represents the updated submatrix $\m{A}_k$ as an implicit
179summation of a set of dense sub-matrices (referred to as {\em elements},
180borrowing a phrase from finite-element methods) that arise when the frontal
181matrices are factorized and their pivot rows and columns eliminated.
182
183Each frontal matrix represents the elimination of one or more columns;
184each column of $\m{A}$ will be eliminated in a specific frontal matrix,
185and which frontal matrix will be used for which column is determined by
186the pre-analysis phase.  The pre-analysis phase also determines the worst-case
187size of each frontal matrix so that they can hold any candidate pivot column
188and any candidate pivot row.  From the perspective of the analysis phase, any
189candidate pivot column in the frontal matrix is identical (in terms of nonzero
190pattern), and so is any row.  However, the numeric factorization phase has
191more information than the analysis phase.  It uses this information to reorder
192the columns within each frontal matrix to reduce fill-in.  Similarly, since
193the number of nonzeros in each row and column are maintained (more precisely,
194COLMMD-style approximate degrees \cite{GilbertMolerSchreiber}), a pivot row can
195be selected based on sparsity-preserving criteria (low degree) as well as
196numerical considerations (relaxed threshold partial pivoting).
197
198When the symmetric strategy are used,
199the column preordering is not refined during numeric factorization.
200Row pivoting for sparsity and numerical accuracy is performed if the
201diagonal entry is too small.
202
203More details of the method, including experimental results, are
204described in \cite{Davis03,Davis03_algo}, available at
205http://www.suitesparse.com.
206
207%-------------------------------------------------------------------------------
208\section{Availability}
209%-------------------------------------------------------------------------------
210
211In addition to appearing as a Collected Algorithm of the ACM,
212UMFPACK is available at \newline http://www.suitesparse.com.
213It is included as a built-in routine in MATLAB.
214Version 4.0 (in MATLAB 6.5)
215does not have the symmetric strategy and it takes
216less advantage of the level-3
217BLAS \cite{DaydeDuff99,ACM679a,ATLAS,GotoVandeGeijn02}.
218Versions 5.x through v4.1 tend to be much faster than Version 4.0,
219particularly on unsymmetric matrices with mostly symmetric
220nonzero pattern (such as finite element and circuit simulation matrices).
221Version 3.0 and following make
222use of a modified version of COLAMD V2.0 by Timothy A.~Davis, Stefan
223Larimore, John Gilbert, and Esmond Ng.  The original COLAMD V2.1 is available in
224as a built-in routine in MATLAB V6.0 (or later), and at
225http://www.suitesparse.com.
226These codes are also available in Netlib \cite{netlib} at
227http://www.netlib.org.
228UMFPACK Versions 2.2.1 and earlier, co-authored with Iain Duff,
229are available as
230MA38 (functionally equivalent to Version 2.2.1) in the Harwell
231Subroutine Library.
232
233{\bf NOTE: you must use the correct version of AMD with UMFPACK; using
234an old version of AMD with a newer version of UMFPACK can fail.}
235
236%-------------------------------------------------------------------------------
237\section{Primary changes from prior versions}
238%-------------------------------------------------------------------------------
239
240A detailed list of changes is in the {\tt ChangeLog} file.
241
242%-------------------------------------------------------------------------------
243\subsection{Version 5.7.3}
244%-------------------------------------------------------------------------------
245
246Renamed the MATLAB interface back to {\tt umfpack} from {\tt umpfack2}, since
247it no longer conflicts with any built-in MATLAB function of the same name.
248See the comments for Version 5.0.3 below.
249
250%-------------------------------------------------------------------------------
251\subsection{Version 5.7.0}
252%-------------------------------------------------------------------------------
253
254Replaced memory manager functions and printf with functions
255in {\tt SuiteSparse\_config}.
256
257%-------------------------------------------------------------------------------
258\subsection{Version 5.6.0}
259%-------------------------------------------------------------------------------
260
261Replaced {\tt UFconfig} with {\tt SuiteSparse\_config}, for {\tt
262SuiteSparse\_timer}.  The {\tt UF\_long} type is replaced with {\tt
263SuiteSparse\_long} to avoid potential name conflicts with the {\tt UF\_}
264prefix, but the former is still available for user programs.
265User programs may safely {\tt \#undef} the {\tt UF\_long} macro and
266use {\tt SuiteSparse\_long} instead.  In future versions, {\tt UF\_long} will
267be removed completely from {\tt SuiteSparse\_config.h}.
268
269%-------------------------------------------------------------------------------
270\subsection{Version 5.5.0}
271%-------------------------------------------------------------------------------
272
273Added more user ordering options (interface to CHOLMOD and METIS orderings,
274and the ability to pass in a user ordering function).
275Minor change to how the 64-bit BLAS is used.
276Added an option to disable the search for singletons.
277
278%-------------------------------------------------------------------------------
279\subsection{Version 5.4.0}
280%-------------------------------------------------------------------------------
281
282Disabled the 2-by-2 ordering strategy.
283Bug fix for \verb'umfpack_make.m' for Windows.
284
285%-------------------------------------------------------------------------------
286\subsection{Version 5.3.0}
287%-------------------------------------------------------------------------------
288
289A bug fix for structurally singular matrices, and a compiler workaround for
290gcc versions 4.2.(3 and 4).
291
292%-------------------------------------------------------------------------------
293\subsection{Version 5.2.0}
294%-------------------------------------------------------------------------------
295
296Change of default license.
297Note that your license may be entirely different;
298refer to UMFPACK/Doc/License.txt for your actual license.
299
300%-------------------------------------------------------------------------------
301\subsection{Version 5.1.0}
302%-------------------------------------------------------------------------------
303
304Port of MATLAB interface to 64-bit MATLAB.
305
306%-------------------------------------------------------------------------------
307\subsection{Version 5.0.3}
308%-------------------------------------------------------------------------------
309
310Renamed the MATLAB function to {\tt umfpack2}, so as not to confict with
311itself (the MATLAB built-in version of UMFPACK).
312
313%-------------------------------------------------------------------------------
314\subsection{Version 5.0}
315%-------------------------------------------------------------------------------
316
317Changed {\tt long} to {\tt UF\_long}, controlled by the {\tt UFconfig.h}
318file.  A {\tt UF\_long} is normally just {\tt long}, except on the Windows 64
319(WIN64) platform.  In that case, it becomes {\tt \_\_int64}.
320
321%-------------------------------------------------------------------------------
322\subsection{Version 4.6}
323%-------------------------------------------------------------------------------
324
325Added additional options to {\tt umf\_solve.c}.
326
327%-------------------------------------------------------------------------------
328\subsection{Version 4.5}
329%-------------------------------------------------------------------------------
330
331Added function pointers for malloc, calloc, realloc, free, printf, hypot,
332and complex divisiion, so that these functions can be redefined at run-time.
333Added a version number so you can determine the
334version of UMFPACK at run time or compile time.  UMFPACK requires AMD v2.0
335or later.
336
337%-------------------------------------------------------------------------------
338\subsection{Version 4.4}
339%-------------------------------------------------------------------------------
340
341Bug fix in strategy selection in {\tt umfpack\_*\_qsymbolic}.
342Added packed complex case for all complex input/output arguments.
343Added {\tt umfpack\_get\_determinant}.
344Added minimal support for Microsoft Visual Studio
345(the {\tt umf\_multicompile.c} file).
346
347%-------------------------------------------------------------------------------
348\subsection{Version 4.3.1}
349%-------------------------------------------------------------------------------
350
351Minor bug fix in the forward/backsolve.  This bug had the effect of turning
352off iterative refinement when solving $\m{A}\tr\m{x}=\m{b}$ after factorizing
353$\m{A}$.  UMFPACK mexFunction now factorizes $\m{A}\tr$ in its forward-slash
354operation.
355
356%-------------------------------------------------------------------------------
357\subsection{Version 4.3}
358%-------------------------------------------------------------------------------
359
360No changes are visible to the C or MATLAB user, except the presence of
361one new control parameter in the {\tt Control} array,
362and three new statistics in the {\tt Info} array.
363The primary change is the addition of an (optional) drop tolerance.
364
365%-------------------------------------------------------------------------------
366\subsection{Version 4.1}
367%-------------------------------------------------------------------------------
368
369The following is a summary of the main changes that are visible to the C
370or MATLAB user:
371
372\begin{enumerate}
373
374\item New ordering strategies added.  No changes are required in user code
375    (either C or MATLAB) to use the new default strategy, which is an automatic
376    selection of the unsymmetric, symmetric, or 2-by-2 strategies.
377
378\item Row scaling added.  This is only visible to the MATLAB caller when using
379    the form {\tt [L,U,P,Q,R] = umfpack (A)}, to retrieve the LU factors.
380    Likewise, it is only visible to the C caller when the LU factors are
381    retrieved, or when solving systems with just $\m{L}$ or $\m{U}$.
382    New C-callable and MATLAB-callable routines are included to get and to
383    apply the scale factors computed by UMFPACK.  Row scaling is enabled by
384    default, but can be disabled.  Row scaling usually leads to a better
385    factorization, particularly when the symmetric strategy is used.
386
387\item Error code {\tt UMFPACK\_ERROR\_problem\_to\_large} removed.
388    Version 4.0 would generate this error when the upper bound memory usage
389    exceeded 2GB (for the {\tt int} version), even when the actual memory
390    usage was less than this.  The new version properly handles this case,
391    and can successfully factorize the matrix if sufficient memory is
392    available.
393
394\item New control parameters and statistics provided.
395
396\item The AMD symmetric approximate minimum degree ordering routine added
397    \cite{AmestoyDavisDuff96,AmestoyDavisDuff03}.
398    It is used by UMFPACK, and can also be called independently from C or
399    MATLAB.
400
401\item The {\tt umfpack} mexFunction now returns permutation matrices, not
402    permutation vectors, when using the form {\tt [L,U,P,Q] = umfpack (A)}
403    or the new form {\tt [L,U,P,Q,R] = umfpack (A)}.
404
405\item New arguments added to the user-callable routines
406    {\tt umfpack\_*\_symbolic},
407    {\tt umfpack\_*\_qsymbolic},
408    {\tt umfpack\_*\_get\_numeric}, and
409    {\tt umfpack\_*\_get\_symbolic}.
410    The symbolic analysis now makes use of the numerical values of the matrix
411    $\m{A}$, to guide the 2-by-2 strategy.  The subsequent matrix passed to
412    the numeric factorization step does not have to have the same numerical
413    values.  All of the new arguments are optional.  If you do not wish to
414    include them, simply pass {\tt NULL} pointers instead.  The 2-by-2 strategy
415    will assume all entries are numerically large, for example.
416
417\item New routines added to save and load the {\tt Numeric} and {\tt Symbolic}
418    objects to and from a binary file.
419
420\item A Fortran interface added.  It provides access to a subset of
421    UMFPACK's features.
422
423\item You can compute an incomplete LU factorization, by dropping small
424    entries from $\m{L}$ and $\m{U}$.  By default, no nonzero entry is
425    dropped, no matter how small in absolute value.  This feature is new
426    to Version 4.3.
427
428\end{enumerate}
429
430%-------------------------------------------------------------------------------
431\section{Using UMFPACK in MATLAB}
432%-------------------------------------------------------------------------------
433
434The easiest way to use UMFPACK is within MATLAB.  Version 4.3 is a built-in
435routine in MATLAB 7.0.4, and is used in {\tt x = A}$\backslash${\tt b} when
436{\tt A} is sparse, square, unsymmetric (or symmetric but not positive definite),
437and with nonzero entries that are not confined in a narrow band.
438It is also used for the {\tt [L,U,P,Q] = lu (A)} usage of {\tt lu}.
439Type {\tt help lu} in MATLAB 6.5 or later for more details.
440
441To compile both the UMFPACK and AMD mexFunctions, just type {\tt umfpack\_make}
442in MATLAB, while in the {\tt UMFPACK/MATLAB} directory.
443
444See Section~\ref{Install} for more details on how to install UMFPACK.
445Once installed, the UMFPACK mexFunction can analyze, factor, and solve linear
446systems.  Table~\ref{matlab} summarizes some of the more common uses
447of the UMFPACK mexFunction within MATLAB.
448
449An optional input argument can be used to modify the control parameters for
450UMFPACK, and an optional output argument provides statistics on the
451factorization.
452
453Refer to the AMD User Guide for more details about the AMD mexFunction.
454
455\begin{table}
456\caption{Using UMFPACK's MATLAB interface}
457\label{matlab}
458\vspace{0.1in}
459{\footnotesize
460\begin{tabular}{l|l|l}
461\hline
462Function & Using UMFPACK & MATLAB 6.0 equivalent \\
463\hline
464 & & \\
465\begin{minipage}[t]{1.5in}
466Solve $\m{Ax}=\m{b}$.
467\end{minipage}
468&
469\begin{minipage}[t]{2.2in}
470\begin{verbatim}
471x = umfpack (A,'\',b) ;
472\end{verbatim}
473\end{minipage}
474&
475\begin{minipage}[t]{2.2in}
476\begin{verbatim}
477x = A \ b ;
478\end{verbatim}
479\end{minipage}
480 \\
481 & & \\
482\hline
483 & & \\
484\begin{minipage}[t]{1.5in}
485Solve $\m{Ax}=\m{b}$ using a different row and column pre-ordering
486(symmetric strategy).
487\end{minipage}
488&
489\begin{minipage}[t]{2.2in}
490\begin{verbatim}
491S = spones (A) ;
492Q = symamd (S+S') ;
493Control = umfpack ;
494Control.strategy = 'symmetric' ;
495x = umfpack (A,Q,'\',b,Control) ;
496\end{verbatim}
497\end{minipage}
498&
499\begin{minipage}[t]{2.2in}
500\begin{verbatim}
501spparms ('autommd',0) ;
502S = spones (A) ;
503Q = symamd (S+S') ;
504x = A (Q,Q) \ b (Q) ;
505x (Q) = x ;
506spparms ('autommd',1) ;
507\end{verbatim}
508\end{minipage}
509 \\
510 & & \\
511\hline
512 & & \\
513\begin{minipage}[t]{1.5in}
514Solve $\m{A}\tr\m{x}\tr = \m{b}\tr$.
515\end{minipage}
516&
517\begin{minipage}[t]{2.2in}
518\begin{verbatim}
519x = umfpack (b,'/',A) ;
520\end{verbatim}
521Note: $\m{A}$ is factorized.
522\end{minipage}
523&
524\begin{minipage}[t]{2.2in}
525\begin{verbatim}
526x = b / A ;
527\end{verbatim}
528Note: $\m{A}\tr$ is factorized.
529\end{minipage}
530 \\
531 & & \\
532\hline
533 & & \\
534\begin{minipage}[t]{1.5in}
535Scale and factorize $\m{A}$, then solve $\m{Ax}=\m{b}$.
536\end{minipage}
537&
538\begin{minipage}[t]{2.2in}
539\begin{verbatim}
540[L,U,P,Q,R] = umfpack (A) ;
541c = P * (R \ b) ;
542x = Q * (U \ (L \ c)) ;
543\end{verbatim}
544\end{minipage}
545&
546\begin{minipage}[t]{2.2in}
547\begin{verbatim}
548[m n] = size (A) ;
549r = full (sum (abs (A), 2)) ;
550r (find (r == 0)) = 1 ;
551R = spdiags (r, 0, m, m) ;
552I = speye (n) ;
553Q = I (:, colamd (A)) ;
554[L,U,P] = lu ((R\A)*Q) ;
555c = P * (R \ b) ;
556x = Q * (U \ (L \ c)) ;
557\end{verbatim}
558\end{minipage}
559 \\
560 & & \\
561\hline
562\end{tabular}
563}
564\end{table}
565
566Note: in MATLAB 6.5 or later, use {\tt spparms ('autoamd',0)} in addition to
567{\tt spparms ('autommd',0)}, in Table~\ref{matlab}, to turn off MATLAB's
568default reordering.
569
570UMFPACK requires
571{\tt b} to be a dense vector (real or complex) of the appropriate dimension.
572This is more restrictive than what you can do with MATLAB's
573backslash or forward slash.  See {\tt umfpack\_solve} for an M-file that
574removes this restriction.
575This restriction does not apply to the built-in backslash operator
576in MATLAB 6.5 or later, which uses UMFPACK to factorize the matrix.
577You can do this yourself in MATLAB:
578
579{\footnotesize
580\begin{verbatim}
581    [L,U,P,Q,R] = umfpack (A) ;
582    x = Q * (U \ (L \ (P * (R \ b)))) ;
583\end{verbatim}
584}
585
586or, with no row scaling:
587
588{\footnotesize
589\begin{verbatim}
590    [L,U,P,Q] = umfpack (A) ;
591    x = Q * (U \ (L \ (P * b))) ;
592\end{verbatim}
593}
594
595The above examples do not make use of the iterative refinement
596that is built into
597{\tt x = }{\tt umfpack (A,'}$\backslash${\tt ',b)}
598however.
599
600MATLAB's {\tt [L,U,P] = lu(A)} returns a lower triangular {\tt L}, an upper
601triangular {\tt U}, and a permutation matrix {\tt P} such that {\tt P*A} is
602equal to {\tt L*U}.  UMFPACK behaves differently.  By default, it scales
603the rows of {\tt A} and reorders the columns of {\tt A} prior to
604factorization, so that {\tt L*U} is equal to {\tt P*(R}$\backslash${\tt A)*Q},
605where {\tt R} is a diagonal sparse matrix of scale factors for the rows
606of {\tt A}.  The scale factors {\tt R} are applied to {\tt A} via the MATLAB
607expression {\tt R}$\backslash${\tt A} to avoid multiplying by
608the reciprocal, which can be numerically inaccurate.
609
610There are more options; you can provide your own column pre-ordering (in which
611case UMFPACK does not call COLAMD or AMD), you can modify other control settings
612(similar to the {\tt spparms} in MATLAB), and you can get various statistics on
613the analysis, factorization, and solution of the linear system.  Type
614{\tt umfpack\_details} and {\tt umfpack\_report} in MATLAB for more
615information.  Two demo M-files are provided.   Just type {\tt umfpack\_simple}
616and {\tt umfpack\_demo} to run them.
617The output of these two programs should be about the same
618as the files {\tt umfpack\_simple.m.out} and {\tt umfpack\_demo.m.out}
619that are provided.
620
621Factorizing {\tt A'} (or {\tt A.'}) and using the transposed factors can
622sometimes be faster than factorizing {\tt A}.  It can also be preferable to
623factorize {\tt A'} if {\tt A} is rectangular.  UMFPACK pre-orders the columns
624to maintain sparsity; the row ordering is not determined until the matrix
625is factorized.  Thus, if {\tt A} is {\tt m} by {\tt n} with structural
626rank {\tt m}
627and {\tt m} $<$ {\tt n}, then {\tt umfpack} might not find a factor
628{\tt U} with a structurally zero-free diagonal.
629Unless the matrix ill-conditioned or
630poorly scaled, factorizing {\tt A'} in this case will guarantee that both
631factors will have zero-free diagonals (in the structural sense; they may
632be numerically zero).  Note that there is no guarantee
633as to the size of the diagonal entries of {\tt U}; UMFPACK does not do a
634rank-revealing factorization.  Here's how you can factorize {\tt A'}
635and get the factors of {\tt A} instead:
636
637\begin{verbatim}
638    [l,u,p,q] = umfpack (A') ;
639    L = u' ;
640    U = l' ;
641    P = q ;
642    Q = p ;
643    clear l u p q
644\end{verbatim}
645
646This is an alternative to {\tt [L,U,P,Q]=umfpack(A)}.
647
648A simple M-file ({\tt umfpack\_btf}) is provided that first permutes the matrix
649to upper block triangular form, using MATLAB's {\tt dmperm} routine, and then
650solves each block.  The LU factors are not returned.  Its usage is simple:
651{\tt x = umfpack\_btf(A,b)}.  Type {\tt help umfpack\_btf} for more options.
652An estimate of the 1-norm of {\tt L*U-P*A*Q} can be computed in MATLAB
653as {\tt lu\_normest(P*A*Q,L,U)}, using the {\tt lu\_normest.m} M-file
654by Hager and Davis \cite{DavisHager99} that is included with the
655UMFPACK distribution.  With row scaling enabled, use
656{\tt lu\_normest(P*(R}$\backslash${\tt A)*Q,L,U)} instead.
657
658One issue you may encounter is how UMFPACK allocates its memory when being used
659in a mexFunction.  One part of its working space is of variable size.   The
660symbolic analysis phase determines an upper bound on the size of this memory,
661but not all of this memory will typically be used in the numerical
662factorization.  UMFPACK tries to allocate a decent amount of working space.
663This is 70\% of the upper bound, by default, for the unsymmetric strategy.
664For the symmetric strategy, the fraction of the upper bound is computed
665automatically (assuming a best-case scenario with no numerical pivoting
666required during numeric factorization).
667If this initial allocation fails, it reduces its request
668and uses less memory.   If the space is not large enough during factorization,
669it is increased via {\tt mxRealloc}.
670
671However, {\tt mxMalloc} and {\tt mxRealloc} abort the {\tt umfpack} mexFunction
672if they fail, so this strategy does not work in MATLAB.
673
674To compute the determinant with UMFPACK:
675
676\begin{verbatim}
677    d = umfpack (A, 'det') ;
678    [d e] = umfpack (A, 'det') ;
679\end{verbatim}
680
681The first case is identical to MATLAB's {\tt det}.
682The second case returns the determinant in the form
683$d \times 10^e$, which avoids overflow if $e$ is large.
684
685%-------------------------------------------------------------------------------
686\section{Using UMFPACK in a C program}
687\label{C}
688%-------------------------------------------------------------------------------
689
690The C-callable UMFPACK library consists of 32 user-callable routines and one
691include file.  All but three of the routines come in four versions, with
692different sizes of integers and for real or complex floating-point numbers:
693\begin{enumerate}
694\item {\tt umfpack\_di\_*}: real double precision, {\tt int} integers.
695\item {\tt umfpack\_dl\_*}: real double precision, {\tt SuiteSparse\_long} integers.
696\item {\tt umfpack\_zi\_*}: complex double precision, {\tt int} integers.
697\item {\tt umfpack\_zl\_*}: complex double precision, {\tt SuiteSparse\_long} integers.
698\end{enumerate}
699where {\tt *} denotes the specific name of one of the routines.
700Routine names beginning with {\tt umf\_} are internal to the package,
701and should not be called by the user.  The include file {\tt umfpack.h}
702must be included in any C program that uses UMFPACK.
703The other three routines are the same for all four versions.
704
705In addition, the C-callable AMD library distributed with UMFPACK
706includes 4 user-callable routines (in two versions with {\tt int} and
707{\tt SuiteSparse\_long} integers) and one include file.  Refer to the AMD documentation
708for more details.
709
710Use only one version for any one problem; do not attempt to use one version
711to analyze the matrix and another version to factorize the matrix, for example.
712
713The notation {\tt umfpack\_di\_*} refers to all user-callable routines
714for the real double precision and {\tt int} integer case.  The notation
715{\tt umfpack\_*\_numeric}, for example, refers all four versions
716(real/complex, int/SuiteSparse\_long) of a single operation
717(in this case numeric factorization).
718
719%-------------------------------------------------------------------------------
720\subsection{The size of an integer}
721%-------------------------------------------------------------------------------
722
723The {\tt umfpack\_di\_*} and {\tt umfpack\_zi\_*} routines use {\tt int} integer
724arguments; those starting with {\tt umfpack\_dl\_} or {\tt umfpack\_zl\_}
725use {\tt SuiteSparse\_long} integer arguments.  If you compile UMFPACK in the standard
726ILP32 mode (32-bit {\tt int}'s, {\tt long}'s, and pointers) then the versions
727are essentially identical.  You will be able to solve problems using up to 2GB
728of memory.  If you compile UMFPACK in the standard LP64 mode, the size of an
729{\tt int} remains 32-bits, but the size of a {\tt long} and a pointer both get
730promoted to 64-bits.  In the LP64 mode, the {\tt umfpack\_dl\_*}
731and {\tt umfpack\_zl\_*} routines can solve huge
732problems (not limited to 2GB), limited of course by the amount of available
733memory.  The only drawback to the 64-bit mode is that not all BLAS libraries
734support 64-bit integers.  This limits the performance you will obtain.
735Those that do support 64-bit integers are specific to particular
736architectures, and are not portable.  UMFPACK and AMD should be compiled
737in the same mode.
738If you compile UMFPACK and AMD in the LP64 mode,
739be sure to add {\tt -DLP64} to the compilation command.  See the examples in
740the {\tt SuiteSparse\_config/SuiteSparse\_config.mk} file.
741
742%-------------------------------------------------------------------------------
743\subsection{Real and complex floating-point}
744%-------------------------------------------------------------------------------
745
746The {\tt umfpack\_di\_*} and {\tt umfpack\_dl\_*} routines take (real) double
747precision arguments, and return double precision arguments.  In the
748{\tt umfpack\_zi\_*} and {\tt umfpack\_zl\_*} routines, these same arguments
749hold the real part of the matrices; and second double precision arrays hold
750the imaginary part of the input and output matrices.  Internally, complex
751numbers are stored in arrays with their real and imaginary parts interleaved,
752as required by the BLAS (``packed'' complex form).
753
754New to Version 4.4 is the option of providing input/output arguments
755in packed complex form.
756
757%-------------------------------------------------------------------------------
758\subsection{Primary routines, and a simple example}
759%-------------------------------------------------------------------------------
760
761Five primary UMFPACK routines are required to factorize $\m{A}$ or
762solve $\m{Ax}=\m{b}$.  They are fully described in Section~\ref{Primary}:
763
764\begin{itemize}
765\item {\tt umfpack\_*\_symbolic}:
766
767    Pre-orders the columns of $\m{A}$ to reduce fill-in.
768    Returns an opaque {\tt Symbolic} object as a {\tt void *}
769    pointer.  The object contains the symbolic analysis and is needed for the
770    numeric factorization.  This routine requires only $O(|\m{A}|)$ space,
771    where $|\m{A}|$ is the number of nonzero entries in the matrix.  It computes
772    upper bounds on the nonzeros in $\m{L}$ and $\m{U}$, the floating-point
773    operations required, and the memory usage of {\tt umfpack\_*\_numeric}.  The
774    {\tt Symbolic} object is small; it contains just the column pre-ordering,
775    the supernodal column elimination tree, and information about each frontal
776    matrix. It is no larger than about $13n$ integers if $\m{A}$ is
777    $n$-by-$n$.
778
779\item {\tt umfpack\_*\_numeric}:
780
781    Numerically scales and then factorizes a sparse matrix into
782    $\m{PAQ}$, $\m{PRAQ}$, or $\m{PR}^{-1}\m{AQ}$ into the product $\m{LU}$,
783    where
784    $\m{P}$ and $\m{Q}$ are permutation matrices, $\m{R}$ is a diagonal
785    matrix of scale factors, $\m{L}$ is lower triangular with unit diagonal,
786    and $\m{U}$ is upper triangular.  Requires the
787    symbolic ordering and analysis computed by {\tt umfpack\_*\_symbolic}
788    or {\tt umfpack\_*\_qsymbolic}.
789    Returns an opaque {\tt Numeric} object as a
790    {\tt void *} pointer.  The object contains the numerical factorization and
791    is used by {\tt umfpack\_*\_solve}.  You can factorize a new matrix with a
792    different values (but identical pattern) as the matrix analyzed by
793    {\tt umfpack\_*\_symbolic} or {\tt umfpack\_*\_qsymbolic} by re-using the
794    {\tt Symbolic} object (this feature is available when using UMFPACK in a
795    C or Fortran program, but not in MATLAB).
796    The matrix
797    $\m{U}$ will have zeros on the diagonal if $\m{A}$ is singular; this
798    produces a warning, but the factorization is still valid.
799
800\item {\tt umfpack\_*\_solve}:
801
802    Solves a sparse linear system ($\m{Ax}=\m{b}$, $\m{A}\tr\m{x}=\m{b}$, or
803    systems involving just $\m{L}$ or $\m{U}$), using the numeric factorization
804    computed by {\tt umfpack\_*\_numeric}.  Iterative refinement with sparse
805    backward error \cite{ardd:89} is used by default.  The matrix $\m{A}$ must
806    be square.  If it is singular, then a divide-by-zero will occur, and your
807    solution with contain IEEE Inf's or NaN's in the appropriate places.
808
809\item {\tt umfpack\_*\_free\_symbolic}:
810
811    Frees the {\tt Symbolic} object created by {\tt umfpack\_*\_symbolic}
812    or {\tt umfpack\_*\_qsymbolic}.
813
814\item {\tt umfpack\_*\_free\_numeric}:
815
816    Frees the {\tt Numeric} object created by {\tt umfpack\_*\_numeric}.
817
818\end{itemize}
819
820Be careful not to free a {\tt Symbolic} object with
821{\tt umfpack\_*\_free\_numeric}.  Nor should you attempt to free a {\tt Numeric}
822object with {\tt umfpack\_*\_free\_symbolic}.
823Failure to free these objects will lead to memory leaks.
824
825The matrix $\m{A}$ is represented in compressed column form, which is
826identical to the sparse matrix representation used by MATLAB.  It consists
827of three or four arrays, where the matrix is {\tt m}-by-{\tt n},
828with {\tt nz} entries.  For the {\tt int} version of UMFPACK:
829
830{\footnotesize
831\begin{verbatim}
832     int Ap [n+1] ;
833     int Ai [nz] ;
834     double Ax [nz] ;
835\end{verbatim}
836}
837
838For the {\tt SuiteSparse\_long} version of UMFPACK:
839
840{\footnotesize
841\begin{verbatim}
842     SuiteSparse_long Ap [n+1] ;
843     SuiteSparse_long Ai [nz] ;
844     double Ax [nz] ;
845\end{verbatim}
846}
847
848The complex versions add another array for the imaginary part:
849
850{\footnotesize
851\begin{verbatim}
852     double Az [nz] ;
853\end{verbatim}
854}
855
856Alternatively, if {\tt Az} is {\tt NULL},
857the real part of the $k$th entry is located in
858{\tt Ax[2*k]} and the imaginary part is located in
859{\tt Ax[2*k+1]}, and the {\tt Ax} array is of size {\tt 2*nz}.
860
861All nonzeros are entries, but an entry may be numerically zero.  The row indices
862of entries in column {\tt j} are stored in
863    {\tt Ai[Ap[j]} \ldots {\tt Ap[j+1]-1]}.
864The corresponding numerical values are stored in
865    {\tt Ax[Ap[j]} \ldots {\tt Ap[j+1]-1]}.
866The imaginary part, for the complex versions, is stored in
867    {\tt Az[Ap[j]} \ldots {\tt Ap[j+1]-1]}
868    (see above for the packed complex case).
869
870No duplicate row indices may be present, and the row indices in any given
871column must be sorted in ascending order.  The first entry {\tt Ap[0]} must be
872zero.  The total number of entries in the matrix is thus {\tt nz = Ap[n]}.
873Except for the fact that extra zero entries can be included, there is thus a
874unique compressed column representation of any given matrix $\m{A}$.
875For a more flexible method for providing an input matrix to UMFPACK,
876see Section~\ref{triplet}.
877
878Here is a simple main program, {\tt umfpack\_simple.c}, that illustrates the
879basic usage of UMFPACK.  See Section~\ref{Synopsis} for a short description
880of each calling sequence, including a list of options for the first
881argument of {\tt umfpack\_di\_solve}.
882
883{\footnotesize
884\begin{verbatim}
885INCLUDE umfpack_simple.c via sed
886\end{verbatim}
887}
888
889The {\tt Ap}, {\tt Ai}, and {\tt Ax} arrays represent the matrix
890\[
891\m{A} = \left[
892\begin{array}{rrrrr}
893 2 &  3 &  0 &  0 &  0 \\
894 3 &  0 &  4 &  0 &  6 \\
895 0 & -1 & -3 &  2 &  0 \\
896 0 &  0 &  1 &  0 &  0 \\
897 0 &  4 &  2 &  0 &  1 \\
898\end{array}
899\right].
900\]
901and the solution to $\m{Ax}=\m{b}$ is $\m{x} = [1 \, 2 \, 3 \, 4 \, 5]\tr$.
902The program uses default control settings and does not return any statistics
903about the ordering, factorization, or solution ({\tt Control} and {\tt Info}
904are both {\tt (double *) NULL}).  It also ignores the status value returned by
905most user-callable UMFPACK routines.
906
907%-------------------------------------------------------------------------------
908\subsection{A note about zero-sized arrays}
909%-------------------------------------------------------------------------------
910
911UMFPACK uses many user-provided arrays of
912size {\tt m} or {\tt n} (the order of the matrix), and of size
913{\tt nz} (the number of nonzeros in a matrix).  UMFPACK does not handle
914zero-dimensioned arrays;
915it returns an error code if {\tt m} or {\tt n}
916are zero.  However, {\tt nz} can be zero, since all singular matrices are
917handled correctly.  If you attempt to {\tt malloc} an array of size {\tt nz}
918= 0, however, {\tt malloc} will return a null pointer which UMFPACK will report
919as a missing argument.  If you {\tt malloc} an array of
920size {\tt nz} to pass to UMFPACK, make sure that you handle the {\tt nz} = 0
921case correctly (use a size equal to the maximum of {\tt nz} and 1, or use a
922size of {\tt nz+1}).
923
924%-------------------------------------------------------------------------------
925\subsection{Alternative routines}
926%-------------------------------------------------------------------------------
927
928Three alternative routines are provided that modify UMFPACK's default
929behavior.  They are fully described in Section~\ref{Alternative}:
930
931\begin{itemize}
932\item {\tt umfpack\_*\_defaults}:
933
934    Sets the default control parameters in the {\tt Control} array.  These can
935    then be modified as desired before passing the array to the other UMFPACK
936    routines.  Control parameters are summarized in Section~\ref{control_param}.
937    Three particular parameters deserve special notice.
938    UMFPACK uses relaxed partial pivoting, where a candidate pivot entry is
939    numerically acceptable if its magnitude is greater than or equal to a
940    tolerance parameter times the magnitude of the largest entry in the same
941    column.  The parameter {\tt Control [UMFPACK\_PIVOT\_TOLERANCE]} has a
942    default value of 0.1, and is used for the unsymmetric strategy.
943    For complex matrices, a cheap approximation of the absolute value is
944    used for the threshold pivoting test
945    ($|a| \approx |a_{\mbox{real}}|+|a_{\mbox{imag}}|$).
946
947    For the symmetric strategy, a second tolerance is used for diagonal
948    entries: \newline {\tt Control [UMFPACK\_SYM\_PIVOT\_TOLERANCE]}, with
949    a default value of 0.001.  The first parameter (with a default of 0.1)
950    is used for any off-diagonal candidate pivot entries.
951
952    These two parameters may be too small for some matrices, particularly for
953    ill-conditioned or poorly scaled ones.  With the default pivot tolerances
954    and default iterative refinement,
955        {\tt x = umfpack (A,'}$\backslash${\tt ',b)}
956    is just as accurate as (or more accurate) than
957        {\tt x = A}$\backslash${\tt b}
958    in MATLAB 6.1 for nearly all matrices.
959
960    If {\tt Control [UMFPACK\_PIVOT\_TOLERANCE]} is zero, than any
961    nonzero entry is acceptable as a pivot (this is changed from Version 4.0,
962    which treated a value of 0.0 the same as 1.0).  If the symmetric strategy is
963    used, and {\tt Control [UMFPACK\_SYM\_PIVOT\_TOLERANCE]} is zero, then any
964    nonzero entry on the diagonal is accepted as a pivot.  Off-diagonal pivoting
965    will still occur if the diagonal entry is exactly zero.  The
966    {\tt Control [UMFPACK\_SYM\_PIVOT\_TOLERANCE]} parameter is new to Version
967    4.1.  It is similar in function to the pivot tolerance for left-looking
968    methods (the MATLAB {\tt THRESH} option in {\tt [L,U,P] = lu (A, THRESH)},
969    and the pivot tolerance parameter in SuperLU).
970
971    The parameter {\tt Control [UMFPACK\_STRATEGY]} can be used to bypass
972    UMFPACK's automatic strategy selection.  The automatic strategy nearly
973    always selects the best method.  When it does not, the different methods
974    nearly always give about the same quality of results.  There may be
975    cases where the automatic strategy fails to pick a good strategy. Also,
976    you can save some computing time if you know the right strategy for your
977    set of matrix problems.  The default is \verb'UMFPACK_STRATEGY_AUTO',
978    in which UMFPACK selects the strategy by itself.
979    \verb'UMFPACK_STRATEGY_UNSYMMETRIC' gives the unsymmetric
980    strategy, which is to use a column pre-ordering (such as COLAMD)
981    and to give no preference to the diagonal during partial pivoting.
982    \verb'UMFPACK_STRATEGY_SYMMETRIC' gives the symmetric strategy,
983    which is to use a symmetric row and column ordering (such as AMD)
984    and to give strong preference to the diagonal during partial pivoting.
985
986    The parameter {\tt Control [UMFPACK\_ORDERING]} defines what ordering
987    method UMFPACK should use.  The options are:
988
989    \begin{itemize}
990    \item \verb'UMFPACK_ORDERING_CHOLMOD' (0).
991        This is the method used by CHOLMOD.  It first tries AMD or COLAMD
992        (depending on what strategy is used).  If that method gives low
993        fill-in, it is used without trying METIS at all.
994        Otherwise METIS is tried (on $\m{A}^T\m{A}$ for the unsymmetric
995        strategy, or $\m{A}+\m{A}^T$ for the symmetric strategy),
996        and the ordering (AMD/COLAMD or METIS) giving the lowest fill-in is
997        used.
998
999    \item \verb'UMFPACK_ORDERING_DEFAULT'.  This is the same as
1000    \verb'UMFPACK_ORDERING_AMD'.
1001
1002    \item \verb'UMFPACK_ORDERING_AMD' (1).  This is the default.
1003        AMD is used for the symmetric strategy, on the pattern of
1004        $\m{A}+\m{A}^T$.  COLAMD is used on $\m{A}$ for the unsymmetric
1005        strategy.
1006
1007    \item \verb'UMFPACK_ORDERING_GIVEN' (2).  This is assumed if
1008        a permutation is provided to \newline
1009        {\tt umfpack\_*\_qsymbolic}.
1010
1011    \item \verb'UMFPACK_ORDERING_METIS' (3).
1012        Use METIS (on $\m{A}^T\m{A}$ for the unsymmetric
1013        strategy, or $\m{A}+\m{A}^T$ for the symmetric strategy).
1014
1015    \item \verb'UMFPACK_ORDERING_BEST' (4).
1016        Try three methods and take the best.  The three methods are
1017        AMD/COLAMD, METIS, and NESDIS (CHOLMOD's nested dissection
1018        ordering, based on METIS and CCAMD/CCOLAMD).  This results the
1019        highest analysis time, but the lowest numerical factorization time.
1020
1021    \item \verb'UMFPACK_ORDERING_NONE' (5).
1022        The matrix is factorized as-is, except that singletons are still
1023        removed.
1024
1025    \item \verb'UMFPACK_ORDERING_USER' (6).  Use the user-ordering function
1026        passed to \newline {\tt umfpack\_*\_fsymbolic}.    Refer to
1027        \verb'UMFPACK/Source/umf_cholmod.c' for an example.
1028
1029    \end{itemize}
1030
1031    To disable the singleton filter, set \verb'Control [UMFPACK_SINGLETONS]' to
1032    0.  Disabling this search for singletons can slow UMFPACK down quite a bit
1033    for some matrices, but it does ensure that $\m{L}$ is well-conditioned and
1034    that any ill-conditioning of $\m{A}$ is captured only in $\m{U}$.
1035
1036\item {\tt umfpack\_*\_qsymbolic}:
1037
1038    An alternative to {\tt umfpack\_*\_symbolic}.  Allows the user to specify
1039    his or her own column pre-ordering, rather than using the default COLAMD
1040    or AMD pre-orderings.  For example, a graph partitioning-based order
1041    of $\m{A}\tr\m{A}$ would be suitable for UMFPACK's unsymmetric strategy.
1042    A partitioning of $\m{A}+\m{A}\tr$ would be suitable for UMFPACK's
1043    symmetric strategy.
1044
1045\item {\tt umfpack\_*\_fsymbolic}:
1046
1047    An alternative to {\tt umfpack\_*\_symbolic}.
1048    Allows the user to pass a pointer to a function,
1049    which is called to compute the ordering on the matrix
1050    (or on a submatrix with singletons removed, if any exist).
1051
1052\item {\tt umfpack\_*\_wsolve}:
1053
1054    An alternative to {\tt umfpack\_*\_solve} which does not dynamically
1055    allocate any memory.  Requires the user to pass two additional work
1056    arrays.
1057
1058\end{itemize}
1059
1060%-------------------------------------------------------------------------------
1061\subsection{Matrix manipulation routines}
1062\label{triplet}
1063%-------------------------------------------------------------------------------
1064
1065The compressed column data structure is compact, and simplifies the UMFPACK
1066routines that operate on the sparse matrix $\m{A}$.  However, it can be
1067inconvenient for the user to generate.  Section~\ref{Manipulate} presents the
1068details of routines for manipulating sparse matrices in {\em triplet} form,
1069compressed column form, and compressed row form (the transpose of the
1070compressed column form).  The triplet form of a matrix consists of three or
1071four arrays.  For the {\tt int} version of UMFPACK:
1072
1073{\footnotesize
1074\begin{verbatim}
1075     int Ti [nz] ;
1076     int Tj [nz] ;
1077     double Tx [nz] ;
1078\end{verbatim}
1079}
1080
1081For the {\tt SuiteSparse\_long} version:
1082
1083{\footnotesize
1084\begin{verbatim}
1085     SuiteSparse_long Ti [nz] ;
1086     SuiteSparse_long Tj [nz] ;
1087     double Tx [nz] ;
1088\end{verbatim}
1089}
1090
1091The complex versions use another array to hold the imaginary part:
1092
1093{\footnotesize
1094\begin{verbatim}
1095     double Tz [nz] ;
1096\end{verbatim}
1097}
1098
1099The {\tt k}-th triplet is $(i,j,a_{ij})$, where $i =$ {\tt Ti[k]},
1100$j =$ {\tt Tj[k]}, and $a_{ij} =$ {\tt Tx[k]}.  For the complex versions,
1101{\tt Tx[k]} is the real part of $a_{ij}$ and
1102{\tt Tz[k]} is the imaginary part.
1103The triplets can be in any
1104order in the {\tt Ti}, {\tt Tj}, and {\tt Tx} arrays (and {\tt Tz} for
1105the complex versions), and duplicate entries may
1106exist.
1107If {\tt Tz} is NULL, then the array {\tt Tx} becomes of size {\tt 2*nz},
1108and the real and imaginary parts of the
1109{\tt k}-th triplet are located in {\tt Tx[2*k]} and {\tt Tx[2*k+1]},
1110respectively.
1111Any duplicate entries are summed when the triplet form is converted to
1112compressed column form.  This is a convenient way to create a matrix arising in
1113finite-element methods, for example.
1114
1115Four routines are provided for manipulating sparse matrices:
1116
1117\begin{itemize}
1118\item {\tt umfpack\_*\_triplet\_to\_col}:
1119
1120    Converts a triplet form of a matrix to compressed column form (ready for
1121    input to \newline
1122    {\tt umfpack\_*\_symbolic}, {\tt umfpack\_*\_qsymbolic}, and
1123    {\tt umfpack\_*\_numeric}).  Identical to {\tt A = spconvert(i,j,x)} in
1124    MATLAB, except that zero entries are not removed, so that the pattern of
1125    entries in the compressed column form of $\m{A}$ are fully under user
1126    control.  This is important if you want to factorize a new matrix with the
1127    {\tt Symbolic} object from a prior matrix with the same pattern as the new
1128    one.
1129
1130\item {\tt umfpack\_*\_col\_to\_triplet}:
1131
1132    The opposite of {\tt umfpack\_*\_triplet\_to\_col}.  Identical to
1133    {\tt [i,j,x] = find(A)} in MATLAB, except that numerically zero entries
1134    may be included.
1135
1136\item {\tt umfpack\_*\_transpose}:
1137
1138    Transposes and optionally permutes a column form matrix \cite{Gustavson78}.
1139    Identical to
1140    {\tt R = A(P,Q)'} (linear algebraic transpose, using the complex conjugate)
1141    or {\tt R = A(P,Q).'} (the array transpose)
1142    in MATLAB, except for the presence of numerically zero entries.
1143
1144    Factorizing $\m{A}\tr$ and then solving $\m{Ax}=\m{b}$ with the transposed
1145    factors can sometimes be much faster or much slower than factorizing
1146    $\m{A}$.  It is highly dependent on your particular matrix.
1147
1148\item {\tt umfpack\_*\_scale}:
1149
1150    Applies the row scale factors to a user-provided vector.  This is not
1151    required to solve the sparse linear system $\m{Ax}=\m{b}$ or
1152    $\m{A}\tr\m{x}=\m{b}$, since {\tt umfpack\_*\_solve} applies the scale
1153    factors for those systems.
1154
1155\end{itemize}
1156
1157It is quite easy to add matrices in triplet form, subtract them, transpose
1158them, permute them, construct a submatrix, and multiply a triplet-form matrix
1159times a vector.  UMFPACK does not provide code for these basic operations,
1160however.  Refer to the discussion of
1161{\tt umfpack\_*\_triplet\_to\_col} in Section~\ref{Manipulate} for more details
1162on how to compute these operations in your own code.
1163The only primary matrix operation not provided by UMFPACK is the
1164multiplication of two sparse matrices \cite{Gustavson78}.
1165The CHOLMOD provides many of these matrix operations, which
1166can then be used in conjunction with UMFPACK.
1167See my web page for details.
1168
1169%-------------------------------------------------------------------------------
1170\subsection{Getting the contents of opaque objects}
1171%-------------------------------------------------------------------------------
1172
1173There are cases where you may wish to do more with the LU factorization
1174of a matrix than solve a linear system.  The opaque {\tt Symbolic} and
1175{\tt Numeric} objects are just that - opaque.  You cannot do anything with them
1176except to pass them back to subsequent calls to UMFPACK.  Three routines
1177are provided for copying their contents into user-provided arrays using simpler
1178data structures.  Four routines are provided for saving and loading the
1179{\tt Numeric} and {\tt Symbolic} objects to/from binary files.
1180An additional routine is provided that computes the determinant.
1181They are fully described in Section~\ref{Get}:
1182
1183\begin{itemize}
1184\item {\tt umfpack\_*\_get\_lunz}:
1185
1186    Returns the number of nonzeros in $\m{L}$ and $\m{U}$.
1187
1188\item {\tt umfpack\_*\_get\_numeric}:
1189
1190    Copies $\m{L}$, $\m{U}$, $\m{P}$, $\m{Q}$, and $\m{R}$
1191    from the {\tt Numeric} object
1192    into arrays provided by the user.  The matrix $\m{L}$ is returned in
1193    compressed row form (with the column indices in each row sorted in ascending
1194    order).  The matrix $\m{U}$ is returned in compressed column form (with
1195    sorted columns).  There are no explicit zero entries in $\m{L}$ and $\m{U}$,
1196    but such entries may exist in the {\tt Numeric} object.  The permutations
1197    $\m{P}$ and $\m{Q}$ are represented as permutation vectors, where
1198    {\tt P[k] = i} means that row {\tt i} of the original matrix is the
1199    the {\tt k}-th row of $\m{PAQ}$, and where
1200    {\tt Q[k] = j} means that column {\tt j} of the original matrix is the
1201    {\tt k}-th column of $\m{PAQ}$.  This is identical to how MATLAB uses
1202    permutation vectors (type {\tt help colamd} in MATLAB 6.1 or later).
1203
1204\item {\tt umfpack\_*\_get\_symbolic}:
1205
1206    Copies the contents of the {\tt Symbolic} object (the initial row and column
1207    preordering, supernodal column elimination tree, and information
1208    about each frontal matrix) into arrays provided by the user.
1209
1210\item {\tt umfpack\_*\_get\_determinant}:
1211
1212    Computes the determinant from the diagonal of $\m{U}$ and the permutations
1213    $\m{P}$ and $\m{Q}$.  This is mostly of theoretical interest.
1214    It is not a good test to determine if your matrix is singular or not.
1215
1216\item {\tt umfpack\_*\_save\_numeric}:
1217
1218    Saves a copy of the {\tt Numeric} object to a file, in binary format.
1219
1220\item {\tt umfpack\_*\_load\_numeric}:
1221
1222    Creates a {\tt Numeric} object by loading it from a file created
1223    by {\tt umfpack\_*\_save\_numeric}.
1224
1225\item {\tt umfpack\_*\_save\_symbolic}:
1226
1227    Saves a copy of the {\tt Symbolic} object to a file, in binary format.
1228
1229\item {\tt umfpack\_*\_load\_symbolic}:
1230
1231    Creates a {\tt Symbolic} object by loading it from a file created
1232    by {\tt umfpack\_*\_save\_symbolic}.
1233
1234\end{itemize}
1235
1236UMFPACK itself does not make use of these routines;
1237they are provided solely for returning the contents of the opaque
1238{\tt Symbolic} and {\tt Numeric} objects to the user, and saving/loading
1239them to/from a binary file.  None of them do any computation, except for
1240{\tt umfpack\_*\_get\_determinant}.
1241
1242%-------------------------------------------------------------------------------
1243\subsection{Reporting routines}
1244\label{Reporting}
1245%-------------------------------------------------------------------------------
1246
1247None of the UMFPACK routines discussed so far prints anything, even when an
1248error occurs.  UMFPACK provides you with nine routines for printing the input
1249and output arguments (including the {\tt Control} settings and {\tt Info}
1250statistics) of UMFPACK routines discussed above.  They are fully described in
1251Section~\ref{Report}:
1252
1253\begin{itemize}
1254\item {\tt umfpack\_*\_report\_status}:
1255
1256    Prints the status (return value) of other {\tt umfpack\_*} routines.
1257
1258\item {\tt umfpack\_*\_report\_info}:
1259
1260    Prints the statistics returned in the {\tt Info} array by
1261    {\tt umfpack\_*\_*symbolic},
1262    {\tt umfpack\_*\_numeric}, and {\tt umfpack\_*\_*solve}.
1263
1264\item {\tt umfpack\_*\_report\_control}:
1265
1266    Prints the {\tt Control} settings.
1267
1268\item {\tt umfpack\_*\_report\_matrix}:
1269
1270    Verifies and prints a compressed column-form or compressed row-form sparse
1271    matrix.
1272
1273\item {\tt umfpack\_*\_report\_triplet}:
1274
1275    Verifies and prints a matrix in triplet form.
1276
1277\item {\tt umfpack\_*\_report\_symbolic}:
1278
1279    Verifies and prints a {\tt Symbolic} object.
1280
1281\item {\tt umfpack\_*\_report\_numeric}:
1282
1283    Verifies and prints a {\tt Numeric} object.
1284
1285\item {\tt umfpack\_*\_report\_perm}:
1286
1287    Verifies and prints a permutation vector.
1288
1289\item {\tt umfpack\_*\_report\_vector}:
1290
1291    Verifies and prints a real or complex vector.
1292
1293\end{itemize}
1294
1295The {\tt umfpack\_*\_report\_*} routines behave slightly differently when
1296compiled
1297into the C-callable UMFPACK library than when used in the MATLAB mexFunction.
1298MATLAB stores its sparse matrices using the same compressed column data
1299structure discussed above, where row and column indices of an $m$-by-$n$
1300matrix are in the range 0 to $m-1$ or $n-1$, respectively\footnote{Complex
1301matrices in MATLAB use the split array form, with one {\tt double} array
1302for the real part and another array for the imaginary part.  UMFPACK
1303supports that format, as well as the packed complex format (new to Version 4.4).}
1304It prints them as if they are in the range 1 to $m$ or $n$.
1305The UMFPACK mexFunction behaves the same way.
1306
1307You can control how much the {\tt umfpack\_*\_report\_*} routines print by
1308modifying the {\tt Control [UMFPACK\_PRL]} parameter.  Its default value is 1.
1309Here is a summary of how the routines use this print level parameter:
1310
1311\begin{itemize}
1312\item {\tt umfpack\_*\_report\_status}:
1313
1314    No output if the print level is 0 or less, even when an error occurs.
1315    If 1, then error messages are printed, and nothing is printed if
1316    the status is {\tt UMFPACK\_OK}.  A warning message is printed if
1317    the matrix is singular.  If 2 or more, then the status is always
1318    printed.  If 4 or more, then the UMFPACK Copyright is printed.
1319    If 6 or more, then the UMFPACK License is printed.  See also the first page
1320    of this User Guide for the Copyright and License.
1321
1322\item {\tt umfpack\_*\_report\_control}:
1323
1324    No output if the print level is 1 or less.  If 2 or more, all of
1325    {\tt Control} is printed.
1326
1327\item {\tt umfpack\_*\_report\_info}:
1328
1329    No output if the print level is 1 or less.  If 2 or more, all of
1330    {\tt Info} is printed.
1331
1332\item all other {\tt umfpack\_*\_report\_*} routines:
1333
1334    If the print level is 2 or less, then these routines return silently without
1335    checking their inputs.  If 3 or more, the inputs are fully verified and a
1336    short status summary is printed.  If 4, then the first few entries of the
1337    input arguments are printed.  If 5, then all of the input arguments are
1338    printed.
1339
1340\end{itemize}
1341
1342This print level parameter has an additional effect on the MATLAB mexFunction.
1343If zero, then no warnings of singular or nearly singular matrices are
1344printed (similar to the MATLAB commands
1345{\tt warning off MATLAB:singularMatrix} and
1346{\tt warning off MATLAB:nearlySingularMatrix}).
1347
1348%-------------------------------------------------------------------------------
1349\subsection{Utility routines}
1350%-------------------------------------------------------------------------------
1351
1352Three timing routines are provided in UMFPACK Version 4.1 and following,
1353{\tt umfpack\_tic}, {\tt umfpack\_toc}, and {\tt umfpack\_timer}.
1354These three routines are the only user-callable
1355routine that is identical in all four {\tt int}/{\tt SuiteSparse\_long}, real/complex
1356versions (there is no {\tt umfpack\_di\_timer} routine, for example).
1357
1358%-------------------------------------------------------------------------------
1359\subsection{Control parameters}
1360\label{control_param}
1361%-------------------------------------------------------------------------------
1362
1363UMFPACK uses an optional {\tt double} array (currently of size 20)
1364to modify its control parameters.  If you pass {\tt (double *) NULL} instead
1365of a {\tt Control} array, then defaults are used.  These defaults provide
1366nearly optimal performance (both speed, memory usage, and numerical accuracy)
1367for a wide range of matrices from real applications.
1368
1369This array will almost certainly grow in size in future releases,
1370so be sure to dimension your {\tt Control} array to be of size
1371{\tt UMFPACK\_CONTROL}.  That constant is currently defined to be 20,
1372but may increase in future versions, since all 20 entries are in use.
1373
1374The contents of this array may be modified by the user
1375(see {\tt umfpack\_*\_defaults}).  Each
1376user-callable routine includes a complete description of how each control
1377setting modifies its behavior.  Table~\ref{control} summarizes the entire
1378contents of the {\tt Control} array.
1379Note that ANSI C uses 0-based indexing, while MATLAB uses 1-based
1380indexing.  Thus, {\tt Control(1)} in MATLAB is the same as
1381{\tt Control[0]} or {\tt Control[UMFPACK\_PRL]} in ANSI C.
1382
1383\begin{table}
1384\caption{UMFPACK Control parameters}
1385\label{control}
1386{\footnotesize
1387\begin{tabular}{llll}
1388\hline
1389
1390MATLAB & ANSI C & default & description \\
1391\verb'struct' & & & \\
1392\hline
1393{\tt prl}           & {\tt Control[UMFPACK\_PRL]} & 1 & printing level \\
1394-                   & {\tt Control[UMFPACK\_DENSE\_ROW]} & 0.2 & dense row parameter \\
1395-                   & {\tt Control[UMFPACK\_DENSE\_COL]} & 0.2 & dense column parameter \\
1396{\tt tol}           & {\tt Control[UMFPACK\_PIVOT\_TOLERANCE]} & 0.1 & partial pivoting tolerance \\
1397-                   & {\tt Control[UMFPACK\_BLOCK\_SIZE]} & 32 & BLAS block size \\
1398{\tt strategy}      & {\tt Control[UMFPACK\_STRATEGY]} & 0 (auto) & select strategy \\
1399{\tt ordering}      & {\tt Control[UMFPACK\_ORDERING]} & 1 (AMD) & select ordering \\
1400-                   & {\tt Control[UMFPACK\_ALLOC\_INIT]} & 0.7 & initial memory allocation  \\
1401{\tt irstep}        & {\tt Control[UMFPACK\_IRSTEP]} & 2 & max iter. refinement steps \\
1402% {\tt Control(13)} & ignored & & \\
1403-                   & {\tt Control[UMFPACK\_FIXQ]} & 0 (auto) & fix or modify Q \\
1404-                   & {\tt Control[UMFPACK\_AMD\_DENSE]} & 10 & AMD dense row/col param. \\
1405{\tt symtol}        & {\tt Control[UMFPACK\_SYM\_PIVOT\_TOLERANCE]} & 0.001 & for diagonal entries \\
1406{\tt scale}         & {\tt Control[UMFPACK\_SCALE]} & 1 (sum) & row scaling (none, sum, max) \\
1407-                   & {\tt Control[UMFPACK\_FRONT\_ALLOC\_INIT]} & 0.5 & frontal matrix allocation ratio \\
1408-                   & {\tt Control[UMFPACK\_DROPTOL]} & 0 & drop tolerance \\
1409-                   & {\tt Control[UMFPACK\_AGGRESSIVE]} & 1 (yes) & aggressive absorption \\
1410{\tt singletons}    & {\tt Control[UMFPACK\_SINGLETONS]} & 1 (enable) & enable singleton filter \\
1411%
1412\hline
1413\end{tabular}
1414}
1415\end{table}
1416
1417Let $\alpha_r = ${\tt Control [UMFPACK\_DENSE\_ROW]},
1418    $\alpha_c = ${\tt Control [UMFPACK\_DENSE\_COL]}, and
1419    $\alpha = ${\tt Control [UMFPACK\_AMD\_DENSE]}.
1420Suppose the submatrix $\m{S}$, obtained after eliminating pivots with
1421zero Markowitz cost, is $m$-by-$n$.
1422Then a row is considered ``dense'' if it has more than
1423$\max (16, 16 \alpha_r \sqrt{n})$ entries.
1424A column is considered ``dense'' if it has more than
1425$\max (16, 16 \alpha_c \sqrt{m})$ entries.
1426These rows and columns are treated different in COLAMD and during numerical
1427factorization.   In COLAMD, dense columns are placed last in their natural
1428order, and dense rows are ignored.  During numerical factorization, dense
1429rows are stored differently.
1430In AMD, a row/column of the square matrix $\m{S}+\m{S}\tr$ is
1431considered ``dense'' if it has more than $\max (16, \alpha \sqrt{n})$ entries.
1432These rows/columns are placed last in AMD's output ordering.
1433For more details on the control parameters, refer to the documentation of
1434{\tt umfpack\_*\_qsymbolic},
1435{\tt umfpack\_*\_fsymbolic},
1436{\tt umfpack\_*\_numeric}, {\tt umfpack\_*\_solve},
1437and the {\tt umfpack\_*\_report\_*} routines,
1438in Sections~\ref{Primary}~through~\ref{Report}, below.
1439
1440%-------------------------------------------------------------------------------
1441\subsection{Error codes}
1442\label{error_codes}
1443%-------------------------------------------------------------------------------
1444
1445Many of the routines return a {\tt status} value.
1446This is also returned as the first entry in the {\tt Info} array, for
1447those routines with that argument.  The following list summarizes
1448all of the error codes in UMFPACK.  Each error code is given a
1449specific name in the {\tt umfpack.h} include file, so you can use
1450those constants instead of hard-coded values in your program.
1451Future versions may report additional error codes.
1452
1453A value of zero means everything was successful, and the matrix is
1454non-singular.  A value greater than zero means the routine was successful,
1455but a warning occurred.
1456A negative value means the routine was not successful.
1457In this case, no {\tt Symbolic} or {\tt Numeric} object was created.
1458
1459\begin{itemize}
1460\item {\tt UMFPACK\_OK},  (0):  UMFPACK was successful.
1461
1462\item {\tt UMFPACK\_WARNING\_singular\_matrix},  (1):  Matrix is singular.
1463    There are exact zeros on the diagonal of $\m{U}$.
1464
1465\item {\tt UMFPACK\_WARNING\_determinant\_underflow}, (2):
1466    The determinant is nonzero, but smaller in magnitude than
1467    the smallest positive floating-point number.
1468
1469\item {\tt UMFPACK\_WARNING\_determinant\_overflow}, (3):
1470    The determinant is larger in magnitude than
1471    the largest positive floating-point number (IEEE Inf).
1472
1473\item {\tt UMFPACK\_ERROR\_out\_of\_memory},  (-1):  Not enough memory.
1474    The ANSI C {\tt malloc} or {\tt realloc} routine failed.
1475
1476\item {\tt UMFPACK\_ERROR\_invalid\_Numeric\_object},  (-3):
1477    Routines that take a {\tt Numeric} object as input (or load it
1478    from a file) check this object and return this error code if it is
1479    invalid.  This can be caused by a memory leak or overrun in your
1480    program, which can overwrite part of the Numeric object.  It can also
1481    be caused by passing a Symbolic object by mistake, or some other pointer.
1482    If you try to factorize a matrix using one version of UMFPACK and
1483    then use the factors in another version, this error code will trigger as
1484    well.  You cannot factor your matrix using
1485    version 4.0 and then solve with version 4.1, for example.\footnote{
1486    Exception: v4.3, v4.3.1, and v4.4 use identical data structures
1487    for the {\tt Numeric} and {\tt Symbolic} objects}.
1488    You cannot use different precisions of the same version
1489    (real and complex, for example).
1490    It is possible for the {\tt Numeric} object to be corrupted by your
1491    program in subtle ways that are not detectable by this quick check.
1492    In this case, you may see an
1493    {\tt UMFPACK\_ERROR\_different\_pattern} error code, or even an
1494    {\tt UMFPACK\_ERROR\_internal\_error}.
1495
1496\item {\tt UMFPACK\_ERROR\_invalid\_Symbolic\_object},  (-4):
1497    Routines that take a {\tt Symbolic} object as input (or load it
1498    from a file) check this object and return this error code if it is
1499    invalid.  The causes of this error are analogous to the
1500    {\tt UMFPACK\_ERROR\_invalid\_Numeric\_object} error described above.
1501
1502\item {\tt UMFPACK\_ERROR\_argument\_missing},  (-5):
1503    Some arguments of some are optional (you can pass a {\tt NULL} pointer
1504    instead of an array).  This error code occurs if you pass a {\tt NULL}
1505    pointer when that argument is required to be present.
1506
1507\item {\tt UMFPACK\_ERROR\_n\_nonpositive}  (-6):
1508    The number of rows or columns of the matrix must be greater than zero.
1509
1510\item {\tt UMFPACK\_ERROR\_invalid\_matrix}  (-8):
1511    The matrix is invalid.  For the column-oriented input, this error
1512    code will occur if the contents of {\tt Ap} and/or {\tt Ai} are invalid.
1513
1514    {\tt Ap} is an integer array of size {\tt n\_col+1}.
1515    On input, it holds the
1516    ``pointers'' for the column form of the sparse matrix $\m{A}$.
1517    Column {\tt j} of
1518    the matrix A is held in {\tt Ai [(Ap [j])} \ldots {\tt (Ap [j+1]-1)]}.
1519    The first entry, {\tt Ap [0]}, must be zero,
1520    and {\tt Ap [j]} $\le$ {\tt Ap [j+1]} must hold for all
1521    {\tt j} in the range 0 to {\tt n\_col-1}.
1522    The value {\tt nz = Ap [n\_col]} is thus the
1523    total number of entries in the pattern of the matrix A.
1524    {\tt nz} must be greater than or equal to zero.
1525
1526    The nonzero pattern (row indices) for column {\tt j} is stored in
1527    {\tt Ai [(Ap [j])} \ldots {\tt (Ap [j+1]-1)]}.  The row indices in a given
1528    column {\tt j}
1529    must be in ascending order, and no duplicate row indices may be present.
1530    Row indices must be in the range 0 to {\tt n\_row-1}
1531    (the matrix is 0-based).
1532
1533    Some routines take a triplet-form input, with arguments
1534    {\tt nz}, {\tt Ti}, and {\tt Tj}.  This error code is returned
1535    if {\tt nz} is less than zero,
1536    if any row    index in {\tt Ti} is outside the range 0 to {\tt n\_col-1}, or
1537    if any column index in {\tt Tj} is outside the range 0 to {\tt n\_row-1}.
1538
1539\item {\tt UMFPACK\_ERROR\_different\_pattern},  (-11):
1540    The most common cause of this error is that the pattern of the
1541    matrix has changed between the symbolic and numeric factorization.
1542    It can also occur if the {\tt Numeric} or {\tt Symbolic} object has
1543    been subtly corrupted by your program.
1544
1545\item {\tt UMFPACK\_ERROR\_invalid\_system},  (-13):
1546    The {\tt sys} argument provided to one of the solve routines is invalid.
1547
1548\item {\tt UMFPACK\_ERROR\_invalid\_permutation},  (-15):
1549    The permutation vector provided as input is invalid.
1550
1551\item {\tt UMFPACK\_ERROR\_file\_IO},  (-17):
1552    This error code is returned by the routines that save and load
1553    the {\tt Numeric} or {\tt Symbolic} objects to/from a file, if a
1554    file I/O error has occurred.  The file may not exist or may not be readable,
1555    you may be trying to create a file that you don't have permission to create,
1556    or you may be out of disk space.  The file you are trying to read might
1557    be the wrong one, and an earlier end-of-file condition would then result
1558    in this error.
1559
1560\item {\tt UMFPACK\_ERROR\_ordering\_failed},  (-18):
1561    The ordering method failed.
1562
1563\item {\tt UMFPACK\_ERROR\_internal\_error},  (-911):
1564    An internal error has occurred, of unknown cause.  This is either a bug
1565    in UMFPACK, or the result of a memory overrun from your program.
1566    Try modifying the file {\tt AMD/Include/amd\_internal.h} and adding
1567    the statement {\tt \#undef NDEBUG}, to enable the debugging mode.
1568    Recompile UMFPACK and rerun your program.
1569    A failed assertion might occur which
1570    can give you a better indication as to what is going wrong.  Be aware that
1571    UMFPACK will be extraordinarily slow when running in debug mode.
1572    If all else fails, contact the developer (DrTimothyAldenDavis@gmail.com)
1573    with as many details as possible.
1574
1575\end{itemize}
1576
1577%-------------------------------------------------------------------------------
1578\subsection{Larger examples}
1579%-------------------------------------------------------------------------------
1580
1581Full examples of all user-callable UMFPACK routines
1582are available in four stand-alone C main programs, {\tt umfpack\_*\_demo.c}.
1583Another example is
1584the UMFPACK mexFunction, {\tt umfpackmex.c}.  The mexFunction accesses only the
1585user-callable C interface to UMFPACK.  The only features that it does not use
1586are the support for the triplet form (MATLAB's sparse arrays are already in the
1587compressed column form) and the ability to reuse the {\tt Symbolic} object to
1588numerically factorize a matrix whose pattern is the same as a prior matrix
1589analyzed by {\tt umfpack\_*\_symbolic},
1590{\tt umfpack\_*\_qsymbolic} or
1591{\tt umfpack\_*\_fsymbolic}.
1592The
1593latter is an important feature, but the mexFunction does not return its opaque
1594{\tt Symbolic} and {\tt Numeric} objects to MATLAB.  Instead, it gets the
1595contents of these objects after extracting them via the {\tt umfpack\_*\_get\_*}
1596routines, and returns them as MATLAB sparse matrices.
1597
1598The {\tt umf4.c} program for reading matrices in Harwell/Boeing format
1599\cite{DuffGrimesLewis87b} is provided.  It requires three Fortran 77 programs
1600({\tt readhb.f}, {\tt readhb\_nozeros.f}, and {\tt readhb\_size.f})
1601for reading in the sample Harwell/Boeing files in the {\tt UMFPACK/Demo/HB}
1602directory.  More matrices are available at
1603http://www.suitesparse.com.
1604Type {\tt make hb} in the {\tt UMFPACK/Demo/HB} directory
1605to compile and run this demo.  This program was used for the experimental
1606results in \cite{Davis03}.
1607
1608%-------------------------------------------------------------------------------
1609\section{Synopsis of C-callable routines}
1610\label{Synopsis}
1611%-------------------------------------------------------------------------------
1612
1613Each subsection, below, summarizes the input variables, output variables, return
1614values, and calling sequences of the routines in one category.  Variables with
1615the same name as those already listed in a prior category have the same size
1616and type.
1617
1618The real, {\tt SuiteSparse\_long} integer {\tt umfpack\_dl\_*} routines are
1619identical to the real, {\tt int} routines, except that {\tt \_di\_} is replaced
1620with {\tt \_dl\_} in the name, and all {\tt int} arguments become
1621{\tt SuiteSparse\_long}.
1622Similarly, the complex, {\tt SuiteSparse\_long} integer {\tt umfpack\_zl\_*} routines are
1623identical to the complex, {\tt int} routines, except that {\tt \_zi\_} is
1624replaced
1625with {\tt \_zl\_} in the name, and all {\tt int} arguments become
1626{\tt SuiteSparse\_long}.
1627Only the real and complex {\tt int} versions are listed in the synopsis below.
1628
1629The matrix $\m{A}$ is {\tt m}-by-{\tt n} with {\tt nz} entries.
1630
1631The {\tt sys} argument of {\tt umfpack\_*\_solve}
1632is an integer in the range 0 to 14 which defines which linear system is
1633to be solved.
1634\footnote{Integer values for {\tt sys} are used instead of strings (as in LINPACK
1635and LAPACK) to avoid C-to-Fortran portability issues.}
1636Valid values are listed in Table~\ref{sys}.
1637The notation $\m{A}\he$ refers to the matrix transpose, which is the
1638complex conjugate transpose for complex matrices ({\tt A'} in MATLAB).
1639The array transpose is $\m{A}\tr$, which is {\tt A.'} in MATLAB.
1640
1641\begin{table}
1642\begin{center}
1643\caption{UMFPACK {\tt sys} parameter}
1644\label{sys}
1645{\footnotesize
1646\begin{tabular}{ll|l}
1647\hline
1648Value & & system \\
1649\hline
1650& & \\
1651{\tt UMFPACK\_A}      &  (0) & $\m{Ax}=\m{b}$ \\
1652{\tt UMFPACK\_At}     &  (1) & $\m{A}\he\m{x}=\m{b}$ \\
1653{\tt UMFPACK\_Aat}    &  (2) & $\m{A}\tr\m{x}=\m{b}$ \\
1654& & \\
1655\hline
1656& & \\
1657{\tt UMFPACK\_Pt\_L}  &  (3) & $\m{P}\tr\m{Lx}=\m{b}$ \\
1658{\tt UMFPACK\_L}      &  (4) & $\m{Lx}=\m{b}$ \\
1659{\tt UMFPACK\_Lt\_P}  &  (5) & $\m{L}\he\m{Px}=\m{b}$ \\
1660{\tt UMFPACK\_Lat\_P} &  (6) & $\m{L}\tr\m{Px}=\m{b}$ \\
1661{\tt UMFPACK\_Lt}     &  (7) & $\m{L}\he\m{x}=\m{b}$ \\
1662{\tt UMFPACK\_Lat}    &  (8) & $\m{L}\tr\m{x}=\m{b}$ \\
1663& & \\
1664\hline
1665& & \\
1666{\tt UMFPACK\_U\_Qt}  &  (9) & $\m{UQ}\tr\m{x}=\m{b}$ \\
1667{\tt UMFPACK\_U}      & (10) & $\m{Ux}=\m{b}$ \\
1668{\tt UMFPACK\_Q\_Ut}  & (11) & $\m{QU}\he\m{x}=\m{b}$ \\
1669{\tt UMFPACK\_Q\_Uat} & (12) & $\m{QU}\tr\m{x}=\m{b}$ \\
1670{\tt UMFPACK\_Ut}     & (13) & $\m{U}\he\m{x}=\m{b}$ \\
1671{\tt UMFPACK\_Uat}    & (14) & $\m{U}\tr\m{x}=\m{b}$ \\
1672& & \\
1673\hline
1674\end{tabular}
1675}
1676\end{center}
1677\end{table}
1678
1679%-------------------------------------------------------------------------------
1680\subsection{Primary routines: real/{\tt int}}
1681%-------------------------------------------------------------------------------
1682
1683{\footnotesize
1684\begin{verbatim}
1685#include "umfpack.h"
1686int status, sys, n, m, nz, Ap [n+1], Ai [nz] ;
1687double Control [UMFPACK_CONTROL], Info [UMFPACK_INFO], Ax [nz], X [n], B [n] ;
1688void *Symbolic, *Numeric ;
1689
1690status = umfpack_di_symbolic (m, n, Ap, Ai, Ax, &Symbolic, Control, Info) ;
1691status = umfpack_di_numeric (Ap, Ai, Ax, Symbolic, &Numeric, Control, Info) ;
1692status = umfpack_di_solve (sys, Ap, Ai, Ax, X, B, Numeric, Control, Info) ;
1693umfpack_di_free_symbolic (&Symbolic) ;
1694umfpack_di_free_numeric (&Numeric) ;
1695\end{verbatim}
1696}
1697
1698%-------------------------------------------------------------------------------
1699\subsection{Alternative routines: real/{\tt int}}
1700%-------------------------------------------------------------------------------
1701
1702{\footnotesize
1703\begin{verbatim}
1704int Qinit [n], Wi [n] ;
1705double W [5*n] ;
1706
1707umfpack_di_defaults (Control) ;
1708status = umfpack_di_qsymbolic (m, n, Ap, Ai, Ax, Qinit, &Symbolic, Control, Info) ;
1709status = umfpack_di_fsymbolic (m, n, Ap, Ai, Ax, &user_ordering, user_params, &Symbolic, Control, Info) ;
1710status = umfpack_di_wsolve (sys, Ap, Ai, Ax, X, B, Numeric, Control, Info, Wi, W) ;
1711\end{verbatim}
1712}
1713
1714%-------------------------------------------------------------------------------
1715\subsection{Matrix manipulation routines: real/{\tt int}}
1716%-------------------------------------------------------------------------------
1717
1718{\footnotesize
1719\begin{verbatim}
1720int Ti [nz], Tj [nz], P [m], Q [n], Rp [m+1], Ri [nz], Map [nz] ;
1721double Tx [nz], Rx [nz], Y [m], Z [m] ;
1722
1723status = umfpack_di_col_to_triplet (n, Ap, Tj) ;
1724status = umfpack_di_triplet_to_col (m, n, nz, Ti, Tj, Tx, Ap, Ai, Ax, Map) ;
1725status = umfpack_di_transpose (m, n, Ap, Ai, Ax, P, Q, Rp, Ri, Rx) ;
1726status = umfpack_di_scale (Y, Z, Numeric) ;
1727\end{verbatim}
1728}
1729
1730%-------------------------------------------------------------------------------
1731\subsection{Getting the contents of opaque objects: real/{\tt int}}
1732%-------------------------------------------------------------------------------
1733
1734The {\tt filename} string should be large enough to hold the name of a file.
1735
1736{\footnotesize
1737\begin{verbatim}
1738int lnz, unz, Lp [m+1], Lj [lnz], Up [n+1], Ui [unz], do_recip ;
1739double Lx [lnz], Ux [unz], D [min (m,n)], Rs [m], Mx [1], Ex [1] ;
1740int nfr, nchains, P1 [m], Q1 [n], Front_npivcol [n+1], Front_parent [n+1], Front_1strow [n+1],
1741    Front_leftmostdesc [n+1], Chain_start [n+1], Chain_maxrows [n+1], Chain_maxcols [n+1] ;
1742char filename [100] ;
1743
1744status = umfpack_di_get_lunz (&lnz, &unz, &m, &n, &nz_udiag, Numeric) ;
1745status = umfpack_di_get_numeric (Lp, Lj, Lx, Up, Ui, Ux, P, Q, D,
1746    &do_recip, Rs, Numeric) ;
1747status = umfpack_di_get_symbolic (&m, &n, &n1, &nz, &nfr, &nchains, P1, Q1,
1748    Front_npivcol, Front_parent, Front_1strow, Front_leftmostdesc,
1749    Chain_start, Chain_maxrows, Chain_maxcols, Symbolic) ;
1750status = umfpack_di_load_numeric (&Numeric, filename) ;
1751status = umfpack_di_save_numeric (Numeric, filename) ;
1752status = umfpack_di_load_symbolic (&Symbolic, filename) ;
1753status = umfpack_di_save_symbolic (Symbolic, filename) ;
1754status = umfapck_di_get_determinant (Mx, Ex, Numeric, Info) ;
1755\end{verbatim}
1756}
1757
1758%-------------------------------------------------------------------------------
1759\subsection{Reporting routines: real/{\tt int}}
1760%-------------------------------------------------------------------------------
1761
1762{\footnotesize
1763\begin{verbatim}
1764
1765umfpack_di_report_status (Control, status) ;
1766umfpack_di_report_control (Control) ;
1767umfpack_di_report_info (Control, Info) ;
1768status = umfpack_di_report_matrix (m, n, Ap, Ai, Ax, 1, Control) ;
1769status = umfpack_di_report_matrix (m, n, Rp, Ri, Rx, 0, Control) ;
1770status = umfpack_di_report_numeric (Numeric, Control) ;
1771status = umfpack_di_report_perm (m, P, Control) ;
1772status = umfpack_di_report_perm (n, Q, Control) ;
1773status = umfpack_di_report_symbolic (Symbolic, Control) ;
1774status = umfpack_di_report_triplet (m, n, nz, Ti, Tj, Tx, Control) ;
1775status = umfpack_di_report_vector (n, X, Control) ;
1776\end{verbatim}
1777}
1778
1779
1780
1781
1782
1783
1784%-------------------------------------------------------------------------------
1785\subsection{Primary routines: complex/{\tt int}}
1786%-------------------------------------------------------------------------------
1787
1788{\footnotesize
1789\begin{verbatim}
1790double Az [nz], Xx [n], Xz [n], Bx [n], Bz [n] ;
1791
1792status = umfpack_zi_symbolic (m, n, Ap, Ai, Ax, Az, &Symbolic, Control, Info) ;
1793status = umfpack_zi_numeric (Ap, Ai, Ax, Az, Symbolic, &Numeric, Control, Info) ;
1794status = umfpack_zi_solve (sys, Ap, Ai, Ax, Az, Xx, Xz, Bx, Bz, Numeric, Control, Info) ;
1795umfpack_zi_free_symbolic (&Symbolic) ;
1796umfpack_zi_free_numeric (&Numeric) ;
1797\end{verbatim}
1798}
1799
1800The arrays {\tt Ax}, {\tt Bx}, and {\tt Xx} double in size if
1801any imaginary argument ({\tt Az}, {\tt Xz}, or {\tt Bz}) is {\tt NULL}.
1802
1803%-------------------------------------------------------------------------------
1804\subsection{Alternative routines: complex/{\tt int}}
1805%-------------------------------------------------------------------------------
1806
1807{\footnotesize
1808\begin{verbatim}
1809double Wz [10*n] ;
1810
1811umfpack_zi_defaults (Control) ;
1812status = umfpack_zi_qsymbolic (m, n, Ap, Ai, Ax, Az, Qinit, &Symbolic, Control, Info) ;
1813status = umfpack_zi_fsymbolic (m, n, Ap, Ai, Ax, Az, &user_ordering, user_params, &Symbolic, Control, Info) ;
1814status = umfpack_zi_wsolve (sys, Ap, Ai, Ax, Az, Xx, Xz, Bx, Bz, Numeric, Control, Info, Wi, Wz) ;
1815\end{verbatim}
1816}
1817
1818%-------------------------------------------------------------------------------
1819\subsection{Matrix manipulation routines: complex/{\tt int}}
1820%-------------------------------------------------------------------------------
1821
1822{\footnotesize
1823\begin{verbatim}
1824double Tz [nz], Rz [nz], Yx [m], Yz [m], Zx [m], Zz [m] ;
1825
1826status = umfpack_zi_col_to_triplet (n, Ap, Tj) ;
1827status = umfpack_zi_triplet_to_col (m, n, nz, Ti, Tj, Tx, Tz, Ap, Ai, Ax, Az, Map) ;
1828status = umfpack_zi_transpose (m, n, Ap, Ai, Ax, Az, P, Q, Rp, Ri, Rx, Rz, 1) ;
1829status = umfpack_zi_transpose (m, n, Ap, Ai, Ax, Az, P, Q, Rp, Ri, Rx, Rz, 0) ;
1830status = umfpack_zi_scale (Yx, Yz, Zx, Zz, Numeric) ;
1831\end{verbatim}
1832}
1833
1834The arrays {\tt Tx}, {\tt Rx}, {\tt Yx}, and {\tt Zx} double in size if
1835any imaginary argument ({\tt Tz}, {\tt Rz}, {\tt Yz}, or {\tt Zz}) is {\tt NULL}.
1836
1837%-------------------------------------------------------------------------------
1838\subsection{Getting the contents of opaque objects: complex/{\tt int}}
1839%-------------------------------------------------------------------------------
1840
1841{\footnotesize
1842\begin{verbatim}
1843double Lz [lnz], Uz [unz], Dx [min (m,n)], Dz [min (m,n)], Mz [1] ;
1844
1845status = umfpack_zi_get_lunz (&lnz, &unz, &m, &n, &nz_udiag, Numeric) ;
1846status = umfpack_zi_get_numeric (Lp, Lj, Lx, Lz, Up, Ui, Ux, Uz, P, Q, Dx, Dz,
1847    &do_recip, Rs, Numeric) ;
1848status = umfpack_zi_get_symbolic (&m, &n, &n1, &nz, &nfr, &nchains, P1, Q1,
1849    Front_npivcol, Front_parent, Front_1strow, Front_leftmostdesc,
1850    Chain_start, Chain_maxrows, Chain_maxcols, Symbolic) ;
1851status = umfpack_zi_load_numeric (&Numeric, filename) ;
1852status = umfpack_zi_save_numeric (Numeric, filename) ;
1853status = umfpack_zi_load_symbolic (&Symbolic, filename) ;
1854status = umfpack_zi_save_symbolic (Symbolic, filename) ;
1855status = umfapck_zi_get_determinant (Mx, Mz, Ex, Numeric, Info) ;
1856\end{verbatim}
1857}
1858
1859The arrays {\tt Lx}, {\tt Ux}, {\tt Dx}, and {\tt Mx} double in size if
1860any imaginary argument ({\tt Lz}, {\tt Uz}, {\tt Dz}, or {\tt Mz}) is {\tt NULL}.
1861
1862%-------------------------------------------------------------------------------
1863\subsection{Reporting routines: complex/{\tt int}}
1864%-------------------------------------------------------------------------------
1865
1866{\footnotesize
1867\begin{verbatim}
1868
1869umfpack_zi_report_status (Control, status) ;
1870umfpack_zi_report_control (Control) ;
1871umfpack_zi_report_info (Control, Info) ;
1872status = umfpack_zi_report_matrix (m, n, Ap, Ai, Ax, Az, 1, Control) ;
1873status = umfpack_zi_report_matrix (m, n, Rp, Ri, Rx, Rz, 0, Control) ;
1874status = umfpack_zi_report_numeric (Numeric, Control) ;
1875status = umfpack_zi_report_perm (m, P, Control) ;
1876status = umfpack_zi_report_perm (n, Q, Control) ;
1877status = umfpack_zi_report_symbolic (Symbolic, Control) ;
1878status = umfpack_zi_report_triplet (m, n, nz, Ti, Tj, Tx, Tz, Control) ;
1879status = umfpack_zi_report_vector (n, Xx, Xz, Control) ;
1880\end{verbatim}
1881}
1882
1883The arrays {\tt Ax}, {\tt Rx}, {\tt Tx}, and {\tt Xx} double in size if
1884any imaginary argument ({\tt Az}, {\tt Rz}, {\tt Tz}, or {\tt Xz}) is {\tt NULL}.
1885
1886
1887
1888
1889%-------------------------------------------------------------------------------
1890\subsection{Utility routines}
1891%-------------------------------------------------------------------------------
1892
1893These routines are the same in all four versions of UMFPACK.
1894
1895{\footnotesize
1896\begin{verbatim}
1897double t, s [2] ;
1898
1899t = umfpack_timer ( ) ;
1900umfpack_tic (s) ;
1901umfpack_toc (s) ;
1902
1903\end{verbatim}
1904}
1905
1906%-------------------------------------------------------------------------------
1907\subsection{AMD ordering routines}
1908%-------------------------------------------------------------------------------
1909
1910UMFPACK makes use of the AMD ordering package for its symmetric ordering
1911strategy.  You may also use these four user-callable routines in your own C
1912programs.  You need to include the {\tt amd.h} file only if you make direct
1913calls to the AMD routines themselves.  The {\tt int} versions are summarized
1914below; {\tt SuiteSparse\_long} versions are also available.  Refer to the AMD User Guide
1915for more information, or to the file {\tt amd.h} which documents these routines.
1916
1917{\footnotesize
1918\begin{verbatim}
1919#include "amd.h"
1920double amd_control [AMD_CONTROL], amd_info [AMD_INFO] ;
1921
1922amd_defaults (amd_control) ;
1923status = amd_order (n, Ap, Ai, P, amd_control, amd_info) ;
1924amd_control (amd_control) ;
1925amd_info (amd_info) ;
1926
1927\end{verbatim}
1928}
1929
1930%-------------------------------------------------------------------------------
1931\section{Using UMFPACK in a Fortran program}
1932%-------------------------------------------------------------------------------
1933
1934UMFPACK includes a basic Fortran 77 interface to some of the C-callable
1935UMFPACK routines.
1936Since interfacing C and Fortran programs is not portable, this interface might
1937not work with all C and Fortran compilers.  Refer to Section~\ref{Install} for
1938more details.  The following Fortran routines are provided.
1939The list includes the C-callable routines that the Fortran interface
1940routine calls.  Refer to the corresponding C routines in Section~\ref{C} for
1941more details on what the Fortran routine does.
1942
1943\begin{itemize}
1944\item {\tt umf4def}: sets the default control parameters
1945    ({\tt umfpack\_di\_defaults}).
1946
1947\item {\tt umf4sym}: pre-ordering and symbolic factorization
1948    ({\tt umfpack\_di\_symbolic}).
1949
1950\item {\tt umf4num}: numeric factorization
1951    ({\tt umfpack\_di\_numeric}).
1952
1953\item {\tt umf4solr}: solve a linear system with iterative refinement
1954    ({\tt umfpack\_di\_solve}).
1955
1956\item {\tt umf4sol}: solve a linear system without iterative refinement
1957    ({\tt umfpack\_di\_solve}).  Sets {\tt Control [UMFPACK\_IRSTEP]}
1958    to zero, and does not require the matrix $\m{A}$.
1959
1960\item {\tt umf4scal}: scales a vector using UMFPACK's scale factors
1961    ({\tt umfpack\_di\_scale}).
1962
1963\item {\tt umf4fnum}: free the {\tt Numeric} object
1964    ({\tt umfpack\_di\_free\_numeric}).
1965
1966\item {\tt umf4fsym}: free the {\tt Symbolic} object
1967    ({\tt umfpack\_di\_free\_symbolic}).
1968
1969\item {\tt umf4pcon}: prints the control parameters
1970    ({\tt umfpack\_di\_report\_control}).
1971
1972\item {\tt umf4pinf}: print statistics
1973    ({\tt umfpack\_di\_report\_info}).
1974
1975\item {\tt umf4snum}: save the {\tt Numeric} object to a file
1976    ({\tt umfpack\_di\_save\_numeric}).
1977
1978\item {\tt umf4ssym}: save the {\tt Symbolic} object to a file
1979    ({\tt umfpack\_di\_save\_symbolic}).
1980
1981\item {\tt umf4lnum}: load the {\tt Numeric} object from a file
1982    ({\tt umfpack\_di\_load\_numeric}).
1983
1984\item {\tt umf4lsym}: load the {\tt Symbolic} object from a file
1985    ({\tt umfpack\_di\_load\_symbolic}).
1986\end{itemize}
1987
1988The matrix $\m{A}$ is passed to UMFPACK in compressed column form, with 0-based
1989indices.  In Fortran, for an {\tt m}-by-{\tt n} matrix $\m{A}$ with {\tt nz}
1990entries, the row indices of the first column (column 1) are in
1991{\tt Ai (Ap(1)+1} \ldots {\tt Ap(2))}, with values in
1992{\tt Ax (Ap(1)+1} \ldots {\tt Ap(2))}.  The last column (column {\tt n}) is in
1993{\tt Ai (Ap(n)+1} \ldots {\tt Ap(n+1))} and
1994{\tt Ax (Ap(n)+1} \ldots {\tt Ap(n+1))}.
1995The number of entries in the matrix is thus {\tt nz = Ap (n+1)}.
1996The row indices in {\tt Ai} are in the range 0 to {\tt m}-1.  They must be
1997sorted, with no duplicate entries allowed.  None of the UMFPACK routines
1998modify the input matrix $\m{A}$.
1999The following definitions apply for the Fortran routines:
2000
2001{\footnotesize
2002\begin{verbatim}
2003    integer m, n, Ap (n+1), Ai (nz), symbolic, numeric, filenum, status
2004    double precision Ax (nz), control (20), info (90), x (n), b (n)
2005\end{verbatim}
2006}
2007
2008UMFPACK's status is returned in either a {\tt status} argument, or in
2009{\tt info (1)}.
2010It is zero if UMFPACK was successful, 1 if the matrix is singular (this is a
2011warning, not an error), and negative if an error occurred.
2012Section~\ref{error_codes} summarizes the possible values of {\tt status}
2013and {\tt info (1)}.
2014See Table~\ref{sys} for a list of the values of the {\tt sys} argument.
2015See Table~\ref{control} for a list of the control parameters (the
2016Fortran usage is the same as the MATLAB usage for this array).
2017
2018For the {\tt Numeric} and {\tt Symbolic} handles, it is probably safe to
2019assume that a Fortran {\tt integer} is sufficient to store a C pointer.  If
2020that does not work, try defining {\tt numeric} and {\tt symbolic} in your
2021Fortran program as integer arrays of size 2.  You will need to define them
2022as {\tt integer*8} if you compile UMFPACK in the 64-bit mode.
2023
2024To avoid passing strings between C and Fortran in the load/save routines,
2025a file number is passed instead, and the C interface constructs a file name
2026(if {\tt filenum} is 42, the {\tt Numeric} file name is {\tt n42.umf}, and
2027the {\tt Symbolic} file name is {\tt s42.umf}).
2028
2029The following is a summary of the calling sequence of each Fortran
2030interface routine.  An example of their use is in the {\tt Demo/umf4hb.f}
2031file.  That routine also includes an example of how to convert a 1-based
2032sparse matrix into 0-based form.  For more details on the arguments of each
2033routine, refer to the arguments of the same name in the corresponding
2034C-callable routine, in Sections~\ref{Primary}~through~\ref{Utility}.
2035The only exception is the {\tt control} argument of {\tt umf4sol},
2036which sets {\tt control (8)} to zero to disable iterative refinement.
2037Note that the solve routines do not overwrite {\tt b} with the solution,
2038but return their solution in a different array, {\tt x}.
2039
2040{\footnotesize
2041\begin{verbatim}
2042    call umf4def (control)
2043    call umf4sym (m, n, Ap, Ai, Ax, symbolic, control, info)
2044    call umf4num (Ap, Ai, Ax, symbolic, numeric, control, info)
2045    call umf4solr (sys, Ap, Ai, Ax, x, b, numeric, control, info)
2046    call umf4sol (sys, x, b, numeric, control, info)
2047    call umf4scal (x, b, numeric, status)
2048    call umf4fnum (numeric)
2049    call umf4fsym (symbolic)
2050    call umf4pcon (control)
2051    call umf4pinf (control)
2052    call umf4snum (numeric, filenum, status)
2053    call umf4ssym (symbolic, filenum, status)
2054    call umf4lnum (numeric, filenum, status)
2055    call umf4lsym (symbolic, filenum, status)
2056\end{verbatim}
2057}
2058
2059Access to the complex routines in UMFPACK is provided by the interface
2060routines in {\tt umf4\_f77zwrapper.c}.  The following is a synopsis
2061of each routine.  All the arguments are the same as the real versions,
2062except {\tt Az}, {\tt xz}, and {\tt bz} are the imaginary parts of
2063the matrix, solution, and right-hand side, respectively.  The
2064{\tt Ax}, {\tt x}, and {\tt b} are the real parts.
2065
2066{\footnotesize
2067\begin{verbatim}
2068    call umf4zdef (control)
2069    call umf4zsym (m, n, Ap, Ai, Ax, Az, symbolic, control, info)
2070    call umf4znum (Ap, Ai, Ax, Az, symbolic, numeric, control, info)
2071    call umf4zsolr (sys, Ap, Ai, Ax, Az, x, xz, b, bz, numeric, control, info)
2072    call umf4zsol (sys, x, xz, b, bz, numeric, control, info)
2073    call umf4zscal (x, xz, b, bz, numeric, status)
2074    call umf4zfnum (numeric)
2075    call umf4zfsym (symbolic)
2076    call umf4zpcon (control)
2077    call umf4zpinf (control)
2078    call umf4zsnum (numeric, filenum, status)
2079    call umf4zssym (symbolic, filenum, status)
2080    call umf4zlnum (numeric, filenum, status)
2081    call umf4zlsym (symbolic, filenum, status)
2082\end{verbatim}
2083}
2084
2085The Fortran interface does not support the packed complex case.
2086
2087%-------------------------------------------------------------------------------
2088\section{Installation}
2089\label{Install}
2090%-------------------------------------------------------------------------------
2091
2092%-------------------------------------------------------------------------------
2093\subsection{Installing the C library}
2094%-------------------------------------------------------------------------------
2095
2096The following discussion assumes you have the {\tt make} program, either in
2097Unix, or in Windows with Cygwin\footnote{www.cygwin.com}.
2098You can skip this section and go to next one if all you want to use is
2099the UMFPACK and AMD mexFunctions in MATLAB.
2100
2101You will need to install both UMFPACK and AMD to use UMFPACK.
2102The {\tt UMFPACK} and {\tt AMD} subdirectories must be placed side-by-side
2103within the same parent directory.  AMD is a stand-alone package that
2104is required by UMFPACK.  UMFPACK can be compiled without the
2105BLAS \cite{DaydeDuff99,ACM679a,ATLAS,GotoVandeGeijn02},
2106but your performance will be much less than what it should be.
2107
2108UMFPACK also requires CHOLMOD, CCAMD, CCOLAMD, COLAMD, and metis-5.1.0
2109by default.  You can remove this dependency by compiling with
2110{\tt -DNCHOLMOD}.  Add this to the {\tt UMFPACK\_CONFIG} definition
2111in {\tt SuiteSparse\_config/SuiteSparse\_config.mk}.
2112
2113System-dependent configurations are in the {\tt SuiteSparse\_config/SuiteSparse\_config.mk}
2114file.  The default
2115settings will work on most systems, except that UMFPACK will be compiled so
2116that it does not use the BLAS.  Sample configurations are provided
2117for Linux, Mac, Sun Solaris, SGI IRIX, IBM AIX, and the DEC/Compaq Alpha.
2118The Makefile file will automatically detect what system you have
2119and will set the compile parameters accordingly.
2120
2121To compile both packages,
2122go to the {\tt UMFPACK} directory and type {\tt make}.  This will compile the
2123static
2124libraries ({\tt AMD/Lib/libamd.a} and {\tt UMFPACK/Lib/libumfpack.a}),
2125and the shared libraries ({\tt *.so} on Linux/Unix or {\tt *.dylib} on the Mac).
2126A demo of the AMD ordering routine will be compiled and tested in
2127the {\tt AMD/Demo} directory, and five demo programs will then be
2128compiled and tested in the {\tt UMFPACK/Demo} directory.
2129The outputs of these demo programs will then be compared with output
2130files in the distribution.  Expect to see a few differences, such as
2131residual norms, compile-time control settings, and perhaps memory usage
2132differences.
2133
2134To install into /usr/local/lib and /usr/local/include,
2135do {\tt make install}.
2136To uninstall from there,
2137do {\tt make uninstall}.
2138For more options, see the {\tt SuiteSparse/README.txt} file.
2139
2140Use the MATLAB command {\tt umfpack\_make} in the MATLAB directory
2141to compile UMFPACK and AMD for use in MATLAB.
2142
2143If you compile UMFPACK and AMD and then later change the
2144{\tt SuiteSparse\_config/SuiteSparse\_config.mk} file
2145then you should type {\tt make purge} and then {\tt make} to recompile.
2146
2147Here are a few of the various parameters that you can control in your
2148{\tt SuiteSparse\_config/SuiteSparse\_config.mk} file.
2149To list them all, do {\tt make config}.
2150
2151\begin{itemize}
2152\item {\tt CC = } your C compiler, such as {\tt cc}.
2153\item {\tt RANLIB = } your system's {\tt ranlib} program, if needed.
2154\item {\tt CFLAGS = } optimization flags, such as {\tt -O}.
2155\item {\tt UMFPACK\_CONFIG = } configuration settings for the BLAS,
2156    memory allocation routines, and timing routines.
2157\item {\tt LIB = } your libraries, such as {\tt -lm} or {\tt -lblas}.
2158\item {\tt RM =} the command to delete a file.
2159\item {\tt MV =} the command to rename a file.
2160\item {\tt F77 =} the command to compile a Fortran program (optional).
2161\item {\tt F77FLAGS =} the Fortran compiler flags (optional).
2162\item {\tt F77LIB =} the Fortran libraries (optional).
2163\end{itemize}
2164
2165The {\tt UMFPACK\_CONFIG} string can include combinations of the following;
2166most deal with how the BLAS are called:
2167\begin{itemize}
2168\item {\tt -DNBLAS} if you do not have any BLAS at all.
2169\item {\tt -DNSUNPERF} if you are on Solaris but do not have the Sun
2170    Performance Library (for the BLAS).
2171\item {\tt -DLONGBLAS} if your BLAS takes non-{\tt int} integer arguments.
2172\item {\tt -DBLAS\_INT = } the integer used by the BLAS.
2173
2174\item {\tt -DBLAS\_NO\_UNDERSCORE}
2175    for controlling how C calls the Fortran BLAS.
2176    This is set automatically for Windows,
2177    Sun Solaris, SGI Irix, Red Hat Linux, Compaq Alpha, and
2178    AIX (the IBM RS 6000).
2179
2180\item {\tt -DNRECIPROCAL} controls a trade-off between speed and accuracy.
2181    If defined (or if the pivot value itself is less than $10^{-12}$),
2182    then the pivot column is divided by the pivot value during numeric
2183    factorization.  Otherwise, it is multiplied by the reciprocal of the
2184    pivot, which is faster but can be less accurate.  The default is
2185    to multiply by the reciprocal unless the pivot value is small.
2186    This option also modifies how the rows of the matrix $\m{A}$ are
2187    scaled.  If {\tt -DNRECIPROCAL} is defined (or if any scale factor is
2188    less than $10^{-12}$), entries in the rows of $\m{A}$ are divided
2189    by the scale factors.  Otherwise, they are multiplied by the reciprocal.
2190    When compiling the complex routines with the {\tt gcc} compiler, the
2191    pivot column is always divided by the pivot entry, because of a
2192    numerical accuracy issue encountered with {\tt gcc} version 3.2 with a
2193    few complex matrices on a Pentium 4M (running Linux).  You can still
2194    use {\tt -DNRECIPROCAL} to control how the scale factors
2195    for the rows of $\m{A}$ are applied.
2196\item {\tt -DNO\_DIVIDE\_BY\_ZERO} controls how UMFPACK treats zeros
2197    on the diagonal of $\m{U}$, for a singular matrix $\m{A}$.
2198    If defined, then no division by
2199    zero is performed (a zero entry on the diagonal of $\m{U}$ is
2200    treated as if it were equal to one).  By default,
2201    UMFPACK will divide by zero.
2202
2203\end{itemize}
2204
2205If a Fortran BLAS package is used you may see compiler warnings.  The BLAS
2206routines
2207{\tt dgemm}, {\tt dgemv}, {\tt dger}, {\tt dtrsm}, {\tt dtrsv}, {\tt dscal}
2208and their corresponding complex versions are used.
2209Header files are not provided for the Fortran
2210BLAS.  You may safely ignore all of these warnings.
2211
2212When you compile your program that uses the C-callable UMFPACK library,
2213you need to link your program with all libraries:
2214-lumfpack -lamd -lcholmod -lcolamd -lccolamd -lcamd -lmetis.
2215If you don't compile UMFPACK to use METIS, then you can
2216All libraries are placed in {\tt SuiteSparse/lib} and all include
2217files are placed in {\tt SuiteSparse/include}.
2218
2219You do not need to directly include any AMD include files in your
2220program, unless you directly call AMD routines.  You only need the
2221\begin{verbatim}
2222#include "umfpack.h"
2223\end{verbatim}
2224statement, as described in Section~\ref{Synopsis}.
2225
2226If you do {\tt make install}, then the compiler will know where to
2227find all of the SuiteSparse libraries.  Otherwise, add
2228{\tt -L/home/myself/SuiteSparse/lib} and
2229{\tt -I/home/myself/SuiteSparse/include} to your compiler flags,
2230if your copy of SuiteSparse is located in {\tt /home/myself/SuiteSparse}.
2231
2232Type {\tt make hb} in the {\tt UMFPACK/Demo/HB} directory
2233to compile and run a C program that reads in and factorizes
2234Harwell/Boeing matrices.  Note that this uses a stand-alone Fortran
2235program to read in the Fortran-formatted Harwell/Boeing matrices and
2236write them to a file which can be read by a C program.
2237
2238The {\tt umf\_multicompile.c} file has been added to assist in the
2239compilation of UMFPACK in Microsoft Visual Studio, for Windows.
2240
2241%-------------------------------------------------------------------------------
2242\subsection{Installing the MATLAB interface}
2243%-------------------------------------------------------------------------------
2244
2245Simply type
2246{\tt umfpack\_make} in MATLAB while in the {\tt UMFPACK/MATLAB} directory.
2247You can also type {\tt amd\_make} in the {\tt AMD/MATLAB} directory
2248to compile the stand-alone AMD mexFunction (this is not required to
2249compile the UMFPACK mexFunction).
2250
2251If you are using Windows and the {\tt lcc} compiler bundled with
2252MATLAB 6.1, then you may need to copy the
2253{\tt UMFPACK}$\backslash${\tt MATLAB}$\backslash${\tt lcc\_lib}$\backslash${\tt libmwlapack.lib}
2254file into the
2255{\tt <matlab>}$\backslash${\tt extern}$\backslash${\tt lib}$\backslash${\tt win32}$\backslash${\tt lcc}$\backslash$
2256directory.
2257Next, type {\tt mex -setup}
2258at the MATLAB prompt, and ask MATLAB to select the {\tt lcc} compiler.
2259MATLAB 6.1 has built-in BLAS, but in that version of MATLAB the BLAS
2260cannot be accessed by a mexFunction compiled by {\tt lcc} without first copying
2261this file to the location listed above.
2262If you have MATLAB 6.5 or later, you can probably skip this step.
2263
2264%-------------------------------------------------------------------------------
2265\subsection{Installing the Fortran interface}
2266%-------------------------------------------------------------------------------
2267
2268Once the 32-bit C-callable UMFPACK library is compiled, you can also compile
2269the Fortran interface, by typing {\tt make fortran}.  This will create
2270the {\tt umf4hb} program, test it, and compare the output with the
2271file {\tt umf4hb.out} in the distribution.
2272If you compiled UMFPACK in 64-bit mode, you need to use {\tt make fortran64}
2273instead, which compiles the {\tt umf4hb64} program and compares its output
2274with the file {\tt umf4hb64.out}.
2275Refer to the comments in the {\tt Demo/umf4\_f77wrapper.c} file
2276for more details.
2277
2278This interface is {\bf highly} non-portable, since it depends
2279on how C and Fortran are interfaced.
2280Because of this issue, the interface is included in the {\tt Demo} directory,
2281and not as a primary part of the UMFPACK library.  The interface routines are
2282not included in the compiled {\tt UMFPACK/Lib/libumfpack.a} library, but left
2283as stand-alone compiled files ({\tt umf4\_f77wrapper.o} and
2284{\tt umf4\_f77wrapper64.o} in the {\tt Demo} directory).
2285You may need to modify the interface routines in the file
2286{\tt umf4\_f77wrapper.c} if you are using compilers for which this interface
2287has not been tested.
2288
2289In particular, I was not able to get C and Fortran to work together on
2290the Mac (Snow Leopard).
2291
2292%-------------------------------------------------------------------------------
2293\subsection{Known Issues}
2294%-------------------------------------------------------------------------------
2295
2296The Microsoft C or C++ compilers on a Pentium badly break the IEEE 754 standard,
2297and do not treat NaN's properly.  According to IEEE 754, the expression
2298{\tt (x != x)} is supposed to be true if and only if {\tt x} is NaN.  For
2299non-compliant compilers in Windows that expression is always false, and another
2300test must be used: {\tt (x < x)} is true if and only if {\tt x}
2301is NaN.  For compliant compilers, {\tt (x < x)} is always false, for any
2302value of {\tt x} (including NaN).
2303To cover both cases, UMFPACK when running under Microsoft Windows
2304defines the following macro, which is true if and only if {\tt x} is NaN,
2305regardless of whether your compiler is compliant or not:
2306
2307\begin{verbatim}
2308#define SCALAR_IS_NAN(x) (((x) != (x)) || ((x) < (x)))
2309\end{verbatim}
2310
2311If your compiler breaks this test, then UMFPACK will fail catastrophically
2312if it encounters a NaN.  You will not just see NaN's in your output; UMFPACK
2313will probably crash with a segmentation fault.  In that case, you might try to
2314see if the common (but non-ANSI C) routine {\tt isnan} is available, and modify
2315the macro {\tt SCALAR\_IS\_NAN} in {\tt umf\_version.h} accordingly.  The
2316simpler (and IEEE 754-compliant) test {\tt (x != x)} is always true with Linux
2317on a PC, and on every Unix compiler I have tested.
2318
2319Some compilers will complain about the Fortran BLAS being defined implicitly.
2320C prototypes for the BLAS are not used, except the C-BLAS.  Some compilers
2321will complain about unrecognized {\tt \#pragma}'s.  You may safely ignore
2322all of these warnings.
2323
2324%-------------------------------------------------------------------------------
2325\section{Future work}
2326\label{Future}
2327%-------------------------------------------------------------------------------
2328
2329Here are a few features that are not in the current version of UMFPACK,
2330in no particular
2331order.  They may appear in a future release of UMFPACK.  If you are interested,
2332let me know and I could consider including them:
2333
2334\begin{enumerate}
2335
2336\item Remove the restriction that the column-oriented form be given with
2337    sorted columns.  This has already been done in AMD Version 2.0.
2338
2339\item Future versions may have different default {\tt Control} parameters.
2340    Future versions may return more statistics in the {\tt Info} array, and
2341    they may use more entries in the {\tt Control} array.
2342    These two arrays will probably become larger, since there are very few
2343    unused entries.  If they change in size, the constants
2344    {\tt UMFPACK\_CONTROL} and {\tt UMFPACK\_INFO} defined in {\tt umfpack.h}
2345    will be changed to reflect their new size.  Your C program should use
2346    these constants when declaring the size of these two arrays.  Do not
2347    define them as {\tt Control [20]} and {\tt Info [90]}.
2348
2349\item Forward/back solvers for the conventional row or column-form data
2350    structure for $\m{L}$ and $\m{U}$ (the output of
2351    {\tt umfpack\_*\_di\_get\_numeric}).  This would enable a separate
2352    solver that could be used to write a MATLAB mexFunction
2353    {\tt x = lu\_refine (A, b, L, U, P, Q, R)} that gives MATLAB access
2354    to the iterative refinement algorithm with sparse backward error
2355    analysis.  It would also be easier to handle sparse right-hand sides
2356    in this data structure, and end up with good asymptotic run-time
2357    in this case
2358    (particularly for $\m{Lx}=\m{b}$; see \cite{GilbertPeierls88}).
2359    See also CSparse and
2360    CXSparse for software for handling sparse right-hand sides.
2361
2362\item Complex absolute value computations could be
2363    based on FDLIBM (see \newline
2364    http://www.netlib.org/fdlibm),
2365    using the {\tt hypot(x,y)} routine.
2366
2367\item When using iterative refinement, the residual $\m{Ax}-\m{b}$ could be
2368    returned by {\tt umfpack\_solve}.
2369
2370\item The solve routines could handle multiple right-hand sides, and sparse
2371    right-hand sides.  See {\tt umfpack\_solve} for the MATLAB version
2372    of this feature.
2373    See also CSparse and
2374    CXSparse for software for handling sparse right-hand sides.
2375
2376\item An option to redirect the error and diagnostic output.
2377
2378\item Permutation to block-triangular-form \cite{Duff78a} for the C-callable
2379    interface.  There are two routines in the ACM Collected
2380    Algorithms (529 and 575) \cite{Duff81b,Duff78b}
2381    that could be translated from Fortran
2382    to C and included in UMFPACK.  This would result in better performance
2383    for matrices from circuit simulation and
2384    chemical process engineering.  See {\tt umfpack\_btf.m} for the MATLAB
2385    version of this feature.  KLU includes this feature.
2386    See also {\tt cs\_dmperm} in CSparse and CXSparse.
2387
2388\item The ability to use user-provided work arrays, so that {\tt malloc},
2389    {\tt free}, and {\tt realloc} realloc are not called.  The
2390    {\tt umfpack\_*\_wsolve} routine is one example.
2391
2392\item A method that takes time proportional to the number of nonzeros in
2393    $\m{A}$ to compute the symbolic factorization \cite{GilbertNgPeyton94}.
2394    This would improve the performance of the symmetric strategy,
2395    and the unsymmetric strategy when dense rows are present.
2396    The current method takes
2397    time proportional to the number of nonzeros in the upper bound of $\m{U}$.
2398    The method used in UMFPACK exploits super-columns, however, so this
2399    bound is rarely reached.
2400    See {\tt cs\_counts} in CSparse and CXSparse,
2401    and {\tt cholmod\_analyze} in CHOLMOD.
2402
2403\item Other basic sparse matrix operations, such as sparse matrix
2404    multiplication, could be included.
2405
2406\item A more complete Fortran interface.
2407
2408\item A C++ interface.
2409
2410\item A parallel version using MPI.  This would require a large amount
2411    of effort.
2412
2413\end{enumerate}
2414
2415
2416%-------------------------------------------------------------------------------
2417\newpage
2418\section{The primary UMFPACK routines}
2419\label{Primary}
2420%-------------------------------------------------------------------------------
2421
2422The include files are the same for all four versions of
2423UMFPACK.  The generic integer type is {\tt Int}, which is an {\tt int} or
2424{\tt SuiteSparse\_long}, depending on which version of UMFPACK you are using.
2425
2426\subsection{umfpack\_*\_symbolic}
2427
2428{\footnotesize
2429\begin{verbatim}
2430INCLUDE umfpack_symbolic.h via sed
2431\end{verbatim}
2432}
2433
2434\newpage
2435\subsection{umfpack\_*\_numeric}
2436
2437{\footnotesize
2438\begin{verbatim}
2439INCLUDE umfpack_numeric.h via sed
2440\end{verbatim}
2441}
2442
2443\newpage
2444\subsection{umfpack\_*\_solve}
2445
2446{\footnotesize
2447\begin{verbatim}
2448INCLUDE umfpack_solve.h via sed
2449\end{verbatim}
2450}
2451
2452\newpage
2453\subsection{umfpack\_*\_free\_symbolic}
2454
2455{\footnotesize
2456\begin{verbatim}
2457INCLUDE umfpack_free_symbolic.h via sed
2458\end{verbatim}
2459}
2460
2461\newpage
2462\subsection{umfpack\_*\_free\_numeric}
2463
2464{\footnotesize
2465\begin{verbatim}
2466INCLUDE umfpack_free_numeric.h via sed
2467\end{verbatim}
2468}
2469
2470%-------------------------------------------------------------------------------
2471\newpage
2472\section{Alternative routines}
2473\label{Alternative}
2474%-------------------------------------------------------------------------------
2475
2476\subsection{umfpack\_*\_defaults}
2477
2478{\footnotesize
2479\begin{verbatim}
2480INCLUDE umfpack_defaults.h via sed
2481\end{verbatim}
2482}
2483
2484\newpage
2485\subsection{umfpack\_*\_qsymbolic and umfpack\_*\_fsymbolic}
2486
2487{\footnotesize
2488\begin{verbatim}
2489INCLUDE umfpack_qsymbolic.h via sed
2490\end{verbatim}
2491}
2492
2493\newpage
2494\subsection{umfpack\_*\_wsolve}
2495
2496{\footnotesize
2497\begin{verbatim}
2498INCLUDE umfpack_wsolve.h via sed
2499\end{verbatim}
2500}
2501
2502%-------------------------------------------------------------------------------
2503\newpage
2504\section{Matrix manipulation routines}
2505\label{Manipulate}
2506%-------------------------------------------------------------------------------
2507
2508\subsection{umfpack\_*\_col\_to\_triplet}
2509
2510{\footnotesize
2511\begin{verbatim}
2512INCLUDE umfpack_col_to_triplet.h via sed
2513\end{verbatim}
2514}
2515
2516\newpage
2517\subsection{umfpack\_*\_triplet\_to\_col}
2518
2519{\footnotesize
2520\begin{verbatim}
2521INCLUDE umfpack_triplet_to_col.h via sed
2522\end{verbatim}
2523}
2524
2525\newpage
2526\subsection{umfpack\_*\_transpose}
2527
2528{\footnotesize
2529\begin{verbatim}
2530INCLUDE umfpack_transpose.h via sed
2531\end{verbatim}
2532}
2533
2534\newpage
2535\subsection{umfpack\_*\_scale}
2536
2537{\footnotesize
2538\begin{verbatim}
2539INCLUDE umfpack_scale.h via sed
2540\end{verbatim}
2541}
2542
2543%-------------------------------------------------------------------------------
2544\newpage
2545\section{Getting the contents of opaque objects}
2546\label{Get}
2547%-------------------------------------------------------------------------------
2548
2549\subsection{umfpack\_*\_get\_lunz}
2550
2551{\footnotesize
2552\begin{verbatim}
2553INCLUDE umfpack_get_lunz.h via sed
2554\end{verbatim}
2555}
2556
2557\newpage
2558\subsection{umfpack\_*\_get\_numeric}
2559
2560{\footnotesize
2561\begin{verbatim}
2562INCLUDE umfpack_get_numeric.h via sed
2563\end{verbatim}
2564}
2565
2566\newpage
2567\subsection{umfpack\_*\_get\_symbolic}
2568
2569{\footnotesize
2570\begin{verbatim}
2571INCLUDE umfpack_get_symbolic.h via sed
2572\end{verbatim}
2573}
2574
2575\newpage
2576\subsection{umfpack\_*\_save\_numeric}
2577
2578{\footnotesize
2579\begin{verbatim}
2580INCLUDE umfpack_save_numeric.h via sed
2581\end{verbatim}
2582}
2583
2584\newpage
2585\subsection{umfpack\_*\_load\_numeric}
2586
2587{\footnotesize
2588\begin{verbatim}
2589INCLUDE umfpack_load_numeric.h via sed
2590\end{verbatim}
2591}
2592
2593\newpage
2594\subsection{umfpack\_*\_save\_symbolic}
2595
2596{\footnotesize
2597\begin{verbatim}
2598INCLUDE umfpack_save_symbolic.h via sed
2599\end{verbatim}
2600}
2601
2602\newpage
2603\subsection{umfpack\_*\_load\_symbolic}
2604
2605{\footnotesize
2606\begin{verbatim}
2607INCLUDE umfpack_load_symbolic.h via sed
2608\end{verbatim}
2609}
2610
2611\newpage
2612\subsection{umfpack\_*\_get\_determinant}
2613
2614{\footnotesize
2615\begin{verbatim}
2616INCLUDE umfpack_get_determinant.h via sed
2617\end{verbatim}
2618}
2619
2620%-------------------------------------------------------------------------------
2621\newpage
2622\section{Reporting routines}
2623\label{Report}
2624%-------------------------------------------------------------------------------
2625
2626\subsection{umfpack\_*\_report\_status}
2627
2628{\footnotesize
2629\begin{verbatim}
2630INCLUDE umfpack_report_status.h via sed
2631\end{verbatim}
2632}
2633
2634\newpage
2635\subsection{umfpack\_*\_report\_control}
2636
2637{\footnotesize
2638\begin{verbatim}
2639INCLUDE umfpack_report_control.h via sed
2640\end{verbatim}
2641}
2642
2643\newpage
2644\subsection{umfpack\_*\_report\_info}
2645
2646{\footnotesize
2647\begin{verbatim}
2648INCLUDE umfpack_report_info.h via sed
2649\end{verbatim}
2650}
2651
2652\newpage
2653\subsection{umfpack\_*\_report\_matrix}
2654
2655{\footnotesize
2656\begin{verbatim}
2657INCLUDE umfpack_report_matrix.h via sed
2658\end{verbatim}
2659}
2660
2661\newpage
2662\subsection{umfpack\_*\_report\_numeric}
2663
2664{\footnotesize
2665\begin{verbatim}
2666INCLUDE umfpack_report_numeric.h via sed
2667\end{verbatim}
2668}
2669
2670\newpage
2671\subsection{umfpack\_*\_report\_perm}
2672
2673{\footnotesize
2674\begin{verbatim}
2675INCLUDE umfpack_report_perm.h via sed
2676\end{verbatim}
2677}
2678
2679\newpage
2680\subsection{umfpack\_*\_report\_symbolic}
2681
2682{\footnotesize
2683\begin{verbatim}
2684INCLUDE umfpack_report_symbolic.h via sed
2685\end{verbatim}
2686}
2687
2688\newpage
2689\subsection{umfpack\_*\_report\_triplet}
2690
2691{\footnotesize
2692\begin{verbatim}
2693INCLUDE umfpack_report_triplet.h via sed
2694\end{verbatim}
2695}
2696
2697\newpage
2698\subsection{umfpack\_*\_report\_vector}
2699
2700{\footnotesize
2701\begin{verbatim}
2702INCLUDE umfpack_report_vector.h via sed
2703\end{verbatim}
2704}
2705
2706%-------------------------------------------------------------------------------
2707\newpage
2708\section{Utility routines}
2709\label{Utility}
2710%-------------------------------------------------------------------------------
2711
2712\subsection{umfpack\_timer}
2713
2714{\footnotesize
2715\begin{verbatim}
2716INCLUDE umfpack_timer.h via sed
2717\end{verbatim}
2718}
2719
2720\newpage
2721\subsection{umfpack\_tic and umfpack\_toc}
2722
2723{\footnotesize
2724\begin{verbatim}
2725INCLUDE umfpack_tictoc.h via sed
2726\end{verbatim}
2727}
2728
2729
2730%-------------------------------------------------------------------------------
2731\newpage
2732% References
2733%-------------------------------------------------------------------------------
2734
2735\bibliographystyle{plain}
2736\bibliography{UserGuide}
2737
2738\end{document}
2739