1be1d678aSKris Buschelman #define PETSCMAT_DLL 21c2a3de1SBarry Smith 3397b6df1SKris Buschelman /* 4c2b5dc30SHong Zhang Provides an interface to the MUMPS sparse solver 5397b6df1SKris Buschelman */ 6397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h" 7397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h" 8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h" 9397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h" 10397b6df1SKris Buschelman 11397b6df1SKris Buschelman EXTERN_C_BEGIN 12397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 13397b6df1SKris Buschelman #include "zmumps_c.h" 14397b6df1SKris Buschelman #else 15397b6df1SKris Buschelman #include "dmumps_c.h" 16397b6df1SKris Buschelman #endif 17397b6df1SKris Buschelman EXTERN_C_END 18397b6df1SKris Buschelman #define JOB_INIT -1 19397b6df1SKris Buschelman #define JOB_END -2 20397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 21397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 22397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 23397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 24a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 25397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 26adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 27397b6df1SKris Buschelman 28397b6df1SKris Buschelman typedef struct { 29397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 30397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 31397b6df1SKris Buschelman #else 32397b6df1SKris Buschelman DMUMPS_STRUC_C id; 33397b6df1SKris Buschelman #endif 34397b6df1SKris Buschelman MatStructure matstruc; 35c1490034SHong Zhang PetscMPIInt myid,size; 36329ec9b3SHong Zhang PetscInt *irn,*jcn,sym,nSolve; 37397b6df1SKris Buschelman PetscScalar *val; 38397b6df1SKris Buschelman MPI_Comm comm_mumps; 39329ec9b3SHong Zhang VecScatter scat_rhs, scat_sol; 40c338a77dSKris Buschelman PetscTruth isAIJ,CleanUpMUMPS; 41329ec9b3SHong Zhang Vec b_seq,x_seq; 42f0c56d0fSKris Buschelman } Mat_MUMPS; 43f0c56d0fSKris Buschelman 44dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 45b24902e0SBarry Smith 46397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 47397b6df1SKris Buschelman /* 48397b6df1SKris Buschelman input: 4975747be1SHong Zhang A - matrix in mpiaij or mpisbaij (bs=1) format 50397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 51397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 52397b6df1SKris Buschelman TRUE: only the values in v array are updated 53397b6df1SKris Buschelman output: 54397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 55397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 56397b6df1SKris Buschelman */ 57b24902e0SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) 58b24902e0SBarry Smith { 59c1490034SHong Zhang PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray; 60dfbe8321SBarry Smith PetscErrorCode ierr; 612a4c71feSBarry Smith PetscInt i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol; 62c1490034SHong Zhang PetscInt *row,*col; 63397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 64*5c9eb25fSBarry Smith PetscTruth isAIJ; 65397b6df1SKris Buschelman 66397b6df1SKris Buschelman PetscFunctionBegin; 67*5c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr); 68*5c9eb25fSBarry Smith if (isAIJ){ 69397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 70397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 71397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 72397b6df1SKris Buschelman nz = aa->nz + bb->nz; 732a4c71feSBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart; 74397b6df1SKris Buschelman garray = mat->garray; 75397b6df1SKris Buschelman av=aa->a; bv=bb->a; 76397b6df1SKris Buschelman 77397b6df1SKris Buschelman } else { 78397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 79397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 80397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 810c0e133fSBarry Smith if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs); 826c6c5352SBarry Smith nz = aa->nz + bb->nz; 832a4c71feSBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart; 84397b6df1SKris Buschelman garray = mat->garray; 85397b6df1SKris Buschelman av=aa->a; bv=bb->a; 86397b6df1SKris Buschelman } 87397b6df1SKris Buschelman 88397b6df1SKris Buschelman if (!valOnly){ 897c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr); 907c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr); 91397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 92397b6df1SKris Buschelman *r = row; *c = col; *v = val; 93397b6df1SKris Buschelman } else { 94397b6df1SKris Buschelman row = *r; col = *c; val = *v; 95397b6df1SKris Buschelman } 96397b6df1SKris Buschelman *nnz = nz; 97397b6df1SKris Buschelman 98028e57e8SHong Zhang jj = 0; irow = rstart; 99397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 100397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 101397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 102397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 103397b6df1SKris Buschelman bjj = bj + bi[i]; 104397b6df1SKris Buschelman 105397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 106397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 10775747be1SHong Zhang j=-1; 10875747be1SHong Zhang do { 10975747be1SHong Zhang j++; 11075747be1SHong Zhang if (j == countB) break; 111397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11275747be1SHong Zhang } while (jcol < colA_start); 11375747be1SHong Zhang jB = j; 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; 12175747be1SHong Zhang 122397b6df1SKris Buschelman } 123397b6df1SKris Buschelman val[jj++] = *bv++; 124397b6df1SKris Buschelman } 125397b6df1SKris Buschelman /* A-part */ 126397b6df1SKris Buschelman for (j=0; j<countA; j++){ 127397b6df1SKris Buschelman if (!valOnly){ 128397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 129397b6df1SKris Buschelman } 130397b6df1SKris Buschelman val[jj++] = *av++; 131397b6df1SKris Buschelman } 132397b6df1SKris Buschelman /* B-part, larger col index */ 133397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 134397b6df1SKris Buschelman if (!valOnly){ 135397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 136397b6df1SKris Buschelman } 137397b6df1SKris Buschelman val[jj++] = *bv++; 138397b6df1SKris Buschelman } 139397b6df1SKris Buschelman irow++; 140397b6df1SKris Buschelman } 141397b6df1SKris Buschelman 142397b6df1SKris Buschelman PetscFunctionReturn(0); 143397b6df1SKris Buschelman } 144397b6df1SKris Buschelman 145397b6df1SKris Buschelman #undef __FUNCT__ 1463924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 147dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 148dfbe8321SBarry Smith { 149f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 150dfbe8321SBarry Smith PetscErrorCode ierr; 151c1490034SHong Zhang PetscMPIInt size=lu->size; 1526849ba73SBarry Smith PetscErrorCode (*specialdestroy)(Mat); 153b24902e0SBarry Smith 154397b6df1SKris Buschelman PetscFunctionBegin; 155397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 156397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 157329ec9b3SHong Zhang if (size > 1){ 158329ec9b3SHong Zhang ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 159329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 160329ec9b3SHong Zhang ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 1612750af12SHong Zhang if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);} 1622750af12SHong Zhang if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);} 163329ec9b3SHong Zhang ierr = PetscFree(lu->val);CHKERRQ(ierr); 164329ec9b3SHong Zhang } 165397b6df1SKris Buschelman lu->id.job=JOB_END; 166397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 167397b6df1SKris Buschelman zmumps_c(&lu->id); 168397b6df1SKris Buschelman #else 169397b6df1SKris Buschelman dmumps_c(&lu->id); 170397b6df1SKris Buschelman #endif 171c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 172c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 173397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 174397b6df1SKris Buschelman } 175c338a77dSKris Buschelman ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 176397b6df1SKris Buschelman PetscFunctionReturn(0); 177397b6df1SKris Buschelman } 178397b6df1SKris Buschelman 179397b6df1SKris Buschelman #undef __FUNCT__ 180f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 181b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 182b24902e0SBarry Smith { 183f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 184d54de34fSKris Buschelman PetscScalar *array; 185397b6df1SKris Buschelman Vec x_seq; 186329ec9b3SHong Zhang IS is_iden,is_petsc; 187dfbe8321SBarry Smith PetscErrorCode ierr; 188329ec9b3SHong Zhang PetscInt i; 189397b6df1SKris Buschelman 190397b6df1SKris Buschelman PetscFunctionBegin; 191329ec9b3SHong Zhang lu->id.nrhs = 1; 192329ec9b3SHong Zhang x_seq = lu->b_seq; 193397b6df1SKris Buschelman if (lu->size > 1){ 194329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 195f6cfb2d1SLisandro Dalcin ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 196f6cfb2d1SLisandro Dalcin ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 197397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 198397b6df1SKris Buschelman } else { /* size == 1 */ 199397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 200397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 201397b6df1SKris Buschelman } 202397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 203397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 204397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 205397b6df1SKris Buschelman #else 206397b6df1SKris Buschelman lu->id.rhs = array; 207397b6df1SKris Buschelman #endif 208397b6df1SKris Buschelman } 209329ec9b3SHong Zhang if (lu->size == 1){ 210329ec9b3SHong Zhang ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 211329ec9b3SHong Zhang } else if (!lu->myid){ 212329ec9b3SHong Zhang ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 213329ec9b3SHong Zhang } 214329ec9b3SHong Zhang 215329ec9b3SHong Zhang if (lu->size > 1){ 216329ec9b3SHong Zhang /* distributed solution */ 217329ec9b3SHong Zhang lu->id.ICNTL(21) = 1; 218329ec9b3SHong Zhang if (!lu->nSolve){ 219329ec9b3SHong Zhang /* Create x_seq=sol_loc for repeated use */ 220329ec9b3SHong Zhang PetscInt lsol_loc; 221329ec9b3SHong Zhang PetscScalar *sol_loc; 222329ec9b3SHong Zhang lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 223329ec9b3SHong Zhang ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr); 224329ec9b3SHong Zhang lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc); 225329ec9b3SHong Zhang lu->id.lsol_loc = lsol_loc; 2266f8312c5SHong Zhang #if defined(PETSC_USE_COMPLEX) 2276f8312c5SHong Zhang lu->id.sol_loc = (ZMUMPS_DOUBLE *)sol_loc; 2286f8312c5SHong Zhang #else 2296f8312c5SHong Zhang lu->id.sol_loc = (DMUMPS_DOUBLE *)sol_loc; 2306f8312c5SHong Zhang #endif 231329ec9b3SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 232329ec9b3SHong Zhang } 233329ec9b3SHong Zhang } 234397b6df1SKris Buschelman 235397b6df1SKris Buschelman /* solve phase */ 236329ec9b3SHong Zhang /*-------------*/ 237397b6df1SKris Buschelman lu->id.job = 3; 238397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 239397b6df1SKris Buschelman zmumps_c(&lu->id); 240397b6df1SKris Buschelman #else 241397b6df1SKris Buschelman dmumps_c(&lu->id); 242397b6df1SKris Buschelman #endif 243397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 24479a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 245397b6df1SKris Buschelman } 246397b6df1SKris Buschelman 247329ec9b3SHong Zhang if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 248329ec9b3SHong Zhang if (!lu->nSolve){ /* create scatter scat_sol */ 249329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 250329ec9b3SHong Zhang for (i=0; i<lu->id.lsol_loc; i++){ 251329ec9b3SHong Zhang lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 252397b6df1SKris Buschelman } 253329ec9b3SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr); /* to */ 254329ec9b3SHong Zhang ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 255329ec9b3SHong Zhang ierr = ISDestroy(is_iden);CHKERRQ(ierr); 256329ec9b3SHong Zhang ierr = ISDestroy(is_petsc);CHKERRQ(ierr); 257397b6df1SKris Buschelman } 258ca9f406cSSatish Balay ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 259ca9f406cSSatish Balay ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 260329ec9b3SHong Zhang } 261329ec9b3SHong Zhang lu->nSolve++; 262397b6df1SKris Buschelman PetscFunctionReturn(0); 263397b6df1SKris Buschelman } 264397b6df1SKris Buschelman 265a58c3f20SHong Zhang /* 266a58c3f20SHong Zhang input: 267a58c3f20SHong Zhang F: numeric factor 268a58c3f20SHong Zhang output: 269a58c3f20SHong Zhang nneg: total number of negative pivots 270a58c3f20SHong Zhang nzero: 0 271a58c3f20SHong Zhang npos: (global dimension of F) - nneg 272a58c3f20SHong Zhang */ 273a58c3f20SHong Zhang 274a58c3f20SHong Zhang #undef __FUNCT__ 275a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 276dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 277a58c3f20SHong Zhang { 278a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 279dfbe8321SBarry Smith PetscErrorCode ierr; 280c1490034SHong Zhang PetscMPIInt size; 281a58c3f20SHong Zhang 282a58c3f20SHong Zhang PetscFunctionBegin; 2837adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 284bcb30aebSHong 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 */ 285bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 28679a5c55eSBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13)); 287bcb30aebSHong Zhang } 288a58c3f20SHong Zhang if (nneg){ 289a58c3f20SHong Zhang if (!lu->myid){ 290a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 291a58c3f20SHong Zhang } 292bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 293a58c3f20SHong Zhang } 294a58c3f20SHong Zhang if (nzero) *nzero = 0; 2952a4c71feSBarry Smith if (npos) *npos = F->rmap.N - (*nneg); 296a58c3f20SHong Zhang PetscFunctionReturn(0); 297a58c3f20SHong Zhang } 298a58c3f20SHong Zhang 299397b6df1SKris Buschelman #undef __FUNCT__ 300f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 301f6c57405SHong Zhang PetscErrorCode MatFactorNumeric_MUMPS(Mat A,MatFactorInfo *info,Mat *F) 302af281ebdSHong Zhang { 303f0c56d0fSKris Buschelman Mat_MUMPS *lu =(Mat_MUMPS*)(*F)->spptr; 3046849ba73SBarry Smith PetscErrorCode ierr; 3052a4c71feSBarry Smith PetscInt rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl; 306397b6df1SKris Buschelman PetscTruth valOnly,flg; 307e09efc27SHong Zhang Mat F_diag; 308c349612cSHong Zhang IS is_iden; 309c349612cSHong Zhang Vec b; 310*5c9eb25fSBarry Smith PetscTruth isSeqAIJ,isSeqSBAIJ; 311397b6df1SKris Buschelman 312397b6df1SKris Buschelman PetscFunctionBegin; 313*5c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 314*5c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 315397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 316f6c57405SHong Zhang (*F)->ops->solve = MatSolve_MUMPS; 317397b6df1SKris Buschelman 318397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 3197adad957SLisandro Dalcin ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 3207adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 321397b6df1SKris Buschelman lu->id.job = JOB_INIT; 3227adad957SLisandro Dalcin ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 3236a1dac61SBarry Smith lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 324397b6df1SKris Buschelman 325397b6df1SKris Buschelman /* Set mumps options */ 3267adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 327397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 328397b6df1SKris Buschelman lu->id.sym=lu->sym; 329397b6df1SKris Buschelman if (lu->sym == 2){ 330397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 331397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 332397b6df1SKris Buschelman } 333397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 334397b6df1SKris Buschelman zmumps_c(&lu->id); 335397b6df1SKris Buschelman #else 336397b6df1SKris Buschelman dmumps_c(&lu->id); 337397b6df1SKris Buschelman #endif 338397b6df1SKris Buschelman 339*5c9eb25fSBarry Smith if (isSeqAIJ || isSeqSBAIJ){ 340397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 341397b6df1SKris Buschelman } else { 342397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 343397b6df1SKris Buschelman } 344397b6df1SKris Buschelman 345397b6df1SKris Buschelman icntl=-1; 34621f4b680SHong Zhang lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 347397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 34819facb7aSBarry Smith if ((flg && icntl > 0) || PetscLogPrintInfo) { 349397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 350397b6df1SKris Buschelman } else { /* no output */ 351397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 352397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 353397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 354397b6df1SKris Buschelman } 355397b6df1SKris 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); 356397b6df1SKris Buschelman icntl=-1; 357397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 358397b6df1SKris Buschelman if (flg) { 359397b6df1SKris Buschelman if (icntl== 1){ 360397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 361397b6df1SKris Buschelman } else { 362397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 363397b6df1SKris Buschelman } 364397b6df1SKris Buschelman } 365397b6df1SKris 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); 366397b6df1SKris 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); 36794b7f48cSBarry 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); 368397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 369397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 370adc1d99fSHong 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); 371397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 372397b6df1SKris Buschelman 373397b6df1SKris 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); 374397b6df1SKris 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); 375397b6df1SKris 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); 37625f9c88cSHong Zhang ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr); 377397b6df1SKris Buschelman PetscOptionsEnd(); 378397b6df1SKris Buschelman } 379397b6df1SKris Buschelman 380397b6df1SKris Buschelman /* define matrix A */ 381397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 382397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 383397b6df1SKris Buschelman if (!lu->myid) { 384*5c9eb25fSBarry Smith if (isSeqAIJ){ 385397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 386397b6df1SKris Buschelman nz = aa->nz; 387397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 388*5c9eb25fSBarry Smith } else if (isSeqSBAIJ) { 389397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 3906c6c5352SBarry Smith nz = aa->nz; 391397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 392*5c9eb25fSBarry Smith } else { 393*5c9eb25fSBarry Smith SETERRQ(PETSC_ERR_SUP,"No mumps factorization for this matrix type"); 394397b6df1SKris Buschelman } 395397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 3967c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 3977c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 398397b6df1SKris Buschelman nz = 0; 399397b6df1SKris Buschelman for (i=0; i<M; i++){ 400397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 401397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 402397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 403397b6df1SKris Buschelman } 404397b6df1SKris Buschelman } 405397b6df1SKris Buschelman } 406397b6df1SKris Buschelman } 407397b6df1SKris Buschelman break; 408397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 409397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 410397b6df1SKris Buschelman valOnly = PETSC_FALSE; 411397b6df1SKris Buschelman } else { 412397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 413397b6df1SKris Buschelman } 414397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 415397b6df1SKris Buschelman break; 416397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 417397b6df1SKris Buschelman } 418397b6df1SKris Buschelman 419397b6df1SKris Buschelman /* analysis phase */ 420329ec9b3SHong Zhang /*----------------*/ 421397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 422329ec9b3SHong Zhang lu->id.job = 1; 423329ec9b3SHong Zhang 424397b6df1SKris Buschelman lu->id.n = M; 425397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 426397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 427397b6df1SKris Buschelman if (!lu->myid) { 428397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 429397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 430397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 431397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 432397b6df1SKris Buschelman #else 433397b6df1SKris Buschelman lu->id.a = lu->val; 434397b6df1SKris Buschelman #endif 435397b6df1SKris Buschelman } 436397b6df1SKris Buschelman } 437397b6df1SKris Buschelman break; 438397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 439397b6df1SKris Buschelman lu->id.nz_loc = nnz; 440397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 441397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 442397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 443397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 444397b6df1SKris Buschelman #else 445397b6df1SKris Buschelman lu->id.a_loc = lu->val; 446397b6df1SKris Buschelman #endif 447397b6df1SKris Buschelman } 448329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 449329ec9b3SHong Zhang if (!lu->myid){ 450329ec9b3SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr); 451329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr); 452329ec9b3SHong Zhang } else { 453329ec9b3SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 454329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 455329ec9b3SHong Zhang } 4567adad957SLisandro Dalcin ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 457329ec9b3SHong Zhang ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 458329ec9b3SHong Zhang ierr = VecSetFromOptions(b);CHKERRQ(ierr); 459329ec9b3SHong Zhang 460329ec9b3SHong Zhang ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 461329ec9b3SHong Zhang ierr = ISDestroy(is_iden);CHKERRQ(ierr); 462329ec9b3SHong Zhang ierr = VecDestroy(b);CHKERRQ(ierr); 463397b6df1SKris Buschelman break; 464397b6df1SKris Buschelman } 465397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 466397b6df1SKris Buschelman zmumps_c(&lu->id); 467397b6df1SKris Buschelman #else 468397b6df1SKris Buschelman dmumps_c(&lu->id); 469397b6df1SKris Buschelman #endif 470397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 47179a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 472397b6df1SKris Buschelman } 473397b6df1SKris Buschelman } 474397b6df1SKris Buschelman 475397b6df1SKris Buschelman /* numerical factorization phase */ 476329ec9b3SHong Zhang /*-------------------------------*/ 477329ec9b3SHong Zhang lu->id.job = 2; 478958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 479a7aca84bSHong Zhang if (!lu->myid) { 480397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 481397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 482397b6df1SKris Buschelman #else 483397b6df1SKris Buschelman lu->id.a = lu->val; 484397b6df1SKris Buschelman #endif 485397b6df1SKris Buschelman } 486397b6df1SKris Buschelman } else { 487397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 488397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 489397b6df1SKris Buschelman #else 490397b6df1SKris Buschelman lu->id.a_loc = lu->val; 491397b6df1SKris Buschelman #endif 492397b6df1SKris Buschelman } 493397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 494397b6df1SKris Buschelman zmumps_c(&lu->id); 495397b6df1SKris Buschelman #else 496397b6df1SKris Buschelman dmumps_c(&lu->id); 497397b6df1SKris Buschelman #endif 498397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 49919facb7aSBarry Smith if (lu->id.INFO(1) == -13) { 50019facb7aSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 50119facb7aSBarry Smith } else { 50279a5c55eSBarry Smith SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2)); 503397b6df1SKris Buschelman } 50419facb7aSBarry Smith } 505397b6df1SKris Buschelman 50619facb7aSBarry Smith if (!lu->myid && lu->id.ICNTL(16) > 0){ 50779a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 508397b6df1SKris Buschelman } 509397b6df1SKris Buschelman 5108ada1bb4SHong Zhang if (lu->size > 1){ 511*5c9eb25fSBarry Smith if ((*F)->factor == MAT_FACTOR_LU){ 512e09efc27SHong Zhang F_diag = ((Mat_MPIAIJ *)(*F)->data)->A; 513e09efc27SHong Zhang } else { 514e09efc27SHong Zhang F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A; 515e09efc27SHong Zhang } 516e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 517329ec9b3SHong Zhang if (lu->nSolve){ 518329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 519329ec9b3SHong Zhang ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 520329ec9b3SHong Zhang ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 521329ec9b3SHong Zhang } 5228ada1bb4SHong Zhang } 523397b6df1SKris Buschelman (*F)->assembled = PETSC_TRUE; 524397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 525ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 526329ec9b3SHong Zhang lu->nSolve = 0; 527397b6df1SKris Buschelman PetscFunctionReturn(0); 528397b6df1SKris Buschelman } 529397b6df1SKris Buschelman 530b24902e0SBarry Smith 531397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 532397b6df1SKris Buschelman #undef __FUNCT__ 533f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 534b24902e0SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) 535b24902e0SBarry Smith { 536f0c56d0fSKris Buschelman Mat_MUMPS *lu; 537dfbe8321SBarry Smith PetscErrorCode ierr; 538397b6df1SKris Buschelman 539397b6df1SKris Buschelman PetscFunctionBegin; 540b24902e0SBarry Smith 541b24902e0SBarry Smith /* Create the factorization matrix */ 542b24902e0SBarry Smith lu = (Mat_MUMPS*)(*F)->spptr; 543b24902e0SBarry Smith lu->sym = 0; 544b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 545b24902e0SBarry Smith PetscFunctionReturn(0); 546b24902e0SBarry Smith } 547b24902e0SBarry Smith 548*5c9eb25fSBarry Smith EXTERN_C_BEGIN 549*5c9eb25fSBarry Smith /* 550*5c9eb25fSBarry Smith The seq and mpi versions of this function are the same 551*5c9eb25fSBarry Smith */ 552b24902e0SBarry Smith #undef __FUNCT__ 553*5c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_mumps" 554*5c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F) 555b24902e0SBarry Smith { 556b24902e0SBarry Smith Mat B; 557b24902e0SBarry Smith PetscErrorCode ierr; 558b24902e0SBarry Smith Mat_MUMPS *mumps; 559b24902e0SBarry Smith 560b24902e0SBarry Smith PetscFunctionBegin; 561*5c9eb25fSBarry Smith if (ftype != MAT_FACTOR_LU) { 562*5c9eb25fSBarry Smith SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix"); 563*5c9eb25fSBarry Smith } 564397b6df1SKris Buschelman /* Create the factorization matrix */ 5657adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 5662a4c71feSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 5677adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 568397b6df1SKris Buschelman ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 569397b6df1SKris Buschelman 570f6c57405SHong Zhang B->ops->lufactornumeric = MatFactorNumeric_MUMPS; 571b24902e0SBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 572*5c9eb25fSBarry Smith B->factor = MAT_FACTOR_LU; 573b24902e0SBarry Smith 574b24902e0SBarry Smith ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 575b24902e0SBarry Smith mumps->CleanUpMUMPS = PETSC_FALSE; 576b24902e0SBarry Smith mumps->isAIJ = PETSC_TRUE; 577b24902e0SBarry Smith mumps->scat_rhs = PETSC_NULL; 578b24902e0SBarry Smith mumps->scat_sol = PETSC_NULL; 579b24902e0SBarry Smith mumps->nSolve = 0; 580b24902e0SBarry Smith 581b24902e0SBarry Smith B->spptr = (void*)mumps; 582397b6df1SKris Buschelman 583397b6df1SKris Buschelman *F = B; 584397b6df1SKris Buschelman PetscFunctionReturn(0); 585397b6df1SKris Buschelman } 586*5c9eb25fSBarry Smith EXTERN_C_END 587397b6df1SKris Buschelman 588*5c9eb25fSBarry Smith EXTERN_C_BEGIN 589*5c9eb25fSBarry Smith #undef __FUNCT__ 590*5c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_mpiaij_mumps" 591*5c9eb25fSBarry Smith PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F) 592*5c9eb25fSBarry Smith { 593*5c9eb25fSBarry Smith Mat B; 594*5c9eb25fSBarry Smith PetscErrorCode ierr; 595*5c9eb25fSBarry Smith Mat_MUMPS *mumps; 596*5c9eb25fSBarry Smith 597*5c9eb25fSBarry Smith PetscFunctionBegin; 598*5c9eb25fSBarry Smith if (ftype != MAT_FACTOR_LU) { 599*5c9eb25fSBarry Smith SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix"); 600*5c9eb25fSBarry Smith } 601*5c9eb25fSBarry Smith /* Create the factorization matrix */ 602*5c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 603*5c9eb25fSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 604*5c9eb25fSBarry Smith ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 605*5c9eb25fSBarry Smith ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 606*5c9eb25fSBarry Smith ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 607*5c9eb25fSBarry Smith 608*5c9eb25fSBarry Smith B->ops->lufactornumeric = MatFactorNumeric_MUMPS; 609*5c9eb25fSBarry Smith B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 610*5c9eb25fSBarry Smith B->factor = MAT_FACTOR_LU; 611*5c9eb25fSBarry Smith 612*5c9eb25fSBarry Smith ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 613*5c9eb25fSBarry Smith mumps->CleanUpMUMPS = PETSC_FALSE; 614*5c9eb25fSBarry Smith mumps->isAIJ = PETSC_TRUE; 615*5c9eb25fSBarry Smith mumps->scat_rhs = PETSC_NULL; 616*5c9eb25fSBarry Smith mumps->scat_sol = PETSC_NULL; 617*5c9eb25fSBarry Smith mumps->nSolve = 0; 618*5c9eb25fSBarry Smith 619*5c9eb25fSBarry Smith B->spptr = (void*)mumps; 620*5c9eb25fSBarry Smith 621*5c9eb25fSBarry Smith *F = B; 622*5c9eb25fSBarry Smith PetscFunctionReturn(0); 623*5c9eb25fSBarry Smith } 624*5c9eb25fSBarry Smith EXTERN_C_END 625b24902e0SBarry Smith 626397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 627397b6df1SKris Buschelman #undef __FUNCT__ 628f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 629b24902e0SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) 630b24902e0SBarry Smith { 631f0c56d0fSKris Buschelman Mat_MUMPS *lu; 632dfbe8321SBarry Smith PetscErrorCode ierr; 633397b6df1SKris Buschelman 634397b6df1SKris Buschelman PetscFunctionBegin; 635b24902e0SBarry Smith lu = (Mat_MUMPS*)(*F)->spptr; 636b24902e0SBarry Smith lu->sym = 2; 637b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 638b24902e0SBarry Smith PetscFunctionReturn(0); 639b24902e0SBarry Smith } 640b24902e0SBarry Smith 641*5c9eb25fSBarry Smith EXTERN_C_BEGIN 642b24902e0SBarry Smith #undef __FUNCT__ 643*5c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_seqsbaij_mumps" 644*5c9eb25fSBarry Smith PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 645b24902e0SBarry Smith { 646b24902e0SBarry Smith Mat B; 647b24902e0SBarry Smith PetscErrorCode ierr; 648b24902e0SBarry Smith Mat_MUMPS *mumps; 649b24902e0SBarry Smith 650b24902e0SBarry Smith PetscFunctionBegin; 651*5c9eb25fSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) { 652*5c9eb25fSBarry Smith SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 653*5c9eb25fSBarry Smith } 654397b6df1SKris Buschelman /* Create the factorization matrix */ 6557adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 6562a4c71feSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 6577adad957SLisandro Dalcin ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 658efc670deSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 659efc670deSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 660397b6df1SKris Buschelman 661b24902e0SBarry Smith B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 662f6c57405SHong Zhang B->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 663a58c3f20SHong Zhang B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 664*5c9eb25fSBarry Smith B->factor = MAT_FACTOR_CHOLESKY; 665397b6df1SKris Buschelman 666b24902e0SBarry Smith ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 667b24902e0SBarry Smith mumps->CleanUpMUMPS = PETSC_FALSE; 668b24902e0SBarry Smith mumps->isAIJ = PETSC_TRUE; 669b24902e0SBarry Smith mumps->scat_rhs = PETSC_NULL; 670b24902e0SBarry Smith mumps->scat_sol = PETSC_NULL; 671b24902e0SBarry Smith mumps->nSolve = 0; 672b24902e0SBarry Smith B->spptr = (void*)mumps; 673397b6df1SKris Buschelman *F = B; 674397b6df1SKris Buschelman PetscFunctionReturn(0); 675397b6df1SKris Buschelman } 676*5c9eb25fSBarry Smith EXTERN_C_END 677*5c9eb25fSBarry Smith 678*5c9eb25fSBarry Smith EXTERN_C_BEGIN 679*5c9eb25fSBarry Smith #undef __FUNCT__ 680*5c9eb25fSBarry Smith #define __FUNCT__ "MatGetFactor_mpisbaij_mumps" 681*5c9eb25fSBarry Smith PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 682*5c9eb25fSBarry Smith { 683*5c9eb25fSBarry Smith Mat B; 684*5c9eb25fSBarry Smith PetscErrorCode ierr; 685*5c9eb25fSBarry Smith Mat_MUMPS *mumps; 686*5c9eb25fSBarry Smith 687*5c9eb25fSBarry Smith PetscFunctionBegin; 688*5c9eb25fSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) { 689*5c9eb25fSBarry Smith SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 690*5c9eb25fSBarry Smith } 691*5c9eb25fSBarry Smith /* Create the factorization matrix */ 692*5c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 693*5c9eb25fSBarry Smith ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr); 694*5c9eb25fSBarry Smith ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 695*5c9eb25fSBarry Smith ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 696*5c9eb25fSBarry Smith ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 697*5c9eb25fSBarry Smith 698*5c9eb25fSBarry Smith B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 699*5c9eb25fSBarry Smith B->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 700*5c9eb25fSBarry Smith B->ops->getinertia = MatGetInertia_SBAIJMUMPS; 701*5c9eb25fSBarry Smith B->factor = MAT_FACTOR_CHOLESKY; 702*5c9eb25fSBarry Smith 703*5c9eb25fSBarry Smith ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 704*5c9eb25fSBarry Smith mumps->CleanUpMUMPS = PETSC_FALSE; 705*5c9eb25fSBarry Smith mumps->isAIJ = PETSC_TRUE; 706*5c9eb25fSBarry Smith mumps->scat_rhs = PETSC_NULL; 707*5c9eb25fSBarry Smith mumps->scat_sol = PETSC_NULL; 708*5c9eb25fSBarry Smith mumps->nSolve = 0; 709*5c9eb25fSBarry Smith B->spptr = (void*)mumps; 710*5c9eb25fSBarry Smith *F = B; 711*5c9eb25fSBarry Smith PetscFunctionReturn(0); 712*5c9eb25fSBarry Smith } 713*5c9eb25fSBarry Smith EXTERN_C_END 714397b6df1SKris Buschelman 715397b6df1SKris Buschelman #undef __FUNCT__ 716f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS" 717f6c57405SHong Zhang PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 718f6c57405SHong Zhang Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 719f6c57405SHong Zhang PetscErrorCode ierr; 720f6c57405SHong Zhang 721f6c57405SHong Zhang PetscFunctionBegin; 722f6c57405SHong Zhang /* check if matrix is mumps type */ 723f6c57405SHong Zhang if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 724f6c57405SHong Zhang 725f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 726f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 727f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 728f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 729f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 730f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 731f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 732f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 733f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 734f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 735f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 736f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 737f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 738f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 739f6c57405SHong Zhang if (!lu->myid && lu->id.ICNTL(11)>0) { 740f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 741f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 742f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 743f6c57405SHong Zhang 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); 744f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 745f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr); 746f6c57405SHong Zhang 747f6c57405SHong Zhang } 748f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 749f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 750f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 751f6c57405SHong Zhang /* ICNTL(15-17) not used */ 752f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 753f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 754f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 755f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 756f6c57405SHong Zhang 757f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 758f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 759f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 760f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 761f6c57405SHong Zhang 762f6c57405SHong Zhang /* infomation local to each processor */ 763f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);} 7647adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 7657adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 766f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);} 7677adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 7687adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 769f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);} 7707adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 7717adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 772f6c57405SHong Zhang /* 773f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);} 7747adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr); 7757adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 776f6c57405SHong Zhang */ 777f6c57405SHong Zhang 778f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);} 7797adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 7807adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 781f6c57405SHong Zhang 782f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);} 7837adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 7847adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 785f6c57405SHong Zhang 786f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);} 7877adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 7887adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 789f6c57405SHong Zhang 790f6c57405SHong Zhang if (!lu->myid){ /* information from the host */ 791f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 792f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 793f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 794f6c57405SHong Zhang 795f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 796f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 797f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 798f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 799f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 800f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 801f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr); 802f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 803f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 804f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 805f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 806f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 807f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 808f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr); 809f6c57405SHong 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); 810f6c57405SHong 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); 811f6c57405SHong 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); 812f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 813f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr); 814f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr); 815f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 816f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 817f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 818f6c57405SHong Zhang } 819f6c57405SHong Zhang 820f6c57405SHong Zhang PetscFunctionReturn(0); 821f6c57405SHong Zhang } 822f6c57405SHong Zhang 823f6c57405SHong Zhang #undef __FUNCT__ 824f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS" 825b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 826b24902e0SBarry Smith { 827f6c57405SHong Zhang PetscErrorCode ierr; 828f6c57405SHong Zhang PetscTruth iascii; 829f6c57405SHong Zhang PetscViewerFormat format; 830f6c57405SHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)(A->spptr); 831f6c57405SHong Zhang 832f6c57405SHong Zhang PetscFunctionBegin; 833f6c57405SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 834f6c57405SHong Zhang if (iascii) { 835f6c57405SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 836f6c57405SHong Zhang if (format == PETSC_VIEWER_ASCII_INFO){ 837f6c57405SHong Zhang ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 838f6c57405SHong Zhang } 839f6c57405SHong Zhang } 840f6c57405SHong Zhang PetscFunctionReturn(0); 841f6c57405SHong Zhang } 842f6c57405SHong Zhang 843397b6df1SKris Buschelman 84424b6179bSKris Buschelman /*MC 845fafad747SKris Buschelman MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 84624b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 84724b6179bSKris Buschelman 84824b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 84924b6179bSKris Buschelman on how to declare the existence of external packages), 85024b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 851175b88e8SBarry Smith After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS), then 852175b88e8SBarry Smith optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT 853175b88e8SBarry Smith call MatCreateSeqAIJ/MPIAIJ() directly or the preallocation information will be LOST! 85424b6179bSKris Buschelman 85524b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 85624b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 8573ec795f1SBarry Smith MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported 85824b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 8593ec795f1SBarry Smith the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 86028b08bd3SKris Buschelman conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 861175b88e8SBarry Smith without data copy AFTER the matrix values are set. 86224b6179bSKris Buschelman 86324b6179bSKris Buschelman Options Database Keys: 8640bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 86524b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 86624b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 86724b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 86824b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 86924b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 87024b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 87194b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 87224b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 87324b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 87424b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 87524b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 87624b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 87724b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 87824b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 87924b6179bSKris Buschelman 88024b6179bSKris Buschelman Level: beginner 88124b6179bSKris Buschelman 88224b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 88324b6179bSKris Buschelman M*/ 88424b6179bSKris Buschelman 885f0c56d0fSKris Buschelman 88624b6179bSKris Buschelman /*MC 887fafad747SKris Buschelman MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 88824b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 88924b6179bSKris Buschelman 89024b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 89124b6179bSKris Buschelman on how to declare the existence of external packages), 89224b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 893175b88e8SBarry Smith After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS), then 894175b88e8SBarry Smith optionally call MatSeqSBAIJSetPreallocation() or MatMPISBAIJSetPreallocation() DO NOT 895175b88e8SBarry Smith call MatCreateSeqSBAIJ/MPISBAIJ() directly or the preallocation information will be LOST! 89624b6179bSKris Buschelman 89724b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 89824b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 899175b88e8SBarry Smith MatSeqSBAIJSetPreallocation() is supported, and similarly MatMPISBAIJSetPreallocation() is supported 90024b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 901175b88e8SBarry Smith the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 90228b08bd3SKris Buschelman conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 903175b88e8SBarry Smith without data copy AFTER the matrix values have been set. 90424b6179bSKris Buschelman 90524b6179bSKris Buschelman Options Database Keys: 9060bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 90724b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 90824b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 90924b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 91024b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 91124b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 91224b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 91394b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 91424b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 91524b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 91624b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 91724b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 91824b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 91924b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 92024b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 92124b6179bSKris Buschelman 92224b6179bSKris Buschelman Level: beginner 92324b6179bSKris Buschelman 92424b6179bSKris Buschelman .seealso: MATAIJMUMPS 92524b6179bSKris Buschelman M*/ 92624b6179bSKris Buschelman 927