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