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