SuperLU 7.0.0
slu_ddefs.h
Go to the documentation of this file.
1
74#ifndef __SUPERLU_dSP_DEFS /* allow multiple inclusions */
75#define __SUPERLU_dSP_DEFS
76
77/*
78 * File name: dsp_defs.h
79 * Purpose: Sparse matrix types and function prototypes
80 * History:
81 */
82
83#ifdef _CRAY
84#include <fortran.h>
85#endif
86
87#include <math.h>
88#include <limits.h>
89#include <stdio.h>
90#include <stdlib.h>
91#include <stdint.h>
92#include <string.h>
93#include "slu_Cnames.h"
94#include "superlu_config.h"
95#include "supermatrix.h"
96#include "slu_util.h"
97
98
99/* -------- Prototypes -------- */
100
101#ifdef __cplusplus
102extern "C" {
103#endif
104
106extern void
107dgssv(superlu_options_t *, SuperMatrix *, int *, int *, SuperMatrix *,
109extern void
110dgssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *,
111 char *, double *, double *, SuperMatrix *, SuperMatrix *,
112 void *, int_t lwork, SuperMatrix *, SuperMatrix *,
113 double *, double *, double *, double *,
115 /* ILU */
116extern void
118 SuperMatrix *, SuperMatrix *, SuperLUStat_t *, int *);
119extern void
120dgsisx(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r,
121 int *etree, char *equed, double *R, double *C,
122 SuperMatrix *L, SuperMatrix *U, void *work, int_t lwork,
123 SuperMatrix *B, SuperMatrix *X, double *recip_pivot_growth, double *rcond,
124 GlobalLU_t *Glu, mem_usage_t *mem_usage, SuperLUStat_t *stat, int_t *info);
125
126
128extern void
129dCreate_CompCol_Matrix(SuperMatrix *, int, int, int_t, double *,
131extern void
132dCreate_CompRow_Matrix(SuperMatrix *, int, int, int_t, double *,
134extern void dCompRow_to_CompCol(int, int, int_t, double*, int_t*, int_t*,
135 double **, int_t **, int_t **);
136extern void
138extern void
139dCreate_Dense_Matrix(SuperMatrix *, int, int, double *, int,
141extern void
142dCreate_SuperNode_Matrix(SuperMatrix *, int, int, int_t, double *,
143 int_t *, int_t *, int_t *, int *, int *,
145extern void
146dCopy_Dense_Matrix(int, int, double *, int, double *, int);
147
148extern void dallocateA (int, int_t, double **, int_t **, int_t **);
149extern void dgstrf (superlu_options_t*, SuperMatrix*,
150 int, int, int*, void *, int_t, int *, int *,
152 SuperLUStat_t*, int_t *info);
153extern int_t dsnode_dfs (const int, const int, const int_t *, const int_t *,
154 const int_t *, int_t *, int *, GlobalLU_t *);
155extern int dsnode_bmod (const int, const int, const int, double *,
156 double *, GlobalLU_t *, SuperLUStat_t*);
157extern void dpanel_dfs (const int, const int, const int, SuperMatrix *,
158 int *, int *, double *, int *, int *, int *,
159 int_t *, int *, int *, int_t *, GlobalLU_t *);
160extern void dpanel_bmod (const int, const int, const int, const int,
161 double *, double *, int *, int *,
163extern int dcolumn_dfs (const int, const int, int *, int *, int *, int *,
164 int *, int_t *, int *, int *, int_t *, GlobalLU_t *);
165extern int dcolumn_bmod (const int, const int, double *,
166 double *, int *, int *, int,
168extern int dcopy_to_ucol (int, int, int *, int *, int *,
169 double *, GlobalLU_t *);
170extern int dpivotL (const int, const double, int *, int *,
171 int *, int *, int *, GlobalLU_t *, SuperLUStat_t*);
172extern void dpruneL (const int, const int *, const int, const int,
173 const int *, const int *, int_t *, GlobalLU_t *);
174extern void dreadmt (int *, int *, int_t *, double **, int_t **, int_t **);
175extern void dGenXtrue (int, int, double *, int);
176extern void dFillRHS (trans_t, int, double *, int, SuperMatrix *,
177 SuperMatrix *);
178extern void dgstrs (trans_t, SuperMatrix *, SuperMatrix *, int *, int *,
179 SuperMatrix *, SuperLUStat_t*, int *);
180/* ILU */
181extern void dgsitrf (superlu_options_t*, SuperMatrix*, int, int, int*,
182 void *, int_t, int *, int *, SuperMatrix *, SuperMatrix *,
183 GlobalLU_t *, SuperLUStat_t*, int_t *info);
184extern int dldperm(int, int, int_t, int_t [], int_t [], double [],
185 int [], double [], double []);
186extern int ilu_dsnode_dfs (const int, const int, const int_t *, const int_t *,
187 const int_t *, int *, GlobalLU_t *);
188extern void ilu_dpanel_dfs (const int, const int, const int, SuperMatrix *,
189 int *, int *, double *, double *, int *, int *,
190 int *, int *, int *, int_t *, GlobalLU_t *);
191extern int ilu_dcolumn_dfs (const int, const int, int *, int *, int *,
192 int *, int *, int *, int *, int_t *, GlobalLU_t *);
193extern int ilu_dcopy_to_ucol (int, int, int *, int *, int *,
194 double *, int, milu_t, double, int,
195 double *, int *, GlobalLU_t *, double *);
196extern int ilu_dpivotL (const int, const double, int *, int *, int, int *,
197 int *, int *, int *, double, milu_t,
198 double, GlobalLU_t *, SuperLUStat_t*);
199extern int ilu_ddrop_row (superlu_options_t *, int, int, double,
200 int, int *, double *, GlobalLU_t *,
201 double *, double *, int);
202
203
206extern void dgsequ (SuperMatrix *, double *, double *, double *,
207 double *, double *, int *);
208extern void dlaqgs (SuperMatrix *, double *, double *, double,
209 double, double, char *);
210extern void dgscon (char *, SuperMatrix *, SuperMatrix *,
211 double, double *, SuperLUStat_t*, int *);
212extern double dPivotGrowth(int, SuperMatrix *, int *,
214extern void dgsrfs (trans_t, SuperMatrix *, SuperMatrix *,
215 SuperMatrix *, int *, int *, char *, double *,
216 double *, SuperMatrix *, SuperMatrix *,
217 double *, double *, SuperLUStat_t*, int *);
218
219extern int sp_dtrsv (char *, char *, char *, SuperMatrix *,
220 SuperMatrix *, double *, SuperLUStat_t*, int *);
221extern int sp_dgemv (char *, double, SuperMatrix *, double *,
222 int, double, double *, int);
223
224extern int sp_dgemm (char *, char *, int, int, int, double,
225 SuperMatrix *, double *, int, double,
226 double *, int);
227extern double dmach(char *); /* from C99 standard, in float.h */
228
230extern int_t dLUMemInit (fact_t, void *, int_t, int, int, int_t, int,
231 double, SuperMatrix *, SuperMatrix *,
232 GlobalLU_t *, int **, double **);
233extern void dSetRWork (int, int, double *, double **, double **);
234extern void dLUWorkFree (int *, double *, GlobalLU_t *);
235extern int_t dLUMemXpand (int, int_t, MemType, int_t *, GlobalLU_t *);
236
237extern double *doubleMalloc(size_t);
238extern double *doubleCalloc(size_t);
239extern int_t dmemory_usage(const int_t, const int_t, const int_t, const int);
242
244extern void dreadhb(FILE *, int *, int *, int_t *, double **, int_t **, int_t **);
245extern void dreadrb(int *, int *, int_t *, double **, int_t **, int_t **);
246extern void dreadtriple(int *, int *, int_t *, double **, int_t **, int_t **);
247extern void dreadtriple_noheader(int *, int *, int_t *, double **, int_t **, int_t **);
248extern void dreadMM(FILE *, int *, int *, int_t *, double **, int_t **, int_t **);
249extern void dfill (double *, int, double);
250extern void dinf_norm_error (int, SuperMatrix *, double *);
251extern double dqselect(int, double *, int);
252
253
255extern void dPrint_CompCol_Matrix(char *, SuperMatrix *);
256extern void dPrint_SuperNode_Matrix(char *, SuperMatrix *);
257extern void dPrint_Dense_Matrix(char *, SuperMatrix *);
258extern void dprint_lu_col(char *, int, int, int_t *, GlobalLU_t *);
259extern int print_double_vec(char *, int, double *);
260extern void dcheck_tempv(int, double *);
261
264extern int dgemm_(const char*, const char*, const int*, const int*, const int*,
265 const double*, const double*, const int*, const double*,
266 const int*, const double*, double*, const int*);
267extern int dtrsv_(char*, char*, char*, int*, double*, int*,
268 double*, int*);
269extern int dtrsm_(char*, char*, char*, char*, int*, int*,
270 double*, double*, int*, double*, int*);
271extern int dgemv_(char *, int *, int *, double *, double *a, int *,
272 double *, int *, double *, double *, int *);
273
274extern void dusolve(int, int, double*, double*);
275extern void dlsolve(int, int, double*, double*);
276extern void dmatvec(int, int, int, double*, double*, double*);
277
278#ifdef __cplusplus
279 }
280#endif
281
282#endif /* __SUPERLU_dSP_DEFS */
283
#define X(I)
#define A(I, J)
#define U(I)
Macros defining how C routines will be called.
void dgsisv(superlu_options_t *, SuperMatrix *, int *, int *, SuperMatrix *, SuperMatrix *, SuperMatrix *, SuperLUStat_t *, int *)
void ilu_dpanel_dfs(const int, const int, const int, SuperMatrix *, int *, int *, double *, double *, int *, int *, int *, int *, int *, int_t *, GlobalLU_t *)
Definition: ilu_dpanel_dfs.c:56
double * doubleMalloc(size_t)
Definition: dmemory.c:689
void dCreate_Dense_Matrix(SuperMatrix *, int, int, double *, int, Stype_t, Dtype_t, Mtype_t)
Definition: dutil.c:103
int dpivotL(const int, const double, int *, int *, int *, int *, int *, GlobalLU_t *, SuperLUStat_t *)
Definition: dpivotL.c:67
void dreadtriple_noheader(int *, int *, int_t *, double **, int_t **, int_t **)
Read matrix in triplet format from stdin.
Definition: dreadtriple_noheader.c:38
int ilu_dpivotL(const int, const double, int *, int *, int, int *, int *, int *, int *, double, milu_t, double, GlobalLU_t *, SuperLUStat_t *)
Definition: ilu_dpivotL.c:57
void dgscon(char *, SuperMatrix *, SuperMatrix *, double, double *, SuperLUStat_t *, int *)
Definition: dgscon.c:85
void dCopy_CompCol_Matrix(SuperMatrix *, SuperMatrix *)
Copy matrix A into matrix B.
Definition: dutil.c:82
void dpanel_bmod(const int, const int, const int, const int, double *, double *, int *, int *, GlobalLU_t *, SuperLUStat_t *)
Definition: dpanel_bmod.c:61
void dreadtriple(int *, int *, int_t *, double **, int_t **, int_t **)
Definition: dreadtriple.c:26
void dCreate_CompCol_Matrix(SuperMatrix *, int, int, int_t, double *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
Supernodal LU factor related.
Definition: dutil.c:39
int dsnode_bmod(const int, const int, const int, double *, double *, GlobalLU_t *, SuperLUStat_t *)
Performs numeric block updates within the relaxed snode.
Definition: dsnode_bmod.c:41
double dqselect(int, double *, int)
double dPivotGrowth(int, SuperMatrix *, int *, SuperMatrix *, SuperMatrix *)
Definition: dpivotgrowth.c:59
void dpanel_dfs(const int, const int, const int, SuperMatrix *, int *, int *, double *, int *, int *, int *, int_t *, int *, int *, int_t *, GlobalLU_t *)
Definition: dpanel_dfs.c:69
void dreadrb(int *, int *, int_t *, double **, int_t **, int_t **)
Definition: dreadrb.c:281
int dcolumn_bmod(const int, const int, double *, double *, int *, int *, int, GlobalLU_t *, SuperLUStat_t *)
Definition: dcolumn_bmod.c:52
void dreadhb(FILE *, int *, int *, int_t *, double **, int_t **, int_t **)
Auxiliary routines.
Definition: dreadhb.c:283
int_t dsnode_dfs(const int, const int, const int_t *, const int_t *, const int_t *, int_t *, int *, GlobalLU_t *)
Definition: dsnode_dfs.c:55
void dgsrfs(trans_t, SuperMatrix *, SuperMatrix *, SuperMatrix *, int *, int *, char *, double *, double *, SuperMatrix *, SuperMatrix *, double *, double *, SuperLUStat_t *, int *)
Definition: dgsrfs.c:141
void dgssvx(superlu_options_t *, SuperMatrix *, int *, int *, int *, char *, double *, double *, SuperMatrix *, SuperMatrix *, void *, int_t lwork, SuperMatrix *, SuperMatrix *, double *, double *, double *, double *, GlobalLU_t *, mem_usage_t *, SuperLUStat_t *, int_t *info)
Definition: dgssvx.c:358
int dldperm(int, int, int_t, int_t[], int_t[], double[], int[], double[], double[])
void dPrint_CompCol_Matrix(char *, SuperMatrix *)
Routines for debugging.
Definition: dutil.c:202
void dallocateA(int, int_t, double **, int_t **, int_t **)
Allocate storage for original matrix A.
Definition: dmemory.c:681
int sp_dgemv(char *, double, SuperMatrix *, double *, int, double, double *, int)
Performs one of the matrix-vector operations y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y,
Definition: dsp_blas2.c:374
double dmach(char *)
Definition: dmach.c:18
void dlsolve(int, int, double *, double *)
Solves a dense UNIT lower triangular system.
Definition: dmyblas2.c:38
void dgssv(superlu_options_t *, SuperMatrix *, int *, int *, SuperMatrix *, SuperMatrix *, SuperMatrix *, SuperLUStat_t *, int_t *info)
Driver routines.
Definition: dgssv.c:144
void dgstrf(superlu_options_t *, SuperMatrix *, int, int, int *, void *, int_t, int *, int *, SuperMatrix *, SuperMatrix *, GlobalLU_t *, SuperLUStat_t *, int_t *info)
Definition: dgstrf.c:201
void dSetRWork(int, int, double *, double **, double **)
Set up pointers for real working arrays.
Definition: dmemory.c:392
void dmatvec(int, int, int, double *, double *, double *)
Performs a dense matrix-vector multiply: Mxvec = Mxvec + M * vec.
Definition: dmyblas2.c:161
int ilu_dcopy_to_ucol(int, int, int *, int *, int *, double *, int, milu_t, double, int, double *, int *, GlobalLU_t *, double *)
Definition: ilu_dcopy_to_ucol.c:44
int dcolumn_dfs(const int, const int, int *, int *, int *, int *, int *, int_t *, int *, int *, int_t *, GlobalLU_t *)
Definition: dcolumn_dfs.c:76
int ilu_dQuerySpace(SuperMatrix *, SuperMatrix *, mem_usage_t *)
Definition: dmemory.c:149
void dlaqgs(SuperMatrix *, double *, double *, double, double, double, char *)
Definition: dlaqgs.c:92
void dreadMM(FILE *, int *, int *, int_t *, double **, int_t **, int_t **)
Definition: dreadMM.c:35
void dgsitrf(superlu_options_t *, SuperMatrix *, int, int, int *, void *, int_t, int *, int *, SuperMatrix *, SuperMatrix *, GlobalLU_t *, SuperLUStat_t *, int_t *info)
Definition: dgsitrf.c:188
int ilu_dsnode_dfs(const int, const int, const int_t *, const int_t *, const int_t *, int *, GlobalLU_t *)
Definition: ilu_dsnode_dfs.c:42
int dtrsv_(char *, char *, char *, int *, double *, int *, double *, int *)
void dgsisx(superlu_options_t *options, SuperMatrix *A, int *perm_c, int *perm_r, int *etree, char *equed, double *R, double *C, SuperMatrix *L, SuperMatrix *U, void *work, int_t lwork, SuperMatrix *B, SuperMatrix *X, double *recip_pivot_growth, double *rcond, GlobalLU_t *Glu, mem_usage_t *mem_usage, SuperLUStat_t *stat, int_t *info)
Definition: dgsisx.c:404
void dgstrs(trans_t, SuperMatrix *, SuperMatrix *, int *, int *, SuperMatrix *, SuperLUStat_t *, int *)
Definition: dgstrs.c:93
int dtrsm_(char *, char *, char *, char *, int *, int *, double *, double *, int *, double *, int *)
void dCopy_Dense_Matrix(int, int, double *, int, double *, int)
Definition: dutil.c:121
int dcopy_to_ucol(int, int, int *, int *, int *, double *, GlobalLU_t *)
Definition: dcopy_to_ucol.c:36
int dQuerySpace(SuperMatrix *, SuperMatrix *, mem_usage_t *)
Definition: dmemory.c:111
void dpruneL(const int, const int *, const int, const int, const int *, const int *, int_t *, GlobalLU_t *)
Definition: dpruneL.c:48
void dFillRHS(trans_t, int, double *, int, SuperMatrix *, SuperMatrix *)
Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's.
Definition: dutil.c:366
int print_double_vec(char *, int, double *)
Definition: dutil.c:473
void dfill(double *, int, double)
Fills a double precision array with a given value.
Definition: dutil.c:393
int sp_dgemm(char *, char *, int, int, int, double, SuperMatrix *, double *, int, double, double *, int)
Definition: dsp_blas3.c:126
int dgemv_(char *, int *, int *, double *, double *a, int *, double *, int *, double *, double *, int *)
int ilu_ddrop_row(superlu_options_t *, int, int, double, int, int *, double *, GlobalLU_t *, double *, double *, int)
void dCompRow_to_CompCol(int, int, int_t, double *, int_t *, int_t *, double **, int_t **, int_t **)
Convert a row compressed storage into a column compressed storage.
Definition: dutil.c:164
void dusolve(int, int, double *, double *)
Solves a dense upper triangular system.
Definition: dmyblas2.c:136
int dgemm_(const char *, const char *, const int *, const int *, const int *, const double *, const double *, const int *, const double *, const int *, const double *, double *, const int *)
BLAS.
double * doubleCalloc(size_t)
Definition: dmemory.c:699
void dcheck_tempv(int, double *)
Check whether tempv[] == 0. This should be true before and after calling any numeric routines,...
Definition: dutil.c:339
void dreadmt(int *, int *, int_t *, double **, int_t **, int_t **)
void dprint_lu_col(char *, int, int, int_t *, GlobalLU_t *)
Diagnostic print of column "jcol" in the U/L factor.
Definition: dutil.c:298
void dPrint_Dense_Matrix(char *, SuperMatrix *)
Definition: dutil.c:276
void dPrint_SuperNode_Matrix(char *, SuperMatrix *)
Definition: dutil.c:226
void dinf_norm_error(int, SuperMatrix *, double *)
Check the inf-norm of the error vector.
Definition: dutil.c:403
int_t dmemory_usage(const int_t, const int_t, const int_t, const int)
Definition: dmemory.c:713
void dLUWorkFree(int *, double *, GlobalLU_t *)
Free the working storage used by factor routines.
Definition: dmemory.c:407
void dCreate_CompRow_Matrix(SuperMatrix *, int, int, int_t, double *, int_t *, int_t *, Stype_t, Dtype_t, Mtype_t)
Definition: dutil.c:60
void dGenXtrue(int, int, double *, int)
Definition: dutil.c:354
int_t dLUMemXpand(int, int_t, MemType, int_t *, GlobalLU_t *)
Expand the data structures for L and U during the factorization.
Definition: dmemory.c:430
int_t dLUMemInit(fact_t, void *, int_t, int, int, int_t, int, double, SuperMatrix *, SuperMatrix *, GlobalLU_t *, int **, double **)
Memory-related.
Definition: dmemory.c:190
int ilu_dcolumn_dfs(const int, const int, int *, int *, int *, int *, int *, int *, int *, int_t *, GlobalLU_t *)
Definition: ilu_dcolumn_dfs.c:61
int sp_dtrsv(char *, char *, char *, SuperMatrix *, SuperMatrix *, double *, SuperLUStat_t *, int *)
Solves one of the systems of equations A*x = b, or A'*x = b.
Definition: dsp_blas2.c:86
void dgsequ(SuperMatrix *, double *, double *, double *, double *, double *, int *)
Driver related.
Definition: dgsequ.c:94
void dCreate_SuperNode_Matrix(SuperMatrix *, int, int, int_t, double *, int_t *, int_t *, int_t *, int *, int *, Stype_t, Dtype_t, Mtype_t)
Definition: dutil.c:134
Utility header file.
Definition: slu_util.h:336
Definition: slu_util.h:321
Definition: supermatrix.h:54
Definition: slu_util.h:330
Definition: slu_util.h:277
int int_t
Definition: superlu_config.h:20
trans_t
Definition: superlu_enum_consts.h:34
milu_t
Definition: superlu_enum_consts.h:46
MemType
Definition: superlu_enum_consts.h:37
fact_t
Definition: superlu_enum_consts.h:30
Matrix type definitions.
Mtype_t
Definition: supermatrix.h:42
Dtype_t
Definition: supermatrix.h:35
Stype_t
Definition: supermatrix.h:22