1 2 /* 3 Provides an interface to the MUMPS sparse solver 4 */ 5 6 #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/ 7 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 8 #include <../src/mat/impls/sell/mpi/mpisell.h> 9 10 EXTERN_C_BEGIN 11 #if defined(PETSC_USE_COMPLEX) 12 #if defined(PETSC_USE_REAL_SINGLE) 13 #include <cmumps_c.h> 14 #else 15 #include <zmumps_c.h> 16 #endif 17 #else 18 #if defined(PETSC_USE_REAL_SINGLE) 19 #include <smumps_c.h> 20 #else 21 #include <dmumps_c.h> 22 #endif 23 #endif 24 EXTERN_C_END 25 #define JOB_INIT -1 26 #define JOB_FACTSYMBOLIC 1 27 #define JOB_FACTNUMERIC 2 28 #define JOB_SOLVE 3 29 #define JOB_END -2 30 31 /* calls to MUMPS */ 32 #if defined(PETSC_USE_COMPLEX) 33 #if defined(PETSC_USE_REAL_SINGLE) 34 #define PetscMUMPS_c cmumps_c 35 #else 36 #define PetscMUMPS_c zmumps_c 37 #endif 38 #else 39 #if defined(PETSC_USE_REAL_SINGLE) 40 #define PetscMUMPS_c smumps_c 41 #else 42 #define PetscMUMPS_c dmumps_c 43 #endif 44 #endif 45 46 /* declare MumpsScalar */ 47 #if defined(PETSC_USE_COMPLEX) 48 #if defined(PETSC_USE_REAL_SINGLE) 49 #define MumpsScalar mumps_complex 50 #else 51 #define MumpsScalar mumps_double_complex 52 #endif 53 #else 54 #define MumpsScalar PetscScalar 55 #endif 56 57 /* macros s.t. indices match MUMPS documentation */ 58 #define ICNTL(I) icntl[(I)-1] 59 #define CNTL(I) cntl[(I)-1] 60 #define INFOG(I) infog[(I)-1] 61 #define INFO(I) info[(I)-1] 62 #define RINFOG(I) rinfog[(I)-1] 63 #define RINFO(I) rinfo[(I)-1] 64 65 typedef struct { 66 #if defined(PETSC_USE_COMPLEX) 67 #if defined(PETSC_USE_REAL_SINGLE) 68 CMUMPS_STRUC_C id; 69 #else 70 ZMUMPS_STRUC_C id; 71 #endif 72 #else 73 #if defined(PETSC_USE_REAL_SINGLE) 74 SMUMPS_STRUC_C id; 75 #else 76 DMUMPS_STRUC_C id; 77 #endif 78 #endif 79 80 MatStructure matstruc; 81 PetscMPIInt myid,size; 82 PetscInt *irn,*jcn,nz,sym; 83 PetscScalar *val; 84 MPI_Comm comm_mumps; 85 PetscInt ICNTL9_pre; /* check if ICNTL(9) is changed from previous MatSolve */ 86 VecScatter scat_rhs, scat_sol; /* used by MatSolve() */ 87 Vec b_seq,x_seq; 88 PetscInt ninfo,*info; /* display INFO */ 89 PetscInt sizeredrhs; 90 PetscScalar *schur_sol; 91 PetscInt schur_sizesol; 92 93 PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**); 94 } Mat_MUMPS; 95 96 extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 97 98 static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps) 99 { 100 PetscErrorCode ierr; 101 102 PetscFunctionBegin; 103 ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr); 104 ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr); 105 ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr); 106 mumps->id.size_schur = 0; 107 mumps->id.schur_lld = 0; 108 mumps->id.ICNTL(19) = 0; 109 PetscFunctionReturn(0); 110 } 111 112 /* solve with rhs in mumps->id.redrhs and return in the same location */ 113 static PetscErrorCode MatMumpsSolveSchur_Private(Mat F) 114 { 115 Mat_MUMPS *mumps=(Mat_MUMPS*)F->data; 116 Mat S,B,X; 117 MatFactorSchurStatus schurstatus; 118 PetscInt sizesol; 119 PetscErrorCode ierr; 120 121 PetscFunctionBegin; 122 ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr); 123 ierr = MatFactorGetSchurComplement(F,&S,&schurstatus);CHKERRQ(ierr); 124 ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B);CHKERRQ(ierr); 125 switch (schurstatus) { 126 case MAT_FACTOR_SCHUR_FACTORED: 127 ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X);CHKERRQ(ierr); 128 if (!mumps->id.ICNTL(9)) { /* transpose solve */ 129 ierr = MatMatSolveTranspose(S,B,X);CHKERRQ(ierr); 130 } else { 131 ierr = MatMatSolve(S,B,X);CHKERRQ(ierr); 132 } 133 break; 134 case MAT_FACTOR_SCHUR_INVERTED: 135 sizesol = mumps->id.nrhs*mumps->id.size_schur; 136 if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) { 137 ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr); 138 ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr); 139 mumps->schur_sizesol = sizesol; 140 } 141 ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X);CHKERRQ(ierr); 142 if (!mumps->id.ICNTL(9)) { /* transpose solve */ 143 ierr = MatTransposeMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 144 } else { 145 ierr = MatMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 146 } 147 ierr = MatCopy(X,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 148 break; 149 default: 150 SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status); 151 break; 152 } 153 ierr = MatFactorRestoreSchurComplement(F,&S,schurstatus);CHKERRQ(ierr); 154 ierr = MatDestroy(&B);CHKERRQ(ierr); 155 ierr = MatDestroy(&X);CHKERRQ(ierr); 156 PetscFunctionReturn(0); 157 } 158 159 static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion) 160 { 161 Mat_MUMPS *mumps=(Mat_MUMPS*)F->data; 162 PetscErrorCode ierr; 163 164 PetscFunctionBegin; 165 if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */ 166 PetscFunctionReturn(0); 167 } 168 if (!expansion) { /* prepare for the condensation step */ 169 PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur; 170 /* allocate MUMPS internal array to store reduced right-hand sides */ 171 if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) { 172 ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr); 173 mumps->id.lredrhs = mumps->id.size_schur; 174 ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr); 175 mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs; 176 } 177 mumps->id.ICNTL(26) = 1; /* condensation phase */ 178 } else { /* prepare for the expansion step */ 179 /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */ 180 ierr = MatMumpsSolveSchur_Private(F);CHKERRQ(ierr); 181 mumps->id.ICNTL(26) = 2; /* expansion phase */ 182 PetscMUMPS_c(&mumps->id); 183 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 184 /* restore defaults */ 185 mumps->id.ICNTL(26) = -1; 186 /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */ 187 if (mumps->id.nrhs > 1) { 188 ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr); 189 mumps->id.lredrhs = 0; 190 mumps->sizeredrhs = 0; 191 } 192 } 193 PetscFunctionReturn(0); 194 } 195 196 /* 197 MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz] 198 199 input: 200 A - matrix in aij,baij or sbaij (bs=1) format 201 shift - 0: C style output triple; 1: Fortran style output triple. 202 reuse - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple 203 MAT_REUSE_MATRIX: only the values in v array are updated 204 output: 205 nnz - dim of r, c, and v (number of local nonzero entries of A) 206 r, c, v - row and col index, matrix values (matrix triples) 207 208 The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is 209 freed with PetscFree(mumps->irn); This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means 210 that the PetscMalloc() cannot easily be replaced with a PetscMalloc3(). 211 212 */ 213 214 PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 215 { 216 const PetscInt *ai,*aj,*ajj,M=A->rmap->n; 217 PetscInt nz,rnz,i,j; 218 PetscErrorCode ierr; 219 PetscInt *row,*col; 220 Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 221 222 PetscFunctionBegin; 223 *v=aa->a; 224 if (reuse == MAT_INITIAL_MATRIX) { 225 nz = aa->nz; 226 ai = aa->i; 227 aj = aa->j; 228 *nnz = nz; 229 ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 230 col = row + nz; 231 232 nz = 0; 233 for (i=0; i<M; i++) { 234 rnz = ai[i+1] - ai[i]; 235 ajj = aj + ai[i]; 236 for (j=0; j<rnz; j++) { 237 row[nz] = i+shift; col[nz++] = ajj[j] + shift; 238 } 239 } 240 *r = row; *c = col; 241 } 242 PetscFunctionReturn(0); 243 } 244 245 PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 246 { 247 Mat_SeqSELL *a=(Mat_SeqSELL*)A->data; 248 PetscInt *ptr; 249 250 PetscFunctionBegin; 251 *v = a->val; 252 if (reuse == MAT_INITIAL_MATRIX) { 253 PetscInt nz,i,j,row; 254 PetscErrorCode ierr; 255 256 nz = a->sliidx[a->totalslices]; 257 *nnz = nz; 258 ierr = PetscMalloc1(2*nz, &ptr);CHKERRQ(ierr); 259 *r = ptr; 260 *c = ptr + nz; 261 262 for (i=0; i<a->totalslices; i++) { 263 for (j=a->sliidx[i],row=0; j<a->sliidx[i+1]; j++,row=((row+1)&0x07)) { 264 *ptr++ = 8*i + row + shift; 265 } 266 } 267 for (i=0;i<nz;i++) *ptr++ = a->colidx[i] + shift; 268 } 269 PetscFunctionReturn(0); 270 } 271 272 PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 273 { 274 Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 275 const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2; 276 PetscInt bs,M,nz,idx=0,rnz,i,j,k,m; 277 PetscErrorCode ierr; 278 PetscInt *row,*col; 279 280 PetscFunctionBegin; 281 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 282 M = A->rmap->N/bs; 283 *v = aa->a; 284 if (reuse == MAT_INITIAL_MATRIX) { 285 ai = aa->i; aj = aa->j; 286 nz = bs2*aa->nz; 287 *nnz = nz; 288 ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 289 col = row + nz; 290 291 for (i=0; i<M; i++) { 292 ajj = aj + ai[i]; 293 rnz = ai[i+1] - ai[i]; 294 for (k=0; k<rnz; k++) { 295 for (j=0; j<bs; j++) { 296 for (m=0; m<bs; m++) { 297 row[idx] = i*bs + m + shift; 298 col[idx++] = bs*(ajj[k]) + j + shift; 299 } 300 } 301 } 302 } 303 *r = row; *c = col; 304 } 305 PetscFunctionReturn(0); 306 } 307 308 PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 309 { 310 const PetscInt *ai, *aj,*ajj,M=A->rmap->n; 311 PetscInt nz,rnz,i,j; 312 PetscErrorCode ierr; 313 PetscInt *row,*col; 314 Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)A->data; 315 316 PetscFunctionBegin; 317 *v = aa->a; 318 if (reuse == MAT_INITIAL_MATRIX) { 319 nz = aa->nz; 320 ai = aa->i; 321 aj = aa->j; 322 *v = aa->a; 323 *nnz = nz; 324 ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr); 325 col = row + nz; 326 327 nz = 0; 328 for (i=0; i<M; i++) { 329 rnz = ai[i+1] - ai[i]; 330 ajj = aj + ai[i]; 331 for (j=0; j<rnz; j++) { 332 row[nz] = i+shift; col[nz++] = ajj[j] + shift; 333 } 334 } 335 *r = row; *c = col; 336 } 337 PetscFunctionReturn(0); 338 } 339 340 PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 341 { 342 const PetscInt *ai,*aj,*ajj,*adiag,M=A->rmap->n; 343 PetscInt nz,rnz,i,j; 344 const PetscScalar *av,*v1; 345 PetscScalar *val; 346 PetscErrorCode ierr; 347 PetscInt *row,*col; 348 Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 349 PetscBool missing; 350 351 PetscFunctionBegin; 352 ai = aa->i; aj = aa->j; av = aa->a; 353 adiag = aa->diag; 354 ierr = MatMissingDiagonal_SeqAIJ(A,&missing,&i);CHKERRQ(ierr); 355 if (reuse == MAT_INITIAL_MATRIX) { 356 /* count nz in the upper triangular part of A */ 357 nz = 0; 358 if (missing) { 359 for (i=0; i<M; i++) { 360 if (PetscUnlikely(adiag[i] >= ai[i+1])) { 361 for (j=ai[i];j<ai[i+1];j++) { 362 if (aj[j] < i) continue; 363 nz++; 364 } 365 } else { 366 nz += ai[i+1] - adiag[i]; 367 } 368 } 369 } else { 370 for (i=0; i<M; i++) nz += ai[i+1] - adiag[i]; 371 } 372 *nnz = nz; 373 374 ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 375 col = row + nz; 376 val = (PetscScalar*)(col + nz); 377 378 nz = 0; 379 if (missing) { 380 for (i=0; i<M; i++) { 381 if (PetscUnlikely(adiag[i] >= ai[i+1])) { 382 for (j=ai[i];j<ai[i+1];j++) { 383 if (aj[j] < i) continue; 384 row[nz] = i+shift; 385 col[nz] = aj[j]+shift; 386 val[nz] = av[j]; 387 nz++; 388 } 389 } else { 390 rnz = ai[i+1] - adiag[i]; 391 ajj = aj + adiag[i]; 392 v1 = av + adiag[i]; 393 for (j=0; j<rnz; j++) { 394 row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 395 } 396 } 397 } 398 } else { 399 for (i=0; i<M; i++) { 400 rnz = ai[i+1] - adiag[i]; 401 ajj = aj + adiag[i]; 402 v1 = av + adiag[i]; 403 for (j=0; j<rnz; j++) { 404 row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j]; 405 } 406 } 407 } 408 *r = row; *c = col; *v = val; 409 } else { 410 nz = 0; val = *v; 411 if (missing) { 412 for (i=0; i <M; i++) { 413 if (PetscUnlikely(adiag[i] >= ai[i+1])) { 414 for (j=ai[i];j<ai[i+1];j++) { 415 if (aj[j] < i) continue; 416 val[nz++] = av[j]; 417 } 418 } else { 419 rnz = ai[i+1] - adiag[i]; 420 v1 = av + adiag[i]; 421 for (j=0; j<rnz; j++) { 422 val[nz++] = v1[j]; 423 } 424 } 425 } 426 } else { 427 for (i=0; i <M; i++) { 428 rnz = ai[i+1] - adiag[i]; 429 v1 = av + adiag[i]; 430 for (j=0; j<rnz; j++) { 431 val[nz++] = v1[j]; 432 } 433 } 434 } 435 } 436 PetscFunctionReturn(0); 437 } 438 439 PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 440 { 441 const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 442 PetscErrorCode ierr; 443 PetscInt rstart,nz,i,j,jj,irow,countA,countB; 444 PetscInt *row,*col; 445 const PetscScalar *av, *bv,*v1,*v2; 446 PetscScalar *val; 447 Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 448 Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)(mat->A)->data; 449 Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 450 451 PetscFunctionBegin; 452 ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 453 av=aa->a; bv=bb->a; 454 455 garray = mat->garray; 456 457 if (reuse == MAT_INITIAL_MATRIX) { 458 nz = aa->nz + bb->nz; 459 *nnz = nz; 460 ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 461 col = row + nz; 462 val = (PetscScalar*)(col + nz); 463 464 *r = row; *c = col; *v = val; 465 } else { 466 row = *r; col = *c; val = *v; 467 } 468 469 jj = 0; irow = rstart; 470 for (i=0; i<m; i++) { 471 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 472 countA = ai[i+1] - ai[i]; 473 countB = bi[i+1] - bi[i]; 474 bjj = bj + bi[i]; 475 v1 = av + ai[i]; 476 v2 = bv + bi[i]; 477 478 /* A-part */ 479 for (j=0; j<countA; j++) { 480 if (reuse == MAT_INITIAL_MATRIX) { 481 row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 482 } 483 val[jj++] = v1[j]; 484 } 485 486 /* B-part */ 487 for (j=0; j < countB; j++) { 488 if (reuse == MAT_INITIAL_MATRIX) { 489 row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 490 } 491 val[jj++] = v2[j]; 492 } 493 irow++; 494 } 495 PetscFunctionReturn(0); 496 } 497 498 PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 499 { 500 const PetscInt *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 501 PetscErrorCode ierr; 502 PetscInt rstart,nz,i,j,jj,irow,countA,countB; 503 PetscInt *row,*col; 504 const PetscScalar *av, *bv,*v1,*v2; 505 PetscScalar *val; 506 Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 507 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(mat->A)->data; 508 Mat_SeqAIJ *bb = (Mat_SeqAIJ*)(mat->B)->data; 509 510 PetscFunctionBegin; 511 ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 512 av=aa->a; bv=bb->a; 513 514 garray = mat->garray; 515 516 if (reuse == MAT_INITIAL_MATRIX) { 517 nz = aa->nz + bb->nz; 518 *nnz = nz; 519 ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 520 col = row + nz; 521 val = (PetscScalar*)(col + nz); 522 523 *r = row; *c = col; *v = val; 524 } else { 525 row = *r; col = *c; val = *v; 526 } 527 528 jj = 0; irow = rstart; 529 for (i=0; i<m; i++) { 530 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 531 countA = ai[i+1] - ai[i]; 532 countB = bi[i+1] - bi[i]; 533 bjj = bj + bi[i]; 534 v1 = av + ai[i]; 535 v2 = bv + bi[i]; 536 537 /* A-part */ 538 for (j=0; j<countA; j++) { 539 if (reuse == MAT_INITIAL_MATRIX) { 540 row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 541 } 542 val[jj++] = v1[j]; 543 } 544 545 /* B-part */ 546 for (j=0; j < countB; j++) { 547 if (reuse == MAT_INITIAL_MATRIX) { 548 row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 549 } 550 val[jj++] = v2[j]; 551 } 552 irow++; 553 } 554 PetscFunctionReturn(0); 555 } 556 557 PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 558 { 559 Mat_MPIBAIJ *mat = (Mat_MPIBAIJ*)A->data; 560 Mat_SeqBAIJ *aa = (Mat_SeqBAIJ*)(mat->A)->data; 561 Mat_SeqBAIJ *bb = (Mat_SeqBAIJ*)(mat->B)->data; 562 const PetscInt *ai = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj; 563 const PetscInt *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart; 564 const PetscInt bs2=mat->bs2; 565 PetscErrorCode ierr; 566 PetscInt bs,nz,i,j,k,n,jj,irow,countA,countB,idx; 567 PetscInt *row,*col; 568 const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2; 569 PetscScalar *val; 570 571 PetscFunctionBegin; 572 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 573 if (reuse == MAT_INITIAL_MATRIX) { 574 nz = bs2*(aa->nz + bb->nz); 575 *nnz = nz; 576 ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 577 col = row + nz; 578 val = (PetscScalar*)(col + nz); 579 580 *r = row; *c = col; *v = val; 581 } else { 582 row = *r; col = *c; val = *v; 583 } 584 585 jj = 0; irow = rstart; 586 for (i=0; i<mbs; i++) { 587 countA = ai[i+1] - ai[i]; 588 countB = bi[i+1] - bi[i]; 589 ajj = aj + ai[i]; 590 bjj = bj + bi[i]; 591 v1 = av + bs2*ai[i]; 592 v2 = bv + bs2*bi[i]; 593 594 idx = 0; 595 /* A-part */ 596 for (k=0; k<countA; k++) { 597 for (j=0; j<bs; j++) { 598 for (n=0; n<bs; n++) { 599 if (reuse == MAT_INITIAL_MATRIX) { 600 row[jj] = irow + n + shift; 601 col[jj] = rstart + bs*ajj[k] + j + shift; 602 } 603 val[jj++] = v1[idx++]; 604 } 605 } 606 } 607 608 idx = 0; 609 /* B-part */ 610 for (k=0; k<countB; k++) { 611 for (j=0; j<bs; j++) { 612 for (n=0; n<bs; n++) { 613 if (reuse == MAT_INITIAL_MATRIX) { 614 row[jj] = irow + n + shift; 615 col[jj] = bs*garray[bjj[k]] + j + shift; 616 } 617 val[jj++] = v2[idx++]; 618 } 619 } 620 } 621 irow += bs; 622 } 623 PetscFunctionReturn(0); 624 } 625 626 PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v) 627 { 628 const PetscInt *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj; 629 PetscErrorCode ierr; 630 PetscInt rstart,nz,nza,nzb,i,j,jj,irow,countA,countB; 631 PetscInt *row,*col; 632 const PetscScalar *av, *bv,*v1,*v2; 633 PetscScalar *val; 634 Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 635 Mat_SeqAIJ *aa =(Mat_SeqAIJ*)(mat->A)->data; 636 Mat_SeqAIJ *bb =(Mat_SeqAIJ*)(mat->B)->data; 637 638 PetscFunctionBegin; 639 ai=aa->i; aj=aa->j; adiag=aa->diag; 640 bi=bb->i; bj=bb->j; garray = mat->garray; 641 av=aa->a; bv=bb->a; 642 643 rstart = A->rmap->rstart; 644 645 if (reuse == MAT_INITIAL_MATRIX) { 646 nza = 0; /* num of upper triangular entries in mat->A, including diagonals */ 647 nzb = 0; /* num of upper triangular entries in mat->B */ 648 for (i=0; i<m; i++) { 649 nza += (ai[i+1] - adiag[i]); 650 countB = bi[i+1] - bi[i]; 651 bjj = bj + bi[i]; 652 for (j=0; j<countB; j++) { 653 if (garray[bjj[j]] > rstart) nzb++; 654 } 655 } 656 657 nz = nza + nzb; /* total nz of upper triangular part of mat */ 658 *nnz = nz; 659 ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr); 660 col = row + nz; 661 val = (PetscScalar*)(col + nz); 662 663 *r = row; *c = col; *v = val; 664 } else { 665 row = *r; col = *c; val = *v; 666 } 667 668 jj = 0; irow = rstart; 669 for (i=0; i<m; i++) { 670 ajj = aj + adiag[i]; /* ptr to the beginning of the diagonal of this row */ 671 v1 = av + adiag[i]; 672 countA = ai[i+1] - adiag[i]; 673 countB = bi[i+1] - bi[i]; 674 bjj = bj + bi[i]; 675 v2 = bv + bi[i]; 676 677 /* A-part */ 678 for (j=0; j<countA; j++) { 679 if (reuse == MAT_INITIAL_MATRIX) { 680 row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 681 } 682 val[jj++] = v1[j]; 683 } 684 685 /* B-part */ 686 for (j=0; j < countB; j++) { 687 if (garray[bjj[j]] > rstart) { 688 if (reuse == MAT_INITIAL_MATRIX) { 689 row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 690 } 691 val[jj++] = v2[j]; 692 } 693 } 694 irow++; 695 } 696 PetscFunctionReturn(0); 697 } 698 699 PetscErrorCode MatDestroy_MUMPS(Mat A) 700 { 701 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 702 PetscErrorCode ierr; 703 704 PetscFunctionBegin; 705 ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 706 ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr); 707 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 708 ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr); 709 ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 710 ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr); 711 ierr = PetscFree(mumps->irn);CHKERRQ(ierr); 712 ierr = PetscFree(mumps->info);CHKERRQ(ierr); 713 ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr); 714 mumps->id.job = JOB_END; 715 PetscMUMPS_c(&mumps->id); 716 ierr = MPI_Comm_free(&mumps->comm_mumps);CHKERRQ(ierr); 717 ierr = PetscFree(A->data);CHKERRQ(ierr); 718 719 /* clear composed functions */ 720 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr); 721 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr); 722 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr); 723 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr); 724 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr); 725 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr); 726 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr); 727 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr); 728 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr); 729 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr); 730 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr); 731 PetscFunctionReturn(0); 732 } 733 734 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 735 { 736 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 737 PetscScalar *array; 738 Vec b_seq; 739 IS is_iden,is_petsc; 740 PetscErrorCode ierr; 741 PetscInt i; 742 PetscBool second_solve = PETSC_FALSE; 743 static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE; 744 745 PetscFunctionBegin; 746 ierr = PetscCitationsRegister("@article{MUMPS01,\n author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n journal = {SIAM Journal on Matrix Analysis and Applications},\n volume = {23},\n number = {1},\n pages = {15--41},\n year = {2001}\n}\n",&cite1);CHKERRQ(ierr); 747 ierr = PetscCitationsRegister("@article{MUMPS02,\n author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n title = {Hybrid scheduling for the parallel solution of linear systems},\n journal = {Parallel Computing},\n volume = {32},\n number = {2},\n pages = {136--156},\n year = {2006}\n}\n",&cite2);CHKERRQ(ierr); 748 749 if (A->factorerrortype) { 750 ierr = PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 751 ierr = VecSetInf(x);CHKERRQ(ierr); 752 PetscFunctionReturn(0); 753 } 754 755 mumps->id.nrhs = 1; 756 b_seq = mumps->b_seq; 757 if (mumps->size > 1) { 758 /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 759 ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 760 ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 761 if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 762 } else { /* size == 1 */ 763 ierr = VecCopy(b,x);CHKERRQ(ierr); 764 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 765 } 766 if (!mumps->myid) { /* define rhs on the host */ 767 mumps->id.nrhs = 1; 768 mumps->id.rhs = (MumpsScalar*)array; 769 } 770 771 /* 772 handle condensation step of Schur complement (if any) 773 We set by default ICNTL(26) == -1 when Schur indices have been provided by the user. 774 According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase 775 Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system. 776 This requires an extra call to PetscMUMPS_c and the computation of the factors for S 777 */ 778 if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) { 779 if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n"); 780 second_solve = PETSC_TRUE; 781 ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr); 782 } 783 /* solve phase */ 784 /*-------------*/ 785 mumps->id.job = JOB_SOLVE; 786 PetscMUMPS_c(&mumps->id); 787 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 788 789 /* handle expansion step of Schur complement (if any) */ 790 if (second_solve) { 791 ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr); 792 } 793 794 if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 795 if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) { 796 /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */ 797 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 798 } 799 if (!mumps->scat_sol) { /* create scatter scat_sol */ 800 ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 801 for (i=0; i<mumps->id.lsol_loc; i++) { 802 mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 803 } 804 ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 805 ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr); 806 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 807 ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 808 809 mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */ 810 } 811 812 ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 813 ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 814 } 815 PetscFunctionReturn(0); 816 } 817 818 PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 819 { 820 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 821 PetscErrorCode ierr; 822 823 PetscFunctionBegin; 824 mumps->id.ICNTL(9) = 0; 825 ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 826 mumps->id.ICNTL(9) = 1; 827 PetscFunctionReturn(0); 828 } 829 830 PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 831 { 832 PetscErrorCode ierr; 833 Mat Bt = NULL; 834 PetscBool flg, flgT; 835 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 836 PetscInt i,nrhs,M; 837 PetscScalar *array,*bray; 838 839 PetscFunctionBegin; 840 ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 841 ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);CHKERRQ(ierr); 842 if (flgT) { 843 if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 844 ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr); 845 } else { 846 if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 847 } 848 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 849 if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 850 if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution"); 851 852 ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr); 853 mumps->id.nrhs = nrhs; 854 mumps->id.lrhs = M; 855 856 if (mumps->size == 1) { 857 PetscScalar *aa; 858 PetscInt spnr,*ia,*ja; 859 PetscBool second_solve = PETSC_FALSE; 860 861 /* copy B to X */ 862 ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 863 mumps->id.rhs = (MumpsScalar*)array; 864 if (!Bt) { 865 ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 866 ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 867 ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 868 } else { 869 PetscBool done; 870 871 ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr); 872 ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 873 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure"); 874 mumps->id.irhs_ptr = ia; 875 mumps->id.irhs_sparse = ja; 876 mumps->id.nz_rhs = ia[spnr] - 1; 877 mumps->id.rhs_sparse = (MumpsScalar*)aa; 878 mumps->id.ICNTL(20) = 1; 879 } 880 /* handle condensation step of Schur complement (if any) */ 881 if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) { 882 second_solve = PETSC_TRUE; 883 ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr); 884 } 885 /* solve phase */ 886 /*-------------*/ 887 mumps->id.job = JOB_SOLVE; 888 PetscMUMPS_c(&mumps->id); 889 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 890 891 /* handle expansion step of Schur complement (if any) */ 892 if (second_solve) { 893 ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr); 894 } 895 if (Bt) { 896 PetscBool done; 897 898 ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr); 899 ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 900 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure"); 901 mumps->id.ICNTL(20) = 0; 902 } 903 ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 904 } else { /*--------- parallel case --------*/ 905 PetscInt lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save; 906 MumpsScalar *sol_loc,*sol_loc_save; 907 IS is_to,is_from; 908 PetscInt k,proc,j,m; 909 const PetscInt *rstart; 910 Vec v_mpi,b_seq,x_seq; 911 VecScatter scat_rhs,scat_sol; 912 913 if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n"); 914 915 /* create x_seq to hold local solution */ 916 isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */ 917 sol_loc_save = mumps->id.sol_loc; 918 919 lsol_loc = mumps->id.INFO(23); 920 nlsol_loc = nrhs*lsol_loc; /* length of sol_loc */ 921 ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr); 922 mumps->id.sol_loc = (MumpsScalar*)sol_loc; 923 mumps->id.isol_loc = isol_loc; 924 925 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr); 926 927 /* copy rhs matrix B into vector v_mpi */ 928 ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr); 929 ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 930 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr); 931 ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 932 933 /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */ 934 /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B; 935 iidx: inverse of idx, will be used by scattering xx_seq -> X */ 936 ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr); 937 ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr); 938 k = 0; 939 for (proc=0; proc<mumps->size; proc++){ 940 for (j=0; j<nrhs; j++){ 941 for (i=rstart[proc]; i<rstart[proc+1]; i++){ 942 iidx[j*M + i] = k; 943 idx[k++] = j*M + i; 944 } 945 } 946 } 947 948 if (!mumps->myid) { 949 ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr); 950 ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 951 ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr); 952 } else { 953 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr); 954 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr); 955 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr); 956 } 957 ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr); 958 ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 959 ierr = ISDestroy(&is_to);CHKERRQ(ierr); 960 ierr = ISDestroy(&is_from);CHKERRQ(ierr); 961 ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 962 963 if (!mumps->myid) { /* define rhs on the host */ 964 ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr); 965 mumps->id.rhs = (MumpsScalar*)bray; 966 ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr); 967 } 968 969 /* solve phase */ 970 /*-------------*/ 971 mumps->id.job = JOB_SOLVE; 972 PetscMUMPS_c(&mumps->id); 973 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 974 975 /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */ 976 ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 977 ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr); 978 979 /* create scatter scat_sol */ 980 ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr); 981 ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr); 982 for (i=0; i<lsol_loc; i++) { 983 isol_loc[i] -= 1; /* change Fortran style to C style */ 984 idxx[i] = iidx[isol_loc[i]]; 985 for (j=1; j<nrhs; j++){ 986 idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M]; 987 } 988 } 989 ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 990 ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr); 991 ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 992 ierr = ISDestroy(&is_from);CHKERRQ(ierr); 993 ierr = ISDestroy(&is_to);CHKERRQ(ierr); 994 ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 995 ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 996 997 /* free spaces */ 998 mumps->id.sol_loc = sol_loc_save; 999 mumps->id.isol_loc = isol_loc_save; 1000 1001 ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr); 1002 ierr = PetscFree2(idx,iidx);CHKERRQ(ierr); 1003 ierr = PetscFree(idxx);CHKERRQ(ierr); 1004 ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 1005 ierr = VecDestroy(&v_mpi);CHKERRQ(ierr); 1006 ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 1007 ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr); 1008 ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr); 1009 } 1010 PetscFunctionReturn(0); 1011 } 1012 1013 #if !defined(PETSC_USE_COMPLEX) 1014 /* 1015 input: 1016 F: numeric factor 1017 output: 1018 nneg: total number of negative pivots 1019 nzero: total number of zero pivots 1020 npos: (global dimension of F) - nneg - nzero 1021 */ 1022 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 1023 { 1024 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1025 PetscErrorCode ierr; 1026 PetscMPIInt size; 1027 1028 PetscFunctionBegin; 1029 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr); 1030 /* 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 */ 1031 if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13)); 1032 1033 if (nneg) *nneg = mumps->id.INFOG(12); 1034 if (nzero || npos) { 1035 if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection"); 1036 if (nzero) *nzero = mumps->id.INFOG(28); 1037 if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28)); 1038 } 1039 PetscFunctionReturn(0); 1040 } 1041 #endif 1042 1043 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 1044 { 1045 Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->data; 1046 PetscErrorCode ierr; 1047 PetscBool isMPIAIJ; 1048 1049 PetscFunctionBegin; 1050 if (mumps->id.INFOG(1) < 0) { 1051 if (mumps->id.INFOG(1) == -6) { 1052 ierr = PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1053 } 1054 ierr = PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1055 PetscFunctionReturn(0); 1056 } 1057 1058 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1059 1060 /* numerical factorization phase */ 1061 /*-------------------------------*/ 1062 mumps->id.job = JOB_FACTNUMERIC; 1063 if (!mumps->id.ICNTL(18)) { /* A is centralized */ 1064 if (!mumps->myid) { 1065 mumps->id.a = (MumpsScalar*)mumps->val; 1066 } 1067 } else { 1068 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1069 } 1070 PetscMUMPS_c(&mumps->id); 1071 if (mumps->id.INFOG(1) < 0) { 1072 if (A->erroriffailure) { 1073 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)); 1074 } else { 1075 if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */ 1076 ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1077 F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1078 } else if (mumps->id.INFOG(1) == -13) { 1079 ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1080 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 1081 } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) { 1082 ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1083 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 1084 } else { 1085 ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1086 F->factorerrortype = MAT_FACTOR_OTHER; 1087 } 1088 } 1089 } 1090 if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB," mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16)); 1091 1092 F->assembled = PETSC_TRUE; 1093 mumps->matstruc = SAME_NONZERO_PATTERN; 1094 if (F->schur) { /* reset Schur status to unfactored */ 1095 if (mumps->id.ICNTL(19) == 1) { /* stored by rows */ 1096 mumps->id.ICNTL(19) = 2; 1097 ierr = MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);CHKERRQ(ierr); 1098 } 1099 ierr = MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1100 } 1101 1102 /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */ 1103 if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3; 1104 1105 if (mumps->size > 1) { 1106 PetscInt lsol_loc; 1107 PetscScalar *sol_loc; 1108 1109 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 1110 1111 /* distributed solution; Create x_seq=sol_loc for repeated use */ 1112 if (mumps->x_seq) { 1113 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 1114 ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 1115 ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 1116 } 1117 lsol_loc = mumps->id.INFO(23); /* length of sol_loc */ 1118 ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr); 1119 mumps->id.lsol_loc = lsol_loc; 1120 mumps->id.sol_loc = (MumpsScalar*)sol_loc; 1121 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr); 1122 } 1123 PetscFunctionReturn(0); 1124 } 1125 1126 /* Sets MUMPS options from the options database */ 1127 PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 1128 { 1129 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1130 PetscErrorCode ierr; 1131 PetscInt icntl,info[40],i,ninfo=40; 1132 PetscBool flg; 1133 1134 PetscFunctionBegin; 1135 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 1136 ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 1137 if (flg) mumps->id.ICNTL(1) = icntl; 1138 ierr = PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);CHKERRQ(ierr); 1139 if (flg) mumps->id.ICNTL(2) = icntl; 1140 ierr = PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);CHKERRQ(ierr); 1141 if (flg) mumps->id.ICNTL(3) = icntl; 1142 1143 ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 1144 if (flg) mumps->id.ICNTL(4) = icntl; 1145 if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 1146 1147 ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr); 1148 if (flg) mumps->id.ICNTL(6) = icntl; 1149 1150 ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 1151 if (flg) { 1152 if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 1153 else mumps->id.ICNTL(7) = icntl; 1154 } 1155 1156 ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr); 1157 /* ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */ 1158 ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr); 1159 ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr); 1160 ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr); 1161 ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr); 1162 ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr); 1163 ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr); 1164 if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */ 1165 ierr = MatDestroy(&F->schur);CHKERRQ(ierr); 1166 ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr); 1167 } 1168 /* ierr = PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */ 1169 /* ierr = PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL);CHKERRQ(ierr); we only use distributed solution vector */ 1170 1171 ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr); 1172 ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);CHKERRQ(ierr); 1173 ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);CHKERRQ(ierr); 1174 if (mumps->id.ICNTL(24)) { 1175 mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 1176 } 1177 1178 ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computes a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr); 1179 ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr); 1180 ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr); 1181 ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);CHKERRQ(ierr); 1182 ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr); 1183 ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);CHKERRQ(ierr); 1184 ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr); 1185 /* ierr = PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);CHKERRQ(ierr); -- not supported by PETSc API */ 1186 ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr); 1187 ierr = PetscOptionsInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Lock Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);CHKERRQ(ierr); 1188 1189 ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr); 1190 ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr); 1191 ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr); 1192 ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr); 1193 ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr); 1194 ierr = PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);CHKERRQ(ierr); 1195 1196 ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);CHKERRQ(ierr); 1197 1198 ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr); 1199 if (ninfo) { 1200 if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo); 1201 ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr); 1202 mumps->ninfo = ninfo; 1203 for (i=0; i<ninfo; i++) { 1204 if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo); 1205 else mumps->info[i] = info[i]; 1206 } 1207 } 1208 1209 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1210 PetscFunctionReturn(0); 1211 } 1212 1213 PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps) 1214 { 1215 PetscErrorCode ierr; 1216 1217 PetscFunctionBegin; 1218 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);CHKERRQ(ierr); 1219 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr); 1220 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr); 1221 1222 mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 1223 1224 mumps->id.job = JOB_INIT; 1225 mumps->id.par = 1; /* host participates factorizaton and solve */ 1226 mumps->id.sym = mumps->sym; 1227 PetscMUMPS_c(&mumps->id); 1228 1229 mumps->scat_rhs = NULL; 1230 mumps->scat_sol = NULL; 1231 1232 /* set PETSc-MUMPS default options - override MUMPS default */ 1233 mumps->id.ICNTL(3) = 0; 1234 mumps->id.ICNTL(4) = 0; 1235 if (mumps->size == 1) { 1236 mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 1237 } else { 1238 mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 1239 mumps->id.ICNTL(20) = 0; /* rhs is in dense format */ 1240 mumps->id.ICNTL(21) = 1; /* distributed solution */ 1241 } 1242 1243 /* schur */ 1244 mumps->id.size_schur = 0; 1245 mumps->id.listvar_schur = NULL; 1246 mumps->id.schur = NULL; 1247 mumps->sizeredrhs = 0; 1248 mumps->schur_sol = NULL; 1249 mumps->schur_sizesol = 0; 1250 PetscFunctionReturn(0); 1251 } 1252 1253 PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps) 1254 { 1255 PetscErrorCode ierr; 1256 1257 PetscFunctionBegin; 1258 if (mumps->id.INFOG(1) < 0) { 1259 if (A->erroriffailure) { 1260 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 1261 } else { 1262 if (mumps->id.INFOG(1) == -6) { 1263 ierr = PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1264 F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT; 1265 } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) { 1266 ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1267 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 1268 } else { 1269 ierr = PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1270 F->factorerrortype = MAT_FACTOR_OTHER; 1271 } 1272 } 1273 } 1274 PetscFunctionReturn(0); 1275 } 1276 1277 /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */ 1278 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1279 { 1280 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1281 PetscErrorCode ierr; 1282 Vec b; 1283 IS is_iden; 1284 const PetscInt M = A->rmap->N; 1285 1286 PetscFunctionBegin; 1287 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1288 1289 /* Set MUMPS options from the options database */ 1290 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1291 1292 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1293 1294 /* analysis phase */ 1295 /*----------------*/ 1296 mumps->id.job = JOB_FACTSYMBOLIC; 1297 mumps->id.n = M; 1298 switch (mumps->id.ICNTL(18)) { 1299 case 0: /* centralized assembled matrix input */ 1300 if (!mumps->myid) { 1301 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1302 if (mumps->id.ICNTL(6)>1) { 1303 mumps->id.a = (MumpsScalar*)mumps->val; 1304 } 1305 if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */ 1306 /* 1307 PetscBool flag; 1308 ierr = ISEqual(r,c,&flag);CHKERRQ(ierr); 1309 if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm"); 1310 ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF); 1311 */ 1312 if (!mumps->myid) { 1313 const PetscInt *idx; 1314 PetscInt i,*perm_in; 1315 1316 ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr); 1317 ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 1318 1319 mumps->id.perm_in = perm_in; 1320 for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 1321 ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 1322 } 1323 } 1324 } 1325 break; 1326 case 3: /* distributed assembled matrix input (size>1) */ 1327 mumps->id.nz_loc = mumps->nz; 1328 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1329 if (mumps->id.ICNTL(6)>1) { 1330 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1331 } 1332 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1333 if (!mumps->myid) { 1334 ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr); 1335 ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr); 1336 } else { 1337 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1338 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1339 } 1340 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1341 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1342 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1343 ierr = VecDestroy(&b);CHKERRQ(ierr); 1344 break; 1345 } 1346 PetscMUMPS_c(&mumps->id); 1347 ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 1348 1349 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1350 F->ops->solve = MatSolve_MUMPS; 1351 F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1352 F->ops->matsolve = MatMatSolve_MUMPS; 1353 PetscFunctionReturn(0); 1354 } 1355 1356 /* Note the Petsc r and c permutations are ignored */ 1357 PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1358 { 1359 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1360 PetscErrorCode ierr; 1361 Vec b; 1362 IS is_iden; 1363 const PetscInt M = A->rmap->N; 1364 1365 PetscFunctionBegin; 1366 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1367 1368 /* Set MUMPS options from the options database */ 1369 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1370 1371 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1372 1373 /* analysis phase */ 1374 /*----------------*/ 1375 mumps->id.job = JOB_FACTSYMBOLIC; 1376 mumps->id.n = M; 1377 switch (mumps->id.ICNTL(18)) { 1378 case 0: /* centralized assembled matrix input */ 1379 if (!mumps->myid) { 1380 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1381 if (mumps->id.ICNTL(6)>1) { 1382 mumps->id.a = (MumpsScalar*)mumps->val; 1383 } 1384 } 1385 break; 1386 case 3: /* distributed assembled matrix input (size>1) */ 1387 mumps->id.nz_loc = mumps->nz; 1388 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1389 if (mumps->id.ICNTL(6)>1) { 1390 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1391 } 1392 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1393 if (!mumps->myid) { 1394 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 1395 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 1396 } else { 1397 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1398 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1399 } 1400 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1401 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1402 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1403 ierr = VecDestroy(&b);CHKERRQ(ierr); 1404 break; 1405 } 1406 PetscMUMPS_c(&mumps->id); 1407 ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 1408 1409 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1410 F->ops->solve = MatSolve_MUMPS; 1411 F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1412 PetscFunctionReturn(0); 1413 } 1414 1415 /* Note the Petsc r permutation and factor info are ignored */ 1416 PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 1417 { 1418 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1419 PetscErrorCode ierr; 1420 Vec b; 1421 IS is_iden; 1422 const PetscInt M = A->rmap->N; 1423 1424 PetscFunctionBegin; 1425 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1426 1427 /* Set MUMPS options from the options database */ 1428 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1429 1430 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1431 1432 /* analysis phase */ 1433 /*----------------*/ 1434 mumps->id.job = JOB_FACTSYMBOLIC; 1435 mumps->id.n = M; 1436 switch (mumps->id.ICNTL(18)) { 1437 case 0: /* centralized assembled matrix input */ 1438 if (!mumps->myid) { 1439 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1440 if (mumps->id.ICNTL(6)>1) { 1441 mumps->id.a = (MumpsScalar*)mumps->val; 1442 } 1443 } 1444 break; 1445 case 3: /* distributed assembled matrix input (size>1) */ 1446 mumps->id.nz_loc = mumps->nz; 1447 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1448 if (mumps->id.ICNTL(6)>1) { 1449 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1450 } 1451 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1452 if (!mumps->myid) { 1453 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 1454 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 1455 } else { 1456 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1457 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1458 } 1459 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1460 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1461 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1462 ierr = VecDestroy(&b);CHKERRQ(ierr); 1463 break; 1464 } 1465 PetscMUMPS_c(&mumps->id); 1466 ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 1467 1468 F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1469 F->ops->solve = MatSolve_MUMPS; 1470 F->ops->solvetranspose = MatSolve_MUMPS; 1471 F->ops->matsolve = MatMatSolve_MUMPS; 1472 #if defined(PETSC_USE_COMPLEX) 1473 F->ops->getinertia = NULL; 1474 #else 1475 F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1476 #endif 1477 PetscFunctionReturn(0); 1478 } 1479 1480 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 1481 { 1482 PetscErrorCode ierr; 1483 PetscBool iascii; 1484 PetscViewerFormat format; 1485 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 1486 1487 PetscFunctionBegin; 1488 /* check if matrix is mumps type */ 1489 if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 1490 1491 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1492 if (iascii) { 1493 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1494 if (format == PETSC_VIEWER_ASCII_INFO) { 1495 ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1496 ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1497 ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1498 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr); 1499 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr); 1500 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr); 1501 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr); 1502 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr); 1503 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr); 1504 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr); 1505 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scaling strategy): %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr); 1506 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr); 1507 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr); 1508 if (mumps->id.ICNTL(11)>0) { 1509 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr); 1510 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr); 1511 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr); 1512 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr); 1513 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr); 1514 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr); 1515 } 1516 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr); 1517 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr); 1518 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr); 1519 /* ICNTL(15-17) not used */ 1520 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr); 1521 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Schur complement info): %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr); 1522 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr); 1523 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr); 1524 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr); 1525 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr); 1526 1527 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr); 1528 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr); 1529 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr); 1530 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr); 1531 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr); 1532 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr); 1533 1534 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr); 1535 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr); 1536 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr); 1537 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(35) (activate BLR based factorization): %d \n",mumps->id.ICNTL(35));CHKERRQ(ierr); 1538 1539 ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));CHKERRQ(ierr); 1540 ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr); 1541 ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",mumps->id.CNTL(3));CHKERRQ(ierr); 1542 ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",mumps->id.CNTL(4));CHKERRQ(ierr); 1543 ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));CHKERRQ(ierr); 1544 ierr = PetscViewerASCIIPrintf(viewer," CNTL(7) (dropping parameter for BLR): %g \n",mumps->id.CNTL(7));CHKERRQ(ierr); 1545 1546 /* infomation local to each processor */ 1547 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 1548 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 1549 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr); 1550 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1551 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 1552 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr); 1553 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1554 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 1555 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr); 1556 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1557 1558 ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 1559 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr); 1560 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1561 1562 ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 1563 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr); 1564 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1565 1566 ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 1567 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr); 1568 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1569 1570 if (mumps->ninfo && mumps->ninfo <= 40){ 1571 PetscInt i; 1572 for (i=0; i<mumps->ninfo; i++){ 1573 ierr = PetscViewerASCIIPrintf(viewer, " INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr); 1574 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr); 1575 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1576 } 1577 } 1578 1579 1580 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 1581 1582 if (!mumps->myid) { /* information from the host */ 1583 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr); 1584 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr); 1585 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr); 1586 ierr = PetscViewerASCIIPrintf(viewer," (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));CHKERRQ(ierr); 1587 1588 ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr); 1589 ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr); 1590 ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr); 1591 ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr); 1592 ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr); 1593 ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr); 1594 ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr); 1595 ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr); 1596 ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr); 1597 ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr); 1598 ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr); 1599 ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr); 1600 ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr); 1601 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",mumps->id.INFOG(16));CHKERRQ(ierr); 1602 ierr = PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));CHKERRQ(ierr); 1603 ierr = PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));CHKERRQ(ierr); 1604 ierr = PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));CHKERRQ(ierr); 1605 ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr); 1606 ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));CHKERRQ(ierr); 1607 ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));CHKERRQ(ierr); 1608 ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr); 1609 ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr); 1610 ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr); 1611 ierr = PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr); 1612 ierr = PetscViewerASCIIPrintf(viewer," INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr); 1613 ierr = PetscViewerASCIIPrintf(viewer," INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));CHKERRQ(ierr); 1614 ierr = PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr); 1615 ierr = PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr); 1616 ierr = PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr); 1617 } 1618 } 1619 } 1620 PetscFunctionReturn(0); 1621 } 1622 1623 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 1624 { 1625 Mat_MUMPS *mumps =(Mat_MUMPS*)A->data; 1626 1627 PetscFunctionBegin; 1628 info->block_size = 1.0; 1629 info->nz_allocated = mumps->id.INFOG(20); 1630 info->nz_used = mumps->id.INFOG(20); 1631 info->nz_unneeded = 0.0; 1632 info->assemblies = 0.0; 1633 info->mallocs = 0.0; 1634 info->memory = 0.0; 1635 info->fill_ratio_given = 0; 1636 info->fill_ratio_needed = 0; 1637 info->factor_mallocs = 0; 1638 PetscFunctionReturn(0); 1639 } 1640 1641 /* -------------------------------------------------------------------------------------------*/ 1642 PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is) 1643 { 1644 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1645 const PetscInt *idxs; 1646 PetscInt size,i; 1647 PetscErrorCode ierr; 1648 1649 PetscFunctionBegin; 1650 ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 1651 if (mumps->size > 1) { 1652 PetscBool ls,gs; 1653 1654 ls = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE; 1655 ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);CHKERRQ(ierr); 1656 if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n"); 1657 } 1658 if (mumps->id.size_schur != size) { 1659 ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr); 1660 mumps->id.size_schur = size; 1661 mumps->id.schur_lld = size; 1662 ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr); 1663 } 1664 1665 /* Schur complement matrix */ 1666 ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);CHKERRQ(ierr); 1667 if (mumps->sym == 1) { 1668 ierr = MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); 1669 } 1670 1671 /* MUMPS expects Fortran style indices */ 1672 ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 1673 ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr); 1674 for (i=0;i<size;i++) mumps->id.listvar_schur[i]++; 1675 ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 1676 if (mumps->size > 1) { 1677 mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */ 1678 } else { 1679 if (F->factortype == MAT_FACTOR_LU) { 1680 mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */ 1681 } else { 1682 mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */ 1683 } 1684 } 1685 /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */ 1686 mumps->id.ICNTL(26) = -1; 1687 PetscFunctionReturn(0); 1688 } 1689 1690 /* -------------------------------------------------------------------------------------------*/ 1691 PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S) 1692 { 1693 Mat St; 1694 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1695 PetscScalar *array; 1696 #if defined(PETSC_USE_COMPLEX) 1697 PetscScalar im = PetscSqrtScalar((PetscScalar)-1.0); 1698 #endif 1699 PetscErrorCode ierr; 1700 1701 PetscFunctionBegin; 1702 if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 1703 ierr = MatCreate(PETSC_COMM_SELF,&St);CHKERRQ(ierr); 1704 ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr); 1705 ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr); 1706 ierr = MatSetUp(St);CHKERRQ(ierr); 1707 ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr); 1708 if (!mumps->sym) { /* MUMPS always return a full matrix */ 1709 if (mumps->id.ICNTL(19) == 1) { /* stored by rows */ 1710 PetscInt i,j,N=mumps->id.size_schur; 1711 for (i=0;i<N;i++) { 1712 for (j=0;j<N;j++) { 1713 #if !defined(PETSC_USE_COMPLEX) 1714 PetscScalar val = mumps->id.schur[i*N+j]; 1715 #else 1716 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1717 #endif 1718 array[j*N+i] = val; 1719 } 1720 } 1721 } else { /* stored by columns */ 1722 ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1723 } 1724 } else { /* either full or lower-triangular (not packed) */ 1725 if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */ 1726 PetscInt i,j,N=mumps->id.size_schur; 1727 for (i=0;i<N;i++) { 1728 for (j=i;j<N;j++) { 1729 #if !defined(PETSC_USE_COMPLEX) 1730 PetscScalar val = mumps->id.schur[i*N+j]; 1731 #else 1732 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1733 #endif 1734 array[i*N+j] = val; 1735 array[j*N+i] = val; 1736 } 1737 } 1738 } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */ 1739 ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1740 } else { /* ICNTL(19) == 1 lower triangular stored by rows */ 1741 PetscInt i,j,N=mumps->id.size_schur; 1742 for (i=0;i<N;i++) { 1743 for (j=0;j<i+1;j++) { 1744 #if !defined(PETSC_USE_COMPLEX) 1745 PetscScalar val = mumps->id.schur[i*N+j]; 1746 #else 1747 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1748 #endif 1749 array[i*N+j] = val; 1750 array[j*N+i] = val; 1751 } 1752 } 1753 } 1754 } 1755 ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr); 1756 *S = St; 1757 PetscFunctionReturn(0); 1758 } 1759 1760 /* -------------------------------------------------------------------------------------------*/ 1761 PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 1762 { 1763 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1764 1765 PetscFunctionBegin; 1766 mumps->id.ICNTL(icntl) = ival; 1767 PetscFunctionReturn(0); 1768 } 1769 1770 PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival) 1771 { 1772 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1773 1774 PetscFunctionBegin; 1775 *ival = mumps->id.ICNTL(icntl); 1776 PetscFunctionReturn(0); 1777 } 1778 1779 /*@ 1780 MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 1781 1782 Logically Collective on Mat 1783 1784 Input Parameters: 1785 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1786 . icntl - index of MUMPS parameter array ICNTL() 1787 - ival - value of MUMPS ICNTL(icntl) 1788 1789 Options Database: 1790 . -mat_mumps_icntl_<icntl> <ival> 1791 1792 Level: beginner 1793 1794 References: 1795 . MUMPS Users' Guide 1796 1797 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1798 @*/ 1799 PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 1800 { 1801 PetscErrorCode ierr; 1802 1803 PetscFunctionBegin; 1804 PetscValidType(F,1); 1805 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1806 PetscValidLogicalCollectiveInt(F,icntl,2); 1807 PetscValidLogicalCollectiveInt(F,ival,3); 1808 ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 1809 PetscFunctionReturn(0); 1810 } 1811 1812 /*@ 1813 MatMumpsGetIcntl - Get MUMPS parameter ICNTL() 1814 1815 Logically Collective on Mat 1816 1817 Input Parameters: 1818 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1819 - icntl - index of MUMPS parameter array ICNTL() 1820 1821 Output Parameter: 1822 . ival - value of MUMPS ICNTL(icntl) 1823 1824 Level: beginner 1825 1826 References: 1827 . MUMPS Users' Guide 1828 1829 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1830 @*/ 1831 PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival) 1832 { 1833 PetscErrorCode ierr; 1834 1835 PetscFunctionBegin; 1836 PetscValidType(F,1); 1837 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1838 PetscValidLogicalCollectiveInt(F,icntl,2); 1839 PetscValidIntPointer(ival,3); 1840 ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1841 PetscFunctionReturn(0); 1842 } 1843 1844 /* -------------------------------------------------------------------------------------------*/ 1845 PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val) 1846 { 1847 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1848 1849 PetscFunctionBegin; 1850 mumps->id.CNTL(icntl) = val; 1851 PetscFunctionReturn(0); 1852 } 1853 1854 PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val) 1855 { 1856 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1857 1858 PetscFunctionBegin; 1859 *val = mumps->id.CNTL(icntl); 1860 PetscFunctionReturn(0); 1861 } 1862 1863 /*@ 1864 MatMumpsSetCntl - Set MUMPS parameter CNTL() 1865 1866 Logically Collective on Mat 1867 1868 Input Parameters: 1869 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1870 . icntl - index of MUMPS parameter array CNTL() 1871 - val - value of MUMPS CNTL(icntl) 1872 1873 Options Database: 1874 . -mat_mumps_cntl_<icntl> <val> 1875 1876 Level: beginner 1877 1878 References: 1879 . MUMPS Users' Guide 1880 1881 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1882 @*/ 1883 PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val) 1884 { 1885 PetscErrorCode ierr; 1886 1887 PetscFunctionBegin; 1888 PetscValidType(F,1); 1889 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1890 PetscValidLogicalCollectiveInt(F,icntl,2); 1891 PetscValidLogicalCollectiveReal(F,val,3); 1892 ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr); 1893 PetscFunctionReturn(0); 1894 } 1895 1896 /*@ 1897 MatMumpsGetCntl - Get MUMPS parameter CNTL() 1898 1899 Logically Collective on Mat 1900 1901 Input Parameters: 1902 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1903 - icntl - index of MUMPS parameter array CNTL() 1904 1905 Output Parameter: 1906 . val - value of MUMPS CNTL(icntl) 1907 1908 Level: beginner 1909 1910 References: 1911 . MUMPS Users' Guide 1912 1913 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1914 @*/ 1915 PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val) 1916 { 1917 PetscErrorCode ierr; 1918 1919 PetscFunctionBegin; 1920 PetscValidType(F,1); 1921 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1922 PetscValidLogicalCollectiveInt(F,icntl,2); 1923 PetscValidRealPointer(val,3); 1924 ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1925 PetscFunctionReturn(0); 1926 } 1927 1928 PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info) 1929 { 1930 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1931 1932 PetscFunctionBegin; 1933 *info = mumps->id.INFO(icntl); 1934 PetscFunctionReturn(0); 1935 } 1936 1937 PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog) 1938 { 1939 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1940 1941 PetscFunctionBegin; 1942 *infog = mumps->id.INFOG(icntl); 1943 PetscFunctionReturn(0); 1944 } 1945 1946 PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo) 1947 { 1948 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1949 1950 PetscFunctionBegin; 1951 *rinfo = mumps->id.RINFO(icntl); 1952 PetscFunctionReturn(0); 1953 } 1954 1955 PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog) 1956 { 1957 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1958 1959 PetscFunctionBegin; 1960 *rinfog = mumps->id.RINFOG(icntl); 1961 PetscFunctionReturn(0); 1962 } 1963 1964 /*@ 1965 MatMumpsGetInfo - Get MUMPS parameter INFO() 1966 1967 Logically Collective on Mat 1968 1969 Input Parameters: 1970 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1971 - icntl - index of MUMPS parameter array INFO() 1972 1973 Output Parameter: 1974 . ival - value of MUMPS INFO(icntl) 1975 1976 Level: beginner 1977 1978 References: 1979 . MUMPS Users' Guide 1980 1981 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1982 @*/ 1983 PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival) 1984 { 1985 PetscErrorCode ierr; 1986 1987 PetscFunctionBegin; 1988 PetscValidType(F,1); 1989 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1990 PetscValidIntPointer(ival,3); 1991 ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1992 PetscFunctionReturn(0); 1993 } 1994 1995 /*@ 1996 MatMumpsGetInfog - Get MUMPS parameter INFOG() 1997 1998 Logically Collective on Mat 1999 2000 Input Parameters: 2001 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2002 - icntl - index of MUMPS parameter array INFOG() 2003 2004 Output Parameter: 2005 . ival - value of MUMPS INFOG(icntl) 2006 2007 Level: beginner 2008 2009 References: 2010 . MUMPS Users' Guide 2011 2012 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2013 @*/ 2014 PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival) 2015 { 2016 PetscErrorCode ierr; 2017 2018 PetscFunctionBegin; 2019 PetscValidType(F,1); 2020 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2021 PetscValidIntPointer(ival,3); 2022 ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2023 PetscFunctionReturn(0); 2024 } 2025 2026 /*@ 2027 MatMumpsGetRinfo - Get MUMPS parameter RINFO() 2028 2029 Logically Collective on Mat 2030 2031 Input Parameters: 2032 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2033 - icntl - index of MUMPS parameter array RINFO() 2034 2035 Output Parameter: 2036 . val - value of MUMPS RINFO(icntl) 2037 2038 Level: beginner 2039 2040 References: 2041 . MUMPS Users' Guide 2042 2043 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2044 @*/ 2045 PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val) 2046 { 2047 PetscErrorCode ierr; 2048 2049 PetscFunctionBegin; 2050 PetscValidType(F,1); 2051 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2052 PetscValidRealPointer(val,3); 2053 ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2054 PetscFunctionReturn(0); 2055 } 2056 2057 /*@ 2058 MatMumpsGetRinfog - Get MUMPS parameter RINFOG() 2059 2060 Logically Collective on Mat 2061 2062 Input Parameters: 2063 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2064 - icntl - index of MUMPS parameter array RINFOG() 2065 2066 Output Parameter: 2067 . val - value of MUMPS RINFOG(icntl) 2068 2069 Level: beginner 2070 2071 References: 2072 . MUMPS Users' Guide 2073 2074 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2075 @*/ 2076 PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val) 2077 { 2078 PetscErrorCode ierr; 2079 2080 PetscFunctionBegin; 2081 PetscValidType(F,1); 2082 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2083 PetscValidRealPointer(val,3); 2084 ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2085 PetscFunctionReturn(0); 2086 } 2087 2088 /*MC 2089 MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 2090 distributed and sequential matrices via the external package MUMPS. 2091 2092 Works with MATAIJ and MATSBAIJ matrices 2093 2094 Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS 2095 2096 Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver 2097 2098 Options Database Keys: 2099 + -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages 2100 . -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning 2101 . -mat_mumps_icntl_3 - ICNTL(3): output stream for global information, collected on the host 2102 . -mat_mumps_icntl_4 - ICNTL(4): level of printing (0 to 4) 2103 . -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) 2104 . -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis 2105 . -mat_mumps_icntl_8 - ICNTL(8): scaling strategy (-2 to 8 or 77) 2106 . -mat_mumps_icntl_10 - ICNTL(10): max num of refinements 2107 . -mat_mumps_icntl_11 - ICNTL(11): statistics related to an error analysis (via -ksp_view) 2108 . -mat_mumps_icntl_12 - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) 2109 . -mat_mumps_icntl_13 - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting 2110 . -mat_mumps_icntl_14 - ICNTL(14): percentage increase in the estimated working space 2111 . -mat_mumps_icntl_19 - ICNTL(19): computes the Schur complement 2112 . -mat_mumps_icntl_22 - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) 2113 . -mat_mumps_icntl_23 - ICNTL(23): max size of the working memory (MB) that can allocate per processor 2114 . -mat_mumps_icntl_24 - ICNTL(24): detection of null pivot rows (0 or 1) 2115 . -mat_mumps_icntl_25 - ICNTL(25): compute a solution of a deficient matrix and a null space basis 2116 . -mat_mumps_icntl_26 - ICNTL(26): drives the solution phase if a Schur complement matrix 2117 . -mat_mumps_icntl_28 - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering 2118 . -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis 2119 . -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A) 2120 . -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization 2121 . -mat_mumps_icntl_33 - ICNTL(33): compute determinant 2122 . -mat_mumps_cntl_1 - CNTL(1): relative pivoting threshold 2123 . -mat_mumps_cntl_2 - CNTL(2): stopping criterion of refinement 2124 . -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold 2125 . -mat_mumps_cntl_4 - CNTL(4): value for static pivoting 2126 - -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots 2127 2128 Level: beginner 2129 2130 Notes: When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PCSETUP_FAILED, one can find the MUMPS information about the failure by calling 2131 $ KSPGetPC(ksp,&pc); 2132 $ PCFactorGetMatrix(pc,&mat); 2133 $ MatMumpsGetInfo(mat,....); 2134 $ MatMumpsGetInfog(mat,....); etc. 2135 Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message. 2136 2137 .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix() 2138 2139 M*/ 2140 2141 static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type) 2142 { 2143 PetscFunctionBegin; 2144 *type = MATSOLVERMUMPS; 2145 PetscFunctionReturn(0); 2146 } 2147 2148 /* MatGetFactor for Seq and MPI AIJ matrices */ 2149 static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 2150 { 2151 Mat B; 2152 PetscErrorCode ierr; 2153 Mat_MUMPS *mumps; 2154 PetscBool isSeqAIJ; 2155 2156 PetscFunctionBegin; 2157 /* Create the factorization matrix */ 2158 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 2159 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2160 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2161 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2162 ierr = MatSetUp(B);CHKERRQ(ierr); 2163 2164 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2165 2166 B->ops->view = MatView_MUMPS; 2167 B->ops->getinfo = MatGetInfo_MUMPS; 2168 2169 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2170 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2171 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2172 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2173 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2174 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2175 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2176 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2177 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2178 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2179 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2180 2181 if (ftype == MAT_FACTOR_LU) { 2182 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 2183 B->factortype = MAT_FACTOR_LU; 2184 if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 2185 else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 2186 mumps->sym = 0; 2187 } else { 2188 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2189 B->factortype = MAT_FACTOR_CHOLESKY; 2190 if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 2191 else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 2192 #if defined(PETSC_USE_COMPLEX) 2193 mumps->sym = 2; 2194 #else 2195 if (A->spd_set && A->spd) mumps->sym = 1; 2196 else mumps->sym = 2; 2197 #endif 2198 } 2199 2200 /* set solvertype */ 2201 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2202 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2203 2204 B->ops->destroy = MatDestroy_MUMPS; 2205 B->data = (void*)mumps; 2206 2207 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2208 2209 *F = B; 2210 PetscFunctionReturn(0); 2211 } 2212 2213 /* MatGetFactor for Seq and MPI SBAIJ matrices */ 2214 static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 2215 { 2216 Mat B; 2217 PetscErrorCode ierr; 2218 Mat_MUMPS *mumps; 2219 PetscBool isSeqSBAIJ; 2220 2221 PetscFunctionBegin; 2222 if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 2223 if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead"); 2224 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 2225 /* Create the factorization matrix */ 2226 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2227 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2228 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2229 ierr = MatSetUp(B);CHKERRQ(ierr); 2230 2231 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2232 if (isSeqSBAIJ) { 2233 mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 2234 } else { 2235 mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 2236 } 2237 2238 B->ops->getinfo = MatGetInfo_External; 2239 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2240 B->ops->view = MatView_MUMPS; 2241 2242 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2243 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2244 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2245 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2246 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2247 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2248 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2249 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2250 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2251 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2252 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2253 2254 B->factortype = MAT_FACTOR_CHOLESKY; 2255 #if defined(PETSC_USE_COMPLEX) 2256 mumps->sym = 2; 2257 #else 2258 if (A->spd_set && A->spd) mumps->sym = 1; 2259 else mumps->sym = 2; 2260 #endif 2261 2262 /* set solvertype */ 2263 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2264 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2265 2266 B->ops->destroy = MatDestroy_MUMPS; 2267 B->data = (void*)mumps; 2268 2269 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2270 2271 *F = B; 2272 PetscFunctionReturn(0); 2273 } 2274 2275 static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 2276 { 2277 Mat B; 2278 PetscErrorCode ierr; 2279 Mat_MUMPS *mumps; 2280 PetscBool isSeqBAIJ; 2281 2282 PetscFunctionBegin; 2283 /* Create the factorization matrix */ 2284 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 2285 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2286 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2287 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2288 ierr = MatSetUp(B);CHKERRQ(ierr); 2289 2290 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2291 if (ftype == MAT_FACTOR_LU) { 2292 B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 2293 B->factortype = MAT_FACTOR_LU; 2294 if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 2295 else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 2296 mumps->sym = 0; 2297 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 2298 2299 B->ops->getinfo = MatGetInfo_External; 2300 B->ops->view = MatView_MUMPS; 2301 2302 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2303 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2304 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2305 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2306 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2307 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2308 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2309 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2310 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2311 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2312 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2313 2314 /* set solvertype */ 2315 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2316 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2317 2318 B->ops->destroy = MatDestroy_MUMPS; 2319 B->data = (void*)mumps; 2320 2321 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2322 2323 *F = B; 2324 PetscFunctionReturn(0); 2325 } 2326 2327 /* MatGetFactor for Seq and MPI SELL matrices */ 2328 static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F) 2329 { 2330 Mat B; 2331 PetscErrorCode ierr; 2332 Mat_MUMPS *mumps; 2333 PetscBool isSeqSELL; 2334 2335 PetscFunctionBegin; 2336 /* Create the factorization matrix */ 2337 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL);CHKERRQ(ierr); 2338 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2339 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2340 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2341 ierr = MatSetUp(B);CHKERRQ(ierr); 2342 2343 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2344 2345 B->ops->view = MatView_MUMPS; 2346 B->ops->getinfo = MatGetInfo_MUMPS; 2347 2348 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2349 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2350 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2351 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2352 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2353 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2354 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2355 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2356 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2357 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2358 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2359 2360 if (ftype == MAT_FACTOR_LU) { 2361 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 2362 B->factortype = MAT_FACTOR_LU; 2363 if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij; 2364 else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented"); 2365 mumps->sym = 0; 2366 } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented"); 2367 2368 /* set solvertype */ 2369 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2370 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2371 2372 B->ops->destroy = MatDestroy_MUMPS; 2373 B->data = (void*)mumps; 2374 2375 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2376 2377 *F = B; 2378 PetscFunctionReturn(0); 2379 } 2380 2381 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void) 2382 { 2383 PetscErrorCode ierr; 2384 2385 PetscFunctionBegin; 2386 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2387 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2388 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2389 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2390 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 2391 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2392 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2393 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2394 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2395 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 2396 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps);CHKERRQ(ierr); 2397 PetscFunctionReturn(0); 2398 } 2399 2400