11c2a3de1SBarry Smith 2397b6df1SKris Buschelman /* 3a58c3f20SHong Zhang Provides an interface to the MUMPS_4.3.1 sparse solver 4397b6df1SKris Buschelman */ 5397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 6397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h" 7397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h" 8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h" 9397b6df1SKris Buschelman 10397b6df1SKris Buschelman EXTERN_C_BEGIN 11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 12397b6df1SKris Buschelman #include "zmumps_c.h" 13397b6df1SKris Buschelman #else 14397b6df1SKris Buschelman #include "dmumps_c.h" 15397b6df1SKris Buschelman #endif 16397b6df1SKris Buschelman EXTERN_C_END 17397b6df1SKris Buschelman #define JOB_INIT -1 18397b6df1SKris Buschelman #define JOB_END -2 19397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 20397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 21397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 22397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 23a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 24397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 25adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 26397b6df1SKris Buschelman 27397b6df1SKris Buschelman typedef struct { 28397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 29397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 30397b6df1SKris Buschelman #else 31397b6df1SKris Buschelman DMUMPS_STRUC_C id; 32397b6df1SKris Buschelman #endif 33397b6df1SKris Buschelman MatStructure matstruc; 34397b6df1SKris Buschelman int myid,size,*irn,*jcn,sym; 35397b6df1SKris Buschelman PetscScalar *val; 36397b6df1SKris Buschelman MPI_Comm comm_mumps; 37397b6df1SKris Buschelman 38c338a77dSKris Buschelman PetscTruth isAIJ,CleanUpMUMPS; 39*6849ba73SBarry Smith PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 40*6849ba73SBarry Smith PetscErrorCode (*MatView)(Mat,PetscViewer); 41*6849ba73SBarry Smith PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType); 42*6849ba73SBarry Smith PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 43*6849ba73SBarry Smith PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*); 44*6849ba73SBarry Smith PetscErrorCode (*MatDestroy)(Mat); 45*6849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 46*6849ba73SBarry Smith PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*); 47f0c56d0fSKris Buschelman } Mat_MUMPS; 48f0c56d0fSKris Buschelman 49dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 50892f6c3fSKris Buschelman EXTERN_C_BEGIN 51dfbe8321SBarry Smith PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,Mat*); 52892f6c3fSKris Buschelman EXTERN_C_END 53397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 54397b6df1SKris Buschelman /* 55397b6df1SKris Buschelman input: 5675747be1SHong Zhang A - matrix in mpiaij or mpisbaij (bs=1) format 57397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 58397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 59397b6df1SKris Buschelman TRUE: only the values in v array are updated 60397b6df1SKris Buschelman output: 61397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 62397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 63397b6df1SKris Buschelman */ 64dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) { 65397b6df1SKris Buschelman int *ai, *aj, *bi, *bj, rstart,nz, *garray; 66dfbe8321SBarry Smith PetscErrorCode ierr; 67dfbe8321SBarry Smith int i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol; 68d54de34fSKris Buschelman int *row,*col; 69397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 70f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 71397b6df1SKris Buschelman 72397b6df1SKris Buschelman PetscFunctionBegin; 73397b6df1SKris Buschelman if (mumps->isAIJ){ 74397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 75397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 76397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 77397b6df1SKris Buschelman nz = aa->nz + bb->nz; 78397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 79397b6df1SKris Buschelman garray = mat->garray; 80397b6df1SKris Buschelman av=aa->a; bv=bb->a; 81397b6df1SKris Buschelman 82397b6df1SKris Buschelman } else { 83397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 84397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 85397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 86847143adSKris Buschelman if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs); 876c6c5352SBarry Smith nz = aa->nz + bb->nz; 88397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 89397b6df1SKris Buschelman garray = mat->garray; 90397b6df1SKris Buschelman av=aa->a; bv=bb->a; 91397b6df1SKris Buschelman } 92397b6df1SKris Buschelman 93397b6df1SKris Buschelman if (!valOnly){ 94397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr); 95397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr); 96397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 97397b6df1SKris Buschelman *r = row; *c = col; *v = val; 98397b6df1SKris Buschelman } else { 99397b6df1SKris Buschelman row = *r; col = *c; val = *v; 100397b6df1SKris Buschelman } 101397b6df1SKris Buschelman *nnz = nz; 102397b6df1SKris Buschelman 103028e57e8SHong Zhang jj = 0; irow = rstart; 104397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 105397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 106397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 107397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 108397b6df1SKris Buschelman bjj = bj + bi[i]; 109397b6df1SKris Buschelman 110397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 111397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 11275747be1SHong Zhang j=-1; 11375747be1SHong Zhang do { 11475747be1SHong Zhang j++; 11575747be1SHong Zhang if (j == countB) break; 116397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11775747be1SHong Zhang } while (jcol < colA_start); 11875747be1SHong Zhang jB = j; 119397b6df1SKris Buschelman 120397b6df1SKris Buschelman /* B-part, smaller col index */ 121397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 122397b6df1SKris Buschelman for (j=0; j<jB; j++){ 123397b6df1SKris Buschelman jcol = garray[bjj[j]]; 124397b6df1SKris Buschelman if (!valOnly){ 125397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 12675747be1SHong Zhang 127397b6df1SKris Buschelman } 128397b6df1SKris Buschelman val[jj++] = *bv++; 129397b6df1SKris Buschelman } 130397b6df1SKris Buschelman /* A-part */ 131397b6df1SKris Buschelman for (j=0; j<countA; j++){ 132397b6df1SKris Buschelman if (!valOnly){ 133397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 134397b6df1SKris Buschelman } 135397b6df1SKris Buschelman val[jj++] = *av++; 136397b6df1SKris Buschelman } 137397b6df1SKris Buschelman /* B-part, larger col index */ 138397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 139397b6df1SKris Buschelman if (!valOnly){ 140397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 141397b6df1SKris Buschelman } 142397b6df1SKris Buschelman val[jj++] = *bv++; 143397b6df1SKris Buschelman } 144397b6df1SKris Buschelman irow++; 145397b6df1SKris Buschelman } 146397b6df1SKris Buschelman 147397b6df1SKris Buschelman PetscFunctionReturn(0); 148397b6df1SKris Buschelman } 149397b6df1SKris Buschelman 150c338a77dSKris Buschelman EXTERN_C_BEGIN 151c338a77dSKris Buschelman #undef __FUNCT__ 152c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base" 153dfbe8321SBarry Smith PetscErrorCode MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) { 154dfbe8321SBarry Smith PetscErrorCode ierr; 155c338a77dSKris Buschelman Mat B=*newmat; 156f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 157c338a77dSKris Buschelman 158c338a77dSKris Buschelman PetscFunctionBegin; 159c338a77dSKris Buschelman if (B != A) { 160c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 161c338a77dSKris Buschelman } 162f0c56d0fSKris Buschelman B->ops->duplicate = mumps->MatDuplicate; 163f0c56d0fSKris Buschelman B->ops->view = mumps->MatView; 164f0c56d0fSKris Buschelman B->ops->assemblyend = mumps->MatAssemblyEnd; 165f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic; 166f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic; 167f0c56d0fSKris Buschelman B->ops->destroy = mumps->MatDestroy; 1683924e44cSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr); 169f0c56d0fSKris Buschelman ierr = PetscFree(mumps);CHKERRQ(ierr); 170c338a77dSKris Buschelman *newmat = B; 171c338a77dSKris Buschelman PetscFunctionReturn(0); 172c338a77dSKris Buschelman } 173c338a77dSKris Buschelman EXTERN_C_END 174c338a77dSKris Buschelman 175397b6df1SKris Buschelman #undef __FUNCT__ 1763924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 177dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 178dfbe8321SBarry Smith { 179f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 180dfbe8321SBarry Smith PetscErrorCode ierr; 181dfbe8321SBarry Smith int size=lu->size; 182*6849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 183397b6df1SKris Buschelman PetscFunctionBegin; 184397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 185397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 186397b6df1SKris Buschelman lu->id.job=JOB_END; 187397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 188397b6df1SKris Buschelman zmumps_c(&lu->id); 189397b6df1SKris Buschelman #else 190397b6df1SKris Buschelman dmumps_c(&lu->id); 191397b6df1SKris Buschelman #endif 192c338a77dSKris Buschelman if (lu->irn) { 193c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 194c338a77dSKris Buschelman } 195c338a77dSKris Buschelman if (lu->jcn) { 196c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 197c338a77dSKris Buschelman } 198c338a77dSKris Buschelman if (size>1 && lu->val) { 199c338a77dSKris Buschelman ierr = PetscFree(lu->val);CHKERRQ(ierr); 200c338a77dSKris Buschelman } 201397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 202397b6df1SKris Buschelman } 203a39386dcSKris Buschelman specialdestroy = lu->specialdestroy; 204a39386dcSKris Buschelman ierr = (*specialdestroy)(A);CHKERRQ(ierr); 205c338a77dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 206397b6df1SKris Buschelman PetscFunctionReturn(0); 207397b6df1SKris Buschelman } 208397b6df1SKris Buschelman 209397b6df1SKris Buschelman #undef __FUNCT__ 210a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS" 211dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A) 212dfbe8321SBarry Smith { 213*6849ba73SBarry Smith PetscErrorCode ierr; 214*6849ba73SBarry Smith int size; 215a39386dcSKris Buschelman 216a39386dcSKris Buschelman PetscFunctionBegin; 217a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 218a39386dcSKris Buschelman if (size==1) { 219a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr); 220a39386dcSKris Buschelman } else { 221a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr); 222a39386dcSKris Buschelman } 223a39386dcSKris Buschelman PetscFunctionReturn(0); 224a39386dcSKris Buschelman } 225a39386dcSKris Buschelman 226a39386dcSKris Buschelman #undef __FUNCT__ 227a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS" 228dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A) 229dfbe8321SBarry Smith { 230*6849ba73SBarry Smith PetscErrorCode ierr; 231*6849ba73SBarry Smith int size; 232a39386dcSKris Buschelman 233a39386dcSKris Buschelman PetscFunctionBegin; 234a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 235a39386dcSKris Buschelman if (size==1) { 236a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr); 237a39386dcSKris Buschelman } else { 238a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr); 239a39386dcSKris Buschelman } 240a39386dcSKris Buschelman PetscFunctionReturn(0); 241a39386dcSKris Buschelman } 242a39386dcSKris Buschelman 243a39386dcSKris Buschelman #undef __FUNCT__ 244c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS" 245dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 246f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 247dfbe8321SBarry Smith PetscErrorCode ierr; 248397b6df1SKris Buschelman 249397b6df1SKris Buschelman PetscFunctionBegin; 250c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 251c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 252c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 253c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 254c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 255c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 256c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 257c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 258c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 259c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 260958c9bccSBarry Smith if (!lu->myid && lu->id.ICNTL(11)>0) { 261c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 262c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 263c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 264c338a77dSKris 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); 265c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 266c338a77dSKris 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); 267c338a77dSKris Buschelman 268c338a77dSKris Buschelman } 269c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 270c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 271adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 272c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(15) (efficiency control): %d \n",lu->id.ICNTL(15));CHKERRQ(ierr); 273c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 274c338a77dSKris Buschelman 275c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 276c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 277c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 27857f0c58bSHong Zhang 27957f0c58bSHong Zhang /* infomation local to each processor */ 280958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 28157f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 28257f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 283958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 28457f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 28557f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 286958c9bccSBarry Smith if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 28757f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 28857f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 289adc1d99fSHong Zhang 290958c9bccSBarry Smith if (!lu->myid){ /* information from the host */ 291adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 292adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 293adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 294adc1d99fSHong Zhang 295adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 296adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 297adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 298adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 299adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 300adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 301adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 302adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 303adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 304adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 305adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 306adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 307adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 308adc1d99fSHong 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); 309adc1d99fSHong 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); 310adc1d99fSHong 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); 311adc1d99fSHong 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); 312adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 313adc1d99fSHong Zhang } 314adc1d99fSHong Zhang 315397b6df1SKris Buschelman PetscFunctionReturn(0); 316397b6df1SKris Buschelman } 317397b6df1SKris Buschelman 318397b6df1SKris Buschelman #undef __FUNCT__ 319f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS" 320dfbe8321SBarry Smith PetscErrorCode MatView_AIJMUMPS(Mat A,PetscViewer viewer) { 321dfbe8321SBarry Smith PetscErrorCode ierr; 32232077d6dSBarry Smith PetscTruth iascii; 323397b6df1SKris Buschelman PetscViewerFormat format; 324f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 325397b6df1SKris Buschelman 326397b6df1SKris Buschelman PetscFunctionBegin; 327397b6df1SKris Buschelman ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr); 328397b6df1SKris Buschelman 32932077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 33032077d6dSBarry Smith if (iascii) { 331397b6df1SKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 332397b6df1SKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 333397b6df1SKris Buschelman ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 334397b6df1SKris Buschelman } 335397b6df1SKris Buschelman } 336397b6df1SKris Buschelman PetscFunctionReturn(0); 337397b6df1SKris Buschelman } 338397b6df1SKris Buschelman 339397b6df1SKris Buschelman #undef __FUNCT__ 340f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS" 341dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) { 342f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 343d54de34fSKris Buschelman PetscScalar *array; 344397b6df1SKris Buschelman Vec x_seq; 345397b6df1SKris Buschelman IS iden; 346397b6df1SKris Buschelman VecScatter scat; 347dfbe8321SBarry Smith PetscErrorCode ierr; 348397b6df1SKris Buschelman 349397b6df1SKris Buschelman PetscFunctionBegin; 350397b6df1SKris Buschelman if (lu->size > 1){ 351397b6df1SKris Buschelman if (!lu->myid){ 352397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr); 353397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr); 354397b6df1SKris Buschelman } else { 355397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr); 356397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr); 357397b6df1SKris Buschelman } 358397b6df1SKris Buschelman ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr); 359397b6df1SKris Buschelman ierr = ISDestroy(iden);CHKERRQ(ierr); 360397b6df1SKris Buschelman 361397b6df1SKris Buschelman ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 362397b6df1SKris Buschelman ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 363397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 364397b6df1SKris Buschelman } else { /* size == 1 */ 365397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 366397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 367397b6df1SKris Buschelman } 368397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 369397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 370397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 371397b6df1SKris Buschelman #else 372397b6df1SKris Buschelman lu->id.rhs = array; 373397b6df1SKris Buschelman #endif 374397b6df1SKris Buschelman } 375397b6df1SKris Buschelman 376397b6df1SKris Buschelman /* solve phase */ 377397b6df1SKris Buschelman lu->id.job=3; 378397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 379397b6df1SKris Buschelman zmumps_c(&lu->id); 380397b6df1SKris Buschelman #else 381397b6df1SKris Buschelman dmumps_c(&lu->id); 382397b6df1SKris Buschelman #endif 383397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 384397b6df1SKris Buschelman SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 385397b6df1SKris Buschelman } 386397b6df1SKris Buschelman 387397b6df1SKris Buschelman /* convert mumps solution x_seq to petsc mpi x */ 388397b6df1SKris Buschelman if (lu->size > 1) { 389397b6df1SKris Buschelman if (!lu->myid){ 390397b6df1SKris Buschelman ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 391397b6df1SKris Buschelman } 392397b6df1SKris Buschelman ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 393397b6df1SKris Buschelman ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 394397b6df1SKris Buschelman ierr = VecScatterDestroy(scat);CHKERRQ(ierr); 395397b6df1SKris Buschelman ierr = VecDestroy(x_seq);CHKERRQ(ierr); 396397b6df1SKris Buschelman } else { 397397b6df1SKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 398397b6df1SKris Buschelman } 399397b6df1SKris Buschelman 400397b6df1SKris Buschelman PetscFunctionReturn(0); 401397b6df1SKris Buschelman } 402397b6df1SKris Buschelman 403a58c3f20SHong Zhang /* 404a58c3f20SHong Zhang input: 405a58c3f20SHong Zhang F: numeric factor 406a58c3f20SHong Zhang output: 407a58c3f20SHong Zhang nneg: total number of negative pivots 408a58c3f20SHong Zhang nzero: 0 409a58c3f20SHong Zhang npos: (global dimension of F) - nneg 410a58c3f20SHong Zhang */ 411a58c3f20SHong Zhang 412a58c3f20SHong Zhang #undef __FUNCT__ 413a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 414dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 415a58c3f20SHong Zhang { 416a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 417dfbe8321SBarry Smith PetscErrorCode ierr; 418dfbe8321SBarry Smith int size; 419a58c3f20SHong Zhang 420a58c3f20SHong Zhang PetscFunctionBegin; 421bcb30aebSHong Zhang ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr); 422bcb30aebSHong 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 */ 423bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 424bcb30aebSHong Zhang SETERRQ1(1,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13)); 425bcb30aebSHong Zhang } 426a58c3f20SHong Zhang if (nneg){ 427a58c3f20SHong Zhang if (!lu->myid){ 428a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 429a58c3f20SHong Zhang } 430bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 431a58c3f20SHong Zhang } 432a58c3f20SHong Zhang if (nzero) *nzero = 0; 433a58c3f20SHong Zhang if (npos) *npos = F->M - (*nneg); 434a58c3f20SHong Zhang PetscFunctionReturn(0); 435a58c3f20SHong Zhang } 436a58c3f20SHong Zhang 437397b6df1SKris Buschelman #undef __FUNCT__ 438f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS" 439dfbe8321SBarry Smith PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) { 440f0c56d0fSKris Buschelman Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 441f0c56d0fSKris Buschelman Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr; 442*6849ba73SBarry Smith PetscErrorCode ierr; 443*6849ba73SBarry Smith int rnz,nnz,nz,i,M=A->M,*ai,*aj,icntl; 444397b6df1SKris Buschelman PetscTruth valOnly,flg; 445397b6df1SKris Buschelman 446397b6df1SKris Buschelman PetscFunctionBegin; 447397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 448f0c56d0fSKris Buschelman (*F)->ops->solve = MatSolve_AIJMUMPS; 449397b6df1SKris Buschelman 450397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 451397b6df1SKris Buschelman ierr = MPI_Comm_rank(A->comm, &lu->myid); 452397b6df1SKris Buschelman ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr); 45375747be1SHong Zhang lua->myid = lu->myid; lua->size = lu->size; 454397b6df1SKris Buschelman lu->id.job = JOB_INIT; 455397b6df1SKris Buschelman ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr); 456397b6df1SKris Buschelman lu->id.comm_fortran = lu->comm_mumps; 457397b6df1SKris Buschelman 458397b6df1SKris Buschelman /* Set mumps options */ 459397b6df1SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 460397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 461397b6df1SKris Buschelman lu->id.sym=lu->sym; 462397b6df1SKris Buschelman if (lu->sym == 2){ 463397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 464397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 465397b6df1SKris Buschelman } 466397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 467397b6df1SKris Buschelman zmumps_c(&lu->id); 468397b6df1SKris Buschelman #else 469397b6df1SKris Buschelman dmumps_c(&lu->id); 470397b6df1SKris Buschelman #endif 471397b6df1SKris Buschelman 472397b6df1SKris Buschelman if (lu->size == 1){ 473397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 474397b6df1SKris Buschelman } else { 475397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 476397b6df1SKris Buschelman } 477397b6df1SKris Buschelman 478397b6df1SKris Buschelman icntl=-1; 479397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 480397b6df1SKris Buschelman if (flg && icntl > 0) { 481397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 482397b6df1SKris Buschelman } else { /* no output */ 483397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 484397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 485397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 486397b6df1SKris Buschelman lu->id.ICNTL(4) = 0; /* level of printing, 0,1,2,3,4, default=2 */ 487397b6df1SKris Buschelman } 488397b6df1SKris 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); 489397b6df1SKris Buschelman icntl=-1; 490397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 491397b6df1SKris Buschelman if (flg) { 492397b6df1SKris Buschelman if (icntl== 1){ 493397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 494397b6df1SKris Buschelman } else { 495397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 496397b6df1SKris Buschelman } 497397b6df1SKris Buschelman } 498397b6df1SKris 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); 499397b6df1SKris 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); 50094b7f48cSBarry 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); 501397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 502397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 503adc1d99fSHong 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); 504397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 505397b6df1SKris Buschelman 506397b6df1SKris Buschelman /* 507397b6df1SKris 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); 508397b6df1SKris Buschelman if (flg){ 509397b6df1SKris Buschelman if (icntl >-1 && icntl <3 ){ 510397b6df1SKris Buschelman if (lu->myid==0) lu->id.ICNTL(16) = icntl; 511397b6df1SKris Buschelman } else { 512397b6df1SKris Buschelman SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl); 513397b6df1SKris Buschelman } 514397b6df1SKris Buschelman } 515397b6df1SKris Buschelman */ 516397b6df1SKris Buschelman 517397b6df1SKris 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); 518397b6df1SKris 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); 519397b6df1SKris 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); 520397b6df1SKris Buschelman PetscOptionsEnd(); 521397b6df1SKris Buschelman } 522397b6df1SKris Buschelman 523397b6df1SKris Buschelman /* define matrix A */ 524397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 525397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 526397b6df1SKris Buschelman if (!lu->myid) { 527c36ead0aSKris Buschelman if (lua->isAIJ){ 528397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 529397b6df1SKris Buschelman nz = aa->nz; 530397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 531397b6df1SKris Buschelman } else { 532397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 5336c6c5352SBarry Smith nz = aa->nz; 534397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 535397b6df1SKris Buschelman } 536397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 537397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr); 538397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr); 539397b6df1SKris Buschelman nz = 0; 540397b6df1SKris Buschelman for (i=0; i<M; i++){ 541397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 542397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 543397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 544397b6df1SKris Buschelman } 545397b6df1SKris Buschelman } 546397b6df1SKris Buschelman } 547397b6df1SKris Buschelman } 548397b6df1SKris Buschelman break; 549397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 550397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 551397b6df1SKris Buschelman valOnly = PETSC_FALSE; 552397b6df1SKris Buschelman } else { 553397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 554397b6df1SKris Buschelman } 555397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 556397b6df1SKris Buschelman break; 557397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 558397b6df1SKris Buschelman } 559397b6df1SKris Buschelman 560397b6df1SKris Buschelman /* analysis phase */ 561397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 562397b6df1SKris Buschelman lu->id.n = M; 563397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 564397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 565397b6df1SKris Buschelman if (!lu->myid) { 566397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 567397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 568397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 569397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 570397b6df1SKris Buschelman #else 571397b6df1SKris Buschelman lu->id.a = lu->val; 572397b6df1SKris Buschelman #endif 573397b6df1SKris Buschelman } 574397b6df1SKris Buschelman } 575397b6df1SKris Buschelman break; 576397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 577397b6df1SKris Buschelman lu->id.nz_loc = nnz; 578397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 579397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 580397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 581397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 582397b6df1SKris Buschelman #else 583397b6df1SKris Buschelman lu->id.a_loc = lu->val; 584397b6df1SKris Buschelman #endif 585397b6df1SKris Buschelman } 586397b6df1SKris Buschelman break; 587397b6df1SKris Buschelman } 588397b6df1SKris Buschelman lu->id.job=1; 589397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 590397b6df1SKris Buschelman zmumps_c(&lu->id); 591397b6df1SKris Buschelman #else 592397b6df1SKris Buschelman dmumps_c(&lu->id); 593397b6df1SKris Buschelman #endif 594397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 595397b6df1SKris Buschelman SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 596397b6df1SKris Buschelman } 597397b6df1SKris Buschelman } 598397b6df1SKris Buschelman 599397b6df1SKris Buschelman /* numerical factorization phase */ 600958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 601a7aca84bSHong Zhang if (!lu->myid) { 602397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 603397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 604397b6df1SKris Buschelman #else 605397b6df1SKris Buschelman lu->id.a = lu->val; 606397b6df1SKris Buschelman #endif 607397b6df1SKris Buschelman } 608397b6df1SKris Buschelman } else { 609397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 610397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 611397b6df1SKris Buschelman #else 612397b6df1SKris Buschelman lu->id.a_loc = lu->val; 613397b6df1SKris Buschelman #endif 614397b6df1SKris Buschelman } 615397b6df1SKris Buschelman lu->id.job=2; 616397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 617397b6df1SKris Buschelman zmumps_c(&lu->id); 618397b6df1SKris Buschelman #else 619397b6df1SKris Buschelman dmumps_c(&lu->id); 620397b6df1SKris Buschelman #endif 621397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 622a7aca84bSHong Zhang SETERRQ2(1,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2)); 623397b6df1SKris Buschelman } 624397b6df1SKris Buschelman 625397b6df1SKris Buschelman if (lu->myid==0 && lu->id.ICNTL(16) > 0){ 626397b6df1SKris Buschelman SETERRQ1(1," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 627397b6df1SKris Buschelman } 628397b6df1SKris Buschelman 629397b6df1SKris Buschelman (*F)->assembled = PETSC_TRUE; 630397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 631ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 632397b6df1SKris Buschelman PetscFunctionReturn(0); 633397b6df1SKris Buschelman } 634397b6df1SKris Buschelman 635397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 636397b6df1SKris Buschelman #undef __FUNCT__ 637f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 638dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 639397b6df1SKris Buschelman Mat B; 640f0c56d0fSKris Buschelman Mat_MUMPS *lu; 641dfbe8321SBarry Smith PetscErrorCode ierr; 642397b6df1SKris Buschelman 643397b6df1SKris Buschelman PetscFunctionBegin; 644397b6df1SKris Buschelman 645397b6df1SKris Buschelman /* Create the factorization matrix */ 646397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 647be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 648397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 649397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 650397b6df1SKris Buschelman 651f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS; 652397b6df1SKris Buschelman B->factor = FACTOR_LU; 653f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 654397b6df1SKris Buschelman lu->sym = 0; 655397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 656397b6df1SKris Buschelman 657397b6df1SKris Buschelman *F = B; 658397b6df1SKris Buschelman PetscFunctionReturn(0); 659397b6df1SKris Buschelman } 660397b6df1SKris Buschelman 661397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 662397b6df1SKris Buschelman #undef __FUNCT__ 663f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 664dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) { 665397b6df1SKris Buschelman Mat B; 666f0c56d0fSKris Buschelman Mat_MUMPS *lu; 667dfbe8321SBarry Smith PetscErrorCode ierr; 668397b6df1SKris Buschelman 669397b6df1SKris Buschelman PetscFunctionBegin; 670397b6df1SKris Buschelman 671397b6df1SKris Buschelman /* Create the factorization matrix */ 672397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 673be5d1d56SKris Buschelman ierr = MatSetType(B,A->type_name);CHKERRQ(ierr); 674efc670deSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 675efc670deSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 676397b6df1SKris Buschelman 677f0c56d0fSKris Buschelman B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS; 678a58c3f20SHong Zhang B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 679397b6df1SKris Buschelman B->factor = FACTOR_CHOLESKY; 680f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 681397b6df1SKris Buschelman lu->sym = 2; 682397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 683397b6df1SKris Buschelman 684397b6df1SKris Buschelman *F = B; 685397b6df1SKris Buschelman PetscFunctionReturn(0); 686397b6df1SKris Buschelman } 687397b6df1SKris Buschelman 688397b6df1SKris Buschelman #undef __FUNCT__ 689f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS" 690dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) { 691dfbe8321SBarry Smith PetscErrorCode ierr; 692f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 693c338a77dSKris Buschelman 694397b6df1SKris Buschelman PetscFunctionBegin; 695c338a77dSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 696f0c56d0fSKris Buschelman 697c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 698c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 699f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 700397b6df1SKris Buschelman PetscFunctionReturn(0); 701397b6df1SKris Buschelman } 702397b6df1SKris Buschelman 703c338a77dSKris Buschelman EXTERN_C_BEGIN 704c338a77dSKris Buschelman #undef __FUNCT__ 705f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS" 706dfbe8321SBarry Smith PetscErrorCode MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat)\ 707dfbe8321SBarry Smith { 708dfbe8321SBarry Smith PetscErrorCode ierr; 709dfbe8321SBarry Smith int size; 710c338a77dSKris Buschelman MPI_Comm comm; 711c338a77dSKris Buschelman Mat B=*newmat; 712f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 713397b6df1SKris Buschelman 714397b6df1SKris Buschelman PetscFunctionBegin; 715c338a77dSKris Buschelman if (B != A) { 716c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 717397b6df1SKris Buschelman } 718397b6df1SKris Buschelman 719c338a77dSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 720f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 721c338a77dSKris Buschelman 722f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 723c338a77dSKris Buschelman mumps->MatView = A->ops->view; 724c338a77dSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 725c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 726c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 727c338a77dSKris Buschelman mumps->MatDestroy = A->ops->destroy; 728a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_AIJMUMPS; 729c338a77dSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 730f579278aSKris Buschelman mumps->isAIJ = PETSC_TRUE; 731c338a77dSKris Buschelman 7324b68dd72SKris Buschelman B->spptr = (void*)mumps; 733422e82a1SHong Zhang B->ops->duplicate = MatDuplicate_MUMPS; 734f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 735f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS; 736f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 7373924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 738c338a77dSKris Buschelman 739c338a77dSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 740c338a77dSKris Buschelman if (size == 1) { 741c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C", 742f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 743c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C", 744c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 745c338a77dSKris Buschelman } else { 746c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C", 747f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 748c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C", 749c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 750c338a77dSKris Buschelman } 751c338a77dSKris Buschelman 752f579278aSKris Buschelman PetscLogInfo(0,"Using MUMPS for LU factorization and solves."); 753c338a77dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 754c338a77dSKris Buschelman *newmat = B; 755397b6df1SKris Buschelman PetscFunctionReturn(0); 756397b6df1SKris Buschelman } 757c338a77dSKris Buschelman EXTERN_C_END 758397b6df1SKris Buschelman 75924b6179bSKris Buschelman /*MC 760fafad747SKris Buschelman MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 76124b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 76224b6179bSKris Buschelman 76324b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 76424b6179bSKris Buschelman on how to declare the existence of external packages), 76524b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 76624b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS). 76724b6179bSKris Buschelman This matrix type is only supported for double precision real. 76824b6179bSKris Buschelman 76924b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 77024b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 77124b6179bSKris Buschelman MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported 77224b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 77328b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 77428b08bd3SKris Buschelman conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 77528b08bd3SKris Buschelman without data copy. 77624b6179bSKris Buschelman 77724b6179bSKris Buschelman Options Database Keys: 7780bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 77924b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 78024b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 78124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 78224b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 78324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 78424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 78594b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 78624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 78724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 78824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 78924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 79024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 79124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 79224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 79324b6179bSKris Buschelman 79424b6179bSKris Buschelman Level: beginner 79524b6179bSKris Buschelman 79624b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 79724b6179bSKris Buschelman M*/ 79824b6179bSKris Buschelman 799397b6df1SKris Buschelman EXTERN_C_BEGIN 800397b6df1SKris Buschelman #undef __FUNCT__ 801f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS" 802dfbe8321SBarry Smith PetscErrorCode MatCreate_AIJMUMPS(Mat A) 803dfbe8321SBarry Smith { 804dfbe8321SBarry Smith PetscErrorCode ierr; 805dfbe8321SBarry Smith int size; 806e2d9671bSKris Buschelman Mat A_diag; 807397b6df1SKris Buschelman MPI_Comm comm; 808397b6df1SKris Buschelman 809397b6df1SKris Buschelman PetscFunctionBegin; 8105441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */ 8115441df8eSKris Buschelman /* and AIJMUMPS types */ 8125441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr); 813397b6df1SKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 814397b6df1SKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 815397b6df1SKris Buschelman if (size == 1) { 816397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 817397b6df1SKris Buschelman } else { 818397b6df1SKris Buschelman ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 819e2d9671bSKris Buschelman A_diag = ((Mat_MPIAIJ *)A->data)->A; 820e2d9671bSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr); 821397b6df1SKris Buschelman } 822f0c56d0fSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr); 823397b6df1SKris Buschelman PetscFunctionReturn(0); 824397b6df1SKris Buschelman } 825397b6df1SKris Buschelman EXTERN_C_END 826397b6df1SKris Buschelman 827f579278aSKris Buschelman #undef __FUNCT__ 828f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS" 829dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) 830dfbe8321SBarry Smith { 831dfbe8321SBarry Smith PetscErrorCode ierr; 832f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 833f579278aSKris Buschelman 834f579278aSKris Buschelman PetscFunctionBegin; 835f579278aSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 836f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 837f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 838f0c56d0fSKris Buschelman A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 839f579278aSKris Buschelman PetscFunctionReturn(0); 840f579278aSKris Buschelman } 841f579278aSKris Buschelman 842f579278aSKris Buschelman EXTERN_C_BEGIN 843f579278aSKris Buschelman #undef __FUNCT__ 8449c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS" 845dfbe8321SBarry Smith PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz) 8469c097c71SKris Buschelman { 8479c097c71SKris Buschelman Mat A; 84802217bfdSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr; 849dfbe8321SBarry Smith PetscErrorCode ierr; 8509c097c71SKris Buschelman 8519c097c71SKris Buschelman PetscFunctionBegin; 8529c097c71SKris Buschelman /* 8539c097c71SKris Buschelman After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix 8549c097c71SKris Buschelman into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS. I would 8559c097c71SKris Buschelman like this to be done in the MatCreate routine, but the creation of this inner matrix requires 8569c097c71SKris Buschelman block size info so that PETSc can determine the local size properly. The block size info is set 8579c097c71SKris Buschelman in the preallocation routine. 8589c097c71SKris Buschelman */ 8599c097c71SKris Buschelman ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz); 8609c097c71SKris Buschelman A = ((Mat_MPISBAIJ *)B->data)->A; 8619c097c71SKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr); 8629c097c71SKris Buschelman PetscFunctionReturn(0); 8639c097c71SKris Buschelman } 8649c097c71SKris Buschelman EXTERN_C_END 8659c097c71SKris Buschelman 8669c097c71SKris Buschelman EXTERN_C_BEGIN 8679c097c71SKris Buschelman #undef __FUNCT__ 868f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS" 869dfbe8321SBarry Smith PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat) 870dfbe8321SBarry Smith { 871dfbe8321SBarry Smith PetscErrorCode ierr; 872dfbe8321SBarry Smith int size; 873f579278aSKris Buschelman MPI_Comm comm; 874f579278aSKris Buschelman Mat B=*newmat; 875422e82a1SHong Zhang Mat_MUMPS *mumps; 8769c097c71SKris Buschelman void (*f)(void); 877f579278aSKris Buschelman 878f579278aSKris Buschelman PetscFunctionBegin; 879f579278aSKris Buschelman if (B != A) { 880f579278aSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 881f579278aSKris Buschelman } 882f579278aSKris Buschelman 883f579278aSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 884f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 885f579278aSKris Buschelman 886f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 887f579278aSKris Buschelman mumps->MatView = A->ops->view; 888f579278aSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 889f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 890f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 891f579278aSKris Buschelman mumps->MatDestroy = A->ops->destroy; 892a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_SBAIJMUMPS; 893f579278aSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 894f579278aSKris Buschelman mumps->isAIJ = PETSC_FALSE; 895f579278aSKris Buschelman 896f579278aSKris Buschelman B->spptr = (void*)mumps; 897422e82a1SHong Zhang B->ops->duplicate = MatDuplicate_MUMPS; 898f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 899f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS; 900f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 9013924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 902f579278aSKris Buschelman 903f579278aSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 904f579278aSKris Buschelman if (size == 1) { 905f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C", 906f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 907f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C", 908f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 909f579278aSKris Buschelman } else { 9109c097c71SKris Buschelman /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */ 9119c097c71SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr); 9129c097c71SKris Buschelman if (f) { 913*6849ba73SBarry Smith mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f; 9149c097c71SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C", 9159c097c71SKris Buschelman "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS", 9169c097c71SKris Buschelman MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr); 9179c097c71SKris Buschelman } 918f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C", 919f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 920f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C", 921f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 922f579278aSKris Buschelman } 923f579278aSKris Buschelman 924f579278aSKris Buschelman PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves."); 925f579278aSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 926f579278aSKris Buschelman *newmat = B; 927f579278aSKris Buschelman PetscFunctionReturn(0); 928f579278aSKris Buschelman } 929f579278aSKris Buschelman EXTERN_C_END 930f579278aSKris Buschelman 931f0c56d0fSKris Buschelman #undef __FUNCT__ 932422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS" 933dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) { 934dfbe8321SBarry Smith PetscErrorCode ierr; 9358e393735SKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr; 9368f340917SKris Buschelman 937f0c56d0fSKris Buschelman PetscFunctionBegin; 9388f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 9398e393735SKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr); 940f0c56d0fSKris Buschelman PetscFunctionReturn(0); 941f0c56d0fSKris Buschelman } 942f0c56d0fSKris Buschelman 94324b6179bSKris Buschelman /*MC 944fafad747SKris Buschelman MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 94524b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 94624b6179bSKris Buschelman 94724b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 94824b6179bSKris Buschelman on how to declare the existence of external packages), 94924b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 95024b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS). 95124b6179bSKris Buschelman This matrix type is only supported for double precision real. 95224b6179bSKris Buschelman 95324b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 95424b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 95524b6179bSKris Buschelman MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported 95624b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 95728b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 95828b08bd3SKris Buschelman conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 95928b08bd3SKris Buschelman without data copy. 96024b6179bSKris Buschelman 96124b6179bSKris Buschelman Options Database Keys: 9620bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 96324b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 96424b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 96524b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 96624b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 96724b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 96824b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 96994b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 97024b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 97124b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 97224b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 97324b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 97424b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 97524b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 97624b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 97724b6179bSKris Buschelman 97824b6179bSKris Buschelman Level: beginner 97924b6179bSKris Buschelman 98024b6179bSKris Buschelman .seealso: MATAIJMUMPS 98124b6179bSKris Buschelman M*/ 98224b6179bSKris Buschelman 983397b6df1SKris Buschelman EXTERN_C_BEGIN 984397b6df1SKris Buschelman #undef __FUNCT__ 985f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS" 986dfbe8321SBarry Smith PetscErrorCode MatCreate_SBAIJMUMPS(Mat A) 987dfbe8321SBarry Smith { 988*6849ba73SBarry Smith PetscErrorCode ierr; 989*6849ba73SBarry Smith int size; 990397b6df1SKris Buschelman 991397b6df1SKris Buschelman PetscFunctionBegin; 9925441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */ 9935441df8eSKris Buschelman /* and SBAIJMUMPS types */ 9945441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr); 9955441df8eSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 996397b6df1SKris Buschelman if (size == 1) { 997397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 998397b6df1SKris Buschelman } else { 999397b6df1SKris Buschelman ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr); 1000397b6df1SKris Buschelman } 1001f0c56d0fSKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr); 1002397b6df1SKris Buschelman PetscFunctionReturn(0); 1003397b6df1SKris Buschelman } 1004397b6df1SKris Buschelman EXTERN_C_END 1005