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 mumps->id.nrhs = 1; 878 b_seq = mumps->b_seq; 879 if (mumps->size > 1) { 880 /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 881 ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 882 ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 883 if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);} 884 } else { /* size == 1 */ 885 ierr = VecCopy(b,x);CHKERRQ(ierr); 886 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 887 } 888 if (!mumps->myid) { /* define rhs on the host */ 889 mumps->id.nrhs = 1; 890 mumps->id.rhs = (MumpsScalar*)array; 891 } 892 893 /* 894 handle condensation step of Schur complement (if any) 895 We set by default ICNTL(26) == -1 when Schur indices have been provided by the user. 896 According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase 897 Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system. 898 This requires an extra call to PetscMUMPS_c and the computation of the factors for S 899 */ 900 if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) { 901 second_solve = PETSC_TRUE; 902 ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 903 } 904 /* solve phase */ 905 /*-------------*/ 906 mumps->id.job = JOB_SOLVE; 907 PetscMUMPS_c(&mumps->id); 908 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)); 909 910 /* handle expansion step of Schur complement (if any) */ 911 if (second_solve) { 912 ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr); 913 } 914 915 if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 916 if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) { 917 /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */ 918 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 919 } 920 if (!mumps->scat_sol) { /* create scatter scat_sol */ 921 ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 922 for (i=0; i<mumps->id.lsol_loc; i++) { 923 mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 924 } 925 ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr); /* to */ 926 ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr); 927 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 928 ierr = ISDestroy(&is_petsc);CHKERRQ(ierr); 929 930 mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */ 931 } 932 933 ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 934 ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 935 } 936 PetscFunctionReturn(0); 937 } 938 939 #undef __FUNCT__ 940 #define __FUNCT__ "MatSolveTranspose_MUMPS" 941 PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x) 942 { 943 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 944 PetscErrorCode ierr; 945 946 PetscFunctionBegin; 947 mumps->id.ICNTL(9) = 0; 948 ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr); 949 mumps->id.ICNTL(9) = 1; 950 PetscFunctionReturn(0); 951 } 952 953 #undef __FUNCT__ 954 #define __FUNCT__ "MatMatSolve_MUMPS" 955 PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X) 956 { 957 PetscErrorCode ierr; 958 PetscBool flg; 959 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 960 PetscInt i,nrhs,M; 961 PetscScalar *array,*bray; 962 963 PetscFunctionBegin; 964 ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 965 if (!flg) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 966 ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 967 if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 968 if (B->rmap->n != X->rmap->n) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution"); 969 970 ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr); 971 mumps->id.nrhs = nrhs; 972 mumps->id.lrhs = M; 973 974 if (mumps->size == 1) { 975 /* copy B to X */ 976 ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 977 ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 978 ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 979 ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 980 mumps->id.rhs = (MumpsScalar*)array; 981 /* handle condensation step of Schur complement (if any) */ 982 ierr = MatMumpsHandleSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 983 984 /* solve phase */ 985 /*-------------*/ 986 mumps->id.job = JOB_SOLVE; 987 PetscMUMPS_c(&mumps->id); 988 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)); 989 990 /* handle expansion step of Schur complement (if any) */ 991 ierr = MatMumpsHandleSchur_Private(mumps,PETSC_TRUE);CHKERRQ(ierr); 992 ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 993 } else { /*--------- parallel case --------*/ 994 PetscInt lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save; 995 MumpsScalar *sol_loc,*sol_loc_save; 996 IS is_to,is_from; 997 PetscInt k,proc,j,m; 998 const PetscInt *rstart; 999 Vec v_mpi,b_seq,x_seq; 1000 VecScatter scat_rhs,scat_sol; 1001 1002 /* create x_seq to hold local solution */ 1003 isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */ 1004 sol_loc_save = mumps->id.sol_loc; 1005 1006 lsol_loc = mumps->id.INFO(23); 1007 nlsol_loc = nrhs*lsol_loc; /* length of sol_loc */ 1008 ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr); 1009 mumps->id.sol_loc = (MumpsScalar*)sol_loc; 1010 mumps->id.isol_loc = isol_loc; 1011 1012 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr); 1013 1014 /* copy rhs matrix B into vector v_mpi */ 1015 ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr); 1016 ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr); 1017 ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr); 1018 ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr); 1019 1020 /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */ 1021 /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B; 1022 iidx: inverse of idx, will be used by scattering xx_seq -> X */ 1023 ierr = PetscMalloc2(nrhs*M,&idx,nrhs*M,&iidx);CHKERRQ(ierr); 1024 ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr); 1025 k = 0; 1026 for (proc=0; proc<mumps->size; proc++){ 1027 for (j=0; j<nrhs; j++){ 1028 for (i=rstart[proc]; i<rstart[proc+1]; i++){ 1029 iidx[j*M + i] = k; 1030 idx[k++] = j*M + i; 1031 } 1032 } 1033 } 1034 1035 if (!mumps->myid) { 1036 ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr); 1037 ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 1038 ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr); 1039 } else { 1040 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr); 1041 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr); 1042 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr); 1043 } 1044 ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr); 1045 ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1046 ierr = ISDestroy(&is_to);CHKERRQ(ierr); 1047 ierr = ISDestroy(&is_from);CHKERRQ(ierr); 1048 ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1049 1050 if (!mumps->myid) { /* define rhs on the host */ 1051 ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr); 1052 mumps->id.rhs = (MumpsScalar*)bray; 1053 ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr); 1054 } 1055 1056 /* solve phase */ 1057 /*-------------*/ 1058 mumps->id.job = JOB_SOLVE; 1059 PetscMUMPS_c(&mumps->id); 1060 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)); 1061 1062 /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */ 1063 ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr); 1064 ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr); 1065 1066 /* create scatter scat_sol */ 1067 ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr); 1068 ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr); 1069 for (i=0; i<lsol_loc; i++) { 1070 isol_loc[i] -= 1; /* change Fortran style to C style */ 1071 idxx[i] = iidx[isol_loc[i]]; 1072 for (j=1; j<nrhs; j++){ 1073 idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M]; 1074 } 1075 } 1076 ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr); 1077 ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr); 1078 ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1079 ierr = ISDestroy(&is_from);CHKERRQ(ierr); 1080 ierr = ISDestroy(&is_to);CHKERRQ(ierr); 1081 ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1082 ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr); 1083 1084 /* free spaces */ 1085 mumps->id.sol_loc = sol_loc_save; 1086 mumps->id.isol_loc = isol_loc_save; 1087 1088 ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr); 1089 ierr = PetscFree2(idx,iidx);CHKERRQ(ierr); 1090 ierr = PetscFree(idxx);CHKERRQ(ierr); 1091 ierr = VecDestroy(&x_seq);CHKERRQ(ierr); 1092 ierr = VecDestroy(&v_mpi);CHKERRQ(ierr); 1093 ierr = VecDestroy(&b_seq);CHKERRQ(ierr); 1094 ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr); 1095 ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr); 1096 } 1097 PetscFunctionReturn(0); 1098 } 1099 1100 #if !defined(PETSC_USE_COMPLEX) 1101 /* 1102 input: 1103 F: numeric factor 1104 output: 1105 nneg: total number of negative pivots 1106 nzero: 0 1107 npos: (global dimension of F) - nneg 1108 */ 1109 1110 #undef __FUNCT__ 1111 #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 1112 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 1113 { 1114 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1115 PetscErrorCode ierr; 1116 PetscMPIInt size; 1117 1118 PetscFunctionBegin; 1119 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr); 1120 /* 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 */ 1121 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)); 1122 1123 if (nneg) *nneg = mumps->id.INFOG(12); 1124 if (nzero || npos) { 1125 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"); 1126 if (nzero) *nzero = mumps->id.INFOG(28); 1127 if (npos) *npos = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28)); 1128 } 1129 PetscFunctionReturn(0); 1130 } 1131 #endif /* !defined(PETSC_USE_COMPLEX) */ 1132 1133 #undef __FUNCT__ 1134 #define __FUNCT__ "MatFactorNumeric_MUMPS" 1135 PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 1136 { 1137 Mat_MUMPS *mumps =(Mat_MUMPS*)(F)->spptr; 1138 PetscErrorCode ierr; 1139 Mat F_diag; 1140 PetscBool isMPIAIJ; 1141 1142 PetscFunctionBegin; 1143 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1144 1145 /* numerical factorization phase */ 1146 /*-------------------------------*/ 1147 mumps->id.job = JOB_FACTNUMERIC; 1148 if (!mumps->id.ICNTL(18)) { /* A is centralized */ 1149 if (!mumps->myid) { 1150 mumps->id.a = (MumpsScalar*)mumps->val; 1151 } 1152 } else { 1153 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1154 } 1155 PetscMUMPS_c(&mumps->id); 1156 if (mumps->id.INFOG(1) < 0) { 1157 if (mumps->id.INFO(1) == -13) { 1158 if (mumps->id.INFO(2) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",-mumps->id.INFO(2)); 1159 else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d bytes\n",mumps->id.INFO(2)); 1160 } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",mumps->id.INFO(1),mumps->id.INFO(2)); 1161 } 1162 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)); 1163 1164 (F)->assembled = PETSC_TRUE; 1165 mumps->matstruc = SAME_NONZERO_PATTERN; 1166 mumps->schur_factored = PETSC_FALSE; 1167 mumps->schur_inverted = PETSC_FALSE; 1168 1169 /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */ 1170 if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3; 1171 1172 if (mumps->size > 1) { 1173 PetscInt lsol_loc; 1174 PetscScalar *sol_loc; 1175 1176 ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 1177 if (isMPIAIJ) F_diag = ((Mat_MPIAIJ*)(F)->data)->A; 1178 else F_diag = ((Mat_MPISBAIJ*)(F)->data)->A; 1179 F_diag->assembled = PETSC_TRUE; 1180 1181 /* distributed solution; Create x_seq=sol_loc for repeated use */ 1182 if (mumps->x_seq) { 1183 ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr); 1184 ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr); 1185 ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr); 1186 } 1187 lsol_loc = mumps->id.INFO(23); /* length of sol_loc */ 1188 ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr); 1189 mumps->id.lsol_loc = lsol_loc; 1190 mumps->id.sol_loc = (MumpsScalar*)sol_loc; 1191 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr); 1192 } 1193 PetscFunctionReturn(0); 1194 } 1195 1196 /* Sets MUMPS options from the options database */ 1197 #undef __FUNCT__ 1198 #define __FUNCT__ "PetscSetMUMPSFromOptions" 1199 PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A) 1200 { 1201 Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1202 PetscErrorCode ierr; 1203 PetscInt icntl,info[40],i,ninfo=40; 1204 PetscBool flg; 1205 1206 PetscFunctionBegin; 1207 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 1208 ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr); 1209 if (flg) mumps->id.ICNTL(1) = icntl; 1210 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); 1211 if (flg) mumps->id.ICNTL(2) = icntl; 1212 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); 1213 if (flg) mumps->id.ICNTL(3) = icntl; 1214 1215 ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 1216 if (flg) mumps->id.ICNTL(4) = icntl; 1217 if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */ 1218 1219 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); 1220 if (flg) mumps->id.ICNTL(6) = icntl; 1221 1222 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); 1223 if (flg) { 1224 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"); 1225 else mumps->id.ICNTL(7) = icntl; 1226 } 1227 1228 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); 1229 /* 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() */ 1230 ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr); 1231 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); 1232 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); 1233 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); 1234 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); 1235 ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr); 1236 if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */ 1237 ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr); 1238 } 1239 /* 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 */ 1240 /* 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 */ 1241 1242 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); 1243 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); 1244 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); 1245 if (mumps->id.ICNTL(24)) { 1246 mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */ 1247 } 1248 1249 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); 1250 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); 1251 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); 1252 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); 1253 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); 1254 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); 1255 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); 1256 /* 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 */ 1257 ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr); 1258 1259 ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr); 1260 ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr); 1261 ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr); 1262 ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr); 1263 ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr); 1264 1265 ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL); 1266 1267 ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr); 1268 if (ninfo) { 1269 if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo); 1270 ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr); 1271 mumps->ninfo = ninfo; 1272 for (i=0; i<ninfo; i++) { 1273 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); 1274 else { 1275 mumps->info[i] = info[i]; 1276 } 1277 } 1278 } 1279 1280 PetscOptionsEnd(); 1281 PetscFunctionReturn(0); 1282 } 1283 1284 #undef __FUNCT__ 1285 #define __FUNCT__ "PetscInitializeMUMPS" 1286 PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps) 1287 { 1288 PetscErrorCode ierr; 1289 1290 PetscFunctionBegin; 1291 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid); 1292 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr); 1293 ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr); 1294 1295 mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps); 1296 1297 mumps->id.job = JOB_INIT; 1298 mumps->id.par = 1; /* host participates factorizaton and solve */ 1299 mumps->id.sym = mumps->sym; 1300 PetscMUMPS_c(&mumps->id); 1301 1302 mumps->scat_rhs = NULL; 1303 mumps->scat_sol = NULL; 1304 1305 /* set PETSc-MUMPS default options - override MUMPS default */ 1306 mumps->id.ICNTL(3) = 0; 1307 mumps->id.ICNTL(4) = 0; 1308 if (mumps->size == 1) { 1309 mumps->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 1310 } else { 1311 mumps->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 1312 mumps->id.ICNTL(20) = 0; /* rhs is in dense format */ 1313 mumps->id.ICNTL(21) = 1; /* distributed solution */ 1314 } 1315 1316 /* schur */ 1317 mumps->id.size_schur = 0; 1318 mumps->id.listvar_schur = NULL; 1319 mumps->id.schur = NULL; 1320 mumps->sizeredrhs = 0; 1321 mumps->schur_pivots = NULL; 1322 mumps->schur_work = NULL; 1323 mumps->schur_sol = NULL; 1324 mumps->schur_sizesol = 0; 1325 mumps->schur_factored = PETSC_FALSE; 1326 mumps->schur_inverted = PETSC_FALSE; 1327 PetscFunctionReturn(0); 1328 } 1329 1330 /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */ 1331 #undef __FUNCT__ 1332 #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 1333 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1334 { 1335 Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1336 PetscErrorCode ierr; 1337 Vec b; 1338 IS is_iden; 1339 const PetscInt M = A->rmap->N; 1340 1341 PetscFunctionBegin; 1342 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1343 1344 /* Set MUMPS options from the options database */ 1345 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1346 1347 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1348 1349 /* analysis phase */ 1350 /*----------------*/ 1351 mumps->id.job = JOB_FACTSYMBOLIC; 1352 mumps->id.n = M; 1353 switch (mumps->id.ICNTL(18)) { 1354 case 0: /* centralized assembled matrix input */ 1355 if (!mumps->myid) { 1356 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1357 if (mumps->id.ICNTL(6)>1) { 1358 mumps->id.a = (MumpsScalar*)mumps->val; 1359 } 1360 if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */ 1361 /* 1362 PetscBool flag; 1363 ierr = ISEqual(r,c,&flag);CHKERRQ(ierr); 1364 if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm"); 1365 ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF); 1366 */ 1367 if (!mumps->myid) { 1368 const PetscInt *idx; 1369 PetscInt i,*perm_in; 1370 1371 ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr); 1372 ierr = ISGetIndices(r,&idx);CHKERRQ(ierr); 1373 1374 mumps->id.perm_in = perm_in; 1375 for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */ 1376 ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr); 1377 } 1378 } 1379 } 1380 break; 1381 case 3: /* distributed assembled matrix input (size>1) */ 1382 mumps->id.nz_loc = mumps->nz; 1383 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1384 if (mumps->id.ICNTL(6)>1) { 1385 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1386 } 1387 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1388 if (!mumps->myid) { 1389 ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr); 1390 ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr); 1391 } else { 1392 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1393 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1394 } 1395 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1396 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1397 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1398 ierr = VecDestroy(&b);CHKERRQ(ierr); 1399 break; 1400 } 1401 PetscMUMPS_c(&mumps->id); 1402 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 1403 1404 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1405 F->ops->solve = MatSolve_MUMPS; 1406 F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1407 F->ops->matsolve = MatMatSolve_MUMPS; 1408 PetscFunctionReturn(0); 1409 } 1410 1411 /* Note the Petsc r and c permutations are ignored */ 1412 #undef __FUNCT__ 1413 #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 1414 PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1415 { 1416 Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1417 PetscErrorCode ierr; 1418 Vec b; 1419 IS is_iden; 1420 const PetscInt M = A->rmap->N; 1421 1422 PetscFunctionBegin; 1423 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1424 1425 /* Set MUMPS options from the options database */ 1426 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1427 1428 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1429 1430 /* analysis phase */ 1431 /*----------------*/ 1432 mumps->id.job = JOB_FACTSYMBOLIC; 1433 mumps->id.n = M; 1434 switch (mumps->id.ICNTL(18)) { 1435 case 0: /* centralized assembled matrix input */ 1436 if (!mumps->myid) { 1437 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1438 if (mumps->id.ICNTL(6)>1) { 1439 mumps->id.a = (MumpsScalar*)mumps->val; 1440 } 1441 } 1442 break; 1443 case 3: /* distributed assembled matrix input (size>1) */ 1444 mumps->id.nz_loc = mumps->nz; 1445 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1446 if (mumps->id.ICNTL(6)>1) { 1447 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1448 } 1449 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1450 if (!mumps->myid) { 1451 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 1452 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 1453 } else { 1454 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1455 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1456 } 1457 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1458 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1459 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1460 ierr = VecDestroy(&b);CHKERRQ(ierr); 1461 break; 1462 } 1463 PetscMUMPS_c(&mumps->id); 1464 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 1465 1466 F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 1467 F->ops->solve = MatSolve_MUMPS; 1468 F->ops->solvetranspose = MatSolveTranspose_MUMPS; 1469 PetscFunctionReturn(0); 1470 } 1471 1472 /* Note the Petsc r permutation and factor info are ignored */ 1473 #undef __FUNCT__ 1474 #define __FUNCT__ "MatCholeskyFactorSymbolic_MUMPS" 1475 PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 1476 { 1477 Mat_MUMPS *mumps = (Mat_MUMPS*)F->spptr; 1478 PetscErrorCode ierr; 1479 Vec b; 1480 IS is_iden; 1481 const PetscInt M = A->rmap->N; 1482 1483 PetscFunctionBegin; 1484 mumps->matstruc = DIFFERENT_NONZERO_PATTERN; 1485 1486 /* Set MUMPS options from the options database */ 1487 ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr); 1488 1489 ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr); 1490 1491 /* analysis phase */ 1492 /*----------------*/ 1493 mumps->id.job = JOB_FACTSYMBOLIC; 1494 mumps->id.n = M; 1495 switch (mumps->id.ICNTL(18)) { 1496 case 0: /* centralized assembled matrix input */ 1497 if (!mumps->myid) { 1498 mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn; 1499 if (mumps->id.ICNTL(6)>1) { 1500 mumps->id.a = (MumpsScalar*)mumps->val; 1501 } 1502 } 1503 break; 1504 case 3: /* distributed assembled matrix input (size>1) */ 1505 mumps->id.nz_loc = mumps->nz; 1506 mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn; 1507 if (mumps->id.ICNTL(6)>1) { 1508 mumps->id.a_loc = (MumpsScalar*)mumps->val; 1509 } 1510 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 1511 if (!mumps->myid) { 1512 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr); 1513 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 1514 } else { 1515 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr); 1516 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 1517 } 1518 ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr); 1519 ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr); 1520 ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 1521 ierr = VecDestroy(&b);CHKERRQ(ierr); 1522 break; 1523 } 1524 PetscMUMPS_c(&mumps->id); 1525 if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1)); 1526 1527 F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 1528 F->ops->solve = MatSolve_MUMPS; 1529 F->ops->solvetranspose = MatSolve_MUMPS; 1530 F->ops->matsolve = MatMatSolve_MUMPS; 1531 #if defined(PETSC_USE_COMPLEX) 1532 F->ops->getinertia = NULL; 1533 #else 1534 F->ops->getinertia = MatGetInertia_SBAIJMUMPS; 1535 #endif 1536 PetscFunctionReturn(0); 1537 } 1538 1539 #undef __FUNCT__ 1540 #define __FUNCT__ "MatView_MUMPS" 1541 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 1542 { 1543 PetscErrorCode ierr; 1544 PetscBool iascii; 1545 PetscViewerFormat format; 1546 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 1547 1548 PetscFunctionBegin; 1549 /* check if matrix is mumps type */ 1550 if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 1551 1552 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1553 if (iascii) { 1554 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1555 if (format == PETSC_VIEWER_ASCII_INFO) { 1556 ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 1557 ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 1558 ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 1559 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr); 1560 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr); 1561 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr); 1562 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr); 1563 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr); 1564 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr); 1565 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (sequentia matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr); 1566 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr); 1567 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr); 1568 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr); 1569 if (mumps->id.ICNTL(11)>0) { 1570 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(4) (inf norm of input mat): %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr); 1571 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(5) (inf norm of solution): %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr); 1572 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(6) (inf norm of residual): %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr); 1573 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr); 1574 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(9) (error estimate): %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr); 1575 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr); 1576 } 1577 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr); 1578 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr); 1579 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr); 1580 /* ICNTL(15-17) not used */ 1581 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr); 1582 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr); 1583 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr); 1584 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr); 1585 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr); 1586 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr); 1587 1588 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr); 1589 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr); 1590 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr); 1591 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr); 1592 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(28) (use parallel or sequential ordering): %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr); 1593 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(29) (parallel ordering): %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr); 1594 1595 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(30) (user-specified set of entries in inv(A)): %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr); 1596 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(31) (factors is discarded in the solve phase): %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr); 1597 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(33) (compute determinant): %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr); 1598 1599 ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",mumps->id.CNTL(1));CHKERRQ(ierr); 1600 ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr); 1601 ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",mumps->id.CNTL(3));CHKERRQ(ierr); 1602 ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",mumps->id.CNTL(4));CHKERRQ(ierr); 1603 ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",mumps->id.CNTL(5));CHKERRQ(ierr); 1604 1605 /* infomation local to each processor */ 1606 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 1607 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 1608 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr); 1609 ierr = PetscViewerFlush(viewer); 1610 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 1611 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr); 1612 ierr = PetscViewerFlush(viewer); 1613 ierr = PetscViewerASCIIPrintf(viewer, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 1614 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr); 1615 ierr = PetscViewerFlush(viewer); 1616 1617 ierr = PetscViewerASCIIPrintf(viewer, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr); 1618 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr); 1619 ierr = PetscViewerFlush(viewer); 1620 1621 ierr = PetscViewerASCIIPrintf(viewer, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr); 1622 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr); 1623 ierr = PetscViewerFlush(viewer); 1624 1625 ierr = PetscViewerASCIIPrintf(viewer, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr); 1626 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr); 1627 ierr = PetscViewerFlush(viewer); 1628 1629 if (mumps->ninfo && mumps->ninfo <= 40){ 1630 PetscInt i; 1631 for (i=0; i<mumps->ninfo; i++){ 1632 ierr = PetscViewerASCIIPrintf(viewer, " INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr); 1633 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr); 1634 ierr = PetscViewerFlush(viewer); 1635 } 1636 } 1637 1638 1639 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 1640 1641 if (!mumps->myid) { /* information from the host */ 1642 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr); 1643 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr); 1644 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr); 1645 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); 1646 1647 ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr); 1648 ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr); 1649 ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr); 1650 ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr); 1651 ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr); 1652 ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr); 1653 ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr); 1654 ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr); 1655 ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr); 1656 ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr); 1657 ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr); 1658 ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr); 1659 ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr); 1660 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); 1661 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); 1662 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); 1663 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); 1664 ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr); 1665 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); 1666 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); 1667 ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr); 1668 ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr); 1669 ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr); 1670 ierr = PetscViewerASCIIPrintf(viewer," INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr); 1671 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); 1672 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); 1673 ierr = PetscViewerASCIIPrintf(viewer," INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr); 1674 ierr = PetscViewerASCIIPrintf(viewer," INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr); 1675 ierr = PetscViewerASCIIPrintf(viewer," INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr); 1676 } 1677 } 1678 } 1679 PetscFunctionReturn(0); 1680 } 1681 1682 #undef __FUNCT__ 1683 #define __FUNCT__ "MatGetInfo_MUMPS" 1684 PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 1685 { 1686 Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 1687 1688 PetscFunctionBegin; 1689 info->block_size = 1.0; 1690 info->nz_allocated = mumps->id.INFOG(20); 1691 info->nz_used = mumps->id.INFOG(20); 1692 info->nz_unneeded = 0.0; 1693 info->assemblies = 0.0; 1694 info->mallocs = 0.0; 1695 info->memory = 0.0; 1696 info->fill_ratio_given = 0; 1697 info->fill_ratio_needed = 0; 1698 info->factor_mallocs = 0; 1699 PetscFunctionReturn(0); 1700 } 1701 1702 /* -------------------------------------------------------------------------------------------*/ 1703 #undef __FUNCT__ 1704 #define __FUNCT__ "MatFactorSetSchurIS_MUMPS" 1705 PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is) 1706 { 1707 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1708 const PetscInt *idxs; 1709 PetscInt size,i; 1710 PetscErrorCode ierr; 1711 1712 PetscFunctionBegin; 1713 if (mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS parallel Schur complements not yet supported from PETSc\n"); 1714 ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 1715 if (mumps->id.size_schur != size) { 1716 ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr); 1717 mumps->id.size_schur = size; 1718 mumps->id.schur_lld = size; 1719 ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr); 1720 } 1721 ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 1722 ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr); 1723 /* MUMPS expects Fortran style indices */ 1724 for (i=0;i<size;i++) mumps->id.listvar_schur[i]++; 1725 ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 1726 if (size) { /* turn on Schur switch if we the set of indices is not empty */ 1727 if (F->factortype == MAT_FACTOR_LU) { 1728 mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */ 1729 } else { 1730 mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */ 1731 } 1732 /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */ 1733 mumps->id.ICNTL(26) = -1; 1734 } 1735 PetscFunctionReturn(0); 1736 } 1737 1738 /* -------------------------------------------------------------------------------------------*/ 1739 #undef __FUNCT__ 1740 #define __FUNCT__ "MatFactorCreateSchurComplement_MUMPS" 1741 PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S) 1742 { 1743 Mat St; 1744 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1745 PetscScalar *array; 1746 #if defined(PETSC_USE_COMPLEX) 1747 PetscScalar im = PetscSqrtScalar((PetscScalar)-1.0); 1748 #endif 1749 PetscErrorCode ierr; 1750 1751 PetscFunctionBegin; 1752 if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 1753 else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 1754 1755 ierr = MatCreate(PetscObjectComm((PetscObject)F),&St);CHKERRQ(ierr); 1756 ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr); 1757 ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr); 1758 ierr = MatSetUp(St);CHKERRQ(ierr); 1759 ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr); 1760 if (!mumps->sym) { /* MUMPS always return a full matrix */ 1761 if (mumps->id.ICNTL(19) == 1) { /* stored by rows */ 1762 PetscInt i,j,N=mumps->id.size_schur; 1763 for (i=0;i<N;i++) { 1764 for (j=0;j<N;j++) { 1765 #if !defined(PETSC_USE_COMPLEX) 1766 PetscScalar val = mumps->id.schur[i*N+j]; 1767 #else 1768 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1769 #endif 1770 array[j*N+i] = val; 1771 } 1772 } 1773 } else { /* stored by columns */ 1774 ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1775 } 1776 } else { /* either full or lower-triangular (not packed) */ 1777 if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */ 1778 PetscInt i,j,N=mumps->id.size_schur; 1779 for (i=0;i<N;i++) { 1780 for (j=i;j<N;j++) { 1781 #if !defined(PETSC_USE_COMPLEX) 1782 PetscScalar val = mumps->id.schur[i*N+j]; 1783 #else 1784 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1785 #endif 1786 array[i*N+j] = val; 1787 array[j*N+i] = val; 1788 } 1789 } 1790 } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */ 1791 ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr); 1792 } else { /* ICNTL(19) == 1 lower triangular stored by rows */ 1793 PetscInt i,j,N=mumps->id.size_schur; 1794 for (i=0;i<N;i++) { 1795 for (j=0;j<i+1;j++) { 1796 #if !defined(PETSC_USE_COMPLEX) 1797 PetscScalar val = mumps->id.schur[i*N+j]; 1798 #else 1799 PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i; 1800 #endif 1801 array[i*N+j] = val; 1802 array[j*N+i] = val; 1803 } 1804 } 1805 } 1806 } 1807 ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr); 1808 *S = St; 1809 PetscFunctionReturn(0); 1810 } 1811 1812 /* -------------------------------------------------------------------------------------------*/ 1813 #undef __FUNCT__ 1814 #define __FUNCT__ "MatFactorGetSchurComplement_MUMPS" 1815 PetscErrorCode MatFactorGetSchurComplement_MUMPS(Mat F,Mat* S) 1816 { 1817 Mat St; 1818 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1819 PetscErrorCode ierr; 1820 1821 PetscFunctionBegin; 1822 if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 1823 else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 1824 1825 /* 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 */ 1826 ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)F),mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&St);CHKERRQ(ierr); 1827 *S = St; 1828 PetscFunctionReturn(0); 1829 } 1830 1831 /* -------------------------------------------------------------------------------------------*/ 1832 #undef __FUNCT__ 1833 #define __FUNCT__ "MatFactorInvertSchurComplement_MUMPS" 1834 PetscErrorCode MatFactorInvertSchurComplement_MUMPS(Mat F) 1835 { 1836 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1837 PetscErrorCode ierr; 1838 1839 PetscFunctionBegin; 1840 if (!mumps->id.ICNTL(19)) { /* do nothing */ 1841 PetscFunctionReturn(0); 1842 } 1843 if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 1844 ierr = MatMumpsInvertSchur_Private(mumps);CHKERRQ(ierr); 1845 PetscFunctionReturn(0); 1846 } 1847 1848 /* -------------------------------------------------------------------------------------------*/ 1849 #undef __FUNCT__ 1850 #define __FUNCT__ "MatFactorSolveSchurComplement_MUMPS" 1851 PetscErrorCode MatFactorSolveSchurComplement_MUMPS(Mat F, Vec rhs, Vec sol) 1852 { 1853 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1854 MumpsScalar *orhs; 1855 PetscScalar *osol,*nrhs,*nsol; 1856 PetscInt orhs_size,osol_size,olrhs_size; 1857 PetscErrorCode ierr; 1858 1859 PetscFunctionBegin; 1860 if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 1861 if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 1862 1863 /* swap pointers */ 1864 orhs = mumps->id.redrhs; 1865 olrhs_size = mumps->id.lredrhs; 1866 orhs_size = mumps->sizeredrhs; 1867 osol = mumps->schur_sol; 1868 osol_size = mumps->schur_sizesol; 1869 ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr); 1870 ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr); 1871 mumps->id.redrhs = (MumpsScalar*)nrhs; 1872 ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr); 1873 mumps->id.lredrhs = mumps->sizeredrhs; 1874 mumps->schur_sol = nsol; 1875 ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr); 1876 1877 /* solve Schur complement */ 1878 mumps->id.nrhs = 1; 1879 ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 1880 /* restore pointers */ 1881 ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr); 1882 ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr); 1883 mumps->id.redrhs = orhs; 1884 mumps->id.lredrhs = olrhs_size; 1885 mumps->sizeredrhs = orhs_size; 1886 mumps->schur_sol = osol; 1887 mumps->schur_sizesol = osol_size; 1888 PetscFunctionReturn(0); 1889 } 1890 1891 /* -------------------------------------------------------------------------------------------*/ 1892 #undef __FUNCT__ 1893 #define __FUNCT__ "MatFactorSolveSchurComplementTranspose_MUMPS" 1894 PetscErrorCode MatFactorSolveSchurComplementTranspose_MUMPS(Mat F, Vec rhs, Vec sol) 1895 { 1896 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1897 MumpsScalar *orhs; 1898 PetscScalar *osol,*nrhs,*nsol; 1899 PetscInt orhs_size,osol_size; 1900 PetscErrorCode ierr; 1901 1902 PetscFunctionBegin; 1903 if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 1904 else if (!mumps->id.size_schur) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 1905 1906 /* swap pointers */ 1907 orhs = mumps->id.redrhs; 1908 orhs_size = mumps->sizeredrhs; 1909 osol = mumps->schur_sol; 1910 osol_size = mumps->schur_sizesol; 1911 ierr = VecGetArray(rhs,&nrhs);CHKERRQ(ierr); 1912 ierr = VecGetArray(sol,&nsol);CHKERRQ(ierr); 1913 mumps->id.redrhs = (MumpsScalar*)nrhs; 1914 ierr = VecGetLocalSize(rhs,&mumps->sizeredrhs);CHKERRQ(ierr); 1915 mumps->schur_sol = nsol; 1916 ierr = VecGetLocalSize(sol,&mumps->schur_sizesol);CHKERRQ(ierr); 1917 1918 /* solve Schur complement */ 1919 mumps->id.nrhs = 1; 1920 mumps->id.ICNTL(9) = 0; 1921 ierr = MatMumpsSolveSchur_Private(mumps,PETSC_FALSE);CHKERRQ(ierr); 1922 mumps->id.ICNTL(9) = 1; 1923 /* restore pointers */ 1924 ierr = VecRestoreArray(rhs,&nrhs);CHKERRQ(ierr); 1925 ierr = VecRestoreArray(sol,&nsol);CHKERRQ(ierr); 1926 mumps->id.redrhs = orhs; 1927 mumps->sizeredrhs = orhs_size; 1928 mumps->schur_sol = osol; 1929 mumps->schur_sizesol = osol_size; 1930 PetscFunctionReturn(0); 1931 } 1932 1933 /* -------------------------------------------------------------------------------------------*/ 1934 #undef __FUNCT__ 1935 #define __FUNCT__ "MatMumpsSetIcntl_MUMPS" 1936 PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival) 1937 { 1938 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1939 1940 PetscFunctionBegin; 1941 mumps->id.ICNTL(icntl) = ival; 1942 PetscFunctionReturn(0); 1943 } 1944 1945 #undef __FUNCT__ 1946 #define __FUNCT__ "MatMumpsGetIcntl_MUMPS" 1947 PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival) 1948 { 1949 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 1950 1951 PetscFunctionBegin; 1952 *ival = mumps->id.ICNTL(icntl); 1953 PetscFunctionReturn(0); 1954 } 1955 1956 #undef __FUNCT__ 1957 #define __FUNCT__ "MatMumpsSetIcntl" 1958 /*@ 1959 MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 1960 1961 Logically Collective on Mat 1962 1963 Input Parameters: 1964 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1965 . icntl - index of MUMPS parameter array ICNTL() 1966 - ival - value of MUMPS ICNTL(icntl) 1967 1968 Options Database: 1969 . -mat_mumps_icntl_<icntl> <ival> 1970 1971 Level: beginner 1972 1973 References: MUMPS Users' Guide 1974 1975 .seealso: MatGetFactor() 1976 @*/ 1977 PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival) 1978 { 1979 PetscErrorCode ierr; 1980 1981 PetscFunctionBegin; 1982 PetscValidLogicalCollectiveInt(F,icntl,2); 1983 PetscValidLogicalCollectiveInt(F,ival,3); 1984 ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 1985 PetscFunctionReturn(0); 1986 } 1987 1988 #undef __FUNCT__ 1989 #define __FUNCT__ "MatMumpsGetIcntl" 1990 /*@ 1991 MatMumpsGetIcntl - Get MUMPS parameter ICNTL() 1992 1993 Logically Collective on Mat 1994 1995 Input Parameters: 1996 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 1997 - icntl - index of MUMPS parameter array ICNTL() 1998 1999 Output Parameter: 2000 . ival - value of MUMPS ICNTL(icntl) 2001 2002 Level: beginner 2003 2004 References: MUMPS Users' Guide 2005 2006 .seealso: MatGetFactor() 2007 @*/ 2008 PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival) 2009 { 2010 PetscErrorCode ierr; 2011 2012 PetscFunctionBegin; 2013 PetscValidLogicalCollectiveInt(F,icntl,2); 2014 PetscValidIntPointer(ival,3); 2015 ierr = PetscTryMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2016 PetscFunctionReturn(0); 2017 } 2018 2019 /* -------------------------------------------------------------------------------------------*/ 2020 #undef __FUNCT__ 2021 #define __FUNCT__ "MatMumpsSetCntl_MUMPS" 2022 PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val) 2023 { 2024 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2025 2026 PetscFunctionBegin; 2027 mumps->id.CNTL(icntl) = val; 2028 PetscFunctionReturn(0); 2029 } 2030 2031 #undef __FUNCT__ 2032 #define __FUNCT__ "MatMumpsGetCntl_MUMPS" 2033 PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val) 2034 { 2035 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2036 2037 PetscFunctionBegin; 2038 *val = mumps->id.CNTL(icntl); 2039 PetscFunctionReturn(0); 2040 } 2041 2042 #undef __FUNCT__ 2043 #define __FUNCT__ "MatMumpsSetCntl" 2044 /*@ 2045 MatMumpsSetCntl - Set MUMPS parameter CNTL() 2046 2047 Logically Collective on Mat 2048 2049 Input Parameters: 2050 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2051 . icntl - index of MUMPS parameter array CNTL() 2052 - val - value of MUMPS CNTL(icntl) 2053 2054 Options Database: 2055 . -mat_mumps_cntl_<icntl> <val> 2056 2057 Level: beginner 2058 2059 References: MUMPS Users' Guide 2060 2061 .seealso: MatGetFactor() 2062 @*/ 2063 PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val) 2064 { 2065 PetscErrorCode ierr; 2066 2067 PetscFunctionBegin; 2068 PetscValidLogicalCollectiveInt(F,icntl,2); 2069 PetscValidLogicalCollectiveReal(F,val,3); 2070 ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr); 2071 PetscFunctionReturn(0); 2072 } 2073 2074 #undef __FUNCT__ 2075 #define __FUNCT__ "MatMumpsGetCntl" 2076 /*@ 2077 MatMumpsGetCntl - Get MUMPS parameter CNTL() 2078 2079 Logically Collective on Mat 2080 2081 Input Parameters: 2082 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2083 - icntl - index of MUMPS parameter array CNTL() 2084 2085 Output Parameter: 2086 . val - value of MUMPS CNTL(icntl) 2087 2088 Level: beginner 2089 2090 References: MUMPS Users' Guide 2091 2092 .seealso: MatGetFactor() 2093 @*/ 2094 PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val) 2095 { 2096 PetscErrorCode ierr; 2097 2098 PetscFunctionBegin; 2099 PetscValidLogicalCollectiveInt(F,icntl,2); 2100 PetscValidRealPointer(val,3); 2101 ierr = PetscTryMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2102 PetscFunctionReturn(0); 2103 } 2104 2105 #undef __FUNCT__ 2106 #define __FUNCT__ "MatMumpsGetInfo_MUMPS" 2107 PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info) 2108 { 2109 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2110 2111 PetscFunctionBegin; 2112 *info = mumps->id.INFO(icntl); 2113 PetscFunctionReturn(0); 2114 } 2115 2116 #undef __FUNCT__ 2117 #define __FUNCT__ "MatMumpsGetInfog_MUMPS" 2118 PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog) 2119 { 2120 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2121 2122 PetscFunctionBegin; 2123 *infog = mumps->id.INFOG(icntl); 2124 PetscFunctionReturn(0); 2125 } 2126 2127 #undef __FUNCT__ 2128 #define __FUNCT__ "MatMumpsGetRinfo_MUMPS" 2129 PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo) 2130 { 2131 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2132 2133 PetscFunctionBegin; 2134 *rinfo = mumps->id.RINFO(icntl); 2135 PetscFunctionReturn(0); 2136 } 2137 2138 #undef __FUNCT__ 2139 #define __FUNCT__ "MatMumpsGetRinfog_MUMPS" 2140 PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog) 2141 { 2142 Mat_MUMPS *mumps =(Mat_MUMPS*)F->spptr; 2143 2144 PetscFunctionBegin; 2145 *rinfog = mumps->id.RINFOG(icntl); 2146 PetscFunctionReturn(0); 2147 } 2148 2149 #undef __FUNCT__ 2150 #define __FUNCT__ "MatMumpsGetInfo" 2151 /*@ 2152 MatMumpsGetInfo - Get MUMPS parameter INFO() 2153 2154 Logically Collective on Mat 2155 2156 Input Parameters: 2157 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2158 - icntl - index of MUMPS parameter array INFO() 2159 2160 Output Parameter: 2161 . ival - value of MUMPS INFO(icntl) 2162 2163 Level: beginner 2164 2165 References: MUMPS Users' Guide 2166 2167 .seealso: MatGetFactor() 2168 @*/ 2169 PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival) 2170 { 2171 PetscErrorCode ierr; 2172 2173 PetscFunctionBegin; 2174 PetscValidIntPointer(ival,3); 2175 ierr = PetscTryMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2176 PetscFunctionReturn(0); 2177 } 2178 2179 #undef __FUNCT__ 2180 #define __FUNCT__ "MatMumpsGetInfog" 2181 /*@ 2182 MatMumpsGetInfog - Get MUMPS parameter INFOG() 2183 2184 Logically Collective on Mat 2185 2186 Input Parameters: 2187 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2188 - icntl - index of MUMPS parameter array INFOG() 2189 2190 Output Parameter: 2191 . ival - value of MUMPS INFOG(icntl) 2192 2193 Level: beginner 2194 2195 References: MUMPS Users' Guide 2196 2197 .seealso: MatGetFactor() 2198 @*/ 2199 PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival) 2200 { 2201 PetscErrorCode ierr; 2202 2203 PetscFunctionBegin; 2204 PetscValidIntPointer(ival,3); 2205 ierr = PetscTryMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr); 2206 PetscFunctionReturn(0); 2207 } 2208 2209 #undef __FUNCT__ 2210 #define __FUNCT__ "MatMumpsGetRinfo" 2211 /*@ 2212 MatMumpsGetRinfo - Get MUMPS parameter RINFO() 2213 2214 Logically Collective on Mat 2215 2216 Input Parameters: 2217 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2218 - icntl - index of MUMPS parameter array RINFO() 2219 2220 Output Parameter: 2221 . val - value of MUMPS RINFO(icntl) 2222 2223 Level: beginner 2224 2225 References: MUMPS Users' Guide 2226 2227 .seealso: MatGetFactor() 2228 @*/ 2229 PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val) 2230 { 2231 PetscErrorCode ierr; 2232 2233 PetscFunctionBegin; 2234 PetscValidRealPointer(val,3); 2235 ierr = PetscTryMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2236 PetscFunctionReturn(0); 2237 } 2238 2239 #undef __FUNCT__ 2240 #define __FUNCT__ "MatMumpsGetRinfog" 2241 /*@ 2242 MatMumpsGetRinfog - Get MUMPS parameter RINFOG() 2243 2244 Logically Collective on Mat 2245 2246 Input Parameters: 2247 + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 2248 - icntl - index of MUMPS parameter array RINFOG() 2249 2250 Output Parameter: 2251 . val - value of MUMPS RINFOG(icntl) 2252 2253 Level: beginner 2254 2255 References: MUMPS Users' Guide 2256 2257 .seealso: MatGetFactor() 2258 @*/ 2259 PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val) 2260 { 2261 PetscErrorCode ierr; 2262 2263 PetscFunctionBegin; 2264 PetscValidRealPointer(val,3); 2265 ierr = PetscTryMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr); 2266 PetscFunctionReturn(0); 2267 } 2268 2269 /*MC 2270 MATSOLVERMUMPS - A matrix type providing direct solvers (LU and Cholesky) for 2271 distributed and sequential matrices via the external package MUMPS. 2272 2273 Works with MATAIJ and MATSBAIJ matrices 2274 2275 Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch to have PETSc installed with MUMPS 2276 2277 Use -pc_type cholesky or lu -pc_factor_mat_solver_package mumps to us this direct solver 2278 2279 Options Database Keys: 2280 + -mat_mumps_icntl_1 <6>: ICNTL(1): output stream for error messages (None) 2281 . -mat_mumps_icntl_2 <0>: ICNTL(2): output stream for diagnostic printing, statistics, and warning (None) 2282 . -mat_mumps_icntl_3 <0>: ICNTL(3): output stream for global information, collected on the host (None) 2283 . -mat_mumps_icntl_4 <0>: ICNTL(4): level of printing (0 to 4) (None) 2284 . -mat_mumps_icntl_6 <7>: ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7) (None) 2285 . -mat_mumps_icntl_7 <7>: ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis (None) 2286 . -mat_mumps_icntl_8 <77>: ICNTL(8): scaling strategy (-2 to 8 or 77) (None) 2287 . -mat_mumps_icntl_10 <0>: ICNTL(10): max num of refinements (None) 2288 . -mat_mumps_icntl_11 <0>: ICNTL(11): statistics related to an error analysis (via -ksp_view) (None) 2289 . -mat_mumps_icntl_12 <1>: ICNTL(12): an ordering strategy for symmetric matrices (0 to 3) (None) 2290 . -mat_mumps_icntl_13 <0>: ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting (None) 2291 . -mat_mumps_icntl_14 <20>: ICNTL(14): percentage increase in the estimated working space (None) 2292 . -mat_mumps_icntl_19 <0>: ICNTL(19): computes the Schur complement (None) 2293 . -mat_mumps_icntl_22 <0>: ICNTL(22): in-core/out-of-core factorization and solve (0 or 1) (None) 2294 . -mat_mumps_icntl_23 <0>: ICNTL(23): max size of the working memory (MB) that can allocate per processor (None) 2295 . -mat_mumps_icntl_24 <0>: ICNTL(24): detection of null pivot rows (0 or 1) (None) 2296 . -mat_mumps_icntl_25 <0>: ICNTL(25): compute a solution of a deficient matrix and a null space basis (None) 2297 . -mat_mumps_icntl_26 <0>: ICNTL(26): drives the solution phase if a Schur complement matrix (None) 2298 . -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) 2299 . -mat_mumps_icntl_29 <0>: ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis (None) 2300 . -mat_mumps_icntl_30 <0>: ICNTL(30): compute user-specified set of entries in inv(A) (None) 2301 . -mat_mumps_icntl_31 <0>: ICNTL(31): indicates which factors may be discarded during factorization (None) 2302 . -mat_mumps_icntl_33 <0>: ICNTL(33): compute determinant (None) 2303 . -mat_mumps_cntl_1 <0.01>: CNTL(1): relative pivoting threshold (None) 2304 . -mat_mumps_cntl_2 <1.49012e-08>: CNTL(2): stopping criterion of refinement (None) 2305 . -mat_mumps_cntl_3 <0>: CNTL(3): absolute pivoting threshold (None) 2306 . -mat_mumps_cntl_4 <-1>: CNTL(4): value for static pivoting (None) 2307 - -mat_mumps_cntl_5 <0>: CNTL(5): fixation for null pivots (None) 2308 2309 Level: beginner 2310 2311 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 2312 2313 M*/ 2314 2315 #undef __FUNCT__ 2316 #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 2317 static PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 2318 { 2319 PetscFunctionBegin; 2320 *type = MATSOLVERMUMPS; 2321 PetscFunctionReturn(0); 2322 } 2323 2324 /* MatGetFactor for Seq and MPI AIJ matrices */ 2325 #undef __FUNCT__ 2326 #define __FUNCT__ "MatGetFactor_aij_mumps" 2327 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F) 2328 { 2329 Mat B; 2330 PetscErrorCode ierr; 2331 Mat_MUMPS *mumps; 2332 PetscBool isSeqAIJ; 2333 2334 PetscFunctionBegin; 2335 /* Create the factorization matrix */ 2336 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 2337 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2338 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2339 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2340 if (isSeqAIJ) { 2341 ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 2342 } else { 2343 ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr); 2344 } 2345 2346 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2347 2348 B->ops->view = MatView_MUMPS; 2349 B->ops->getinfo = MatGetInfo_MUMPS; 2350 B->ops->getdiagonal = MatGetDiagonal_MUMPS; 2351 2352 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 2353 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2354 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr); 2355 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2356 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr); 2357 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr); 2358 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr); 2359 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2360 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2361 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2362 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2363 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2364 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2365 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2366 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2367 2368 if (ftype == MAT_FACTOR_LU) { 2369 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 2370 B->factortype = MAT_FACTOR_LU; 2371 if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij; 2372 else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij; 2373 mumps->sym = 0; 2374 } else { 2375 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2376 B->factortype = MAT_FACTOR_CHOLESKY; 2377 if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij; 2378 else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij; 2379 #if defined(PETSC_USE_COMPLEX) 2380 mumps->sym = 2; 2381 #else 2382 if (A->spd_set && A->spd) mumps->sym = 1; 2383 else mumps->sym = 2; 2384 #endif 2385 } 2386 2387 mumps->isAIJ = PETSC_TRUE; 2388 mumps->Destroy = B->ops->destroy; 2389 B->ops->destroy = MatDestroy_MUMPS; 2390 B->spptr = (void*)mumps; 2391 2392 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2393 2394 *F = B; 2395 PetscFunctionReturn(0); 2396 } 2397 2398 /* MatGetFactor for Seq and MPI SBAIJ matrices */ 2399 #undef __FUNCT__ 2400 #define __FUNCT__ "MatGetFactor_sbaij_mumps" 2401 PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 2402 { 2403 Mat B; 2404 PetscErrorCode ierr; 2405 Mat_MUMPS *mumps; 2406 PetscBool isSeqSBAIJ; 2407 2408 PetscFunctionBegin; 2409 if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 2410 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"); 2411 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 2412 /* Create the factorization matrix */ 2413 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2414 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2415 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2416 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2417 if (isSeqSBAIJ) { 2418 ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr); 2419 2420 mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij; 2421 } else { 2422 ierr = MatMPISBAIJSetPreallocation(B,1,0,NULL,0,NULL);CHKERRQ(ierr); 2423 2424 mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij; 2425 } 2426 2427 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS; 2428 B->ops->view = MatView_MUMPS; 2429 B->ops->getdiagonal = MatGetDiagonal_MUMPS; 2430 2431 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 2432 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2433 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr); 2434 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2435 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr); 2436 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr); 2437 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr); 2438 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2439 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2440 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2441 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2442 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2443 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2444 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2445 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2446 2447 B->factortype = MAT_FACTOR_CHOLESKY; 2448 #if defined(PETSC_USE_COMPLEX) 2449 mumps->sym = 2; 2450 #else 2451 if (A->spd_set && A->spd) mumps->sym = 1; 2452 else mumps->sym = 2; 2453 #endif 2454 2455 mumps->isAIJ = PETSC_FALSE; 2456 mumps->Destroy = B->ops->destroy; 2457 B->ops->destroy = MatDestroy_MUMPS; 2458 B->spptr = (void*)mumps; 2459 2460 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2461 2462 *F = B; 2463 PetscFunctionReturn(0); 2464 } 2465 2466 #undef __FUNCT__ 2467 #define __FUNCT__ "MatGetFactor_baij_mumps" 2468 PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F) 2469 { 2470 Mat B; 2471 PetscErrorCode ierr; 2472 Mat_MUMPS *mumps; 2473 PetscBool isSeqBAIJ; 2474 2475 PetscFunctionBegin; 2476 /* Create the factorization matrix */ 2477 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 2478 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 2479 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 2480 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 2481 if (isSeqBAIJ) { 2482 ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 2483 } else { 2484 ierr = MatMPIBAIJSetPreallocation(B,A->rmap->bs,0,NULL,0,NULL);CHKERRQ(ierr); 2485 } 2486 2487 ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr); 2488 if (ftype == MAT_FACTOR_LU) { 2489 B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 2490 B->factortype = MAT_FACTOR_LU; 2491 if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij; 2492 else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij; 2493 mumps->sym = 0; 2494 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n"); 2495 2496 B->ops->view = MatView_MUMPS; 2497 B->ops->getdiagonal = MatGetDiagonal_MUMPS; 2498 2499 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 2500 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr); 2501 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MUMPS);CHKERRQ(ierr); 2502 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr); 2503 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MUMPS);CHKERRQ(ierr); 2504 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MUMPS);CHKERRQ(ierr); 2505 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MUMPS);CHKERRQ(ierr); 2506 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr); 2507 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr); 2508 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr); 2509 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr); 2510 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr); 2511 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr); 2512 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr); 2513 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr); 2514 2515 mumps->isAIJ = PETSC_TRUE; 2516 mumps->Destroy = B->ops->destroy; 2517 B->ops->destroy = MatDestroy_MUMPS; 2518 B->spptr = (void*)mumps; 2519 2520 ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr); 2521 2522 *F = B; 2523 PetscFunctionReturn(0); 2524 } 2525 2526 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mumps(Mat,MatFactorType,Mat*); 2527 PETSC_EXTERN PetscErrorCode MatGetFactor_baij_mumps(Mat,MatFactorType,Mat*); 2528 PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mumps(Mat,MatFactorType,Mat*); 2529 2530 #undef __FUNCT__ 2531 #define __FUNCT__ "MatSolverPackageRegister_MUMPS" 2532 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MUMPS(void) 2533 { 2534 PetscErrorCode ierr; 2535 2536 PetscFunctionBegin; 2537 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2538 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2539 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2540 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2541 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 2542 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2543 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr); 2544 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2545 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr); 2546 ierr = MatSolverPackageRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr); 2547 PetscFunctionReturn(0); 2548 } 2549 2550