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 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL);CHKERRQ(ierr); 732 PetscFunctionReturn(0); 733 } 734 735 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 736 { 737 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 738 PetscScalar *array; 739 Vec b_seq; 740 IS is_iden,is_petsc; 741 PetscErrorCode ierr; 742 PetscInt i; 743 PetscBool second_solve = PETSC_FALSE; 744 static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE; 745 746 PetscFunctionBegin; 747 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); 748 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); 749 750 if (A->factorerrortype) { 751 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); 752 ierr = VecSetInf(x);CHKERRQ(ierr); 753 PetscFunctionReturn(0); 754 } 755 756 mumps->id.nrhs = 1; 757 b_seq = mumps->b_seq; 758 if (mumps->size > 1) { 759 /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 760 ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 761 ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 762 if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 763 } else { /* size == 1 */ 764 ierr = VecCopy(b,x);CHKERRQ(ierr); 765 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 766 } 767 if (!mumps->myid) { /* define rhs on the host */ 768 mumps->id.nrhs = 1; 769 mumps->id.rhs = (MumpsScalar*)array; 770 } 771 772 /* 773 handle condensation step of Schur complement (if any) 774 We set by default ICNTL(26) == -1 when Schur indices have been provided by the user. 775 According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase 776 Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system. 777 This requires an extra call to PetscMUMPS_c and the computation of the factors for S 778 */ 779 if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) { 780 if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n"); 781 second_solve = PETSC_TRUE; 782 ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr); 783 } 784 /* solve phase */ 785 /*-------------*/ 786 mumps->id.job = JOB_SOLVE; 787 PetscMUMPS_c(&mumps->id); 788 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)); 789 790 /* handle expansion step of Schur complement (if any) */ 791 if (second_solve) { 792 ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr); 793 } 794 795 if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 796 if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) { 797 /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */ 798 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 799 } 800 if (!mumps->scat_sol) { /* create scatter scat_sol */ 801 ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 802 for (i=0; i<mumps->id.lsol_loc; i++) { 803 mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 804 } 805 ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 806 ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr); 807 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 808 ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 809 810 mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */ 811 } 812 813 ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 814 ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 815 } 816 PetscFunctionReturn(0); 817 } 818 819 PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 820 { 821 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 822 PetscErrorCode ierr; 823 824 PetscFunctionBegin; 825 mumps->id.ICNTL(9) = 0; 826 ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 827 mumps->id.ICNTL(9) = 1; 828 PetscFunctionReturn(0); 829 } 830 831 PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 832 { 833 PetscErrorCode ierr; 834 Mat Bt = NULL; 835 PetscBool flg, flgT; 836 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 837 PetscInt i,nrhs,M; 838 PetscScalar *array,*bray; 839 840 PetscFunctionBegin; 841 ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 842 ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);CHKERRQ(ierr); 843 if (flgT) { 844 /* sparse B */ 845 ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr); 846 } else { 847 /* dense B */ 848 if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 849 if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution"); 850 } 851 852 ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr); 853 mumps->id.nrhs = nrhs; 854 mumps->id.lrhs = M; 855 mumps->id.rhs = NULL; 856 857 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 858 if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 859 860 if (mumps->size == 1) { 861 PetscScalar *aa; 862 PetscInt spnr,*ia,*ja; 863 PetscBool second_solve = PETSC_FALSE; 864 865 ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 866 mumps->id.rhs = (MumpsScalar*)array; 867 868 if (!Bt) { /* dense B */ 869 /* copy B to X */ 870 ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 871 ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 872 ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 873 } else { /* sparse B */ 874 PetscBool done; 875 876 ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr); 877 ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 878 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure"); 879 /* mumps requires ia and ja start at 1! */ 880 mumps->id.irhs_ptr = ia; 881 mumps->id.irhs_sparse = ja; 882 mumps->id.nz_rhs = ia[spnr] - 1; 883 mumps->id.rhs_sparse = (MumpsScalar*)aa; 884 mumps->id.ICNTL(20) = 1; 885 } 886 /* handle condensation step of Schur complement (if any) */ 887 if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) { 888 second_solve = PETSC_TRUE; 889 ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr); 890 } 891 /* solve phase */ 892 /*-------------*/ 893 mumps->id.job = JOB_SOLVE; 894 PetscMUMPS_c(&mumps->id); 895 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)); 896 897 /* handle expansion step of Schur complement (if any) */ 898 if (second_solve) { 899 ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr); 900 } 901 if (Bt) { /* sparse B */ 902 PetscBool done; 903 904 ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr); 905 ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 906 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure"); 907 mumps->id.ICNTL(20) = 0; 908 } 909 ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 910 } else { /*--------- parallel case --------*/ 911 PetscInt lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save; 912 MumpsScalar *sol_loc,*sol_loc_save; 913 IS is_to,is_from; 914 PetscInt k,proc,j,m; 915 const PetscInt *rstart; 916 Vec v_mpi,b_seq,x_seq; 917 VecScatter scat_rhs,scat_sol; 918 PetscScalar *aa; 919 PetscInt spnr,*ia,*ja; 920 PetscBool done; 921 Mat_MPIAIJ *b; 922 923 if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n"); 924 925 /* create x_seq to hold local solution */ 926 isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */ 927 sol_loc_save = mumps->id.sol_loc; 928 929 lsol_loc = mumps->id.INFO(23); 930 nlsol_loc = nrhs*lsol_loc; /* length of sol_loc */ 931 ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr); 932 mumps->id.sol_loc = (MumpsScalar*)sol_loc; 933 mumps->id.isol_loc = isol_loc; 934 935 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr); 936 937 /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */ 938 /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B; 939 iidx: inverse of idx, will be used by scattering mumps x_seq -> petsc X */ 940 ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr); 941 942 ierr = MatGetOwnershipRanges(X,&rstart);CHKERRQ(ierr); 943 k = 0; 944 for (proc=0; proc<mumps->size; proc++){ 945 for (j=0; j<nrhs; j++){ 946 for (i=rstart[proc]; i<rstart[proc+1]; i++){ 947 iidx[j*M + i] = k; 948 idx[k++] = j*M + i; 949 } 950 } 951 } 952 953 if (!Bt) { /* dense B */ 954 /* copy dense rhs matrix B into vector v_mpi */ 955 ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr); 956 ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 957 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr); 958 ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 959 960 if (!mumps->myid) { 961 ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr); 962 ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 963 ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr); 964 } else { 965 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr); 966 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr); 967 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr); 968 } 969 ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr); 970 ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 971 ierr = ISDestroy(&is_to);CHKERRQ(ierr); 972 ierr = ISDestroy(&is_from);CHKERRQ(ierr); 973 ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 974 975 if (!mumps->myid) { /* define rhs on the host */ 976 ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr); 977 mumps->id.rhs = (MumpsScalar*)bray; 978 ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr); 979 } 980 981 } else { /* sparse B */ 982 b = (Mat_MPIAIJ*)Bt->data; 983 984 /* scatter v_mpi to shard local arrays with X */ 985 ierr = MatGetLocalSize(X,&m,NULL);CHKERRQ(ierr); 986 ierr = MatDenseGetArray(X,&bray);CHKERRQ(ierr); 987 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr); 988 ierr = MatDenseRestoreArray(X,&bray);CHKERRQ(ierr); 989 990 if (!mumps->myid) { 991 ierr = MatSeqAIJGetArray(b->A,&aa);CHKERRQ(ierr); 992 ierr = MatGetRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 993 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure"); 994 /* mumps requires ia and ja start at 1! */ 995 mumps->id.irhs_ptr = ia; 996 mumps->id.irhs_sparse = ja; 997 mumps->id.nz_rhs = ia[spnr] - 1; 998 mumps->id.rhs_sparse = (MumpsScalar*)aa; 999 } else { 1000 mumps->id.irhs_ptr = NULL; 1001 mumps->id.irhs_sparse = NULL; 1002 mumps->id.nz_rhs = 0; 1003 mumps->id.rhs_sparse = NULL; 1004 } 1005 mumps->id.ICNTL(20) = 1; /* rhs is sparse */ 1006 } 1007 1008 /* solve phase */ 1009 /*-------------*/ 1010 mumps->id.job = JOB_SOLVE; 1011 PetscMUMPS_c(&mumps->id); 1012 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)); 1013 1014 /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */ 1015 ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 1016 ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr); 1017 1018 /* create scatter scat_sol */ 1019 ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr); 1020 ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr); 1021 for (i=0; i<lsol_loc; i++) { 1022 isol_loc[i] -= 1; /* change Fortran style to C style */ 1023 idxx[i] = iidx[isol_loc[i]]; 1024 for (j=1; j<nrhs; j++){ 1025 idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M]; 1026 } 1027 } 1028 ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); //is_to is gabage!!! 1029 ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr); 1030 ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1031 ierr = ISDestroy(&is_from);CHKERRQ(ierr); 1032 ierr = ISDestroy(&is_to);CHKERRQ(ierr); 1033 ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1034 ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 1035 1036 /* free spaces */ 1037 mumps->id.sol_loc = sol_loc_save; 1038 mumps->id.isol_loc = isol_loc_save; 1039 1040 ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr); 1041 ierr = PetscFree2(idx,iidx);CHKERRQ(ierr); 1042 ierr = PetscFree(idxx);CHKERRQ(ierr); 1043 ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 1044 ierr = VecDestroy(&v_mpi);CHKERRQ(ierr); 1045 if (Bt) { 1046 if (!mumps->myid) { 1047 ierr = MatSeqAIJRestoreArray(b->A,&aa);CHKERRQ(ierr); 1048 ierr = MatRestoreRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 1049 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure"); 1050 } 1051 } else { 1052 ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 1053 ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr); 1054 } 1055 ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr); 1056 } 1057 PetscFunctionReturn(0); 1058 } 1059 1060 #if !defined(PETSC_USE_COMPLEX) 1061 /* 1062 input: 1063 F: numeric factor 1064 output: 1065 nneg: total number of negative pivots 1066 nzero: total number of zero pivots 1067 npos: (global dimension of F) - nneg - nzero 1068 */ 1069 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 1070 { 1071 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1072 PetscErrorCode ierr; 1073 PetscMPIInt size; 1074 1075 PetscFunctionBegin; 1076 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr); 1077 /* 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 */ 1078 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)); 1079 1080 if (nneg) *nneg = mumps->id.INFOG(12); 1081 if (nzero || npos) { 1082 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"); 1083 if (nzero) *nzero = mumps->id.INFOG(28); 1084 if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28)); 1085 } 1086 PetscFunctionReturn(0); 1087 } 1088 #endif 1089 1090 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 1091 { 1092 Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->data; 1093 PetscErrorCode ierr; 1094 PetscBool isMPIAIJ; 1095 1096 PetscFunctionBegin; 1097 if (mumps->id.INFOG(1) < 0) { 1098 if (mumps->id.INFOG(1) == -6) { 1099 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); 1100 } 1101 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); 1102 PetscFunctionReturn(0); 1103 } 1104 1105 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1106 1107 /* numerical factorization phase */ 1108 /*-------------------------------*/ 1109 mumps->id.job = JOB_FACTNUMERIC; 1110 if (!mumps->id.ICNTL(18)) { /* A is centralized */ 1111 if (!mumps->myid) { 1112 mumps->id.a = (MumpsScalar*)mumps->val; 1113 } 1114 } else { 1115 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1116 } 1117 PetscMUMPS_c(&mumps->id); 1118 if (mumps->id.INFOG(1) < 0) { 1119 if (A->erroriffailure) { 1120 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)); 1121 } else { 1122 if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */ 1123 ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1124 F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1125 } else if (mumps->id.INFOG(1) == -13) { 1126 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); 1127 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 1128 } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) { 1129 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); 1130 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 1131 } else { 1132 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); 1133 F->factorerrortype = MAT_FACTOR_OTHER; 1134 } 1135 } 1136 } 1137 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)); 1138 1139 F->assembled = PETSC_TRUE; 1140 mumps->matstruc = SAME_NONZERO_PATTERN; 1141 if (F->schur) { /* reset Schur status to unfactored */ 1142 if (mumps->id.ICNTL(19) == 1) { /* stored by rows */ 1143 mumps->id.ICNTL(19) = 2; 1144 ierr = MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);CHKERRQ(ierr); 1145 } 1146 ierr = MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr); 1147 } 1148 1149 /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */ 1150 if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3; 1151 1152 if (mumps->size > 1) { 1153 PetscInt lsol_loc; 1154 PetscScalar *sol_loc; 1155 1156 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 1157 1158 /* distributed solution; Create x_seq=sol_loc for repeated use */ 1159 if (mumps->x_seq) { 1160 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 1161 ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 1162 ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 1163 } 1164 lsol_loc = mumps->id.INFO(23); /* length of sol_loc */ 1165 ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr); 1166 mumps->id.lsol_loc = lsol_loc; 1167 mumps->id.sol_loc = (MumpsScalar*)sol_loc; 1168 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr); 1169 } 1170 PetscFunctionReturn(0); 1171 } 1172 1173 /* Sets MUMPS options from the options database */ 1174 PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 1175 { 1176 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1177 PetscErrorCode ierr; 1178 PetscInt icntl,info[40],i,ninfo=40; 1179 PetscBool flg; 1180 1181 PetscFunctionBegin; 1182 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 1183 ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 1184 if (flg) mumps->id.ICNTL(1) = icntl; 1185 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); 1186 if (flg) mumps->id.ICNTL(2) = icntl; 1187 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); 1188 if (flg) mumps->id.ICNTL(3) = icntl; 1189 1190 ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 1191 if (flg) mumps->id.ICNTL(4) = icntl; 1192 if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 1193 1194 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); 1195 if (flg) mumps->id.ICNTL(6) = icntl; 1196 1197 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); 1198 if (flg) { 1199 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"); 1200 else mumps->id.ICNTL(7) = icntl; 1201 } 1202 1203 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); 1204 /* 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() */ 1205 ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr); 1206 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); 1207 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); 1208 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); 1209 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); 1210 ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr); 1211 if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */ 1212 ierr = MatDestroy(&F->schur);CHKERRQ(ierr); 1213 ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr); 1214 } 1215 /* 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 */ 1216 /* 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 */ 1217 1218 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); 1219 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); 1220 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); 1221 if (mumps->id.ICNTL(24)) { 1222 mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 1223 } 1224 1225 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); 1226 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); 1227 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); 1228 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); 1229 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); 1230 /* 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); */ /* call MatMumpsGetInverse() directly */ 1231 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); 1232 /* 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 */ 1233 ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr); 1234 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); 1235 1236 ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr); 1237 ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr); 1238 ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr); 1239 ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr); 1240 ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr); 1241 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); 1242 1243 ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);CHKERRQ(ierr); 1244 1245 ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr); 1246 if (ninfo) { 1247 if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo); 1248 ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr); 1249 mumps->ninfo = ninfo; 1250 for (i=0; i<ninfo; i++) { 1251 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); 1252 else mumps->info[i] = info[i]; 1253 } 1254 } 1255 1256 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1257 PetscFunctionReturn(0); 1258 } 1259 1260 PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps) 1261 { 1262 PetscErrorCode ierr; 1263 1264 PetscFunctionBegin; 1265 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);CHKERRQ(ierr); 1266 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr); 1267 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr); 1268 1269 mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 1270 1271 mumps->id.job = JOB_INIT; 1272 mumps->id.par = 1; /* host participates factorizaton and solve */ 1273 mumps->id.sym = mumps->sym; 1274 PetscMUMPS_c(&mumps->id); 1275 1276 mumps->scat_rhs = NULL; 1277 mumps->scat_sol = NULL; 1278 1279 /* set PETSc-MUMPS default options - override MUMPS default */ 1280 mumps->id.ICNTL(3) = 0; 1281 mumps->id.ICNTL(4) = 0; 1282 if (mumps->size == 1) { 1283 mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 1284 } else { 1285 mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 1286 mumps->id.ICNTL(20) = 0; /* rhs is in dense format */ 1287 mumps->id.ICNTL(21) = 1; /* distributed solution */ 1288 } 1289 1290 /* schur */ 1291 mumps->id.size_schur = 0; 1292 mumps->id.listvar_schur = NULL; 1293 mumps->id.schur = NULL; 1294 mumps->sizeredrhs = 0; 1295 mumps->schur_sol = NULL; 1296 mumps->schur_sizesol = 0; 1297 PetscFunctionReturn(0); 1298 } 1299 1300 PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps) 1301 { 1302 PetscErrorCode ierr; 1303 1304 PetscFunctionBegin; 1305 if (mumps->id.INFOG(1) < 0) { 1306 if (A->erroriffailure) { 1307 SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 1308 } else { 1309 if (mumps->id.INFOG(1) == -6) { 1310 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); 1311 F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT; 1312 } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) { 1313 ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr); 1314 F->factorerrortype = MAT_FACTOR_OUTMEMORY; 1315 } else { 1316 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); 1317 F->factorerrortype = MAT_FACTOR_OTHER; 1318 } 1319 } 1320 } 1321 PetscFunctionReturn(0); 1322 } 1323 1324 /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */ 1325 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1326 { 1327 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1328 PetscErrorCode ierr; 1329 Vec b; 1330 IS is_iden; 1331 const PetscInt M = A->rmap->N; 1332 1333 PetscFunctionBegin; 1334 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1335 1336 /* Set MUMPS options from the options database */ 1337 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1338 1339 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1340 1341 /* analysis phase */ 1342 /*----------------*/ 1343 mumps->id.job = JOB_FACTSYMBOLIC; 1344 mumps->id.n = M; 1345 switch (mumps->id.ICNTL(18)) { 1346 case 0: /* centralized assembled matrix input */ 1347 if (!mumps->myid) { 1348 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1349 if (mumps->id.ICNTL(6)>1) { 1350 mumps->id.a = (MumpsScalar*)mumps->val; 1351 } 1352 if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */ 1353 /* 1354 PetscBool flag; 1355 ierr = ISEqual(r,c,&flag);CHKERRQ(ierr); 1356 if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm"); 1357 ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF); 1358 */ 1359 if (!mumps->myid) { 1360 const PetscInt *idx; 1361 PetscInt i,*perm_in; 1362 1363 ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr); 1364 ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 1365 1366 mumps->id.perm_in = perm_in; 1367 for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 1368 ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 1369 } 1370 } 1371 } 1372 break; 1373 case 3: /* distributed assembled matrix input (size>1) */ 1374 mumps->id.nz_loc = mumps->nz; 1375 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1376 if (mumps->id.ICNTL(6)>1) { 1377 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1378 } 1379 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1380 if (!mumps->myid) { 1381 ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr); 1382 ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr); 1383 } else { 1384 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1385 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1386 } 1387 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1388 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1389 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1390 ierr = VecDestroy(&b);CHKERRQ(ierr); 1391 break; 1392 } 1393 PetscMUMPS_c(&mumps->id); 1394 ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 1395 1396 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1397 F->ops->solve = MatSolve_MUMPS; 1398 F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1399 F->ops->matsolve = MatMatSolve_MUMPS; 1400 PetscFunctionReturn(0); 1401 } 1402 1403 /* Note the Petsc r and c permutations are ignored */ 1404 PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1405 { 1406 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1407 PetscErrorCode ierr; 1408 Vec b; 1409 IS is_iden; 1410 const PetscInt M = A->rmap->N; 1411 1412 PetscFunctionBegin; 1413 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1414 1415 /* Set MUMPS options from the options database */ 1416 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1417 1418 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1419 1420 /* analysis phase */ 1421 /*----------------*/ 1422 mumps->id.job = JOB_FACTSYMBOLIC; 1423 mumps->id.n = M; 1424 switch (mumps->id.ICNTL(18)) { 1425 case 0: /* centralized assembled matrix input */ 1426 if (!mumps->myid) { 1427 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1428 if (mumps->id.ICNTL(6)>1) { 1429 mumps->id.a = (MumpsScalar*)mumps->val; 1430 } 1431 } 1432 break; 1433 case 3: /* distributed assembled matrix input (size>1) */ 1434 mumps->id.nz_loc = mumps->nz; 1435 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1436 if (mumps->id.ICNTL(6)>1) { 1437 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1438 } 1439 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1440 if (!mumps->myid) { 1441 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 1442 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 1443 } else { 1444 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1445 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1446 } 1447 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1448 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1449 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1450 ierr = VecDestroy(&b);CHKERRQ(ierr); 1451 break; 1452 } 1453 PetscMUMPS_c(&mumps->id); 1454 ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 1455 1456 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1457 F->ops->solve = MatSolve_MUMPS; 1458 F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1459 PetscFunctionReturn(0); 1460 } 1461 1462 /* Note the Petsc r permutation and factor info are ignored */ 1463 PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 1464 { 1465 Mat_MUMPS *mumps = (Mat_MUMPS*)F->data; 1466 PetscErrorCode ierr; 1467 Vec b; 1468 IS is_iden; 1469 const PetscInt M = A->rmap->N; 1470 1471 PetscFunctionBegin; 1472 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1473 1474 /* Set MUMPS options from the options database */ 1475 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1476 1477 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1478 1479 /* analysis phase */ 1480 /*----------------*/ 1481 mumps->id.job = JOB_FACTSYMBOLIC; 1482 mumps->id.n = M; 1483 switch (mumps->id.ICNTL(18)) { 1484 case 0: /* centralized assembled matrix input */ 1485 if (!mumps->myid) { 1486 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1487 if (mumps->id.ICNTL(6)>1) { 1488 mumps->id.a = (MumpsScalar*)mumps->val; 1489 } 1490 } 1491 break; 1492 case 3: /* distributed assembled matrix input (size>1) */ 1493 mumps->id.nz_loc = mumps->nz; 1494 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1495 if (mumps->id.ICNTL(6)>1) { 1496 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1497 } 1498 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1499 if (!mumps->myid) { 1500 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 1501 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 1502 } else { 1503 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1504 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1505 } 1506 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1507 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1508 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1509 ierr = VecDestroy(&b);CHKERRQ(ierr); 1510 break; 1511 } 1512 PetscMUMPS_c(&mumps->id); 1513 ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr); 1514 1515 F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1516 F->ops->solve = MatSolve_MUMPS; 1517 F->ops->solvetranspose = MatSolve_MUMPS; 1518 F->ops->matsolve = MatMatSolve_MUMPS; 1519 #if defined(PETSC_USE_COMPLEX) 1520 F->ops->getinertia = NULL; 1521 #else 1522 F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1523 #endif 1524 PetscFunctionReturn(0); 1525 } 1526 1527 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 1528 { 1529 PetscErrorCode ierr; 1530 PetscBool iascii; 1531 PetscViewerFormat format; 1532 Mat_MUMPS *mumps=(Mat_MUMPS*)A->data; 1533 1534 PetscFunctionBegin; 1535 /* check if matrix is mumps type */ 1536 if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 1537 1538 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1539 if (iascii) { 1540 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1541 if (format == PETSC_VIEWER_ASCII_INFO) { 1542 ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1543 ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1544 ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1545 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr); 1546 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr); 1547 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr); 1548 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr); 1549 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr); 1550 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr); 1551 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr); 1552 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scaling strategy): %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr); 1553 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr); 1554 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr); 1555 if (mumps->id.ICNTL(11)>0) { 1556 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr); 1557 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr); 1558 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr); 1559 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr); 1560 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr); 1561 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr); 1562 } 1563 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr); 1564 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr); 1565 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr); 1566 /* ICNTL(15-17) not used */ 1567 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr); 1568 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Schur complement info): %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr); 1569 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr); 1570 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr); 1571 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr); 1572 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr); 1573 1574 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr); 1575 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr); 1576 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr); 1577 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr); 1578 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr); 1579 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr); 1580 1581 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr); 1582 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr); 1583 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr); 1584 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(35) (activate BLR based factorization): %d \n",mumps->id.ICNTL(35));CHKERRQ(ierr); 1585 1586 ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));CHKERRQ(ierr); 1587 ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr); 1588 ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",mumps->id.CNTL(3));CHKERRQ(ierr); 1589 ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",mumps->id.CNTL(4));CHKERRQ(ierr); 1590 ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));CHKERRQ(ierr); 1591 ierr = PetscViewerASCIIPrintf(viewer," CNTL(7) (dropping parameter for BLR): %g \n",mumps->id.CNTL(7));CHKERRQ(ierr); 1592 1593 /* infomation local to each processor */ 1594 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 1595 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 1596 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr); 1597 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1598 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 1599 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr); 1600 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1601 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 1602 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr); 1603 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1604 1605 ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 1606 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr); 1607 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1608 1609 ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 1610 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr); 1611 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1612 1613 ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 1614 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr); 1615 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1616 1617 if (mumps->ninfo && mumps->ninfo <= 40){ 1618 PetscInt i; 1619 for (i=0; i<mumps->ninfo; i++){ 1620 ierr = PetscViewerASCIIPrintf(viewer, " INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr); 1621 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr); 1622 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1623 } 1624 } 1625 1626 1627 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 1628 1629 if (!mumps->myid) { /* information from the host */ 1630 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr); 1631 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr); 1632 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr); 1633 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); 1634 1635 ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr); 1636 ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr); 1637 ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr); 1638 ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr); 1639 ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr); 1640 ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr); 1641 ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr); 1642 ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr); 1643 ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr); 1644 ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr); 1645 ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr); 1646 ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr); 1647 ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr); 1648 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); 1649 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); 1650 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); 1651 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); 1652 ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr); 1653 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); 1654 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); 1655 ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr); 1656 ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr); 1657 ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr); 1658 ierr = PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr); 1659 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); 1660 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); 1661 ierr = PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr); 1662 ierr = PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr); 1663 ierr = PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr); 1664 } 1665 } 1666 } 1667 PetscFunctionReturn(0); 1668 } 1669 1670 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 1671 { 1672 Mat_MUMPS *mumps =(Mat_MUMPS*)A->data; 1673 1674 PetscFunctionBegin; 1675 info->block_size = 1.0; 1676 info->nz_allocated = mumps->id.INFOG(20); 1677 info->nz_used = mumps->id.INFOG(20); 1678 info->nz_unneeded = 0.0; 1679 info->assemblies = 0.0; 1680 info->mallocs = 0.0; 1681 info->memory = 0.0; 1682 info->fill_ratio_given = 0; 1683 info->fill_ratio_needed = 0; 1684 info->factor_mallocs = 0; 1685 PetscFunctionReturn(0); 1686 } 1687 1688 /* -------------------------------------------------------------------------------------------*/ 1689 PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is) 1690 { 1691 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1692 const PetscInt *idxs; 1693 PetscInt size,i; 1694 PetscErrorCode ierr; 1695 1696 PetscFunctionBegin; 1697 ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 1698 if (mumps->size > 1) { 1699 PetscBool ls,gs; 1700 1701 ls = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE; 1702 ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);CHKERRQ(ierr); 1703 if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n"); 1704 } 1705 if (mumps->id.size_schur != size) { 1706 ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr); 1707 mumps->id.size_schur = size; 1708 mumps->id.schur_lld = size; 1709 ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr); 1710 } 1711 1712 /* Schur complement matrix */ 1713 ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);CHKERRQ(ierr); 1714 if (mumps->sym == 1) { 1715 ierr = MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); 1716 } 1717 1718 /* MUMPS expects Fortran style indices */ 1719 ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 1720 ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr); 1721 for (i=0;i<size;i++) mumps->id.listvar_schur[i]++; 1722 ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 1723 if (mumps->size > 1) { 1724 mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */ 1725 } else { 1726 if (F->factortype == MAT_FACTOR_LU) { 1727 mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */ 1728 } else { 1729 mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */ 1730 } 1731 } 1732 /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */ 1733 mumps->id.ICNTL(26) = -1; 1734 PetscFunctionReturn(0); 1735 } 1736 1737 /* -------------------------------------------------------------------------------------------*/ 1738 PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S) 1739 { 1740 Mat St; 1741 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1742 PetscScalar *array; 1743 #if defined(PETSC_USE_COMPLEX) 1744 PetscScalar im = PetscSqrtScalar((PetscScalar)-1.0); 1745 #endif 1746 PetscErrorCode ierr; 1747 1748 PetscFunctionBegin; 1749 if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 1750 ierr = MatCreate(PETSC_COMM_SELF,&St);CHKERRQ(ierr); 1751 ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr); 1752 ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr); 1753 ierr = MatSetUp(St);CHKERRQ(ierr); 1754 ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr); 1755 if (!mumps->sym) { /* MUMPS always return a full matrix */ 1756 if (mumps->id.ICNTL(19) == 1) { /* stored by rows */ 1757 PetscInt i,j,N=mumps->id.size_schur; 1758 for (i=0;i<N;i++) { 1759 for (j=0;j<N;j++) { 1760 #if !defined(PETSC_USE_COMPLEX) 1761 PetscScalar val = mumps->id.schur[i*N+j]; 1762 #else 1763 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1764 #endif 1765 array[j*N+i] = val; 1766 } 1767 } 1768 } else { /* stored by columns */ 1769 ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1770 } 1771 } else { /* either full or lower-triangular (not packed) */ 1772 if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */ 1773 PetscInt i,j,N=mumps->id.size_schur; 1774 for (i=0;i<N;i++) { 1775 for (j=i;j<N;j++) { 1776 #if !defined(PETSC_USE_COMPLEX) 1777 PetscScalar val = mumps->id.schur[i*N+j]; 1778 #else 1779 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1780 #endif 1781 array[i*N+j] = val; 1782 array[j*N+i] = val; 1783 } 1784 } 1785 } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */ 1786 ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1787 } else { /* ICNTL(19) == 1 lower triangular stored by rows */ 1788 PetscInt i,j,N=mumps->id.size_schur; 1789 for (i=0;i<N;i++) { 1790 for (j=0;j<i+1;j++) { 1791 #if !defined(PETSC_USE_COMPLEX) 1792 PetscScalar val = mumps->id.schur[i*N+j]; 1793 #else 1794 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1795 #endif 1796 array[i*N+j] = val; 1797 array[j*N+i] = val; 1798 } 1799 } 1800 } 1801 } 1802 ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr); 1803 *S = St; 1804 PetscFunctionReturn(0); 1805 } 1806 1807 /* -------------------------------------------------------------------------------------------*/ 1808 PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 1809 { 1810 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1811 1812 PetscFunctionBegin; 1813 mumps->id.ICNTL(icntl) = ival; 1814 PetscFunctionReturn(0); 1815 } 1816 1817 PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival) 1818 { 1819 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1820 1821 PetscFunctionBegin; 1822 *ival = mumps->id.ICNTL(icntl); 1823 PetscFunctionReturn(0); 1824 } 1825 1826 /*@ 1827 MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 1828 1829 Logically Collective on Mat 1830 1831 Input Parameters: 1832 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1833 . icntl - index of MUMPS parameter array ICNTL() 1834 - ival - value of MUMPS ICNTL(icntl) 1835 1836 Options Database: 1837 . -mat_mumps_icntl_<icntl> <ival> 1838 1839 Level: beginner 1840 1841 References: 1842 . MUMPS Users' Guide 1843 1844 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1845 @*/ 1846 PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 1847 { 1848 PetscErrorCode ierr; 1849 1850 PetscFunctionBegin; 1851 PetscValidType(F,1); 1852 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1853 PetscValidLogicalCollectiveInt(F,icntl,2); 1854 PetscValidLogicalCollectiveInt(F,ival,3); 1855 ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 1856 PetscFunctionReturn(0); 1857 } 1858 1859 /*@ 1860 MatMumpsGetIcntl - Get MUMPS parameter ICNTL() 1861 1862 Logically Collective on Mat 1863 1864 Input Parameters: 1865 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1866 - icntl - index of MUMPS parameter array ICNTL() 1867 1868 Output Parameter: 1869 . ival - value of MUMPS ICNTL(icntl) 1870 1871 Level: beginner 1872 1873 References: 1874 . MUMPS Users' Guide 1875 1876 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1877 @*/ 1878 PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival) 1879 { 1880 PetscErrorCode ierr; 1881 1882 PetscFunctionBegin; 1883 PetscValidType(F,1); 1884 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1885 PetscValidLogicalCollectiveInt(F,icntl,2); 1886 PetscValidIntPointer(ival,3); 1887 ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 1888 PetscFunctionReturn(0); 1889 } 1890 1891 /* -------------------------------------------------------------------------------------------*/ 1892 PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val) 1893 { 1894 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1895 1896 PetscFunctionBegin; 1897 mumps->id.CNTL(icntl) = val; 1898 PetscFunctionReturn(0); 1899 } 1900 1901 PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val) 1902 { 1903 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1904 1905 PetscFunctionBegin; 1906 *val = mumps->id.CNTL(icntl); 1907 PetscFunctionReturn(0); 1908 } 1909 1910 /*@ 1911 MatMumpsSetCntl - Set MUMPS parameter CNTL() 1912 1913 Logically Collective on Mat 1914 1915 Input Parameters: 1916 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1917 . icntl - index of MUMPS parameter array CNTL() 1918 - val - value of MUMPS CNTL(icntl) 1919 1920 Options Database: 1921 . -mat_mumps_cntl_<icntl> <val> 1922 1923 Level: beginner 1924 1925 References: 1926 . MUMPS Users' Guide 1927 1928 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1929 @*/ 1930 PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val) 1931 { 1932 PetscErrorCode ierr; 1933 1934 PetscFunctionBegin; 1935 PetscValidType(F,1); 1936 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1937 PetscValidLogicalCollectiveInt(F,icntl,2); 1938 PetscValidLogicalCollectiveReal(F,val,3); 1939 ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr); 1940 PetscFunctionReturn(0); 1941 } 1942 1943 /*@ 1944 MatMumpsGetCntl - Get MUMPS parameter CNTL() 1945 1946 Logically Collective on Mat 1947 1948 Input Parameters: 1949 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1950 - icntl - index of MUMPS parameter array CNTL() 1951 1952 Output Parameter: 1953 . val - value of MUMPS CNTL(icntl) 1954 1955 Level: beginner 1956 1957 References: 1958 . MUMPS Users' Guide 1959 1960 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 1961 @*/ 1962 PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val) 1963 { 1964 PetscErrorCode ierr; 1965 1966 PetscFunctionBegin; 1967 PetscValidType(F,1); 1968 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 1969 PetscValidLogicalCollectiveInt(F,icntl,2); 1970 PetscValidRealPointer(val,3); 1971 ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 1972 PetscFunctionReturn(0); 1973 } 1974 1975 PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info) 1976 { 1977 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1978 1979 PetscFunctionBegin; 1980 *info = mumps->id.INFO(icntl); 1981 PetscFunctionReturn(0); 1982 } 1983 1984 PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog) 1985 { 1986 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1987 1988 PetscFunctionBegin; 1989 *infog = mumps->id.INFOG(icntl); 1990 PetscFunctionReturn(0); 1991 } 1992 1993 PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo) 1994 { 1995 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 1996 1997 PetscFunctionBegin; 1998 *rinfo = mumps->id.RINFO(icntl); 1999 PetscFunctionReturn(0); 2000 } 2001 2002 PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog) 2003 { 2004 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 2005 2006 PetscFunctionBegin; 2007 *rinfog = mumps->id.RINFOG(icntl); 2008 PetscFunctionReturn(0); 2009 } 2010 2011 PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS) 2012 { 2013 PetscErrorCode ierr; 2014 Mat Bt = NULL; 2015 PetscBool flgT; 2016 Mat_MUMPS *mumps =(Mat_MUMPS*)F->data; 2017 PetscBool done; 2018 PetscScalar *aa; 2019 PetscInt spnr,*ia,*ja; 2020 2021 PetscFunctionBegin; 2022 if (!mumps->myid) { 2023 PetscValidIntPointer(spRHS,2); 2024 ierr = PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flgT);CHKERRQ(ierr); 2025 if (flgT) { 2026 ierr = MatTransposeGetMat(spRHS,&Bt);CHKERRQ(ierr); 2027 } else { 2028 SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix"); 2029 } 2030 } 2031 2032 ierr = MatMumpsSetIcntl(F,30,1);CHKERRQ(ierr); 2033 2034 if (!mumps->myid) { 2035 ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr); 2036 ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 2037 if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure"); 2038 2039 mumps->id.irhs_ptr = ia; 2040 mumps->id.irhs_sparse = ja; 2041 mumps->id.nz_rhs = ia[spnr] - 1; 2042 mumps->id.rhs_sparse = (MumpsScalar*)aa; 2043 } else { 2044 mumps->id.irhs_ptr = NULL; 2045 mumps->id.irhs_sparse = NULL; 2046 mumps->id.nz_rhs = 0; 2047 mumps->id.rhs_sparse = NULL; 2048 } 2049 mumps->id.ICNTL(20) = 1; /* rhs is sparse */ 2050 mumps->id.ICNTL(21) = 0; /* solution is in assembled centralized format */ 2051 2052 /* solve phase */ 2053 /*-------------*/ 2054 mumps->id.job = JOB_SOLVE; 2055 PetscMUMPS_c(&mumps->id); 2056 if (mumps->id.INFOG(1) < 0) 2057 SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2)); 2058 2059 if (!mumps->myid) { 2060 ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr); 2061 ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr); 2062 } 2063 PetscFunctionReturn(0); 2064 } 2065 2066 /*@ 2067 MatMumpsGetInverse - Get user-specified set of entries in inverse of A 2068 2069 Logically Collective on Mat 2070 2071 Input Parameters: 2072 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2073 - spRHS - sequential sparse matrix in MATTRANSPOSEMAT format holding specified indices in processor[0] 2074 2075 Output Parameter: 2076 . spRHS - requested entries of inverse of A 2077 2078 Level: beginner 2079 2080 References: 2081 . MUMPS Users' Guide 2082 2083 .seealso: MatGetFactor(), MatCreateTranspose() 2084 @*/ 2085 PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS) 2086 { 2087 PetscErrorCode ierr; 2088 2089 PetscFunctionBegin; 2090 PetscValidType(F,1); 2091 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2092 ierr = PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));CHKERRQ(ierr); 2093 PetscFunctionReturn(0); 2094 } 2095 2096 /*@ 2097 MatMumpsGetInfo - Get MUMPS parameter INFO() 2098 2099 Logically Collective on Mat 2100 2101 Input Parameters: 2102 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2103 - icntl - index of MUMPS parameter array INFO() 2104 2105 Output Parameter: 2106 . ival - value of MUMPS INFO(icntl) 2107 2108 Level: beginner 2109 2110 References: 2111 . MUMPS Users' Guide 2112 2113 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2114 @*/ 2115 PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival) 2116 { 2117 PetscErrorCode ierr; 2118 2119 PetscFunctionBegin; 2120 PetscValidType(F,1); 2121 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2122 PetscValidIntPointer(ival,3); 2123 ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2124 PetscFunctionReturn(0); 2125 } 2126 2127 /*@ 2128 MatMumpsGetInfog - Get MUMPS parameter INFOG() 2129 2130 Logically Collective on Mat 2131 2132 Input Parameters: 2133 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2134 - icntl - index of MUMPS parameter array INFOG() 2135 2136 Output Parameter: 2137 . ival - value of MUMPS INFOG(icntl) 2138 2139 Level: beginner 2140 2141 References: 2142 . MUMPS Users' Guide 2143 2144 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2145 @*/ 2146 PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival) 2147 { 2148 PetscErrorCode ierr; 2149 2150 PetscFunctionBegin; 2151 PetscValidType(F,1); 2152 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2153 PetscValidIntPointer(ival,3); 2154 ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2155 PetscFunctionReturn(0); 2156 } 2157 2158 /*@ 2159 MatMumpsGetRinfo - Get MUMPS parameter RINFO() 2160 2161 Logically Collective on Mat 2162 2163 Input Parameters: 2164 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2165 - icntl - index of MUMPS parameter array RINFO() 2166 2167 Output Parameter: 2168 . val - value of MUMPS RINFO(icntl) 2169 2170 Level: beginner 2171 2172 References: 2173 . MUMPS Users' Guide 2174 2175 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2176 @*/ 2177 PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val) 2178 { 2179 PetscErrorCode ierr; 2180 2181 PetscFunctionBegin; 2182 PetscValidType(F,1); 2183 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2184 PetscValidRealPointer(val,3); 2185 ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2186 PetscFunctionReturn(0); 2187 } 2188 2189 /*@ 2190 MatMumpsGetRinfog - Get MUMPS parameter RINFOG() 2191 2192 Logically Collective on Mat 2193 2194 Input Parameters: 2195 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2196 - icntl - index of MUMPS parameter array RINFOG() 2197 2198 Output Parameter: 2199 . val - value of MUMPS RINFOG(icntl) 2200 2201 Level: beginner 2202 2203 References: 2204 . MUMPS Users' Guide 2205 2206 .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog() 2207 @*/ 2208 PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val) 2209 { 2210 PetscErrorCode ierr; 2211 2212 PetscFunctionBegin; 2213 PetscValidType(F,1); 2214 if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix"); 2215 PetscValidRealPointer(val,3); 2216 ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2217 PetscFunctionReturn(0); 2218 } 2219 2220 /*MC 2221 MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 2222 distributed and sequential matrices via the external package MUMPS. 2223 2224 Works with MATAIJ and MATSBAIJ matrices 2225 2226 Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS 2227 2228 Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver 2229 2230 Options Database Keys: 2231 + -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages 2232 . -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning 2233 . -mat_mumps_icntl_3 - ICNTL(3): output stream for global information, collected on the host 2234 . -mat_mumps_icntl_4 - ICNTL(4): level of printing (0 to 4) 2235 . -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) 2236 . -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis 2237 . -mat_mumps_icntl_8 - ICNTL(8): scaling strategy (-2 to 8 or 77) 2238 . -mat_mumps_icntl_10 - ICNTL(10): max num of refinements 2239 . -mat_mumps_icntl_11 - ICNTL(11): statistics related to an error analysis (via -ksp_view) 2240 . -mat_mumps_icntl_12 - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) 2241 . -mat_mumps_icntl_13 - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting 2242 . -mat_mumps_icntl_14 - ICNTL(14): percentage increase in the estimated working space 2243 . -mat_mumps_icntl_19 - ICNTL(19): computes the Schur complement 2244 . -mat_mumps_icntl_22 - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) 2245 . -mat_mumps_icntl_23 - ICNTL(23): max size of the working memory (MB) that can allocate per processor 2246 . -mat_mumps_icntl_24 - ICNTL(24): detection of null pivot rows (0 or 1) 2247 . -mat_mumps_icntl_25 - ICNTL(25): compute a solution of a deficient matrix and a null space basis 2248 . -mat_mumps_icntl_26 - ICNTL(26): drives the solution phase if a Schur complement matrix 2249 . -mat_mumps_icntl_28 - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering 2250 . -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis 2251 . -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A) 2252 . -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization 2253 . -mat_mumps_icntl_33 - ICNTL(33): compute determinant 2254 . -mat_mumps_cntl_1 - CNTL(1): relative pivoting threshold 2255 . -mat_mumps_cntl_2 - CNTL(2): stopping criterion of refinement 2256 . -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold 2257 . -mat_mumps_cntl_4 - CNTL(4): value for static pivoting 2258 - -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots 2259 2260 Level: beginner 2261 2262 Notes: 2263 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 2264 $ KSPGetPC(ksp,&pc); 2265 $ PCFactorGetMatrix(pc,&mat); 2266 $ MatMumpsGetInfo(mat,....); 2267 $ MatMumpsGetInfog(mat,....); etc. 2268 Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message. 2269 2270 .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix() 2271 2272 M*/ 2273 2274 static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type) 2275 { 2276 PetscFunctionBegin; 2277 *type = MATSOLVERMUMPS; 2278 PetscFunctionReturn(0); 2279 } 2280 2281 /* MatGetFactor for Seq and MPI AIJ matrices */ 2282 static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 2283 { 2284 Mat B; 2285 PetscErrorCode ierr; 2286 Mat_MUMPS *mumps; 2287 PetscBool isSeqAIJ; 2288 2289 PetscFunctionBegin; 2290 /* Create the factorization matrix */ 2291 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 2292 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2293 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2294 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2295 ierr = MatSetUp(B);CHKERRQ(ierr); 2296 2297 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2298 2299 B->ops->view = MatView_MUMPS; 2300 B->ops->getinfo = MatGetInfo_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 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);CHKERRQ(ierr); 2314 2315 if (ftype == MAT_FACTOR_LU) { 2316 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 2317 B->factortype = MAT_FACTOR_LU; 2318 if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 2319 else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 2320 mumps->sym = 0; 2321 } else { 2322 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2323 B->factortype = MAT_FACTOR_CHOLESKY; 2324 if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 2325 else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 2326 #if defined(PETSC_USE_COMPLEX) 2327 mumps->sym = 2; 2328 #else 2329 if (A->spd_set && A->spd) mumps->sym = 1; 2330 else mumps->sym = 2; 2331 #endif 2332 } 2333 2334 /* set solvertype */ 2335 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2336 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2337 2338 B->ops->destroy = MatDestroy_MUMPS; 2339 B->data = (void*)mumps; 2340 2341 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2342 2343 *F = B; 2344 PetscFunctionReturn(0); 2345 } 2346 2347 /* MatGetFactor for Seq and MPI SBAIJ matrices */ 2348 static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 2349 { 2350 Mat B; 2351 PetscErrorCode ierr; 2352 Mat_MUMPS *mumps; 2353 PetscBool isSeqSBAIJ; 2354 2355 PetscFunctionBegin; 2356 if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 2357 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"); 2358 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 2359 /* Create the factorization matrix */ 2360 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2361 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2362 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2363 ierr = MatSetUp(B);CHKERRQ(ierr); 2364 2365 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2366 if (isSeqSBAIJ) { 2367 mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 2368 } else { 2369 mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 2370 } 2371 2372 B->ops->getinfo = MatGetInfo_External; 2373 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2374 B->ops->view = MatView_MUMPS; 2375 2376 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2377 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2378 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2379 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2380 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2381 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2382 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2383 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2384 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2385 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2386 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2387 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);CHKERRQ(ierr); 2388 2389 B->factortype = MAT_FACTOR_CHOLESKY; 2390 #if defined(PETSC_USE_COMPLEX) 2391 mumps->sym = 2; 2392 #else 2393 if (A->spd_set && A->spd) mumps->sym = 1; 2394 else mumps->sym = 2; 2395 #endif 2396 2397 /* set solvertype */ 2398 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2399 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2400 2401 B->ops->destroy = MatDestroy_MUMPS; 2402 B->data = (void*)mumps; 2403 2404 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2405 2406 *F = B; 2407 PetscFunctionReturn(0); 2408 } 2409 2410 static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 2411 { 2412 Mat B; 2413 PetscErrorCode ierr; 2414 Mat_MUMPS *mumps; 2415 PetscBool isSeqBAIJ; 2416 2417 PetscFunctionBegin; 2418 /* Create the factorization matrix */ 2419 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 2420 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2421 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2422 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2423 ierr = MatSetUp(B);CHKERRQ(ierr); 2424 2425 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2426 if (ftype == MAT_FACTOR_LU) { 2427 B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 2428 B->factortype = MAT_FACTOR_LU; 2429 if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 2430 else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 2431 mumps->sym = 0; 2432 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 2433 2434 B->ops->getinfo = MatGetInfo_External; 2435 B->ops->view = MatView_MUMPS; 2436 2437 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2438 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2439 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2440 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2441 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2442 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2443 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2444 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2445 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2446 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2447 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2448 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);CHKERRQ(ierr); 2449 2450 /* set solvertype */ 2451 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2452 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2453 2454 B->ops->destroy = MatDestroy_MUMPS; 2455 B->data = (void*)mumps; 2456 2457 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2458 2459 *F = B; 2460 PetscFunctionReturn(0); 2461 } 2462 2463 /* MatGetFactor for Seq and MPI SELL matrices */ 2464 static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F) 2465 { 2466 Mat B; 2467 PetscErrorCode ierr; 2468 Mat_MUMPS *mumps; 2469 PetscBool isSeqSELL; 2470 2471 PetscFunctionBegin; 2472 /* Create the factorization matrix */ 2473 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL);CHKERRQ(ierr); 2474 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2475 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2476 ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr); 2477 ierr = MatSetUp(B);CHKERRQ(ierr); 2478 2479 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2480 2481 B->ops->view = MatView_MUMPS; 2482 B->ops->getinfo = MatGetInfo_MUMPS; 2483 2484 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr); 2485 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2486 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2487 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2488 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2489 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2490 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2491 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2492 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2493 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2494 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2495 2496 if (ftype == MAT_FACTOR_LU) { 2497 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 2498 B->factortype = MAT_FACTOR_LU; 2499 if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij; 2500 else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented"); 2501 mumps->sym = 0; 2502 } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented"); 2503 2504 /* set solvertype */ 2505 ierr = PetscFree(B->solvertype);CHKERRQ(ierr); 2506 ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr); 2507 2508 B->ops->destroy = MatDestroy_MUMPS; 2509 B->data = (void*)mumps; 2510 2511 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2512 2513 *F = B; 2514 PetscFunctionReturn(0); 2515 } 2516 2517 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void) 2518 { 2519 PetscErrorCode ierr; 2520 2521 PetscFunctionBegin; 2522 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2523 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2524 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2525 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2526 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 2527 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2528 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2529 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2530 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2531 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 2532 ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps);CHKERRQ(ierr); 2533 PetscFunctionReturn(0); 2534 } 2535 2536