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