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