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