1be1d678aSKris Buschelman #define PETSCMAT_DLL 21c2a3de1SBarry Smith 3397b6df1SKris Buschelman /* 4c2b5dc30SHong Zhang Provides an interface to the MUMPS 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; 35c1490034SHong Zhang PetscMPIInt myid,size; 36c1490034SHong Zhang int *irn,*jcn,sym; 37397b6df1SKris Buschelman PetscScalar *val; 38397b6df1SKris Buschelman MPI_Comm comm_mumps; 39397b6df1SKris Buschelman 40c338a77dSKris Buschelman PetscTruth isAIJ,CleanUpMUMPS; 416849ba73SBarry Smith PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 426849ba73SBarry Smith PetscErrorCode (*MatView)(Mat,PetscViewer); 436849ba73SBarry Smith PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 446849ba73SBarry Smith PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 456849ba73SBarry Smith PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*); 466849ba73SBarry Smith PetscErrorCode (*MatDestroy)(Mat); 476849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 486849ba73SBarry Smith PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*); 49f0c56d0fSKris Buschelman } Mat_MUMPS; 50f0c56d0fSKris Buschelman 51dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 52892f6c3fSKris Buschelman EXTERN_C_BEGIN 5375179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*); 54892f6c3fSKris Buschelman EXTERN_C_END 55397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 56397b6df1SKris Buschelman /* 57397b6df1SKris Buschelman input: 5875747be1SHong Zhang A - matrix in mpiaij or mpisbaij (bs=1) format 59397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 60397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 61397b6df1SKris Buschelman TRUE: only the values in v array are updated 62397b6df1SKris Buschelman output: 63397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 64397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 65397b6df1SKris Buschelman */ 66dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) { 67c1490034SHong Zhang PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray; 68dfbe8321SBarry Smith PetscErrorCode ierr; 692a4c71feSBarry Smith PetscInt i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol; 70c1490034SHong Zhang PetscInt *row,*col; 71397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 72f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 73397b6df1SKris Buschelman 74397b6df1SKris Buschelman PetscFunctionBegin; 75397b6df1SKris Buschelman if (mumps->isAIJ){ 76397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 77397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 78397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 79397b6df1SKris Buschelman nz = aa->nz + bb->nz; 802a4c71feSBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart; 81397b6df1SKris Buschelman garray = mat->garray; 82397b6df1SKris Buschelman av=aa->a; bv=bb->a; 83397b6df1SKris Buschelman 84397b6df1SKris Buschelman } else { 85397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 86397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 87397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 880c0e133fSBarry Smith if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs); 896c6c5352SBarry Smith nz = aa->nz + bb->nz; 902a4c71feSBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart; 91397b6df1SKris Buschelman garray = mat->garray; 92397b6df1SKris Buschelman av=aa->a; bv=bb->a; 93397b6df1SKris Buschelman } 94397b6df1SKris Buschelman 95397b6df1SKris Buschelman if (!valOnly){ 967c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr); 977c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr); 98397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 99397b6df1SKris Buschelman *r = row; *c = col; *v = val; 100397b6df1SKris Buschelman } else { 101397b6df1SKris Buschelman row = *r; col = *c; val = *v; 102397b6df1SKris Buschelman } 103397b6df1SKris Buschelman *nnz = nz; 104397b6df1SKris Buschelman 105028e57e8SHong Zhang jj = 0; irow = rstart; 106397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 107397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 108397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 109397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 110397b6df1SKris Buschelman bjj = bj + bi[i]; 111397b6df1SKris Buschelman 112397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 113397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 11475747be1SHong Zhang j=-1; 11575747be1SHong Zhang do { 11675747be1SHong Zhang j++; 11775747be1SHong Zhang if (j == countB) break; 118397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11975747be1SHong Zhang } while (jcol < colA_start); 12075747be1SHong Zhang jB = j; 121397b6df1SKris Buschelman 122397b6df1SKris Buschelman /* B-part, smaller col index */ 123397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 124397b6df1SKris Buschelman for (j=0; j<jB; j++){ 125397b6df1SKris Buschelman jcol = garray[bjj[j]]; 126397b6df1SKris Buschelman if (!valOnly){ 127397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 12875747be1SHong Zhang 129397b6df1SKris Buschelman } 130397b6df1SKris Buschelman val[jj++] = *bv++; 131397b6df1SKris Buschelman } 132397b6df1SKris Buschelman /* A-part */ 133397b6df1SKris Buschelman for (j=0; j<countA; j++){ 134397b6df1SKris Buschelman if (!valOnly){ 135397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 136397b6df1SKris Buschelman } 137397b6df1SKris Buschelman val[jj++] = *av++; 138397b6df1SKris Buschelman } 139397b6df1SKris Buschelman /* B-part, larger col index */ 140397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 141397b6df1SKris Buschelman if (!valOnly){ 142397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 143397b6df1SKris Buschelman } 144397b6df1SKris Buschelman val[jj++] = *bv++; 145397b6df1SKris Buschelman } 146397b6df1SKris Buschelman irow++; 147397b6df1SKris Buschelman } 148397b6df1SKris Buschelman 149397b6df1SKris Buschelman PetscFunctionReturn(0); 150397b6df1SKris Buschelman } 151397b6df1SKris Buschelman 152c338a77dSKris Buschelman EXTERN_C_BEGIN 153c338a77dSKris Buschelman #undef __FUNCT__ 154c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base" 15575179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat) \ 156521d7252SBarry Smith { 157dfbe8321SBarry Smith PetscErrorCode ierr; 158c338a77dSKris Buschelman Mat B=*newmat; 159f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 160901853e0SKris Buschelman void (*f)(void); 161c338a77dSKris Buschelman 162c338a77dSKris Buschelman PetscFunctionBegin; 163ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 164c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 165c338a77dSKris Buschelman } 166f0c56d0fSKris Buschelman B->ops->duplicate = mumps->MatDuplicate; 167f0c56d0fSKris Buschelman B->ops->view = mumps->MatView; 168f0c56d0fSKris Buschelman B->ops->assemblyend = mumps->MatAssemblyEnd; 169f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic; 170f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic; 171f0c56d0fSKris Buschelman B->ops->destroy = mumps->MatDestroy; 172901853e0SKris Buschelman 173f68b968cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 174901853e0SKris Buschelman if (f) { 175f68b968cSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr); 176901853e0SKris Buschelman } 177f0c56d0fSKris Buschelman ierr = PetscFree(mumps);CHKERRQ(ierr); 178901853e0SKris Buschelman 179901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 180901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr); 181901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 182901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr); 1832895b8caSSatish Balay ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 1842895b8caSSatish Balay ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr); 185901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr); 186901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr); 187901853e0SKris Buschelman 188901853e0SKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr); 189c338a77dSKris Buschelman *newmat = B; 190c338a77dSKris Buschelman PetscFunctionReturn(0); 191c338a77dSKris Buschelman } 192c338a77dSKris Buschelman EXTERN_C_END 193c338a77dSKris Buschelman 194397b6df1SKris Buschelman #undef __FUNCT__ 1953924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 196dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 197dfbe8321SBarry Smith { 198f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 199dfbe8321SBarry Smith PetscErrorCode ierr; 200c1490034SHong Zhang PetscMPIInt size=lu->size; 2016849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 202397b6df1SKris Buschelman PetscFunctionBegin; 203397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 204397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 205397b6df1SKris Buschelman lu->id.job=JOB_END; 206397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 207397b6df1SKris Buschelman zmumps_c(&lu->id); 208397b6df1SKris Buschelman #else 209397b6df1SKris Buschelman dmumps_c(&lu->id); 210397b6df1SKris Buschelman #endif 211c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 212c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 213c338a77dSKris Buschelman if (size>1 && lu->val) { 214c338a77dSKris Buschelman ierr = PetscFree(lu->val);CHKERRQ(ierr); 215c338a77dSKris Buschelman } 216397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 217397b6df1SKris Buschelman } 218a39386dcSKris Buschelman specialdestroy = lu->specialdestroy; 219a39386dcSKris Buschelman ierr = (*specialdestroy)(A);CHKERRQ(ierr); 220c338a77dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 221397b6df1SKris Buschelman PetscFunctionReturn(0); 222397b6df1SKris Buschelman } 223397b6df1SKris Buschelman 224397b6df1SKris Buschelman #undef __FUNCT__ 225a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS" 226dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A) 227dfbe8321SBarry Smith { 2286849ba73SBarry Smith PetscErrorCode ierr; 2296849ba73SBarry Smith int size; 230a39386dcSKris Buschelman 231a39386dcSKris Buschelman PetscFunctionBegin; 232a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 233a39386dcSKris Buschelman if (size==1) { 234ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 235a39386dcSKris Buschelman } else { 236ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 237a39386dcSKris Buschelman } 238a39386dcSKris Buschelman PetscFunctionReturn(0); 239a39386dcSKris Buschelman } 240a39386dcSKris Buschelman 241a39386dcSKris Buschelman #undef __FUNCT__ 242a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS" 243dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A) 244dfbe8321SBarry Smith { 2456849ba73SBarry Smith PetscErrorCode ierr; 2466849ba73SBarry Smith int size; 247a39386dcSKris Buschelman 248a39386dcSKris Buschelman PetscFunctionBegin; 249a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 250a39386dcSKris Buschelman if (size==1) { 251ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 252a39386dcSKris Buschelman } else { 253ceb03754SKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 254a39386dcSKris Buschelman } 255a39386dcSKris Buschelman PetscFunctionReturn(0); 256a39386dcSKris Buschelman } 257a39386dcSKris Buschelman 258a39386dcSKris Buschelman #undef __FUNCT__ 259c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS" 260dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 261f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 262dfbe8321SBarry Smith PetscErrorCode ierr; 263397b6df1SKris Buschelman 264397b6df1SKris Buschelman PetscFunctionBegin; 265c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 266c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 267c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 268c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 269c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 270c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 271c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 272c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 273c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 274c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 275958c9bccSBarry Smith if (!lu->myid && lu->id.ICNTL(11)>0) { 276c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 277c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 278c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 279c338a77dSKris 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); 280c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 281c338a77dSKris 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); 282c338a77dSKris Buschelman 283c338a77dSKris Buschelman } 284c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 285c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 286adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 287c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(15) (efficiency control): %d \n",lu->id.ICNTL(15));CHKERRQ(ierr); 288c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 289c338a77dSKris Buschelman 290c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 291c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 292c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 29357f0c58bSHong Zhang 29457f0c58bSHong Zhang /* infomation local to each processor */ 295958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 29657f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 29757f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 298958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 29957f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 30057f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 301958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 30257f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 30357f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 304adc1d99fSHong Zhang 305958c9bccSBarry Smith if (!lu->myid){ /* information from the host */ 306adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 307adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 308adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 309adc1d99fSHong Zhang 310adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 311adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 312adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 313adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 314adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 315adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 316adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 317adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 318adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 319adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 320adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 321adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 322adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 323adc1d99fSHong 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); 324adc1d99fSHong 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); 325adc1d99fSHong 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); 326adc1d99fSHong 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); 327adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 328adc1d99fSHong Zhang } 329adc1d99fSHong Zhang 330397b6df1SKris Buschelman PetscFunctionReturn(0); 331397b6df1SKris Buschelman } 332397b6df1SKris Buschelman 333397b6df1SKris Buschelman #undef __FUNCT__ 334c1490034SHong Zhang #define __FUNCT__ "MatView_MUMPS" 335c1490034SHong Zhang PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) { 336dfbe8321SBarry Smith PetscErrorCode ierr; 33732077d6dSBarry Smith PetscTruth iascii; 338397b6df1SKris Buschelman PetscViewerFormat format; 339f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 340397b6df1SKris Buschelman 341397b6df1SKris Buschelman PetscFunctionBegin; 342397b6df1SKris Buschelman ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr); 343397b6df1SKris Buschelman 34432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 34532077d6dSBarry Smith if (iascii) { 346397b6df1SKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 347397b6df1SKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 348397b6df1SKris Buschelman ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 349397b6df1SKris Buschelman } 350397b6df1SKris Buschelman } 351397b6df1SKris Buschelman PetscFunctionReturn(0); 352397b6df1SKris Buschelman } 353397b6df1SKris Buschelman 354397b6df1SKris Buschelman #undef __FUNCT__ 355f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS" 356dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) { 357f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 358d54de34fSKris Buschelman PetscScalar *array; 359397b6df1SKris Buschelman Vec x_seq; 360397b6df1SKris Buschelman IS iden; 361397b6df1SKris Buschelman VecScatter scat; 362dfbe8321SBarry Smith PetscErrorCode ierr; 363397b6df1SKris Buschelman 364397b6df1SKris Buschelman PetscFunctionBegin; 365397b6df1SKris Buschelman if (lu->size > 1){ 366397b6df1SKris Buschelman if (!lu->myid){ 3672a4c71feSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap.N,&x_seq);CHKERRQ(ierr); 3682a4c71feSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap.N,0,1,&iden);CHKERRQ(ierr); 369397b6df1SKris Buschelman } else { 370397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr); 371397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr); 372397b6df1SKris Buschelman } 373397b6df1SKris Buschelman ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr); 374397b6df1SKris Buschelman ierr = ISDestroy(iden);CHKERRQ(ierr); 375397b6df1SKris Buschelman 376397b6df1SKris Buschelman ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 377397b6df1SKris Buschelman ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 378397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 379397b6df1SKris Buschelman } else { /* size == 1 */ 380397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 381397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 382397b6df1SKris Buschelman } 383397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 384397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 385397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 386397b6df1SKris Buschelman #else 387397b6df1SKris Buschelman lu->id.rhs = array; 388397b6df1SKris Buschelman #endif 389397b6df1SKris Buschelman } 390397b6df1SKris Buschelman 391397b6df1SKris Buschelman /* solve phase */ 392397b6df1SKris Buschelman lu->id.job=3; 393397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 394397b6df1SKris Buschelman zmumps_c(&lu->id); 395397b6df1SKris Buschelman #else 396397b6df1SKris Buschelman dmumps_c(&lu->id); 397397b6df1SKris Buschelman #endif 398397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 39979a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 400397b6df1SKris Buschelman } 401397b6df1SKris Buschelman 402397b6df1SKris Buschelman /* convert mumps solution x_seq to petsc mpi x */ 403397b6df1SKris Buschelman if (lu->size > 1) { 404397b6df1SKris Buschelman if (!lu->myid){ 405397b6df1SKris Buschelman ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 406397b6df1SKris Buschelman } 407397b6df1SKris Buschelman ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 408397b6df1SKris Buschelman ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 409397b6df1SKris Buschelman ierr = VecScatterDestroy(scat);CHKERRQ(ierr); 410397b6df1SKris Buschelman ierr = VecDestroy(x_seq);CHKERRQ(ierr); 411397b6df1SKris Buschelman } else { 412397b6df1SKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 413397b6df1SKris Buschelman } 414397b6df1SKris Buschelman 415397b6df1SKris Buschelman PetscFunctionReturn(0); 416397b6df1SKris Buschelman } 417397b6df1SKris Buschelman 418a58c3f20SHong Zhang /* 419a58c3f20SHong Zhang input: 420a58c3f20SHong Zhang F: numeric factor 421a58c3f20SHong Zhang output: 422a58c3f20SHong Zhang nneg: total number of negative pivots 423a58c3f20SHong Zhang nzero: 0 424a58c3f20SHong Zhang npos: (global dimension of F) - nneg 425a58c3f20SHong Zhang */ 426a58c3f20SHong Zhang 427a58c3f20SHong Zhang #undef __FUNCT__ 428a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 429dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 430a58c3f20SHong Zhang { 431a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 432dfbe8321SBarry Smith PetscErrorCode ierr; 433c1490034SHong Zhang PetscMPIInt size; 434a58c3f20SHong Zhang 435a58c3f20SHong Zhang PetscFunctionBegin; 436bcb30aebSHong Zhang ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr); 437bcb30aebSHong 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 */ 438bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 43979a5c55eSBarry 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)); 440bcb30aebSHong Zhang } 441a58c3f20SHong Zhang if (nneg){ 442a58c3f20SHong Zhang if (!lu->myid){ 443a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 444a58c3f20SHong Zhang } 445bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 446a58c3f20SHong Zhang } 447a58c3f20SHong Zhang if (nzero) *nzero = 0; 4482a4c71feSBarry Smith if (npos) *npos = F->rmap.N - (*nneg); 449a58c3f20SHong Zhang PetscFunctionReturn(0); 450a58c3f20SHong Zhang } 451a58c3f20SHong Zhang 452397b6df1SKris Buschelman #undef __FUNCT__ 45319facb7aSBarry Smith #define __FUNCT__ "MatFactorNumeric_AIJMUMPS" 454af281ebdSHong Zhang PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F) 455af281ebdSHong Zhang { 456f0c56d0fSKris Buschelman Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 457f0c56d0fSKris Buschelman Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr; 4586849ba73SBarry Smith PetscErrorCode ierr; 4592a4c71feSBarry Smith PetscInt rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl; 460397b6df1SKris Buschelman PetscTruth valOnly,flg; 461e09efc27SHong Zhang Mat F_diag; 462397b6df1SKris Buschelman 463397b6df1SKris Buschelman PetscFunctionBegin; 464397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 465f0c56d0fSKris Buschelman (*F)->ops->solve = MatSolve_AIJMUMPS; 466397b6df1SKris Buschelman 467397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 468397b6df1SKris Buschelman ierr = MPI_Comm_rank(A->comm, &lu->myid); 469397b6df1SKris Buschelman ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr); 47075747be1SHong Zhang lua->myid = lu->myid; lua->size = lu->size; 471397b6df1SKris Buschelman lu->id.job = JOB_INIT; 472397b6df1SKris Buschelman ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr); 473a0e2756fSSatish Balay ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr); 474397b6df1SKris Buschelman 475397b6df1SKris Buschelman /* Set mumps options */ 476397b6df1SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 477397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 478397b6df1SKris Buschelman lu->id.sym=lu->sym; 479397b6df1SKris Buschelman if (lu->sym == 2){ 480397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 481397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 482397b6df1SKris Buschelman } 483397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 484397b6df1SKris Buschelman zmumps_c(&lu->id); 485397b6df1SKris Buschelman #else 486397b6df1SKris Buschelman dmumps_c(&lu->id); 487397b6df1SKris Buschelman #endif 488397b6df1SKris Buschelman 489397b6df1SKris Buschelman if (lu->size == 1){ 490397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 491397b6df1SKris Buschelman } else { 492397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 493397b6df1SKris Buschelman } 494397b6df1SKris Buschelman 495397b6df1SKris Buschelman icntl=-1; 496397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 49719facb7aSBarry Smith if ((flg && icntl > 0) || PetscLogPrintInfo) { 498397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 499397b6df1SKris Buschelman } else { /* no output */ 500397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 501397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 502397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 503397b6df1SKris Buschelman lu->id.ICNTL(4) = 0; /* level of printing, 0,1,2,3,4, default=2 */ 504397b6df1SKris Buschelman } 505397b6df1SKris 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); 506397b6df1SKris Buschelman icntl=-1; 507397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 508397b6df1SKris Buschelman if (flg) { 509397b6df1SKris Buschelman if (icntl== 1){ 510397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 511397b6df1SKris Buschelman } else { 512397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 513397b6df1SKris Buschelman } 514397b6df1SKris Buschelman } 515397b6df1SKris 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); 516397b6df1SKris 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); 51794b7f48cSBarry 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); 518397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 519397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 520adc1d99fSHong 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); 521397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 522397b6df1SKris Buschelman 523397b6df1SKris 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); 524397b6df1SKris 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); 525397b6df1SKris 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); 526397b6df1SKris Buschelman PetscOptionsEnd(); 527397b6df1SKris Buschelman } 528397b6df1SKris Buschelman 529397b6df1SKris Buschelman /* define matrix A */ 530397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 531397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 532397b6df1SKris Buschelman if (!lu->myid) { 533c36ead0aSKris Buschelman if (lua->isAIJ){ 534397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 535397b6df1SKris Buschelman nz = aa->nz; 536397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 537397b6df1SKris Buschelman } else { 538397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 5396c6c5352SBarry Smith nz = aa->nz; 540397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 541397b6df1SKris Buschelman } 542397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 5437c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 5447c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 545397b6df1SKris Buschelman nz = 0; 546397b6df1SKris Buschelman for (i=0; i<M; i++){ 547397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 548397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 549397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 550397b6df1SKris Buschelman } 551397b6df1SKris Buschelman } 552397b6df1SKris Buschelman } 553397b6df1SKris Buschelman } 554397b6df1SKris Buschelman break; 555397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 556397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 557397b6df1SKris Buschelman valOnly = PETSC_FALSE; 558397b6df1SKris Buschelman } else { 559397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 560397b6df1SKris Buschelman } 561397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 562397b6df1SKris Buschelman break; 563397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 564397b6df1SKris Buschelman } 565397b6df1SKris Buschelman 566397b6df1SKris Buschelman /* analysis phase */ 567397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 568397b6df1SKris Buschelman lu->id.n = M; 569397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 570397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 571397b6df1SKris Buschelman if (!lu->myid) { 572397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 573397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 574397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 575397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 576397b6df1SKris Buschelman #else 577397b6df1SKris Buschelman lu->id.a = lu->val; 578397b6df1SKris Buschelman #endif 579397b6df1SKris Buschelman } 580397b6df1SKris Buschelman } 581397b6df1SKris Buschelman break; 582397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 583397b6df1SKris Buschelman lu->id.nz_loc = nnz; 584397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 585397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 586397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 587397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 588397b6df1SKris Buschelman #else 589397b6df1SKris Buschelman lu->id.a_loc = lu->val; 590397b6df1SKris Buschelman #endif 591397b6df1SKris Buschelman } 592397b6df1SKris Buschelman break; 593397b6df1SKris Buschelman } 594397b6df1SKris Buschelman lu->id.job=1; 595397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 596397b6df1SKris Buschelman zmumps_c(&lu->id); 597397b6df1SKris Buschelman #else 598397b6df1SKris Buschelman dmumps_c(&lu->id); 599397b6df1SKris Buschelman #endif 600397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 60179a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 602397b6df1SKris Buschelman } 603397b6df1SKris Buschelman } 604397b6df1SKris Buschelman 605397b6df1SKris Buschelman /* numerical factorization phase */ 606958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 607a7aca84bSHong Zhang if (!lu->myid) { 608397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 609397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 610397b6df1SKris Buschelman #else 611397b6df1SKris Buschelman lu->id.a = lu->val; 612397b6df1SKris Buschelman #endif 613397b6df1SKris Buschelman } 614397b6df1SKris Buschelman } else { 615397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 616397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 617397b6df1SKris Buschelman #else 618397b6df1SKris Buschelman lu->id.a_loc = lu->val; 619397b6df1SKris Buschelman #endif 620397b6df1SKris Buschelman } 621397b6df1SKris Buschelman lu->id.job=2; 622397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 623397b6df1SKris Buschelman zmumps_c(&lu->id); 624397b6df1SKris Buschelman #else 625397b6df1SKris Buschelman dmumps_c(&lu->id); 626397b6df1SKris Buschelman #endif 627397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 62819facb7aSBarry Smith if (lu->id.INFO(1) == -13) { 62919facb7aSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 63019facb7aSBarry Smith } else { 63179a5c55eSBarry 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)); 632397b6df1SKris Buschelman } 63319facb7aSBarry Smith } 634397b6df1SKris Buschelman 63519facb7aSBarry Smith if (!lu->myid && lu->id.ICNTL(16) > 0){ 63679a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 637397b6df1SKris Buschelman } 638397b6df1SKris Buschelman 6398ada1bb4SHong Zhang if (lu->size > 1){ 640e09efc27SHong Zhang if ((*F)->factor == FACTOR_LU){ 641e09efc27SHong Zhang F_diag = ((Mat_MPIAIJ *)(*F)->data)->A; 642e09efc27SHong Zhang } else { 643e09efc27SHong Zhang F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A; 644e09efc27SHong Zhang } 645e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 6468ada1bb4SHong Zhang } 647397b6df1SKris Buschelman (*F)->assembled = PETSC_TRUE; 648397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 649ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 650397b6df1SKris Buschelman PetscFunctionReturn(0); 651397b6df1SKris Buschelman } 652397b6df1SKris Buschelman 653397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 654397b6df1SKris Buschelman #undef __FUNCT__ 655f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 656dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 657397b6df1SKris Buschelman Mat B; 658f0c56d0fSKris Buschelman Mat_MUMPS *lu; 659dfbe8321SBarry Smith PetscErrorCode ierr; 660397b6df1SKris Buschelman 661397b6df1SKris Buschelman PetscFunctionBegin; 662397b6df1SKris Buschelman /* Create the factorization matrix */ 663f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 6642a4c71feSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 665be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 666397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 667397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 668397b6df1SKris Buschelman 669f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS; 670397b6df1SKris Buschelman B->factor = FACTOR_LU; 671f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 672397b6df1SKris Buschelman lu->sym = 0; 673397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 674397b6df1SKris Buschelman 675397b6df1SKris Buschelman *F = B; 676397b6df1SKris Buschelman PetscFunctionReturn(0); 677397b6df1SKris Buschelman } 678397b6df1SKris Buschelman 679397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 680397b6df1SKris Buschelman #undef __FUNCT__ 681f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 682dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) { 683397b6df1SKris Buschelman Mat B; 684f0c56d0fSKris Buschelman Mat_MUMPS *lu; 685dfbe8321SBarry Smith PetscErrorCode ierr; 686397b6df1SKris Buschelman 687397b6df1SKris Buschelman PetscFunctionBegin; 688397b6df1SKris Buschelman /* Create the factorization matrix */ 689f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&B);CHKERRQ(ierr); 6902a4c71feSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 691be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 692efc670deSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 693efc670deSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 694397b6df1SKris Buschelman 695f0c56d0fSKris Buschelman B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS; 696a58c3f20SHong Zhang B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 697397b6df1SKris Buschelman B->factor = FACTOR_CHOLESKY; 698f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 699397b6df1SKris Buschelman lu->sym = 2; 700397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 701397b6df1SKris Buschelman 702397b6df1SKris Buschelman *F = B; 703397b6df1SKris Buschelman PetscFunctionReturn(0); 704397b6df1SKris Buschelman } 705397b6df1SKris Buschelman 706397b6df1SKris Buschelman #undef __FUNCT__ 707f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS" 708dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) { 709dfbe8321SBarry Smith PetscErrorCode ierr; 710f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 711c338a77dSKris Buschelman 712397b6df1SKris Buschelman PetscFunctionBegin; 713c338a77dSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 714f0c56d0fSKris Buschelman 715c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 716c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 717f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 718397b6df1SKris Buschelman PetscFunctionReturn(0); 719397b6df1SKris Buschelman } 720397b6df1SKris Buschelman 721c338a77dSKris Buschelman EXTERN_C_BEGIN 722c338a77dSKris Buschelman #undef __FUNCT__ 723f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS" 72475179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 725dfbe8321SBarry Smith { 726dfbe8321SBarry Smith PetscErrorCode ierr; 727521d7252SBarry Smith PetscMPIInt size; 728c338a77dSKris Buschelman MPI_Comm comm; 729c338a77dSKris Buschelman Mat B=*newmat; 730f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 731397b6df1SKris Buschelman 732397b6df1SKris Buschelman PetscFunctionBegin; 733ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 734c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 735397b6df1SKris Buschelman } 736397b6df1SKris Buschelman 737c338a77dSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 738f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 739c338a77dSKris Buschelman 740f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 741c338a77dSKris Buschelman mumps->MatView = A->ops->view; 742c338a77dSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 743c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 744c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 745c338a77dSKris Buschelman mumps->MatDestroy = A->ops->destroy; 746a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_AIJMUMPS; 747c338a77dSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 748f579278aSKris Buschelman mumps->isAIJ = PETSC_TRUE; 749c338a77dSKris Buschelman 7504b68dd72SKris Buschelman B->spptr = (void*)mumps; 751422e82a1SHong Zhang B->ops->duplicate = MatDuplicate_MUMPS; 752c1490034SHong Zhang B->ops->view = MatView_MUMPS; 753f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS; 754f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 7553924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 756c338a77dSKris Buschelman 757c338a77dSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 758c338a77dSKris Buschelman if (size == 1) { 759c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C", 760f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 761c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C", 762c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 763c338a77dSKris Buschelman } else { 764c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C", 765f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 766c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C", 767c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 768c338a77dSKris Buschelman } 769c338a77dSKris Buschelman 770ae15b995SBarry Smith ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr); 771c338a77dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 772c338a77dSKris Buschelman *newmat = B; 773397b6df1SKris Buschelman PetscFunctionReturn(0); 774397b6df1SKris Buschelman } 775c338a77dSKris Buschelman EXTERN_C_END 776397b6df1SKris Buschelman 77724b6179bSKris Buschelman /*MC 778fafad747SKris Buschelman MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 77924b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 78024b6179bSKris Buschelman 78124b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 78224b6179bSKris Buschelman on how to declare the existence of external packages), 78324b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 78424b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS). 78524b6179bSKris Buschelman 78624b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 78724b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 78824b6179bSKris Buschelman MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported 78924b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 79028b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 79128b08bd3SKris Buschelman conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 79228b08bd3SKris Buschelman without data copy. 79324b6179bSKris Buschelman 79424b6179bSKris Buschelman Options Database Keys: 7950bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 79624b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 79724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 79824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 79924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 80024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 80124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 80294b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 80324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 80424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 80524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 80624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 80724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 80824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 80924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 81024b6179bSKris Buschelman 81124b6179bSKris Buschelman Level: beginner 81224b6179bSKris Buschelman 81324b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 81424b6179bSKris Buschelman M*/ 81524b6179bSKris Buschelman 816397b6df1SKris Buschelman EXTERN_C_BEGIN 817397b6df1SKris Buschelman #undef __FUNCT__ 818f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS" 819be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A) 820dfbe8321SBarry Smith { 821dfbe8321SBarry Smith PetscErrorCode ierr; 822c1490034SHong Zhang PetscMPIInt size; 823397b6df1SKris Buschelman MPI_Comm comm; 824397b6df1SKris Buschelman 825397b6df1SKris Buschelman PetscFunctionBegin; 8265441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */ 8275441df8eSKris Buschelman /* and AIJMUMPS types */ 8285441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr); 829397b6df1SKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 830397b6df1SKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 831397b6df1SKris Buschelman if (size == 1) { 832397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 833397b6df1SKris Buschelman } else { 834397b6df1SKris Buschelman ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 835*00ff2a26SHong Zhang /* 836*00ff2a26SHong Zhang Mat A_diag = ((Mat_MPIAIJ *)A->data)->A; 837ceb03754SKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr); 838*00ff2a26SHong Zhang */ 839397b6df1SKris Buschelman } 840ceb03754SKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 841397b6df1SKris Buschelman PetscFunctionReturn(0); 842397b6df1SKris Buschelman } 843397b6df1SKris Buschelman EXTERN_C_END 844397b6df1SKris Buschelman 845f579278aSKris Buschelman #undef __FUNCT__ 846f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS" 847dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) 848dfbe8321SBarry Smith { 849dfbe8321SBarry Smith PetscErrorCode ierr; 850f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 851f579278aSKris Buschelman 852f579278aSKris Buschelman PetscFunctionBegin; 853f579278aSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 854f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 855f0c56d0fSKris Buschelman A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 856f579278aSKris Buschelman PetscFunctionReturn(0); 857f579278aSKris Buschelman } 858f579278aSKris Buschelman 859f579278aSKris Buschelman EXTERN_C_BEGIN 860f579278aSKris Buschelman #undef __FUNCT__ 8619c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS" 862c1490034SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz) 8639c097c71SKris Buschelman { 8649c097c71SKris Buschelman Mat A; 86502217bfdSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr; 866dfbe8321SBarry Smith PetscErrorCode ierr; 8679c097c71SKris Buschelman 8689c097c71SKris Buschelman PetscFunctionBegin; 8699c097c71SKris Buschelman /* 8709c097c71SKris Buschelman After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix 8719c097c71SKris Buschelman into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS. I would 8729c097c71SKris Buschelman like this to be done in the MatCreate routine, but the creation of this inner matrix requires 8739c097c71SKris Buschelman block size info so that PETSc can determine the local size properly. The block size info is set 8749c097c71SKris Buschelman in the preallocation routine. 8759c097c71SKris Buschelman */ 8769c097c71SKris Buschelman ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz); 8779c097c71SKris Buschelman A = ((Mat_MPISBAIJ *)B->data)->A; 878ceb03754SKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 8799c097c71SKris Buschelman PetscFunctionReturn(0); 8809c097c71SKris Buschelman } 8819c097c71SKris Buschelman EXTERN_C_END 8829c097c71SKris Buschelman 8839c097c71SKris Buschelman EXTERN_C_BEGIN 8849c097c71SKris Buschelman #undef __FUNCT__ 885f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS" 88675179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 887dfbe8321SBarry Smith { 888dfbe8321SBarry Smith PetscErrorCode ierr; 889521d7252SBarry Smith PetscMPIInt size; 890f579278aSKris Buschelman MPI_Comm comm; 891f579278aSKris Buschelman Mat B=*newmat; 892422e82a1SHong Zhang Mat_MUMPS *mumps; 8939c097c71SKris Buschelman void (*f)(void); 894f579278aSKris Buschelman 895f579278aSKris Buschelman PetscFunctionBegin; 896ceb03754SKris Buschelman if (reuse == MAT_INITIAL_MATRIX) { 897f579278aSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 898f579278aSKris Buschelman } 899f579278aSKris Buschelman 900f579278aSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 901f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 902f579278aSKris Buschelman 903f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 904f579278aSKris Buschelman mumps->MatView = A->ops->view; 905f579278aSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 906f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 907f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 908f579278aSKris Buschelman mumps->MatDestroy = A->ops->destroy; 909a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_SBAIJMUMPS; 910f579278aSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 911f579278aSKris Buschelman mumps->isAIJ = PETSC_FALSE; 912f579278aSKris Buschelman 913f579278aSKris Buschelman B->spptr = (void*)mumps; 914422e82a1SHong Zhang B->ops->duplicate = MatDuplicate_MUMPS; 915c1490034SHong Zhang B->ops->view = MatView_MUMPS; 916f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS; 917f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 9183924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 919f579278aSKris Buschelman 920f579278aSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 921f579278aSKris Buschelman if (size == 1) { 922f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C", 923f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 924f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C", 925f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 926f579278aSKris Buschelman } else { 9279c097c71SKris Buschelman /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */ 928f68b968cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr); 929901853e0SKris Buschelman if (f) { /* This case should always be true when this routine is called */ 9306849ba73SBarry Smith mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f; 9319c097c71SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C", 9329c097c71SKris Buschelman "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS", 9339c097c71SKris Buschelman MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr); 9349c097c71SKris Buschelman } 935f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C", 936f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 937f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C", 938f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 939f579278aSKris Buschelman } 940f579278aSKris Buschelman 941ae15b995SBarry Smith ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr); 942f579278aSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 943f579278aSKris Buschelman *newmat = B; 944f579278aSKris Buschelman PetscFunctionReturn(0); 945f579278aSKris Buschelman } 946f579278aSKris Buschelman EXTERN_C_END 947f579278aSKris Buschelman 948f0c56d0fSKris Buschelman #undef __FUNCT__ 949422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS" 950dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) { 951dfbe8321SBarry Smith PetscErrorCode ierr; 9528e393735SKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr; 9538f340917SKris Buschelman 954f0c56d0fSKris Buschelman PetscFunctionBegin; 9558f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 9568e393735SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr); 957f0c56d0fSKris Buschelman PetscFunctionReturn(0); 958f0c56d0fSKris Buschelman } 959f0c56d0fSKris Buschelman 96024b6179bSKris Buschelman /*MC 961fafad747SKris Buschelman MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 96224b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 96324b6179bSKris Buschelman 96424b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 96524b6179bSKris Buschelman on how to declare the existence of external packages), 96624b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 96724b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS). 96824b6179bSKris Buschelman 96924b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 97024b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 97124b6179bSKris Buschelman MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported 97224b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 97328b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 97428b08bd3SKris Buschelman conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 97528b08bd3SKris Buschelman without data copy. 97624b6179bSKris Buschelman 97724b6179bSKris Buschelman Options Database Keys: 9780bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 97924b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 98024b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 98124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 98224b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 98324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 98424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 98594b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 98624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 98724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 98824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 98924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 99024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 99124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 99224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 99324b6179bSKris Buschelman 99424b6179bSKris Buschelman Level: beginner 99524b6179bSKris Buschelman 99624b6179bSKris Buschelman .seealso: MATAIJMUMPS 99724b6179bSKris Buschelman M*/ 99824b6179bSKris Buschelman 999397b6df1SKris Buschelman EXTERN_C_BEGIN 1000397b6df1SKris Buschelman #undef __FUNCT__ 1001f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS" 1002be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A) 1003dfbe8321SBarry Smith { 10046849ba73SBarry Smith PetscErrorCode ierr; 10058ada1bb4SHong Zhang PetscMPIInt size; 1006397b6df1SKris Buschelman 1007397b6df1SKris Buschelman PetscFunctionBegin; 10085441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */ 10095441df8eSKris Buschelman /* and SBAIJMUMPS types */ 10105441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr); 10115441df8eSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 1012397b6df1SKris Buschelman if (size == 1) { 1013397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 1014397b6df1SKris Buschelman } else { 1015397b6df1SKris Buschelman ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr); 1016397b6df1SKris Buschelman } 1017ceb03754SKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr); 1018397b6df1SKris Buschelman PetscFunctionReturn(0); 1019397b6df1SKris Buschelman } 1020397b6df1SKris Buschelman EXTERN_C_END 1021