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