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