|
| for (j=0;j< N;++j) for(i=0 |
|
| if (! L->Store) |
|
void | sCopy_CompRowLoc_Matrix_dist (SuperMatrix *A, SuperMatrix *B) |
|
void | sZero_CompRowLoc_Matrix_dist (SuperMatrix *A) |
| Sets all entries of a matrix to zero, A_{i,j}=0, for i,j=1,..,n. More...
|
|
void | sScaleAddId_CompRowLoc_Matrix_dist (SuperMatrix *A, float c) |
| Scale and add I: scales a matrix and adds an identity. A_{i,j} = c * A_{i,j} + \delta_{i,j} for i,j=1,...,n and \delta_{i,j} is the Kronecker delta. More...
|
|
void | sScaleAdd_CompRowLoc_Matrix_dist (SuperMatrix *A, SuperMatrix *B, float c) |
| Scale and add: adds a scalar multiple of one matrix to another. A_{i,j} = c * A_{i,j} + B_{i,j}$ for i,j=1,...,n. More...
|
|
void | sScalePermstructInit (const int_t m, const int_t n, sScalePermstruct_t *ScalePermstruct) |
| Allocate storage in ScalePermstruct. More...
|
|
void | sScalePermstructFree (sScalePermstruct_t *ScalePermstruct) |
| Deallocate ScalePermstruct. More...
|
|
int | sAllocGlu_3d (int_t n, int_t nsupers, sLUstruct_t *LUstruct) |
|
int | sDeAllocGlu_3d (sLUstruct_t *LUstruct) |
|
int | sDeAllocLlu_3d (int_t n, sLUstruct_t *LUstruct, gridinfo3d_t *grid3d) |
|
void | sGenXtrue_dist (int_t n, int_t nrhs, float *x, int_t ldx) |
|
void | sFillRHS_dist (char *trans, int_t nrhs, float *x, int_t ldx, SuperMatrix *A, float *rhs, int_t ldb) |
| Let rhs[i] = sum of i-th row of A, so the solution vector is all 1's. More...
|
|
void | sfill_dist (float *a, int_t alen, float dval) |
| Fills a float precision array with a given value. More...
|
|
void | sinf_norm_error_dist (int_t n, int_t nrhs, float *x, int_t ldx, float *xtrue, int_t ldxtrue, gridinfo_t *grid) |
| Check the inf-norm of the error vector. More...
|
|
void | Printfloat5 (char *name, int_t len, float *x) |
|
int | file_Printfloat5 (FILE *fp, char *name, int_t len, float *x) |
|
void | sPrintLblocks (int iam, int_t nsupers, gridinfo_t *grid, Glu_persist_t *Glu_persist, sLocalLU_t *Llu) |
| Print the blocks in the factored matrix L. More...
|
|
void | sZeroLblocks (int iam, int n, gridinfo_t *grid, sLUstruct_t *LUstruct) |
| Sets all entries of matrix L to zero. More...
|
|
void | sDumpLblocks (int iam, int_t nsupers, gridinfo_t *grid, Glu_persist_t *Glu_persist, sLocalLU_t *Llu) |
| Dump the factored matrix L using matlab triple-let format. More...
|
|
void | sPrintUblocks (int iam, int_t nsupers, gridinfo_t *grid, Glu_persist_t *Glu_persist, sLocalLU_t *Llu) |
| Print the blocks in the factored matrix U. More...
|
|
void | sZeroUblocks (int iam, int n, gridinfo_t *grid, sLUstruct_t *LUstruct) |
| Sets all entries of matrix U to zero. More...
|
|
int | sprint_gsmv_comm (FILE *fp, int_t m_loc, psgsmv_comm_t *gsmv_comm, gridinfo_t *grid) |
|
void | sGenXtrueRHS (int nrhs, SuperMatrix *A, Glu_persist_t *Glu_persist, gridinfo_t *grid, float **xact, int *ldx, float **b, int *ldb) |
|
Several matrix utilities.
Copyright (c) 2003, The Regents of the University of California, through Lawrence Berkeley National Laboratory (subject to receipt of any required approvals from U.S. Dept. of Energy)
All rights reserved.
The source code is distributed under BSD license, see the file License.txt at the top-level directory.
-- Distributed SuperLU routine (version 7.1.0) --
Lawrence Berkeley National Lab, Univ. of California Berkeley.
March 15, 2003
October 5, 2021
*/
void
sCreate_CompCol_Matrix_dist(SuperMatrix *A, int_t m, int_t n, int_t nnz,
float *nzval, int_t *rowind, int_t *colptr,
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
{
NCformat *Astore;
A->Stype = stype;
A->Dtype = dtype;
A->Mtype = mtype;
A->nrow = m;
A->ncol = n;
A->Store = (void *) SUPERLU_MALLOC( sizeof(NCformat) );
if ( !(A->Store) ) ABORT("SUPERLU_MALLOC fails for A->Store");
Astore = (NCformat *) A->Store;
Astore->nnz = nnz;
Astore->nzval = nzval;
Astore->rowind = rowind;
Astore->colptr = colptr;
}
void
sCreate_CompRowLoc_Matrix_dist(SuperMatrix *A, int_t m, int_t n,
int_t nnz_loc, int_t m_loc, int_t fst_row,
float *nzval, int_t *colind, int_t *rowptr,
Stype_t stype, Dtype_t dtype, Mtype_t mtype)
{
NRformat_loc *Astore;
A->Stype = stype;
A->Dtype = dtype;
A->Mtype = mtype;
A->nrow = m;
A->ncol = n;
A->Store = (void *) SUPERLU_MALLOC( sizeof(NRformat_loc) );
if ( !(A->Store) ) ABORT("SUPERLU_MALLOC fails for A->Store");
Astore = (NRformat_loc *) A->Store;
Astore->nnz_loc = nnz_loc;
Astore->fst_row = fst_row;
Astore->m_loc = m_loc;
Astore->nzval = nzval;
Astore->colind = colind;
Astore->rowptr = rowptr;
}
/*!
Convert a row compressed storage into a column compressed storage.
*/
void
sCompRow_to_CompCol_dist(int_t m, int_t n, int_t nnz,
float *a, int_t *colind, int_t *rowptr,
float **at, int_t **rowind, int_t **colptr)
{
register int i, j, col, relpos;
int_t *marker;
/* Allocate storage for another copy of the matrix. */
*at = (float *) floatMalloc_dist(nnz);
*rowind = intMalloc_dist(nnz);
*colptr = intMalloc_dist(n+1);
marker = intCalloc_dist(n);
/* Get counts of each column of A, and set up column pointers */
for (i = 0; i < m; ++i)
for (j = rowptr[i]; j < rowptr[i+1]; ++j) ++marker[colind[j]];
(*colptr)[0] = 0;
for (j = 0; j < n; ++j) {
(*colptr)[j+1] = (*colptr)[j] + marker[j];
marker[j] = (*colptr)[j];
}
/* Transfer the matrix into the compressed column storage. */
for (i = 0; i < m; ++i) {
for (j = rowptr[i]; j < rowptr[i+1]; ++j) {
col = colind[j];
relpos = marker[col];
(*rowind)[relpos] = i;
(*at)[relpos] = a[j];
++marker[col];
}
}
SUPERLU_FREE(marker);
}
/*!
Copy matrix A into matrix B. */
void
sCopy_CompCol_Matrix_dist(SuperMatrix *A, SuperMatrix *B)
{
NCformat *Astore, *Bstore;
int ncol, nnz, i;
B->Stype = A->Stype;
B->Dtype = A->Dtype;
B->Mtype = A->Mtype;
B->nrow = A->nrow;;
B->ncol = ncol = A->ncol;
Astore = (NCformat *) A->Store;
Bstore = (NCformat *) B->Store;
Bstore->nnz = nnz = Astore->nnz;
for (i = 0; i < nnz; ++i)
((float *)Bstore->nzval)[i] = ((float *)Astore->nzval)[i];
for (i = 0; i < nnz; ++i) Bstore->rowind[i] = Astore->rowind[i];
for (i = 0; i <= ncol; ++i) Bstore->colptr[i] = Astore->colptr[i];
}
void sPrint_CompCol_Matrix_dist(SuperMatrix *A)
{
NCformat *Astore;
register int i;
float *dp;
printf("\nCompCol matrix: ");
printf("Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
Astore = (NCformat *) A->Store;
printf("nrow %lld, ncol %lld, nnz %lld\n", (long long) A->nrow,
(long long) A->ncol, (long long) Astore->nnz);
if ( (dp = (float *) Astore->nzval) != NULL ) {
printf("nzval:\n");
for (i = 0; i < Astore->nnz; ++i) printf("%f ", dp[i]);
}
printf("\nrowind:\n");
for (i = 0; i < Astore->nnz; ++i)
printf("%lld ", (long long) Astore->rowind[i]);
printf("\ncolptr:\n");
for (i = 0; i <= A->ncol; ++i)
printf("%lld ", (long long) Astore->colptr[i]);
printf("\nend CompCol matrix.\n");
}
void sPrint_Dense_Matrix_dist(SuperMatrix *A)
{
DNformat *Astore;
register int i;
float *dp;
printf("\nDense matrix: ");
printf("Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
Astore = (DNformat *) A->Store;
dp = (float *) Astore->nzval;
printf("nrow %lld, ncol %lld, lda %lld\n",
(long long) A->nrow, (long long) A->ncol, (long long) Astore->lda);
printf("\nnzval: ");
for (i = 0; i < A->nrow; ++i) printf("%f ", dp[i]);
printf("\nend Dense matrix.\n");
}
int sPrint_CompRowLoc_Matrix_dist(SuperMatrix *A)
{
NRformat_loc *Astore;
int_t nnz_loc, m_loc;
float *dp;
printf("\n==== CompRowLoc matrix: ");
printf("Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
Astore = (NRformat_loc *) A->Store;
printf("nrow %ld, ncol %ld\n",
(long int) A->nrow, (long int) A->ncol);
nnz_loc = Astore->nnz_loc; m_loc = Astore->m_loc;
printf("nnz_loc %ld, m_loc %ld, fst_row %ld\n", (long int) nnz_loc,
(long int) m_loc, (long int) Astore->fst_row);
PrintInt10("rowptr", m_loc+1, Astore->rowptr);
PrintInt10("colind", nnz_loc, Astore->colind);
if ( (dp = (float *) Astore->nzval) != NULL )
Printfloat5("nzval", nnz_loc, dp);
printf("==== end CompRowLoc matrix\n");
return 0;
}
int file_sPrint_CompRowLoc_Matrix_dist(FILE *fp, SuperMatrix *A)
{
NRformat_loc *Astore;
int_t nnz_loc, m_loc;
float *dp;
fprintf(fp, "\n==== CompRowLoc matrix: ");
fprintf(fp, "Stype %d, Dtype %d, Mtype %d\n", A->Stype,A->Dtype,A->Mtype);
Astore = (NRformat_loc *) A->Store;
fprintf(fp, "nrow %ld, ncol %ld\n", (long int) A->nrow, (long int) A->ncol);
nnz_loc = Astore->nnz_loc; m_loc = Astore->m_loc;
fprintf(fp, "nnz_loc %ld, m_loc %ld, fst_row %ld\n", (long int) nnz_loc,
(long int) m_loc, (long int) Astore->fst_row);
file_PrintInt10(fp, "rowptr", m_loc+1, Astore->rowptr);
file_PrintInt10(fp, "colind", nnz_loc, Astore->colind);
if ( (dp = (float *) Astore->nzval) != NULL )
file_Printfloat5(fp, "nzval", nnz_loc, dp);
fprintf(fp, "==== end CompRowLoc matrix\n");
return 0;
}
void
sCreate_Dense_Matrix_dist(SuperMatrix *X, int_t m, int_t n, float *x,
int_t ldx, Stype_t stype, Dtype_t dtype,
Mtype_t mtype)
{
DNformat *Xstore;
X->Stype = stype;
X->Dtype = dtype;
X->Mtype = mtype;
X->nrow = m;
X->ncol = n;
X->Store = (void *) SUPERLU_MALLOC( sizeof(DNformat) );
if ( !(X->Store) ) ABORT("SUPERLU_MALLOC fails for X->Store");
Xstore = (DNformat *) X->Store;
Xstore->lda = ldx;
Xstore->nzval = (float *) x;
}
void
sCopy_Dense_Matrix_dist(int_t M, int_t N, float *X, int_t ldx,
float *Y, int_t ldy)
{
/*!
Purpose
=======
Copies a two-dimensional matrix X to another matrix Y.