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