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