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 "zsp_defs.h" 12 #else 13 #include "dsp_defs.h" 14 #endif 15 #include "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,const MatType,MatReuse,Mat*); 49 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat,const 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,numRows,numNullCols,numRows,numNullCols,nullMat);CHKERRQ(ierr); 141 ierr = MatSetType(*nullMat,MATSEQDENSE);CHKERRQ(ierr); 142 ierr = MatSeqDenseSetPreallocation(*nullMat,PETSC_NULL);CHKERRQ(ierr); 143 if (!numNullCols) { 144 ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 145 ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146 PetscFunctionReturn(0); 147 } 148 #if defined(PETSC_USE_COMPLEX) 149 nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v; 150 #else 151 nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v; 152 #endif 153 /* Copy in the columns */ 154 Lstore = (SCformat*)lu->L.Store; 155 for(block = 0; block <= Lstore->nsuper; block++) { 156 newRow = Lstore->sup_to_col[block]; 157 size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block]; 158 for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) { 159 newCol = Lstore->rowind[col]; 160 if (newCol >= numRows) { 161 for(b = 0; b < size; b++) 162 #if defined(PETSC_USE_COMPLEX) 163 nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 164 #else 165 nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col]; 166 #endif 167 } 168 } 169 } 170 /* Permute rhs to form P^T_c B */ 171 ierr = PetscMalloc(numRows*sizeof(PetscReal),&workVals);CHKERRQ(ierr); 172 for(b = 0; b < numNullCols; b++) { 173 for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row]; 174 for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row]; 175 } 176 /* Backward solve the upper triangle A x = b */ 177 for(b = 0; b < numNullCols; b++) { 178 #if defined(PETSC_USE_COMPLEX) 179 sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 180 #else 181 sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&stat,&ierr); 182 #endif 183 if (ierr < 0) 184 SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %D was invalid",-ierr); 185 } 186 ierr = PetscFree(workVals);CHKERRQ(ierr); 187 188 ierr = MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 189 ierr = MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 190 PetscFunctionReturn(0); 191 } 192 #endif 193 194 #undef __FUNCT__ 195 #define __FUNCT__ "MatSolve_SuperLU" 196 PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x) 197 { 198 Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr; 199 PetscScalar *barray,*xarray; 200 PetscErrorCode ierr; 201 PetscInt info,i; 202 SuperLUStat_t stat; 203 PetscReal ferr,berr; 204 205 PetscFunctionBegin; 206 if ( lu->lwork == -1 ) { 207 PetscFunctionReturn(0); 208 } 209 lu->B.ncol = 1; /* Set the number of right-hand side */ 210 ierr = VecGetArray(b,&barray);CHKERRQ(ierr); 211 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 212 213 #if defined(PETSC_USE_COMPLEX) 214 ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray; 215 ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray; 216 #else 217 ((DNformat*)lu->B.Store)->nzval = barray; 218 ((DNformat*)lu->X.Store)->nzval = xarray; 219 #endif 220 221 /* Initialize the statistics variables. */ 222 StatInit(&stat); 223 224 lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */ 225 lu->options.Trans = TRANS; 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__ "MatLUFactorNumeric_SuperLU" 265 PetscErrorCode MatLUFactorNumeric_SuperLU(Mat A,MatFactorInfo *info,Mat *F) 266 { 267 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data; 268 Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr; 269 PetscErrorCode ierr; 270 PetscInt sinfo; 271 SuperLUStat_t stat; 272 PetscReal ferr, berr; 273 NCformat *Ustore; 274 SCformat *Lstore; 275 276 PetscFunctionBegin; 277 if (lu->flg == SAME_NONZERO_PATTERN){ /* successing numerical factorization */ 278 lu->options.Fact = SamePattern; 279 /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */ 280 Destroy_SuperMatrix_Store(&lu->A); 281 if ( lu->lwork >= 0 ) { 282 Destroy_SuperNode_Matrix(&lu->L); 283 Destroy_CompCol_Matrix(&lu->U); 284 lu->options.Fact = SamePattern; 285 } 286 } 287 288 /* Create the SuperMatrix for lu->A=A^T: 289 Since SuperLU likes column-oriented matrices,we pass it the transpose, 290 and then solve A^T X = B in MatSolve(). */ 291 #if defined(PETSC_USE_COMPLEX) 292 zCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i, 293 SLU_NC,SLU_Z,SLU_GE); 294 #else 295 dCreate_CompCol_Matrix(&lu->A,A->n,A->m,aa->nz,aa->a,aa->j,aa->i, 296 SLU_NC,SLU_D,SLU_GE); 297 #endif 298 299 /* Initialize the statistics variables. */ 300 StatInit(&stat); 301 302 /* Numerical factorization */ 303 lu->B.ncol = 0; /* Indicate not to solve the system */ 304 #if defined(PETSC_USE_COMPLEX) 305 zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 306 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 307 &lu->mem_usage, &stat, &sinfo); 308 #else 309 dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C, 310 &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr, 311 &lu->mem_usage, &stat, &sinfo); 312 #endif 313 if ( !sinfo || sinfo == lu->A.ncol+1 ) { 314 if ( lu->options.PivotGrowth ) 315 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. pivot growth = %e\n", lu->rpg); 316 if ( lu->options.ConditionNumber ) 317 ierr = PetscPrintf(PETSC_COMM_SELF," Recip. condition number = %e\n", lu->rcond); 318 } else if ( sinfo > 0 ){ 319 if ( lu->lwork == -1 ) { 320 ierr = PetscPrintf(PETSC_COMM_SELF," ** Estimated memory: %D bytes\n", sinfo - lu->A.ncol); 321 } else { 322 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: gssvx() returns info %D\n",sinfo); 323 } 324 } else { /* sinfo < 0 */ 325 SETERRQ2(PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo); 326 } 327 328 if ( lu->options.PrintStat ) { 329 ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n"); 330 StatPrint(&stat); 331 Lstore = (SCformat *) lu->L.Store; 332 Ustore = (NCformat *) lu->U.Store; 333 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor L = %D\n", Lstore->nnz); 334 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in factor U = %D\n", Ustore->nnz); 335 ierr = PetscPrintf(PETSC_COMM_SELF," No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol); 336 ierr = PetscPrintf(PETSC_COMM_SELF," L\\U MB %.3f\ttotal MB needed %.3f\texpansions %D\n", 337 lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6, 338 lu->mem_usage.expansions); 339 } 340 StatFree(&stat); 341 342 lu->flg = SAME_NONZERO_PATTERN; 343 PetscFunctionReturn(0); 344 } 345 346 /* 347 Note the r permutation is ignored 348 */ 349 #undef __FUNCT__ 350 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU" 351 PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 352 { 353 Mat B; 354 Mat_SuperLU *lu; 355 PetscErrorCode ierr; 356 PetscInt m=A->m,n=A->n,indx; 357 PetscTruth flg; 358 const char *colperm[]={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */ 359 const char *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"}; 360 const char *rowperm[]={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */ 361 362 PetscFunctionBegin; 363 ierr = MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);CHKERRQ(ierr); 364 ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 365 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 366 367 B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU; 368 B->ops->solve = MatSolve_SuperLU; 369 B->factor = FACTOR_LU; 370 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 371 372 lu = (Mat_SuperLU*)(B->spptr); 373 374 /* Set SuperLU options */ 375 /* the default values for options argument: 376 options.Fact = DOFACT; 377 options.Equil = YES; 378 options.ColPerm = COLAMD; 379 options.DiagPivotThresh = 1.0; 380 options.Trans = NOTRANS; 381 options.IterRefine = NOREFINE; 382 options.SymmetricMode = NO; 383 options.PivotGrowth = NO; 384 options.ConditionNumber = NO; 385 options.PrintStat = YES; 386 */ 387 set_default_options(&lu->options); 388 /* equilibration causes error in solve(), thus not supported here. See dgssvx.c for possible reason. */ 389 lu->options.Equil = NO; 390 lu->options.PrintStat = NO; 391 lu->lwork = 0; /* allocate space internally by system malloc */ 392 393 ierr = PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");CHKERRQ(ierr); 394 /* 395 ierr = PetscOptionsLogical("-mat_superlu_equil","Equil","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 396 if (flg) lu->options.Equil = YES; -- not supported by the interface !!! 397 */ 398 ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr); 399 if (flg) {lu->options.ColPerm = (colperm_t)indx;} 400 ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr); 401 if (flg) { lu->options.IterRefine = (IterRefine_t)indx;} 402 ierr = PetscOptionsLogical("-mat_superlu_symmetricmode","SymmetricMode","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 403 if (flg) lu->options.SymmetricMode = YES; 404 ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&lu->options.DiagPivotThresh,PETSC_NULL);CHKERRQ(ierr); 405 ierr = PetscOptionsLogical("-mat_superlu_pivotgrowth","PivotGrowth","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 406 if (flg) lu->options.PivotGrowth = YES; 407 ierr = PetscOptionsLogical("-mat_superlu_conditionnumber","ConditionNumber","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 408 if (flg) lu->options.ConditionNumber = YES; 409 ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr); 410 if (flg) {lu->options.RowPerm = (rowperm_t)indx;} 411 ierr = PetscOptionsLogical("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 412 if (flg) lu->options.ReplaceTinyPivot = YES; 413 ierr = PetscOptionsLogical("-mat_superlu_printstat","PrintStat","None",PETSC_FALSE,&flg,0);CHKERRQ(ierr); 414 if (flg) lu->options.PrintStat = YES; 415 ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,PETSC_NULL);CHKERRQ(ierr); 416 if (lu->lwork > 0 ){ 417 ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr); 418 } else if (lu->lwork != 0 && lu->lwork != -1){ 419 ierr = PetscPrintf(PETSC_COMM_SELF," Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork); 420 lu->lwork = 0; 421 } 422 PetscOptionsEnd(); 423 424 #ifdef SUPERLU2 425 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU", 426 (void(*)(void))MatCreateNull_SuperLU);CHKERRQ(ierr); 427 #endif 428 429 /* Allocate spaces (notice sizes are for the transpose) */ 430 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->etree);CHKERRQ(ierr); 431 ierr = PetscMalloc(n*sizeof(PetscInt),&lu->perm_r);CHKERRQ(ierr); 432 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->perm_c);CHKERRQ(ierr); 433 ierr = PetscMalloc(n*sizeof(PetscInt),&lu->R);CHKERRQ(ierr); 434 ierr = PetscMalloc(m*sizeof(PetscInt),&lu->C);CHKERRQ(ierr); 435 436 /* create rhs and solution x without allocate space for .Store */ 437 #if defined(PETSC_USE_COMPLEX) 438 zCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 439 zCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_Z, SLU_GE); 440 #else 441 dCreate_Dense_Matrix(&lu->B, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 442 dCreate_Dense_Matrix(&lu->X, m, 1, PETSC_NULL, m, SLU_DN, SLU_D, SLU_GE); 443 #endif 444 445 lu->flg = DIFFERENT_NONZERO_PATTERN; 446 lu->CleanUpSuperLU = PETSC_TRUE; 447 448 *F = B; 449 ierr = PetscLogObjectMemory(B,(A->m+A->n)*sizeof(PetscInt)+sizeof(Mat_SuperLU));CHKERRQ(ierr); 450 PetscFunctionReturn(0); 451 } 452 453 /* used by -ksp_view */ 454 #undef __FUNCT__ 455 #define __FUNCT__ "MatFactorInfo_SuperLU" 456 PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer) 457 { 458 Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr; 459 PetscErrorCode ierr; 460 superlu_options_t options; 461 462 PetscFunctionBegin; 463 /* check if matrix is superlu_dist type */ 464 if (A->ops->solve != MatSolve_SuperLU) PetscFunctionReturn(0); 465 466 options = lu->options; 467 ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr); 468 ierr = PetscViewerASCIIPrintf(viewer," Equil: %s\n",(options.Equil != NO) ? "YES": "NO");CHKERRQ(ierr); 469 ierr = PetscViewerASCIIPrintf(viewer," ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr); 470 ierr = PetscViewerASCIIPrintf(viewer," IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr); 471 ierr = PetscViewerASCIIPrintf(viewer," SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES": "NO");CHKERRQ(ierr); 472 ierr = PetscViewerASCIIPrintf(viewer," DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr); 473 ierr = PetscViewerASCIIPrintf(viewer," PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES": "NO");CHKERRQ(ierr); 474 ierr = PetscViewerASCIIPrintf(viewer," ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES": "NO");CHKERRQ(ierr); 475 ierr = PetscViewerASCIIPrintf(viewer," RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr); 476 ierr = PetscViewerASCIIPrintf(viewer," ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES": "NO");CHKERRQ(ierr); 477 ierr = PetscViewerASCIIPrintf(viewer," PrintStat: %s\n",(options.PrintStat != NO) ? "YES": "NO");CHKERRQ(ierr); 478 ierr = PetscViewerASCIIPrintf(viewer," lwork: %D\n",lu->lwork);CHKERRQ(ierr); 479 480 PetscFunctionReturn(0); 481 } 482 483 #undef __FUNCT__ 484 #define __FUNCT__ "MatDuplicate_SuperLU" 485 PetscErrorCode MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) { 486 PetscErrorCode ierr; 487 Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 488 489 PetscFunctionBegin; 490 ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 491 ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_SuperLU));CHKERRQ(ierr); 492 PetscFunctionReturn(0); 493 } 494 495 EXTERN_C_BEGIN 496 #undef __FUNCT__ 497 #define __FUNCT__ "MatConvert_SuperLU_SeqAIJ" 498 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SuperLU_SeqAIJ(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 499 { 500 /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */ 501 /* to its base PETSc type, so we will ignore 'MatType type'. */ 502 PetscErrorCode ierr; 503 Mat B=*newmat; 504 Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr; 505 506 PetscFunctionBegin; 507 if (reuse == MAT_INITIAL_MATRIX) { 508 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 509 } 510 /* Reset the original function pointers */ 511 B->ops->duplicate = lu->MatDuplicate; 512 B->ops->view = lu->MatView; 513 B->ops->assemblyend = lu->MatAssemblyEnd; 514 B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic; 515 B->ops->destroy = lu->MatDestroy; 516 /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */ 517 ierr = PetscFree(lu);CHKERRQ(ierr); 518 519 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_superlu_C","",PETSC_NULL);CHKERRQ(ierr); 520 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_superlu_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 521 522 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); 523 *newmat = B; 524 PetscFunctionReturn(0); 525 } 526 EXTERN_C_END 527 528 EXTERN_C_BEGIN 529 #undef __FUNCT__ 530 #define __FUNCT__ "MatConvert_SeqAIJ_SuperLU" 531 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SuperLU(Mat A,const MatType type,MatReuse reuse,Mat *newmat) 532 { 533 /* This routine is only called to convert to MATSUPERLU */ 534 /* from MATSEQAIJ, so we will ignore 'MatType type'. */ 535 PetscErrorCode ierr; 536 Mat B=*newmat; 537 Mat_SuperLU *lu; 538 539 PetscFunctionBegin; 540 if (reuse == MAT_INITIAL_MATRIX) { 541 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 542 } 543 544 ierr = PetscNew(Mat_SuperLU,&lu);CHKERRQ(ierr); 545 lu->MatDuplicate = A->ops->duplicate; 546 lu->MatView = A->ops->view; 547 lu->MatAssemblyEnd = A->ops->assemblyend; 548 lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 549 lu->MatDestroy = A->ops->destroy; 550 lu->CleanUpSuperLU = PETSC_FALSE; 551 552 B->spptr = (void*)lu; 553 B->ops->duplicate = MatDuplicate_SuperLU; 554 B->ops->view = MatView_SuperLU; 555 B->ops->assemblyend = MatAssemblyEnd_SuperLU; 556 B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU; 557 B->ops->choleskyfactorsymbolic = 0; 558 B->ops->destroy = MatDestroy_SuperLU; 559 560 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C", 561 "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);CHKERRQ(ierr); 562 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C", 563 "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);CHKERRQ(ierr); 564 ierr = PetscLogInfo((0,"MatConvert_SeqAIJ_SuperLU:Using SuperLU for SeqAIJ LU factorization and solves.\n"));CHKERRQ(ierr); 565 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);CHKERRQ(ierr); 566 *newmat = B; 567 PetscFunctionReturn(0); 568 } 569 EXTERN_C_END 570 571 /*MC 572 MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices 573 via the external package SuperLU. 574 575 If SuperLU is installed (see the manual for 576 instructions on how to declare the existence of external packages), 577 a matrix type can be constructed which invokes SuperLU solvers. 578 After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU). 579 580 This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is 581 supported for this matrix type. One can also call MatConvert for an inplace conversion to or from 582 the MATSEQAIJ type without data copy. 583 584 Options Database Keys: 585 + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions() 586 - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering, 587 1: MMD applied to A'*A, 588 2: MMD applied to A'+A, 589 3: COLAMD, approximate minimum degree column ordering 590 591 Level: beginner 592 593 .seealso: PCLU 594 M*/ 595 596 EXTERN_C_BEGIN 597 #undef __FUNCT__ 598 #define __FUNCT__ "MatCreate_SuperLU" 599 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SuperLU(Mat A) 600 { 601 PetscErrorCode ierr; 602 603 PetscFunctionBegin; 604 /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */ 605 ierr = PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);CHKERRQ(ierr); 606 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 607 ierr = MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 608 PetscFunctionReturn(0); 609 } 610 EXTERN_C_END 611