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