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