1 | /* |
---|
2 | * ----------------------------------------------------------------- |
---|
3 | * $Revision$ |
---|
4 | * $Date$ |
---|
5 | * ----------------------------------------------------------------- |
---|
6 | * Programmer: Radu Serban @ LLNL |
---|
7 | * ----------------------------------------------------------------- |
---|
8 | * LLNS Copyright Start |
---|
9 | * Copyright (c) 2014, Lawrence Livermore National Security |
---|
10 | * This work was performed under the auspices of the U.S. Department |
---|
11 | * of Energy by Lawrence Livermore National Laboratory in part under |
---|
12 | * Contract W-7405-Eng-48 and in part under Contract DE-AC52-07NA27344. |
---|
13 | * Produced at the Lawrence Livermore National Laboratory. |
---|
14 | * All rights reserved. |
---|
15 | * For details, see the LICENSE file. |
---|
16 | * LLNS Copyright End |
---|
17 | * ----------------------------------------------------------------- |
---|
18 | * This is the header file for a generic package of DENSE matrix |
---|
19 | * operations, based on the DlsMat type defined in sundials_direct.h. |
---|
20 | * |
---|
21 | * There are two sets of dense solver routines listed in |
---|
22 | * this file: one set uses type DlsMat defined below and the |
---|
23 | * other set uses the type realtype ** for dense matrix arguments. |
---|
24 | * Routines that work with the type DlsMat begin with "Dense". |
---|
25 | * Routines that work with realtype** begin with "dense". |
---|
26 | * ----------------------------------------------------------------- |
---|
27 | */ |
---|
28 | |
---|
29 | #ifndef _SUNDIALS_DENSE_H |
---|
30 | #define _SUNDIALS_DENSE_H |
---|
31 | |
---|
32 | #include <sundials/sundials_direct.h> |
---|
33 | |
---|
34 | #ifdef __cplusplus /* wrapper to enable C++ usage */ |
---|
35 | extern "C" { |
---|
36 | #endif |
---|
37 | |
---|
38 | /* |
---|
39 | * ----------------------------------------------------------------- |
---|
40 | * Functions: DenseGETRF and DenseGETRS |
---|
41 | * ----------------------------------------------------------------- |
---|
42 | * DenseGETRF performs the LU factorization of the M by N dense |
---|
43 | * matrix A. This is done using standard Gaussian elimination |
---|
44 | * with partial (row) pivoting. Note that this applies only |
---|
45 | * to matrices with M >= N and full column rank. |
---|
46 | * |
---|
47 | * A successful LU factorization leaves the matrix A and the |
---|
48 | * pivot array p with the following information: |
---|
49 | * |
---|
50 | * (1) p[k] contains the row number of the pivot element chosen |
---|
51 | * at the beginning of elimination step k, k=0, 1, ..., N-1. |
---|
52 | * |
---|
53 | * (2) If the unique LU factorization of A is given by PA = LU, |
---|
54 | * where P is a permutation matrix, L is a lower trapezoidal |
---|
55 | * matrix with all 1's on the diagonal, and U is an upper |
---|
56 | * triangular matrix, then the upper triangular part of A |
---|
57 | * (including its diagonal) contains U and the strictly lower |
---|
58 | * trapezoidal part of A contains the multipliers, I-L. |
---|
59 | * |
---|
60 | * For square matrices (M = N), L is unit lower triangular. |
---|
61 | * |
---|
62 | * DenseGETRF returns 0 if successful. Otherwise it encountered |
---|
63 | * a zero diagonal element during the factorization. In this case |
---|
64 | * it returns the column index (numbered from one) at which |
---|
65 | * it encountered the zero. |
---|
66 | * |
---|
67 | * DenseGETRS solves the N-dimensional system A x = b using |
---|
68 | * the LU factorization in A and the pivot information in p |
---|
69 | * computed in DenseGETRF. The solution x is returned in b. This |
---|
70 | * routine cannot fail if the corresponding call to DenseGETRF |
---|
71 | * did not fail. |
---|
72 | * DenseGETRS does NOT check for a square matrix! |
---|
73 | * |
---|
74 | * ----------------------------------------------------------------- |
---|
75 | * DenseGETRF and DenseGETRS are simply wrappers around denseGETRF |
---|
76 | * and denseGETRS, respectively, which perform all the work by |
---|
77 | * directly accessing the data in the DlsMat A (i.e. in A->cols). |
---|
78 | * ----------------------------------------------------------------- |
---|
79 | */ |
---|
80 | |
---|
81 | SUNDIALS_EXPORT sunindextype DenseGETRF(DlsMat A, sunindextype *p); |
---|
82 | SUNDIALS_EXPORT void DenseGETRS(DlsMat A, sunindextype *p, realtype *b); |
---|
83 | |
---|
84 | SUNDIALS_EXPORT sunindextype denseGETRF(realtype **a, sunindextype m, sunindextype n, sunindextype *p); |
---|
85 | SUNDIALS_EXPORT void denseGETRS(realtype **a, sunindextype n, sunindextype *p, realtype *b); |
---|
86 | |
---|
87 | /* |
---|
88 | * ----------------------------------------------------------------- |
---|
89 | * Functions : DensePOTRF and DensePOTRS |
---|
90 | * ----------------------------------------------------------------- |
---|
91 | * DensePOTRF computes the Cholesky factorization of a real symmetric |
---|
92 | * positive definite matrix A. |
---|
93 | * ----------------------------------------------------------------- |
---|
94 | * DensePOTRS solves a system of linear equations A*X = B with a |
---|
95 | * symmetric positive definite matrix A using the Cholesky factorization |
---|
96 | * A = L*L**T computed by DensePOTRF. |
---|
97 | * |
---|
98 | * ----------------------------------------------------------------- |
---|
99 | * DensePOTRF and DensePOTRS are simply wrappers around densePOTRF |
---|
100 | * and densePOTRS, respectively, which perform all the work by |
---|
101 | * directly accessing the data in the DlsMat A (i.e. the field cols) |
---|
102 | * ----------------------------------------------------------------- |
---|
103 | */ |
---|
104 | |
---|
105 | SUNDIALS_EXPORT sunindextype DensePOTRF(DlsMat A); |
---|
106 | SUNDIALS_EXPORT void DensePOTRS(DlsMat A, realtype *b); |
---|
107 | |
---|
108 | SUNDIALS_EXPORT sunindextype densePOTRF(realtype **a, sunindextype m); |
---|
109 | SUNDIALS_EXPORT void densePOTRS(realtype **a, sunindextype m, realtype *b); |
---|
110 | |
---|
111 | /* |
---|
112 | * ----------------------------------------------------------------- |
---|
113 | * Functions : DenseGEQRF and DenseORMQR |
---|
114 | * ----------------------------------------------------------------- |
---|
115 | * DenseGEQRF computes a QR factorization of a real M-by-N matrix A: |
---|
116 | * A = Q * R (with M>= N). |
---|
117 | * |
---|
118 | * DenseGEQRF requires a temporary work vector wrk of length M. |
---|
119 | * ----------------------------------------------------------------- |
---|
120 | * DenseORMQR computes the product w = Q * v where Q is a real |
---|
121 | * orthogonal matrix defined as the product of k elementary reflectors |
---|
122 | * |
---|
123 | * Q = H(1) H(2) . . . H(k) |
---|
124 | * |
---|
125 | * as returned by DenseGEQRF. Q is an M-by-N matrix, v is a vector |
---|
126 | * of length N and w is a vector of length M (with M >= N). |
---|
127 | * |
---|
128 | * DenseORMQR requires a temporary work vector wrk of length M. |
---|
129 | * |
---|
130 | * ----------------------------------------------------------------- |
---|
131 | * DenseGEQRF and DenseORMQR are simply wrappers around denseGEQRF |
---|
132 | * and denseORMQR, respectively, which perform all the work by |
---|
133 | * directly accessing the data in the DlsMat A (i.e. the field cols) |
---|
134 | * ----------------------------------------------------------------- |
---|
135 | */ |
---|
136 | |
---|
137 | SUNDIALS_EXPORT int DenseGEQRF(DlsMat A, realtype *beta, realtype *wrk); |
---|
138 | SUNDIALS_EXPORT int DenseORMQR(DlsMat A, realtype *beta, realtype *vn, realtype *vm, |
---|
139 | realtype *wrk); |
---|
140 | |
---|
141 | SUNDIALS_EXPORT int denseGEQRF(realtype **a, sunindextype m, sunindextype n, realtype *beta, |
---|
142 | realtype *wrk); |
---|
143 | SUNDIALS_EXPORT int denseORMQR(realtype **a, sunindextype m, sunindextype n, realtype *beta, |
---|
144 | realtype *v, realtype *w, realtype *wrk); |
---|
145 | |
---|
146 | /* |
---|
147 | * ----------------------------------------------------------------- |
---|
148 | * Function : DenseCopy |
---|
149 | * ----------------------------------------------------------------- |
---|
150 | * DenseCopy copies the contents of the M-by-N matrix A into the |
---|
151 | * M-by-N matrix B. |
---|
152 | * |
---|
153 | * DenseCopy is a wrapper around denseCopy which accesses the data |
---|
154 | * in the DlsMat A and DlsMat B (i.e. the fields cols) |
---|
155 | * ----------------------------------------------------------------- |
---|
156 | */ |
---|
157 | |
---|
158 | SUNDIALS_EXPORT void DenseCopy(DlsMat A, DlsMat B); |
---|
159 | SUNDIALS_EXPORT void denseCopy(realtype **a, realtype **b, sunindextype m, sunindextype n); |
---|
160 | |
---|
161 | /* |
---|
162 | * ----------------------------------------------------------------- |
---|
163 | * Function: DenseScale |
---|
164 | * ----------------------------------------------------------------- |
---|
165 | * DenseScale scales the elements of the M-by-N matrix A by the |
---|
166 | * constant c and stores the result back in A. |
---|
167 | * |
---|
168 | * DenseScale is a wrapper around denseScale which performs the actual |
---|
169 | * scaling by accessing the data in the DlsMat A (i.e. in A->cols). |
---|
170 | * ----------------------------------------------------------------- |
---|
171 | */ |
---|
172 | |
---|
173 | SUNDIALS_EXPORT void DenseScale(realtype c, DlsMat A); |
---|
174 | SUNDIALS_EXPORT void denseScale(realtype c, realtype **a, sunindextype m, sunindextype n); |
---|
175 | |
---|
176 | |
---|
177 | /* |
---|
178 | * ----------------------------------------------------------------- |
---|
179 | * Function: denseAddIdentity |
---|
180 | * ----------------------------------------------------------------- |
---|
181 | * denseAddIdentity adds the identity matrix to the n-by-n matrix |
---|
182 | * stored in a realtype** array. |
---|
183 | * ----------------------------------------------------------------- |
---|
184 | */ |
---|
185 | |
---|
186 | SUNDIALS_EXPORT void denseAddIdentity(realtype **a, sunindextype n); |
---|
187 | |
---|
188 | |
---|
189 | /* |
---|
190 | * ----------------------------------------------------------------- |
---|
191 | * Function: DenseMatvec |
---|
192 | * ----------------------------------------------------------------- |
---|
193 | * DenseMatvec computes the matrix-vector product y = A*x, where A |
---|
194 | * is an M-by-N matrix, x is a vector of length N, and y is a vector |
---|
195 | * of length M. No error checking is performed on the length of the |
---|
196 | * arrays x and y. Only y is modified in this routine. |
---|
197 | * |
---|
198 | * DenseMatvec is a wrapper around denseMatvec which performs the |
---|
199 | * actual product by accessing the data in the DlsMat A. |
---|
200 | * ----------------------------------------------------------------- |
---|
201 | */ |
---|
202 | |
---|
203 | SUNDIALS_EXPORT void DenseMatvec(DlsMat A, realtype *x, realtype *y); |
---|
204 | SUNDIALS_EXPORT void denseMatvec(realtype **a, realtype *x, realtype *y, sunindextype m, sunindextype n); |
---|
205 | |
---|
206 | |
---|
207 | #ifdef __cplusplus |
---|
208 | } |
---|
209 | #endif |
---|
210 | |
---|
211 | #endif |
---|