1 #define PETSCMAT_DLL 2 3 /* 4 Provides an interface to the MUMPS sparse solver 5 */ 6 #include "../src/mat/impls/aij/seq/aij.h" 7 #include "../src/mat/impls/aij/mpi/mpiaij.h" 8 #include "../src/mat/impls/sbaij/seq/sbaij.h" 9 #include "../src/mat/impls/sbaij/mpi/mpisbaij.h" 10 11 EXTERN_C_BEGIN 12 #if defined(PETSC_USE_COMPLEX) 13 #include "zmumps_c.h" 14 #else 15 #include "dmumps_c.h" 16 #endif 17 EXTERN_C_END 18 #define JOB_INIT -1 19 #define JOB_END -2 20 /* macros s.t. indices match MUMPS documentation */ 21 #define ICNTL(I) icntl[(I)-1] 22 #define CNTL(I) cntl[(I)-1] 23 #define INFOG(I) infog[(I)-1] 24 #define INFO(I) info[(I)-1] 25 #define RINFOG(I) rinfog[(I)-1] 26 #define RINFO(I) rinfo[(I)-1] 27 28 typedef struct { 29 #if defined(PETSC_USE_COMPLEX) 30 ZMUMPS_STRUC_C id; 31 #else 32 DMUMPS_STRUC_C id; 33 #endif 34 MatStructure matstruc; 35 PetscMPIInt myid,size; 36 PetscInt *irn,*jcn,sym,nSolve; 37 PetscScalar *val; 38 MPI_Comm comm_mumps; 39 VecScatter scat_rhs, scat_sol; 40 PetscTruth isAIJ,CleanUpMUMPS; 41 Vec b_seq,x_seq; 42 PetscErrorCode (*MatDestroy)(Mat); 43 } Mat_MUMPS; 44 45 EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 46 47 /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 48 /* 49 input: 50 A - matrix in mpiaij or mpisbaij (bs=1) format 51 shift - 0: C style output triple; 1: Fortran style output triple. 52 valOnly - FALSE: spaces are allocated and values are set for the triple 53 TRUE: only the values in v array are updated 54 output: 55 nnz - dim of r, c, and v (number of local nonzero entries of A) 56 r, c, v - row and col index, matrix values (matrix triples) 57 */ 58 PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) 59 { 60 PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray; 61 PetscErrorCode ierr; 62 PetscInt i,j,jj,jB,irow,m=A->rmap->n,*ajj,*bjj,countA,countB,colA_start,jcol; 63 PetscInt *row,*col; 64 PetscScalar *av, *bv,*val; 65 PetscTruth isAIJ; 66 67 PetscFunctionBegin; 68 ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr); 69 if (isAIJ){ 70 Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 71 Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 72 Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 73 nz = aa->nz + bb->nz; 74 ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 75 garray = mat->garray; 76 av=aa->a; bv=bb->a; 77 78 } else { 79 Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 80 Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 81 Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 82 if (A->rmap->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap->bs); 83 nz = aa->nz + bb->nz; 84 ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 85 garray = mat->garray; 86 av=aa->a; bv=bb->a; 87 } 88 89 if (!valOnly){ 90 ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr); 91 ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr); 92 ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 93 *r = row; *c = col; *v = val; 94 } else { 95 row = *r; col = *c; val = *v; 96 } 97 *nnz = nz; 98 99 jj = 0; irow = rstart; 100 for ( i=0; i<m; i++ ) { 101 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 102 countA = ai[i+1] - ai[i]; 103 countB = bi[i+1] - bi[i]; 104 bjj = bj + bi[i]; 105 106 /* get jB, the starting local col index for the 2nd B-part */ 107 colA_start = rstart + ajj[0]; /* the smallest col index for A */ 108 j=-1; 109 do { 110 j++; 111 if (j == countB) break; 112 jcol = garray[bjj[j]]; 113 } while (jcol < colA_start); 114 jB = j; 115 116 /* B-part, smaller col index */ 117 colA_start = rstart + ajj[0]; /* the smallest col index for A */ 118 for (j=0; j<jB; j++){ 119 jcol = garray[bjj[j]]; 120 if (!valOnly){ 121 row[jj] = irow + shift; col[jj] = jcol + shift; 122 123 } 124 val[jj++] = *bv++; 125 } 126 /* A-part */ 127 for (j=0; j<countA; j++){ 128 if (!valOnly){ 129 row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 130 } 131 val[jj++] = *av++; 132 } 133 /* B-part, larger col index */ 134 for (j=jB; j<countB; j++){ 135 if (!valOnly){ 136 row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 137 } 138 val[jj++] = *bv++; 139 } 140 irow++; 141 } 142 143 PetscFunctionReturn(0); 144 } 145 146 #undef __FUNCT__ 147 #define __FUNCT__ "MatDestroy_MUMPS" 148 PetscErrorCode MatDestroy_MUMPS(Mat A) 149 { 150 Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 151 PetscErrorCode ierr; 152 PetscMPIInt size=lu->size; 153 154 PetscFunctionBegin; 155 if (lu->CleanUpMUMPS) { 156 /* Terminate instance, deallocate memories */ 157 if (size > 1){ 158 ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 159 ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 160 ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 161 if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);} 162 if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);} 163 ierr = PetscFree(lu->val);CHKERRQ(ierr); 164 } 165 lu->id.job=JOB_END; 166 #if defined(PETSC_USE_COMPLEX) 167 zmumps_c(&lu->id); 168 #else 169 dmumps_c(&lu->id); 170 #endif 171 ierr = PetscFree(lu->irn);CHKERRQ(ierr); 172 ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 173 ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 174 } 175 ierr = (lu->MatDestroy)(A);CHKERRQ(ierr); 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "MatSolve_MUMPS" 181 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 182 { 183 Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 184 PetscScalar *array; 185 Vec x_seq; 186 IS is_iden,is_petsc; 187 PetscErrorCode ierr; 188 PetscInt i; 189 190 PetscFunctionBegin; 191 lu->id.nrhs = 1; 192 x_seq = lu->b_seq; 193 if (lu->size > 1){ 194 /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 195 ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 196 ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 197 if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 198 } else { /* size == 1 */ 199 ierr = VecCopy(b,x);CHKERRQ(ierr); 200 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 201 } 202 if (!lu->myid) { /* define rhs on the host */ 203 #if defined(PETSC_USE_COMPLEX) 204 lu->id.rhs = (mumps_double_complex*)array; 205 #else 206 lu->id.rhs = array; 207 #endif 208 } 209 if (lu->size == 1){ 210 ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 211 } else if (!lu->myid){ 212 ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 213 } 214 215 if (lu->size > 1){ 216 /* distributed solution */ 217 lu->id.ICNTL(21) = 1; 218 if (!lu->nSolve){ 219 /* Create x_seq=sol_loc for repeated use */ 220 PetscInt lsol_loc; 221 PetscScalar *sol_loc; 222 lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 223 ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr); 224 lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc); 225 lu->id.lsol_loc = lsol_loc; 226 lu->id.sol_loc = (PetscScalar *)sol_loc; 227 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 228 } 229 } 230 231 /* solve phase */ 232 /*-------------*/ 233 lu->id.job = 3; 234 #if defined(PETSC_USE_COMPLEX) 235 zmumps_c(&lu->id); 236 #else 237 dmumps_c(&lu->id); 238 #endif 239 if (lu->id.INFOG(1) < 0) { 240 SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 241 } 242 243 if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 244 if (!lu->nSolve){ /* create scatter scat_sol */ 245 ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 246 for (i=0; i<lu->id.lsol_loc; i++){ 247 lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 248 } 249 ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr); /* to */ 250 ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 251 ierr = ISDestroy(is_iden);CHKERRQ(ierr); 252 ierr = ISDestroy(is_petsc);CHKERRQ(ierr); 253 } 254 ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 255 ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 256 } 257 lu->nSolve++; 258 PetscFunctionReturn(0); 259 } 260 261 #if !defined(PETSC_USE_COMPLEX) 262 /* 263 input: 264 F: numeric factor 265 output: 266 nneg: total number of negative pivots 267 nzero: 0 268 npos: (global dimension of F) - nneg 269 */ 270 271 #undef __FUNCT__ 272 #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 273 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 274 { 275 Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 276 PetscErrorCode ierr; 277 PetscMPIInt size; 278 279 PetscFunctionBegin; 280 ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 281 /* MUMPS 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK */ 282 if (size > 1 && lu->id.ICNTL(13) != 1){ 283 SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13)); 284 } 285 if (nneg){ 286 if (!lu->myid){ 287 *nneg = lu->id.INFOG(12); 288 } 289 ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 290 } 291 if (nzero) *nzero = 0; 292 if (npos) *npos = F->rmap->N - (*nneg); 293 PetscFunctionReturn(0); 294 } 295 #endif /* !defined(PETSC_USE_COMPLEX) */ 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "MatFactorNumeric_MUMPS" 299 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 300 { 301 Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 302 PetscErrorCode ierr; 303 PetscInt rnz,nnz,nz=0,i,M=A->rmap->N,*ai,*aj,icntl; 304 PetscTruth valOnly,flg; 305 Mat F_diag; 306 IS is_iden; 307 Vec b; 308 PetscTruth isSeqAIJ,isSeqSBAIJ; 309 310 PetscFunctionBegin; 311 ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 312 ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 313 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 314 (F)->ops->solve = MatSolve_MUMPS; 315 316 /* Initialize a MUMPS instance */ 317 ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 318 ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 319 lu->id.job = JOB_INIT; 320 ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 321 lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 322 323 /* Set mumps options */ 324 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 325 lu->id.par=1; /* host participates factorizaton and solve */ 326 lu->id.sym=lu->sym; 327 if (lu->sym == 2){ 328 ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 329 if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 330 } 331 #if defined(PETSC_USE_COMPLEX) 332 zmumps_c(&lu->id); 333 #else 334 dmumps_c(&lu->id); 335 #endif 336 337 if (isSeqAIJ || isSeqSBAIJ){ 338 lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 339 } else { 340 lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 341 } 342 343 icntl=-1; 344 lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 345 ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 346 if ((flg && icntl > 0) || PetscLogPrintInfo) { 347 lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 348 } else { /* no output */ 349 lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 350 lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 351 lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 352 } 353 ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): matrix prescaling (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr); 354 icntl=-1; 355 ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 356 if (flg) { 357 if (icntl== 1){ 358 SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 359 } else { 360 lu->id.ICNTL(7) = icntl; 361 } 362 } 363 ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): A or A^T x=b to be solved. 1: A; otherwise: A^T","None",lu->id.ICNTL(9),&lu->id.ICNTL(9),PETSC_NULL);CHKERRQ(ierr); 364 ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr); 365 ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr); 366 ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 367 ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 368 ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr); 369 ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 370 371 ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr); 372 ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr); 373 ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr); 374 ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr); 375 PetscOptionsEnd(); 376 } 377 378 /* define matrix A */ 379 switch (lu->id.ICNTL(18)){ 380 case 0: /* centralized assembled matrix input (size=1) */ 381 if (!lu->myid) { 382 if (isSeqAIJ){ 383 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 384 nz = aa->nz; 385 ai = aa->i; aj = aa->j; lu->val = aa->a; 386 } else if (isSeqSBAIJ) { 387 Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 388 nz = aa->nz; 389 ai = aa->i; aj = aa->j; lu->val = aa->a; 390 } else { 391 SETERRQ(PETSC_ERR_SUP,"No mumps factorization for this matrix type"); 392 } 393 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 394 ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 395 ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 396 nz = 0; 397 for (i=0; i<M; i++){ 398 rnz = ai[i+1] - ai[i]; 399 while (rnz--) { /* Fortran row/col index! */ 400 lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 401 } 402 } 403 } 404 } 405 break; 406 case 3: /* distributed assembled matrix input (size>1) */ 407 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 408 valOnly = PETSC_FALSE; 409 } else { 410 valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 411 } 412 ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 413 break; 414 default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 415 } 416 417 /* analysis phase */ 418 /*----------------*/ 419 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 420 lu->id.job = 1; 421 422 lu->id.n = M; 423 switch (lu->id.ICNTL(18)){ 424 case 0: /* centralized assembled matrix input */ 425 if (!lu->myid) { 426 lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 427 if (lu->id.ICNTL(6)>1){ 428 #if defined(PETSC_USE_COMPLEX) 429 lu->id.a = (mumps_double_complex*)lu->val; 430 #else 431 lu->id.a = lu->val; 432 #endif 433 } 434 } 435 break; 436 case 3: /* distributed assembled matrix input (size>1) */ 437 lu->id.nz_loc = nnz; 438 lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 439 if (lu->id.ICNTL(6)>1) { 440 #if defined(PETSC_USE_COMPLEX) 441 lu->id.a_loc = (mumps_double_complex*)lu->val; 442 #else 443 lu->id.a_loc = lu->val; 444 #endif 445 } 446 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 447 if (!lu->myid){ 448 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 449 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 450 } else { 451 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 452 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 453 } 454 ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 455 ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 456 ierr = VecSetFromOptions(b);CHKERRQ(ierr); 457 458 ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 459 ierr = ISDestroy(is_iden);CHKERRQ(ierr); 460 ierr = VecDestroy(b);CHKERRQ(ierr); 461 break; 462 } 463 #if defined(PETSC_USE_COMPLEX) 464 zmumps_c(&lu->id); 465 #else 466 dmumps_c(&lu->id); 467 #endif 468 if (lu->id.INFOG(1) < 0) { 469 SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 470 } 471 } 472 473 /* numerical factorization phase */ 474 /*-------------------------------*/ 475 lu->id.job = 2; 476 if(!lu->id.ICNTL(18)) { 477 if (!lu->myid) { 478 #if defined(PETSC_USE_COMPLEX) 479 lu->id.a = (mumps_double_complex*)lu->val; 480 #else 481 lu->id.a = lu->val; 482 #endif 483 } 484 } else { 485 #if defined(PETSC_USE_COMPLEX) 486 lu->id.a_loc = (mumps_double_complex*)lu->val; 487 #else 488 lu->id.a_loc = lu->val; 489 #endif 490 } 491 #if defined(PETSC_USE_COMPLEX) 492 zmumps_c(&lu->id); 493 #else 494 dmumps_c(&lu->id); 495 #endif 496 if (lu->id.INFOG(1) < 0) { 497 if (lu->id.INFO(1) == -13) { 498 SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 499 } else { 500 SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2)); 501 } 502 } 503 504 if (!lu->myid && lu->id.ICNTL(16) > 0){ 505 SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 506 } 507 508 if (lu->size > 1){ 509 if ((F)->factor == MAT_FACTOR_LU){ 510 F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 511 } else { 512 F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 513 } 514 F_diag->assembled = PETSC_TRUE; 515 if (lu->nSolve){ 516 ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 517 ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 518 ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 519 } 520 } 521 (F)->assembled = PETSC_TRUE; 522 lu->matstruc = SAME_NONZERO_PATTERN; 523 lu->CleanUpMUMPS = PETSC_TRUE; 524 lu->nSolve = 0; 525 PetscFunctionReturn(0); 526 } 527 528 529 /* Note the Petsc r and c permutations are ignored */ 530 #undef __FUNCT__ 531 #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 532 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 533 { 534 Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 535 536 PetscFunctionBegin; 537 lu->sym = 0; 538 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 539 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 540 PetscFunctionReturn(0); 541 } 542 543 544 /* Note the Petsc r permutation is ignored */ 545 #undef __FUNCT__ 546 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 547 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 548 { 549 Mat_MUMPS *lu = (Mat_MUMPS*)(F)->spptr; 550 551 PetscFunctionBegin; 552 lu->sym = 2; 553 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 554 (F)->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 555 #if !defined(PETSC_USE_COMPLEX) 556 (F)->ops->getinertia = MatGetInertia_SBAIJMUMPS; 557 #endif 558 PetscFunctionReturn(0); 559 } 560 561 #undef __FUNCT__ 562 #define __FUNCT__ "MatFactorInfo_MUMPS" 563 PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 564 Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 565 PetscErrorCode ierr; 566 567 PetscFunctionBegin; 568 /* check if matrix is mumps type */ 569 if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 570 571 ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 572 ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 573 ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 574 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 575 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 576 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 577 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 578 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 579 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 580 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 581 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 582 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 583 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 584 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 585 if (!lu->myid && lu->id.ICNTL(11)>0) { 586 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 587 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 588 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 589 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr); 590 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 591 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 592 593 } 594 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 595 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 596 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 597 /* ICNTL(15-17) not used */ 598 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 599 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 600 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 601 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 602 603 ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 604 ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 605 ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 606 ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 607 608 /* infomation local to each processor */ 609 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);} 610 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 611 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 612 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);} 613 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 614 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 615 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);} 616 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 617 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 618 /* 619 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);} 620 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr); 621 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 622 */ 623 624 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);} 625 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 626 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 627 628 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);} 629 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 630 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 631 632 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);} 633 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 634 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 635 636 if (!lu->myid){ /* information from the host */ 637 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 638 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 639 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 640 641 ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 642 ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 643 ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 644 ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 645 ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 646 ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 647 ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 648 ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 649 ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 650 ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 651 ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 652 ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 653 ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 654 ierr = PetscViewerASCIIPrintf(viewer," INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr); 655 ierr = PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));CHKERRQ(ierr); 656 ierr = PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));CHKERRQ(ierr); 657 ierr = PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));CHKERRQ(ierr); 658 ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 659 ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr); 660 ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr); 661 ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 662 ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 663 ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 664 } 665 666 PetscFunctionReturn(0); 667 } 668 669 #undef __FUNCT__ 670 #define __FUNCT__ "MatView_MUMPS" 671 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 672 { 673 PetscErrorCode ierr; 674 PetscTruth iascii; 675 PetscViewerFormat format; 676 677 PetscFunctionBegin; 678 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 679 if (iascii) { 680 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 681 if (format == PETSC_VIEWER_ASCII_INFO){ 682 ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 683 } 684 } 685 PetscFunctionReturn(0); 686 } 687 688 689 /*MC 690 MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 691 and sequential matrices via the external package MUMPS. 692 693 If MUMPS is installed (see the manual for instructions 694 on how to declare the existence of external packages), 695 a matrix type can be constructed which invokes MUMPS solvers. 696 After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS), then 697 optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT 698 call MatCreateSeqAIJ/MPIAIJ() directly or the preallocation information will be LOST! 699 700 If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 701 Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 702 MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported 703 for communicators controlling multiple processes. It is recommended that you call both of 704 the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 705 conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 706 without data copy AFTER the matrix values are set. 707 708 Options Database Keys: 709 + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 710 . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 711 . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 712 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 713 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 714 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 715 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 716 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 717 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 718 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 719 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 720 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 721 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 722 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 723 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 724 725 Level: beginner 726 727 .seealso: MATSBAIJMUMPS 728 M*/ 729 730 731 #undef __FUNCT__ 732 #define __FUNCT__ "MatGetInfo_MUMPS" 733 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 734 { 735 Mat_MUMPS *lu =(Mat_MUMPS*)A->spptr; 736 737 PetscFunctionBegin; 738 info->block_size = 1.0; 739 info->nz_allocated = lu->id.INFOG(20); 740 info->nz_used = lu->id.INFOG(20); 741 info->nz_unneeded = 0.0; 742 info->assemblies = 0.0; 743 info->mallocs = 0.0; 744 info->memory = 0.0; 745 info->fill_ratio_given = 0; 746 info->fill_ratio_needed = 0; 747 info->factor_mallocs = 0; 748 PetscFunctionReturn(0); 749 } 750 751 /*MC 752 MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 753 distributed and sequential matrices via the external package MUMPS. 754 755 If MUMPS is installed (see the manual for instructions 756 on how to declare the existence of external packages), 757 a matrix type can be constructed which invokes MUMPS solvers. 758 After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS), then 759 optionally call MatSeqSBAIJSetPreallocation() or MatMPISBAIJSetPreallocation() DO NOT 760 call MatCreateSeqSBAIJ/MPISBAIJ() directly or the preallocation information will be LOST! 761 762 If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 763 Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 764 MatSeqSBAIJSetPreallocation() is supported, and similarly MatMPISBAIJSetPreallocation() is supported 765 for communicators controlling multiple processes. It is recommended that you call both of 766 the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 767 conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 768 without data copy AFTER the matrix values have been set. 769 770 Options Database Keys: 771 + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 772 . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 773 . -mat_mumps_icntl_4 <0,...,4> - print level 774 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 775 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 776 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 777 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 778 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 779 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 780 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 781 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 782 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 783 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 784 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 785 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 786 787 Level: beginner 788 789 .seealso: MATAIJMUMPS 790 M*/ 791 792 EXTERN_C_BEGIN 793 #undef __FUNCT__ 794 #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 795 PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 796 { 797 PetscFunctionBegin; 798 *type = MAT_SOLVER_MUMPS; 799 PetscFunctionReturn(0); 800 } 801 EXTERN_C_END 802 803 EXTERN_C_BEGIN 804 /* 805 The seq and mpi versions of this function are the same 806 */ 807 #undef __FUNCT__ 808 #define __FUNCT__ "MatGetFactor_seqaij_mumps" 809 PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F) 810 { 811 Mat B; 812 PetscErrorCode ierr; 813 Mat_MUMPS *mumps; 814 815 PetscFunctionBegin; 816 if (ftype != MAT_FACTOR_LU) { 817 SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix"); 818 } 819 /* Create the factorization matrix */ 820 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 821 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 822 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 823 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 824 825 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 826 B->ops->view = MatView_MUMPS; 827 B->ops->getinfo = MatGetInfo_MUMPS; 828 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 829 B->factor = MAT_FACTOR_LU; 830 831 ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 832 mumps->CleanUpMUMPS = PETSC_FALSE; 833 mumps->isAIJ = PETSC_TRUE; 834 mumps->scat_rhs = PETSC_NULL; 835 mumps->scat_sol = PETSC_NULL; 836 mumps->nSolve = 0; 837 mumps->MatDestroy = B->ops->destroy; 838 B->ops->destroy = MatDestroy_MUMPS; 839 B->spptr = (void*)mumps; 840 841 *F = B; 842 PetscFunctionReturn(0); 843 } 844 EXTERN_C_END 845 846 EXTERN_C_BEGIN 847 #undef __FUNCT__ 848 #define __FUNCT__ "MatGetFactor_mpiaij_mumps" 849 PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F) 850 { 851 Mat B; 852 PetscErrorCode ierr; 853 Mat_MUMPS *mumps; 854 855 PetscFunctionBegin; 856 if (ftype != MAT_FACTOR_LU) { 857 SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix"); 858 } 859 /* Create the factorization matrix */ 860 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 861 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 862 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 863 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 864 ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 865 866 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 867 B->ops->view = MatView_MUMPS; 868 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 869 B->factor = MAT_FACTOR_LU; 870 871 ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 872 mumps->CleanUpMUMPS = PETSC_FALSE; 873 mumps->isAIJ = PETSC_TRUE; 874 mumps->scat_rhs = PETSC_NULL; 875 mumps->scat_sol = PETSC_NULL; 876 mumps->nSolve = 0; 877 mumps->MatDestroy = B->ops->destroy; 878 B->ops->destroy = MatDestroy_MUMPS; 879 B->spptr = (void*)mumps; 880 881 *F = B; 882 PetscFunctionReturn(0); 883 } 884 EXTERN_C_END 885 886 EXTERN_C_BEGIN 887 #undef __FUNCT__ 888 #define __FUNCT__ "MatGetFactor_seqsbaij_mumps" 889 PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 890 { 891 Mat B; 892 PetscErrorCode ierr; 893 Mat_MUMPS *mumps; 894 895 PetscFunctionBegin; 896 if (ftype != MAT_FACTOR_CHOLESKY) { 897 SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 898 } 899 /* Create the factorization matrix */ 900 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 901 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 902 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 903 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 904 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 905 906 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 907 B->ops->view = MatView_MUMPS; 908 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 909 910 B->factor = MAT_FACTOR_CHOLESKY; 911 912 ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 913 mumps->CleanUpMUMPS = PETSC_FALSE; 914 mumps->isAIJ = PETSC_TRUE; 915 mumps->scat_rhs = PETSC_NULL; 916 mumps->scat_sol = PETSC_NULL; 917 mumps->nSolve = 0; 918 mumps->MatDestroy = B->ops->destroy; 919 B->ops->destroy = MatDestroy_MUMPS; 920 B->spptr = (void*)mumps; 921 922 *F = B; 923 PetscFunctionReturn(0); 924 } 925 EXTERN_C_END 926 927 EXTERN_C_BEGIN 928 #undef __FUNCT__ 929 #define __FUNCT__ "MatGetFactor_mpisbaij_mumps" 930 PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 931 { 932 Mat B; 933 PetscErrorCode ierr; 934 Mat_MUMPS *mumps; 935 936 PetscFunctionBegin; 937 if (ftype != MAT_FACTOR_CHOLESKY) { 938 SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 939 } 940 /* Create the factorization matrix */ 941 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 942 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 943 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 944 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 945 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 946 947 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 948 B->ops->view = MatView_MUMPS; 949 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 950 B->factor = MAT_FACTOR_CHOLESKY; 951 952 ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 953 mumps->CleanUpMUMPS = PETSC_FALSE; 954 mumps->isAIJ = PETSC_TRUE; 955 mumps->scat_rhs = PETSC_NULL; 956 mumps->scat_sol = PETSC_NULL; 957 mumps->nSolve = 0; 958 mumps->MatDestroy = B->ops->destroy; 959 B->ops->destroy = MatDestroy_MUMPS; 960 B->spptr = (void*)mumps; 961 962 *F = B; 963 PetscFunctionReturn(0); 964 } 965 EXTERN_C_END 966