1397b6df1SKris Buschelman /*$Id: mumps.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/ 2397b6df1SKris Buschelman /* 3397b6df1SKris Buschelman Provides an interface to the MUMPS_4.2_beta sparse solver 4397b6df1SKris Buschelman */ 5397b6df1SKris Buschelman 6397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 7397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h" 8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h" 9397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h" 10397b6df1SKris Buschelman 11397b6df1SKris Buschelman EXTERN_C_BEGIN 12397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 13397b6df1SKris Buschelman #include "zmumps_c.h" 14397b6df1SKris Buschelman #else 15397b6df1SKris Buschelman #include "dmumps_c.h" 16397b6df1SKris Buschelman #endif 17397b6df1SKris Buschelman EXTERN_C_END 18397b6df1SKris Buschelman #define JOB_INIT -1 19397b6df1SKris Buschelman #define JOB_END -2 20397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 21397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 22397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 23397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 24397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 25397b6df1SKris Buschelman 26397b6df1SKris Buschelman typedef struct { 27397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 28397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 29397b6df1SKris Buschelman #else 30397b6df1SKris Buschelman DMUMPS_STRUC_C id; 31397b6df1SKris Buschelman #endif 32397b6df1SKris Buschelman MatStructure matstruc; 33397b6df1SKris Buschelman int myid,size,*irn,*jcn,sym; 34397b6df1SKris Buschelman PetscScalar *val; 35397b6df1SKris Buschelman MPI_Comm comm_mumps; 36397b6df1SKris Buschelman 37c338a77dSKris Buschelman MatType basetype; 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); 45f0c56d0fSKris Buschelman } Mat_MUMPS; 46f0c56d0fSKris Buschelman 47f0c56d0fSKris Buschelman EXTERN int MatDuplicate_AIJMUMPS(Mat,MatDuplicateOption,Mat*); 48f0c56d0fSKris Buschelman EXTERN int MatDuplicate_SBAIJMUMPS(Mat,MatDuplicateOption,Mat*); 490e3434eeSKris Buschelman 50397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 51397b6df1SKris Buschelman /* 52397b6df1SKris Buschelman input: 53397b6df1SKris Buschelman A - matrix in mpiaij format 54397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 55397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 56397b6df1SKris Buschelman TRUE: only the values in v array are updated 57397b6df1SKris Buschelman output: 58397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 59397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 60397b6df1SKris Buschelman */ 61f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) { 62397b6df1SKris Buschelman int *ai, *aj, *bi, *bj, rstart,nz, *garray; 63397b6df1SKris Buschelman int ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol; 64d54de34fSKris Buschelman int *row,*col; 65397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 66f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 67397b6df1SKris Buschelman 68397b6df1SKris Buschelman PetscFunctionBegin; 69397b6df1SKris Buschelman 70397b6df1SKris Buschelman if (mumps->isAIJ){ 71397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 72397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 73397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 74397b6df1SKris Buschelman nz = aa->nz + bb->nz; 75397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 76397b6df1SKris Buschelman garray = mat->garray; 77397b6df1SKris Buschelman av=aa->a; bv=bb->a; 78397b6df1SKris Buschelman 79397b6df1SKris Buschelman } else { 80397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 81397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 82397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 83847143adSKris Buschelman if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs); 84397b6df1SKris Buschelman nz = aa->s_nz + bb->nz; 85397b6df1SKris Buschelman ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart; 86397b6df1SKris Buschelman garray = mat->garray; 87397b6df1SKris Buschelman av=aa->a; bv=bb->a; 88397b6df1SKris Buschelman } 89397b6df1SKris Buschelman 90397b6df1SKris Buschelman if (!valOnly){ 91397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr); 92397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr); 93397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 94397b6df1SKris Buschelman *r = row; *c = col; *v = val; 95397b6df1SKris Buschelman } else { 96397b6df1SKris Buschelman row = *r; col = *c; val = *v; 97397b6df1SKris Buschelman } 98397b6df1SKris Buschelman *nnz = nz; 99397b6df1SKris Buschelman 100397b6df1SKris Buschelman jj = 0; jB = 0; irow = rstart; 101397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 102397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 103397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 104397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 105397b6df1SKris Buschelman bjj = bj + bi[i]; 106397b6df1SKris Buschelman 107397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 108397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 109397b6df1SKris Buschelman for (j=0; j<countB; j++){ 110397b6df1SKris Buschelman jcol = garray[bjj[j]]; 111397b6df1SKris Buschelman if (jcol > colA_start) { jB = j; break; } 112397b6df1SKris Buschelman if (j==countB-1) jB = countB; 113397b6df1SKris Buschelman } 114397b6df1SKris Buschelman 115397b6df1SKris Buschelman /* B-part, smaller col index */ 116397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 117397b6df1SKris Buschelman for (j=0; j<jB; j++){ 118397b6df1SKris Buschelman jcol = garray[bjj[j]]; 119397b6df1SKris Buschelman if (!valOnly){ 120397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 121397b6df1SKris Buschelman } 122397b6df1SKris Buschelman val[jj++] = *bv++; 123397b6df1SKris Buschelman } 124397b6df1SKris Buschelman /* A-part */ 125397b6df1SKris Buschelman for (j=0; j<countA; j++){ 126397b6df1SKris Buschelman if (!valOnly){ 127397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 128397b6df1SKris Buschelman } 129397b6df1SKris Buschelman val[jj++] = *av++; 130397b6df1SKris Buschelman } 131397b6df1SKris Buschelman /* B-part, larger col index */ 132397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 133397b6df1SKris Buschelman if (!valOnly){ 134397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 135397b6df1SKris Buschelman } 136397b6df1SKris Buschelman val[jj++] = *bv++; 137397b6df1SKris Buschelman } 138397b6df1SKris Buschelman irow++; 139397b6df1SKris Buschelman } 140397b6df1SKris Buschelman 141397b6df1SKris Buschelman PetscFunctionReturn(0); 142397b6df1SKris Buschelman } 143397b6df1SKris Buschelman 144c338a77dSKris Buschelman EXTERN_C_BEGIN 145c338a77dSKris Buschelman #undef __FUNCT__ 146c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base" 147c338a77dSKris Buschelman int MatConvert_MUMPS_Base(Mat A,MatType type,Mat *newmat) { 148c338a77dSKris Buschelman /* This routine is only called to convert an unfactored PETSc-MUMPS matrix */ 149c338a77dSKris Buschelman /* to its base PETSc type, so we will ignore 'MatType type'. */ 150c338a77dSKris Buschelman int ierr; 151c338a77dSKris Buschelman Mat B=*newmat; 152f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 153c338a77dSKris Buschelman 154c338a77dSKris Buschelman PetscFunctionBegin; 155c338a77dSKris Buschelman if (B != A) { 156c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 157c338a77dSKris Buschelman } 158f0c56d0fSKris Buschelman B->ops->duplicate = mumps->MatDuplicate; 159f0c56d0fSKris Buschelman B->ops->view = mumps->MatView; 160f0c56d0fSKris Buschelman B->ops->assemblyend = mumps->MatAssemblyEnd; 161f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = mumps->MatLUFactorSymbolic; 162f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic; 163f0c56d0fSKris Buschelman B->ops->destroy = mumps->MatDestroy; 164f0c56d0fSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,mumps->basetype);CHKERRQ(ierr); 165f0c56d0fSKris Buschelman ierr = PetscFree(mumps);CHKERRQ(ierr); 166c338a77dSKris Buschelman *newmat = B; 167c338a77dSKris Buschelman PetscFunctionReturn(0); 168c338a77dSKris Buschelman } 169c338a77dSKris Buschelman EXTERN_C_END 170c338a77dSKris Buschelman 171397b6df1SKris Buschelman #undef __FUNCT__ 172f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS" 173f0c56d0fSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) { 174f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 175c338a77dSKris Buschelman int ierr,size=lu->size; 176397b6df1SKris Buschelman 177397b6df1SKris Buschelman PetscFunctionBegin; 178397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 179397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 180397b6df1SKris Buschelman lu->id.job=JOB_END; 181397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 182397b6df1SKris Buschelman zmumps_c(&lu->id); 183397b6df1SKris Buschelman #else 184397b6df1SKris Buschelman dmumps_c(&lu->id); 185397b6df1SKris Buschelman #endif 186c338a77dSKris Buschelman if (lu->irn) { 187c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 188c338a77dSKris Buschelman } 189c338a77dSKris Buschelman if (lu->jcn) { 190c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 191c338a77dSKris Buschelman } 192c338a77dSKris Buschelman if (size>1 && lu->val) { 193c338a77dSKris Buschelman ierr = PetscFree(lu->val);CHKERRQ(ierr); 194c338a77dSKris Buschelman } 195397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 196397b6df1SKris Buschelman } 197c338a77dSKris Buschelman ierr = MatConvert_MUMPS_Base(A,lu->basetype,&A);CHKERRQ(ierr); 198c338a77dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 199397b6df1SKris Buschelman PetscFunctionReturn(0); 200397b6df1SKris Buschelman } 201397b6df1SKris Buschelman 202397b6df1SKris Buschelman #undef __FUNCT__ 203c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS" 204f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 205f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 206397b6df1SKris Buschelman int ierr; 207397b6df1SKris Buschelman 208397b6df1SKris Buschelman PetscFunctionBegin; 209c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 210c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 211c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 212c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 213c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 214c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 215c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 216c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 217c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 218c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 219c338a77dSKris Buschelman if (lu->myid == 0 && lu->id.ICNTL(11)>0) { 220c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 221c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 222c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 223c338a77dSKris 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); 224c338a77dSKris Buschelman ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 225c338a77dSKris 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); 226c338a77dSKris Buschelman 227c338a77dSKris Buschelman } 228c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 229c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 230c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (efficiency control): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 231c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(15) (efficiency control): %d \n",lu->id.ICNTL(15));CHKERRQ(ierr); 232c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 233c338a77dSKris Buschelman 234c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 235c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 236c338a77dSKris Buschelman ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 237397b6df1SKris Buschelman PetscFunctionReturn(0); 238397b6df1SKris Buschelman } 239397b6df1SKris Buschelman 240397b6df1SKris Buschelman #undef __FUNCT__ 241f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS" 242f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) { 243397b6df1SKris Buschelman int ierr; 244397b6df1SKris Buschelman PetscTruth isascii; 245397b6df1SKris Buschelman PetscViewerFormat format; 246f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 247397b6df1SKris Buschelman 248397b6df1SKris Buschelman PetscFunctionBegin; 249397b6df1SKris Buschelman ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr); 250397b6df1SKris Buschelman 251397b6df1SKris Buschelman ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 252397b6df1SKris Buschelman if (isascii) { 253397b6df1SKris Buschelman ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 254397b6df1SKris Buschelman if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) { 255397b6df1SKris Buschelman ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 256397b6df1SKris Buschelman } 257397b6df1SKris Buschelman } 258397b6df1SKris Buschelman PetscFunctionReturn(0); 259397b6df1SKris Buschelman } 260397b6df1SKris Buschelman 261397b6df1SKris Buschelman #undef __FUNCT__ 262f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS" 263f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) { 264f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 265d54de34fSKris Buschelman PetscScalar *array; 266397b6df1SKris Buschelman Vec x_seq; 267397b6df1SKris Buschelman IS iden; 268397b6df1SKris Buschelman VecScatter scat; 269397b6df1SKris Buschelman int ierr; 270397b6df1SKris Buschelman 271397b6df1SKris Buschelman PetscFunctionBegin; 272397b6df1SKris Buschelman if (lu->size > 1){ 273397b6df1SKris Buschelman if (!lu->myid){ 274397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr); 275397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr); 276397b6df1SKris Buschelman } else { 277397b6df1SKris Buschelman ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr); 278397b6df1SKris Buschelman ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr); 279397b6df1SKris Buschelman } 280397b6df1SKris Buschelman ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr); 281397b6df1SKris Buschelman ierr = ISDestroy(iden);CHKERRQ(ierr); 282397b6df1SKris Buschelman 283397b6df1SKris Buschelman ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 284397b6df1SKris Buschelman ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr); 285397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 286397b6df1SKris Buschelman } else { /* size == 1 */ 287397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 288397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 289397b6df1SKris Buschelman } 290397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 291397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 292397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 293397b6df1SKris Buschelman #else 294397b6df1SKris Buschelman lu->id.rhs = array; 295397b6df1SKris Buschelman #endif 296397b6df1SKris Buschelman } 297397b6df1SKris Buschelman 298397b6df1SKris Buschelman /* solve phase */ 299397b6df1SKris Buschelman lu->id.job=3; 300397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 301397b6df1SKris Buschelman zmumps_c(&lu->id); 302397b6df1SKris Buschelman #else 303397b6df1SKris Buschelman dmumps_c(&lu->id); 304397b6df1SKris Buschelman #endif 305397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 306397b6df1SKris Buschelman SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 307397b6df1SKris Buschelman } 308397b6df1SKris Buschelman 309397b6df1SKris Buschelman /* convert mumps solution x_seq to petsc mpi x */ 310397b6df1SKris Buschelman if (lu->size > 1) { 311397b6df1SKris Buschelman if (!lu->myid){ 312397b6df1SKris Buschelman ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 313397b6df1SKris Buschelman } 314397b6df1SKris Buschelman ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 315397b6df1SKris Buschelman ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr); 316397b6df1SKris Buschelman ierr = VecScatterDestroy(scat);CHKERRQ(ierr); 317397b6df1SKris Buschelman ierr = VecDestroy(x_seq);CHKERRQ(ierr); 318397b6df1SKris Buschelman } else { 319397b6df1SKris Buschelman ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 320397b6df1SKris Buschelman } 321397b6df1SKris Buschelman 322397b6df1SKris Buschelman PetscFunctionReturn(0); 323397b6df1SKris Buschelman } 324397b6df1SKris Buschelman 325397b6df1SKris Buschelman #undef __FUNCT__ 326f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS" 327f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) { 328f0c56d0fSKris Buschelman Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 329f0c56d0fSKris Buschelman Mat_MUMPS *lua=(Mat_MUMPS*)(A)->spptr; 330397b6df1SKris Buschelman int rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl; 331397b6df1SKris Buschelman PetscTruth valOnly,flg; 332397b6df1SKris Buschelman 333397b6df1SKris Buschelman PetscFunctionBegin; 334397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 335f0c56d0fSKris Buschelman (*F)->ops->solve = MatSolve_AIJMUMPS; 336397b6df1SKris Buschelman 337397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 338397b6df1SKris Buschelman ierr = MPI_Comm_rank(A->comm, &lu->myid); 339397b6df1SKris Buschelman ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr); 340397b6df1SKris Buschelman lu->id.job = JOB_INIT; 341397b6df1SKris Buschelman ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr); 342397b6df1SKris Buschelman lu->id.comm_fortran = lu->comm_mumps; 343397b6df1SKris Buschelman 344397b6df1SKris Buschelman /* Set mumps options */ 345397b6df1SKris Buschelman ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 346397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 347397b6df1SKris Buschelman lu->id.sym=lu->sym; 348397b6df1SKris Buschelman if (lu->sym == 2){ 349397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 350397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 351397b6df1SKris Buschelman } 352397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 353397b6df1SKris Buschelman zmumps_c(&lu->id); 354397b6df1SKris Buschelman #else 355397b6df1SKris Buschelman dmumps_c(&lu->id); 356397b6df1SKris Buschelman #endif 357397b6df1SKris Buschelman 358397b6df1SKris Buschelman if (lu->size == 1){ 359397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 360397b6df1SKris Buschelman } else { 361397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 362397b6df1SKris Buschelman } 363397b6df1SKris Buschelman 364397b6df1SKris Buschelman icntl=-1; 365397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 366397b6df1SKris Buschelman if (flg && icntl > 0) { 367397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 368397b6df1SKris Buschelman } else { /* no output */ 369397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 370397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 371397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 372397b6df1SKris Buschelman lu->id.ICNTL(4) = 0; /* level of printing, 0,1,2,3,4, default=2 */ 373397b6df1SKris Buschelman } 374397b6df1SKris 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); 375397b6df1SKris Buschelman icntl=-1; 376397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 377397b6df1SKris Buschelman if (flg) { 378397b6df1SKris Buschelman if (icntl== 1){ 379397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 380397b6df1SKris Buschelman } else { 381397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 382397b6df1SKris Buschelman } 383397b6df1SKris Buschelman } 384397b6df1SKris 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); 385397b6df1SKris 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); 386397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -sles_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr); 387397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 388397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 389397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): efficiency control","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr); 390397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 391397b6df1SKris Buschelman 392397b6df1SKris Buschelman /* 393397b6df1SKris 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); 394397b6df1SKris Buschelman if (flg){ 395397b6df1SKris Buschelman if (icntl >-1 && icntl <3 ){ 396397b6df1SKris Buschelman if (lu->myid==0) lu->id.ICNTL(16) = icntl; 397397b6df1SKris Buschelman } else { 398397b6df1SKris Buschelman SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl); 399397b6df1SKris Buschelman } 400397b6df1SKris Buschelman } 401397b6df1SKris Buschelman */ 402397b6df1SKris Buschelman 403397b6df1SKris 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); 404397b6df1SKris 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); 405397b6df1SKris 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); 406397b6df1SKris Buschelman PetscOptionsEnd(); 407397b6df1SKris Buschelman } 408397b6df1SKris Buschelman 409397b6df1SKris Buschelman /* define matrix A */ 410397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 411397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 412397b6df1SKris Buschelman if (!lu->myid) { 413c36ead0aSKris Buschelman if (lua->isAIJ){ 414397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 415397b6df1SKris Buschelman nz = aa->nz; 416397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 417397b6df1SKris Buschelman } else { 418397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 419397b6df1SKris Buschelman nz = aa->s_nz; 420397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 421397b6df1SKris Buschelman } 422397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 423397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr); 424397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr); 425397b6df1SKris Buschelman nz = 0; 426397b6df1SKris Buschelman for (i=0; i<M; i++){ 427397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 428397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 429397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 430397b6df1SKris Buschelman } 431397b6df1SKris Buschelman } 432397b6df1SKris Buschelman } 433397b6df1SKris Buschelman } 434397b6df1SKris Buschelman break; 435397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 436397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 437397b6df1SKris Buschelman valOnly = PETSC_FALSE; 438397b6df1SKris Buschelman } else { 439397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 440397b6df1SKris Buschelman } 441397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 442397b6df1SKris Buschelman break; 443397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 444397b6df1SKris Buschelman } 445397b6df1SKris Buschelman 446397b6df1SKris Buschelman /* analysis phase */ 447397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 448397b6df1SKris Buschelman lu->id.n = M; 449397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 450397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 451397b6df1SKris Buschelman if (!lu->myid) { 452397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 453397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 454397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 455397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 456397b6df1SKris Buschelman #else 457397b6df1SKris Buschelman lu->id.a = lu->val; 458397b6df1SKris Buschelman #endif 459397b6df1SKris Buschelman } 460397b6df1SKris Buschelman } 461397b6df1SKris Buschelman break; 462397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 463397b6df1SKris Buschelman lu->id.nz_loc = nnz; 464397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 465397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 466397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 467397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 468397b6df1SKris Buschelman #else 469397b6df1SKris Buschelman lu->id.a_loc = lu->val; 470397b6df1SKris Buschelman #endif 471397b6df1SKris Buschelman } 472397b6df1SKris Buschelman break; 473397b6df1SKris Buschelman } 474397b6df1SKris Buschelman lu->id.job=1; 475397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 476397b6df1SKris Buschelman zmumps_c(&lu->id); 477397b6df1SKris Buschelman #else 478397b6df1SKris Buschelman dmumps_c(&lu->id); 479397b6df1SKris Buschelman #endif 480397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 481397b6df1SKris Buschelman SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 482397b6df1SKris Buschelman } 483397b6df1SKris Buschelman } 484397b6df1SKris Buschelman 485397b6df1SKris Buschelman /* numerical factorization phase */ 486397b6df1SKris Buschelman if(lu->id.ICNTL(18) == 0) { 487397b6df1SKris Buschelman if (lu->myid == 0) { 488397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 489397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 490397b6df1SKris Buschelman #else 491397b6df1SKris Buschelman lu->id.a = lu->val; 492397b6df1SKris Buschelman #endif 493397b6df1SKris Buschelman } 494397b6df1SKris Buschelman } else { 495397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 496397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 497397b6df1SKris Buschelman #else 498397b6df1SKris Buschelman lu->id.a_loc = lu->val; 499397b6df1SKris Buschelman #endif 500397b6df1SKris Buschelman } 501397b6df1SKris Buschelman lu->id.job=2; 502397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 503397b6df1SKris Buschelman zmumps_c(&lu->id); 504397b6df1SKris Buschelman #else 505397b6df1SKris Buschelman dmumps_c(&lu->id); 506397b6df1SKris Buschelman #endif 507397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 508397b6df1SKris Buschelman SETERRQ1(1,"1, Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 509397b6df1SKris Buschelman } 510397b6df1SKris Buschelman 511397b6df1SKris Buschelman if (lu->myid==0 && lu->id.ICNTL(16) > 0){ 512397b6df1SKris Buschelman SETERRQ1(1," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 513397b6df1SKris Buschelman } 514397b6df1SKris Buschelman 515397b6df1SKris Buschelman (*F)->assembled = PETSC_TRUE; 516397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 517ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 518397b6df1SKris Buschelman PetscFunctionReturn(0); 519397b6df1SKris Buschelman } 520397b6df1SKris Buschelman 521397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 522397b6df1SKris Buschelman #undef __FUNCT__ 523f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 524f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) { 525397b6df1SKris Buschelman Mat B; 526f0c56d0fSKris Buschelman Mat_MUMPS *lu; 527397b6df1SKris Buschelman int ierr; 528397b6df1SKris Buschelman 529397b6df1SKris Buschelman PetscFunctionBegin; 530397b6df1SKris Buschelman 531397b6df1SKris Buschelman /* Create the factorization matrix */ 532397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 533397b6df1SKris Buschelman ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr); 534397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 535397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 536397b6df1SKris Buschelman 537f0c56d0fSKris Buschelman B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS; 538397b6df1SKris Buschelman B->factor = FACTOR_LU; 539f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 540397b6df1SKris Buschelman lu->sym = 0; 541397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 542397b6df1SKris Buschelman 543397b6df1SKris Buschelman *F = B; 544397b6df1SKris Buschelman PetscFunctionReturn(0); 545397b6df1SKris Buschelman } 546397b6df1SKris Buschelman 547397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 548397b6df1SKris Buschelman #undef __FUNCT__ 549f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 550f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) { 551397b6df1SKris Buschelman Mat B; 552f0c56d0fSKris Buschelman Mat_MUMPS *lu; 553397b6df1SKris Buschelman int ierr; 554397b6df1SKris Buschelman 555397b6df1SKris Buschelman PetscFunctionBegin; 556397b6df1SKris Buschelman 557397b6df1SKris Buschelman /* Create the factorization matrix */ 558397b6df1SKris Buschelman ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr); 559397b6df1SKris Buschelman ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr); 560397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 561397b6df1SKris Buschelman ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 562397b6df1SKris Buschelman 563f0c56d0fSKris Buschelman B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS; 564397b6df1SKris Buschelman B->factor = FACTOR_CHOLESKY; 565f0c56d0fSKris Buschelman lu = (Mat_MUMPS*)B->spptr; 566397b6df1SKris Buschelman lu->sym = 2; 567397b6df1SKris Buschelman lu->matstruc = DIFFERENT_NONZERO_PATTERN; 568397b6df1SKris Buschelman 569397b6df1SKris Buschelman *F = B; 570397b6df1SKris Buschelman PetscFunctionReturn(0); 571397b6df1SKris Buschelman } 572397b6df1SKris Buschelman 573397b6df1SKris Buschelman #undef __FUNCT__ 574f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS" 575f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) { 576c338a77dSKris Buschelman int ierr; 577f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 578c338a77dSKris Buschelman 579397b6df1SKris Buschelman PetscFunctionBegin; 580c338a77dSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 581f0c56d0fSKris Buschelman 582c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 583c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 584f0c56d0fSKris Buschelman A->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 585397b6df1SKris Buschelman PetscFunctionReturn(0); 586397b6df1SKris Buschelman } 587397b6df1SKris Buschelman 588c338a77dSKris Buschelman EXTERN_C_BEGIN 589c338a77dSKris Buschelman #undef __FUNCT__ 590f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS" 591f0c56d0fSKris Buschelman int MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,Mat *newmat) { 592c338a77dSKris Buschelman int ierr,size; 593c338a77dSKris Buschelman MPI_Comm comm; 594c338a77dSKris Buschelman Mat B=*newmat; 595f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 596397b6df1SKris Buschelman 597397b6df1SKris Buschelman PetscFunctionBegin; 598c338a77dSKris Buschelman if (B != A) { 599c338a77dSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 600397b6df1SKris Buschelman } 601397b6df1SKris Buschelman 602c338a77dSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 603f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 604c338a77dSKris Buschelman 605f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 606c338a77dSKris Buschelman mumps->MatView = A->ops->view; 607c338a77dSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 608c338a77dSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 609c338a77dSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 610c338a77dSKris Buschelman mumps->MatDestroy = A->ops->destroy; 611c338a77dSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 612f579278aSKris Buschelman mumps->isAIJ = PETSC_TRUE; 613c338a77dSKris Buschelman 6144b68dd72SKris Buschelman B->spptr = (void *)mumps; 615f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_AIJMUMPS; 616f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 617f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_AIJMUMPS; 618f0c56d0fSKris Buschelman B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 619f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_AIJMUMPS; 620c338a77dSKris Buschelman 621c338a77dSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 622c338a77dSKris Buschelman if (size == 1) { 623c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C", 624f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 625c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C", 626c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 627c338a77dSKris Buschelman } else { 628c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C", 629f0c56d0fSKris Buschelman "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr); 630c338a77dSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C", 631c338a77dSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 632c338a77dSKris Buschelman } 633c338a77dSKris Buschelman 634f579278aSKris Buschelman PetscLogInfo(0,"Using MUMPS for LU factorization and solves."); 635c338a77dSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 636c338a77dSKris Buschelman *newmat = B; 637397b6df1SKris Buschelman PetscFunctionReturn(0); 638397b6df1SKris Buschelman } 639c338a77dSKris Buschelman EXTERN_C_END 640397b6df1SKris Buschelman 641f0c56d0fSKris Buschelman #undef __FUNCT__ 642f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS" 643f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) { 644f0c56d0fSKris Buschelman int ierr; 645f0c56d0fSKris Buschelman PetscFunctionBegin; 646f0c56d0fSKris Buschelman ierr = (*A->ops->duplicate)(A,op,M);CHKERRQ(ierr); 647f0c56d0fSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr); 648f0c56d0fSKris Buschelman PetscFunctionReturn(0); 649f0c56d0fSKris Buschelman } 650f0c56d0fSKris Buschelman 65124b6179bSKris Buschelman /*MC 652f579278aSKris Buschelman MATAIJMUMPS - a matrix type providing direct solvers (LU) for distributed 65324b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 65424b6179bSKris Buschelman 65524b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 65624b6179bSKris Buschelman on how to declare the existence of external packages), 65724b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 65824b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS). 65924b6179bSKris Buschelman This matrix type is only supported for double precision real. 66024b6179bSKris Buschelman 66124b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 66224b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 66324b6179bSKris Buschelman MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported 66424b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 66524b6179bSKris Buschelman the above preallocation routines for simplicity. 66624b6179bSKris Buschelman 66724b6179bSKris Buschelman Options Database Keys: 66824b6179bSKris Buschelman + -mat_type aijmumps 66924b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 67024b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 67124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 67224b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 67324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 67424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 67524b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view 67624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 67724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 67824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 67924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 68024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 68124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 68224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 68324b6179bSKris Buschelman 68424b6179bSKris Buschelman Level: beginner 68524b6179bSKris Buschelman 68624b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 68724b6179bSKris Buschelman M*/ 68824b6179bSKris Buschelman 689397b6df1SKris Buschelman EXTERN_C_BEGIN 690397b6df1SKris Buschelman #undef __FUNCT__ 691f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS" 692f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) { 693397b6df1SKris Buschelman int ierr,size; 694397b6df1SKris Buschelman MPI_Comm comm; 695397b6df1SKris Buschelman 696397b6df1SKris Buschelman PetscFunctionBegin; 697*5441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */ 698*5441df8eSKris Buschelman /* and AIJMUMPS types */ 699*5441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr); 700397b6df1SKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 701397b6df1SKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 702397b6df1SKris Buschelman if (size == 1) { 703397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 704397b6df1SKris Buschelman } else { 705397b6df1SKris Buschelman ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 706397b6df1SKris Buschelman } 707f0c56d0fSKris Buschelman ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr); 708397b6df1SKris Buschelman PetscFunctionReturn(0); 709397b6df1SKris Buschelman } 710397b6df1SKris Buschelman EXTERN_C_END 711397b6df1SKris Buschelman 712f579278aSKris Buschelman EXTERN_C_BEGIN 713f579278aSKris Buschelman #undef __FUNCT__ 714f0c56d0fSKris Buschelman #define __FUNCT__ "MatLoad_AIJMUMPS" 715f0c56d0fSKris Buschelman int MatLoad_AIJMUMPS(PetscViewer viewer,MatType type,Mat *A) { 716f579278aSKris Buschelman int ierr,size,(*r)(PetscViewer,MatType,Mat*); 717f579278aSKris Buschelman MPI_Comm comm; 718f579278aSKris Buschelman 719f579278aSKris Buschelman PetscFunctionBegin; 720f579278aSKris Buschelman ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 721f579278aSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 722f579278aSKris Buschelman if (size == 1) { 723f579278aSKris Buschelman ierr = PetscFListFind(comm,MatLoadList,MATSEQAIJ,(void(**)(void))&r);CHKERRQ(ierr); 724f579278aSKris Buschelman } else { 725f579278aSKris Buschelman ierr = PetscFListFind(comm,MatLoadList,MATMPIAIJ,(void(**)(void))&r);CHKERRQ(ierr); 726f579278aSKris Buschelman } 727f579278aSKris Buschelman ierr = (*r)(viewer,type,A);CHKERRQ(ierr); 728f579278aSKris Buschelman PetscFunctionReturn(0); 729f579278aSKris Buschelman } 730f579278aSKris Buschelman EXTERN_C_END 731f579278aSKris Buschelman 732f579278aSKris Buschelman #undef __FUNCT__ 733f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS" 734f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) { 735f579278aSKris Buschelman int ierr; 736f0c56d0fSKris Buschelman Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 737f579278aSKris Buschelman 738f579278aSKris Buschelman PetscFunctionBegin; 739f579278aSKris Buschelman ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr); 740f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 741f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 742f0c56d0fSKris Buschelman A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 743f579278aSKris Buschelman PetscFunctionReturn(0); 744f579278aSKris Buschelman } 745f579278aSKris Buschelman 746f579278aSKris Buschelman EXTERN_C_BEGIN 747f579278aSKris Buschelman #undef __FUNCT__ 748f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS" 749f0c56d0fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,Mat *newmat) { 750f579278aSKris Buschelman int ierr,size; 751f579278aSKris Buschelman MPI_Comm comm; 752f579278aSKris Buschelman Mat B=*newmat; 753f0c56d0fSKris Buschelman Mat_MUMPS *mumps; 754f579278aSKris Buschelman 755f579278aSKris Buschelman PetscFunctionBegin; 756f579278aSKris Buschelman if (B != A) { 757f579278aSKris Buschelman ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 758f579278aSKris Buschelman } 759f579278aSKris Buschelman 760f579278aSKris Buschelman ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 761f0c56d0fSKris Buschelman ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr); 762f579278aSKris Buschelman 763f0c56d0fSKris Buschelman mumps->MatDuplicate = A->ops->duplicate; 764f579278aSKris Buschelman mumps->MatView = A->ops->view; 765f579278aSKris Buschelman mumps->MatAssemblyEnd = A->ops->assemblyend; 766f579278aSKris Buschelman mumps->MatLUFactorSymbolic = A->ops->lufactorsymbolic; 767f579278aSKris Buschelman mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic; 768f579278aSKris Buschelman mumps->MatDestroy = A->ops->destroy; 769f579278aSKris Buschelman mumps->CleanUpMUMPS = PETSC_FALSE; 770f579278aSKris Buschelman mumps->isAIJ = PETSC_FALSE; 771f579278aSKris Buschelman 772f579278aSKris Buschelman B->spptr = (void *)mumps; 773f0c56d0fSKris Buschelman B->ops->duplicate = MatDuplicate_SBAIJMUMPS; 774f0c56d0fSKris Buschelman B->ops->view = MatView_AIJMUMPS; 775f0c56d0fSKris Buschelman B->ops->assemblyend = MatAssemblyEnd_SBAIJMUMPS; 776f0c56d0fSKris Buschelman B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 777f0c56d0fSKris Buschelman B->ops->destroy = MatDestroy_AIJMUMPS; 778f579278aSKris Buschelman 779f579278aSKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 780f579278aSKris Buschelman if (size == 1) { 781f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C", 782f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 783f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C", 784f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 785f579278aSKris Buschelman } else { 786f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C", 787f0c56d0fSKris Buschelman "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr); 788f0c56d0fSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C", 789f579278aSKris Buschelman "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr); 790f579278aSKris Buschelman } 791f579278aSKris Buschelman 792f579278aSKris Buschelman PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves."); 793f579278aSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr); 794f579278aSKris Buschelman *newmat = B; 795f579278aSKris Buschelman PetscFunctionReturn(0); 796f579278aSKris Buschelman } 797f579278aSKris Buschelman EXTERN_C_END 798f579278aSKris Buschelman 799f0c56d0fSKris Buschelman #undef __FUNCT__ 800f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS" 801f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) { 802f0c56d0fSKris Buschelman int ierr; 803f0c56d0fSKris Buschelman PetscFunctionBegin; 804f0c56d0fSKris Buschelman ierr = (*A->ops->duplicate)(A,op,M);CHKERRQ(ierr); 805f0c56d0fSKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr); 806f0c56d0fSKris Buschelman PetscFunctionReturn(0); 807f0c56d0fSKris Buschelman } 808f0c56d0fSKris Buschelman 80924b6179bSKris Buschelman /*MC 810f579278aSKris Buschelman MATSBAIJMUMPS - a symmetric matrix type providing direct solvers (Cholesky) for 81124b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 81224b6179bSKris Buschelman 81324b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 81424b6179bSKris Buschelman on how to declare the existence of external packages), 81524b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 81624b6179bSKris Buschelman After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS). 81724b6179bSKris Buschelman This matrix type is only supported for double precision real. 81824b6179bSKris Buschelman 81924b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 82024b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 82124b6179bSKris Buschelman MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported 82224b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 82324b6179bSKris Buschelman the above preallocation routines for simplicity. 82424b6179bSKris Buschelman 82524b6179bSKris Buschelman Options Database Keys: 82624b6179bSKris Buschelman + -mat_type aijmumps 82724b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 82824b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 82924b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 83024b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 83124b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 83224b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 83324b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view 83424b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 83524b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 83624b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 83724b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 83824b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 83924b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 84024b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 84124b6179bSKris Buschelman 84224b6179bSKris Buschelman Level: beginner 84324b6179bSKris Buschelman 84424b6179bSKris Buschelman .seealso: MATAIJMUMPS 84524b6179bSKris Buschelman M*/ 84624b6179bSKris Buschelman 847397b6df1SKris Buschelman EXTERN_C_BEGIN 848397b6df1SKris Buschelman #undef __FUNCT__ 849f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS" 850f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) { 851397b6df1SKris Buschelman int ierr,size; 852397b6df1SKris Buschelman 853397b6df1SKris Buschelman PetscFunctionBegin; 854*5441df8eSKris Buschelman /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */ 855*5441df8eSKris Buschelman /* and SBAIJMUMPS types */ 856*5441df8eSKris Buschelman ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr); 857*5441df8eSKris Buschelman ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 858397b6df1SKris Buschelman if (size == 1) { 859397b6df1SKris Buschelman ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr); 860397b6df1SKris Buschelman } else { 861397b6df1SKris Buschelman ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr); 862397b6df1SKris Buschelman } 863f0c56d0fSKris Buschelman ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr); 864397b6df1SKris Buschelman PetscFunctionReturn(0); 865397b6df1SKris Buschelman } 866397b6df1SKris Buschelman EXTERN_C_END 867397b6df1SKris Buschelman 868397b6df1SKris Buschelman EXTERN_C_BEGIN 869397b6df1SKris Buschelman #undef __FUNCT__ 870f0c56d0fSKris Buschelman #define __FUNCT__ "MatLoad_SBAIJMUMPS" 871f0c56d0fSKris Buschelman int MatLoad_SBAIJMUMPS(PetscViewer viewer,MatType type,Mat *A) { 872397b6df1SKris Buschelman int ierr,size,(*r)(PetscViewer,MatType,Mat*); 873397b6df1SKris Buschelman MPI_Comm comm; 874397b6df1SKris Buschelman 875397b6df1SKris Buschelman PetscFunctionBegin; 876397b6df1SKris Buschelman ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 877397b6df1SKris Buschelman ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 878397b6df1SKris Buschelman if (size == 1) { 879397b6df1SKris Buschelman ierr = PetscFListFind(comm,MatLoadList,MATSEQSBAIJ,(void(**)(void))&r);CHKERRQ(ierr); 880397b6df1SKris Buschelman } else { 881397b6df1SKris Buschelman ierr = PetscFListFind(comm,MatLoadList,MATMPISBAIJ,(void(**)(void))&r);CHKERRQ(ierr); 882397b6df1SKris Buschelman } 883397b6df1SKris Buschelman ierr = (*r)(viewer,type,A);CHKERRQ(ierr); 884397b6df1SKris Buschelman PetscFunctionReturn(0); 885397b6df1SKris Buschelman } 886397b6df1SKris Buschelman EXTERN_C_END 887