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