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