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; 4267334b25SHong Zhang PetscErrorCode (*MatDestroy)(Mat); 43f0c56d0fSKris Buschelman } Mat_MUMPS; 44f0c56d0fSKris Buschelman 45dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*); 46b24902e0SBarry Smith 47397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */ 48397b6df1SKris Buschelman /* 49397b6df1SKris Buschelman input: 5075747be1SHong Zhang A - matrix in mpiaij or mpisbaij (bs=1) format 51397b6df1SKris Buschelman shift - 0: C style output triple; 1: Fortran style output triple. 52397b6df1SKris Buschelman valOnly - FALSE: spaces are allocated and values are set for the triple 53397b6df1SKris Buschelman TRUE: only the values in v array are updated 54397b6df1SKris Buschelman output: 55397b6df1SKris Buschelman nnz - dim of r, c, and v (number of local nonzero entries of A) 56397b6df1SKris Buschelman r, c, v - row and col index, matrix values (matrix triples) 57397b6df1SKris Buschelman */ 58b24902e0SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) 59b24902e0SBarry Smith { 60c1490034SHong Zhang PetscInt *ai, *aj, *bi, *bj, rstart,nz, *garray; 61dfbe8321SBarry Smith PetscErrorCode ierr; 62d0f46423SBarry Smith PetscInt i,j,jj,jB,irow,m=A->rmap->n,*ajj,*bjj,countA,countB,colA_start,jcol; 63c1490034SHong Zhang PetscInt *row,*col; 64397b6df1SKris Buschelman PetscScalar *av, *bv,*val; 655c9eb25fSBarry Smith PetscTruth isAIJ; 66397b6df1SKris Buschelman 67397b6df1SKris Buschelman PetscFunctionBegin; 685c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr); 695c9eb25fSBarry Smith if (isAIJ){ 70397b6df1SKris Buschelman Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data; 71397b6df1SKris Buschelman Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data; 72397b6df1SKris Buschelman Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data; 73397b6df1SKris Buschelman nz = aa->nz + bb->nz; 74d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 75397b6df1SKris Buschelman garray = mat->garray; 76397b6df1SKris Buschelman av=aa->a; bv=bb->a; 77397b6df1SKris Buschelman 78397b6df1SKris Buschelman } else { 79397b6df1SKris Buschelman Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data; 80397b6df1SKris Buschelman Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data; 81397b6df1SKris Buschelman Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data; 82d0f46423SBarry Smith if (A->rmap->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap->bs); 836c6c5352SBarry Smith nz = aa->nz + bb->nz; 84d0f46423SBarry Smith ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart; 85397b6df1SKris Buschelman garray = mat->garray; 86397b6df1SKris Buschelman av=aa->a; bv=bb->a; 87397b6df1SKris Buschelman } 88397b6df1SKris Buschelman 89397b6df1SKris Buschelman if (!valOnly){ 907c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr); 917c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr); 92397b6df1SKris Buschelman ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr); 93397b6df1SKris Buschelman *r = row; *c = col; *v = val; 94397b6df1SKris Buschelman } else { 95397b6df1SKris Buschelman row = *r; col = *c; val = *v; 96397b6df1SKris Buschelman } 97397b6df1SKris Buschelman *nnz = nz; 98397b6df1SKris Buschelman 99028e57e8SHong Zhang jj = 0; irow = rstart; 100397b6df1SKris Buschelman for ( i=0; i<m; i++ ) { 101397b6df1SKris Buschelman ajj = aj + ai[i]; /* ptr to the beginning of this row */ 102397b6df1SKris Buschelman countA = ai[i+1] - ai[i]; 103397b6df1SKris Buschelman countB = bi[i+1] - bi[i]; 104397b6df1SKris Buschelman bjj = bj + bi[i]; 105397b6df1SKris Buschelman 106397b6df1SKris Buschelman /* get jB, the starting local col index for the 2nd B-part */ 107397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 10875747be1SHong Zhang j=-1; 10975747be1SHong Zhang do { 11075747be1SHong Zhang j++; 11175747be1SHong Zhang if (j == countB) break; 112397b6df1SKris Buschelman jcol = garray[bjj[j]]; 11375747be1SHong Zhang } while (jcol < colA_start); 11475747be1SHong Zhang jB = j; 115397b6df1SKris Buschelman 116397b6df1SKris Buschelman /* B-part, smaller col index */ 117397b6df1SKris Buschelman colA_start = rstart + ajj[0]; /* the smallest col index for A */ 118397b6df1SKris Buschelman for (j=0; j<jB; j++){ 119397b6df1SKris Buschelman jcol = garray[bjj[j]]; 120397b6df1SKris Buschelman if (!valOnly){ 121397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = jcol + shift; 12275747be1SHong Zhang 123397b6df1SKris Buschelman } 124397b6df1SKris Buschelman val[jj++] = *bv++; 125397b6df1SKris Buschelman } 126397b6df1SKris Buschelman /* A-part */ 127397b6df1SKris Buschelman for (j=0; j<countA; j++){ 128397b6df1SKris Buschelman if (!valOnly){ 129397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift; 130397b6df1SKris Buschelman } 131397b6df1SKris Buschelman val[jj++] = *av++; 132397b6df1SKris Buschelman } 133397b6df1SKris Buschelman /* B-part, larger col index */ 134397b6df1SKris Buschelman for (j=jB; j<countB; j++){ 135397b6df1SKris Buschelman if (!valOnly){ 136397b6df1SKris Buschelman row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift; 137397b6df1SKris Buschelman } 138397b6df1SKris Buschelman val[jj++] = *bv++; 139397b6df1SKris Buschelman } 140397b6df1SKris Buschelman irow++; 141397b6df1SKris Buschelman } 142397b6df1SKris Buschelman 143397b6df1SKris Buschelman PetscFunctionReturn(0); 144397b6df1SKris Buschelman } 145397b6df1SKris Buschelman 146397b6df1SKris Buschelman #undef __FUNCT__ 1473924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS" 148dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A) 149dfbe8321SBarry Smith { 150f0c56d0fSKris Buschelman Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 151dfbe8321SBarry Smith PetscErrorCode ierr; 152c1490034SHong Zhang PetscMPIInt size=lu->size; 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 } 17567334b25SHong Zhang ierr = (lu->MatDestroy)(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 265ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX) 266a58c3f20SHong Zhang /* 267a58c3f20SHong Zhang input: 268a58c3f20SHong Zhang F: numeric factor 269a58c3f20SHong Zhang output: 270a58c3f20SHong Zhang nneg: total number of negative pivots 271a58c3f20SHong Zhang nzero: 0 272a58c3f20SHong Zhang npos: (global dimension of F) - nneg 273a58c3f20SHong Zhang */ 274a58c3f20SHong Zhang 275a58c3f20SHong Zhang #undef __FUNCT__ 276a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS" 277dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos) 278a58c3f20SHong Zhang { 279a58c3f20SHong Zhang Mat_MUMPS *lu =(Mat_MUMPS*)F->spptr; 280dfbe8321SBarry Smith PetscErrorCode ierr; 281c1490034SHong Zhang PetscMPIInt size; 282a58c3f20SHong Zhang 283a58c3f20SHong Zhang PetscFunctionBegin; 2847adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); 285bcb30aebSHong 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 */ 286bcb30aebSHong Zhang if (size > 1 && lu->id.ICNTL(13) != 1){ 28779a5c55eSBarry 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)); 288bcb30aebSHong Zhang } 289a58c3f20SHong Zhang if (nneg){ 290a58c3f20SHong Zhang if (!lu->myid){ 291a58c3f20SHong Zhang *nneg = lu->id.INFOG(12); 292a58c3f20SHong Zhang } 293bcb30aebSHong Zhang ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr); 294a58c3f20SHong Zhang } 295a58c3f20SHong Zhang if (nzero) *nzero = 0; 296d0f46423SBarry Smith if (npos) *npos = F->rmap->N - (*nneg); 297a58c3f20SHong Zhang PetscFunctionReturn(0); 298a58c3f20SHong Zhang } 299ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */ 300a58c3f20SHong Zhang 301397b6df1SKris Buschelman #undef __FUNCT__ 302f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS" 3030481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info) 304af281ebdSHong Zhang { 305719d5645SBarry Smith Mat_MUMPS *lu =(Mat_MUMPS*)(F)->spptr; 3066849ba73SBarry Smith PetscErrorCode ierr; 307d0f46423SBarry Smith PetscInt rnz,nnz,nz=0,i,M=A->rmap->N,*ai,*aj,icntl; 308397b6df1SKris Buschelman PetscTruth valOnly,flg; 309e09efc27SHong Zhang Mat F_diag; 310c349612cSHong Zhang IS is_iden; 311c349612cSHong Zhang Vec b; 3125c9eb25fSBarry Smith PetscTruth isSeqAIJ,isSeqSBAIJ; 313397b6df1SKris Buschelman 314397b6df1SKris Buschelman PetscFunctionBegin; 3155c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 3165c9eb25fSBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 317397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 318719d5645SBarry Smith (F)->ops->solve = MatSolve_MUMPS; 319397b6df1SKris Buschelman 320397b6df1SKris Buschelman /* Initialize a MUMPS instance */ 3217adad957SLisandro Dalcin ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid); 3227adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr); 323397b6df1SKris Buschelman lu->id.job = JOB_INIT; 3247adad957SLisandro Dalcin ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr); 3256a1dac61SBarry Smith lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps); 326397b6df1SKris Buschelman 327397b6df1SKris Buschelman /* Set mumps options */ 3287adad957SLisandro Dalcin ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr); 329397b6df1SKris Buschelman lu->id.par=1; /* host participates factorizaton and solve */ 330397b6df1SKris Buschelman lu->id.sym=lu->sym; 331397b6df1SKris Buschelman if (lu->sym == 2){ 332397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr); 333397b6df1SKris Buschelman if (flg && icntl == 1) lu->id.sym=icntl; /* matrix is spd */ 334397b6df1SKris Buschelman } 335397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 336397b6df1SKris Buschelman zmumps_c(&lu->id); 337397b6df1SKris Buschelman #else 338397b6df1SKris Buschelman dmumps_c(&lu->id); 339397b6df1SKris Buschelman #endif 340397b6df1SKris Buschelman 3415c9eb25fSBarry Smith if (isSeqAIJ || isSeqSBAIJ){ 342397b6df1SKris Buschelman lu->id.ICNTL(18) = 0; /* centralized assembled matrix input */ 343397b6df1SKris Buschelman } else { 344397b6df1SKris Buschelman lu->id.ICNTL(18) = 3; /* distributed assembled matrix input */ 345397b6df1SKris Buschelman } 346397b6df1SKris Buschelman 347397b6df1SKris Buschelman icntl=-1; 34821f4b680SHong Zhang lu->id.ICNTL(4) = 0; /* level of printing; overwrite mumps default ICNTL(4)=2 */ 349397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr); 35019facb7aSBarry Smith if ((flg && icntl > 0) || PetscLogPrintInfo) { 351397b6df1SKris Buschelman lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */ 352397b6df1SKris Buschelman } else { /* no output */ 353397b6df1SKris Buschelman lu->id.ICNTL(1) = 0; /* error message, default= 6 */ 354397b6df1SKris Buschelman lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */ 355397b6df1SKris Buschelman lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */ 356397b6df1SKris Buschelman } 357397b6df1SKris 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); 358397b6df1SKris Buschelman icntl=-1; 359397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr); 360397b6df1SKris Buschelman if (flg) { 361397b6df1SKris Buschelman if (icntl== 1){ 362397b6df1SKris Buschelman SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n"); 363397b6df1SKris Buschelman } else { 364397b6df1SKris Buschelman lu->id.ICNTL(7) = icntl; 365397b6df1SKris Buschelman } 366397b6df1SKris Buschelman } 367397b6df1SKris 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); 368397b6df1SKris 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); 36994b7f48cSBarry 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); 370397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr); 371397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr); 372adc1d99fSHong 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); 373397b6df1SKris Buschelman ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr); 374397b6df1SKris Buschelman 375397b6df1SKris 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); 376397b6df1SKris 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); 377397b6df1SKris 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); 37825f9c88cSHong 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); 379397b6df1SKris Buschelman PetscOptionsEnd(); 380397b6df1SKris Buschelman } 381397b6df1SKris Buschelman 382397b6df1SKris Buschelman /* define matrix A */ 383397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 384397b6df1SKris Buschelman case 0: /* centralized assembled matrix input (size=1) */ 385397b6df1SKris Buschelman if (!lu->myid) { 3865c9eb25fSBarry Smith if (isSeqAIJ){ 387397b6df1SKris Buschelman Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 388397b6df1SKris Buschelman nz = aa->nz; 389397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 3905c9eb25fSBarry Smith } else if (isSeqSBAIJ) { 391397b6df1SKris Buschelman Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 3926c6c5352SBarry Smith nz = aa->nz; 393397b6df1SKris Buschelman ai = aa->i; aj = aa->j; lu->val = aa->a; 3945c9eb25fSBarry Smith } else { 3955c9eb25fSBarry Smith SETERRQ(PETSC_ERR_SUP,"No mumps factorization for this matrix type"); 396397b6df1SKris Buschelman } 397397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */ 3987c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr); 3997c307921SBarry Smith ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr); 400397b6df1SKris Buschelman nz = 0; 401397b6df1SKris Buschelman for (i=0; i<M; i++){ 402397b6df1SKris Buschelman rnz = ai[i+1] - ai[i]; 403397b6df1SKris Buschelman while (rnz--) { /* Fortran row/col index! */ 404397b6df1SKris Buschelman lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++; 405397b6df1SKris Buschelman } 406397b6df1SKris Buschelman } 407397b6df1SKris Buschelman } 408397b6df1SKris Buschelman } 409397b6df1SKris Buschelman break; 410397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 411397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 412397b6df1SKris Buschelman valOnly = PETSC_FALSE; 413397b6df1SKris Buschelman } else { 414397b6df1SKris Buschelman valOnly = PETSC_TRUE; /* only update mat values, not row and col index */ 415397b6df1SKris Buschelman } 416397b6df1SKris Buschelman ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr); 417397b6df1SKris Buschelman break; 418397b6df1SKris Buschelman default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS."); 419397b6df1SKris Buschelman } 420397b6df1SKris Buschelman 421397b6df1SKris Buschelman /* analysis phase */ 422329ec9b3SHong Zhang /*----------------*/ 423397b6df1SKris Buschelman if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 424329ec9b3SHong Zhang lu->id.job = 1; 425329ec9b3SHong Zhang 426397b6df1SKris Buschelman lu->id.n = M; 427397b6df1SKris Buschelman switch (lu->id.ICNTL(18)){ 428397b6df1SKris Buschelman case 0: /* centralized assembled matrix input */ 429397b6df1SKris Buschelman if (!lu->myid) { 430397b6df1SKris Buschelman lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn; 431397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1){ 432397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 433397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 434397b6df1SKris Buschelman #else 435397b6df1SKris Buschelman lu->id.a = lu->val; 436397b6df1SKris Buschelman #endif 437397b6df1SKris Buschelman } 438397b6df1SKris Buschelman } 439397b6df1SKris Buschelman break; 440397b6df1SKris Buschelman case 3: /* distributed assembled matrix input (size>1) */ 441397b6df1SKris Buschelman lu->id.nz_loc = nnz; 442397b6df1SKris Buschelman lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn; 443397b6df1SKris Buschelman if (lu->id.ICNTL(6)>1) { 444397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 445397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 446397b6df1SKris Buschelman #else 447397b6df1SKris Buschelman lu->id.a_loc = lu->val; 448397b6df1SKris Buschelman #endif 449397b6df1SKris Buschelman } 450329ec9b3SHong Zhang /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 451329ec9b3SHong Zhang if (!lu->myid){ 452d0f46423SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 453d0f46423SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 454329ec9b3SHong Zhang } else { 455329ec9b3SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr); 456329ec9b3SHong Zhang ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr); 457329ec9b3SHong Zhang } 4587adad957SLisandro Dalcin ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 459d0f46423SBarry Smith ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 460329ec9b3SHong Zhang ierr = VecSetFromOptions(b);CHKERRQ(ierr); 461329ec9b3SHong Zhang 462329ec9b3SHong Zhang ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 463329ec9b3SHong Zhang ierr = ISDestroy(is_iden);CHKERRQ(ierr); 464329ec9b3SHong Zhang ierr = VecDestroy(b);CHKERRQ(ierr); 465397b6df1SKris Buschelman break; 466397b6df1SKris Buschelman } 467397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 468397b6df1SKris Buschelman zmumps_c(&lu->id); 469397b6df1SKris Buschelman #else 470397b6df1SKris Buschelman dmumps_c(&lu->id); 471397b6df1SKris Buschelman #endif 472397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 47379a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1)); 474397b6df1SKris Buschelman } 475397b6df1SKris Buschelman } 476397b6df1SKris Buschelman 477397b6df1SKris Buschelman /* numerical factorization phase */ 478329ec9b3SHong Zhang /*-------------------------------*/ 479329ec9b3SHong Zhang lu->id.job = 2; 480958c9bccSBarry Smith if(!lu->id.ICNTL(18)) { 481a7aca84bSHong Zhang if (!lu->myid) { 482397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 483397b6df1SKris Buschelman lu->id.a = (mumps_double_complex*)lu->val; 484397b6df1SKris Buschelman #else 485397b6df1SKris Buschelman lu->id.a = lu->val; 486397b6df1SKris Buschelman #endif 487397b6df1SKris Buschelman } 488397b6df1SKris Buschelman } else { 489397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 490397b6df1SKris Buschelman lu->id.a_loc = (mumps_double_complex*)lu->val; 491397b6df1SKris Buschelman #else 492397b6df1SKris Buschelman lu->id.a_loc = lu->val; 493397b6df1SKris Buschelman #endif 494397b6df1SKris Buschelman } 495397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX) 496397b6df1SKris Buschelman zmumps_c(&lu->id); 497397b6df1SKris Buschelman #else 498397b6df1SKris Buschelman dmumps_c(&lu->id); 499397b6df1SKris Buschelman #endif 500397b6df1SKris Buschelman if (lu->id.INFOG(1) < 0) { 50119facb7aSBarry Smith if (lu->id.INFO(1) == -13) { 50219facb7aSBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2)); 50319facb7aSBarry Smith } else { 50479a5c55eSBarry 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)); 505397b6df1SKris Buschelman } 50619facb7aSBarry Smith } 507397b6df1SKris Buschelman 50819facb7aSBarry Smith if (!lu->myid && lu->id.ICNTL(16) > 0){ 50979a5c55eSBarry Smith SETERRQ1(PETSC_ERR_LIB," lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16)); 510397b6df1SKris Buschelman } 511397b6df1SKris Buschelman 5128ada1bb4SHong Zhang if (lu->size > 1){ 513719d5645SBarry Smith if ((F)->factor == MAT_FACTOR_LU){ 514719d5645SBarry Smith F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 515e09efc27SHong Zhang } else { 516719d5645SBarry Smith F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 517e09efc27SHong Zhang } 518e09efc27SHong Zhang F_diag->assembled = PETSC_TRUE; 519329ec9b3SHong Zhang if (lu->nSolve){ 520329ec9b3SHong Zhang ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 521329ec9b3SHong Zhang ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr); 522329ec9b3SHong Zhang ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 523329ec9b3SHong Zhang } 5248ada1bb4SHong Zhang } 525719d5645SBarry Smith (F)->assembled = PETSC_TRUE; 526397b6df1SKris Buschelman lu->matstruc = SAME_NONZERO_PATTERN; 527ace87b0dSHong Zhang lu->CleanUpMUMPS = PETSC_TRUE; 528329ec9b3SHong Zhang lu->nSolve = 0; 529397b6df1SKris Buschelman PetscFunctionReturn(0); 530397b6df1SKris Buschelman } 531397b6df1SKris Buschelman 532b24902e0SBarry Smith 533397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */ 534397b6df1SKris Buschelman #undef __FUNCT__ 535f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS" 5360481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 537b24902e0SBarry Smith { 538719d5645SBarry Smith Mat_MUMPS *lu = (Mat_MUMPS*)F->spptr; 539397b6df1SKris Buschelman 540397b6df1SKris Buschelman PetscFunctionBegin; 541b24902e0SBarry Smith lu->sym = 0; 542b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 543719d5645SBarry Smith F->ops->lufactornumeric = MatFactorNumeric_MUMPS; 544b24902e0SBarry Smith PetscFunctionReturn(0); 545b24902e0SBarry Smith } 546b24902e0SBarry Smith 547b24902e0SBarry Smith 548397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */ 549397b6df1SKris Buschelman #undef __FUNCT__ 550f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS" 5510481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info) 552b24902e0SBarry Smith { 553719d5645SBarry Smith Mat_MUMPS *lu = (Mat_MUMPS*)(F)->spptr; 554397b6df1SKris Buschelman 555397b6df1SKris Buschelman PetscFunctionBegin; 556b24902e0SBarry Smith lu->sym = 2; 557b24902e0SBarry Smith lu->matstruc = DIFFERENT_NONZERO_PATTERN; 558719d5645SBarry Smith (F)->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS; 559db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX) 560719d5645SBarry Smith (F)->ops->getinertia = MatGetInertia_SBAIJMUMPS; 561db4efbfdSBarry Smith #endif 562b24902e0SBarry Smith PetscFunctionReturn(0); 563b24902e0SBarry Smith } 564b24902e0SBarry Smith 565397b6df1SKris Buschelman #undef __FUNCT__ 566f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS" 567f6c57405SHong Zhang PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) { 568f6c57405SHong Zhang Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr; 569f6c57405SHong Zhang PetscErrorCode ierr; 570f6c57405SHong Zhang 571f6c57405SHong Zhang PetscFunctionBegin; 572f6c57405SHong Zhang /* check if matrix is mumps type */ 573f6c57405SHong Zhang if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0); 574f6c57405SHong Zhang 575f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr); 576f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," SYM (matrix type): %d \n",lu->id.sym);CHKERRQ(ierr); 577f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," PAR (host participation): %d \n",lu->id.par);CHKERRQ(ierr); 578f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(1) (output for error): %d \n",lu->id.ICNTL(1));CHKERRQ(ierr); 579f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr); 580f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(3) (output for global info): %d \n",lu->id.ICNTL(3));CHKERRQ(ierr); 581f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(4) (level of printing): %d \n",lu->id.ICNTL(4));CHKERRQ(ierr); 582f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(5) (input mat struct): %d \n",lu->id.ICNTL(5));CHKERRQ(ierr); 583f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(6) (matrix prescaling): %d \n",lu->id.ICNTL(6));CHKERRQ(ierr); 584f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(7) (matrix ordering): %d \n",lu->id.ICNTL(7));CHKERRQ(ierr); 585f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(8) (scalling strategy): %d \n",lu->id.ICNTL(8));CHKERRQ(ierr); 586f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(9) (A/A^T x=b is solved): %d \n",lu->id.ICNTL(9));CHKERRQ(ierr); 587f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr); 588f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(11) (error analysis): %d \n",lu->id.ICNTL(11));CHKERRQ(ierr); 589f6c57405SHong Zhang if (!lu->myid && lu->id.ICNTL(11)>0) { 590f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(4) (inf norm of input mat): %g\n",lu->id.RINFOG(4));CHKERRQ(ierr); 591f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(5) (inf norm of solution): %g\n",lu->id.RINFOG(5));CHKERRQ(ierr); 592f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(6) (inf norm of residual): %g\n",lu->id.RINFOG(6));CHKERRQ(ierr); 593f6c57405SHong 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); 594f6c57405SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," RINFOG(9) (error estimate): %g \n",lu->id.RINFOG(9));CHKERRQ(ierr); 595f6c57405SHong 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); 596f6c57405SHong Zhang 597f6c57405SHong Zhang } 598f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(12) (efficiency control): %d \n",lu->id.ICNTL(12));CHKERRQ(ierr); 599f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(13) (efficiency control): %d \n",lu->id.ICNTL(13));CHKERRQ(ierr); 600f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr); 601f6c57405SHong Zhang /* ICNTL(15-17) not used */ 602f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(18) (input mat struct): %d \n",lu->id.ICNTL(18));CHKERRQ(ierr); 603f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(19) (Shur complement info): %d \n",lu->id.ICNTL(19));CHKERRQ(ierr); 604f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(20) (rhs sparse pattern): %d \n",lu->id.ICNTL(20));CHKERRQ(ierr); 605f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ICNTL(21) (solution struct): %d \n",lu->id.ICNTL(21));CHKERRQ(ierr); 606f6c57405SHong Zhang 607f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(1) (relative pivoting threshold): %g \n",lu->id.CNTL(1));CHKERRQ(ierr); 608f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr); 609f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(3) (absolute pivoting threshold): %g \n",lu->id.CNTL(3));CHKERRQ(ierr); 610f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," CNTL(4) (value of static pivoting): %g \n",lu->id.CNTL(4));CHKERRQ(ierr); 611f6c57405SHong Zhang 612f6c57405SHong Zhang /* infomation local to each processor */ 613f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);} 6147adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr); 6157adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 616f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);} 6177adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr); 6187adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 619f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);} 6207adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr); 6217adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 622f6c57405SHong Zhang /* 623f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);} 6247adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr); 6257adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 626f6c57405SHong Zhang */ 627f6c57405SHong Zhang 628f6c57405SHong 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);} 6297adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr); 6307adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 631f6c57405SHong Zhang 632f6c57405SHong 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);} 6337adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr); 6347adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 635f6c57405SHong Zhang 636f6c57405SHong Zhang if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, " INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);} 6377adad957SLisandro Dalcin ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm," [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr); 6387adad957SLisandro Dalcin ierr = PetscSynchronizedFlush(((PetscObject)A)->comm); 639f6c57405SHong Zhang 640f6c57405SHong Zhang if (!lu->myid){ /* information from the host */ 641f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr); 642f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr); 643f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr); 644f6c57405SHong Zhang 645f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr); 646f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr); 647f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr); 648f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr); 649f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr); 650f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr); 651f6c57405SHong 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); 652f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr); 653f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr); 654f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr); 655f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr); 656f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr); 657f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr); 658f6c57405SHong 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); 659f6c57405SHong 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); 660f6c57405SHong 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); 661f6c57405SHong 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); 662f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr); 663f6c57405SHong 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); 664f6c57405SHong 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); 665f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr); 666f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr); 667f6c57405SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr); 668f6c57405SHong Zhang } 669f6c57405SHong Zhang 670f6c57405SHong Zhang PetscFunctionReturn(0); 671f6c57405SHong Zhang } 672f6c57405SHong Zhang 673f6c57405SHong Zhang #undef __FUNCT__ 674f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS" 675b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) 676b24902e0SBarry Smith { 677f6c57405SHong Zhang PetscErrorCode ierr; 678f6c57405SHong Zhang PetscTruth iascii; 679f6c57405SHong Zhang PetscViewerFormat format; 680f6c57405SHong Zhang 681f6c57405SHong Zhang PetscFunctionBegin; 682f6c57405SHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 683f6c57405SHong Zhang if (iascii) { 684f6c57405SHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 685f6c57405SHong Zhang if (format == PETSC_VIEWER_ASCII_INFO){ 686f6c57405SHong Zhang ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr); 687f6c57405SHong Zhang } 688f6c57405SHong Zhang } 689f6c57405SHong Zhang PetscFunctionReturn(0); 690f6c57405SHong Zhang } 691f6c57405SHong Zhang 692397b6df1SKris Buschelman 69324b6179bSKris Buschelman /*MC 694fafad747SKris Buschelman MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed 69524b6179bSKris Buschelman and sequential matrices via the external package MUMPS. 69624b6179bSKris Buschelman 69724b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 69824b6179bSKris Buschelman on how to declare the existence of external packages), 69924b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 700175b88e8SBarry Smith After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS), then 701175b88e8SBarry Smith optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT 702175b88e8SBarry Smith call MatCreateSeqAIJ/MPIAIJ() directly or the preallocation information will be LOST! 70324b6179bSKris Buschelman 70424b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQAIJ. 70524b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPIAIJ. Hence for single process communicators, 7063ec795f1SBarry Smith MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported 70724b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 7083ec795f1SBarry Smith the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 70928b08bd3SKris Buschelman conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size) 710175b88e8SBarry Smith without data copy AFTER the matrix values are set. 71124b6179bSKris Buschelman 71224b6179bSKris Buschelman Options Database Keys: 7130bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions() 71424b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 71524b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level 71624b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 71724b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 71824b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 71924b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 72094b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 72124b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 72224b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 72324b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 72424b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 72524b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 72624b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 72724b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 72824b6179bSKris Buschelman 72924b6179bSKris Buschelman Level: beginner 73024b6179bSKris Buschelman 73124b6179bSKris Buschelman .seealso: MATSBAIJMUMPS 73224b6179bSKris Buschelman M*/ 73324b6179bSKris Buschelman 734f0c56d0fSKris Buschelman 73524b6179bSKris Buschelman /*MC 736fafad747SKris Buschelman MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for 73724b6179bSKris Buschelman distributed and sequential matrices via the external package MUMPS. 73824b6179bSKris Buschelman 73924b6179bSKris Buschelman If MUMPS is installed (see the manual for instructions 74024b6179bSKris Buschelman on how to declare the existence of external packages), 74124b6179bSKris Buschelman a matrix type can be constructed which invokes MUMPS solvers. 742175b88e8SBarry Smith After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS), then 743175b88e8SBarry Smith optionally call MatSeqSBAIJSetPreallocation() or MatMPISBAIJSetPreallocation() DO NOT 744175b88e8SBarry Smith call MatCreateSeqSBAIJ/MPISBAIJ() directly or the preallocation information will be LOST! 74524b6179bSKris Buschelman 74624b6179bSKris Buschelman If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ. 74724b6179bSKris Buschelman Otherwise, this matrix type inherits from MATMPISBAIJ. Hence for single process communicators, 748175b88e8SBarry Smith MatSeqSBAIJSetPreallocation() is supported, and similarly MatMPISBAIJSetPreallocation() is supported 74924b6179bSKris Buschelman for communicators controlling multiple processes. It is recommended that you call both of 750175b88e8SBarry Smith the above preallocation routines for simplicity. One can also call MatConvert() for an inplace 75128b08bd3SKris Buschelman conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size) 752175b88e8SBarry Smith without data copy AFTER the matrix values have been set. 75324b6179bSKris Buschelman 75424b6179bSKris Buschelman Options Database Keys: 7550bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions() 75624b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric 75724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level 75824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide) 75924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide) 76024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T 76124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements 76294b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view 76324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide) 76424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide) 76524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide) 76624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide) 76724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold 76824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement 76924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold 77024b6179bSKris Buschelman 77124b6179bSKris Buschelman Level: beginner 77224b6179bSKris Buschelman 77324b6179bSKris Buschelman .seealso: MATAIJMUMPS 77424b6179bSKris Buschelman M*/ 77524b6179bSKris Buschelman 7762877fffaSHong Zhang EXTERN_C_BEGIN 7772877fffaSHong Zhang /* 7782877fffaSHong Zhang The seq and mpi versions of this function are the same 7792877fffaSHong Zhang */ 7802877fffaSHong Zhang #undef __FUNCT__ 7812877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_mumps" 7822877fffaSHong Zhang PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F) 7832877fffaSHong Zhang { 7842877fffaSHong Zhang Mat B; 7852877fffaSHong Zhang PetscErrorCode ierr; 7862877fffaSHong Zhang Mat_MUMPS *mumps; 7872877fffaSHong Zhang 7882877fffaSHong Zhang PetscFunctionBegin; 7892877fffaSHong Zhang if (ftype != MAT_FACTOR_LU) { 7902877fffaSHong Zhang SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix"); 7912877fffaSHong Zhang } 7922877fffaSHong Zhang /* Create the factorization matrix */ 7932877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 7942877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 7952877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 7962877fffaSHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 7972877fffaSHong Zhang 7982877fffaSHong Zhang B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 7992877fffaSHong Zhang B->ops->view = MatView_MUMPS; 8002877fffaSHong Zhang B->factor = MAT_FACTOR_LU; 8012877fffaSHong Zhang 8022877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 8032877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 8042877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 8052877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 8062877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 8072877fffaSHong Zhang mumps->nSolve = 0; 8082877fffaSHong Zhang mumps->MatDestroy = B->ops->destroy; 8092877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 8102877fffaSHong Zhang B->spptr = (void*)mumps; 8112877fffaSHong Zhang 8122877fffaSHong Zhang *F = B; 8132877fffaSHong Zhang PetscFunctionReturn(0); 8142877fffaSHong Zhang } 8152877fffaSHong Zhang EXTERN_C_END 8162877fffaSHong Zhang 8172877fffaSHong Zhang EXTERN_C_BEGIN 8182877fffaSHong Zhang #undef __FUNCT__ 8192877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_mpiaij_mumps" 8202877fffaSHong Zhang PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F) 8212877fffaSHong Zhang { 8222877fffaSHong Zhang Mat B; 8232877fffaSHong Zhang PetscErrorCode ierr; 8242877fffaSHong Zhang Mat_MUMPS *mumps; 8252877fffaSHong Zhang 8262877fffaSHong Zhang PetscFunctionBegin; 8272877fffaSHong Zhang if (ftype != MAT_FACTOR_LU) { 8282877fffaSHong Zhang SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix"); 8292877fffaSHong Zhang } 8302877fffaSHong Zhang /* Create the factorization matrix */ 8312877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 8322877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 8332877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 8342877fffaSHong Zhang ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 8352877fffaSHong Zhang ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 8362877fffaSHong Zhang 8372877fffaSHong Zhang B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS; 8382877fffaSHong Zhang B->ops->view = MatView_MUMPS; 8392877fffaSHong Zhang B->factor = MAT_FACTOR_LU; 8402877fffaSHong Zhang 8412877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 8422877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 8432877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 8442877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 8452877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 8462877fffaSHong Zhang mumps->nSolve = 0; 847*f3c0ef26SHong Zhang mumps->MatDestroy = B->ops->destroy; 848*f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 8492877fffaSHong Zhang B->spptr = (void*)mumps; 8502877fffaSHong Zhang 8512877fffaSHong Zhang *F = B; 8522877fffaSHong Zhang PetscFunctionReturn(0); 8532877fffaSHong Zhang } 8542877fffaSHong Zhang EXTERN_C_END 8552877fffaSHong Zhang 8562877fffaSHong Zhang EXTERN_C_BEGIN 8572877fffaSHong Zhang #undef __FUNCT__ 8582877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_seqsbaij_mumps" 8592877fffaSHong Zhang PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 8602877fffaSHong Zhang { 8612877fffaSHong Zhang Mat B; 8622877fffaSHong Zhang PetscErrorCode ierr; 8632877fffaSHong Zhang Mat_MUMPS *mumps; 8642877fffaSHong Zhang 8652877fffaSHong Zhang PetscFunctionBegin; 8662877fffaSHong Zhang if (ftype != MAT_FACTOR_CHOLESKY) { 8672877fffaSHong Zhang SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 8682877fffaSHong Zhang } 8692877fffaSHong Zhang /* Create the factorization matrix */ 8702877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 8712877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 8722877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 8732877fffaSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 8742877fffaSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 8752877fffaSHong Zhang 8762877fffaSHong Zhang B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 8772877fffaSHong Zhang B->ops->view = MatView_MUMPS; 8782877fffaSHong Zhang B->factor = MAT_FACTOR_CHOLESKY; 8792877fffaSHong Zhang 8802877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 8812877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 8822877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 8832877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 8842877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 8852877fffaSHong Zhang mumps->nSolve = 0; 8862877fffaSHong Zhang mumps->MatDestroy = B->ops->destroy; 8872877fffaSHong Zhang B->ops->destroy = MatDestroy_MUMPS; 8882877fffaSHong Zhang B->spptr = (void*)mumps; 889*f3c0ef26SHong Zhang 8902877fffaSHong Zhang *F = B; 8912877fffaSHong Zhang PetscFunctionReturn(0); 8922877fffaSHong Zhang } 8932877fffaSHong Zhang EXTERN_C_END 8942877fffaSHong Zhang 8952877fffaSHong Zhang EXTERN_C_BEGIN 8962877fffaSHong Zhang #undef __FUNCT__ 8972877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_mpisbaij_mumps" 8982877fffaSHong Zhang PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F) 8992877fffaSHong Zhang { 9002877fffaSHong Zhang Mat B; 9012877fffaSHong Zhang PetscErrorCode ierr; 9022877fffaSHong Zhang Mat_MUMPS *mumps; 9032877fffaSHong Zhang 9042877fffaSHong Zhang PetscFunctionBegin; 9052877fffaSHong Zhang if (ftype != MAT_FACTOR_CHOLESKY) { 9062877fffaSHong Zhang SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix"); 9072877fffaSHong Zhang } 9082877fffaSHong Zhang /* Create the factorization matrix */ 9092877fffaSHong Zhang ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 9102877fffaSHong Zhang ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 9112877fffaSHong Zhang ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 9122877fffaSHong Zhang ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 9132877fffaSHong Zhang ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 9142877fffaSHong Zhang 9152877fffaSHong Zhang B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS; 9162877fffaSHong Zhang B->ops->view = MatView_MUMPS; 9172877fffaSHong Zhang B->factor = MAT_FACTOR_CHOLESKY; 9182877fffaSHong Zhang 9192877fffaSHong Zhang ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr); 9202877fffaSHong Zhang mumps->CleanUpMUMPS = PETSC_FALSE; 9212877fffaSHong Zhang mumps->isAIJ = PETSC_TRUE; 9222877fffaSHong Zhang mumps->scat_rhs = PETSC_NULL; 9232877fffaSHong Zhang mumps->scat_sol = PETSC_NULL; 9242877fffaSHong Zhang mumps->nSolve = 0; 925*f3c0ef26SHong Zhang mumps->MatDestroy = B->ops->destroy; 926*f3c0ef26SHong Zhang B->ops->destroy = MatDestroy_MUMPS; 9272877fffaSHong Zhang B->spptr = (void*)mumps; 928*f3c0ef26SHong Zhang 9292877fffaSHong Zhang *F = B; 9302877fffaSHong Zhang PetscFunctionReturn(0); 9312877fffaSHong Zhang } 9322877fffaSHong Zhang EXTERN_C_END 933