1 #define PETSCMAT_DLL 2 3 /* 4 Provides an interface to the MUMPS sparse solver 5 */ 6 #include "src/mat/impls/aij/seq/aij.h" 7 #include "src/mat/impls/aij/mpi/mpiaij.h" 8 #include "src/mat/impls/sbaij/seq/sbaij.h" 9 #include "src/mat/impls/sbaij/mpi/mpisbaij.h" 10 11 EXTERN_C_BEGIN 12 #if defined(PETSC_USE_COMPLEX) 13 #include "zmumps_c.h" 14 #else 15 #include "dmumps_c.h" 16 #endif 17 EXTERN_C_END 18 #define JOB_INIT -1 19 #define JOB_END -2 20 /* macros s.t. indices match MUMPS documentation */ 21 #define ICNTL(I) icntl[(I)-1] 22 #define CNTL(I) cntl[(I)-1] 23 #define INFOG(I) infog[(I)-1] 24 #define INFO(I) info[(I)-1] 25 #define RINFOG(I) rinfog[(I)-1] 26 #define RINFO(I) rinfo[(I)-1] 27 28 typedef struct { 29 #if defined(PETSC_USE_COMPLEX) 30 ZMUMPS_STRUC_C id; 31 #else 32 DMUMPS_STRUC_C id; 33 #endif 34 MatStructure matstruc; 35 PetscMPIInt myid,size; 36 PetscInt *irn,*jcn,sym,nSolve; 37 PetscScalar *val; 38 MPI_Comm comm_mumps; 39 VecScatter scat_rhs, scat_sol; 40 PetscTruth isAIJ,CleanUpMUMPS; 41 Vec b_seq,x_seq; 42 PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 43 PetscErrorCode (*MatView)(Mat,PetscViewer); 44 PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 45 PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 46 PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*); 47 PetscErrorCode (*MatDestroy)(Mat); 48 PetscErrorCode (*specialdestroy)(Mat); 49 PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*); 50 } Mat_MUMPS; 51 52 EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 53 EXTERN_C_BEGIN 54 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*); 55 EXTERN_C_END 56 /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 57 /* 58 input: 59 A - matrix in mpiaij or mpisbaij (bs=1) format 60 shift - 0: C style output triple; 1: Fortran style output triple. 61 valOnly - FALSE: spaces are allocated and values are set for the triple 62 TRUE: only the values in v array are updated 63 output: 64 nnz - dim of r, c, and v (number of local nonzero entries of A) 65 r, c, v - row and col index, matrix values (matrix triples) 66 */ 67 PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) { 68 PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray; 69 PetscErrorCode ierr; 70 PetscInt i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol; 71 PetscInt *row,*col; 72 PetscScalar *av, *bv,*val; 73 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 74 75 PetscFunctionBegin; 76 if (mumps->isAIJ){ 77 Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 78 Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 79 Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 80 nz = aa->nz + bb->nz; 81 ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart; 82 garray = mat->garray; 83 av=aa->a; bv=bb->a; 84 85 } else { 86 Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 87 Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 88 Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 89 if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs); 90 nz = aa->nz + bb->nz; 91 ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart; 92 garray = mat->garray; 93 av=aa->a; bv=bb->a; 94 } 95 96 if (!valOnly){ 97 ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr); 98 ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr); 99 ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 100 *r = row; *c = col; *v = val; 101 } else { 102 row = *r; col = *c; val = *v; 103 } 104 *nnz = nz; 105 106 jj = 0; irow = rstart; 107 for ( i=0; i<m; i++ ) { 108 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 109 countA = ai[i+1] - ai[i]; 110 countB = bi[i+1] - bi[i]; 111 bjj = bj + bi[i]; 112 113 /* get jB, the starting local col index for the 2nd B-part */ 114 colA_start = rstart + ajj[0]; /* the smallest col index for A */ 115 j=-1; 116 do { 117 j++; 118 if (j == countB) break; 119 jcol = garray[bjj[j]]; 120 } while (jcol < colA_start); 121 jB = j; 122 123 /* B-part, smaller col index */ 124 colA_start = rstart + ajj[0]; /* the smallest col index for A */ 125 for (j=0; j<jB; j++){ 126 jcol = garray[bjj[j]]; 127 if (!valOnly){ 128 row[jj] = irow + shift; col[jj] = jcol + shift; 129 130 } 131 val[jj++] = *bv++; 132 } 133 /* A-part */ 134 for (j=0; j<countA; j++){ 135 if (!valOnly){ 136 row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 137 } 138 val[jj++] = *av++; 139 } 140 /* B-part, larger col index */ 141 for (j=jB; j<countB; j++){ 142 if (!valOnly){ 143 row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 144 } 145 val[jj++] = *bv++; 146 } 147 irow++; 148 } 149 150 PetscFunctionReturn(0); 151 } 152 153 EXTERN_C_BEGIN 154 #undef __FUNCT__ 155 #define __FUNCT__ "MatConvert_MUMPS_Base" 156 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat) 157 { 158 PetscErrorCode ierr; 159 Mat B=*newmat; 160 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 161 void (*f)(void); 162 163 PetscFunctionBegin; 164 if (reuse == MAT_INITIAL_MATRIX) { 165 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 166 } 167 B->ops->duplicate = mumps->MatDuplicate; 168 B->ops->view = mumps->MatView; 169 B->ops->assemblyend = mumps->MatAssemblyEnd; 170 B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic; 171 B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic; 172 B->ops->destroy = mumps->MatDestroy; 173 174 /* put back original composed preallocation function */ 175 ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 176 if (f) { 177 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr); 178 } 179 ierr = PetscFree(mumps);CHKERRQ(ierr); 180 A->spptr = PETSC_NULL; 181 182 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 183 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 184 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 185 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr); 186 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 187 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr); 188 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 189 ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr); 190 191 ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr); 192 *newmat = B; 193 PetscFunctionReturn(0); 194 } 195 EXTERN_C_END 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "MatDestroy_MUMPS" 199 PetscErrorCode MatDestroy_MUMPS(Mat A) 200 { 201 Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 202 PetscErrorCode ierr; 203 PetscMPIInt size=lu->size; 204 PetscErrorCode (*specialdestroy)(Mat); 205 PetscFunctionBegin; 206 if (lu->CleanUpMUMPS) { 207 /* Terminate instance, deallocate memories */ 208 if (size > 1){ 209 ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 210 ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 211 ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 212 if (lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);} 213 if (lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);} 214 ierr = PetscFree(lu->val);CHKERRQ(ierr); 215 } 216 lu->id.job=JOB_END; 217 #if defined(PETSC_USE_COMPLEX) 218 zmumps_c(&lu->id); 219 #else 220 dmumps_c(&lu->id); 221 #endif 222 ierr = PetscFree(lu->irn);CHKERRQ(ierr); 223 ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 224 ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 225 } 226 specialdestroy = lu->specialdestroy; 227 ierr = (*specialdestroy)(A);CHKERRQ(ierr); 228 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 229 PetscFunctionReturn(0); 230 } 231 232 #undef __FUNCT__ 233 #define __FUNCT__ "MatDestroy_AIJMUMPS" 234 PetscErrorCode MatDestroy_AIJMUMPS(Mat A) 235 { 236 PetscErrorCode ierr; 237 PetscMPIInt size; 238 239 PetscFunctionBegin; 240 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 241 if (size==1) { 242 ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 243 } else { 244 ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 245 } 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "MatDestroy_SBAIJMUMPS" 251 PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A) 252 { 253 PetscErrorCode ierr; 254 PetscMPIInt size; 255 256 PetscFunctionBegin; 257 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 258 if (size==1) { 259 ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 260 } else { 261 ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 262 } 263 PetscFunctionReturn(0); 264 } 265 266 #undef __FUNCT__ 267 #define __FUNCT__ "MatSolve_MUMPS" 268 PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) { 269 Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 270 PetscScalar *array; 271 Vec x_seq; 272 IS is_iden,is_petsc; 273 PetscErrorCode ierr; 274 PetscInt i; 275 276 PetscFunctionBegin; 277 lu->id.nrhs = 1; 278 x_seq = lu->b_seq; 279 if (lu->size > 1){ 280 /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 281 ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 282 ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 283 if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 284 } else { /* size == 1 */ 285 ierr = VecCopy(b,x);CHKERRQ(ierr); 286 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 287 } 288 if (!lu->myid) { /* define rhs on the host */ 289 #if defined(PETSC_USE_COMPLEX) 290 lu->id.rhs = (mumps_double_complex*)array; 291 #else 292 lu->id.rhs = array; 293 #endif 294 } 295 if (lu->size == 1){ 296 ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 297 } else if (!lu->myid){ 298 ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 299 } 300 301 if (lu->size > 1){ 302 /* distributed solution */ 303 lu->id.ICNTL(21) = 1; 304 if (!lu->nSolve){ 305 /* Create x_seq=sol_loc for repeated use */ 306 PetscInt lsol_loc; 307 PetscScalar *sol_loc; 308 lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 309 ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr); 310 lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc); 311 lu->id.lsol_loc = lsol_loc; 312 #if defined(PETSC_USE_COMPLEX) 313 lu->id.sol_loc = (ZMUMPS_DOUBLE *)sol_loc; 314 #else 315 lu->id.sol_loc = (DMUMPS_DOUBLE *)sol_loc; 316 #endif 317 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 318 } 319 } 320 321 /* solve phase */ 322 /*-------------*/ 323 lu->id.job = 3; 324 #if defined(PETSC_USE_COMPLEX) 325 zmumps_c(&lu->id); 326 #else 327 dmumps_c(&lu->id); 328 #endif 329 if (lu->id.INFOG(1) < 0) { 330 SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 331 } 332 333 if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 334 if (!lu->nSolve){ /* create scatter scat_sol */ 335 ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 336 for (i=0; i<lu->id.lsol_loc; i++){ 337 lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 338 } 339 ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr); /* to */ 340 ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 341 ierr = ISDestroy(is_iden);CHKERRQ(ierr); 342 ierr = ISDestroy(is_petsc);CHKERRQ(ierr); 343 } 344 ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 345 ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 346 } 347 lu->nSolve++; 348 PetscFunctionReturn(0); 349 } 350 351 /* 352 input: 353 F: numeric factor 354 output: 355 nneg: total number of negative pivots 356 nzero: 0 357 npos: (global dimension of F) - nneg 358 */ 359 360 #undef __FUNCT__ 361 #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 362 PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 363 { 364 Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 365 PetscErrorCode ierr; 366 PetscMPIInt size; 367 368 PetscFunctionBegin; 369 ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 370 /* 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 */ 371 if (size > 1 && lu->id.ICNTL(13) != 1){ 372 SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13)); 373 } 374 if (nneg){ 375 if (!lu->myid){ 376 *nneg = lu->id.INFOG(12); 377 } 378 ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 379 } 380 if (nzero) *nzero = 0; 381 if (npos) *npos = F->rmap.N - (*nneg); 382 PetscFunctionReturn(0); 383 } 384 385 #undef __FUNCT__ 386 #define __FUNCT__ "MatFactorNumeric_MUMPS" 387 PetscErrorCode MatFactorNumeric_MUMPS(Mat A,MatFactorInfo *info,Mat *F) 388 { 389 Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 390 Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr; 391 PetscErrorCode ierr; 392 PetscInt rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl; 393 PetscTruth valOnly,flg; 394 Mat F_diag; 395 IS is_iden; 396 Vec b; 397 398 PetscFunctionBegin; 399 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 400 (*F)->ops->solve = MatSolve_MUMPS; 401 402 /* Initialize a MUMPS instance */ 403 ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 404 ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 405 lua->myid = lu->myid; lua->size = lu->size; 406 lu->id.job = JOB_INIT; 407 ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 408 ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr); 409 410 /* Set mumps options */ 411 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 412 lu->id.par=1; /* host participates factorizaton and solve */ 413 lu->id.sym=lu->sym; 414 if (lu->sym == 2){ 415 ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 416 if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 417 } 418 #if defined(PETSC_USE_COMPLEX) 419 zmumps_c(&lu->id); 420 #else 421 dmumps_c(&lu->id); 422 #endif 423 424 if (lu->size == 1){ 425 lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 426 } else { 427 lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 428 } 429 430 icntl=-1; 431 lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 432 ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 433 if ((flg && icntl > 0) || PetscLogPrintInfo) { 434 lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 435 } else { /* no output */ 436 lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 437 lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 438 lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 439 } 440 ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): matrix prescaling (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr); 441 icntl=-1; 442 ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 443 if (flg) { 444 if (icntl== 1){ 445 SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 446 } else { 447 lu->id.ICNTL(7) = icntl; 448 } 449 } 450 ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): A or A^T x=b to be solved. 1: A; otherwise: A^T","None",lu->id.ICNTL(9),&lu->id.ICNTL(9),PETSC_NULL);CHKERRQ(ierr); 451 ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",lu->id.ICNTL(10),&lu->id.ICNTL(10),PETSC_NULL);CHKERRQ(ierr); 452 ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr); 453 ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 454 ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 455 ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr); 456 ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 457 458 ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr); 459 ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",lu->id.CNTL(2),&lu->id.CNTL(2),PETSC_NULL);CHKERRQ(ierr); 460 ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr); 461 ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr); 462 PetscOptionsEnd(); 463 } 464 465 /* define matrix A */ 466 switch (lu->id.ICNTL(18)){ 467 case 0: /* centralized assembled matrix input (size=1) */ 468 if (!lu->myid) { 469 if (lua->isAIJ){ 470 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 471 nz = aa->nz; 472 ai = aa->i; aj = aa->j; lu->val = aa->a; 473 } else { 474 Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 475 nz = aa->nz; 476 ai = aa->i; aj = aa->j; lu->val = aa->a; 477 } 478 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 479 ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 480 ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 481 nz = 0; 482 for (i=0; i<M; i++){ 483 rnz = ai[i+1] - ai[i]; 484 while (rnz--) { /* Fortran row/col index! */ 485 lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 486 } 487 } 488 } 489 } 490 break; 491 case 3: /* distributed assembled matrix input (size>1) */ 492 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 493 valOnly = PETSC_FALSE; 494 } else { 495 valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 496 } 497 ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 498 break; 499 default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 500 } 501 502 /* analysis phase */ 503 /*----------------*/ 504 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 505 lu->id.job = 1; 506 507 lu->id.n = M; 508 switch (lu->id.ICNTL(18)){ 509 case 0: /* centralized assembled matrix input */ 510 if (!lu->myid) { 511 lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 512 if (lu->id.ICNTL(6)>1){ 513 #if defined(PETSC_USE_COMPLEX) 514 lu->id.a = (mumps_double_complex*)lu->val; 515 #else 516 lu->id.a = lu->val; 517 #endif 518 } 519 } 520 break; 521 case 3: /* distributed assembled matrix input (size>1) */ 522 lu->id.nz_loc = nnz; 523 lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 524 if (lu->id.ICNTL(6)>1) { 525 #if defined(PETSC_USE_COMPLEX) 526 lu->id.a_loc = (mumps_double_complex*)lu->val; 527 #else 528 lu->id.a_loc = lu->val; 529 #endif 530 } 531 /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 532 if (!lu->myid){ 533 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr); 534 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr); 535 } else { 536 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 537 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 538 } 539 ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 540 ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 541 ierr = VecSetFromOptions(b);CHKERRQ(ierr); 542 543 ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 544 ierr = ISDestroy(is_iden);CHKERRQ(ierr); 545 ierr = VecDestroy(b);CHKERRQ(ierr); 546 break; 547 } 548 #if defined(PETSC_USE_COMPLEX) 549 zmumps_c(&lu->id); 550 #else 551 dmumps_c(&lu->id); 552 #endif 553 if (lu->id.INFOG(1) < 0) { 554 SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 555 } 556 } 557 558 /* numerical factorization phase */ 559 /*-------------------------------*/ 560 lu->id.job = 2; 561 if(!lu->id.ICNTL(18)) { 562 if (!lu->myid) { 563 #if defined(PETSC_USE_COMPLEX) 564 lu->id.a = (mumps_double_complex*)lu->val; 565 #else 566 lu->id.a = lu->val; 567 #endif 568 } 569 } else { 570 #if defined(PETSC_USE_COMPLEX) 571 lu->id.a_loc = (mumps_double_complex*)lu->val; 572 #else 573 lu->id.a_loc = lu->val; 574 #endif 575 } 576 #if defined(PETSC_USE_COMPLEX) 577 zmumps_c(&lu->id); 578 #else 579 dmumps_c(&lu->id); 580 #endif 581 if (lu->id.INFOG(1) < 0) { 582 if (lu->id.INFO(1) == -13) { 583 SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 584 } else { 585 SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2)); 586 } 587 } 588 589 if (!lu->myid && lu->id.ICNTL(16) > 0){ 590 SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 591 } 592 593 if (lu->size > 1){ 594 if ((*F)->factor == FACTOR_LU){ 595 F_diag = ((Mat_MPIAIJ *)(*F)->data)->A; 596 } else { 597 F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A; 598 } 599 F_diag->assembled = PETSC_TRUE; 600 if (lu->nSolve){ 601 ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 602 ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 603 ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 604 } 605 } 606 (*F)->assembled = PETSC_TRUE; 607 lu->matstruc = SAME_NONZERO_PATTERN; 608 lu->CleanUpMUMPS = PETSC_TRUE; 609 lu->nSolve = 0; 610 PetscFunctionReturn(0); 611 } 612 613 /* Note the Petsc r and c permutations are ignored */ 614 #undef __FUNCT__ 615 #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 616 PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 617 Mat B; 618 Mat_MUMPS *lu; 619 PetscErrorCode ierr; 620 621 PetscFunctionBegin; 622 /* Create the factorization matrix */ 623 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 624 ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 625 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 626 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 627 ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 628 629 B->ops->lufactornumeric = MatFactorNumeric_MUMPS; 630 B->factor = FACTOR_LU; 631 lu = (Mat_MUMPS*)B->spptr; 632 lu->sym = 0; 633 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 634 635 *F = B; 636 PetscFunctionReturn(0); 637 } 638 639 /* Note the Petsc r permutation is ignored */ 640 #undef __FUNCT__ 641 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 642 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) { 643 Mat B; 644 Mat_MUMPS *lu; 645 PetscErrorCode ierr; 646 647 PetscFunctionBegin; 648 /* Create the factorization matrix */ 649 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 650 ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 651 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 652 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 653 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 654 655 B->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 656 B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 657 B->factor = FACTOR_CHOLESKY; 658 lu = (Mat_MUMPS*)B->spptr; 659 lu->sym = 2; 660 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 661 662 *F = B; 663 PetscFunctionReturn(0); 664 } 665 666 #undef __FUNCT__ 667 #define __FUNCT__ "MatFactorInfo_MUMPS" 668 PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 669 Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 670 PetscErrorCode ierr; 671 672 PetscFunctionBegin; 673 /* check if matrix is mumps type */ 674 if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 675 676 ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 677 ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 678 ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 679 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 680 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 681 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 682 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 683 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 684 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 685 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 686 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 687 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 688 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 689 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 690 if (!lu->myid && lu->id.ICNTL(11)>0) { 691 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 692 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 693 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 694 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",lu->id.RINFOG(7),lu->id.RINFOG(8));CHKERRQ(ierr); 695 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 696 ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 697 698 } 699 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 700 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 701 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 702 /* ICNTL(15-17) not used */ 703 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 704 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 705 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 706 ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 707 708 ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 709 ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 710 ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 711 ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 712 713 /* infomation local to each processor */ 714 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);} 715 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 716 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 717 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);} 718 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 719 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 720 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);} 721 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 722 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 723 /* 724 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);} 725 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr); 726 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 727 */ 728 729 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);} 730 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 731 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 732 733 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);} 734 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 735 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 736 737 if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);} 738 ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 739 ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 740 741 if (!lu->myid){ /* information from the host */ 742 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 743 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 744 ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 745 746 ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 747 ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 748 ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 749 ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 750 ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 751 ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 752 ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 753 ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 754 ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 755 ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 756 ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 757 ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 758 ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 759 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",lu->id.INFOG(16));CHKERRQ(ierr); 760 ierr = PetscViewerASCIIPrintf(viewer," INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));CHKERRQ(ierr); 761 ierr = PetscViewerASCIIPrintf(viewer," INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));CHKERRQ(ierr); 762 ierr = PetscViewerASCIIPrintf(viewer," INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));CHKERRQ(ierr); 763 ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 764 ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr); 765 ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr); 766 ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 767 ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 768 ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 769 } 770 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "MatView_MUMPS" 776 PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) { 777 PetscErrorCode ierr; 778 PetscTruth iascii; 779 PetscViewerFormat format; 780 Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 781 782 PetscFunctionBegin; 783 ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr); 784 785 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 786 if (iascii) { 787 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 788 if (format == PETSC_VIEWER_ASCII_INFO){ 789 ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 790 } 791 } 792 PetscFunctionReturn(0); 793 } 794 795 #undef __FUNCT__ 796 #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS" 797 PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) { 798 PetscErrorCode ierr; 799 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 800 801 PetscFunctionBegin; 802 ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 803 804 mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 805 mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 806 A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 807 PetscFunctionReturn(0); 808 } 809 810 EXTERN_C_BEGIN 811 #undef __FUNCT__ 812 #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS" 813 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 814 { 815 PetscErrorCode ierr; 816 PetscMPIInt size; 817 MPI_Comm comm; 818 Mat B=*newmat; 819 Mat_MUMPS *mumps; 820 821 PetscFunctionBegin; 822 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 823 ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 824 825 if (reuse == MAT_INITIAL_MATRIX) { 826 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 827 /* A may have special container that is not duplicated, 828 e.g., A is obtainted from MatMatMult(,&A). Save B->ops instead */ 829 mumps->MatDuplicate = B->ops->duplicate; 830 mumps->MatView = B->ops->view; 831 mumps->MatAssemblyEnd = B->ops->assemblyend; 832 mumps->MatLUFactorSymbolic = B->ops->lufactorsymbolic; 833 mumps->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic; 834 mumps->MatDestroy = B->ops->destroy; 835 } else { 836 mumps->MatDuplicate = A->ops->duplicate; 837 mumps->MatView = A->ops->view; 838 mumps->MatAssemblyEnd = A->ops->assemblyend; 839 mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 840 mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 841 mumps->MatDestroy = A->ops->destroy; 842 } 843 mumps->specialdestroy = MatDestroy_AIJMUMPS; 844 mumps->CleanUpMUMPS = PETSC_FALSE; 845 mumps->isAIJ = PETSC_TRUE; 846 847 B->spptr = (void*)mumps; 848 B->ops->duplicate = MatDuplicate_MUMPS; 849 B->ops->view = MatView_MUMPS; 850 B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS; 851 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 852 B->ops->destroy = MatDestroy_MUMPS; 853 854 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 855 if (size == 1) { 856 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C", 857 "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 858 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C", 859 "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 860 } else { 861 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C", 862 "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 863 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C", 864 "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 865 } 866 867 ierr = PetscInfo(A,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr); 868 ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 869 *newmat = B; 870 PetscFunctionReturn(0); 871 } 872 EXTERN_C_END 873 874 /*MC 875 MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 876 and sequential matrices via the external package MUMPS. 877 878 If MUMPS is installed (see the manual for instructions 879 on how to declare the existence of external packages), 880 a matrix type can be constructed which invokes MUMPS solvers. 881 After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS), then 882 optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT 883 call MatCreateSeqAIJ/MPIAIJ() directly or the preallocation information will be LOST! 884 885 If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 886 Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 887 MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported 888 for communicators controlling multiple processes. It is recommended that you call both of 889 the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 890 conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 891 without data copy AFTER the matrix values are set. 892 893 Options Database Keys: 894 + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 895 . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 896 . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 897 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 898 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 899 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 900 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 901 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 902 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 903 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 904 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 905 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 906 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 907 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 908 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 909 910 Level: beginner 911 912 .seealso: MATSBAIJMUMPS 913 M*/ 914 915 EXTERN_C_BEGIN 916 #undef __FUNCT__ 917 #define __FUNCT__ "MatCreate_AIJMUMPS" 918 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A) 919 { 920 PetscErrorCode ierr; 921 PetscMPIInt size; 922 923 PetscFunctionBegin; 924 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 925 ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 926 ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 927 PetscFunctionReturn(0); 928 } 929 EXTERN_C_END 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS" 933 PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) 934 { 935 PetscErrorCode ierr; 936 Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 937 938 PetscFunctionBegin; 939 ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 940 mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 941 A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 942 PetscFunctionReturn(0); 943 } 944 945 EXTERN_C_BEGIN 946 #undef __FUNCT__ 947 #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS" 948 PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz) 949 { 950 Mat A; 951 Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr; 952 PetscErrorCode ierr; 953 954 PetscFunctionBegin; 955 /* 956 After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix 957 into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS. I would 958 like this to be done in the MatCreate routine, but the creation of this inner matrix requires 959 block size info so that PETSc can determine the local size properly. The block size info is set 960 in the preallocation routine. 961 */ 962 ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz); 963 A = ((Mat_MPISBAIJ *)B->data)->A; 964 ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 965 PetscFunctionReturn(0); 966 } 967 EXTERN_C_END 968 969 EXTERN_C_BEGIN 970 #undef __FUNCT__ 971 #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS" 972 PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 973 { 974 PetscErrorCode ierr; 975 PetscMPIInt size; 976 MPI_Comm comm; 977 Mat B=*newmat; 978 Mat_MUMPS *mumps; 979 void (*f)(void); 980 981 PetscFunctionBegin; 982 if (reuse == MAT_INITIAL_MATRIX) { 983 ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 984 } 985 986 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 987 ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 988 989 mumps->MatDuplicate = A->ops->duplicate; 990 mumps->MatView = A->ops->view; 991 mumps->MatAssemblyEnd = A->ops->assemblyend; 992 mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 993 mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 994 mumps->MatDestroy = A->ops->destroy; 995 mumps->specialdestroy = MatDestroy_SBAIJMUMPS; 996 mumps->CleanUpMUMPS = PETSC_FALSE; 997 mumps->isAIJ = PETSC_FALSE; 998 999 B->spptr = (void*)mumps; 1000 B->ops->duplicate = MatDuplicate_MUMPS; 1001 B->ops->view = MatView_MUMPS; 1002 B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS; 1003 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 1004 B->ops->destroy = MatDestroy_MUMPS; 1005 1006 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 1007 if (size == 1) { 1008 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C", 1009 "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 1010 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C", 1011 "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 1012 } else { 1013 /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */ 1014 ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 1015 if (f) { /* This case should always be true when this routine is called */ 1016 mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f; 1017 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C", 1018 "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS", 1019 MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr); 1020 } 1021 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C", 1022 "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 1023 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C", 1024 "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 1025 } 1026 1027 ierr = PetscInfo(A,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr); 1028 ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 1029 *newmat = B; 1030 PetscFunctionReturn(0); 1031 } 1032 EXTERN_C_END 1033 1034 #undef __FUNCT__ 1035 #define __FUNCT__ "MatDuplicate_MUMPS" 1036 PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) { 1037 PetscErrorCode ierr; 1038 Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr; 1039 1040 PetscFunctionBegin; 1041 ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 1042 ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr); 1043 PetscFunctionReturn(0); 1044 } 1045 1046 /*MC 1047 MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 1048 distributed and sequential matrices via the external package MUMPS. 1049 1050 If MUMPS is installed (see the manual for instructions 1051 on how to declare the existence of external packages), 1052 a matrix type can be constructed which invokes MUMPS solvers. 1053 After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS), then 1054 optionally call MatSeqSBAIJSetPreallocation() or MatMPISBAIJSetPreallocation() DO NOT 1055 call MatCreateSeqSBAIJ/MPISBAIJ() directly or the preallocation information will be LOST! 1056 1057 If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 1058 Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 1059 MatSeqSBAIJSetPreallocation() is supported, and similarly MatMPISBAIJSetPreallocation() is supported 1060 for communicators controlling multiple processes. It is recommended that you call both of 1061 the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 1062 conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 1063 without data copy AFTER the matrix values have been set. 1064 1065 Options Database Keys: 1066 + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 1067 . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 1068 . -mat_mumps_icntl_4 <0,...,4> - print level 1069 . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 1070 . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 1071 . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 1072 . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 1073 . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 1074 . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 1075 . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 1076 . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 1077 . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 1078 . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 1079 . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 1080 - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 1081 1082 Level: beginner 1083 1084 .seealso: MATAIJMUMPS 1085 M*/ 1086 1087 EXTERN_C_BEGIN 1088 #undef __FUNCT__ 1089 #define __FUNCT__ "MatCreate_SBAIJMUMPS" 1090 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A) 1091 { 1092 PetscErrorCode ierr; 1093 PetscMPIInt size; 1094 1095 PetscFunctionBegin; 1096 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 1097 ierr = MatSetType(A,MATSBAIJ);CHKERRQ(ierr); 1098 ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 1099 PetscFunctionReturn(0); 1100 } 1101 EXTERN_C_END 1102