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