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