|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! @file sutil.c
|
|
kusano |
7d535a |
* \brief Matrix utility functions
|
|
kusano |
7d535a |
*
|
|
kusano |
7d535a |
*
|
|
kusano |
7d535a |
* -- SuperLU routine (version 3.1) --
|
|
kusano |
7d535a |
* Univ. of California Berkeley, Xerox Palo Alto Research Center,
|
|
kusano |
7d535a |
* and Lawrence Berkeley National Lab.
|
|
kusano |
7d535a |
* August 1, 2008
|
|
kusano |
7d535a |
*
|
|
kusano |
7d535a |
* Copyright (c) 1994 by Xerox Corporation. All rights reserved.
|
|
kusano |
7d535a |
*
|
|
kusano |
7d535a |
* THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY
|
|
kusano |
7d535a |
* EXPRESSED OR IMPLIED. ANY USE IS AT YOUR OWN RISK.
|
|
kusano |
7d535a |
*
|
|
kusano |
7d535a |
* Permission is hereby granted to use or copy this program for any
|
|
kusano |
7d535a |
* purpose, provided the above notices are retained on all copies.
|
|
kusano |
7d535a |
* Permission to modify the code and to distribute modified code is
|
|
kusano |
7d535a |
* granted, provided the above notices are retained, and a notice that
|
|
kusano |
7d535a |
* the code was modified is included with the above copyright notice.
|
|
kusano |
7d535a |
*
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
#include <math.h></math.h>
|
|
kusano |
7d535a |
#include "slu_sdefs.h"
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCreate_CompCol_Matrix(SuperMatrix *A, int m, int n, int nnz,
|
|
kusano |
7d535a |
float *nzval, int *rowind, int *colptr,
|
|
kusano |
7d535a |
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
NCformat *Astore;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
A->Stype = stype;
|
|
kusano |
7d535a |
A->Dtype = dtype;
|
|
kusano |
7d535a |
A->Mtype = mtype;
|
|
kusano |
7d535a |
A->nrow = m;
|
|
kusano |
7d535a |
A->ncol = n;
|
|
kusano |
7d535a |
A->Store = (void *) SUPERLU_MALLOC( sizeof(NCformat) );
|
|
kusano |
7d535a |
if ( !(A->Store) ) ABORT("SUPERLU_MALLOC fails for A->Store");
|
|
kusano |
7d535a |
Astore = A->Store;
|
|
kusano |
7d535a |
Astore->nnz = nnz;
|
|
kusano |
7d535a |
Astore->nzval = nzval;
|
|
kusano |
7d535a |
Astore->rowind = rowind;
|
|
kusano |
7d535a |
Astore->colptr = colptr;
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCreate_CompRow_Matrix(SuperMatrix *A, int m, int n, int nnz,
|
|
kusano |
7d535a |
float *nzval, int *colind, int *rowptr,
|
|
kusano |
7d535a |
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
NRformat *Astore;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
A->Stype = stype;
|
|
kusano |
7d535a |
A->Dtype = dtype;
|
|
kusano |
7d535a |
A->Mtype = mtype;
|
|
kusano |
7d535a |
A->nrow = m;
|
|
kusano |
7d535a |
A->ncol = n;
|
|
kusano |
7d535a |
A->Store = (void *) SUPERLU_MALLOC( sizeof(NRformat) );
|
|
kusano |
7d535a |
if ( !(A->Store) ) ABORT("SUPERLU_MALLOC fails for A->Store");
|
|
kusano |
7d535a |
Astore = A->Store;
|
|
kusano |
7d535a |
Astore->nnz = nnz;
|
|
kusano |
7d535a |
Astore->nzval = nzval;
|
|
kusano |
7d535a |
Astore->colind = colind;
|
|
kusano |
7d535a |
Astore->rowptr = rowptr;
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Copy matrix A into matrix B. */
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCopy_CompCol_Matrix(SuperMatrix *A, SuperMatrix *B)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
NCformat *Astore, *Bstore;
|
|
kusano |
7d535a |
int ncol, nnz, i;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
B->Stype = A->Stype;
|
|
kusano |
7d535a |
B->Dtype = A->Dtype;
|
|
kusano |
7d535a |
B->Mtype = A->Mtype;
|
|
kusano |
7d535a |
B->nrow = A->nrow;;
|
|
kusano |
7d535a |
B->ncol = ncol = A->ncol;
|
|
kusano |
7d535a |
Astore = (NCformat *) A->Store;
|
|
kusano |
7d535a |
Bstore = (NCformat *) B->Store;
|
|
kusano |
7d535a |
Bstore->nnz = nnz = Astore->nnz;
|
|
kusano |
7d535a |
for (i = 0; i < nnz; ++i)
|
|
kusano |
7d535a |
((float *)Bstore->nzval)[i] = ((float *)Astore->nzval)[i];
|
|
kusano |
7d535a |
for (i = 0; i < nnz; ++i) Bstore->rowind[i] = Astore->rowind[i];
|
|
kusano |
7d535a |
for (i = 0; i <= ncol; ++i) Bstore->colptr[i] = Astore->colptr[i];
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCreate_Dense_Matrix(SuperMatrix *X, int m, int n, float *x, int ldx,
|
|
kusano |
7d535a |
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
DNformat *Xstore;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
X->Stype = stype;
|
|
kusano |
7d535a |
X->Dtype = dtype;
|
|
kusano |
7d535a |
X->Mtype = mtype;
|
|
kusano |
7d535a |
X->nrow = m;
|
|
kusano |
7d535a |
X->ncol = n;
|
|
kusano |
7d535a |
X->Store = (void *) SUPERLU_MALLOC( sizeof(DNformat) );
|
|
kusano |
7d535a |
if ( !(X->Store) ) ABORT("SUPERLU_MALLOC fails for X->Store");
|
|
kusano |
7d535a |
Xstore = (DNformat *) X->Store;
|
|
kusano |
7d535a |
Xstore->lda = ldx;
|
|
kusano |
7d535a |
Xstore->nzval = (float *) x;
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCopy_Dense_Matrix(int M, int N, float *X, int ldx,
|
|
kusano |
7d535a |
float *Y, int ldy)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
/*! \brief Copies a two-dimensional matrix X to another matrix Y.
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
int i, j;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
for (j = 0; j < N; ++j)
|
|
kusano |
7d535a |
for (i = 0; i < M; ++i)
|
|
kusano |
7d535a |
Y[i + j*ldy] = X[i + j*ldx];
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCreate_SuperNode_Matrix(SuperMatrix *L, int m, int n, int nnz,
|
|
kusano |
7d535a |
float *nzval, int *nzval_colptr, int *rowind,
|
|
kusano |
7d535a |
int *rowind_colptr, int *col_to_sup, int *sup_to_col,
|
|
kusano |
7d535a |
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
SCformat *Lstore;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
L->Stype = stype;
|
|
kusano |
7d535a |
L->Dtype = dtype;
|
|
kusano |
7d535a |
L->Mtype = mtype;
|
|
kusano |
7d535a |
L->nrow = m;
|
|
kusano |
7d535a |
L->ncol = n;
|
|
kusano |
7d535a |
L->Store = (void *) SUPERLU_MALLOC( sizeof(SCformat) );
|
|
kusano |
7d535a |
if ( !(L->Store) ) ABORT("SUPERLU_MALLOC fails for L->Store");
|
|
kusano |
7d535a |
Lstore = L->Store;
|
|
kusano |
7d535a |
Lstore->nnz = nnz;
|
|
kusano |
7d535a |
Lstore->nsuper = col_to_sup[n];
|
|
kusano |
7d535a |
Lstore->nzval = nzval;
|
|
kusano |
7d535a |
Lstore->nzval_colptr = nzval_colptr;
|
|
kusano |
7d535a |
Lstore->rowind = rowind;
|
|
kusano |
7d535a |
Lstore->rowind_colptr = rowind_colptr;
|
|
kusano |
7d535a |
Lstore->col_to_sup = col_to_sup;
|
|
kusano |
7d535a |
Lstore->sup_to_col = sup_to_col;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Convert a row compressed storage into a column compressed storage.
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sCompRow_to_CompCol(int m, int n, int nnz,
|
|
kusano |
7d535a |
float *a, int *colind, int *rowptr,
|
|
kusano |
7d535a |
float **at, int **rowind, int **colptr)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
register int i, j, col, relpos;
|
|
kusano |
7d535a |
int *marker;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/* Allocate storage for another copy of the matrix. */
|
|
kusano |
7d535a |
*at = (float *) floatMalloc(nnz);
|
|
kusano |
7d535a |
*rowind = (int *) intMalloc(nnz);
|
|
kusano |
7d535a |
*colptr = (int *) intMalloc(n+1);
|
|
kusano |
7d535a |
marker = (int *) intCalloc(n);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/* Get counts of each column of A, and set up column pointers */
|
|
kusano |
7d535a |
for (i = 0; i < m; ++i)
|
|
kusano |
7d535a |
for (j = rowptr[i]; j < rowptr[i+1]; ++j) ++marker[colind[j]];
|
|
kusano |
7d535a |
(*colptr)[0] = 0;
|
|
kusano |
7d535a |
for (j = 0; j < n; ++j) {
|
|
kusano |
7d535a |
(*colptr)[j+1] = (*colptr)[j] + marker[j];
|
|
kusano |
7d535a |
marker[j] = (*colptr)[j];
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/* Transfer the matrix into the compressed column storage. */
|
|
kusano |
7d535a |
for (i = 0; i < m; ++i) {
|
|
kusano |
7d535a |
for (j = rowptr[i]; j < rowptr[i+1]; ++j) {
|
|
kusano |
7d535a |
col = colind[j];
|
|
kusano |
7d535a |
relpos = marker[col];
|
|
kusano |
7d535a |
(*rowind)[relpos] = i;
|
|
kusano |
7d535a |
(*at)[relpos] = a[j];
|
|
kusano |
7d535a |
++marker[col];
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
SUPERLU_FREE(marker);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sPrint_CompCol_Matrix(char *what, SuperMatrix *A)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
NCformat *Astore;
|
|
kusano |
7d535a |
register int i,n;
|
|
kusano |
7d535a |
float *dp;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("\nCompCol matrix %s:\n", what);
|
|
kusano |
7d535a |
printf("Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
|
|
kusano |
7d535a |
n = A->ncol;
|
|
kusano |
7d535a |
Astore = (NCformat *) A->Store;
|
|
kusano |
7d535a |
dp = (float *) Astore->nzval;
|
|
kusano |
7d535a |
printf("nrow %d, ncol %d, nnz %d\n", A->nrow,A->ncol,Astore->nnz);
|
|
kusano |
7d535a |
printf("nzval: ");
|
|
kusano |
7d535a |
for (i = 0; i < Astore->colptr[n]; ++i) printf("%f ", dp[i]);
|
|
kusano |
7d535a |
printf("\nrowind: ");
|
|
kusano |
7d535a |
for (i = 0; i < Astore->colptr[n]; ++i) printf("%d ", Astore->rowind[i]);
|
|
kusano |
7d535a |
printf("\ncolptr: ");
|
|
kusano |
7d535a |
for (i = 0; i <= n; ++i) printf("%d ", Astore->colptr[i]);
|
|
kusano |
7d535a |
printf("\n");
|
|
kusano |
7d535a |
fflush(stdout);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sPrint_SuperNode_Matrix(char *what, SuperMatrix *A)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
SCformat *Astore;
|
|
kusano |
7d535a |
register int i, j, k, c, d, n, nsup;
|
|
kusano |
7d535a |
float *dp;
|
|
kusano |
7d535a |
int *col_to_sup, *sup_to_col, *rowind, *rowind_colptr;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("\nSuperNode matrix %s:\n", what);
|
|
kusano |
7d535a |
printf("Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
|
|
kusano |
7d535a |
n = A->ncol;
|
|
kusano |
7d535a |
Astore = (SCformat *) A->Store;
|
|
kusano |
7d535a |
dp = (float *) Astore->nzval;
|
|
kusano |
7d535a |
col_to_sup = Astore->col_to_sup;
|
|
kusano |
7d535a |
sup_to_col = Astore->sup_to_col;
|
|
kusano |
7d535a |
rowind_colptr = Astore->rowind_colptr;
|
|
kusano |
7d535a |
rowind = Astore->rowind;
|
|
kusano |
7d535a |
printf("nrow %d, ncol %d, nnz %d, nsuper %d\n",
|
|
kusano |
7d535a |
A->nrow,A->ncol,Astore->nnz,Astore->nsuper);
|
|
kusano |
7d535a |
printf("nzval:\n");
|
|
kusano |
7d535a |
for (k = 0; k <= Astore->nsuper; ++k) {
|
|
kusano |
7d535a |
c = sup_to_col[k];
|
|
kusano |
7d535a |
nsup = sup_to_col[k+1] - c;
|
|
kusano |
7d535a |
for (j = c; j < c + nsup; ++j) {
|
|
kusano |
7d535a |
d = Astore->nzval_colptr[j];
|
|
kusano |
7d535a |
for (i = rowind_colptr[c]; i < rowind_colptr[c+1]; ++i) {
|
|
kusano |
7d535a |
printf("%d\t%d\t%e\n", rowind[i], j, dp[d++]);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
#if 0
|
|
kusano |
7d535a |
for (i = 0; i < Astore->nzval_colptr[n]; ++i) printf("%f ", dp[i]);
|
|
kusano |
7d535a |
#endif
|
|
kusano |
7d535a |
printf("\nnzval_colptr: ");
|
|
kusano |
7d535a |
for (i = 0; i <= n; ++i) printf("%d ", Astore->nzval_colptr[i]);
|
|
kusano |
7d535a |
printf("\nrowind: ");
|
|
kusano |
7d535a |
for (i = 0; i < Astore->rowind_colptr[n]; ++i)
|
|
kusano |
7d535a |
printf("%d ", Astore->rowind[i]);
|
|
kusano |
7d535a |
printf("\nrowind_colptr: ");
|
|
kusano |
7d535a |
for (i = 0; i <= n; ++i) printf("%d ", Astore->rowind_colptr[i]);
|
|
kusano |
7d535a |
printf("\ncol_to_sup: ");
|
|
kusano |
7d535a |
for (i = 0; i < n; ++i) printf("%d ", col_to_sup[i]);
|
|
kusano |
7d535a |
printf("\nsup_to_col: ");
|
|
kusano |
7d535a |
for (i = 0; i <= Astore->nsuper+1; ++i)
|
|
kusano |
7d535a |
printf("%d ", sup_to_col[i]);
|
|
kusano |
7d535a |
printf("\n");
|
|
kusano |
7d535a |
fflush(stdout);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sPrint_Dense_Matrix(char *what, SuperMatrix *A)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
DNformat *Astore = (DNformat *) A->Store;
|
|
kusano |
7d535a |
register int i, j, lda = Astore->lda;
|
|
kusano |
7d535a |
float *dp;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("\nDense matrix %s:\n", what);
|
|
kusano |
7d535a |
printf("Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
|
|
kusano |
7d535a |
dp = (float *) Astore->nzval;
|
|
kusano |
7d535a |
printf("nrow %d, ncol %d, lda %d\n", A->nrow,A->ncol,lda);
|
|
kusano |
7d535a |
printf("\nnzval: ");
|
|
kusano |
7d535a |
for (j = 0; j < A->ncol; ++j) {
|
|
kusano |
7d535a |
for (i = 0; i < A->nrow; ++i) printf("%f ", dp[i + j*lda]);
|
|
kusano |
7d535a |
printf("\n");
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
printf("\n");
|
|
kusano |
7d535a |
fflush(stdout);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Diagnostic print of column "jcol" in the U/L factor.
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sprint_lu_col(char *msg, int jcol, int pivrow, int *xprune, GlobalLU_t *Glu)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
int i, k, fsupc;
|
|
kusano |
7d535a |
int *xsup, *supno;
|
|
kusano |
7d535a |
int *xlsub, *lsub;
|
|
kusano |
7d535a |
float *lusup;
|
|
kusano |
7d535a |
int *xlusup;
|
|
kusano |
7d535a |
float *ucol;
|
|
kusano |
7d535a |
int *usub, *xusub;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
xsup = Glu->xsup;
|
|
kusano |
7d535a |
supno = Glu->supno;
|
|
kusano |
7d535a |
lsub = Glu->lsub;
|
|
kusano |
7d535a |
xlsub = Glu->xlsub;
|
|
kusano |
7d535a |
lusup = Glu->lusup;
|
|
kusano |
7d535a |
xlusup = Glu->xlusup;
|
|
kusano |
7d535a |
ucol = Glu->ucol;
|
|
kusano |
7d535a |
usub = Glu->usub;
|
|
kusano |
7d535a |
xusub = Glu->xusub;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("%s", msg);
|
|
kusano |
7d535a |
printf("col %d: pivrow %d, supno %d, xprune %d\n",
|
|
kusano |
7d535a |
jcol, pivrow, supno[jcol], xprune[jcol]);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("\tU-col:\n");
|
|
kusano |
7d535a |
for (i = xusub[jcol]; i < xusub[jcol+1]; i++)
|
|
kusano |
7d535a |
printf("\t%d%10.4f\n", usub[i], ucol[i]);
|
|
kusano |
7d535a |
printf("\tL-col in rectangular snode:\n");
|
|
kusano |
7d535a |
fsupc = xsup[supno[jcol]]; /* first col of the snode */
|
|
kusano |
7d535a |
i = xlsub[fsupc];
|
|
kusano |
7d535a |
k = xlusup[jcol];
|
|
kusano |
7d535a |
while ( i < xlsub[fsupc+1] && k < xlusup[jcol+1] ) {
|
|
kusano |
7d535a |
printf("\t%d\t%10.4f\n", lsub[i], lusup[k]);
|
|
kusano |
7d535a |
i++; k++;
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
fflush(stdout);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Check whether tempv[] == 0. This should be true before and after calling any numeric routines, i.e., "panel_bmod" and "column_bmod".
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
void scheck_tempv(int n, float *tempv)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
int i;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
for (i = 0; i < n; i++) {
|
|
kusano |
7d535a |
if (tempv[i] != 0.0)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
fprintf(stderr,"tempv[%d] = %f\n", i,tempv[i]);
|
|
kusano |
7d535a |
ABORT("scheck_tempv");
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sGenXtrue(int n, int nrhs, float *x, int ldx)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
int i, j;
|
|
kusano |
7d535a |
for (j = 0; j < nrhs; ++j)
|
|
kusano |
7d535a |
for (i = 0; i < n; ++i) {
|
|
kusano |
7d535a |
x[i + j*ldx] = 1.0;/* + (float)(i+1.)/n;*/
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sFillRHS(trans_t trans, int nrhs, float *x, int ldx,
|
|
kusano |
7d535a |
SuperMatrix *A, SuperMatrix *B)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
NCformat *Astore;
|
|
kusano |
7d535a |
float *Aval;
|
|
kusano |
7d535a |
DNformat *Bstore;
|
|
kusano |
7d535a |
float *rhs;
|
|
kusano |
7d535a |
float one = 1.0;
|
|
kusano |
7d535a |
float zero = 0.0;
|
|
kusano |
7d535a |
int ldc;
|
|
kusano |
7d535a |
char transc[1];
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
Astore = A->Store;
|
|
kusano |
7d535a |
Aval = (float *) Astore->nzval;
|
|
kusano |
7d535a |
Bstore = B->Store;
|
|
kusano |
7d535a |
rhs = Bstore->nzval;
|
|
kusano |
7d535a |
ldc = Bstore->lda;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
if ( trans == NOTRANS ) *(unsigned char *)transc = 'N';
|
|
kusano |
7d535a |
else *(unsigned char *)transc = 'T';
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
sp_sgemm(transc, "N", A->nrow, nrhs, A->ncol, one, A,
|
|
kusano |
7d535a |
x, ldx, zero, rhs, ldc);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Fills a float precision array with a given value.
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sfill(float *a, int alen, float dval)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
register int i;
|
|
kusano |
7d535a |
for (i = 0; i < alen; i++) a[i] = dval;
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Check the inf-norm of the error vector
|
|
kusano |
7d535a |
*/
|
|
kusano |
7d535a |
void sinf_norm_error(int nrhs, SuperMatrix *X, float *xtrue)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
DNformat *Xstore;
|
|
kusano |
7d535a |
float err, xnorm;
|
|
kusano |
7d535a |
float *Xmat, *soln_work;
|
|
kusano |
7d535a |
int i, j;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
Xstore = X->Store;
|
|
kusano |
7d535a |
Xmat = Xstore->nzval;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
for (j = 0; j < nrhs; j++) {
|
|
kusano |
7d535a |
soln_work = &Xmat[j*Xstore->lda];
|
|
kusano |
7d535a |
err = xnorm = 0.0;
|
|
kusano |
7d535a |
for (i = 0; i < X->nrow; i++) {
|
|
kusano |
7d535a |
err = SUPERLU_MAX(err, fabs(soln_work[i] - xtrue[i]));
|
|
kusano |
7d535a |
xnorm = SUPERLU_MAX(xnorm, fabs(soln_work[i]));
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
err = err / xnorm;
|
|
kusano |
7d535a |
printf("||X - Xtrue||/||X|| = %e\n", err);
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
/*! \brief Print performance of the code. */
|
|
kusano |
7d535a |
void
|
|
kusano |
7d535a |
sPrintPerf(SuperMatrix *L, SuperMatrix *U, mem_usage_t *mem_usage,
|
|
kusano |
7d535a |
float rpg, float rcond, float *ferr,
|
|
kusano |
7d535a |
float *berr, char *equed, SuperLUStat_t *stat)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
SCformat *Lstore;
|
|
kusano |
7d535a |
NCformat *Ustore;
|
|
kusano |
7d535a |
double *utime;
|
|
kusano |
7d535a |
flops_t *ops;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
utime = stat->utime;
|
|
kusano |
7d535a |
ops = stat->ops;
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
if ( utime[FACT] != 0. )
|
|
kusano |
7d535a |
printf("Factor flops = %e\tMflops = %8.2f\n", ops[FACT],
|
|
kusano |
7d535a |
ops[FACT]*1e-6/utime[FACT]);
|
|
kusano |
7d535a |
printf("Identify relaxed snodes = %8.2f\n", utime[RELAX]);
|
|
kusano |
7d535a |
if ( utime[SOLVE] != 0. )
|
|
kusano |
7d535a |
printf("Solve flops = %.0f, Mflops = %8.2f\n", ops[SOLVE],
|
|
kusano |
7d535a |
ops[SOLVE]*1e-6/utime[SOLVE]);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
Lstore = (SCformat *) L->Store;
|
|
kusano |
7d535a |
Ustore = (NCformat *) U->Store;
|
|
kusano |
7d535a |
printf("\tNo of nonzeros in factor L = %d\n", Lstore->nnz);
|
|
kusano |
7d535a |
printf("\tNo of nonzeros in factor U = %d\n", Ustore->nnz);
|
|
kusano |
7d535a |
printf("\tNo of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("L\\U MB %.3f\ttotal MB needed %.3f\n",
|
|
kusano |
7d535a |
mem_usage->for_lu/1e6, mem_usage->total_needed/1e6);
|
|
kusano |
7d535a |
printf("Number of memory expansions: %d\n", stat->expansions);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("\tFactor\tMflops\tSolve\tMflops\tEtree\tEquil\tRcond\tRefine\n");
|
|
kusano |
7d535a |
printf("PERF:%8.2f%8.2f%8.2f%8.2f%8.2f%8.2f%8.2f%8.2f\n",
|
|
kusano |
7d535a |
utime[FACT], ops[FACT]*1e-6/utime[FACT],
|
|
kusano |
7d535a |
utime[SOLVE], ops[SOLVE]*1e-6/utime[SOLVE],
|
|
kusano |
7d535a |
utime[ETREE], utime[EQUIL], utime[RCOND], utime[REFINE]);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
printf("\tRpg\t\tRcond\t\tFerr\t\tBerr\t\tEquil?\n");
|
|
kusano |
7d535a |
printf("NUM:\t%e\t%e\t%e\t%e\t%s\n",
|
|
kusano |
7d535a |
rpg, rcond, ferr[0], berr[0], equed);
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
|
|
kusano |
7d535a |
print_float_vec(char *what, int n, float *vec)
|
|
kusano |
7d535a |
{
|
|
kusano |
7d535a |
int i;
|
|
kusano |
7d535a |
printf("%s: n %d\n", what, n);
|
|
kusano |
7d535a |
for (i = 0; i < n; ++i) printf("%d\t%f\n", i, vec[i]);
|
|
kusano |
7d535a |
return 0;
|
|
kusano |
7d535a |
}
|
|
kusano |
7d535a |
|