1be1d678aSKris Buschelman #define PETSCMAT_DLL 21c2a3de1SBarry Smith 3397b6df1SKris Buschelman /* 4c2b5dc30SHong Zhang Provides an interface to the MUMPS sparse solver 5397b6df1SKris Buschelman */ 64e1dd20eSHong Zhang #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 77c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 87c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 97c4f633dSBarry Smith #include "../src/mat/impls/sbaij/mpi/mpisbaij.h" 10450b117fSShri Abhyankar #include "../src/mat/impls/baij/seq/baij.h" 11450b117fSShri Abhyankar #include "../src/mat/impls/baij/mpi/mpibaij.h" 12397b6df1SKris Buschelman 13397b6df1SKris Buschelman EXTERN_C_BEGIN 14397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 15397b6df1SKris Buschelman #include "zmumps_c.h" 16397b6df1SKris Buschelman #else 17397b6df1SKris Buschelman #include "dmumps_c.h" 18397b6df1SKris Buschelman #endif 19397b6df1SKris Buschelman EXTERN_C_END 20397b6df1SKris Buschelman #define JOB_INIT -1 21397b6df1SKris Buschelman #define JOB_END -2 22397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */ 23397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1] 24397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1] 25397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1] 26a7aca84bSHong Zhang #define INFO(I) info[(I)-1] 27397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1] 28adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1] 29397b6df1SKris Buschelman 30397b6df1SKris Buschelman typedef struct { 31397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 32397b6df1SKris Buschelman ZMUMPS_STRUC_C id; 33397b6df1SKris Buschelman #else 34397b6df1SKris Buschelman DMUMPS_STRUC_C id; 35397b6df1SKris Buschelman #endif 36397b6df1SKris Buschelman MatStructure matstruc; 37c1490034SHong Zhang PetscMPIInt myid,size; 38329ec9b3SHong Zhang PetscInt *irn,*jcn,sym,nSolve; 39397b6df1SKris Buschelman PetscScalar *val; 40397b6df1SKris Buschelman MPI_Comm comm_mumps; 41329ec9b3SHong Zhang VecScatter scat_rhs, scat_sol; 42cb828f0fSHong Zhang PetscTruth isAIJ,CleanUpMUMPS,mumpsview; 43329ec9b3SHong Zhang Vec b_seq,x_seq; 4467334b25SHong Zhang PetscErrorCode (*MatDestroy)(Mat); 45f0c56d0fSKris Buschelman } Mat_MUMPS; 46f0c56d0fSKris Buschelman 47dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 48b24902e0SBarry Smith 49397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 50397b6df1SKris Buschelman /* 51397b6df1SKris Buschelman input: 5275747be1SHong Zhang A - matrix in mpiaij or mpisbaij (bs=1) format 53397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 54397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 55397b6df1SKris Buschelman TRUE: only the values in v array are updated 56397b6df1SKris Buschelman output: 57397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 58397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 59397b6df1SKris Buschelman */ 60b24902e0SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) 61b24902e0SBarry Smith { 62c1490034SHong Zhang PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray; 63dfbe8321SBarry Smith PetscErrorCode ierr; 64d0f46423SBarry Smith PetscInt i,j,jj,jB,irow,m=A->rmap->n,*ajj,*bjj,countA,countB,colA_start,jcol; 65c1490034SHong Zhang PetscInt *row,*col; 66397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 675c9eb25fSBarry Smith PetscTruth isAIJ; 68397b6df1SKris Buschelman 69397b6df1SKris Buschelman PetscFunctionBegin; 705c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr); 715c9eb25fSBarry Smith if (isAIJ){ 72397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 73397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 74397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 75397b6df1SKris Buschelman nz = aa->nz + bb->nz; 76d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 77397b6df1SKris Buschelman garray = mat->garray; 78397b6df1SKris Buschelman av=aa->a; bv=bb->a; 79397b6df1SKris Buschelman 80397b6df1SKris Buschelman } else { 81397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 82397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 83397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 84e32f2f54SBarry Smith if (A->rmap->bs > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap->bs); 856c6c5352SBarry Smith nz = aa->nz + bb->nz; 86d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 87397b6df1SKris Buschelman garray = mat->garray; 88397b6df1SKris Buschelman av=aa->a; bv=bb->a; 89397b6df1SKris Buschelman } 90397b6df1SKris Buschelman 91397b6df1SKris Buschelman if (!valOnly){ 927c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr); 937c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr); 94397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 95397b6df1SKris Buschelman *r = row; *c = col; *v = val; 96397b6df1SKris Buschelman } else { 97397b6df1SKris Buschelman row = *r; col = *c; val = *v; 98397b6df1SKris Buschelman } 99397b6df1SKris Buschelman *nnz = nz; 100397b6df1SKris Buschelman 101028e57e8SHong Zhang jj = 0; irow = rstart; 102397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 103397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 104397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 105397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 106397b6df1SKris Buschelman bjj = bj + bi[i]; 107397b6df1SKris Buschelman 108397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 109397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 11075747be1SHong Zhang j=-1; 11175747be1SHong Zhang do { 11275747be1SHong Zhang j++; 11375747be1SHong Zhang if (j == countB) break; 114397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11575747be1SHong Zhang } while (jcol < colA_start); 11675747be1SHong Zhang jB = j; 117397b6df1SKris Buschelman 118397b6df1SKris Buschelman /* B-part, smaller col index */ 119397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 120397b6df1SKris Buschelman for (j=0; j<jB; j++){ 121397b6df1SKris Buschelman jcol = garray[bjj[j]]; 122397b6df1SKris Buschelman if (!valOnly){ 123397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 12475747be1SHong Zhang 125397b6df1SKris Buschelman } 126397b6df1SKris Buschelman val[jj++] = *bv++; 127397b6df1SKris Buschelman } 128397b6df1SKris Buschelman /* A-part */ 129397b6df1SKris Buschelman for (j=0; j<countA; j++){ 130397b6df1SKris Buschelman if (!valOnly){ 131397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 132397b6df1SKris Buschelman } 133397b6df1SKris Buschelman val[jj++] = *av++; 134397b6df1SKris Buschelman } 135397b6df1SKris Buschelman /* B-part, larger col index */ 136397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 137397b6df1SKris Buschelman if (!valOnly){ 138397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 139397b6df1SKris Buschelman } 140397b6df1SKris Buschelman val[jj++] = *bv++; 141397b6df1SKris Buschelman } 142397b6df1SKris Buschelman irow++; 143397b6df1SKris Buschelman } 144397b6df1SKris Buschelman PetscFunctionReturn(0); 145397b6df1SKris Buschelman } 146397b6df1SKris Buschelman 147397b6df1SKris Buschelman #undef __FUNCT__ 1483924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 149dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 150dfbe8321SBarry Smith { 151f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 152dfbe8321SBarry Smith PetscErrorCode ierr; 153c1490034SHong Zhang PetscMPIInt size=lu->size; 154b24902e0SBarry Smith 155397b6df1SKris Buschelman PetscFunctionBegin; 156397b6df1SKris Buschelman if (lu->CleanUpMUMPS) { 157397b6df1SKris Buschelman /* Terminate instance, deallocate memories */ 158329ec9b3SHong Zhang if (size > 1){ 15968653410SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 160329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 161329ec9b3SHong Zhang ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 1622750af12SHong Zhang if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);} 1632750af12SHong Zhang if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);} 164329ec9b3SHong Zhang ierr = PetscFree(lu->val);CHKERRQ(ierr); 165329ec9b3SHong Zhang } 166450b117fSShri Abhyankar if( size == 1 && (A)->factortype == MAT_FACTOR_CHOLESKY && lu->isAIJ) { 167450b117fSShri Abhyankar ierr = PetscFree(lu->val);CHKERRQ(ierr); 168450b117fSShri Abhyankar } 169397b6df1SKris Buschelman lu->id.job=JOB_END; 170397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 171397b6df1SKris Buschelman zmumps_c(&lu->id); 172397b6df1SKris Buschelman #else 173397b6df1SKris Buschelman dmumps_c(&lu->id); 174397b6df1SKris Buschelman #endif 175c338a77dSKris Buschelman ierr = PetscFree(lu->irn);CHKERRQ(ierr); 176c338a77dSKris Buschelman ierr = PetscFree(lu->jcn);CHKERRQ(ierr); 177397b6df1SKris Buschelman ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr); 178397b6df1SKris Buschelman } 17997969023SHong Zhang /* clear composed functions */ 18097969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr); 18197969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr); 18267334b25SHong Zhang ierr = (lu->MatDestroy)(A);CHKERRQ(ierr); 183397b6df1SKris Buschelman PetscFunctionReturn(0); 184397b6df1SKris Buschelman } 185397b6df1SKris Buschelman 186397b6df1SKris Buschelman #undef __FUNCT__ 187f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS" 188b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) 189b24902e0SBarry Smith { 190f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 191d54de34fSKris Buschelman PetscScalar *array; 192397b6df1SKris Buschelman Vec x_seq; 193329ec9b3SHong Zhang IS is_iden,is_petsc; 194dfbe8321SBarry Smith PetscErrorCode ierr; 195329ec9b3SHong Zhang PetscInt i; 196397b6df1SKris Buschelman 197397b6df1SKris Buschelman PetscFunctionBegin; 198329ec9b3SHong Zhang lu->id.nrhs = 1; 199329ec9b3SHong Zhang x_seq = lu->b_seq; 200397b6df1SKris Buschelman if (lu->size > 1){ 201329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */ 202f6cfb2d1SLisandro Dalcin ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 203f6cfb2d1SLisandro Dalcin ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 204397b6df1SKris Buschelman if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);} 205397b6df1SKris Buschelman } else { /* size == 1 */ 206397b6df1SKris Buschelman ierr = VecCopy(b,x);CHKERRQ(ierr); 207397b6df1SKris Buschelman ierr = VecGetArray(x,&array);CHKERRQ(ierr); 208397b6df1SKris Buschelman } 209397b6df1SKris Buschelman if (!lu->myid) { /* define rhs on the host */ 2108278f211SHong Zhang lu->id.nrhs = 1; 211397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 212397b6df1SKris Buschelman lu->id.rhs = (mumps_double_complex*)array; 213397b6df1SKris Buschelman #else 214397b6df1SKris Buschelman lu->id.rhs = array; 215397b6df1SKris Buschelman #endif 216397b6df1SKris Buschelman } 217329ec9b3SHong Zhang if (lu->size == 1){ 218329ec9b3SHong Zhang ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 219329ec9b3SHong Zhang } else if (!lu->myid){ 220329ec9b3SHong Zhang ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 221329ec9b3SHong Zhang } 222329ec9b3SHong Zhang 223329ec9b3SHong Zhang if (lu->size > 1){ 224329ec9b3SHong Zhang /* distributed solution */ 225329ec9b3SHong Zhang lu->id.ICNTL(21) = 1; 226329ec9b3SHong Zhang if (!lu->nSolve){ 227329ec9b3SHong Zhang /* Create x_seq=sol_loc for repeated use */ 228329ec9b3SHong Zhang PetscInt lsol_loc; 229329ec9b3SHong Zhang PetscScalar *sol_loc; 230329ec9b3SHong Zhang lsol_loc = lu->id.INFO(23); /* length of sol_loc */ 2310e83c824SBarry Smith ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr); 232329ec9b3SHong Zhang lu->id.lsol_loc = lsol_loc; 23344ea04b1SSatish Balay #if defined(PETSC_USE_COMPLEX) 23444ea04b1SSatish Balay lu->id.sol_loc = (mumps_double_complex*)sol_loc; 23544ea04b1SSatish Balay #else 23644ea04b1SSatish Balay lu->id.sol_loc = sol_loc; 23744ea04b1SSatish Balay #endif 238329ec9b3SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr); 239329ec9b3SHong Zhang } 240329ec9b3SHong Zhang } 241397b6df1SKris Buschelman 242397b6df1SKris Buschelman /* solve phase */ 243329ec9b3SHong Zhang /*-------------*/ 244397b6df1SKris Buschelman lu->id.job = 3; 245397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 246397b6df1SKris Buschelman zmumps_c(&lu->id); 247397b6df1SKris Buschelman #else 248397b6df1SKris Buschelman dmumps_c(&lu->id); 249397b6df1SKris Buschelman #endif 250397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 251e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 252397b6df1SKris Buschelman } 253397b6df1SKris Buschelman 254329ec9b3SHong Zhang if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */ 255329ec9b3SHong Zhang if (!lu->nSolve){ /* create scatter scat_sol */ 256329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */ 257329ec9b3SHong Zhang for (i=0; i<lu->id.lsol_loc; i++){ 258329ec9b3SHong Zhang lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */ 259397b6df1SKris Buschelman } 260329ec9b3SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr); /* to */ 261329ec9b3SHong Zhang ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr); 262329ec9b3SHong Zhang ierr = ISDestroy(is_iden);CHKERRQ(ierr); 263329ec9b3SHong Zhang ierr = ISDestroy(is_petsc);CHKERRQ(ierr); 264397b6df1SKris Buschelman } 265ca9f406cSSatish Balay ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 266ca9f406cSSatish Balay ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 267329ec9b3SHong Zhang } 268329ec9b3SHong Zhang lu->nSolve++; 269397b6df1SKris Buschelman PetscFunctionReturn(0); 270397b6df1SKris Buschelman } 271397b6df1SKris Buschelman 272ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 273a58c3f20SHong Zhang /* 274a58c3f20SHong Zhang input: 275a58c3f20SHong Zhang F: numeric factor 276a58c3f20SHong Zhang output: 277a58c3f20SHong Zhang nneg: total number of negative pivots 278a58c3f20SHong Zhang nzero: 0 279a58c3f20SHong Zhang npos: (global dimension of F) - nneg 280a58c3f20SHong Zhang */ 281a58c3f20SHong Zhang 282a58c3f20SHong Zhang #undef __FUNCT__ 283a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 284dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 285a58c3f20SHong Zhang { 286a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 287dfbe8321SBarry Smith PetscErrorCode ierr; 288c1490034SHong Zhang PetscMPIInt size; 289a58c3f20SHong Zhang 290a58c3f20SHong Zhang PetscFunctionBegin; 2917adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 292bcb30aebSHong 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 */ 293bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 294e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,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)); 295bcb30aebSHong Zhang } 296a58c3f20SHong Zhang if (nneg){ 297a58c3f20SHong Zhang if (!lu->myid){ 298a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 299a58c3f20SHong Zhang } 300bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 301a58c3f20SHong Zhang } 302a58c3f20SHong Zhang if (nzero) *nzero = 0; 303d0f46423SBarry Smith if (npos) *npos = F->rmap->N - (*nneg); 304a58c3f20SHong Zhang PetscFunctionReturn(0); 305a58c3f20SHong Zhang } 306ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 307a58c3f20SHong Zhang 308397b6df1SKris Buschelman #undef __FUNCT__ 309f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 3100481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 311af281ebdSHong Zhang { 312dcd589f8SShri Abhyankar Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 313450b117fSShri Abhyankar Mat newMat; 3146849ba73SBarry Smith PetscErrorCode ierr; 315dcd589f8SShri Abhyankar PetscInt rnz,nnz,nz=0,i,M=A->rmap->N,*ai,*aj,*adiag,jidx; 316450b117fSShri Abhyankar PetscScalar *av; 317dcd589f8SShri Abhyankar PetscTruth valOnly; 318e09efc27SHong Zhang Mat F_diag; 319c349612cSHong Zhang IS is_iden; 320c349612cSHong Zhang Vec b; 321450b117fSShri Abhyankar PetscTruth isSeqAIJ,isSeqSBAIJ,isMPIAIJ; 322397b6df1SKris Buschelman 323397b6df1SKris Buschelman PetscFunctionBegin; 3245c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 3255c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 326a214ac2aSShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 327a214ac2aSShri Abhyankar 328397b6df1SKris Buschelman /* define matrix A */ 329397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 330397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 331397b6df1SKris Buschelman if (!lu->myid) { 3325c9eb25fSBarry Smith if (isSeqAIJ){ 333397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 334397b6df1SKris Buschelman nz = aa->nz; 335450b117fSShri Abhyankar ai = aa->i; aj = aa->j; adiag = aa->diag; av = aa->a; 3365c9eb25fSBarry Smith } else if (isSeqSBAIJ) { 337397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 3386c6c5352SBarry Smith nz = aa->nz; 339450b117fSShri Abhyankar ai = aa->i; aj = aa->j; av = aa->a; 340*e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No mumps factorization for this matrix type"); 341*e7e72b3dSBarry Smith 342dcd589f8SShri Abhyankar if ((F)->factortype == MAT_FACTOR_CHOLESKY && isSeqAIJ) { 343450b117fSShri Abhyankar nz = 0; 344450b117fSShri Abhyankar for (i=0; i<M; i++){ 345450b117fSShri Abhyankar rnz = ai[i+1] - adiag[i]; 346450b117fSShri Abhyankar jidx = adiag[i]; 347450b117fSShri Abhyankar while (rnz--) { /* Fortran row/col index! */ 348450b117fSShri Abhyankar lu->val[nz] = av[jidx]; jidx++; nz++; 349450b117fSShri Abhyankar } 350450b117fSShri Abhyankar } 351450b117fSShri Abhyankar } else { 352450b117fSShri Abhyankar lu->val = av; 353397b6df1SKris Buschelman } 354450b117fSShri Abhyankar } 355397b6df1SKris Buschelman break; 356397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 357397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 358dcd589f8SShri Abhyankar 359dcd589f8SShri Abhyankar if(((F)->factortype == MAT_FACTOR_CHOLESKY) && isMPIAIJ) { 360a214ac2aSShri Abhyankar /* Create an SBAIJ matrix and use this matrix to set the lu values */ 361a214ac2aSShri Abhyankar ierr = MatConvert(A,MATMPISBAIJ,MAT_INITIAL_MATRIX,&newMat);CHKERRQ(ierr); 362a214ac2aSShri Abhyankar ierr = MatConvertToTriples(newMat,1,valOnly,&nnz,&lu->irn , &lu->jcn, &lu->val);CHKERRQ(ierr); 363a214ac2aSShri Abhyankar ierr = MatDestroy(newMat);CHKERRQ(ierr); 364a214ac2aSShri Abhyankar } 365a214ac2aSShri Abhyankar else { 366397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 367a214ac2aSShri Abhyankar } 368397b6df1SKris Buschelman break; 369e32f2f54SBarry Smith default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 370397b6df1SKris Buschelman } 371397b6df1SKris Buschelman 372397b6df1SKris Buschelman /* analysis phase */ 373329ec9b3SHong Zhang /*----------------*/ 374397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 375329ec9b3SHong Zhang lu->id.job = 1; 376329ec9b3SHong Zhang 377397b6df1SKris Buschelman lu->id.n = M; 378397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 379397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 380397b6df1SKris Buschelman if (!lu->myid) { 381397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 382397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 383397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 384397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 385397b6df1SKris Buschelman #else 386397b6df1SKris Buschelman lu->id.a = lu->val; 387397b6df1SKris Buschelman #endif 388397b6df1SKris Buschelman } 389397b6df1SKris Buschelman } 390397b6df1SKris Buschelman break; 391397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 392397b6df1SKris Buschelman lu->id.nz_loc = nnz; 393397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 394397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 395397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 396397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 397397b6df1SKris Buschelman #else 398397b6df1SKris Buschelman lu->id.a_loc = lu->val; 399397b6df1SKris Buschelman #endif 400397b6df1SKris Buschelman } 401329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 402329ec9b3SHong Zhang if (!lu->myid){ 403d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 404d0f46423SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 405329ec9b3SHong Zhang } else { 406329ec9b3SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 407329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 408329ec9b3SHong Zhang } 4097adad957SLisandro Dalcin ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 410d0f46423SBarry Smith ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 411329ec9b3SHong Zhang ierr = VecSetFromOptions(b);CHKERRQ(ierr); 412329ec9b3SHong Zhang 413329ec9b3SHong Zhang ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 414329ec9b3SHong Zhang ierr = ISDestroy(is_iden);CHKERRQ(ierr); 415329ec9b3SHong Zhang ierr = VecDestroy(b);CHKERRQ(ierr); 416397b6df1SKris Buschelman break; 417397b6df1SKris Buschelman } 418397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 419397b6df1SKris Buschelman zmumps_c(&lu->id); 420397b6df1SKris Buschelman #else 421397b6df1SKris Buschelman dmumps_c(&lu->id); 422397b6df1SKris Buschelman #endif 423397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 424e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 425397b6df1SKris Buschelman } 426397b6df1SKris Buschelman } 427397b6df1SKris Buschelman 428397b6df1SKris Buschelman /* numerical factorization phase */ 429329ec9b3SHong Zhang /*-------------------------------*/ 430329ec9b3SHong Zhang lu->id.job = 2; 431958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 432a7aca84bSHong Zhang if (!lu->myid) { 433397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 434397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 435397b6df1SKris Buschelman #else 436397b6df1SKris Buschelman lu->id.a = lu->val; 437397b6df1SKris Buschelman #endif 438397b6df1SKris Buschelman } 439397b6df1SKris Buschelman } else { 440397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 441397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 442397b6df1SKris Buschelman #else 443397b6df1SKris Buschelman lu->id.a_loc = lu->val; 444397b6df1SKris Buschelman #endif 445397b6df1SKris Buschelman } 446397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 447397b6df1SKris Buschelman zmumps_c(&lu->id); 448397b6df1SKris Buschelman #else 449397b6df1SKris Buschelman dmumps_c(&lu->id); 450397b6df1SKris Buschelman #endif 451397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 45219facb7aSBarry Smith if (lu->id.INFO(1) == -13) { 453e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 45419facb7aSBarry Smith } else { 455e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,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)); 456397b6df1SKris Buschelman } 45719facb7aSBarry Smith } 458397b6df1SKris Buschelman 45919facb7aSBarry Smith if (!lu->myid && lu->id.ICNTL(16) > 0){ 460e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 461397b6df1SKris Buschelman } 462397b6df1SKris Buschelman 4638ada1bb4SHong Zhang if (lu->size > 1){ 464a214ac2aSShri Abhyankar if(isMPIAIJ) { 465dcd589f8SShri Abhyankar F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 466e09efc27SHong Zhang } else { 467dcd589f8SShri Abhyankar F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 468e09efc27SHong Zhang } 469e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 470329ec9b3SHong Zhang if (lu->nSolve){ 471329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 4720e83c824SBarry Smith ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr); 473329ec9b3SHong Zhang ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 474329ec9b3SHong Zhang } 4758ada1bb4SHong Zhang } 476dcd589f8SShri Abhyankar (F)->assembled = PETSC_TRUE; 477397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 478ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 479329ec9b3SHong Zhang lu->nSolve = 0; 480397b6df1SKris Buschelman PetscFunctionReturn(0); 481397b6df1SKris Buschelman } 482397b6df1SKris Buschelman 483dcd589f8SShri Abhyankar #undef __FUNCT__ 484dcd589f8SShri Abhyankar #define __FUNCT__ "PetscSetMUMPSOptions" 485dcd589f8SShri Abhyankar PetscErrorCode PetscSetMUMPSOptions(Mat F, Mat A) 486dcd589f8SShri Abhyankar { 487dcd589f8SShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 488dcd589f8SShri Abhyankar PetscErrorCode ierr; 489dcd589f8SShri Abhyankar PetscInt icntl; 490dcd589f8SShri Abhyankar PetscTruth flg; 491dcd589f8SShri Abhyankar 492dcd589f8SShri Abhyankar PetscFunctionBegin; 493dcd589f8SShri Abhyankar ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 494cb828f0fSHong Zhang ierr = PetscOptionsTruth("-mat_mumps_view","View MUMPS parameters","None",lu->mumpsview,&lu->mumpsview,PETSC_NULL);CHKERRQ(ierr); 495dcd589f8SShri Abhyankar if (lu->size == 1){ 496dcd589f8SShri Abhyankar lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 497dcd589f8SShri Abhyankar } else { 498dcd589f8SShri Abhyankar lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 499dcd589f8SShri Abhyankar } 500dcd589f8SShri Abhyankar 501dcd589f8SShri Abhyankar icntl=-1; 502dcd589f8SShri Abhyankar lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 503dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 504dcd589f8SShri Abhyankar if ((flg && icntl > 0) || PetscLogPrintInfo) { 505dcd589f8SShri Abhyankar lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 506dcd589f8SShri Abhyankar } else { /* no output */ 507dcd589f8SShri Abhyankar lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 508dcd589f8SShri Abhyankar lu->id.ICNTL(2) = 0; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 509dcd589f8SShri Abhyankar lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */ 510dcd589f8SShri Abhyankar } 511dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): column permutation and/or scaling to get a zero-free diagonal (0 to 7)","None",lu->id.ICNTL(6),&lu->id.ICNTL(6),PETSC_NULL);CHKERRQ(ierr); 512dcd589f8SShri Abhyankar icntl=-1; 513dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 514dcd589f8SShri Abhyankar if (flg) { 515dcd589f8SShri Abhyankar if (icntl== 1){ 516e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 517dcd589f8SShri Abhyankar } else { 518dcd589f8SShri Abhyankar lu->id.ICNTL(7) = icntl; 519dcd589f8SShri Abhyankar } 520dcd589f8SShri Abhyankar } 521dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 7 or 77)","None",lu->id.ICNTL(8),&lu->id.ICNTL(8),PETSC_NULL);CHKERRQ(ierr); 522dcd589f8SShri Abhyankar 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); 523dcd589f8SShri Abhyankar 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); 524dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to the linear system solved (via -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr); 525dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control: defines the ordering strategy with scaling constraints (0 to 3","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 526dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control: with or without ScaLAPACK","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 527dcd589f8SShri Abhyankar 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); 528dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr); 529dcd589f8SShri Abhyankar 530dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core facility (0 or 1)","None",lu->id.ICNTL(22),&lu->id.ICNTL(22),PETSC_NULL);CHKERRQ(ierr); 531dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",lu->id.ICNTL(23),&lu->id.ICNTL(23),PETSC_NULL);CHKERRQ(ierr); 532dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",lu->id.ICNTL(24),&lu->id.ICNTL(24),PETSC_NULL);CHKERRQ(ierr); 533dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computation of a null space basis","None",lu->id.ICNTL(25),&lu->id.ICNTL(25),PETSC_NULL);CHKERRQ(ierr); 534dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): Schur options for right-hand side or solution vector","None",lu->id.ICNTL(26),&lu->id.ICNTL(26),PETSC_NULL);CHKERRQ(ierr); 535dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr); 536dcd589f8SShri Abhyankar 537dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",lu->id.CNTL(1),&lu->id.CNTL(1),PETSC_NULL);CHKERRQ(ierr); 538dcd589f8SShri Abhyankar 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); 539dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",lu->id.CNTL(3),&lu->id.CNTL(3),PETSC_NULL);CHKERRQ(ierr); 540dcd589f8SShri Abhyankar 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); 541dcd589f8SShri Abhyankar ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",lu->id.CNTL(5),&lu->id.CNTL(5),PETSC_NULL);CHKERRQ(ierr); 542dcd589f8SShri Abhyankar PetscOptionsEnd(); 543dcd589f8SShri Abhyankar PetscFunctionReturn(0); 544dcd589f8SShri Abhyankar } 545dcd589f8SShri Abhyankar 546dcd589f8SShri Abhyankar #undef __FUNCT__ 547dcd589f8SShri Abhyankar #define __FUNCT__ "PetscInitializeMUMPS" 548dcd589f8SShri Abhyankar PetscErrorCode PetscInitializeMUMPS(Mat F) 549dcd589f8SShri Abhyankar { 550dcd589f8SShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 551dcd589f8SShri Abhyankar PetscErrorCode ierr; 552dcd589f8SShri Abhyankar PetscInt icntl; 553dcd589f8SShri Abhyankar PetscTruth flg; 554dcd589f8SShri Abhyankar 555dcd589f8SShri Abhyankar PetscFunctionBegin; 556dcd589f8SShri Abhyankar lu->id.job = JOB_INIT; 557dcd589f8SShri Abhyankar lu->id.par=1; /* host participates factorizaton and solve */ 558dcd589f8SShri Abhyankar lu->id.sym=lu->sym; 559dcd589f8SShri Abhyankar if (lu->sym == 2){ 560dcd589f8SShri Abhyankar ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 561dcd589f8SShri Abhyankar if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 562dcd589f8SShri Abhyankar } 563dcd589f8SShri Abhyankar #if defined(PETSC_USE_COMPLEX) 564dcd589f8SShri Abhyankar zmumps_c(&lu->id); 565dcd589f8SShri Abhyankar #else 566dcd589f8SShri Abhyankar dmumps_c(&lu->id); 567dcd589f8SShri Abhyankar #endif 568dcd589f8SShri Abhyankar PetscFunctionReturn(0); 569dcd589f8SShri Abhyankar } 570dcd589f8SShri Abhyankar 571397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 572397b6df1SKris Buschelman #undef __FUNCT__ 573f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 5740481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 575b24902e0SBarry Smith { 576719d5645SBarry Smith Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 577dcd589f8SShri Abhyankar PetscErrorCode ierr; 578dcd589f8SShri Abhyankar PetscTruth isSeqAIJ,isMPIAIJ; 579dcd589f8SShri Abhyankar PetscInt nz=0,M=A->rmap->N,rnz,i,nnz; 580dcd589f8SShri Abhyankar const PetscInt *ai,*aj; 581dcd589f8SShri Abhyankar PetscTruth valOnly; 582397b6df1SKris Buschelman 583397b6df1SKris Buschelman PetscFunctionBegin; 584b24902e0SBarry Smith lu->sym = 0; 585b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 586dcd589f8SShri Abhyankar 587dcd589f8SShri Abhyankar ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 588dcd589f8SShri Abhyankar ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 589dcd589f8SShri Abhyankar ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 590dcd589f8SShri Abhyankar lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 591dcd589f8SShri Abhyankar 592dcd589f8SShri Abhyankar /* Initialize a MUMPS instance */ 593dcd589f8SShri Abhyankar ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr); 594dcd589f8SShri Abhyankar /* Set MUMPS options */ 595dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 596dcd589f8SShri Abhyankar 597dcd589f8SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 598dcd589f8SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 599dcd589f8SShri Abhyankar switch (lu->id.ICNTL(18)){ 600dcd589f8SShri Abhyankar case 0: /* centralized assembled matrix input (size=1) */ 601dcd589f8SShri Abhyankar if (!lu->myid) { 602dcd589f8SShri Abhyankar Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 603dcd589f8SShri Abhyankar nz = aa->nz; 604dcd589f8SShri Abhyankar ai = aa->i; aj = aa->j; lu->val = aa->a; 605dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 606dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 607dcd589f8SShri Abhyankar nz = 0; 608dcd589f8SShri Abhyankar for (i=0; i<M; i++){ 609dcd589f8SShri Abhyankar rnz = ai[i+1] - ai[i]; 610dcd589f8SShri Abhyankar while (rnz--) { /* Fortran row/col index! */ 611dcd589f8SShri Abhyankar lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 612dcd589f8SShri Abhyankar } 613dcd589f8SShri Abhyankar } 614dcd589f8SShri Abhyankar } 615dcd589f8SShri Abhyankar break; 616dcd589f8SShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 617dcd589f8SShri Abhyankar valOnly = PETSC_FALSE; 618dcd589f8SShri Abhyankar ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 619dcd589f8SShri Abhyankar break; 620e32f2f54SBarry Smith default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 621dcd589f8SShri Abhyankar } 622dcd589f8SShri Abhyankar 623719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 624dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 625b24902e0SBarry Smith PetscFunctionReturn(0); 626b24902e0SBarry Smith } 627b24902e0SBarry Smith 628450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */ 629450b117fSShri Abhyankar #undef __FUNCT__ 630450b117fSShri Abhyankar #define __FUNCT__ "MatLUFactorSymbolic_BAIJMUMPS" 631450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 632450b117fSShri Abhyankar { 633dcd589f8SShri Abhyankar 634450b117fSShri Abhyankar Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 635dcd589f8SShri Abhyankar PetscErrorCode ierr; 636450b117fSShri Abhyankar 637450b117fSShri Abhyankar PetscFunctionBegin; 638450b117fSShri Abhyankar lu->sym = 0; 639450b117fSShri Abhyankar lu->matstruc = DIFFERENT_NONZERO_PATTERN; 640dcd589f8SShri Abhyankar ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 641dcd589f8SShri Abhyankar ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 642dcd589f8SShri Abhyankar ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 643dcd589f8SShri Abhyankar lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 644dcd589f8SShri Abhyankar 645dcd589f8SShri Abhyankar /* Initialize a MUMPS instance */ 646dcd589f8SShri Abhyankar ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr); 647dcd589f8SShri Abhyankar /* Set MUMPS options */ 648dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 649dcd589f8SShri Abhyankar 650450b117fSShri Abhyankar F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 651dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 652450b117fSShri Abhyankar PetscFunctionReturn(0); 653450b117fSShri Abhyankar } 654b24902e0SBarry Smith 655397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 656397b6df1SKris Buschelman #undef __FUNCT__ 657f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 6580481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 659b24902e0SBarry Smith { 6602792810eSHong Zhang Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 661dcd589f8SShri Abhyankar Mat newMat; 662dcd589f8SShri Abhyankar PetscErrorCode ierr; 663dcd589f8SShri Abhyankar PetscInt nz=0,M=A->rmap->N,rnz,i,nnz,jidx; 664dcd589f8SShri Abhyankar const PetscInt *ai,*aj,*adiag; 665dcd589f8SShri Abhyankar PetscScalar *av; 666dcd589f8SShri Abhyankar PetscTruth valOnly,isSeqAIJ,isMPIAIJ; 667dcd589f8SShri Abhyankar 668397b6df1SKris Buschelman 669397b6df1SKris Buschelman PetscFunctionBegin; 670b24902e0SBarry Smith lu->sym = 2; 671b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 672dcd589f8SShri Abhyankar ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 673dcd589f8SShri Abhyankar ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 674dcd589f8SShri Abhyankar ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 675dcd589f8SShri Abhyankar lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 676dcd589f8SShri Abhyankar 677dcd589f8SShri Abhyankar /* Initialize a MUMPS instance */ 678dcd589f8SShri Abhyankar ierr = PetscInitializeMUMPS(F);CHKERRQ(ierr); 679dcd589f8SShri Abhyankar /* Set MUMPS options */ 680dcd589f8SShri Abhyankar ierr = PetscSetMUMPSOptions(F,A);CHKERRQ(ierr); 681dcd589f8SShri Abhyankar 682dcd589f8SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 683dcd589f8SShri Abhyankar ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 684dcd589f8SShri Abhyankar switch (lu->id.ICNTL(18)){ 685dcd589f8SShri Abhyankar case 0: /* centralized assembled matrix input (size=1) */ 686dcd589f8SShri Abhyankar if (!lu->myid) { 687dcd589f8SShri Abhyankar if(isSeqAIJ) { 688dcd589f8SShri Abhyankar Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 689dcd589f8SShri Abhyankar nz = aa->nz; 690dcd589f8SShri Abhyankar ai = aa->i; aj = aa->j; adiag = aa->diag; av = aa->a; 691dcd589f8SShri Abhyankar } else { 692dcd589f8SShri Abhyankar Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 693dcd589f8SShri Abhyankar nz = aa->nz; 694dcd589f8SShri Abhyankar ai = aa->i; aj = aa->j; adiag = aa->diag; av = aa->a; 695dcd589f8SShri Abhyankar } 696dcd589f8SShri Abhyankar if ((F)->factortype == MAT_FACTOR_CHOLESKY && isSeqAIJ) { 697dcd589f8SShri Abhyankar nz = M + (nz - M)/2; /* nz for upper triangle */ 698dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 699dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 700dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscScalar),&lu->val);CHKERRQ(ierr); 701dcd589f8SShri Abhyankar nz = 0; 702dcd589f8SShri Abhyankar for (i=0; i<M; i++){ 703dcd589f8SShri Abhyankar rnz = ai[i+1] - adiag[i]; 704dcd589f8SShri Abhyankar jidx = adiag[i]; 705dcd589f8SShri Abhyankar while (rnz--) { /* Fortran row/col index! */ 706dcd589f8SShri Abhyankar lu->irn[nz] = i+1; lu->jcn[nz] = aj[jidx]+1; 707dcd589f8SShri Abhyankar lu->val[nz] = av[jidx]; jidx++; nz++; 708dcd589f8SShri Abhyankar } 709dcd589f8SShri Abhyankar } 710dcd589f8SShri Abhyankar } else { 711dcd589f8SShri Abhyankar lu->val = av; 712dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 713dcd589f8SShri Abhyankar ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 714dcd589f8SShri Abhyankar nz = 0; 715dcd589f8SShri Abhyankar for (i=0; i<M; i++){ 716dcd589f8SShri Abhyankar rnz = ai[i+1] - ai[i]; 717dcd589f8SShri Abhyankar while (rnz--) { /* Fortran row/col index! */ 718dcd589f8SShri Abhyankar lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 719dcd589f8SShri Abhyankar } 720dcd589f8SShri Abhyankar } 721dcd589f8SShri Abhyankar } 722dcd589f8SShri Abhyankar } 723dcd589f8SShri Abhyankar break; 724dcd589f8SShri Abhyankar case 3: /* distributed assembled matrix input (size>1) */ 725dcd589f8SShri Abhyankar valOnly = PETSC_FALSE; 726dcd589f8SShri Abhyankar if(((F)->factortype == MAT_FACTOR_CHOLESKY) && isMPIAIJ) { 727dcd589f8SShri Abhyankar /* Create an SBAIJ matrix and use this matrix to set the lu values */ 728dcd589f8SShri Abhyankar ierr = MatConvert(A,MATMPISBAIJ,MAT_INITIAL_MATRIX,&newMat);CHKERRQ(ierr); 729dcd589f8SShri Abhyankar ierr = MatConvertToTriples(newMat,1,valOnly,&nnz,&lu->irn , &lu->jcn, &lu->val);CHKERRQ(ierr); 730dcd589f8SShri Abhyankar ierr = MatDestroy(newMat);CHKERRQ(ierr); 731dcd589f8SShri Abhyankar } else { 732dcd589f8SShri Abhyankar ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 733dcd589f8SShri Abhyankar } 734dcd589f8SShri Abhyankar break; 735e32f2f54SBarry Smith default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 736dcd589f8SShri Abhyankar } 737dcd589f8SShri Abhyankar 7382792810eSHong Zhang F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 739dcd589f8SShri Abhyankar F->ops->solve = MatSolve_MUMPS; 740db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX) 741dcd589f8SShri Abhyankar (F)->ops->getinertia = MatGetInertia_SBAIJMUMPS; 742db4efbfdSBarry Smith #endif 743b24902e0SBarry Smith PetscFunctionReturn(0); 744b24902e0SBarry Smith } 745b24902e0SBarry Smith 746397b6df1SKris Buschelman #undef __FUNCT__ 747f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS" 74874ed9c26SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) 74974ed9c26SBarry Smith { 750f6c57405SHong Zhang Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 751f6c57405SHong Zhang PetscErrorCode ierr; 752f6c57405SHong Zhang 753f6c57405SHong Zhang PetscFunctionBegin; 754f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 755f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 756f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 757f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 758f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 759f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 760f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 761f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 762f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 763f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 764f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 765f6c57405SHong Zhang if (!lu->myid && lu->id.ICNTL(11)>0) { 766f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 767f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 768f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 769f6c57405SHong 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); 770f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 771f6c57405SHong 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); 772f6c57405SHong Zhang 773f6c57405SHong Zhang } 774f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 775f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 776f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 777f6c57405SHong Zhang /* ICNTL(15-17) not used */ 778f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 779f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 780f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 781f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 782c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(22) (in-core/out-of-core facility): %d \n",lu->id.ICNTL(22));CHKERRQ(ierr); 783c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr); 784c0165424SHong Zhang 785c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(24) (detection of null pivot rows): %d \n",lu->id.ICNTL(24));CHKERRQ(ierr); 786c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(25) (computation of a null space basis): %d \n",lu->id.ICNTL(25));CHKERRQ(ierr); 787c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(26) (Schur options for rhs or solution): %d \n",lu->id.ICNTL(26));CHKERRQ(ierr); 788c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(27) (experimental parameter): %d \n",lu->id.ICNTL(27));CHKERRQ(ierr); 789f6c57405SHong Zhang 790f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 791f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 792f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 793f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 794c0165424SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(5) (fixation for null pivots): %g \n",lu->id.CNTL(5));CHKERRQ(ierr); 795f6c57405SHong Zhang 796f6c57405SHong Zhang /* infomation local to each processor */ 797f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);} 7987adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 7997adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 800f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);} 8017adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 8027adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 803f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);} 8047adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 8057adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 806f6c57405SHong Zhang 807f6c57405SHong 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);} 8087adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 8097adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 810f6c57405SHong Zhang 811f6c57405SHong 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);} 8127adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 8137adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 814f6c57405SHong Zhang 815f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);} 8167adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 8177adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 818f6c57405SHong Zhang 819f6c57405SHong Zhang if (!lu->myid){ /* information from the host */ 820f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 821f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 822f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 823f6c57405SHong Zhang 824f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 825f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 826f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 827f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 828f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 829f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 830f6c57405SHong 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); 831f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 832f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 833f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 834f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 835f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 836f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 837f6c57405SHong 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); 838f6c57405SHong 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); 839f6c57405SHong 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); 840f6c57405SHong 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); 841f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 842f6c57405SHong 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); 843f6c57405SHong 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); 844f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 845f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 846f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 847f6c57405SHong Zhang } 848f6c57405SHong Zhang PetscFunctionReturn(0); 849f6c57405SHong Zhang } 850f6c57405SHong Zhang 851f6c57405SHong Zhang #undef __FUNCT__ 852f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS" 853b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 854b24902e0SBarry Smith { 855f6c57405SHong Zhang PetscErrorCode ierr; 856f6c57405SHong Zhang PetscTruth iascii; 857f6c57405SHong Zhang PetscViewerFormat format; 858cb828f0fSHong Zhang Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr; 859f6c57405SHong Zhang 860f6c57405SHong Zhang PetscFunctionBegin; 861cb828f0fSHong Zhang /* check if matrix is mumps type */ 862cb828f0fSHong Zhang if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 863cb828f0fSHong Zhang 864f6c57405SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 865f6c57405SHong Zhang if (iascii) { 866f6c57405SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 867cb828f0fSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO || mumps->mumpsview){ 868cb828f0fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 869cb828f0fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",mumps->id.sym);CHKERRQ(ierr); 870cb828f0fSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",mumps->id.par);CHKERRQ(ierr); 871cb828f0fSHong Zhang if (mumps->mumpsview){ /* View all MUMPS parameters */ 872f6c57405SHong Zhang ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 873f6c57405SHong Zhang } 874f6c57405SHong Zhang } 875cb828f0fSHong Zhang } 876f6c57405SHong Zhang PetscFunctionReturn(0); 877f6c57405SHong Zhang } 878f6c57405SHong Zhang 87935bd34faSBarry Smith #undef __FUNCT__ 88035bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS" 88135bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info) 88235bd34faSBarry Smith { 883cb828f0fSHong Zhang Mat_MUMPS *mumps =(Mat_MUMPS*)A->spptr; 88435bd34faSBarry Smith 88535bd34faSBarry Smith PetscFunctionBegin; 88635bd34faSBarry Smith info->block_size = 1.0; 887cb828f0fSHong Zhang info->nz_allocated = mumps->id.INFOG(20); 888cb828f0fSHong Zhang info->nz_used = mumps->id.INFOG(20); 88935bd34faSBarry Smith info->nz_unneeded = 0.0; 89035bd34faSBarry Smith info->assemblies = 0.0; 89135bd34faSBarry Smith info->mallocs = 0.0; 89235bd34faSBarry Smith info->memory = 0.0; 89335bd34faSBarry Smith info->fill_ratio_given = 0; 89435bd34faSBarry Smith info->fill_ratio_needed = 0; 89535bd34faSBarry Smith info->factor_mallocs = 0; 89635bd34faSBarry Smith PetscFunctionReturn(0); 89735bd34faSBarry Smith } 89835bd34faSBarry Smith 89924b6179bSKris Buschelman /*MC 90041c8de11SBarry Smith MAT_SOLVER_MUMPS - A matrix type providing direct solvers (LU and Cholesky) for 90124b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 90224b6179bSKris Buschelman 90341c8de11SBarry Smith Works with MATAIJ and MATSBAIJ matrices 90424b6179bSKris Buschelman 90524b6179bSKris Buschelman Options Database Keys: 90641c8de11SBarry Smith + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 90724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 90824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 90924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 91024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 91124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 91294b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 91324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 91424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 91524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 91624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 91724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 91824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 91924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 92024b6179bSKris Buschelman 92124b6179bSKris Buschelman Level: beginner 92224b6179bSKris Buschelman 92341c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 92441c8de11SBarry Smith 92524b6179bSKris Buschelman M*/ 92624b6179bSKris Buschelman 9272877fffaSHong Zhang EXTERN_C_BEGIN 92835bd34faSBarry Smith #undef __FUNCT__ 92935bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps" 93035bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type) 93135bd34faSBarry Smith { 93235bd34faSBarry Smith PetscFunctionBegin; 93335bd34faSBarry Smith *type = MAT_SOLVER_MUMPS; 93435bd34faSBarry Smith PetscFunctionReturn(0); 93535bd34faSBarry Smith } 93635bd34faSBarry Smith EXTERN_C_END 93735bd34faSBarry Smith 93835bd34faSBarry Smith EXTERN_C_BEGIN 9392877fffaSHong Zhang /* 9402877fffaSHong Zhang The seq and mpi versions of this function are the same 9412877fffaSHong Zhang */ 9422877fffaSHong Zhang #undef __FUNCT__ 9432877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_mumps" 9442877fffaSHong Zhang PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F) 9452877fffaSHong Zhang { 9462877fffaSHong Zhang Mat B; 9472877fffaSHong Zhang PetscErrorCode ierr; 9482877fffaSHong Zhang Mat_MUMPS *mumps; 9492877fffaSHong Zhang 9502877fffaSHong Zhang PetscFunctionBegin; 9512877fffaSHong Zhang /* Create the factorization matrix */ 9522877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 9532877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 9542877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9552877fffaSHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 9562877fffaSHong Zhang 9572877fffaSHong Zhang B->ops->view = MatView_MUMPS; 95835bd34faSBarry Smith B->ops->getinfo = MatGetInfo_MUMPS; 95935bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 96097969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 961450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 962450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 963d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 964dcd589f8SShri Abhyankar } else { 965450b117fSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 966450b117fSShri Abhyankar B->factortype = MAT_FACTOR_CHOLESKY; 967450b117fSShri Abhyankar } 9682877fffaSHong Zhang 9692877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 9702877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 9712877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 9722877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 9732877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 9742877fffaSHong Zhang mumps->nSolve = 0; 9752877fffaSHong Zhang mumps->MatDestroy = B->ops->destroy; 9762877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 9772877fffaSHong Zhang B->spptr = (void*)mumps; 9782877fffaSHong Zhang 9792877fffaSHong Zhang *F = B; 9802877fffaSHong Zhang PetscFunctionReturn(0); 9812877fffaSHong Zhang } 9822877fffaSHong Zhang EXTERN_C_END 9832877fffaSHong Zhang 9842877fffaSHong Zhang EXTERN_C_BEGIN 9852877fffaSHong Zhang #undef __FUNCT__ 9862877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_seqsbaij_mumps" 9872877fffaSHong Zhang PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 9882877fffaSHong Zhang { 9892877fffaSHong Zhang Mat B; 9902877fffaSHong Zhang PetscErrorCode ierr; 9912877fffaSHong Zhang Mat_MUMPS *mumps; 9922877fffaSHong Zhang 9932877fffaSHong Zhang PetscFunctionBegin; 994*e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 9952877fffaSHong Zhang /* Create the factorization matrix */ 9962877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 9972877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 9982877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9992877fffaSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 10002877fffaSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 10012877fffaSHong Zhang 10022877fffaSHong Zhang B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 10032877fffaSHong Zhang B->ops->view = MatView_MUMPS; 100435bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 100597969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1006d5f3da31SBarry Smith B->factortype = MAT_FACTOR_CHOLESKY; 10072877fffaSHong Zhang 10082877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 10092877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 1010450b117fSShri Abhyankar mumps->isAIJ = PETSC_FALSE; 10112877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 10122877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 10132877fffaSHong Zhang mumps->nSolve = 0; 10142877fffaSHong Zhang mumps->MatDestroy = B->ops->destroy; 10152877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 10162877fffaSHong Zhang B->spptr = (void*)mumps; 1017f3c0ef26SHong Zhang 10182877fffaSHong Zhang *F = B; 10192877fffaSHong Zhang PetscFunctionReturn(0); 10202877fffaSHong Zhang } 10212877fffaSHong Zhang EXTERN_C_END 10222877fffaSHong Zhang 10232877fffaSHong Zhang EXTERN_C_BEGIN 10242877fffaSHong Zhang #undef __FUNCT__ 10252877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_mpisbaij_mumps" 10262877fffaSHong Zhang PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 10272877fffaSHong Zhang { 10282877fffaSHong Zhang Mat B; 10292877fffaSHong Zhang PetscErrorCode ierr; 10302877fffaSHong Zhang Mat_MUMPS *mumps; 10312877fffaSHong Zhang 10322877fffaSHong Zhang PetscFunctionBegin; 1033*e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 10342877fffaSHong Zhang /* Create the factorization matrix */ 10352877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 10362877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 10372877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 10382877fffaSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 10392877fffaSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 10402877fffaSHong Zhang 10412877fffaSHong Zhang B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 10422877fffaSHong Zhang B->ops->view = MatView_MUMPS; 104335bd34faSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 104497969023SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1045d5f3da31SBarry Smith B->factortype = MAT_FACTOR_CHOLESKY; 10462877fffaSHong Zhang 10472877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 10482877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 1049a214ac2aSShri Abhyankar mumps->isAIJ = PETSC_FALSE; 1050a214ac2aSShri Abhyankar mumps->scat_rhs = PETSC_NULL; 1051a214ac2aSShri Abhyankar mumps->scat_sol = PETSC_NULL; 1052a214ac2aSShri Abhyankar mumps->nSolve = 0; 1053a214ac2aSShri Abhyankar mumps->MatDestroy = B->ops->destroy; 1054a214ac2aSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 1055a214ac2aSShri Abhyankar B->spptr = (void*)mumps; 1056a214ac2aSShri Abhyankar 1057a214ac2aSShri Abhyankar *F = B; 1058a214ac2aSShri Abhyankar PetscFunctionReturn(0); 1059a214ac2aSShri Abhyankar } 1060a214ac2aSShri Abhyankar EXTERN_C_END 1061a214ac2aSShri Abhyankar 1062a214ac2aSShri Abhyankar EXTERN_C_BEGIN 1063a214ac2aSShri Abhyankar #undef __FUNCT__ 1064a214ac2aSShri Abhyankar #define __FUNCT__ "MatGetFactor_mpiaij_mumps" 1065a214ac2aSShri Abhyankar PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F) 1066a214ac2aSShri Abhyankar { 1067a214ac2aSShri Abhyankar Mat B; 1068a214ac2aSShri Abhyankar PetscErrorCode ierr; 1069a214ac2aSShri Abhyankar Mat_MUMPS *mumps; 1070a214ac2aSShri Abhyankar 1071a214ac2aSShri Abhyankar PetscFunctionBegin; 1072a214ac2aSShri Abhyankar /* Create the factorization matrix */ 1073a214ac2aSShri Abhyankar ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 1074a214ac2aSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1075a214ac2aSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1076a214ac2aSShri Abhyankar ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1077a214ac2aSShri Abhyankar ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1078a214ac2aSShri Abhyankar 1079a214ac2aSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1080a214ac2aSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 1081f4762488SHong Zhang B->factortype = MAT_FACTOR_LU; 1082dcd589f8SShri Abhyankar } else { 1083a214ac2aSShri Abhyankar B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 1084f4762488SHong Zhang B->factortype = MAT_FACTOR_CHOLESKY; 1085450b117fSShri Abhyankar } 1086a214ac2aSShri Abhyankar 1087a214ac2aSShri Abhyankar B->ops->view = MatView_MUMPS; 1088a214ac2aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1089a214ac2aSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1090a214ac2aSShri Abhyankar 1091a214ac2aSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1092a214ac2aSShri Abhyankar mumps->CleanUpMUMPS = PETSC_FALSE; 10932877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 10942877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 10952877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 10962877fffaSHong Zhang mumps->nSolve = 0; 1097f3c0ef26SHong Zhang mumps->MatDestroy = B->ops->destroy; 1098f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 10992877fffaSHong Zhang B->spptr = (void*)mumps; 1100f3c0ef26SHong Zhang 11012877fffaSHong Zhang *F = B; 11022877fffaSHong Zhang PetscFunctionReturn(0); 11032877fffaSHong Zhang } 11042877fffaSHong Zhang EXTERN_C_END 110597969023SHong Zhang 1106450b117fSShri Abhyankar EXTERN_C_BEGIN 1107450b117fSShri Abhyankar #undef __FUNCT__ 1108450b117fSShri Abhyankar #define __FUNCT__ "MatGetFactor_mpibaij_mumps" 1109450b117fSShri Abhyankar PetscErrorCode MatGetFactor_mpibaij_mumps(Mat A,MatFactorType ftype,Mat *F) 1110450b117fSShri Abhyankar { 1111450b117fSShri Abhyankar Mat B; 1112450b117fSShri Abhyankar PetscErrorCode ierr; 1113450b117fSShri Abhyankar Mat_MUMPS *mumps; 1114450b117fSShri Abhyankar 1115450b117fSShri Abhyankar PetscFunctionBegin; 1116450b117fSShri Abhyankar /* Create the factorization matrix */ 1117450b117fSShri Abhyankar ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 1118450b117fSShri Abhyankar ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1119450b117fSShri Abhyankar ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1120450b117fSShri Abhyankar ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 1121450b117fSShri Abhyankar ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 1122450b117fSShri Abhyankar 1123450b117fSShri Abhyankar if (ftype == MAT_FACTOR_LU) { 1124450b117fSShri Abhyankar B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS; 1125450b117fSShri Abhyankar B->factortype = MAT_FACTOR_LU; 1126450b117fSShri Abhyankar } 1127*e7e72b3dSBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cholesky factorization for Petsc BAIJ matrices not supported yet\n"); 1128450b117fSShri Abhyankar B->ops->view = MatView_MUMPS; 1129450b117fSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr); 1130450b117fSShri Abhyankar ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr); 1131450b117fSShri Abhyankar 1132450b117fSShri Abhyankar ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 1133450b117fSShri Abhyankar mumps->CleanUpMUMPS = PETSC_FALSE; 1134450b117fSShri Abhyankar mumps->isAIJ = PETSC_TRUE; 1135450b117fSShri Abhyankar mumps->scat_rhs = PETSC_NULL; 1136450b117fSShri Abhyankar mumps->scat_sol = PETSC_NULL; 1137450b117fSShri Abhyankar mumps->nSolve = 0; 1138450b117fSShri Abhyankar mumps->MatDestroy = B->ops->destroy; 1139450b117fSShri Abhyankar B->ops->destroy = MatDestroy_MUMPS; 1140450b117fSShri Abhyankar B->spptr = (void*)mumps; 1141450b117fSShri Abhyankar 1142450b117fSShri Abhyankar *F = B; 1143450b117fSShri Abhyankar PetscFunctionReturn(0); 1144450b117fSShri Abhyankar } 1145450b117fSShri Abhyankar EXTERN_C_END 1146a214ac2aSShri Abhyankar 114761288eaaSHong Zhang /* -------------------------------------------------------------------------------------------*/ 114861288eaaSHong Zhang /*@ 114961288eaaSHong Zhang MatMumpsSetIcntl - Set MUMPS parameter ICNTL() 115061288eaaSHong Zhang 115161288eaaSHong Zhang Collective on Mat 115261288eaaSHong Zhang 115361288eaaSHong Zhang Input Parameters: 115461288eaaSHong Zhang + F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface 115561288eaaSHong Zhang . idx - index of MUMPS parameter array ICNTL() 115661288eaaSHong Zhang - icntl - value of MUMPS ICNTL(imumps) 115761288eaaSHong Zhang 115861288eaaSHong Zhang Options Database: 115961288eaaSHong Zhang . -mat_mumps_icntl_<idx> <icntl> 116061288eaaSHong Zhang 116161288eaaSHong Zhang Level: beginner 116261288eaaSHong Zhang 116361288eaaSHong Zhang References: MUMPS Users' Guide 116461288eaaSHong Zhang 116561288eaaSHong Zhang .seealso: MatGetFactor() 116661288eaaSHong Zhang @*/ 116797969023SHong Zhang #undef __FUNCT__ 116897969023SHong Zhang #define __FUNCT__ "MatMumpsSetIcntl" 116986e5884dSMatthew Knepley PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt idx,PetscInt icntl) 117097969023SHong Zhang { 1171dcd589f8SShri Abhyankar Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 117297969023SHong Zhang 117397969023SHong Zhang PetscFunctionBegin; 117461288eaaSHong Zhang lu->id.ICNTL(idx) = icntl; 117597969023SHong Zhang PetscFunctionReturn(0); 117697969023SHong Zhang } 117797969023SHong Zhang 1178