1 2 /* 3 Provides an interface to the SuperLU_DIST sparse solver 4 */ 5 6 #include <../src/mat/impls/aij/seq/aij.h> 7 #include <../src/mat/impls/aij/mpi/mpiaij.h> 8 #if defined(PETSC_HAVE_STDLIB_H) /* This is to get around weird problem with SuperLU on cray */ 9 #include <stdlib.h> 10 #endif 11 12 EXTERN_C_BEGIN 13 #if defined(PETSC_USE_COMPLEX) 14 #include <superlu_zdefs.h> 15 #else 16 #include <superlu_ddefs.h> 17 #endif 18 EXTERN_C_END 19 20 /* 21 GLOBAL - The sparse matrix and right hand side are all stored initially on process 0. Should be called centralized 22 DISTRIBUTED - The sparse matrix and right hand size are initially stored across the entire MPI communicator. 23 */ 24 typedef enum {GLOBAL,DISTRIBUTED} SuperLU_MatInputMode; 25 const char *SuperLU_MatInputModes[] = {"GLOBAL","DISTRIBUTED","SuperLU_MatInputMode","PETSC_",0}; 26 27 typedef struct { 28 int_t nprow,npcol,*row,*col; 29 gridinfo_t grid; 30 superlu_dist_options_t options; 31 SuperMatrix A_sup; 32 ScalePermstruct_t ScalePermstruct; 33 LUstruct_t LUstruct; 34 int StatPrint; 35 SuperLU_MatInputMode MatInputMode; 36 SOLVEstruct_t SOLVEstruct; 37 fact_t FactPattern; 38 MPI_Comm comm_superlu; 39 #if defined(PETSC_USE_COMPLEX) 40 doublecomplex *val; 41 #else 42 double *val; 43 #endif 44 PetscBool matsolve_iscalled,matmatsolve_iscalled; 45 PetscBool CleanUpSuperLU_Dist; /* Flag to clean up (non-global) SuperLU objects during Destroy */ 46 } Mat_SuperLU_DIST; 47 48 49 PetscErrorCode MatSuperluDistGetDiagU_SuperLU_DIST(Mat F,PetscScalar *diagU) 50 { 51 Mat_SuperLU_DIST *lu= (Mat_SuperLU_DIST*)F->data; 52 53 PetscFunctionBegin; 54 #if defined(PETSC_USE_COMPLEX) 55 PetscStackCall("SuperLU_DIST:pzGetDiagU",pzGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,(doublecomplex*)diagU)); 56 #else 57 PetscStackCall("SuperLU_DIST:pdGetDiagU",pdGetDiagU(F->rmap->N,&lu->LUstruct,&lu->grid,diagU)); 58 #endif 59 PetscFunctionReturn(0); 60 } 61 62 PetscErrorCode MatSuperluDistGetDiagU(Mat F,PetscScalar *diagU) 63 { 64 PetscErrorCode ierr; 65 66 PetscFunctionBegin; 67 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 68 ierr = PetscTryMethod(F,"MatSuperluDistGetDiagU_C",(Mat,PetscScalar*),(F,diagU));CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 static PetscErrorCode MatDestroy_SuperLU_DIST(Mat A) 73 { 74 PetscErrorCode ierr; 75 Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data; 76 77 PetscFunctionBegin; 78 if (lu->CleanUpSuperLU_Dist) { 79 /* Deallocate SuperLU_DIST storage */ 80 if (lu->MatInputMode == GLOBAL) { 81 PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup)); 82 } else { 83 PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup)); 84 if (lu->options.SolveInitialized) { 85 #if defined(PETSC_USE_COMPLEX) 86 PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct)); 87 #else 88 PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct)); 89 #endif 90 } 91 } 92 PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct)); 93 PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct)); 94 PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct)); 95 96 /* Release the SuperLU_DIST process grid. */ 97 PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid)); 98 ierr = MPI_Comm_free(&(lu->comm_superlu));CHKERRQ(ierr); 99 } 100 ierr = PetscFree(A->data);CHKERRQ(ierr); 101 /* clear composed functions */ 102 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr); 103 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluDistGetDiagU_C",NULL);CHKERRQ(ierr); 104 105 PetscFunctionReturn(0); 106 } 107 108 static PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x) 109 { 110 Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data; 111 PetscErrorCode ierr; 112 PetscMPIInt size; 113 PetscInt m=A->rmap->n,M=A->rmap->N,N=A->cmap->N; 114 SuperLUStat_t stat; 115 double berr[1]; 116 PetscScalar *bptr=NULL; 117 PetscInt nrhs=1; 118 Vec x_seq; 119 IS iden; 120 VecScatter scat; 121 int info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */ 122 static PetscBool cite = PETSC_FALSE; 123 124 PetscFunctionBegin; 125 if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED"); 126 ierr = PetscCitationsRegister("@article{lidemmel03,\n author = {Xiaoye S. Li and James W. Demmel},\n title = {{SuperLU_DIST}: A Scalable Distributed-Memory Sparse Direct\n Solver for Unsymmetric Linear Systems},\n journal = {ACM Trans. Mathematical Software},\n volume = {29},\n number = {2},\n pages = {110-140},\n year = 2003\n}\n",&cite);CHKERRQ(ierr); 127 128 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 129 if (size > 1 && lu->MatInputMode == GLOBAL) { 130 /* global mat input, convert b to x_seq */ 131 ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);CHKERRQ(ierr); 132 ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);CHKERRQ(ierr); 133 ierr = VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);CHKERRQ(ierr); 134 ierr = ISDestroy(&iden);CHKERRQ(ierr); 135 136 ierr = VecScatterBegin(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 137 ierr = VecScatterEnd(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 138 ierr = VecGetArray(x_seq,&bptr);CHKERRQ(ierr); 139 } else { /* size==1 || distributed mat input */ 140 if (lu->options.SolveInitialized && !lu->matsolve_iscalled) { 141 /* see comments in MatMatSolve() */ 142 #if defined(PETSC_USE_COMPLEX) 143 PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct)); 144 #else 145 PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct)); 146 #endif 147 lu->options.SolveInitialized = NO; 148 } 149 ierr = VecCopy(b_mpi,x);CHKERRQ(ierr); 150 ierr = VecGetArray(x,&bptr);CHKERRQ(ierr); 151 } 152 153 PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */ 154 if (lu->MatInputMode == GLOBAL) { 155 #if defined(PETSC_USE_COMPLEX) 156 PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info)); 157 #else 158 PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info)); 159 #endif 160 } else { /* distributed mat input */ 161 #if defined(PETSC_USE_COMPLEX) 162 PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info)); 163 #else 164 PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info)); 165 #endif 166 } 167 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info); 168 169 if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */ 170 PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat)); 171 172 if (size > 1 && lu->MatInputMode == GLOBAL) { 173 /* convert seq x to mpi x */ 174 ierr = VecRestoreArray(x_seq,&bptr);CHKERRQ(ierr); 175 ierr = VecScatterBegin(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 176 ierr = VecScatterEnd(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 177 ierr = VecScatterDestroy(&scat);CHKERRQ(ierr); 178 ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 179 } else { 180 ierr = VecRestoreArray(x,&bptr);CHKERRQ(ierr); 181 182 lu->matsolve_iscalled = PETSC_TRUE; 183 lu->matmatsolve_iscalled = PETSC_FALSE; 184 } 185 PetscFunctionReturn(0); 186 } 187 188 static PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X) 189 { 190 Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->data; 191 PetscErrorCode ierr; 192 PetscMPIInt size; 193 PetscInt M=A->rmap->N,m=A->rmap->n,nrhs; 194 SuperLUStat_t stat; 195 double berr[1]; 196 PetscScalar *bptr; 197 int info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */ 198 PetscBool flg; 199 200 PetscFunctionBegin; 201 if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED"); 202 ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 203 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 204 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 205 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 206 207 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 208 if (size > 1 && lu->MatInputMode == GLOBAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatInputMode=GLOBAL for nproc %d>1 is not supported",size); 209 /* size==1 or distributed mat input */ 210 if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) { 211 /* communication pattern of SOLVEstruct is unlikely created for matmatsolve, 212 thus destroy it and create a new SOLVEstruct. 213 Otherwise it may result in memory corruption or incorrect solution 214 See src/mat/examples/tests/ex125.c */ 215 #if defined(PETSC_USE_COMPLEX) 216 PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct)); 217 #else 218 PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct)); 219 #endif 220 lu->options.SolveInitialized = NO; 221 } 222 ierr = MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 223 224 ierr = MatGetSize(B_mpi,NULL,&nrhs);CHKERRQ(ierr); 225 226 PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */ 227 ierr = MatDenseGetArray(X,&bptr);CHKERRQ(ierr); 228 if (lu->MatInputMode == GLOBAL) { /* size == 1 */ 229 #if defined(PETSC_USE_COMPLEX) 230 PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,(doublecomplex*)bptr, M, nrhs,&lu->grid, &lu->LUstruct, berr, &stat, &info)); 231 #else 232 PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr, M, nrhs, &lu->grid, &lu->LUstruct, berr, &stat, &info)); 233 #endif 234 } else { /* distributed mat input */ 235 #if defined(PETSC_USE_COMPLEX) 236 PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid, &lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info)); 237 #else 238 PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info)); 239 #endif 240 } 241 if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info); 242 ierr = MatDenseRestoreArray(X,&bptr);CHKERRQ(ierr); 243 244 if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */ 245 PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat)); 246 lu->matsolve_iscalled = PETSC_FALSE; 247 lu->matmatsolve_iscalled = PETSC_TRUE; 248 PetscFunctionReturn(0); 249 } 250 251 #if !defined(PETSC_USE_COMPLEX) 252 /* 253 input: 254 F: numeric Cholesky factor 255 output: 256 nneg: total number of negative pivots 257 nzero: total number of zero pivots 258 npos: (global dimension of F) - nneg - nzero 259 */ 260 static PetscErrorCode MatGetInertia_SuperLU_DIST(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos) 261 { 262 PetscErrorCode ierr; 263 Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data; 264 PetscScalar *diagU=NULL; 265 PetscInt M,i,neg=0,zero=0,pos=0; 266 267 PetscFunctionBegin; 268 if (!F->assembled) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix factor F is not assembled"); 269 if (lu->options.RowPerm != NOROWPERM) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must set NOROWPERM"); 270 ierr = MatGetSize(F,&M,NULL);CHKERRQ(ierr); 271 ierr = PetscMalloc1(M,&diagU);CHKERRQ(ierr); 272 ierr = MatSuperluDistGetDiagU(F,diagU);CHKERRQ(ierr); 273 for (i=0; i<M; i++) { 274 if (diagU[i] > 0) { 275 pos++; 276 } else if (diagU[i] < 0) { 277 neg++; 278 } else zero++; 279 } 280 ierr = PetscFree(diagU);CHKERRQ(ierr); 281 *nneg = neg; 282 *nzero = zero; 283 *npos = pos; 284 PetscFunctionReturn(0); 285 } 286 #endif 287 288 static PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info) 289 { 290 Mat *tseq,A_seq = NULL; 291 Mat_SeqAIJ *aa,*bb; 292 Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data; 293 PetscErrorCode ierr; 294 PetscInt M=A->rmap->N,N=A->cmap->N,i,*ai,*aj,*bi,*bj,nz,rstart,*garray, 295 m=A->rmap->n, colA_start,j,jcol,jB,countA,countB,*bjj,*ajj=NULL; 296 int sinfo; /* SuperLU_Dist info flag is always an int even with long long indices */ 297 PetscMPIInt size; 298 SuperLUStat_t stat; 299 double *berr=0; 300 IS isrow; 301 #if defined(PETSC_USE_COMPLEX) 302 doublecomplex *av, *bv; 303 #else 304 double *av, *bv; 305 #endif 306 307 PetscFunctionBegin; 308 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 309 310 if (lu->MatInputMode == GLOBAL) { /* global mat input */ 311 if (size > 1) { /* convert mpi A to seq mat A */ 312 ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 313 ierr = MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 314 ierr = ISDestroy(&isrow);CHKERRQ(ierr); 315 316 A_seq = *tseq; 317 ierr = PetscFree(tseq);CHKERRQ(ierr); 318 aa = (Mat_SeqAIJ*)A_seq->data; 319 } else { 320 PetscBool flg; 321 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 322 if (flg) { 323 Mat_MPIAIJ *At = (Mat_MPIAIJ*)A->data; 324 A = At->A; 325 } 326 aa = (Mat_SeqAIJ*)A->data; 327 } 328 329 /* Convert Petsc NR matrix to SuperLU_DIST NC. 330 Note: memories of lu->val, col and row are allocated by CompRow_to_CompCol_dist()! */ 331 if (lu->options.Fact != DOFACT) {/* successive numeric factorization, sparsity pattern is reused. */ 332 PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup)); 333 if (lu->FactPattern == SamePattern_SameRowPerm) { 334 lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */ 335 } else { /* lu->FactPattern == SamePattern */ 336 PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct)); 337 lu->options.Fact = SamePattern; 338 } 339 } 340 #if defined(PETSC_USE_COMPLEX) 341 PetscStackCall("SuperLU_DIST:zCompRow_to_CompCol_dist",zCompRow_to_CompCol_dist(M,N,aa->nz,(doublecomplex*)aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val,&lu->col, &lu->row)); 342 #else 343 PetscStackCall("SuperLU_DIST:dCompRow_to_CompCol_dist",dCompRow_to_CompCol_dist(M,N,aa->nz,aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val, &lu->col, &lu->row)); 344 #endif 345 346 /* Create compressed column matrix A_sup. */ 347 #if defined(PETSC_USE_COMPLEX) 348 PetscStackCall("SuperLU_DIST:zCreate_CompCol_Matrix_dist",zCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_Z, SLU_GE)); 349 #else 350 PetscStackCall("SuperLU_DIST:dCreate_CompCol_Matrix_dist",dCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_D, SLU_GE)); 351 #endif 352 } else { /* distributed mat input */ 353 Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 354 aa=(Mat_SeqAIJ*)(mat->A)->data; 355 bb=(Mat_SeqAIJ*)(mat->B)->data; 356 ai=aa->i; aj=aa->j; 357 bi=bb->i; bj=bb->j; 358 #if defined(PETSC_USE_COMPLEX) 359 av=(doublecomplex*)aa->a; 360 bv=(doublecomplex*)bb->a; 361 #else 362 av=aa->a; 363 bv=bb->a; 364 #endif 365 rstart = A->rmap->rstart; 366 nz = aa->nz + bb->nz; 367 garray = mat->garray; 368 369 if (lu->options.Fact == DOFACT) { /* first numeric factorization */ 370 #if defined(PETSC_USE_COMPLEX) 371 PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row)); 372 #else 373 PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row)); 374 #endif 375 } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */ 376 if (lu->FactPattern == SamePattern_SameRowPerm) { 377 lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */ 378 } else if (lu->FactPattern == SamePattern) { 379 PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct)); /* Deallocate L and U matrices. */ 380 lu->options.Fact = SamePattern; 381 } else if (lu->FactPattern == DOFACT) { 382 PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup)); 383 PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct)); 384 lu->options.Fact = DOFACT; 385 386 #if defined(PETSC_USE_COMPLEX) 387 PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row)); 388 #else 389 PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row)); 390 #endif 391 } else { 392 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"options.Fact must be one of SamePattern SamePattern_SameRowPerm DOFACT"); 393 } 394 } 395 nz = 0; 396 for (i=0; i<m; i++) { 397 lu->row[i] = nz; 398 countA = ai[i+1] - ai[i]; 399 countB = bi[i+1] - bi[i]; 400 if (aj) { 401 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 402 } else { 403 ajj = NULL; 404 } 405 bjj = bj + bi[i]; 406 407 /* B part, smaller col index */ 408 if (aj) { 409 colA_start = rstart + ajj[0]; /* the smallest global col index of A */ 410 } else { /* superlu_dist does not require matrix has diagonal entries, thus aj=NULL would work */ 411 colA_start = rstart; 412 } 413 jB = 0; 414 for (j=0; j<countB; j++) { 415 jcol = garray[bjj[j]]; 416 if (jcol > colA_start) { 417 jB = j; 418 break; 419 } 420 lu->col[nz] = jcol; 421 lu->val[nz++] = *bv++; 422 if (j==countB-1) jB = countB; 423 } 424 425 /* A part */ 426 for (j=0; j<countA; j++) { 427 lu->col[nz] = rstart + ajj[j]; 428 lu->val[nz++] = *av++; 429 } 430 431 /* B part, larger col index */ 432 for (j=jB; j<countB; j++) { 433 lu->col[nz] = garray[bjj[j]]; 434 lu->val[nz++] = *bv++; 435 } 436 } 437 lu->row[m] = nz; 438 439 if (lu->options.Fact == DOFACT) { 440 #if defined(PETSC_USE_COMPLEX) 441 PetscStackCall("SuperLU_DIST:zCreate_CompRowLoc_Matrix_dist",zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE)); 442 #else 443 PetscStackCall("SuperLU_DIST:dCreate_CompRowLoc_Matrix_dist",dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE)); 444 #endif 445 } 446 } 447 448 /* Factor the matrix. */ 449 PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat)); /* Initialize the statistics variables. */ 450 if (lu->MatInputMode == GLOBAL) { /* global mat input */ 451 #if defined(PETSC_USE_COMPLEX) 452 PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo)); 453 #else 454 PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo)); 455 #endif 456 } else { /* distributed mat input */ 457 #if defined(PETSC_USE_COMPLEX) 458 PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo)); 459 #else 460 PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo)); 461 #endif 462 } 463 464 if (sinfo > 0) { 465 if (A->erroriffailure) { 466 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo); 467 } else { 468 if (sinfo <= lu->A_sup.ncol) { 469 F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 470 ierr = PetscInfo1(F,"U(i,i) is exactly zero, i= %D\n",sinfo);CHKERRQ(ierr); 471 } else if (sinfo > lu->A_sup.ncol) { 472 /* 473 number of bytes allocated when memory allocation 474 failure occurred, plus A->ncol. 475 */ 476 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 477 ierr = PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);CHKERRQ(ierr); 478 } 479 } 480 } else if (sinfo < 0) { 481 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, argument in p*gssvx() had an illegal value", sinfo); 482 } 483 484 if (lu->MatInputMode == GLOBAL && size > 1) { 485 ierr = MatDestroy(&A_seq);CHKERRQ(ierr); 486 } 487 488 if (lu->options.PrintStat) { 489 PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */ 490 } 491 PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat)); 492 F->assembled = PETSC_TRUE; 493 F->preallocated = PETSC_TRUE; 494 lu->options.Fact = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */ 495 PetscFunctionReturn(0); 496 } 497 498 /* Note the Petsc r and c permutations are ignored */ 499 static PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 500 { 501 Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->data; 502 PetscInt M = A->rmap->N,N=A->cmap->N; 503 504 PetscFunctionBegin; 505 /* Initialize the SuperLU process grid. */ 506 PetscStackCall("SuperLU_DIST:superlu_gridinit",superlu_gridinit(lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid)); 507 508 /* Initialize ScalePermstruct and LUstruct. */ 509 PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct)); 510 PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct)); 511 F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST; 512 F->ops->solve = MatSolve_SuperLU_DIST; 513 F->ops->matsolve = MatMatSolve_SuperLU_DIST; 514 F->ops->getinertia = NULL; 515 #if !defined(PETSC_USE_COMPLEX) 516 if (A->symmetric) { 517 F->ops->getinertia = MatGetInertia_SuperLU_DIST; 518 } 519 #endif 520 lu->CleanUpSuperLU_Dist = PETSC_TRUE; 521 PetscFunctionReturn(0); 522 } 523 524 static PetscErrorCode MatCholeskyFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,const MatFactorInfo *info) 525 { 526 PetscErrorCode ierr; 527 528 PetscFunctionBegin; 529 if (!A->symmetric) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Input matrix must be symmetric\n"); 530 ierr = MatLUFactorSymbolic_SuperLU_DIST(F,A,r,r,info);CHKERRQ(ierr); 531 F->ops->choleskyfactornumeric = MatLUFactorNumeric_SuperLU_DIST; 532 PetscFunctionReturn(0); 533 } 534 535 static PetscErrorCode MatFactorGetSolverType_aij_superlu_dist(Mat A,MatSolverType *type) 536 { 537 PetscFunctionBegin; 538 *type = MATSOLVERSUPERLU_DIST; 539 PetscFunctionReturn(0); 540 } 541 542 static PetscErrorCode MatView_Info_SuperLU_DIST(Mat A,PetscViewer viewer) 543 { 544 Mat_SuperLU_DIST *lu=(Mat_SuperLU_DIST*)A->data; 545 superlu_dist_options_t options; 546 PetscErrorCode ierr; 547 548 PetscFunctionBegin; 549 /* check if matrix is superlu_dist type */ 550 if (A->ops->solve != MatSolve_SuperLU_DIST) PetscFunctionReturn(0); 551 552 options = lu->options; 553 ierr = PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");CHKERRQ(ierr); 554 ierr = PetscViewerASCIIPrintf(viewer," Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);CHKERRQ(ierr); 555 ierr = PetscViewerASCIIPrintf(viewer," Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);CHKERRQ(ierr); 556 ierr = PetscViewerASCIIPrintf(viewer," Matrix input mode %d \n",lu->MatInputMode);CHKERRQ(ierr); 557 ierr = PetscViewerASCIIPrintf(viewer," Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);CHKERRQ(ierr); 558 ierr = PetscViewerASCIIPrintf(viewer," Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);CHKERRQ(ierr); 559 ierr = PetscViewerASCIIPrintf(viewer," Processors in row %d col partition %d \n",lu->nprow,lu->npcol);CHKERRQ(ierr); 560 ierr = PetscViewerASCIIPrintf(viewer," Row permutation %s \n",(options.RowPerm == NOROWPERM) ? "NATURAL" : "LargeDiag");CHKERRQ(ierr); 561 562 switch (options.ColPerm) { 563 case NATURAL: 564 ierr = PetscViewerASCIIPrintf(viewer," Column permutation NATURAL\n");CHKERRQ(ierr); 565 break; 566 case MMD_AT_PLUS_A: 567 ierr = PetscViewerASCIIPrintf(viewer," Column permutation MMD_AT_PLUS_A\n");CHKERRQ(ierr); 568 break; 569 case MMD_ATA: 570 ierr = PetscViewerASCIIPrintf(viewer," Column permutation MMD_ATA\n");CHKERRQ(ierr); 571 break; 572 case METIS_AT_PLUS_A: 573 ierr = PetscViewerASCIIPrintf(viewer," Column permutation METIS_AT_PLUS_A\n");CHKERRQ(ierr); 574 break; 575 case PARMETIS: 576 ierr = PetscViewerASCIIPrintf(viewer," Column permutation PARMETIS\n");CHKERRQ(ierr); 577 break; 578 default: 579 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation"); 580 } 581 582 ierr = PetscViewerASCIIPrintf(viewer," Parallel symbolic factorization %s \n",PetscBools[options.ParSymbFact != NO]);CHKERRQ(ierr); 583 584 if (lu->FactPattern == SamePattern) { 585 ierr = PetscViewerASCIIPrintf(viewer," Repeated factorization SamePattern\n");CHKERRQ(ierr); 586 } else if (lu->FactPattern == SamePattern_SameRowPerm) { 587 ierr = PetscViewerASCIIPrintf(viewer," Repeated factorization SamePattern_SameRowPerm\n");CHKERRQ(ierr); 588 } else if (lu->FactPattern == DOFACT) { 589 ierr = PetscViewerASCIIPrintf(viewer," Repeated factorization DOFACT\n");CHKERRQ(ierr); 590 } else { 591 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown factorization pattern"); 592 } 593 PetscFunctionReturn(0); 594 } 595 596 static PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer) 597 { 598 PetscErrorCode ierr; 599 PetscBool iascii; 600 PetscViewerFormat format; 601 602 PetscFunctionBegin; 603 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 604 if (iascii) { 605 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 606 if (format == PETSC_VIEWER_ASCII_INFO) { 607 ierr = MatView_Info_SuperLU_DIST(A,viewer);CHKERRQ(ierr); 608 } 609 } 610 PetscFunctionReturn(0); 611 } 612 613 static PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F) 614 { 615 Mat B; 616 Mat_SuperLU_DIST *lu; 617 PetscErrorCode ierr; 618 PetscInt M=A->rmap->N,N=A->cmap->N,indx; 619 PetscMPIInt size; 620 superlu_dist_options_t options; 621 PetscBool flg; 622 const char *colperm[] = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"}; 623 const char *rowperm[] = {"LargeDiag","NATURAL"}; 624 const char *factPattern[] = {"SamePattern","SamePattern_SameRowPerm","DOFACT"}; 625 PetscBool set; 626 627 PetscFunctionBegin; 628 /* Create the factorization matrix */ 629 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 630 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 631 ierr = PetscStrallocpy("superlu_dist",&((PetscObject)B)->type_name);CHKERRQ(ierr); 632 ierr = MatSetUp(B);CHKERRQ(ierr); 633 B->ops->getinfo = MatGetInfo_External; 634 B->ops->view = MatView_SuperLU_DIST; 635 B->ops->destroy = MatDestroy_SuperLU_DIST; 636 637 if (ftype == MAT_FACTOR_LU) { 638 B->factortype = MAT_FACTOR_LU; 639 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST; 640 } else { 641 B->factortype = MAT_FACTOR_CHOLESKY; 642 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SuperLU_DIST; 643 } 644 645 /* set solvertype */ 646 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 647 ierr = PetscStrallocpy(MATSOLVERSUPERLU_DIST,&B->solvertype);CHKERRQ(ierr); 648 649 ierr = PetscNewLog(B,&lu);CHKERRQ(ierr); 650 B->data = lu; 651 652 /* Set the default input options: 653 options.Fact = DOFACT; 654 options.Equil = YES; 655 options.ParSymbFact = NO; 656 options.ColPerm = METIS_AT_PLUS_A; 657 options.RowPerm = LargeDiag; 658 options.ReplaceTinyPivot = YES; 659 options.IterRefine = DOUBLE; 660 options.Trans = NOTRANS; 661 options.SolveInitialized = NO; -hold the communication pattern used MatSolve() and MatMatSolve() 662 options.RefineInitialized = NO; 663 options.PrintStat = YES; 664 */ 665 set_default_options_dist(&options); 666 667 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(lu->comm_superlu));CHKERRQ(ierr); 668 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 669 /* Default num of process columns and rows */ 670 lu->nprow = (int_t) (0.5 + PetscSqrtReal((PetscReal)size)); 671 if (!lu->nprow) lu->nprow = 1; 672 while (lu->nprow > 0) { 673 lu->npcol = (int_t) (size/lu->nprow); 674 if (size == lu->nprow * lu->npcol) break; 675 lu->nprow--; 676 } 677 678 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");CHKERRQ(ierr); 679 ierr = PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,(PetscInt*)&lu->nprow,NULL);CHKERRQ(ierr); 680 ierr = PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,(PetscInt*)&lu->npcol,NULL);CHKERRQ(ierr); 681 if (size != lu->nprow * lu->npcol) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of processes %d must equal to nprow %d * npcol %d",size,lu->nprow,lu->npcol); 682 683 lu->MatInputMode = DISTRIBUTED; 684 685 ierr = PetscOptionsEnum("-mat_superlu_dist_matinput","Matrix input mode (global or distributed)","None",SuperLU_MatInputModes,(PetscEnum)lu->MatInputMode,(PetscEnum*)&lu->MatInputMode,NULL);CHKERRQ(ierr); 686 if (lu->MatInputMode == DISTRIBUTED && size == 1) lu->MatInputMode = GLOBAL; 687 688 ierr = PetscOptionsBool("-mat_superlu_dist_equil","Equilibrate matrix","None",options.Equil ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 689 if (set && !flg) options.Equil = NO; 690 691 ierr = PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 692 if (flg) { 693 switch (indx) { 694 case 0: 695 options.RowPerm = LargeDiag; 696 break; 697 case 1: 698 options.RowPerm = NOROWPERM; 699 break; 700 } 701 } 702 703 ierr = PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",colperm,5,colperm[3],&indx,&flg);CHKERRQ(ierr); 704 if (flg) { 705 switch (indx) { 706 case 0: 707 options.ColPerm = NATURAL; 708 break; 709 case 1: 710 options.ColPerm = MMD_AT_PLUS_A; 711 break; 712 case 2: 713 options.ColPerm = MMD_ATA; 714 break; 715 case 3: 716 options.ColPerm = METIS_AT_PLUS_A; 717 break; 718 case 4: 719 options.ColPerm = PARMETIS; /* only works for np>1 */ 720 break; 721 default: 722 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation"); 723 } 724 } 725 726 options.ReplaceTinyPivot = NO; 727 ierr = PetscOptionsBool("-mat_superlu_dist_replacetinypivot","Replace tiny pivots","None",options.ReplaceTinyPivot ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 728 if (set && flg) options.ReplaceTinyPivot = YES; 729 730 options.ParSymbFact = NO; 731 ierr = PetscOptionsBool("-mat_superlu_dist_parsymbfact","Parallel symbolic factorization","None",PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 732 if (set && flg && size>1) { 733 if (lu->MatInputMode == GLOBAL) { 734 #if defined(PETSC_USE_INFO) 735 ierr = PetscInfo(A,"Warning: '-mat_superlu_dist_parsymbfact' is ignored because MatInputMode=GLOBAL\n");CHKERRQ(ierr); 736 #endif 737 } else { 738 #if defined(PETSC_HAVE_PARMETIS) 739 options.ParSymbFact = YES; 740 options.ColPerm = PARMETIS; /* in v2.2, PARMETIS is forced for ParSymbFact regardless of user ordering setting */ 741 #else 742 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"parsymbfact needs PARMETIS"); 743 #endif 744 } 745 } 746 747 lu->FactPattern = SamePattern; 748 ierr = PetscOptionsEList("-mat_superlu_dist_fact","Sparsity pattern for repeated matrix factorization","None",factPattern,3,factPattern[0],&indx,&flg);CHKERRQ(ierr); 749 if (flg) { 750 switch (indx) { 751 case 0: 752 lu->FactPattern = SamePattern; 753 break; 754 case 1: 755 lu->FactPattern = SamePattern_SameRowPerm; 756 break; 757 case 2: 758 lu->FactPattern = DOFACT; 759 break; 760 } 761 } 762 763 options.IterRefine = NOREFINE; 764 ierr = PetscOptionsBool("-mat_superlu_dist_iterrefine","Use iterative refinement","None",options.IterRefine == NOREFINE ? PETSC_FALSE : PETSC_TRUE ,&flg,&set);CHKERRQ(ierr); 765 if (set) { 766 if (flg) options.IterRefine = SLU_DOUBLE; 767 else options.IterRefine = NOREFINE; 768 } 769 770 if (PetscLogPrintInfo) options.PrintStat = YES; 771 else options.PrintStat = NO; 772 ierr = PetscOptionsBool("-mat_superlu_dist_statprint","Print factorization information","None",(PetscBool)options.PrintStat,(PetscBool*)&options.PrintStat,NULL);CHKERRQ(ierr); 773 ierr = PetscOptionsEnd();CHKERRQ(ierr); 774 775 lu->options = options; 776 lu->options.Fact = DOFACT; 777 lu->matsolve_iscalled = PETSC_FALSE; 778 lu->matmatsolve_iscalled = PETSC_FALSE; 779 780 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_aij_superlu_dist);CHKERRQ(ierr); 781 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSuperluDistGetDiagU_C",MatSuperluDistGetDiagU_SuperLU_DIST);CHKERRQ(ierr); 782 783 *F = B; 784 PetscFunctionReturn(0); 785 } 786 787 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_SuperLU_DIST(void) 788 { 789 PetscErrorCode ierr; 790 PetscFunctionBegin; 791 ierr = MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ, MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr); 792 ierr = MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ, MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr); 793 ierr = MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr); 794 ierr = MatSolverTypeRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr); 795 PetscFunctionReturn(0); 796 } 797 798 /*MC 799 MATSOLVERSUPERLU_DIST - Parallel direct solver package for LU factorization 800 801 Use ./configure --download-superlu_dist --download-parmetis --download-metis --download-ptscotch to have PETSc installed with SuperLU_DIST 802 803 Use -pc_type lu -pc_factor_mat_solver_type superlu_dist to use this direct solver 804 805 Works with AIJ matrices 806 807 Options Database Keys: 808 + -mat_superlu_dist_r <n> - number of rows in processor partition 809 . -mat_superlu_dist_c <n> - number of columns in processor partition 810 . -mat_superlu_dist_matinput <0,1> - matrix input mode; 0=global, 1=distributed 811 . -mat_superlu_dist_equil - equilibrate the matrix 812 . -mat_superlu_dist_rowperm <LargeDiag,NATURAL> - row permutation 813 . -mat_superlu_dist_colperm <MMD_AT_PLUS_A,MMD_ATA,NATURAL> - column permutation 814 . -mat_superlu_dist_replacetinypivot - replace tiny pivots 815 . -mat_superlu_dist_fact <SamePattern> - (choose one of) SamePattern SamePattern_SameRowPerm DOFACT 816 . -mat_superlu_dist_iterrefine - use iterative refinement 817 - -mat_superlu_dist_statprint - print factorization information 818 819 Level: beginner 820 821 .seealso: PCLU 822 823 .seealso: PCFactorSetMatSolverType(), MatSolverType 824 825 M*/ 826