1 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 2 #include <../src/mat/impls/aij/mpi/mpiaij.h> 3 #include <StrumpackSparseSolver.h> 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "MatGetDiagonal_STRUMPACK" 7 static PetscErrorCode MatGetDiagonal_STRUMPACK(Mat A,Vec v) 8 { 9 PetscFunctionBegin; 10 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: STRUMPACK factor"); 11 PetscFunctionReturn(0); 12 } 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "MatDestroy_STRUMPACK" 16 static PetscErrorCode MatDestroy_STRUMPACK(Mat A) 17 { 18 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 19 PetscErrorCode ierr; 20 PetscBool flg; 21 22 PetscFunctionBegin; 23 /* Deallocate STRUMPACK storage */ 24 PetscStackCall("STRUMPACK_destroy",STRUMPACK_destroy(S)); 25 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 26 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 27 if (flg) { 28 ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr); 29 } else { 30 ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 31 } 32 33 /* clear composed functions */ 34 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 35 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSRelCompTol_C",NULL);CHKERRQ(ierr); 36 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinFrontSize_C",NULL);CHKERRQ(ierr); 37 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetHSSMinSepSize_C",NULL);CHKERRQ(ierr); 38 ierr = PetscObjectComposeFunction((PetscObject)A,"MatSTRUMPACKSetColPerm_C",NULL);CHKERRQ(ierr); 39 40 PetscFunctionReturn(0); 41 } 42 43 #undef __FUNCT__ 44 #define __FUNCT__ "MatSTRUMPACKSetHSSRelCompTol_STRUMPACK" 45 static PetscErrorCode MatSTRUMPACKSetHSSRelCompTol_STRUMPACK(Mat F,PetscReal rctol) 46 { 47 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 48 49 PetscFunctionBegin; 50 PetscStackCall("STRUMPACK_set_hss_rel_tol", STRUMPACK_set_hss_rel_tol(*S,rctol)); 51 PetscFunctionReturn(0); 52 } 53 54 #undef __FUNCT__ 55 #define __FUNCT__ "MatSTRUMPACKSetHSSRelCompTol" 56 /*@ 57 MatSTRUMPACKSetHSSRelCompTol - Set STRUMPACK relative tolerance for HSS compression 58 Logically Collective on Mat 59 60 Input Parameters: 61 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 62 - rctol - relative compression tolerance 63 64 Options Database: 65 . -mat_strumpack_rctol <rctol> 66 67 Level: beginner 68 69 References: 70 . STRUMPACK manual 71 72 .seealso: MatGetFactor() 73 @*/ 74 PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat F,PetscReal rctol) 75 { 76 PetscErrorCode ierr; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 80 PetscValidLogicalCollectiveReal(F,rctol,2); 81 ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSRelCompTol_C",(Mat,PetscReal),(F,rctol));CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK" 87 static PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK(Mat F,PetscInt hssminsize) 88 { 89 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 90 91 PetscFunctionBegin; 92 PetscStackCall("STRUMPACK_set_HSS_min_front_size", STRUMPACK_set_HSS_min_front_size(*S,hssminsize)); 93 PetscFunctionReturn(0); 94 } 95 #undef __FUNCT__ 96 #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize_STRUMPACK" 97 static PetscErrorCode MatSTRUMPACKSetHSSMinSepSize_STRUMPACK(Mat F,PetscInt hssminsize) 98 { 99 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 100 101 PetscFunctionBegin; 102 PetscStackCall("STRUMPACK_set_HSS_min_sep_size", STRUMPACK_set_HSS_min_sep_size(*S,hssminsize)); 103 PetscFunctionReturn(0); 104 } 105 106 #undef __FUNCT__ 107 #define __FUNCT__ "MatSTRUMPACKSetHSSMinFrontSize" 108 /*@ 109 MatSTRUMPACKSetHSSMinFrontSize - Set STRUMPACK minimum dense matrix size for low-rank approximation 110 Logically Collective on Mat 111 112 Input Parameters: 113 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 114 - hssminsize - minimum dense matrix size for low-rank approximation 115 116 Options Database: 117 . -mat_strumpack_hss_min_front_size <hssminsize> 118 119 Level: beginner 120 121 References: 122 . STRUMPACK manual 123 124 .seealso: MatGetFactor() 125 @*/ 126 PetscErrorCode MatSTRUMPACKSetHSSMinFrontSize(Mat F,PetscInt hssminsize) 127 { 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 132 PetscValidLogicalCollectiveInt(F,hssminsize,2); 133 ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinFrontSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "MatSTRUMPACKSetHSSMinSepSize" 139 /*@ 140 MatSTRUMPACKSetHSSMinSepSize - Set STRUMPACK minimum separator size for low-rank approximation 141 Logically Collective on Mat 142 143 Input Parameters: 144 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 145 - hssminsize - minimum dense matrix size for low-rank approximation 146 147 Options Database: 148 . -mat_strumpack_hss_min_sep_size <hssminsize> 149 150 Level: beginner 151 152 References: 153 . STRUMPACK manual 154 155 .seealso: MatGetFactor() 156 @*/ 157 PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat F,PetscInt hssminsize) 158 { 159 PetscErrorCode ierr; 160 161 PetscFunctionBegin; 162 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 163 PetscValidLogicalCollectiveInt(F,hssminsize,2); 164 ierr = PetscTryMethod(F,"MatSTRUMPACKSetHSSMinSepSize_C",(Mat,PetscInt),(F,hssminsize));CHKERRQ(ierr); 165 PetscFunctionReturn(0); 166 } 167 168 169 #undef __FUNCT__ 170 #define __FUNCT__ "MatSTRUMPACKSetColPerm_STRUMPACK" 171 static PetscErrorCode MatSTRUMPACKSetColPerm_STRUMPACK(Mat F,PetscBool cperm) 172 { 173 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 174 175 PetscFunctionBegin; 176 PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,cperm ? 5 : 0)); 177 PetscFunctionReturn(0); 178 } 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "MatSTRUMPACKSetColPerm" 182 /*@ 183 MatSTRUMPACKSetColPerm - Set whether STRUMPACK should try to permute the columns of the matrix in order to get a nonzero diagonal 184 Logically Collective on Mat 185 186 Input Parameters: 187 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-STRUMPACK interface 188 - cperm - whether or not to permute (internally) the columns of the matrix 189 190 Options Database: 191 . -mat_strumpack_colperm <cperm> 192 193 Level: beginner 194 195 References: 196 . STRUMPACK manual 197 198 .seealso: MatGetFactor() 199 @*/ 200 PetscErrorCode MatSTRUMPACKSetColPerm(Mat F,PetscBool cperm) 201 { 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 PetscValidHeaderSpecific(F,MAT_CLASSID,1); 206 PetscValidLogicalCollectiveBool(F,cperm,2); 207 ierr = PetscTryMethod(F,"MatSTRUMPACKSetColPerm_C",(Mat,PetscBool),(F,cperm));CHKERRQ(ierr); 208 PetscFunctionReturn(0); 209 } 210 211 #undef __FUNCT__ 212 #define __FUNCT__ "MatSolve_STRUMPACK" 213 static PetscErrorCode MatSolve_STRUMPACK(Mat A,Vec b_mpi,Vec x) 214 { 215 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)A->spptr; 216 STRUMPACK_RETURN_CODE sp_err; 217 PetscErrorCode ierr; 218 const PetscScalar *bptr; 219 PetscScalar *xptr; 220 221 PetscFunctionBegin; 222 ierr = VecGetArray(x,&xptr);CHKERRQ(ierr); 223 ierr = VecGetArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 224 225 PetscStackCall("STRUMPACK_solve",sp_err = STRUMPACK_solve(*S,(PetscScalar*)bptr,xptr,0)); 226 switch (sp_err) { 227 case STRUMPACK_SUCCESS: break; 228 case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 229 case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 230 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: solve failed"); 231 } 232 ierr = VecRestoreArray(x,&xptr);CHKERRQ(ierr); 233 ierr = VecRestoreArrayRead(b_mpi,&bptr);CHKERRQ(ierr); 234 PetscFunctionReturn(0); 235 } 236 237 #undef __FUNCT__ 238 #define __FUNCT__ "MatMatSolve_STRUMPACK" 239 static PetscErrorCode MatMatSolve_STRUMPACK(Mat A,Mat B_mpi,Mat X) 240 { 241 PetscErrorCode ierr; 242 PetscBool flg; 243 244 PetscFunctionBegin; 245 ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 246 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 247 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 248 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 249 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_STRUMPACK() is not implemented yet"); 250 PetscFunctionReturn(0); 251 } 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "MatFactorInfo_STRUMPACK" 255 static PetscErrorCode MatFactorInfo_STRUMPACK(Mat A,PetscViewer viewer) 256 { 257 PetscErrorCode ierr; 258 259 PetscFunctionBegin; 260 /* check if matrix is strumpack type */ 261 if (A->ops->solve != MatSolve_STRUMPACK) PetscFunctionReturn(0); 262 ierr = PetscViewerASCIIPrintf(viewer,"STRUMPACK sparse solver!\n");CHKERRQ(ierr); 263 PetscFunctionReturn(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "MatView_STRUMPACK" 268 static PetscErrorCode MatView_STRUMPACK(Mat A,PetscViewer viewer) 269 { 270 PetscErrorCode ierr; 271 PetscBool iascii; 272 PetscViewerFormat format; 273 274 PetscFunctionBegin; 275 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 276 if (iascii) { 277 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 278 if (format == PETSC_VIEWER_ASCII_INFO) { 279 ierr = MatFactorInfo_STRUMPACK(A,viewer);CHKERRQ(ierr); 280 } 281 } 282 PetscFunctionReturn(0); 283 } 284 285 #undef __FUNCT__ 286 #define __FUNCT__ "MatLUFactorNumeric_STRUMPACK" 287 static PetscErrorCode MatLUFactorNumeric_STRUMPACK(Mat F,Mat A,const MatFactorInfo *info) 288 { 289 STRUMPACK_SparseSolver *S = (STRUMPACK_SparseSolver*)F->spptr; 290 STRUMPACK_RETURN_CODE sp_err; 291 Mat_SeqAIJ *A_d,*A_o; 292 Mat_MPIAIJ *mat; 293 PetscErrorCode ierr; 294 PetscInt M=A->rmap->N,m=A->rmap->n; 295 PetscBool flg; 296 297 PetscFunctionBegin; 298 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 299 if (flg) { /* A is MATMPIAIJ */ 300 mat = (Mat_MPIAIJ*)A->data; 301 A_d = (Mat_SeqAIJ*)(mat->A)->data; 302 A_o = (Mat_SeqAIJ*)(mat->B)->data; 303 PetscStackCall("STRUMPACK_set_MPIAIJ_matrix",STRUMPACK_set_MPIAIJ_matrix(*S,&m,A_d->i,A_d->j,A_d->a,A_o->i,A_o->j,A_o->a,mat->garray)); 304 } else { /* A is MATSEQAIJ */ 305 A_d = (Mat_SeqAIJ*)A->data; 306 PetscStackCall("STRUMPACK_set_csr_matrix",STRUMPACK_set_csr_matrix(*S,&M,A_d->i,A_d->j,A_d->a,0)); 307 } 308 309 /* Reorder and Factor the matrix. */ 310 /* TODO figure out how to avoid reorder if the matrix values changed, but the pattern remains the same. */ 311 PetscStackCall("STRUMPACK_reorder",sp_err = STRUMPACK_reorder(*S)); 312 PetscStackCall("STRUMPACK_factor",sp_err = STRUMPACK_factor(*S)); 313 switch (sp_err) { 314 case STRUMPACK_SUCCESS: break; 315 case STRUMPACK_MATRIX_NOT_SET: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix was not set"); break; } 316 case STRUMPACK_REORDERING_ERROR: { SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: matrix reordering failed"); break; } 317 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"STRUMPACK error: factorization failed"); 318 } 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNCT__ 323 #define __FUNCT__ "MatLUFactorSymbolic_STRUMPACK" 324 static PetscErrorCode MatLUFactorSymbolic_STRUMPACK(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 325 { 326 PetscFunctionBegin; 327 F->ops->lufactornumeric = MatLUFactorNumeric_STRUMPACK; 328 F->ops->solve = MatSolve_STRUMPACK; 329 F->ops->matsolve = MatMatSolve_STRUMPACK; 330 PetscFunctionReturn(0); 331 } 332 333 #undef __FUNCT__ 334 #define __FUNCT__ "MatFactorGetSolverPackage_aij_strumpack" 335 static PetscErrorCode MatFactorGetSolverPackage_aij_strumpack(Mat A,const MatSolverPackage *type) 336 { 337 PetscFunctionBegin; 338 *type = MATSOLVERSTRUMPACK; 339 PetscFunctionReturn(0); 340 } 341 342 /*MC 343 MATSOLVERSSTRUMPACK = "strumpack" - A solver package providing a direct sparse solver (PCLU) 344 and a preconditioner (PCILU) using low-rank compression via the external package STRUMPACK. 345 346 Consult the STRUMPACK-sparse manual for more info. 347 348 Use 349 ./configure --download-strumpack 350 to have PETSc installed with STRUMPACK 351 352 Use 353 -pc_type lu -pc_factor_mat_solver_package strumpack 354 to use this as an exact (direct) solver, use 355 -pc_type ilu -pc_factor_mat_solver_package strumpack 356 to enable low-rank compression (i.e, use as a preconditioner). 357 358 Works with AIJ matrices 359 360 Options Database Keys: 361 + -mat_strumpack_rctol <1e-4> - Relative compression tolerance (None) 362 . -mat_strumpack_hss_min_front_size - Minimum size of dense block for HSS compression (None) 363 . -mat_strumpack_hss_min_sep_size - Minimum size of separator for HSS compression (None) 364 . -mat_strumpack_colperm <TRUE> - Permute matrix to make diagonal nonzeros (None) 365 366 Level: beginner 367 368 .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage 369 M*/ 370 #undef __FUNCT__ 371 #define __FUNCT__ "MatGetFactor_aij_strumpack" 372 static PetscErrorCode MatGetFactor_aij_strumpack(Mat A,MatFactorType ftype,Mat *F) 373 { 374 Mat B; 375 STRUMPACK_SparseSolver *S; 376 PetscErrorCode ierr; 377 PetscInt M=A->rmap->N,N=A->cmap->N; 378 STRUMPACK_INTERFACE iface; 379 PetscBool verb,flg,set; 380 PetscReal rctol; 381 PetscInt hssminsize; 382 int argc; 383 char **args,*copts,*pname; 384 size_t len; 385 const STRUMPACK_PRECISION table[2][2][2] = {{{STRUMPACK_FLOATCOMPLEX_64,STRUMPACK_DOUBLECOMPLEX_64}, 386 {STRUMPACK_FLOAT_64,STRUMPACK_DOUBLE_64}}, 387 {{STRUMPACK_FLOATCOMPLEX,STRUMPACK_DOUBLECOMPLEX}, 388 {STRUMPACK_FLOAT,STRUMPACK_DOUBLE}}}; 389 const STRUMPACK_PRECISION prec = table[(sizeof(PetscInt)==8)?0:1][(PETSC_SCALAR==PETSC_COMPLEX)?0:1][(PETSC_REAL==PETSC_FLOAT)?0:1]; 390 391 PetscFunctionBegin; 392 /* Create the factorization matrix */ 393 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 394 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr); 395 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 396 ierr = MatSeqAIJSetPreallocation(B,0,NULL); 397 ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 398 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 399 B->ops->lufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 400 B->ops->ilufactorsymbolic = MatLUFactorSymbolic_STRUMPACK; 401 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 402 B->ops->view = MatView_STRUMPACK; 403 B->ops->destroy = MatDestroy_STRUMPACK; 404 B->ops->getdiagonal = MatGetDiagonal_STRUMPACK; 405 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_strumpack);CHKERRQ(ierr); 406 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSRelCompTol_C",MatSTRUMPACKSetHSSRelCompTol_STRUMPACK);CHKERRQ(ierr); 407 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinFrontSize_C",MatSTRUMPACKSetHSSMinFrontSize_STRUMPACK);CHKERRQ(ierr); 408 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetHSSMinSepSize_C",MatSTRUMPACKSetHSSMinSepSize_STRUMPACK);CHKERRQ(ierr); 409 ierr = PetscObjectComposeFunction((PetscObject)B,"MatSTRUMPACKSetColPerm_C",MatSTRUMPACKSetColPerm_STRUMPACK);CHKERRQ(ierr); 410 B->factortype = ftype; 411 ierr = PetscNewLog(B,&S);CHKERRQ(ierr); 412 B->spptr = S; 413 414 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg); 415 iface = flg ? STRUMPACK_MT : STRUMPACK_MPI_DIST; 416 417 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"STRUMPACK Options","Mat");CHKERRQ(ierr); 418 419 verb = PetscLogPrintInfo ? PETSC_TRUE : PETSC_FALSE; 420 ierr = PetscOptionsBool("-mat_strumpack_verbose","Print STRUMPACK information","None",verb,&verb,NULL);CHKERRQ(ierr); 421 422 ierr = PetscOptionsGetAll(NULL,&copts);CHKERRQ(ierr); 423 ierr = PetscStrlen(copts,&len);CHKERRQ(ierr); 424 len += PETSC_MAX_PATH_LEN+1; 425 ierr = PetscMalloc1(len,&pname);CHKERRQ(ierr); 426 /* first string is assumed to be the program name, so add program name to options string */ 427 ierr = PetscGetProgramName(pname,len);CHKERRQ(ierr); 428 ierr = PetscStrcat(pname," ");CHKERRQ(ierr); 429 ierr = PetscStrcat(pname,copts);CHKERRQ(ierr); 430 ierr = PetscStrToArray(pname,' ',&argc,&args);CHKERRQ(ierr); 431 ierr = PetscFree(copts);CHKERRQ(ierr); 432 ierr = PetscFree(pname);CHKERRQ(ierr); 433 434 PetscStackCall("STRUMPACK_init",STRUMPACK_init(S,PetscObjectComm((PetscObject)A),prec,iface,argc,args,verb)); 435 436 PetscStackCall("STRUMPACK_hss_rel_tol",rctol = (PetscReal)STRUMPACK_hss_rel_tol(*S)); 437 ierr = PetscOptionsReal("-mat_strumpack_hss_rel_tol","Relative compression tolerance","None",rctol,&rctol,&set);CHKERRQ(ierr); 438 if (set) PetscStackCall("STRUMPACK_set_hss_rel_tol",STRUMPACK_set_hss_rel_tol(*S,(double)rctol)); 439 440 PetscStackCall("STRUMPACK_mc64job",flg = (STRUMPACK_mc64job(*S) == 0) ? PETSC_FALSE : PETSC_TRUE); 441 ierr = PetscOptionsBool("-mat_strumpack_colperm","Find a col perm to get nonzero diagonal","None",flg,&flg,&set);CHKERRQ(ierr); 442 if (set) PetscStackCall("STRUMPACK_set_mc64job",STRUMPACK_set_mc64job(*S,flg ? 5 : 0)); 443 444 PetscStackCall("STRUMPACK_HSS_min_front_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_front_size(*S)); 445 ierr = PetscOptionsInt("-mat_strumpack_hss_min_front_size","Minimum size of dense block for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 446 if (set) PetscStackCall("STRUMPACK_set_HSS_min_front_size",STRUMPACK_set_HSS_min_front_size(*S,(int)hssminsize)); 447 448 PetscStackCall("STRUMPACK_HSS_min_sep_size",hssminsize = (PetscInt)STRUMPACK_HSS_min_sep_size(*S)); 449 ierr = PetscOptionsInt("-mat_strumpack_hss_min_sep_size","Minimum size of separator for HSS compression","None",hssminsize,&hssminsize,&set);CHKERRQ(ierr); 450 if (set) PetscStackCall("STRUMPACK_set_HSS_min_sep_size",STRUMPACK_set_HSS_min_sep_size(*S,(int)hssminsize)); 451 452 PetscOptionsEnd(); 453 454 if (ftype == MAT_FACTOR_ILU) { 455 /* When enabling HSS compression, the STRUMPACK solver becomes an incomplete */ 456 /* (or approximate) LU factorization. */ 457 PetscStackCall("STRUMPACK_enable_HSS",STRUMPACK_enable_HSS(*S)); 458 /* Disable the outer iterative solver from STRUMPACK. */ 459 /* When STRUMPACK is used as a direct solver, it will by default do iterative refinement. */ 460 /* When STRUMPACK is used with as an approximate factorization preconditioner (by enabling */ 461 /* low-rank compression), it will use it's own GMRES. Here we can disable the */ 462 /* outer iterative solver, as PETSc uses STRUMPACK from within a KSP. */ 463 PetscStackCall("STRUMPACK_set_Krylov_solver", STRUMPACK_set_Krylov_solver(*S, STRUMPACK_DIRECT)); 464 } 465 466 /* Allow the user to set or overwrite the above options from the command line */ 467 PetscStackCall("STRUMPACK_set_from_options",STRUMPACK_set_from_options(*S)); 468 ierr = PetscStrToArrayDestroy(argc,args);CHKERRQ(ierr); 469 470 *F = B; 471 PetscFunctionReturn(0); 472 } 473 474 #undef __FUNCT__ 475 #define __FUNCT__ "MatSolverPackageRegister_STRUMPACK" 476 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_STRUMPACK(void) 477 { 478 PetscErrorCode ierr; 479 480 PetscFunctionBegin; 481 ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 482 ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 483 ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATMPIAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 484 ierr = MatSolverPackageRegister(MATSOLVERSTRUMPACK,MATSEQAIJ,MAT_FACTOR_ILU,MatGetFactor_aij_strumpack);CHKERRQ(ierr); 485 PetscFunctionReturn(0); 486 } 487