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