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