1 2 /* 3 Provides an interface to the SuperLU 3.0 sparse solver 4 */ 5 6 #include "src/mat/impls/aij/seq/aij.h" 7 8 EXTERN_C_BEGIN 9 #if defined(PETSC_USE_COMPLEX) 10 #include "zsp_defs.h" 11 #else 12 #include "dsp_defs.h" 13 #endif 14 #include "util.h" 15 EXTERN_C_END 16 17 typedef struct { 18 SuperMatrix A,L,U,B,X; 19 superlu_options_t options; 20 int *perm_c; /* column permutation vector */ 21 int *perm_r; /* row permutations from partial pivoting */ 22 int *etree; 23 double *R, *C; 24 char equed[1]; 25 int lwork; 26 void *work; 27 double rpg, rcond; 28 mem_usage_t mem_usage; 29 MatStructure flg; 30 31 /* A few function pointers for inheritance */ 32 int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 33 int (*MatView)(Mat,PetscViewer); 34 int (*MatAssemblyEnd)(Mat,MatAssemblyType); 35 int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 36 int (*MatDestroy)(Mat); 37 38 /* Flag to clean up (non-global) SuperLU objects during Destroy */ 39 PetscTruth CleanUpSuperLU; 40 } Mat_SuperLU; 41 42 43 EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer); 44 EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*); 45 46 EXTERN_C_BEGIN 47 EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,const MatType,Mat*); 48 EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,const MatType,Mat*); 49 EXTERN_C_END 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "MatDestroy_SuperLU" 53 int MatDestroy_SuperLU(Mat A) 54 { 55 int ierr; 56 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 57 58 PetscFunctionBegin; 59 if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */ 60 Destroy_SuperMatrix_Store(&lu->A); 61 Destroy_SuperMatrix_Store(&lu->B); 62 Destroy_SuperMatrix_Store(&lu->X); 63 64 ierr = PetscFree(lu->etree);CHKERRQ(ierr); 65 ierr = PetscFree(lu->perm_r);CHKERRQ(ierr); 66 ierr = PetscFree(lu->perm_c);CHKERRQ(ierr); 67 ierr = PetscFree(lu->R);CHKERRQ(ierr); 68 ierr = PetscFree(lu->C);CHKERRQ(ierr); 69 if ( lu->lwork >= 0 ) { 70 Destroy_SuperNode_Matrix(&lu->L); 71 Destroy_CompCol_Matrix(&lu->U); 72 } 73 } 74 ierr = MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);CHKERRQ(ierr); 75 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "MatView_SuperLU" 81 int MatView_SuperLU(Mat A,PetscViewer viewer) 82 { 83 int ierr; 84 PetscTruth iascii; 85 PetscViewerFormat format; 86 Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 87 88 PetscFunctionBegin; 89 ierr = (*lu->MatView)(A,viewer);CHKERRQ(ierr); 90 91 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 92 if (iascii) { 93 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 94 if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 95 ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr); 96 } 97 } 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "MatAssemblyEnd_SuperLU" 103 int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) { 104 int ierr; 105 Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr); 106 107 PetscFunctionBegin; 108 ierr = (*lu->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 109 110 lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 111 A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 112 PetscFunctionReturn(0); 113 } 114 115 /* This function was written for SuperLU 2.0 by Matthew Knepley. Not tested for SuperLU 3.0! */ 116 #ifdef SuperLU2 117 #include "src/mat/impls/dense/seq/dense.h" 118 #undef __FUNCT__ 119 #define __FUNCT__ "MatCreateNull_SuperLU" 120 int MatCreateNull_SuperLU(Mat A,Mat *nullMat) 121 { 122 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 123 int numRows = A->m,numCols = A->n; 124 SCformat *Lstore; 125 int numNullCols,size; 126 SuperLUStat_t stat; 127 #if defined(PETSC_USE_COMPLEX) 128 doublecomplex *nullVals,*workVals; 129 #else 130 PetscScalar *nullVals,*workVals; 131 #endif 132 int row,newRow,col,newCol,block,b,ierr; 133 134 PetscFunctionBegin; 135 if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix"); 136 numNullCols = numCols - numRows; 137 if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems"); 138 /* Create the null matrix using MATSEQDENSE explicitly */ 139 ierr = MatCreate(A->comm,numRows,numNullCols,numRows,numNullCols,nullMat);CHKERRQ(ierr); 140 ierr = MatSetType(*nullMat,MATSEQDENSE);CHKERRQ(ierr); 141 ierr = MatSeqDenseSetPreallocation(*nullMat,PETSC_NULL);CHKERRQ(ierr); 142 if (!numNullCols) { 143 ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 144 ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 #if defined(PETSC_USE_COMPLEX) 148 nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 149 #else 150 nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 151 #endif 152 /* Copy in the columns */ 153 Lstore = (SCformat*)lu->L.Store; 154 for(block = 0; block <= Lstore->nsuper; block++) { 155 newRow = Lstore->sup_to_col[block]; 156 size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 157 for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 158 newCol = Lstore->rowind[col]; 159 if (newCol >= numRows) { 160 for(b = 0; b < size; b++) 161 #if defined(PETSC_USE_COMPLEX) 162 nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 163 #else 164 nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 165 #endif 166 } 167 } 168 } 169 /* Permute rhs to form P^T_c B */ 170 ierr = PetscMalloc(numRows*sizeof(double),&workVals);CHKERRQ(ierr); 171 for(b = 0; b < numNullCols; b++) { 172 for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 173 for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 174 } 175 /* Backward solve the upper triangle A x = b */ 176 for(b = 0; b < numNullCols; b++) { 177 #if defined(PETSC_USE_COMPLEX) 178 sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 179 #else 180 sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 181 #endif 182 if (ierr < 0) 183 SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr); 184 } 185 ierr = PetscFree(workVals);CHKERRQ(ierr); 186 187 ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 188 ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 189 PetscFunctionReturn(0); 190 } 191 #endif 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "MatSolve_SuperLU" 195 int MatSolve_SuperLU(Mat A,Vec b,Vec x) 196 { 197 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 198 PetscScalar *barray,*xarray; 199 int ierr,info,i; 200 SuperLUStat_t stat; 201 double ferr,berr; 202 203 PetscFunctionBegin; 204 if ( lu->lwork == -1 ) { 205 PetscFunctionReturn(0); 206 } 207 lu->B.ncol = 1; /* Set the number of right-hand side */ 208 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 209 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 210 211 #if defined(PETSC_USE_COMPLEX) 212 ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 213 ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 214 #else 215 ((DNformat*)lu->B.Store)->nzval = barray; 216 ((DNformat*)lu->X.Store)->nzval = xarray; 217 #endif 218 219 /* Initialize the statistics variables. */ 220 StatInit(&stat); 221 222 lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 223 lu->options.Trans = TRANS; 224 #if defined(PETSC_USE_COMPLEX) 225 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 226 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 227 &lu->mem_usage, &stat, &info); 228 #else 229 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 230 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 231 &lu->mem_usage, &stat, &info); 232 #endif 233 ierr = VecRestoreArray(b,&barray);CHKERRQ(ierr); 234 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 235 236 if ( !info || info == lu->A.ncol+1 ) { 237 if ( lu->options.IterRefine ) { 238 ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n"); 239 ierr = PetscPrintf(PETSC_COMM_SELF," %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR"); 240 for (i = 0; i < 1; ++i) 241 ierr = PetscPrintf(PETSC_COMM_SELF," %8d%8d%16e%16e\n", i+1, stat.RefineSteps, ferr, berr); 242 } 243 } else if ( info > 0 ){ 244 if ( lu->lwork == -1 ) { 245 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %d bytes\n", info - lu->A.ncol); 246 } else { 247 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %d\n",info); 248 } 249 } else if (info < 0){ 250 SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info); 251 } 252 253 if ( lu->options.PrintStat ) { 254 ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n"); 255 StatPrint(&stat); 256 } 257 StatFree(&stat); 258 PetscFunctionReturn(0); 259 } 260 261 #undef __FUNCT__ 262 #define __FUNCT__ "MatLUFactorNumeric_SuperLU" 263 int MatLUFactorNumeric_SuperLU(Mat A,Mat *F) 264 { 265 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 266 Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 267 int ierr,info; 268 SuperLUStat_t stat; 269 double ferr, berr; 270 NCformat *Ustore; 271 SCformat *Lstore; 272 273 PetscFunctionBegin; 274 if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 275 lu->options.Fact = SamePattern; 276 /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 277 Destroy_SuperMatrix_Store(&lu->A); 278 if ( lu->lwork >= 0 ) { 279 Destroy_SuperNode_Matrix(&lu->L); 280 Destroy_CompCol_Matrix(&lu->U); 281 lu->options.Fact = SamePattern; 282 } 283 } 284 285 /* Create the SuperMatrix for lu->A=A^T: 286 Since SuperLU likes column-oriented matrices,we pass it the transpose, 287 and then solve A^T X = B in MatSolve(). */ 288 #if defined(PETSC_USE_COMPLEX) 289 zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 290 SLU_NC,SLU_Z,SLU_GE); 291 #else 292 dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 293 SLU_NC,SLU_D,SLU_GE); 294 #endif 295 296 /* Initialize the statistics variables. */ 297 StatInit(&stat); 298 299 /* Numerical factorization */ 300 lu->B.ncol = 0; /* Indicate not to solve the system */ 301 #if defined(PETSC_USE_COMPLEX) 302 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 303 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 304 &lu->mem_usage, &stat, &info); 305 #else 306 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 307 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 308 &lu->mem_usage, &stat, &info); 309 #endif 310 if ( !info || info == lu->A.ncol+1 ) { 311 if ( lu->options.PivotGrowth ) 312 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 313 if ( lu->options.ConditionNumber ) 314 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 315 } else if ( info > 0 ){ 316 if ( lu->lwork == -1 ) { 317 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %d bytes\n", info - lu->A.ncol); 318 } else { 319 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %d\n",info); 320 } 321 } else { /* info < 0 */ 322 SETERRQ2(1, "info = %d, the %d-th argument in gssvx() had an illegal value", info,-info); 323 } 324 325 if ( lu->options.PrintStat ) { 326 ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 327 StatPrint(&stat); 328 Lstore = (SCformat *) lu->L.Store; 329 Ustore = (NCformat *) lu->U.Store; 330 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %d\n", Lstore->nnz); 331 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %d\n", Ustore->nnz); 332 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 333 ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n", 334 lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 335 lu->mem_usage.expansions); 336 } 337 StatFree(&stat); 338 339 lu->flg = SAME_NONZERO_PATTERN; 340 PetscFunctionReturn(0); 341 } 342 343 /* 344 Note the r permutation is ignored 345 */ 346 #undef __FUNCT__ 347 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 348 int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 349 { 350 Mat B; 351 Mat_SuperLU *lu; 352 int ierr,m=A->m,n=A->n,indx; 353 PetscTruth flg; 354 const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 355 const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 356 const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 357 358 PetscFunctionBegin; 359 360 ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 361 ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 362 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 363 364 B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 365 B->ops->solve = MatSolve_SuperLU; 366 B->factor = FACTOR_LU; 367 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 368 369 lu = (Mat_SuperLU*)(B->spptr); 370 371 /* Set SuperLU options */ 372 /* the default values for options argument: 373 options.Fact = DOFACT; 374 options.Equil = YES; 375 options.ColPerm = COLAMD; 376 options.DiagPivotThresh = 1.0; 377 options.Trans = NOTRANS; 378 options.IterRefine = NOREFINE; 379 options.SymmetricMode = NO; 380 options.PivotGrowth = NO; 381 options.ConditionNumber = NO; 382 options.PrintStat = YES; 383 */ 384 set_default_options(&lu->options); 385 /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 386 lu->options.Equil = NO; 387 lu->options.PrintStat = NO; 388 lu->lwork = 0; /* allocate space internally by system malloc */ 389 390 ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 391 /* 392 ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 393 if (flg) lu->options.Equil = YES; -- not supported by the interface !!! 394 */ 395 ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 396 if (flg) {lu->options.ColPerm = (colperm_t)indx;} 397 ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 398 if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 399 ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 400 if (flg) lu->options.SymmetricMode = YES; 401 ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 402 ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 403 if (flg) lu->options.PivotGrowth = YES; 404 ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 405 if (flg) lu->options.ConditionNumber = YES; 406 ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 407 if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 408 ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 409 if (flg) lu->options.ReplaceTinyPivot = YES; 410 ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 411 if (flg) lu->options.PrintStat = YES; 412 ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 413 if (lu->lwork > 0 ){ 414 ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 415 } else if (lu->lwork != 0 && lu->lwork != -1){ 416 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %d is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 417 lu->lwork = 0; 418 } 419 PetscOptionsEnd(); 420 421 #ifdef SUPERLU2 422 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 423 (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 424 #endif 425 426 /* Allocate spaces (notice sizes are for the transpose) */ 427 ierr = PetscMalloc(m*sizeof(int),&lu->etree);CHKERRQ(ierr); 428 ierr = PetscMalloc(n*sizeof(int),&lu->perm_r);CHKERRQ(ierr); 429 ierr = PetscMalloc(m*sizeof(int),&lu->perm_c);CHKERRQ(ierr); 430 ierr = PetscMalloc(n*sizeof(int),&lu->R);CHKERRQ(ierr); 431 ierr = PetscMalloc(m*sizeof(int),&lu->C);CHKERRQ(ierr); 432 433 /* create rhs and solution x without allocate space for .Store */ 434 #if defined(PETSC_USE_COMPLEX) 435 zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 436 zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 437 #else 438 dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 439 dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 440 #endif 441 442 lu->flg = DIFFERENT_NONZERO_PATTERN; 443 lu->CleanUpSuperLU = PETSC_TRUE; 444 445 *F = B; 446 PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU)); 447 PetscFunctionReturn(0); 448 } 449 450 /* used by -ksp_view */ 451 #undef __FUNCT__ 452 #define __FUNCT__ "MatFactorInfo_SuperLU" 453 int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 454 { 455 Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 456 int ierr; 457 superlu_options_t options; 458 459 PetscFunctionBegin; 460 /* check if matrix is superlu_dist type */ 461 if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 462 463 options = lu->options; 464 ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 465 ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 466 ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %d\n",options.ColPerm);CHKERRQ(ierr); 467 ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %d\n",options.IterRefine);CHKERRQ(ierr); 468 ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 469 ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 470 ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 471 ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 472 ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %d\n",options.RowPerm);CHKERRQ(ierr); 473 ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 474 ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 475 ierr = PetscViewerASCIIPrintf(viewer," lwork: %d\n",lu->lwork);CHKERRQ(ierr); 476 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "MatDuplicate_SuperLU" 482 int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 483 int ierr; 484 Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 485 486 PetscFunctionBegin; 487 ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 488 ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 EXTERN_C_BEGIN 493 #undef __FUNCT__ 494 #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 495 int MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,Mat *newmat) { 496 /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 497 /* to its base PETSc type, so we will ignore 'MatType type'. */ 498 int ierr; 499 Mat B=*newmat; 500 Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 501 502 PetscFunctionBegin; 503 if (B != A) { 504 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 505 } 506 /* Reset the original function pointers */ 507 B->ops->duplicate = lu->MatDuplicate; 508 B->ops->view = lu->MatView; 509 B->ops->assemblyend = lu->MatAssemblyEnd; 510 B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 511 B->ops->destroy = lu->MatDestroy; 512 /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 513 ierr = PetscFree(lu);CHKERRQ(ierr); 514 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 515 *newmat = B; 516 PetscFunctionReturn(0); 517 } 518 EXTERN_C_END 519 520 EXTERN_C_BEGIN 521 #undef __FUNCT__ 522 #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 523 int MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,Mat *newmat) { 524 /* This routine is only called to convert to MATSUPERLU */ 525 /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 526 int ierr; 527 Mat B=*newmat; 528 Mat_SuperLU *lu; 529 530 PetscFunctionBegin; 531 if (B != A) { 532 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 533 } 534 535 ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 536 lu->MatDuplicate = A->ops->duplicate; 537 lu->MatView = A->ops->view; 538 lu->MatAssemblyEnd = A->ops->assemblyend; 539 lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 540 lu->MatDestroy = A->ops->destroy; 541 lu->CleanUpSuperLU = PETSC_FALSE; 542 543 B->spptr = (void*)lu; 544 B->ops->duplicate = MatDuplicate_SuperLU; 545 B->ops->view = MatView_SuperLU; 546 B->ops->assemblyend = MatAssemblyEnd_SuperLU; 547 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 548 B->ops->choleskyfactorsymbolic = 0; 549 B->ops->destroy = MatDestroy_SuperLU; 550 551 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 552 "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 553 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 554 "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 555 PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves."); 556 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 557 *newmat = B; 558 PetscFunctionReturn(0); 559 } 560 EXTERN_C_END 561 562 /*MC 563 MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 564 via the external package SuperLU. 565 566 If SuperLU is installed (see the manual for 567 instructions on how to declare the existence of external packages), 568 a matrix type can be constructed which invokes SuperLU solvers. 569 After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 570 This matrix type is only supported for double precision real. 571 572 This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 573 supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 574 the MATSEQAIJ type without data copy. 575 576 Options Database Keys: 577 + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 578 - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 579 1: MMD applied to A'*A, 580 2: MMD applied to A'+A, 581 3: COLAMD, approximate minimum degree column ordering 582 583 Level: beginner 584 585 .seealso: PCLU 586 M*/ 587 588 EXTERN_C_BEGIN 589 #undef __FUNCT__ 590 #define __FUNCT__ "MatCreate_SuperLU" 591 int MatCreate_SuperLU(Mat A) { 592 int ierr; 593 594 PetscFunctionBegin; 595 /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 596 ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 597 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 598 ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);CHKERRQ(ierr); 599 PetscFunctionReturn(0); 600 } 601 EXTERN_C_END 602