1*be1d678aSKris Buschelman #define PETSCMAT_DLL 21c2a3de1SBarry Smith 3397b6df1SKris Buschelman /* 4a58c3f20SHong Zhang Provides an interface to the MUMPS_4.3.1 sparse solver 5397b6df1SKris Buschelman */ 6397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 7397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h" 8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h" 9397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h" 10397b6df1SKris Buschelman 11397b6df1SKris Buschelman EXTERN_C_BEGIN 12397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 13397b6df1SKris Buschelman #include "zmumps_c.h" 14397b6df1SKris Buschelman #else 15397b6df1SKris Buschelman #include "dmumps_c.h" 16397b6df1SKris Buschelman #endif 17397b6df1SKris Buschelman EXTERN_C_END 18397b6df1SKris Buschelman #define JOB_INIT -1 19397b6df1SKris Buschelman #define JOB_END -2 20397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 21397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 22397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 23397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 24a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 25397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 26adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 27397b6df1SKris Buschelman 28397b6df1SKris Buschelman typedef struct { 29397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 30397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 31397b6df1SKris Buschelman #else 32397b6df1SKris Buschelman DMUMPS_STRUC_C id; 33397b6df1SKris Buschelman #endif 34397b6df1SKris Buschelman MatStructure matstruc; 35397b6df1SKris Buschelman int myid,size,*irn,*jcn,sym; 36397b6df1SKris Buschelman PetscScalar *val; 37397b6df1SKris Buschelman MPI_Comm comm_mumps; 38397b6df1SKris Buschelman 39c338a77dSKris Buschelman PetscTruth isAIJ,CleanUpMUMPS; 406849ba73SBarry Smith PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 416849ba73SBarry Smith PetscErrorCode (*MatView)(Mat,PetscViewer); 426849ba73SBarry Smith PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 436849ba73SBarry Smith PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 446849ba73SBarry Smith PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*); 456849ba73SBarry Smith PetscErrorCode (*MatDestroy)(Mat); 466849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 476849ba73SBarry Smith PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*); 48f0c56d0fSKris Buschelman } Mat_MUMPS; 49f0c56d0fSKris Buschelman 50dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 51892f6c3fSKris Buschelman EXTERN_C_BEGIN 52*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,MatReuse,Mat*); 53892f6c3fSKris Buschelman EXTERN_C_END 54397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 55397b6df1SKris Buschelman /* 56397b6df1SKris Buschelman input: 5775747be1SHong Zhang A - matrix in mpiaij or mpisbaij (bs=1) format 58397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 59397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 60397b6df1SKris Buschelman TRUE: only the values in v array are updated 61397b6df1SKris Buschelman output: 62397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 63397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 64397b6df1SKris Buschelman */ 65dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) { 66397b6df1SKris Buschelman int *ai, *aj, *bi, *bj, rstart,nz, *garray; 67dfbe8321SBarry Smith PetscErrorCode ierr; 68dfbe8321SBarry Smith int i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol; 69d54de34fSKris Buschelman int *row,*col; 70397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 71f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 72397b6df1SKris Buschelman 73397b6df1SKris Buschelman PetscFunctionBegin; 74397b6df1SKris Buschelman if (mumps->isAIJ){ 75397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 76397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 77397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 78397b6df1SKris Buschelman nz = aa->nz + bb->nz; 79397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 80397b6df1SKris Buschelman garray = mat->garray; 81397b6df1SKris Buschelman av=aa->a; bv=bb->a; 82397b6df1SKris Buschelman 83397b6df1SKris Buschelman } else { 84397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 85397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 86397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 87521d7252SBarry Smith if (A->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->bs); 886c6c5352SBarry Smith nz = aa->nz + bb->nz; 89397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 90397b6df1SKris Buschelman garray = mat->garray; 91397b6df1SKris Buschelman av=aa->a; bv=bb->a; 92397b6df1SKris Buschelman } 93397b6df1SKris Buschelman 94397b6df1SKris Buschelman if (!valOnly){ 95397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr); 96397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr); 97397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 98397b6df1SKris Buschelman *r = row; *c = col; *v = val; 99397b6df1SKris Buschelman } else { 100397b6df1SKris Buschelman row = *r; col = *c; val = *v; 101397b6df1SKris Buschelman } 102397b6df1SKris Buschelman *nnz = nz; 103397b6df1SKris Buschelman 104028e57e8SHong Zhang jj = 0; irow = rstart; 105397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 106397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 107397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 108397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 109397b6df1SKris Buschelman bjj = bj + bi[i]; 110397b6df1SKris Buschelman 111397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 112397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 11375747be1SHong Zhang j=-1; 11475747be1SHong Zhang do { 11575747be1SHong Zhang j++; 11675747be1SHong Zhang if (j == countB) break; 117397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11875747be1SHong Zhang } while (jcol < colA_start); 11975747be1SHong Zhang jB = j; 120397b6df1SKris Buschelman 121397b6df1SKris Buschelman /* B-part, smaller col index */ 122397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 123397b6df1SKris Buschelman for (j=0; j<jB; j++){ 124397b6df1SKris Buschelman jcol = garray[bjj[j]]; 125397b6df1SKris Buschelman if (!valOnly){ 126397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 12775747be1SHong Zhang 128397b6df1SKris Buschelman } 129397b6df1SKris Buschelman val[jj++] = *bv++; 130397b6df1SKris Buschelman } 131397b6df1SKris Buschelman /* A-part */ 132397b6df1SKris Buschelman for (j=0; j<countA; j++){ 133397b6df1SKris Buschelman if (!valOnly){ 134397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 135397b6df1SKris Buschelman } 136397b6df1SKris Buschelman val[jj++] = *av++; 137397b6df1SKris Buschelman } 138397b6df1SKris Buschelman /* B-part, larger col index */ 139397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 140397b6df1SKris Buschelman if (!valOnly){ 141397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 142397b6df1SKris Buschelman } 143397b6df1SKris Buschelman val[jj++] = *bv++; 144397b6df1SKris Buschelman } 145397b6df1SKris Buschelman irow++; 146397b6df1SKris Buschelman } 147397b6df1SKris Buschelman 148397b6df1SKris Buschelman PetscFunctionReturn(0); 149397b6df1SKris Buschelman } 150397b6df1SKris Buschelman 151c338a77dSKris Buschelman EXTERN_C_BEGIN 152c338a77dSKris Buschelman #undef __FUNCT__ 153c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base" 154*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,const MatType type,MatReuse reuse,Mat *newmat) \ 155521d7252SBarry Smith { 156dfbe8321SBarry Smith PetscErrorCode ierr; 157c338a77dSKris Buschelman Mat B=*newmat; 158f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 159901853e0SKris Buschelman void (*f)(void); 160c338a77dSKris Buschelman 161c338a77dSKris Buschelman PetscFunctionBegin; 162ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 163c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 164c338a77dSKris Buschelman } 165f0c56d0fSKris Buschelman B->ops->duplicate = mumps->MatDuplicate; 166f0c56d0fSKris Buschelman B->ops->view = mumps->MatView; 167f0c56d0fSKris Buschelman B->ops->assemblyend = mumps->MatAssemblyEnd; 168f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic; 169f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic; 170f0c56d0fSKris Buschelman B->ops->destroy = mumps->MatDestroy; 171901853e0SKris Buschelman 172901853e0SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr); 173901853e0SKris Buschelman if (f) { 174d1390fe9SSatish Balay ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(FCNVOID)mumps->MatPreallocate);CHKERRQ(ierr); 175901853e0SKris Buschelman } 176f0c56d0fSKris Buschelman ierr = PetscFree(mumps);CHKERRQ(ierr); 177901853e0SKris Buschelman 178901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 179901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 180901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 181901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr); 1822895b8caSSatish Balay ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 1832895b8caSSatish Balay ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr); 184901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 185901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr); 186901853e0SKris Buschelman 187901853e0SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr); 188c338a77dSKris Buschelman *newmat = B; 189c338a77dSKris Buschelman PetscFunctionReturn(0); 190c338a77dSKris Buschelman } 191c338a77dSKris Buschelman EXTERN_C_END 192c338a77dSKris Buschelman 193397b6df1SKris Buschelman #undef __FUNCT__ 1943924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 195dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 196dfbe8321SBarry Smith { 197f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 198dfbe8321SBarry Smith PetscErrorCode ierr; 199dfbe8321SBarry Smith int size=lu->size; 2006849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 201397b6df1SKris Buschelman PetscFunctionBegin; 202397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 203397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 204397b6df1SKris Buschelman lu->id.job=JOB_END; 205397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 206397b6df1SKris Buschelman zmumps_c(&lu->id); 207397b6df1SKris Buschelman #else 208397b6df1SKris Buschelman dmumps_c(&lu->id); 209397b6df1SKris Buschelman #endif 210c338a77dSKris Buschelman if (lu->irn) { 211c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 212c338a77dSKris Buschelman } 213c338a77dSKris Buschelman if (lu->jcn) { 214c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 215c338a77dSKris Buschelman } 216c338a77dSKris Buschelman if (size>1 && lu->val) { 217c338a77dSKris Buschelman ierr = PetscFree(lu->val);CHKERRQ(ierr); 218c338a77dSKris Buschelman } 219397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 220397b6df1SKris Buschelman } 221a39386dcSKris Buschelman specialdestroy = lu->specialdestroy; 222a39386dcSKris Buschelman ierr = (*specialdestroy)(A);CHKERRQ(ierr); 223c338a77dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 224397b6df1SKris Buschelman PetscFunctionReturn(0); 225397b6df1SKris Buschelman } 226397b6df1SKris Buschelman 227397b6df1SKris Buschelman #undef __FUNCT__ 228a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS" 229dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A) 230dfbe8321SBarry Smith { 2316849ba73SBarry Smith PetscErrorCode ierr; 2326849ba73SBarry Smith int size; 233a39386dcSKris Buschelman 234a39386dcSKris Buschelman PetscFunctionBegin; 235a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 236a39386dcSKris Buschelman if (size==1) { 237ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 238a39386dcSKris Buschelman } else { 239ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 240a39386dcSKris Buschelman } 241a39386dcSKris Buschelman PetscFunctionReturn(0); 242a39386dcSKris Buschelman } 243a39386dcSKris Buschelman 244a39386dcSKris Buschelman #undef __FUNCT__ 245a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS" 246dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A) 247dfbe8321SBarry Smith { 2486849ba73SBarry Smith PetscErrorCode ierr; 2496849ba73SBarry Smith int size; 250a39386dcSKris Buschelman 251a39386dcSKris Buschelman PetscFunctionBegin; 252a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 253a39386dcSKris Buschelman if (size==1) { 254ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 255a39386dcSKris Buschelman } else { 256ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 257a39386dcSKris Buschelman } 258a39386dcSKris Buschelman PetscFunctionReturn(0); 259a39386dcSKris Buschelman } 260a39386dcSKris Buschelman 261a39386dcSKris Buschelman #undef __FUNCT__ 262c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS" 263dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 264f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 265dfbe8321SBarry Smith PetscErrorCode ierr; 266397b6df1SKris Buschelman 267397b6df1SKris Buschelman PetscFunctionBegin; 268c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 269c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 270c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 271c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 272c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 273c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 274c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 275c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 276c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 277c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 278958c9bccSBarry Smith if (!lu->myid && lu->id.ICNTL(11)>0) { 279c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 280c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 281c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 282c338a77dSKris Buschelman 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); 283c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 284c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 285c338a77dSKris Buschelman 286c338a77dSKris Buschelman } 287c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 288c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 289adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 290c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(15) (efficiency control): %d \n",lu->id.ICNTL(15));CHKERRQ(ierr); 291c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 292c338a77dSKris Buschelman 293c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 294c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 295c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 29657f0c58bSHong Zhang 29757f0c58bSHong Zhang /* infomation local to each processor */ 298958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 29957f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 30057f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 301958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 30257f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 30357f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 304958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 30557f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 30657f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 307adc1d99fSHong Zhang 308958c9bccSBarry Smith if (!lu->myid){ /* information from the host */ 309adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 310adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 311adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 312adc1d99fSHong Zhang 313adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 314adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 315adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 316adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 317adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 318adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 319adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 320adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 321adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 322adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 323adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 324adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 325adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 326adc1d99fSHong Zhang 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); 327adc1d99fSHong Zhang 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); 328adc1d99fSHong Zhang 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); 329adc1d99fSHong Zhang 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); 330adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 331adc1d99fSHong Zhang } 332adc1d99fSHong Zhang 333397b6df1SKris Buschelman PetscFunctionReturn(0); 334397b6df1SKris Buschelman } 335397b6df1SKris Buschelman 336397b6df1SKris Buschelman #undef __FUNCT__ 337f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS" 338dfbe8321SBarry Smith PetscErrorCode MatView_AIJMUMPS(Mat A,PetscViewer viewer) { 339dfbe8321SBarry Smith PetscErrorCode ierr; 34032077d6dSBarry Smith PetscTruth iascii; 341397b6df1SKris Buschelman PetscViewerFormat format; 342f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 343397b6df1SKris Buschelman 344397b6df1SKris Buschelman PetscFunctionBegin; 345397b6df1SKris Buschelman ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr); 346397b6df1SKris Buschelman 34732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 34832077d6dSBarry Smith if (iascii) { 349397b6df1SKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 350397b6df1SKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 351397b6df1SKris Buschelman ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 352397b6df1SKris Buschelman } 353397b6df1SKris Buschelman } 354397b6df1SKris Buschelman PetscFunctionReturn(0); 355397b6df1SKris Buschelman } 356397b6df1SKris Buschelman 357397b6df1SKris Buschelman #undef __FUNCT__ 358f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS" 359dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) { 360f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 361d54de34fSKris Buschelman PetscScalar *array; 362397b6df1SKris Buschelman Vec x_seq; 363397b6df1SKris Buschelman IS iden; 364397b6df1SKris Buschelman VecScatter scat; 365dfbe8321SBarry Smith PetscErrorCode ierr; 366397b6df1SKris Buschelman 367397b6df1SKris Buschelman PetscFunctionBegin; 368397b6df1SKris Buschelman if (lu->size > 1){ 369397b6df1SKris Buschelman if (!lu->myid){ 370397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr); 371397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr); 372397b6df1SKris Buschelman } else { 373397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr); 374397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr); 375397b6df1SKris Buschelman } 376397b6df1SKris Buschelman ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr); 377397b6df1SKris Buschelman ierr = ISDestroy(iden);CHKERRQ(ierr); 378397b6df1SKris Buschelman 379397b6df1SKris Buschelman ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 380397b6df1SKris Buschelman ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 381397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 382397b6df1SKris Buschelman } else { /* size == 1 */ 383397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 384397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 385397b6df1SKris Buschelman } 386397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 387397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 388397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 389397b6df1SKris Buschelman #else 390397b6df1SKris Buschelman lu->id.rhs = array; 391397b6df1SKris Buschelman #endif 392397b6df1SKris Buschelman } 393397b6df1SKris Buschelman 394397b6df1SKris Buschelman /* solve phase */ 395397b6df1SKris Buschelman lu->id.job=3; 396397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 397397b6df1SKris Buschelman zmumps_c(&lu->id); 398397b6df1SKris Buschelman #else 399397b6df1SKris Buschelman dmumps_c(&lu->id); 400397b6df1SKris Buschelman #endif 401397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 40279a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 403397b6df1SKris Buschelman } 404397b6df1SKris Buschelman 405397b6df1SKris Buschelman /* convert mumps solution x_seq to petsc mpi x */ 406397b6df1SKris Buschelman if (lu->size > 1) { 407397b6df1SKris Buschelman if (!lu->myid){ 408397b6df1SKris Buschelman ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 409397b6df1SKris Buschelman } 410397b6df1SKris Buschelman ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 411397b6df1SKris Buschelman ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 412397b6df1SKris Buschelman ierr = VecScatterDestroy(scat);CHKERRQ(ierr); 413397b6df1SKris Buschelman ierr = VecDestroy(x_seq);CHKERRQ(ierr); 414397b6df1SKris Buschelman } else { 415397b6df1SKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 416397b6df1SKris Buschelman } 417397b6df1SKris Buschelman 418397b6df1SKris Buschelman PetscFunctionReturn(0); 419397b6df1SKris Buschelman } 420397b6df1SKris Buschelman 421a58c3f20SHong Zhang /* 422a58c3f20SHong Zhang input: 423a58c3f20SHong Zhang F: numeric factor 424a58c3f20SHong Zhang output: 425a58c3f20SHong Zhang nneg: total number of negative pivots 426a58c3f20SHong Zhang nzero: 0 427a58c3f20SHong Zhang npos: (global dimension of F) - nneg 428a58c3f20SHong Zhang */ 429a58c3f20SHong Zhang 430a58c3f20SHong Zhang #undef __FUNCT__ 431a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 432dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 433a58c3f20SHong Zhang { 434a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 435dfbe8321SBarry Smith PetscErrorCode ierr; 436dfbe8321SBarry Smith int size; 437a58c3f20SHong Zhang 438a58c3f20SHong Zhang PetscFunctionBegin; 439bcb30aebSHong Zhang ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr); 440bcb30aebSHong Zhang /* 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 */ 441bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 44279a5c55eSBarry Smith 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)); 443bcb30aebSHong Zhang } 444a58c3f20SHong Zhang if (nneg){ 445a58c3f20SHong Zhang if (!lu->myid){ 446a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 447a58c3f20SHong Zhang } 448bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 449a58c3f20SHong Zhang } 450a58c3f20SHong Zhang if (nzero) *nzero = 0; 451a58c3f20SHong Zhang if (npos) *npos = F->M - (*nneg); 452a58c3f20SHong Zhang PetscFunctionReturn(0); 453a58c3f20SHong Zhang } 454a58c3f20SHong Zhang 455397b6df1SKris Buschelman #undef __FUNCT__ 45619facb7aSBarry Smith #define __FUNCT__ "MatFactorNumeric_AIJMUMPS" 457af281ebdSHong Zhang PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F) 458af281ebdSHong Zhang { 459f0c56d0fSKris Buschelman Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 460f0c56d0fSKris Buschelman Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr; 4616849ba73SBarry Smith PetscErrorCode ierr; 462af281ebdSHong Zhang PetscInt rnz,nnz,nz,i,M=A->M,*ai,*aj,icntl; 463397b6df1SKris Buschelman PetscTruth valOnly,flg; 464397b6df1SKris Buschelman 465397b6df1SKris Buschelman PetscFunctionBegin; 466397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 467f0c56d0fSKris Buschelman (*F)->ops->solve = MatSolve_AIJMUMPS; 468397b6df1SKris Buschelman 469397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 470397b6df1SKris Buschelman ierr = MPI_Comm_rank(A->comm, &lu->myid); 471397b6df1SKris Buschelman ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr); 47275747be1SHong Zhang lua->myid = lu->myid; lua->size = lu->size; 473397b6df1SKris Buschelman lu->id.job = JOB_INIT; 474397b6df1SKris Buschelman ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr); 475397b6df1SKris Buschelman lu->id.comm_fortran = lu->comm_mumps; 476397b6df1SKris Buschelman 477397b6df1SKris Buschelman /* Set mumps options */ 478397b6df1SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 479397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 480397b6df1SKris Buschelman lu->id.sym=lu->sym; 481397b6df1SKris Buschelman if (lu->sym == 2){ 482397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 483397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 484397b6df1SKris Buschelman } 485397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 486397b6df1SKris Buschelman zmumps_c(&lu->id); 487397b6df1SKris Buschelman #else 488397b6df1SKris Buschelman dmumps_c(&lu->id); 489397b6df1SKris Buschelman #endif 490397b6df1SKris Buschelman 491397b6df1SKris Buschelman if (lu->size == 1){ 492397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 493397b6df1SKris Buschelman } else { 494397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 495397b6df1SKris Buschelman } 496397b6df1SKris Buschelman 497397b6df1SKris Buschelman icntl=-1; 498397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 49919facb7aSBarry Smith if ((flg && icntl > 0) || PetscLogPrintInfo) { 500397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 501397b6df1SKris Buschelman } else { /* no output */ 502397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 503397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 504397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 505397b6df1SKris Buschelman lu->id.ICNTL(4) = 0; /* level of printing, 0,1,2,3,4, default=2 */ 506397b6df1SKris Buschelman } 507397b6df1SKris Buschelman 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); 508397b6df1SKris Buschelman icntl=-1; 509397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 510397b6df1SKris Buschelman if (flg) { 511397b6df1SKris Buschelman if (icntl== 1){ 512397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 513397b6df1SKris Buschelman } else { 514397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 515397b6df1SKris Buschelman } 516397b6df1SKris Buschelman } 517397b6df1SKris Buschelman 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); 518397b6df1SKris Buschelman 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); 51994b7f48cSBarry Smith 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); 520397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 521397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 522adc1d99fSHong Zhang 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); 523397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 524397b6df1SKris Buschelman 525397b6df1SKris Buschelman /* 526397b6df1SKris Buschelman 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); 527397b6df1SKris Buschelman if (flg){ 528397b6df1SKris Buschelman if (icntl >-1 && icntl <3 ){ 529397b6df1SKris Buschelman if (lu->myid==0) lu->id.ICNTL(16) = icntl; 530397b6df1SKris Buschelman } else { 531397b6df1SKris Buschelman SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl); 532397b6df1SKris Buschelman } 533397b6df1SKris Buschelman } 534397b6df1SKris Buschelman */ 535397b6df1SKris Buschelman 536397b6df1SKris Buschelman ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr); 537397b6df1SKris Buschelman 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); 538397b6df1SKris Buschelman ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr); 539397b6df1SKris Buschelman PetscOptionsEnd(); 540397b6df1SKris Buschelman } 541397b6df1SKris Buschelman 542397b6df1SKris Buschelman /* define matrix A */ 543397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 544397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 545397b6df1SKris Buschelman if (!lu->myid) { 546c36ead0aSKris Buschelman if (lua->isAIJ){ 547397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 548397b6df1SKris Buschelman nz = aa->nz; 549397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 550397b6df1SKris Buschelman } else { 551397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 5526c6c5352SBarry Smith nz = aa->nz; 553397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 554397b6df1SKris Buschelman } 555397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 556397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr); 557397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr); 558397b6df1SKris Buschelman nz = 0; 559397b6df1SKris Buschelman for (i=0; i<M; i++){ 560397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 561397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 562397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 563397b6df1SKris Buschelman } 564397b6df1SKris Buschelman } 565397b6df1SKris Buschelman } 566397b6df1SKris Buschelman } 567397b6df1SKris Buschelman break; 568397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 569397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 570397b6df1SKris Buschelman valOnly = PETSC_FALSE; 571397b6df1SKris Buschelman } else { 572397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 573397b6df1SKris Buschelman } 574397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 575397b6df1SKris Buschelman break; 576397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 577397b6df1SKris Buschelman } 578397b6df1SKris Buschelman 579397b6df1SKris Buschelman /* analysis phase */ 580397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 581397b6df1SKris Buschelman lu->id.n = M; 582397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 583397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 584397b6df1SKris Buschelman if (!lu->myid) { 585397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 586397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 587397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 588397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 589397b6df1SKris Buschelman #else 590397b6df1SKris Buschelman lu->id.a = lu->val; 591397b6df1SKris Buschelman #endif 592397b6df1SKris Buschelman } 593397b6df1SKris Buschelman } 594397b6df1SKris Buschelman break; 595397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 596397b6df1SKris Buschelman lu->id.nz_loc = nnz; 597397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 598397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 599397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 600397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 601397b6df1SKris Buschelman #else 602397b6df1SKris Buschelman lu->id.a_loc = lu->val; 603397b6df1SKris Buschelman #endif 604397b6df1SKris Buschelman } 605397b6df1SKris Buschelman break; 606397b6df1SKris Buschelman } 607397b6df1SKris Buschelman lu->id.job=1; 608397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 609397b6df1SKris Buschelman zmumps_c(&lu->id); 610397b6df1SKris Buschelman #else 611397b6df1SKris Buschelman dmumps_c(&lu->id); 612397b6df1SKris Buschelman #endif 613397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 61479a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 615397b6df1SKris Buschelman } 616397b6df1SKris Buschelman } 617397b6df1SKris Buschelman 618397b6df1SKris Buschelman /* numerical factorization phase */ 619958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 620a7aca84bSHong Zhang if (!lu->myid) { 621397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 622397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 623397b6df1SKris Buschelman #else 624397b6df1SKris Buschelman lu->id.a = lu->val; 625397b6df1SKris Buschelman #endif 626397b6df1SKris Buschelman } 627397b6df1SKris Buschelman } else { 628397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 629397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 630397b6df1SKris Buschelman #else 631397b6df1SKris Buschelman lu->id.a_loc = lu->val; 632397b6df1SKris Buschelman #endif 633397b6df1SKris Buschelman } 634397b6df1SKris Buschelman lu->id.job=2; 635397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 636397b6df1SKris Buschelman zmumps_c(&lu->id); 637397b6df1SKris Buschelman #else 638397b6df1SKris Buschelman dmumps_c(&lu->id); 639397b6df1SKris Buschelman #endif 640397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 64119facb7aSBarry Smith if (lu->id.INFO(1) == -13) { 64219facb7aSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 64319facb7aSBarry Smith } else { 64479a5c55eSBarry Smith 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)); 645397b6df1SKris Buschelman } 64619facb7aSBarry Smith } 647397b6df1SKris Buschelman 64819facb7aSBarry Smith if (!lu->myid && lu->id.ICNTL(16) > 0){ 64979a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 650397b6df1SKris Buschelman } 651397b6df1SKris Buschelman 652397b6df1SKris Buschelman (*F)->assembled = PETSC_TRUE; 653397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 654ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 655397b6df1SKris Buschelman PetscFunctionReturn(0); 656397b6df1SKris Buschelman } 657397b6df1SKris Buschelman 658397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 659397b6df1SKris Buschelman #undef __FUNCT__ 660f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 661dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 662397b6df1SKris Buschelman Mat B; 663f0c56d0fSKris Buschelman Mat_MUMPS *lu; 664dfbe8321SBarry Smith PetscErrorCode ierr; 665397b6df1SKris Buschelman 666397b6df1SKris Buschelman PetscFunctionBegin; 667397b6df1SKris Buschelman 668397b6df1SKris Buschelman /* Create the factorization matrix */ 669397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 670be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 671397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 672397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 673397b6df1SKris Buschelman 674f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS; 675397b6df1SKris Buschelman B->factor = FACTOR_LU; 676f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 677397b6df1SKris Buschelman lu->sym = 0; 678397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 679397b6df1SKris Buschelman 680397b6df1SKris Buschelman *F = B; 681397b6df1SKris Buschelman PetscFunctionReturn(0); 682397b6df1SKris Buschelman } 683397b6df1SKris Buschelman 684397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 685397b6df1SKris Buschelman #undef __FUNCT__ 686f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 687dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) { 688397b6df1SKris Buschelman Mat B; 689f0c56d0fSKris Buschelman Mat_MUMPS *lu; 690dfbe8321SBarry Smith PetscErrorCode ierr; 691397b6df1SKris Buschelman 692397b6df1SKris Buschelman PetscFunctionBegin; 693397b6df1SKris Buschelman /* Create the factorization matrix */ 694397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 695be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 696efc670deSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 697efc670deSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 698397b6df1SKris Buschelman 699f0c56d0fSKris Buschelman B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS; 700a58c3f20SHong Zhang B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 701397b6df1SKris Buschelman B->factor = FACTOR_CHOLESKY; 702f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 703397b6df1SKris Buschelman lu->sym = 2; 704397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 705397b6df1SKris Buschelman 706397b6df1SKris Buschelman *F = B; 707397b6df1SKris Buschelman PetscFunctionReturn(0); 708397b6df1SKris Buschelman } 709397b6df1SKris Buschelman 710397b6df1SKris Buschelman #undef __FUNCT__ 711f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS" 712dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) { 713dfbe8321SBarry Smith PetscErrorCode ierr; 714f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 715c338a77dSKris Buschelman 716397b6df1SKris Buschelman PetscFunctionBegin; 717c338a77dSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 718f0c56d0fSKris Buschelman 719c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 720c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 721f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 722397b6df1SKris Buschelman PetscFunctionReturn(0); 723397b6df1SKris Buschelman } 724397b6df1SKris Buschelman 725c338a77dSKris Buschelman EXTERN_C_BEGIN 726c338a77dSKris Buschelman #undef __FUNCT__ 727f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS" 728*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) 729dfbe8321SBarry Smith { 730dfbe8321SBarry Smith PetscErrorCode ierr; 731521d7252SBarry Smith PetscMPIInt size; 732c338a77dSKris Buschelman MPI_Comm comm; 733c338a77dSKris Buschelman Mat B=*newmat; 734f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 735397b6df1SKris Buschelman 736397b6df1SKris Buschelman PetscFunctionBegin; 737ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 738c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 739397b6df1SKris Buschelman } 740397b6df1SKris Buschelman 741c338a77dSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 742f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 743c338a77dSKris Buschelman 744f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 745c338a77dSKris Buschelman mumps->MatView = A->ops->view; 746c338a77dSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 747c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 748c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 749c338a77dSKris Buschelman mumps->MatDestroy = A->ops->destroy; 750a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_AIJMUMPS; 751c338a77dSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 752f579278aSKris Buschelman mumps->isAIJ = PETSC_TRUE; 753c338a77dSKris Buschelman 7544b68dd72SKris Buschelman B->spptr = (void*)mumps; 755422e82a1SHong Zhang B->ops->duplicate = MatDuplicate_MUMPS; 756f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 757f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS; 758f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 7593924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 760c338a77dSKris Buschelman 761c338a77dSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 762c338a77dSKris Buschelman if (size == 1) { 763c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C", 764f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 765c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C", 766c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 767c338a77dSKris Buschelman } else { 768c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C", 769f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 770c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C", 771c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 772c338a77dSKris Buschelman } 773c338a77dSKris Buschelman 77463ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatConvert_AIJ_AIJMUMPS:Using MUMPS for LU factorization and solves.\n"));CHKERRQ(ierr); 775c338a77dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 776c338a77dSKris Buschelman *newmat = B; 777397b6df1SKris Buschelman PetscFunctionReturn(0); 778397b6df1SKris Buschelman } 779c338a77dSKris Buschelman EXTERN_C_END 780397b6df1SKris Buschelman 78124b6179bSKris Buschelman /*MC 782fafad747SKris Buschelman MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 78324b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 78424b6179bSKris Buschelman 78524b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 78624b6179bSKris Buschelman on how to declare the existence of external packages), 78724b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 78824b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS). 78924b6179bSKris Buschelman 79024b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 79124b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 79224b6179bSKris Buschelman MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported 79324b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 79428b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 79528b08bd3SKris Buschelman conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 79628b08bd3SKris Buschelman without data copy. 79724b6179bSKris Buschelman 79824b6179bSKris Buschelman Options Database Keys: 7990bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 80024b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 80124b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 80224b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 80324b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 80424b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 80524b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 80694b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 80724b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 80824b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 80924b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 81024b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 81124b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 81224b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 81324b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 81424b6179bSKris Buschelman 81524b6179bSKris Buschelman Level: beginner 81624b6179bSKris Buschelman 81724b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 81824b6179bSKris Buschelman M*/ 81924b6179bSKris Buschelman 820397b6df1SKris Buschelman EXTERN_C_BEGIN 821397b6df1SKris Buschelman #undef __FUNCT__ 822f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS" 823*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A) 824dfbe8321SBarry Smith { 825dfbe8321SBarry Smith PetscErrorCode ierr; 826dfbe8321SBarry Smith int size; 827e2d9671bSKris Buschelman Mat A_diag; 828397b6df1SKris Buschelman MPI_Comm comm; 829397b6df1SKris Buschelman 830397b6df1SKris Buschelman PetscFunctionBegin; 8315441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */ 8325441df8eSKris Buschelman /* and AIJMUMPS types */ 8335441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr); 834397b6df1SKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 835397b6df1SKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 836397b6df1SKris Buschelman if (size == 1) { 837397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 838397b6df1SKris Buschelman } else { 839397b6df1SKris Buschelman ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 840e2d9671bSKris Buschelman A_diag = ((Mat_MPIAIJ *)A->data)->A; 841ceb03754SKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr); 842397b6df1SKris Buschelman } 843ceb03754SKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 844397b6df1SKris Buschelman PetscFunctionReturn(0); 845397b6df1SKris Buschelman } 846397b6df1SKris Buschelman EXTERN_C_END 847397b6df1SKris Buschelman 848f579278aSKris Buschelman #undef __FUNCT__ 849f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS" 850dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) 851dfbe8321SBarry Smith { 852dfbe8321SBarry Smith PetscErrorCode ierr; 853f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 854f579278aSKris Buschelman 855f579278aSKris Buschelman PetscFunctionBegin; 856f579278aSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 857f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 858f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 859f0c56d0fSKris Buschelman A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 860f579278aSKris Buschelman PetscFunctionReturn(0); 861f579278aSKris Buschelman } 862f579278aSKris Buschelman 863f579278aSKris Buschelman EXTERN_C_BEGIN 864f579278aSKris Buschelman #undef __FUNCT__ 8659c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS" 866*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz) 8679c097c71SKris Buschelman { 8689c097c71SKris Buschelman Mat A; 86902217bfdSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr; 870dfbe8321SBarry Smith PetscErrorCode ierr; 8719c097c71SKris Buschelman 8729c097c71SKris Buschelman PetscFunctionBegin; 8739c097c71SKris Buschelman /* 8749c097c71SKris Buschelman After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix 8759c097c71SKris Buschelman into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS. I would 8769c097c71SKris Buschelman like this to be done in the MatCreate routine, but the creation of this inner matrix requires 8779c097c71SKris Buschelman block size info so that PETSc can determine the local size properly. The block size info is set 8789c097c71SKris Buschelman in the preallocation routine. 8799c097c71SKris Buschelman */ 8809c097c71SKris Buschelman ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz); 8819c097c71SKris Buschelman A = ((Mat_MPISBAIJ *)B->data)->A; 882ceb03754SKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 8839c097c71SKris Buschelman PetscFunctionReturn(0); 8849c097c71SKris Buschelman } 8859c097c71SKris Buschelman EXTERN_C_END 8869c097c71SKris Buschelman 8879c097c71SKris Buschelman EXTERN_C_BEGIN 8889c097c71SKris Buschelman #undef __FUNCT__ 889f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS" 890*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat) 891dfbe8321SBarry Smith { 892dfbe8321SBarry Smith PetscErrorCode ierr; 893521d7252SBarry Smith PetscMPIInt size; 894f579278aSKris Buschelman MPI_Comm comm; 895f579278aSKris Buschelman Mat B=*newmat; 896422e82a1SHong Zhang Mat_MUMPS *mumps; 8979c097c71SKris Buschelman void (*f)(void); 898f579278aSKris Buschelman 899f579278aSKris Buschelman PetscFunctionBegin; 900ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 901f579278aSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 902f579278aSKris Buschelman } 903f579278aSKris Buschelman 904f579278aSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 905f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 906f579278aSKris Buschelman 907f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 908f579278aSKris Buschelman mumps->MatView = A->ops->view; 909f579278aSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 910f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 911f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 912f579278aSKris Buschelman mumps->MatDestroy = A->ops->destroy; 913a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_SBAIJMUMPS; 914f579278aSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 915f579278aSKris Buschelman mumps->isAIJ = PETSC_FALSE; 916f579278aSKris Buschelman 917f579278aSKris Buschelman B->spptr = (void*)mumps; 918422e82a1SHong Zhang B->ops->duplicate = MatDuplicate_MUMPS; 919f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 920f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS; 921f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 9223924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 923f579278aSKris Buschelman 924f579278aSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 925f579278aSKris Buschelman if (size == 1) { 926f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C", 927f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 928f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C", 929f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 930f579278aSKris Buschelman } else { 9319c097c71SKris Buschelman /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */ 9329c097c71SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr); 933901853e0SKris Buschelman if (f) { /* This case should always be true when this routine is called */ 9346849ba73SBarry Smith mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f; 9359c097c71SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C", 9369c097c71SKris Buschelman "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS", 9379c097c71SKris Buschelman MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr); 9389c097c71SKris Buschelman } 939f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C", 940f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 941f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C", 942f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 943f579278aSKris Buschelman } 944f579278aSKris Buschelman 94563ba0a88SBarry Smith ierr = PetscLogInfo((0,"MatConvert_AIJ_AIJMUMPS:Using MUMPS for Cholesky factorization and solves.\n"));CHKERRQ(ierr); 946f579278aSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 947f579278aSKris Buschelman *newmat = B; 948f579278aSKris Buschelman PetscFunctionReturn(0); 949f579278aSKris Buschelman } 950f579278aSKris Buschelman EXTERN_C_END 951f579278aSKris Buschelman 952f0c56d0fSKris Buschelman #undef __FUNCT__ 953422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS" 954dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) { 955dfbe8321SBarry Smith PetscErrorCode ierr; 9568e393735SKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr; 9578f340917SKris Buschelman 958f0c56d0fSKris Buschelman PetscFunctionBegin; 9598f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 9608e393735SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr); 961f0c56d0fSKris Buschelman PetscFunctionReturn(0); 962f0c56d0fSKris Buschelman } 963f0c56d0fSKris Buschelman 96424b6179bSKris Buschelman /*MC 965fafad747SKris Buschelman MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 96624b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 96724b6179bSKris Buschelman 96824b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 96924b6179bSKris Buschelman on how to declare the existence of external packages), 97024b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 97124b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS). 97224b6179bSKris Buschelman 97324b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 97424b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 97524b6179bSKris Buschelman MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported 97624b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 97728b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 97828b08bd3SKris Buschelman conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 97928b08bd3SKris Buschelman without data copy. 98024b6179bSKris Buschelman 98124b6179bSKris Buschelman Options Database Keys: 9820bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 98324b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 98424b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 98524b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 98624b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 98724b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 98824b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 98994b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 99024b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 99124b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 99224b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 99324b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 99424b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 99524b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 99624b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 99724b6179bSKris Buschelman 99824b6179bSKris Buschelman Level: beginner 99924b6179bSKris Buschelman 100024b6179bSKris Buschelman .seealso: MATAIJMUMPS 100124b6179bSKris Buschelman M*/ 100224b6179bSKris Buschelman 1003397b6df1SKris Buschelman EXTERN_C_BEGIN 1004397b6df1SKris Buschelman #undef __FUNCT__ 1005f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS" 1006*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A) 1007dfbe8321SBarry Smith { 10086849ba73SBarry Smith PetscErrorCode ierr; 10096849ba73SBarry Smith int size; 1010397b6df1SKris Buschelman 1011397b6df1SKris Buschelman PetscFunctionBegin; 10125441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */ 10135441df8eSKris Buschelman /* and SBAIJMUMPS types */ 10145441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr); 10155441df8eSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 1016397b6df1SKris Buschelman if (size == 1) { 1017397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 1018397b6df1SKris Buschelman } else { 1019397b6df1SKris Buschelman ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr); 1020397b6df1SKris Buschelman } 1021ceb03754SKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 1022397b6df1SKris Buschelman PetscFunctionReturn(0); 1023397b6df1SKris Buschelman } 1024397b6df1SKris Buschelman EXTERN_C_END 1025