sundials_dense.h

Go to the documentation of this file.
00001 /*
00002  * -----------------------------------------------------------------
00003  * $Revision: 1.6 $
00004  * $Date: 2006/11/29 00:05:07 $
00005  * -----------------------------------------------------------------
00006  * Programmer: Radu Serban @ LLNL
00007  * -----------------------------------------------------------------
00008  * Copyright (c) 2006, The Regents of the University of California.
00009  * Produced at the Lawrence Livermore National Laboratory.
00010  * All rights reserved.
00011  * For details, see the LICENSE file.
00012  * -----------------------------------------------------------------
00013  * This is the header file for a generic package of DENSE matrix
00014  * operations, based on the DlsMat type defined in sundials_direct.h.
00015  *
00016  * There are two sets of dense solver routines listed in
00017  * this file: one set uses type DlsMat defined below and the
00018  * other set uses the type realtype ** for dense matrix arguments.
00019  * Routines that work with the type DlsMat begin with "Dense".
00020  * Routines that work with realtype** begin with "dense". 
00021  * -----------------------------------------------------------------
00022  */
00023 
00024 #ifndef _SUNDIALS_DENSE_H
00025 #define _SUNDIALS_DENSE_H
00026 
00027 #ifdef __cplusplus  /* wrapper to enable C++ usage */
00028 extern "C" {
00029 #endif
00030 
00031 #include <sundials/sundials_direct.h>
00032 
00033 /*
00034  * -----------------------------------------------------------------
00035  * Functions: DenseGETRF and DenseGETRS
00036  * -----------------------------------------------------------------
00037  * DenseGETRF performs the LU factorization of the M by N dense
00038  * matrix A. This is done using standard Gaussian elimination
00039  * with partial (row) pivoting. Note that this applies only
00040  * to matrices with M >= N and full column rank.
00041  *
00042  * A successful LU factorization leaves the matrix A and the
00043  * pivot array p with the following information:
00044  *
00045  * (1) p[k] contains the row number of the pivot element chosen
00046  *     at the beginning of elimination step k, k=0, 1, ..., N-1.
00047  *
00048  * (2) If the unique LU factorization of A is given by PA = LU,
00049  *     where P is a permutation matrix, L is a lower trapezoidal
00050  *     matrix with all 1's on the diagonal, and U is an upper
00051  *     triangular matrix, then the upper triangular part of A
00052  *     (including its diagonal) contains U and the strictly lower
00053  *     trapezoidal part of A contains the multipliers, I-L.
00054  *
00055  * For square matrices (M=N), L is unit lower triangular.
00056  *
00057  * DenseGETRF returns 0 if successful. Otherwise it encountered
00058  * a zero diagonal element during the factorization. In this case
00059  * it returns the column index (numbered from one) at which
00060  * it encountered the zero.
00061  *
00062  * DenseGETRS solves the N-dimensional system A x = b using
00063  * the LU factorization in A and the pivot information in p
00064  * computed in DenseGETRF. The solution x is returned in b. This
00065  * routine cannot fail if the corresponding call to DenseGETRF
00066  * did not fail.
00067  * DenseGETRS does NOT check for a square matrix!
00068  *
00069  * -----------------------------------------------------------------
00070  * DenseGETRF and DenseGETRS are simply wrappers around denseGETRF
00071  * and denseGETRS, respectively, which perform all the work by
00072  * directly accessing the data in the DlsMat A (i.e. the field cols)
00073  * -----------------------------------------------------------------
00074  */
00075 
00076 SUNDIALS_EXPORT int DenseGETRF(DlsMat A, int *p);
00077 SUNDIALS_EXPORT void DenseGETRS(DlsMat A, int *p, realtype *b);
00078 
00079 SUNDIALS_EXPORT int denseGETRF(realtype **a, int m, int n, int *p);
00080 SUNDIALS_EXPORT void denseGETRS(realtype **a, int n, int *p, realtype *b);
00081 
00082 /*
00083  * -----------------------------------------------------------------
00084  * Functions : DensePOTRF and DensePOTRS
00085  * -----------------------------------------------------------------
00086  * DensePOTRF computes the Cholesky factorization of a real symmetric
00087  * positive definite matrix A.
00088  * -----------------------------------------------------------------
00089  * DensePOTRS solves a system of linear equations A*X = B with a 
00090  * symmetric positive definite matrix A using the Cholesky factorization
00091  * A = L*L**T computed by DensePOTRF.
00092  *
00093  * -----------------------------------------------------------------
00094  * DensePOTRF and DensePOTRS are simply wrappers around densePOTRF
00095  * and densePOTRS, respectively, which perform all the work by
00096  * directly accessing the data in the DlsMat A (i.e. the field cols)
00097  * -----------------------------------------------------------------
00098  */
00099 
00100 SUNDIALS_EXPORT int DensePOTRF(DlsMat A);
00101 SUNDIALS_EXPORT void DensePOTRS(DlsMat A, realtype *b);
00102 
00103 SUNDIALS_EXPORT int densePOTRF(realtype **a, int m);
00104 SUNDIALS_EXPORT void densePOTRS(realtype **a, int m, realtype *b);
00105 
00106 /*
00107  * -----------------------------------------------------------------
00108  * Functions : DenseGEQRF and DenseORMQR
00109  * -----------------------------------------------------------------
00110  * DenseGEQRF computes a QR factorization of a real M-by-N matrix A:
00111  * A = Q * R (with M>= N).
00112  * 
00113  * DenseGEQRF requires a temporary work vector wrk of length M.
00114  * -----------------------------------------------------------------
00115  * DenseORMQR computes the product w = Q * v where Q is a real 
00116  * orthogonal matrix defined as the product of k elementary reflectors
00117  *
00118  *        Q = H(1) H(2) . . . H(k)
00119  *
00120  * as returned by DenseGEQRF. Q is an M-by-N matrix, v is a vector
00121  * of length N and w is a vector of length M (with M>=N).
00122  *
00123  * DenseORMQR requires a temporary work vector wrk of length M.
00124  *
00125  * -----------------------------------------------------------------
00126  * DenseGEQRF and DenseORMQR are simply wrappers around denseGEQRF
00127  * and denseORMQR, respectively, which perform all the work by
00128  * directly accessing the data in the DlsMat A (i.e. the field cols)
00129  * -----------------------------------------------------------------
00130  */
00131 
00132 SUNDIALS_EXPORT int DenseGEQRF(DlsMat A, realtype *beta, realtype *wrk);
00133 SUNDIALS_EXPORT int DenseORMQR(DlsMat A, realtype *beta, realtype *vn, realtype *vm, 
00134                    realtype *wrk);
00135 
00136 SUNDIALS_EXPORT int denseGEQRF(realtype **a, int m, int n, realtype *beta, realtype *v);
00137 SUNDIALS_EXPORT int denseORMQR(realtype **a, int m, int n, realtype *beta,
00138                    realtype *v, realtype *w, realtype *wrk);
00139 
00140 /*
00141  * -----------------------------------------------------------------
00142  * Function : DenseZero
00143  * -----------------------------------------------------------------
00144  * DenseZero sets all the elements of the M-by-N matrix A to 0.0.
00145  *
00146  * DenseZero is a wrapper around denseZero which accesses the data of
00147  * the DlsMat A (i.e. the field cols)
00148  * -----------------------------------------------------------------
00149  */
00150 
00151 SUNDIALS_EXPORT void DenseZero(DlsMat A);
00152 SUNDIALS_EXPORT void denseZero(realtype **a, int m, int n);
00153 
00154 /*
00155  * -----------------------------------------------------------------
00156  * Function : DenseCopy
00157  * -----------------------------------------------------------------
00158  * DenseCopy copies the contents of the M-by-N matrix A into the
00159  * M-by-N matrix B.
00160  * 
00161  * DenseCopy is a wrapper around denseCopy which accesses the data
00162  * in the DlsMat A and B (i.e. the fields cols)
00163  * -----------------------------------------------------------------
00164  */
00165 
00166 SUNDIALS_EXPORT void DenseCopy(DlsMat A, DlsMat B);
00167 SUNDIALS_EXPORT void denseCopy(realtype **a, realtype **b, int m, int n);
00168 
00169 /*
00170  * -----------------------------------------------------------------
00171  * Function: DenseScale
00172  * -----------------------------------------------------------------
00173  * DenseScale scales the elements of the M-by-N matrix A by the
00174  * constant c and stores the result back in A.
00175  *
00176  * DenseScale is a wrapper around denseScale which performs the actual
00177  * scaling by accessing the data in the DlsMat A (i.e. the field
00178  * cols).
00179  * -----------------------------------------------------------------
00180  */
00181 
00182 SUNDIALS_EXPORT void DenseScale(realtype c, DlsMat A);
00183 SUNDIALS_EXPORT void denseScale(realtype c, realtype **a, int m, int n);
00184 
00185 /*
00186  * -----------------------------------------------------------------
00187  * Function : DenseAddI
00188  * -----------------------------------------------------------------
00189  * DenseAddI adds 1.0 to the main diagonal (A_ii, i=1,2,...,N-1) of
00190  * the M-by-N matrix A (M>= N) and stores the result back in A.
00191  * DenseAddI is typically used with square matrices.
00192  * DenseAddI does not check for M >= N and therefore a segmentation
00193  * fault will occur if M < N!
00194  *
00195  * DenseAddI is a wrapper around denseAddI which performs the actual
00196  * work by accessing the data in the DlsMat A (i.e. the field cols)
00197  * -----------------------------------------------------------------
00198  */
00199 
00200 SUNDIALS_EXPORT void DenseAddI(DlsMat A);
00201 SUNDIALS_EXPORT void denseAddI(realtype **a, int n);
00202 
00203 #ifdef __cplusplus
00204 }
00205 #endif
00206 
00207 #endif

Generated on Fri Sep 26 07:44:17 2008 for SimTKcore by  doxygen 1.5.6