1397b6df1SKris Buschelman /*$Id: mumps.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 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; 39f0c56d0fSKris Buschelman int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*); 40c338a77dSKris Buschelman int (*MatView)(Mat,PetscViewer); 41c338a77dSKris Buschelman int (*MatAssemblyEnd)(Mat,MatAssemblyType); 42c338a77dSKris Buschelman int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*); 43c338a77dSKris Buschelman int (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*); 44c338a77dSKris Buschelman int (*MatDestroy)(Mat); 45a39386dcSKris Buschelman int (*specialdestroy)(Mat); 469c097c71SKris Buschelman int (*MatPreallocate)(Mat,int,int,int*,int,int*); 47f0c56d0fSKris Buschelman } Mat_MUMPS; 48f0c56d0fSKris Buschelman 49f0c56d0fSKris Buschelman EXTERN int MatDuplicate_AIJMUMPS(Mat,MatDuplicateOption,Mat*); 50f0c56d0fSKris Buschelman EXTERN int MatDuplicate_SBAIJMUMPS(Mat,MatDuplicateOption,Mat*); 51*892f6c3fSKris Buschelman EXTERN_C_BEGIN 52*892f6c3fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,Mat*); 53*892f6c3fSKris Buschelman EXTERN_C_END 54*892f6c3fSKris Buschelman 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 */ 66f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) { 67397b6df1SKris Buschelman int *ai, *aj, *bi, *bj, rstart,nz, *garray; 68397b6df1SKris Buschelman int ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol; 69d54de34fSKris Buschelman int *row,*col; 70397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 71f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 72397b6df1SKris Buschelman 73397b6df1SKris Buschelman PetscFunctionBegin; 74397b6df1SKris Buschelman if (mumps->isAIJ){ 75397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 76397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 77397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 78397b6df1SKris Buschelman nz = aa->nz + bb->nz; 79397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 80397b6df1SKris Buschelman garray = mat->garray; 81397b6df1SKris Buschelman av=aa->a; bv=bb->a; 82397b6df1SKris Buschelman 83397b6df1SKris Buschelman } else { 84397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 85397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 86397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 87847143adSKris Buschelman if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs); 886c6c5352SBarry Smith nz = aa->nz + bb->nz; 89397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 90397b6df1SKris Buschelman garray = mat->garray; 91397b6df1SKris Buschelman av=aa->a; bv=bb->a; 92397b6df1SKris Buschelman } 93397b6df1SKris Buschelman 94397b6df1SKris Buschelman if (!valOnly){ 95397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr); 96397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr); 97397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 98397b6df1SKris Buschelman *r = row; *c = col; *v = val; 99397b6df1SKris Buschelman } else { 100397b6df1SKris Buschelman row = *r; col = *c; val = *v; 101397b6df1SKris Buschelman } 102397b6df1SKris Buschelman *nnz = nz; 103397b6df1SKris Buschelman 104028e57e8SHong Zhang jj = 0; irow = rstart; 105397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 106397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 107397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 108397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 109397b6df1SKris Buschelman bjj = bj + bi[i]; 110397b6df1SKris Buschelman 111397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 112397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 11375747be1SHong Zhang j=-1; 11475747be1SHong Zhang do { 11575747be1SHong Zhang j++; 11675747be1SHong Zhang if (j == countB) break; 117397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11875747be1SHong Zhang } while (jcol < colA_start); 11975747be1SHong Zhang jB = j; 120397b6df1SKris Buschelman 121397b6df1SKris Buschelman /* B-part, smaller col index */ 122397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 123397b6df1SKris Buschelman for (j=0; j<jB; j++){ 124397b6df1SKris Buschelman jcol = garray[bjj[j]]; 125397b6df1SKris Buschelman if (!valOnly){ 126397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 12775747be1SHong Zhang 128397b6df1SKris Buschelman } 129397b6df1SKris Buschelman val[jj++] = *bv++; 130397b6df1SKris Buschelman } 131397b6df1SKris Buschelman /* A-part */ 132397b6df1SKris Buschelman for (j=0; j<countA; j++){ 133397b6df1SKris Buschelman if (!valOnly){ 134397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 135397b6df1SKris Buschelman } 136397b6df1SKris Buschelman val[jj++] = *av++; 137397b6df1SKris Buschelman } 138397b6df1SKris Buschelman /* B-part, larger col index */ 139397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 140397b6df1SKris Buschelman if (!valOnly){ 141397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 142397b6df1SKris Buschelman } 143397b6df1SKris Buschelman val[jj++] = *bv++; 144397b6df1SKris Buschelman } 145397b6df1SKris Buschelman irow++; 146397b6df1SKris Buschelman } 147397b6df1SKris Buschelman 148397b6df1SKris Buschelman PetscFunctionReturn(0); 149397b6df1SKris Buschelman } 150397b6df1SKris Buschelman 151c338a77dSKris Buschelman EXTERN_C_BEGIN 152c338a77dSKris Buschelman #undef __FUNCT__ 153c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base" 1548e9aea5cSBarry Smith int MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) { 155c338a77dSKris Buschelman int ierr; 156c338a77dSKris Buschelman Mat B=*newmat; 157f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 158c338a77dSKris Buschelman 159c338a77dSKris Buschelman PetscFunctionBegin; 160c338a77dSKris Buschelman if (B != A) { 161c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 162c338a77dSKris Buschelman } 163f0c56d0fSKris Buschelman B->ops->duplicate = mumps->MatDuplicate; 164f0c56d0fSKris Buschelman B->ops->view = mumps->MatView; 165f0c56d0fSKris Buschelman B->ops->assemblyend = mumps->MatAssemblyEnd; 166f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic; 167f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic; 168f0c56d0fSKris Buschelman B->ops->destroy = mumps->MatDestroy; 1693924e44cSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr); 170f0c56d0fSKris Buschelman ierr = PetscFree(mumps);CHKERRQ(ierr); 171c338a77dSKris Buschelman *newmat = B; 172c338a77dSKris Buschelman PetscFunctionReturn(0); 173c338a77dSKris Buschelman } 174c338a77dSKris Buschelman EXTERN_C_END 175c338a77dSKris Buschelman 176397b6df1SKris Buschelman #undef __FUNCT__ 1773924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 1783924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) { 179f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 180c338a77dSKris Buschelman int ierr,size=lu->size; 181a39386dcSKris Buschelman int (*specialdestroy)(Mat); 182397b6df1SKris Buschelman PetscFunctionBegin; 183397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 184397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 185397b6df1SKris Buschelman lu->id.job=JOB_END; 186397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 187397b6df1SKris Buschelman zmumps_c(&lu->id); 188397b6df1SKris Buschelman #else 189397b6df1SKris Buschelman dmumps_c(&lu->id); 190397b6df1SKris Buschelman #endif 191c338a77dSKris Buschelman if (lu->irn) { 192c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 193c338a77dSKris Buschelman } 194c338a77dSKris Buschelman if (lu->jcn) { 195c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 196c338a77dSKris Buschelman } 197c338a77dSKris Buschelman if (size>1 && lu->val) { 198c338a77dSKris Buschelman ierr = PetscFree(lu->val);CHKERRQ(ierr); 199c338a77dSKris Buschelman } 200397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 201397b6df1SKris Buschelman } 202a39386dcSKris Buschelman specialdestroy = lu->specialdestroy; 203a39386dcSKris Buschelman ierr = (*specialdestroy)(A);CHKERRQ(ierr); 204c338a77dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 205397b6df1SKris Buschelman PetscFunctionReturn(0); 206397b6df1SKris Buschelman } 207397b6df1SKris Buschelman 208397b6df1SKris Buschelman #undef __FUNCT__ 209a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS" 210a39386dcSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) { 211a39386dcSKris Buschelman int ierr, size; 212a39386dcSKris Buschelman 213a39386dcSKris Buschelman PetscFunctionBegin; 214a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 215a39386dcSKris Buschelman if (size==1) { 216a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr); 217a39386dcSKris Buschelman } else { 218a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr); 219a39386dcSKris Buschelman } 220a39386dcSKris Buschelman PetscFunctionReturn(0); 221a39386dcSKris Buschelman } 222a39386dcSKris Buschelman 223a39386dcSKris Buschelman #undef __FUNCT__ 224a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS" 225a39386dcSKris Buschelman int MatDestroy_SBAIJMUMPS(Mat A) { 226a39386dcSKris Buschelman int ierr, size; 227a39386dcSKris Buschelman 228a39386dcSKris Buschelman PetscFunctionBegin; 229a39386dcSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 230a39386dcSKris Buschelman if (size==1) { 231a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr); 232a39386dcSKris Buschelman } else { 233a39386dcSKris Buschelman ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr); 234a39386dcSKris Buschelman } 235a39386dcSKris Buschelman PetscFunctionReturn(0); 236a39386dcSKris Buschelman } 237a39386dcSKris Buschelman 238a39386dcSKris Buschelman #undef __FUNCT__ 239c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS" 240f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 241f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 242397b6df1SKris Buschelman int ierr; 243397b6df1SKris Buschelman 244397b6df1SKris Buschelman PetscFunctionBegin; 245c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 246c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 247c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 248c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 249c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 250c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 251c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 252c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 253c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 254c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 255c338a77dSKris Buschelman if (lu->myid == 0 && lu->id.ICNTL(11)>0) { 256c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 257c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 258c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 259c338a77dSKris 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); 260c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 261c338a77dSKris 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); 262c338a77dSKris Buschelman 263c338a77dSKris Buschelman } 264c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 265c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 266adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 267c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(15) (efficiency control): %d \n",lu->id.ICNTL(15));CHKERRQ(ierr); 268c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 269c338a77dSKris Buschelman 270c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 271c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 272c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 27357f0c58bSHong Zhang 27457f0c58bSHong Zhang /* infomation local to each processor */ 27557f0c58bSHong Zhang if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr); 27657f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 27757f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 27857f0c58bSHong Zhang if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr); 27957f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 28057f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 28157f0c58bSHong Zhang if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr); 28257f0c58bSHong Zhang ierr = PetscSynchronizedPrintf(A->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 28357f0c58bSHong Zhang ierr = PetscSynchronizedFlush(A->comm); 284adc1d99fSHong Zhang 285adc1d99fSHong Zhang if (lu->myid == 0){ /* information from the host */ 286adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 287adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 288adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 289adc1d99fSHong Zhang 290adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 291adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 292adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 293adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 294adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 295adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 296adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 297adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 298adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 299adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 300adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 301adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 302adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 303adc1d99fSHong 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); 304adc1d99fSHong 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); 305adc1d99fSHong 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); 306adc1d99fSHong 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); 307adc1d99fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 308adc1d99fSHong Zhang } 309adc1d99fSHong Zhang 310397b6df1SKris Buschelman PetscFunctionReturn(0); 311397b6df1SKris Buschelman } 312397b6df1SKris Buschelman 313397b6df1SKris Buschelman #undef __FUNCT__ 314f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS" 315f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) { 316397b6df1SKris Buschelman int ierr; 317397b6df1SKris Buschelman PetscTruth isascii; 318397b6df1SKris Buschelman PetscViewerFormat format; 319f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 320397b6df1SKris Buschelman 321397b6df1SKris Buschelman PetscFunctionBegin; 322397b6df1SKris Buschelman ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr); 323397b6df1SKris Buschelman 324397b6df1SKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 325397b6df1SKris Buschelman if (isascii) { 326397b6df1SKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 327397b6df1SKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 328397b6df1SKris Buschelman ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 329397b6df1SKris Buschelman } 330397b6df1SKris Buschelman } 331397b6df1SKris Buschelman PetscFunctionReturn(0); 332397b6df1SKris Buschelman } 333397b6df1SKris Buschelman 334397b6df1SKris Buschelman #undef __FUNCT__ 335f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS" 336f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) { 337f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 338d54de34fSKris Buschelman PetscScalar *array; 339397b6df1SKris Buschelman Vec x_seq; 340397b6df1SKris Buschelman IS iden; 341397b6df1SKris Buschelman VecScatter scat; 342397b6df1SKris Buschelman int ierr; 343397b6df1SKris Buschelman 344397b6df1SKris Buschelman PetscFunctionBegin; 345397b6df1SKris Buschelman if (lu->size > 1){ 346397b6df1SKris Buschelman if (!lu->myid){ 347397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr); 348397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr); 349397b6df1SKris Buschelman } else { 350397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr); 351397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr); 352397b6df1SKris Buschelman } 353397b6df1SKris Buschelman ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr); 354397b6df1SKris Buschelman ierr = ISDestroy(iden);CHKERRQ(ierr); 355397b6df1SKris Buschelman 356397b6df1SKris Buschelman ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 357397b6df1SKris Buschelman ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 358397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 359397b6df1SKris Buschelman } else { /* size == 1 */ 360397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 361397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 362397b6df1SKris Buschelman } 363397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 364397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 365397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 366397b6df1SKris Buschelman #else 367397b6df1SKris Buschelman lu->id.rhs = array; 368397b6df1SKris Buschelman #endif 369397b6df1SKris Buschelman } 370397b6df1SKris Buschelman 371397b6df1SKris Buschelman /* solve phase */ 372397b6df1SKris Buschelman lu->id.job=3; 373397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 374397b6df1SKris Buschelman zmumps_c(&lu->id); 375397b6df1SKris Buschelman #else 376397b6df1SKris Buschelman dmumps_c(&lu->id); 377397b6df1SKris Buschelman #endif 378397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 379397b6df1SKris Buschelman SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 380397b6df1SKris Buschelman } 381397b6df1SKris Buschelman 382397b6df1SKris Buschelman /* convert mumps solution x_seq to petsc mpi x */ 383397b6df1SKris Buschelman if (lu->size > 1) { 384397b6df1SKris Buschelman if (!lu->myid){ 385397b6df1SKris Buschelman ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 386397b6df1SKris Buschelman } 387397b6df1SKris Buschelman ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 388397b6df1SKris Buschelman ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 389397b6df1SKris Buschelman ierr = VecScatterDestroy(scat);CHKERRQ(ierr); 390397b6df1SKris Buschelman ierr = VecDestroy(x_seq);CHKERRQ(ierr); 391397b6df1SKris Buschelman } else { 392397b6df1SKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 393397b6df1SKris Buschelman } 394397b6df1SKris Buschelman 395397b6df1SKris Buschelman PetscFunctionReturn(0); 396397b6df1SKris Buschelman } 397397b6df1SKris Buschelman 398a58c3f20SHong Zhang /* 399a58c3f20SHong Zhang input: 400a58c3f20SHong Zhang F: numeric factor 401a58c3f20SHong Zhang output: 402a58c3f20SHong Zhang nneg: total number of negative pivots 403a58c3f20SHong Zhang nzero: 0 404a58c3f20SHong Zhang npos: (global dimension of F) - nneg 405a58c3f20SHong Zhang */ 406a58c3f20SHong Zhang 407a58c3f20SHong Zhang #undef __FUNCT__ 408a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 409a58c3f20SHong Zhang int MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 410a58c3f20SHong Zhang { 411a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 412bcb30aebSHong Zhang int ierr,neg,zero,pos,size; 413a58c3f20SHong Zhang 414a58c3f20SHong Zhang PetscFunctionBegin; 415bcb30aebSHong Zhang ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr); 416bcb30aebSHong 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 */ 417bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 418bcb30aebSHong 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)); 419bcb30aebSHong Zhang } 420a58c3f20SHong Zhang if (nneg){ 421a58c3f20SHong Zhang if (!lu->myid){ 422a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 423a58c3f20SHong Zhang } 424bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 425a58c3f20SHong Zhang } 426a58c3f20SHong Zhang if (nzero) *nzero = 0; 427a58c3f20SHong Zhang if (npos) *npos = F->M - (*nneg); 428a58c3f20SHong Zhang PetscFunctionReturn(0); 429a58c3f20SHong Zhang } 430a58c3f20SHong Zhang 431397b6df1SKris Buschelman #undef __FUNCT__ 432f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS" 433f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) { 434f0c56d0fSKris Buschelman Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 435f0c56d0fSKris Buschelman Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr; 436397b6df1SKris Buschelman int rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl; 437397b6df1SKris Buschelman PetscTruth valOnly,flg; 438397b6df1SKris Buschelman 439397b6df1SKris Buschelman PetscFunctionBegin; 440397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 441f0c56d0fSKris Buschelman (*F)->ops->solve = MatSolve_AIJMUMPS; 442397b6df1SKris Buschelman 443397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 444397b6df1SKris Buschelman ierr = MPI_Comm_rank(A->comm, &lu->myid); 445397b6df1SKris Buschelman ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr); 44675747be1SHong Zhang lua->myid = lu->myid; lua->size = lu->size; 447397b6df1SKris Buschelman lu->id.job = JOB_INIT; 448397b6df1SKris Buschelman ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr); 449397b6df1SKris Buschelman lu->id.comm_fortran = lu->comm_mumps; 450397b6df1SKris Buschelman 451397b6df1SKris Buschelman /* Set mumps options */ 452397b6df1SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 453397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 454397b6df1SKris Buschelman lu->id.sym=lu->sym; 455397b6df1SKris Buschelman if (lu->sym == 2){ 456397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 457397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 458397b6df1SKris Buschelman } 459397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 460397b6df1SKris Buschelman zmumps_c(&lu->id); 461397b6df1SKris Buschelman #else 462397b6df1SKris Buschelman dmumps_c(&lu->id); 463397b6df1SKris Buschelman #endif 464397b6df1SKris Buschelman 465397b6df1SKris Buschelman if (lu->size == 1){ 466397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 467397b6df1SKris Buschelman } else { 468397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 469397b6df1SKris Buschelman } 470397b6df1SKris Buschelman 471397b6df1SKris Buschelman icntl=-1; 472397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 473397b6df1SKris Buschelman if (flg && icntl > 0) { 474397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 475397b6df1SKris Buschelman } else { /* no output */ 476397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 477397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 478397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 479397b6df1SKris Buschelman lu->id.ICNTL(4) = 0; /* level of printing, 0,1,2,3,4, default=2 */ 480397b6df1SKris Buschelman } 481397b6df1SKris 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); 482397b6df1SKris Buschelman icntl=-1; 483397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 484397b6df1SKris Buschelman if (flg) { 485397b6df1SKris Buschelman if (icntl== 1){ 486397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 487397b6df1SKris Buschelman } else { 488397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 489397b6df1SKris Buschelman } 490397b6df1SKris Buschelman } 491397b6df1SKris 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); 492397b6df1SKris 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); 49394b7f48cSBarry 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); 494397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 495397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 496adc1d99fSHong 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); 497397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 498397b6df1SKris Buschelman 499397b6df1SKris Buschelman /* 500397b6df1SKris 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); 501397b6df1SKris Buschelman if (flg){ 502397b6df1SKris Buschelman if (icntl >-1 && icntl <3 ){ 503397b6df1SKris Buschelman if (lu->myid==0) lu->id.ICNTL(16) = icntl; 504397b6df1SKris Buschelman } else { 505397b6df1SKris Buschelman SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl); 506397b6df1SKris Buschelman } 507397b6df1SKris Buschelman } 508397b6df1SKris Buschelman */ 509397b6df1SKris Buschelman 510397b6df1SKris 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); 511397b6df1SKris 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); 512397b6df1SKris 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); 513397b6df1SKris Buschelman PetscOptionsEnd(); 514397b6df1SKris Buschelman } 515397b6df1SKris Buschelman 516397b6df1SKris Buschelman /* define matrix A */ 517397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 518397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 519397b6df1SKris Buschelman if (!lu->myid) { 520c36ead0aSKris Buschelman if (lua->isAIJ){ 521397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 522397b6df1SKris Buschelman nz = aa->nz; 523397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 524397b6df1SKris Buschelman } else { 525397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 5266c6c5352SBarry Smith nz = aa->nz; 527397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 528397b6df1SKris Buschelman } 529397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 530397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr); 531397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr); 532397b6df1SKris Buschelman nz = 0; 533397b6df1SKris Buschelman for (i=0; i<M; i++){ 534397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 535397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 536397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 537397b6df1SKris Buschelman } 538397b6df1SKris Buschelman } 539397b6df1SKris Buschelman } 540397b6df1SKris Buschelman } 541397b6df1SKris Buschelman break; 542397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 543397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 544397b6df1SKris Buschelman valOnly = PETSC_FALSE; 545397b6df1SKris Buschelman } else { 546397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 547397b6df1SKris Buschelman } 548397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 549397b6df1SKris Buschelman break; 550397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 551397b6df1SKris Buschelman } 552397b6df1SKris Buschelman 553397b6df1SKris Buschelman /* analysis phase */ 554397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 555397b6df1SKris Buschelman lu->id.n = M; 556397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 557397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 558397b6df1SKris Buschelman if (!lu->myid) { 559397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 560397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 561397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 562397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 563397b6df1SKris Buschelman #else 564397b6df1SKris Buschelman lu->id.a = lu->val; 565397b6df1SKris Buschelman #endif 566397b6df1SKris Buschelman } 567397b6df1SKris Buschelman } 568397b6df1SKris Buschelman break; 569397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 570397b6df1SKris Buschelman lu->id.nz_loc = nnz; 571397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 572397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 573397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 574397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 575397b6df1SKris Buschelman #else 576397b6df1SKris Buschelman lu->id.a_loc = lu->val; 577397b6df1SKris Buschelman #endif 578397b6df1SKris Buschelman } 579397b6df1SKris Buschelman break; 580397b6df1SKris Buschelman } 581397b6df1SKris Buschelman lu->id.job=1; 582397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 583397b6df1SKris Buschelman zmumps_c(&lu->id); 584397b6df1SKris Buschelman #else 585397b6df1SKris Buschelman dmumps_c(&lu->id); 586397b6df1SKris Buschelman #endif 587397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 588397b6df1SKris Buschelman SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 589397b6df1SKris Buschelman } 590397b6df1SKris Buschelman } 591397b6df1SKris Buschelman 592397b6df1SKris Buschelman /* numerical factorization phase */ 593397b6df1SKris Buschelman if(lu->id.ICNTL(18) == 0) { 594a7aca84bSHong Zhang if (!lu->myid) { 595397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 596397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 597397b6df1SKris Buschelman #else 598397b6df1SKris Buschelman lu->id.a = lu->val; 599397b6df1SKris Buschelman #endif 600397b6df1SKris Buschelman } 601397b6df1SKris Buschelman } else { 602397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 603397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 604397b6df1SKris Buschelman #else 605397b6df1SKris Buschelman lu->id.a_loc = lu->val; 606397b6df1SKris Buschelman #endif 607397b6df1SKris Buschelman } 608397b6df1SKris Buschelman lu->id.job=2; 609397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 610397b6df1SKris Buschelman zmumps_c(&lu->id); 611397b6df1SKris Buschelman #else 612397b6df1SKris Buschelman dmumps_c(&lu->id); 613397b6df1SKris Buschelman #endif 614397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 615a7aca84bSHong 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)); 616397b6df1SKris Buschelman } 617397b6df1SKris Buschelman 618397b6df1SKris Buschelman if (lu->myid==0 && lu->id.ICNTL(16) > 0){ 619397b6df1SKris Buschelman SETERRQ1(1," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 620397b6df1SKris Buschelman } 621397b6df1SKris Buschelman 622397b6df1SKris Buschelman (*F)->assembled = PETSC_TRUE; 623397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 624ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 625397b6df1SKris Buschelman PetscFunctionReturn(0); 626397b6df1SKris Buschelman } 627397b6df1SKris Buschelman 628397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 629397b6df1SKris Buschelman #undef __FUNCT__ 630f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 631f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 632397b6df1SKris Buschelman Mat B; 633f0c56d0fSKris Buschelman Mat_MUMPS *lu; 634397b6df1SKris Buschelman int ierr; 635397b6df1SKris Buschelman 636397b6df1SKris Buschelman PetscFunctionBegin; 637397b6df1SKris Buschelman 638397b6df1SKris Buschelman /* Create the factorization matrix */ 639397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 640397b6df1SKris Buschelman ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr); 641397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 642397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 643397b6df1SKris Buschelman 644f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS; 645397b6df1SKris Buschelman B->factor = FACTOR_LU; 646f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 647397b6df1SKris Buschelman lu->sym = 0; 648397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 649397b6df1SKris Buschelman 650397b6df1SKris Buschelman *F = B; 651397b6df1SKris Buschelman PetscFunctionReturn(0); 652397b6df1SKris Buschelman } 653397b6df1SKris Buschelman 654397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 655397b6df1SKris Buschelman #undef __FUNCT__ 656f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 657f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) { 658397b6df1SKris Buschelman Mat B; 659f0c56d0fSKris Buschelman Mat_MUMPS *lu; 660397b6df1SKris Buschelman int ierr; 661397b6df1SKris Buschelman 662397b6df1SKris Buschelman PetscFunctionBegin; 663397b6df1SKris Buschelman 664397b6df1SKris Buschelman /* Create the factorization matrix */ 665397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 666397b6df1SKris Buschelman ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr); 667397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 668397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 669397b6df1SKris Buschelman 670f0c56d0fSKris Buschelman B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS; 671a58c3f20SHong Zhang B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 672397b6df1SKris Buschelman B->factor = FACTOR_CHOLESKY; 673f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 674397b6df1SKris Buschelman lu->sym = 2; 675397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 676397b6df1SKris Buschelman 677397b6df1SKris Buschelman *F = B; 678397b6df1SKris Buschelman PetscFunctionReturn(0); 679397b6df1SKris Buschelman } 680397b6df1SKris Buschelman 681397b6df1SKris Buschelman #undef __FUNCT__ 682f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS" 683f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) { 684c338a77dSKris Buschelman int ierr; 685f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 686c338a77dSKris Buschelman 687397b6df1SKris Buschelman PetscFunctionBegin; 688c338a77dSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 689f0c56d0fSKris Buschelman 690c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 691c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 692f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 693397b6df1SKris Buschelman PetscFunctionReturn(0); 694397b6df1SKris Buschelman } 695397b6df1SKris Buschelman 696c338a77dSKris Buschelman EXTERN_C_BEGIN 697c338a77dSKris Buschelman #undef __FUNCT__ 698f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS" 6998e9aea5cSBarry Smith int MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat) { 700c338a77dSKris Buschelman int ierr,size; 701c338a77dSKris Buschelman MPI_Comm comm; 702c338a77dSKris Buschelman Mat B=*newmat; 703f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 704397b6df1SKris Buschelman 705397b6df1SKris Buschelman PetscFunctionBegin; 706c338a77dSKris Buschelman if (B != A) { 707c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 708397b6df1SKris Buschelman } 709397b6df1SKris Buschelman 710c338a77dSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 711f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 712c338a77dSKris Buschelman 713f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 714c338a77dSKris Buschelman mumps->MatView = A->ops->view; 715c338a77dSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 716c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 717c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 718c338a77dSKris Buschelman mumps->MatDestroy = A->ops->destroy; 719a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_AIJMUMPS; 720c338a77dSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 721f579278aSKris Buschelman mumps->isAIJ = PETSC_TRUE; 722c338a77dSKris Buschelman 7234b68dd72SKris Buschelman B->spptr = (void *)mumps; 724f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_AIJMUMPS; 725f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 726f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS; 727f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 7283924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 729c338a77dSKris Buschelman 730c338a77dSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 731c338a77dSKris Buschelman if (size == 1) { 732c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C", 733f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 734c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C", 735c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 736c338a77dSKris Buschelman } else { 737c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C", 738f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 739c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C", 740c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 741c338a77dSKris Buschelman } 742c338a77dSKris Buschelman 743f579278aSKris Buschelman PetscLogInfo(0,"Using MUMPS for LU factorization and solves."); 744c338a77dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 745c338a77dSKris Buschelman *newmat = B; 746397b6df1SKris Buschelman PetscFunctionReturn(0); 747397b6df1SKris Buschelman } 748c338a77dSKris Buschelman EXTERN_C_END 749397b6df1SKris Buschelman 750f0c56d0fSKris Buschelman #undef __FUNCT__ 751f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS" 752f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) { 753f0c56d0fSKris Buschelman int ierr; 7548f340917SKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr; 7558f340917SKris Buschelman 756f0c56d0fSKris Buschelman PetscFunctionBegin; 7578f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 758f0c56d0fSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr); 759a39386dcSKris Buschelman ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr); 760f0c56d0fSKris Buschelman PetscFunctionReturn(0); 761f0c56d0fSKris Buschelman } 762f0c56d0fSKris Buschelman 76324b6179bSKris Buschelman /*MC 764fafad747SKris Buschelman MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 76524b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 76624b6179bSKris Buschelman 76724b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 76824b6179bSKris Buschelman on how to declare the existence of external packages), 76924b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 77024b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS). 77124b6179bSKris Buschelman This matrix type is only supported for double precision real. 77224b6179bSKris Buschelman 77324b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 77424b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 77524b6179bSKris Buschelman MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported 77624b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 77728b08bd3SKris Buschelman the above preallocation routines for simplicity. One can also call MatConvert for an inplace 77828b08bd3SKris Buschelman conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 77928b08bd3SKris Buschelman without data copy. 78024b6179bSKris Buschelman 78124b6179bSKris Buschelman Options Database Keys: 7820bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 78324b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 78424b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 78524b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 78624b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 78724b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 78824b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 78994b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 79024b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 79124b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 79224b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 79324b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 79424b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 79524b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 79624b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 79724b6179bSKris Buschelman 79824b6179bSKris Buschelman Level: beginner 79924b6179bSKris Buschelman 80024b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 80124b6179bSKris Buschelman M*/ 80224b6179bSKris Buschelman 803397b6df1SKris Buschelman EXTERN_C_BEGIN 804397b6df1SKris Buschelman #undef __FUNCT__ 805f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS" 806f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) { 807397b6df1SKris Buschelman int ierr,size; 808e2d9671bSKris Buschelman Mat A_diag; 809397b6df1SKris Buschelman MPI_Comm comm; 810397b6df1SKris Buschelman 811397b6df1SKris Buschelman PetscFunctionBegin; 8125441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */ 8135441df8eSKris Buschelman /* and AIJMUMPS types */ 8145441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr); 815397b6df1SKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 816397b6df1SKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 817397b6df1SKris Buschelman if (size == 1) { 818397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 819397b6df1SKris Buschelman } else { 820397b6df1SKris Buschelman ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 821e2d9671bSKris Buschelman A_diag = ((Mat_MPIAIJ *)A->data)->A; 822e2d9671bSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr); 823397b6df1SKris Buschelman } 824f0c56d0fSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr); 825397b6df1SKris Buschelman PetscFunctionReturn(0); 826397b6df1SKris Buschelman } 827397b6df1SKris Buschelman EXTERN_C_END 828397b6df1SKris Buschelman 829f579278aSKris Buschelman #undef __FUNCT__ 830f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS" 831f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) { 832f579278aSKris Buschelman int ierr; 833f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 834f579278aSKris Buschelman 835f579278aSKris Buschelman PetscFunctionBegin; 836f579278aSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 837f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 838f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 839f0c56d0fSKris Buschelman A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 840f579278aSKris Buschelman PetscFunctionReturn(0); 841f579278aSKris Buschelman } 842f579278aSKris Buschelman 843f579278aSKris Buschelman EXTERN_C_BEGIN 844f579278aSKris Buschelman #undef __FUNCT__ 8459c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS" 8469c097c71SKris Buschelman int MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz) 8479c097c71SKris Buschelman { 8489c097c71SKris Buschelman Mat A; 8499c097c71SKris Buschelman Mat_MUMPS *mumps; 8509c097c71SKris Buschelman int ierr; 8519c097c71SKris Buschelman 8529c097c71SKris Buschelman PetscFunctionBegin; 8539c097c71SKris Buschelman /* 8549c097c71SKris Buschelman After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix 8559c097c71SKris Buschelman into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS. I would 8569c097c71SKris Buschelman like this to be done in the MatCreate routine, but the creation of this inner matrix requires 8579c097c71SKris Buschelman block size info so that PETSc can determine the local size properly. The block size info is set 8589c097c71SKris Buschelman in the preallocation routine. 8599c097c71SKris Buschelman */ 8609c097c71SKris Buschelman ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz); 8619c097c71SKris Buschelman A = ((Mat_MPISBAIJ *)B->data)->A; 8629c097c71SKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr); 8639c097c71SKris Buschelman PetscFunctionReturn(0); 8649c097c71SKris Buschelman } 8659c097c71SKris Buschelman EXTERN_C_END 8669c097c71SKris Buschelman 8679c097c71SKris Buschelman EXTERN_C_BEGIN 8689c097c71SKris Buschelman #undef __FUNCT__ 869f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS" 8708e9aea5cSBarry Smith int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat) { 871f579278aSKris Buschelman int ierr,size; 872f579278aSKris Buschelman MPI_Comm comm; 873f579278aSKris Buschelman Mat B=*newmat; 874f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 8759c097c71SKris Buschelman void (*f)(void); 876f579278aSKris Buschelman 877f579278aSKris Buschelman PetscFunctionBegin; 878f579278aSKris Buschelman if (B != A) { 879f579278aSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 880f579278aSKris Buschelman } 881f579278aSKris Buschelman 882f579278aSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 883f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 884f579278aSKris Buschelman 885f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 886f579278aSKris Buschelman mumps->MatView = A->ops->view; 887f579278aSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 888f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 889f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 890f579278aSKris Buschelman mumps->MatDestroy = A->ops->destroy; 891a39386dcSKris Buschelman mumps->specialdestroy = MatDestroy_SBAIJMUMPS; 892f579278aSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 893f579278aSKris Buschelman mumps->isAIJ = PETSC_FALSE; 894f579278aSKris Buschelman 895f579278aSKris Buschelman B->spptr = (void *)mumps; 896f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SBAIJMUMPS; 897f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 898f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS; 899f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 9003924e44cSKris Buschelman B->ops->destroy = MatDestroy_MUMPS; 901f579278aSKris Buschelman 902f579278aSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 903f579278aSKris Buschelman if (size == 1) { 904f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C", 905f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 906f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C", 907f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 908f579278aSKris Buschelman } else { 9099c097c71SKris Buschelman /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */ 9109c097c71SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr); 9119c097c71SKris Buschelman if (f) { 9129c097c71SKris Buschelman mumps->MatPreallocate = (int (*)(Mat,int,int,int*,int,int*))f; 9139c097c71SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C", 9149c097c71SKris Buschelman "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS", 9159c097c71SKris Buschelman MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr); 9169c097c71SKris Buschelman } 917f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C", 918f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 919f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C", 920f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 921f579278aSKris Buschelman } 922f579278aSKris Buschelman 923f579278aSKris Buschelman PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves."); 924f579278aSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 925f579278aSKris Buschelman *newmat = B; 926f579278aSKris Buschelman PetscFunctionReturn(0); 927f579278aSKris Buschelman } 928f579278aSKris Buschelman EXTERN_C_END 929f579278aSKris Buschelman 930f0c56d0fSKris Buschelman #undef __FUNCT__ 931f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS" 932f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) { 933f0c56d0fSKris Buschelman int ierr; 9348f340917SKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr; 9358f340917SKris Buschelman 936f0c56d0fSKris Buschelman PetscFunctionBegin; 9378f340917SKris Buschelman ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr); 938f0c56d0fSKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr); 9393f953163SKris 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" 986f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) { 987397b6df1SKris Buschelman int ierr,size; 988397b6df1SKris Buschelman 989397b6df1SKris Buschelman PetscFunctionBegin; 9905441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */ 9915441df8eSKris Buschelman /* and SBAIJMUMPS types */ 9925441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr); 9935441df8eSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 994397b6df1SKris Buschelman if (size == 1) { 995397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 996397b6df1SKris Buschelman } else { 997397b6df1SKris Buschelman ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr); 998397b6df1SKris Buschelman } 999f0c56d0fSKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr); 1000397b6df1SKris Buschelman PetscFunctionReturn(0); 1001397b6df1SKris Buschelman } 1002397b6df1SKris Buschelman EXTERN_C_END 1003