xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 4e1dd20e4075d61f5b1040c3e2971bc2cfb8b48e)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
3397b6df1SKris Buschelman /*
4c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
5397b6df1SKris Buschelman */
6*4e1dd20eSHong 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"
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   PetscFunctionReturn(0);
143397b6df1SKris Buschelman }
144397b6df1SKris Buschelman 
145397b6df1SKris Buschelman #undef __FUNCT__
1463924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
147dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
148dfbe8321SBarry Smith {
149f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
150dfbe8321SBarry Smith   PetscErrorCode ierr;
151c1490034SHong Zhang   PetscMPIInt    size=lu->size;
152b24902e0SBarry Smith 
153397b6df1SKris Buschelman   PetscFunctionBegin;
154397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
155397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
156329ec9b3SHong Zhang     if (size > 1){
15768653410SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
158329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
159329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
1602750af12SHong Zhang       if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);}
1612750af12SHong Zhang       if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);}
162329ec9b3SHong Zhang       ierr = PetscFree(lu->val);CHKERRQ(ierr);
163329ec9b3SHong Zhang     }
164397b6df1SKris Buschelman     lu->id.job=JOB_END;
165397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
166397b6df1SKris Buschelman     zmumps_c(&lu->id);
167397b6df1SKris Buschelman #else
168397b6df1SKris Buschelman     dmumps_c(&lu->id);
169397b6df1SKris Buschelman #endif
170c338a77dSKris Buschelman     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
171c338a77dSKris Buschelman     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
172397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
173397b6df1SKris Buschelman   }
17497969023SHong Zhang   /* clear composed functions */
17597969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatFactorGetSolverPackage_C","",PETSC_NULL);CHKERRQ(ierr);
17697969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMumpsSetIcntl_C","",PETSC_NULL);CHKERRQ(ierr);
17767334b25SHong Zhang   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
178397b6df1SKris Buschelman   PetscFunctionReturn(0);
179397b6df1SKris Buschelman }
180397b6df1SKris Buschelman 
181397b6df1SKris Buschelman #undef __FUNCT__
182f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
183b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
184b24902e0SBarry Smith {
185f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
186d54de34fSKris Buschelman   PetscScalar    *array;
187397b6df1SKris Buschelman   Vec            x_seq;
188329ec9b3SHong Zhang   IS             is_iden,is_petsc;
189dfbe8321SBarry Smith   PetscErrorCode ierr;
190329ec9b3SHong Zhang   PetscInt       i;
191397b6df1SKris Buschelman 
192397b6df1SKris Buschelman   PetscFunctionBegin;
193329ec9b3SHong Zhang   lu->id.nrhs = 1;
194329ec9b3SHong Zhang   x_seq = lu->b_seq;
195397b6df1SKris Buschelman   if (lu->size > 1){
196329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
197f6cfb2d1SLisandro Dalcin     ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
198f6cfb2d1SLisandro Dalcin     ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
199397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
200397b6df1SKris Buschelman   } else {  /* size == 1 */
201397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
202397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
203397b6df1SKris Buschelman   }
204397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
2058278f211SHong Zhang     lu->id.nrhs = 1;
206397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
207397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
208397b6df1SKris Buschelman #else
209397b6df1SKris Buschelman     lu->id.rhs = array;
210397b6df1SKris Buschelman #endif
211397b6df1SKris Buschelman   }
212329ec9b3SHong Zhang   if (lu->size == 1){
213329ec9b3SHong Zhang     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
214329ec9b3SHong Zhang   } else if (!lu->myid){
215329ec9b3SHong Zhang     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
216329ec9b3SHong Zhang   }
217329ec9b3SHong Zhang 
218329ec9b3SHong Zhang   if (lu->size > 1){
219329ec9b3SHong Zhang     /* distributed solution */
220329ec9b3SHong Zhang     lu->id.ICNTL(21) = 1;
221329ec9b3SHong Zhang     if (!lu->nSolve){
222329ec9b3SHong Zhang       /* Create x_seq=sol_loc for repeated use */
223329ec9b3SHong Zhang       PetscInt    lsol_loc;
224329ec9b3SHong Zhang       PetscScalar *sol_loc;
225329ec9b3SHong Zhang       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
2260e83c824SBarry Smith       ierr = PetscMalloc2(lsol_loc,PetscScalar,&sol_loc,lsol_loc,PetscInt,&lu->id.isol_loc);CHKERRQ(ierr);
227329ec9b3SHong Zhang       lu->id.lsol_loc = lsol_loc;
22844ea04b1SSatish Balay #if defined(PETSC_USE_COMPLEX)
22944ea04b1SSatish Balay       lu->id.sol_loc  = (mumps_double_complex*)sol_loc;
23044ea04b1SSatish Balay #else
23144ea04b1SSatish Balay       lu->id.sol_loc  = sol_loc;
23244ea04b1SSatish Balay #endif
233329ec9b3SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
234329ec9b3SHong Zhang     }
235329ec9b3SHong Zhang   }
236397b6df1SKris Buschelman 
237397b6df1SKris Buschelman   /* solve phase */
238329ec9b3SHong Zhang   /*-------------*/
239397b6df1SKris Buschelman   lu->id.job = 3;
240397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
241397b6df1SKris Buschelman   zmumps_c(&lu->id);
242397b6df1SKris Buschelman #else
243397b6df1SKris Buschelman   dmumps_c(&lu->id);
244397b6df1SKris Buschelman #endif
245397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
24679a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
247397b6df1SKris Buschelman   }
248397b6df1SKris Buschelman 
249329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
250329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
251329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
252329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
253329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
254397b6df1SKris Buschelman       }
255329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
256329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
257329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
258329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
259397b6df1SKris Buschelman     }
260ca9f406cSSatish Balay     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
261ca9f406cSSatish Balay     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
262329ec9b3SHong Zhang   }
263329ec9b3SHong Zhang   lu->nSolve++;
264397b6df1SKris Buschelman   PetscFunctionReturn(0);
265397b6df1SKris Buschelman }
266397b6df1SKris Buschelman 
267ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
268a58c3f20SHong Zhang /*
269a58c3f20SHong Zhang   input:
270a58c3f20SHong Zhang    F:        numeric factor
271a58c3f20SHong Zhang   output:
272a58c3f20SHong Zhang    nneg:     total number of negative pivots
273a58c3f20SHong Zhang    nzero:    0
274a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
275a58c3f20SHong Zhang */
276a58c3f20SHong Zhang 
277a58c3f20SHong Zhang #undef __FUNCT__
278a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
279dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
280a58c3f20SHong Zhang {
281a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
282dfbe8321SBarry Smith   PetscErrorCode ierr;
283c1490034SHong Zhang   PetscMPIInt    size;
284a58c3f20SHong Zhang 
285a58c3f20SHong Zhang   PetscFunctionBegin;
2867adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr);
287bcb30aebSHong 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 */
288bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
28979a5c55eSBarry 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));
290bcb30aebSHong Zhang   }
291a58c3f20SHong Zhang   if (nneg){
292a58c3f20SHong Zhang     if (!lu->myid){
293a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
294a58c3f20SHong Zhang     }
295bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
296a58c3f20SHong Zhang   }
297a58c3f20SHong Zhang   if (nzero) *nzero = 0;
298d0f46423SBarry Smith   if (npos)  *npos  = F->rmap->N - (*nneg);
299a58c3f20SHong Zhang   PetscFunctionReturn(0);
300a58c3f20SHong Zhang }
301ace3df97SHong Zhang #endif /* !defined(PETSC_USE_COMPLEX) */
302a58c3f20SHong Zhang 
303397b6df1SKris Buschelman #undef __FUNCT__
304f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
3050481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
306af281ebdSHong Zhang {
307719d5645SBarry Smith   Mat_MUMPS      *lu =(Mat_MUMPS*)(F)->spptr;
3086849ba73SBarry Smith   PetscErrorCode ierr;
309d0f46423SBarry Smith   PetscInt       rnz,nnz,nz=0,i,M=A->rmap->N,*ai,*aj,icntl;
310397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
311e09efc27SHong Zhang   Mat            F_diag;
312c349612cSHong Zhang   IS             is_iden;
313c349612cSHong Zhang   Vec            b;
3145c9eb25fSBarry Smith   PetscTruth     isSeqAIJ,isSeqSBAIJ;
315397b6df1SKris Buschelman 
316397b6df1SKris Buschelman   PetscFunctionBegin;
3175c9eb25fSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
3185c9eb25fSBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
319397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
320719d5645SBarry Smith     (F)->ops->solve   = MatSolve_MUMPS;
321397b6df1SKris Buschelman 
322397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
3237adad957SLisandro Dalcin     ierr = MPI_Comm_rank(((PetscObject)A)->comm, &lu->myid);
3247adad957SLisandro Dalcin     ierr = MPI_Comm_size(((PetscObject)A)->comm,&lu->size);CHKERRQ(ierr);
325397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
3267adad957SLisandro Dalcin     ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->comm_mumps));CHKERRQ(ierr);
3276a1dac61SBarry Smith     lu->id.comm_fortran = MPI_Comm_c2f(lu->comm_mumps);
328397b6df1SKris Buschelman 
329397b6df1SKris Buschelman     /* Set mumps options */
3307adad957SLisandro Dalcin     ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
331397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
332397b6df1SKris Buschelman     lu->id.sym=lu->sym;
333397b6df1SKris Buschelman     if (lu->sym == 2){
334397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
335397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
336397b6df1SKris Buschelman     }
337397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
338397b6df1SKris Buschelman     zmumps_c(&lu->id);
339397b6df1SKris Buschelman #else
340397b6df1SKris Buschelman     dmumps_c(&lu->id);
341397b6df1SKris Buschelman #endif
342397b6df1SKris Buschelman 
3435c9eb25fSBarry Smith     if (isSeqAIJ || isSeqSBAIJ){
344397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
345397b6df1SKris Buschelman     } else {
346397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
347397b6df1SKris Buschelman     }
348397b6df1SKris Buschelman 
349397b6df1SKris Buschelman     icntl=-1;
350c0165424SHong Zhang     lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
351397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
35219facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
353397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
354397b6df1SKris Buschelman     } else { /* no output */
355397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
3563823f358SHong Zhang       lu->id.ICNTL(2) = 0;  /* output stream for diagnostic printing, statistics, and warning. default=0 */
3573823f358SHong Zhang       lu->id.ICNTL(3) = 0; /* output stream for global information, default=6 */
358397b6df1SKris Buschelman     }
3593823f358SHong Zhang     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);
360397b6df1SKris Buschelman     icntl=-1;
361397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
362397b6df1SKris Buschelman     if (flg) {
363397b6df1SKris Buschelman       if (icntl== 1){
364397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
365397b6df1SKris Buschelman       } else {
366397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
367397b6df1SKris Buschelman       }
368397b6df1SKris Buschelman     }
3693823f358SHong Zhang     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);
370397b6df1SKris 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);
371397b6df1SKris 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);
372c0165424SHong Zhang     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);
3733823f358SHong Zhang     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);
3743823f358SHong Zhang     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);
375adc1d99fSHong 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);
3763823f358SHong Zhang     ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): Schur complement","None",lu->id.ICNTL(19),&lu->id.ICNTL(19),PETSC_NULL);CHKERRQ(ierr);
3773823f358SHong Zhang 
378c0165424SHong Zhang     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);
379c0165424SHong Zhang     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);
380c0165424SHong Zhang     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);
381c0165424SHong Zhang     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);
382c0165424SHong Zhang     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);
383c0165424SHong Zhang     ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): experimental parameter","None",lu->id.ICNTL(27),&lu->id.ICNTL(27),PETSC_NULL);CHKERRQ(ierr);
384397b6df1SKris Buschelman 
385397b6df1SKris 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);
386397b6df1SKris 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);
387397b6df1SKris 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);
38825f9c88cSHong 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);
389c0165424SHong Zhang     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);
390397b6df1SKris Buschelman     PetscOptionsEnd();
391397b6df1SKris Buschelman   }
392397b6df1SKris Buschelman 
393397b6df1SKris Buschelman   /* define matrix A */
394397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
395397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
396397b6df1SKris Buschelman     if (!lu->myid) {
3975c9eb25fSBarry Smith       if (isSeqAIJ){
398397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
399397b6df1SKris Buschelman         nz               = aa->nz;
400397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
4015c9eb25fSBarry Smith       } else if (isSeqSBAIJ) {
402397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
4036c6c5352SBarry Smith         nz                  =  aa->nz;
404397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
4055c9eb25fSBarry Smith       } else {
4065c9eb25fSBarry Smith         SETERRQ(PETSC_ERR_SUP,"No mumps factorization for this matrix type");
407397b6df1SKris Buschelman       }
408397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
4097c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
4107c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
411397b6df1SKris Buschelman         nz = 0;
412397b6df1SKris Buschelman         for (i=0; i<M; i++){
413397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
414397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
415397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
416397b6df1SKris Buschelman           }
417397b6df1SKris Buschelman         }
418397b6df1SKris Buschelman       }
419397b6df1SKris Buschelman     }
420397b6df1SKris Buschelman     break;
421397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
422397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
423397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
424397b6df1SKris Buschelman     } else {
425397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
426397b6df1SKris Buschelman     }
427397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
428397b6df1SKris Buschelman     break;
429397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
430397b6df1SKris Buschelman   }
431397b6df1SKris Buschelman 
432397b6df1SKris Buschelman   /* analysis phase */
433329ec9b3SHong Zhang   /*----------------*/
434397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
435329ec9b3SHong Zhang     lu->id.job = 1;
436329ec9b3SHong Zhang 
437397b6df1SKris Buschelman     lu->id.n = M;
438397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
439397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
440397b6df1SKris Buschelman       if (!lu->myid) {
441397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
442397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
443397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
444397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
445397b6df1SKris Buschelman #else
446397b6df1SKris Buschelman           lu->id.a = lu->val;
447397b6df1SKris Buschelman #endif
448397b6df1SKris Buschelman         }
449397b6df1SKris Buschelman       }
450397b6df1SKris Buschelman       break;
451397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
452397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
453397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
454397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
455397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
456397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
457397b6df1SKris Buschelman #else
458397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
459397b6df1SKris Buschelman #endif
460397b6df1SKris Buschelman       }
461329ec9b3SHong Zhang       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
462329ec9b3SHong Zhang       if (!lu->myid){
463d0f46423SBarry Smith         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
464d0f46423SBarry Smith         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
465329ec9b3SHong Zhang       } else {
466329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
467329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
468329ec9b3SHong Zhang       }
4697adad957SLisandro Dalcin       ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
470d0f46423SBarry Smith       ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
471329ec9b3SHong Zhang       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
472329ec9b3SHong Zhang 
473329ec9b3SHong Zhang       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
474329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
475329ec9b3SHong Zhang       ierr = VecDestroy(b);CHKERRQ(ierr);
476397b6df1SKris Buschelman       break;
477397b6df1SKris Buschelman     }
478397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
479397b6df1SKris Buschelman     zmumps_c(&lu->id);
480397b6df1SKris Buschelman #else
481397b6df1SKris Buschelman     dmumps_c(&lu->id);
482397b6df1SKris Buschelman #endif
483397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
48479a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
485397b6df1SKris Buschelman     }
486397b6df1SKris Buschelman   }
487397b6df1SKris Buschelman 
488397b6df1SKris Buschelman   /* numerical factorization phase */
489329ec9b3SHong Zhang   /*-------------------------------*/
490329ec9b3SHong Zhang   lu->id.job = 2;
491958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
492a7aca84bSHong Zhang     if (!lu->myid) {
493397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
494397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
495397b6df1SKris Buschelman #else
496397b6df1SKris Buschelman       lu->id.a = lu->val;
497397b6df1SKris Buschelman #endif
498397b6df1SKris Buschelman     }
499397b6df1SKris Buschelman   } else {
500397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
501397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
502397b6df1SKris Buschelman #else
503397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
504397b6df1SKris Buschelman #endif
505397b6df1SKris Buschelman   }
506397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
507397b6df1SKris Buschelman   zmumps_c(&lu->id);
508397b6df1SKris Buschelman #else
509397b6df1SKris Buschelman   dmumps_c(&lu->id);
510397b6df1SKris Buschelman #endif
511397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
51219facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
51319facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
51419facb7aSBarry Smith     } else {
51579a5c55eSBarry 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));
516397b6df1SKris Buschelman     }
51719facb7aSBarry Smith   }
518397b6df1SKris Buschelman 
51919facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
52079a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
521397b6df1SKris Buschelman   }
522397b6df1SKris Buschelman 
5238ada1bb4SHong Zhang   if (lu->size > 1){
524719d5645SBarry Smith     if ((F)->factor == MAT_FACTOR_LU){
525719d5645SBarry Smith       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
526e09efc27SHong Zhang     } else {
527719d5645SBarry Smith       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
528e09efc27SHong Zhang     }
529e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
530329ec9b3SHong Zhang     if (lu->nSolve){
531329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
5320e83c824SBarry Smith       ierr = PetscFree2(lu->id.sol_loc,lu->id.isol_loc);CHKERRQ(ierr);
533329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
534329ec9b3SHong Zhang     }
5358ada1bb4SHong Zhang   }
536719d5645SBarry Smith   (F)->assembled   = PETSC_TRUE;
537397b6df1SKris Buschelman   lu->matstruc      = SAME_NONZERO_PATTERN;
538ace87b0dSHong Zhang   lu->CleanUpMUMPS  = PETSC_TRUE;
539329ec9b3SHong Zhang   lu->nSolve        = 0;
540397b6df1SKris Buschelman   PetscFunctionReturn(0);
541397b6df1SKris Buschelman }
542397b6df1SKris Buschelman 
543397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
544397b6df1SKris Buschelman #undef __FUNCT__
545f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
5460481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
547b24902e0SBarry Smith {
548719d5645SBarry Smith   Mat_MUMPS      *lu = (Mat_MUMPS*)F->spptr;
549397b6df1SKris Buschelman 
550397b6df1SKris Buschelman   PetscFunctionBegin;
551b24902e0SBarry Smith   lu->sym                  = 0;
552b24902e0SBarry Smith   lu->matstruc             = DIFFERENT_NONZERO_PATTERN;
553719d5645SBarry Smith   F->ops->lufactornumeric  = MatFactorNumeric_MUMPS;
554b24902e0SBarry Smith   PetscFunctionReturn(0);
555b24902e0SBarry Smith }
556b24902e0SBarry Smith 
557b24902e0SBarry Smith 
558397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
559397b6df1SKris Buschelman #undef __FUNCT__
560f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
5610481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
562b24902e0SBarry Smith {
563719d5645SBarry Smith   Mat_MUMPS      *lu = (Mat_MUMPS*)(F)->spptr;
564397b6df1SKris Buschelman 
565397b6df1SKris Buschelman   PetscFunctionBegin;
566b24902e0SBarry Smith   lu->sym                          = 2;
567b24902e0SBarry Smith   lu->matstruc                     = DIFFERENT_NONZERO_PATTERN;
568719d5645SBarry Smith   (F)->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
569db4efbfdSBarry Smith #if !defined(PETSC_USE_COMPLEX)
570719d5645SBarry Smith   (F)->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
571db4efbfdSBarry Smith #endif
572b24902e0SBarry Smith   PetscFunctionReturn(0);
573b24902e0SBarry Smith }
574b24902e0SBarry Smith 
575397b6df1SKris Buschelman #undef __FUNCT__
576f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS"
57774ed9c26SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer)
57874ed9c26SBarry Smith {
579f6c57405SHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
580f6c57405SHong Zhang   PetscErrorCode ierr;
581f6c57405SHong Zhang 
582f6c57405SHong Zhang   PetscFunctionBegin;
583f6c57405SHong Zhang   /* check if matrix is mumps type */
584f6c57405SHong Zhang   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
585f6c57405SHong Zhang 
586f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
587f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
588f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
589f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
590f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
591f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
592f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
593f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
594f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
595f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
596f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
597f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
598f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
599f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
600f6c57405SHong Zhang   if (!lu->myid && lu->id.ICNTL(11)>0) {
601f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
602f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
603f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
604f6c57405SHong 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);
605f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
606f6c57405SHong 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);
607f6c57405SHong Zhang 
608f6c57405SHong Zhang   }
609f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
610f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
611f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
612f6c57405SHong Zhang   /* ICNTL(15-17) not used */
613f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
614f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
615f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
616f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
617c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",lu->id.ICNTL(22));CHKERRQ(ierr);
618c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",lu->id.ICNTL(23));CHKERRQ(ierr);
619c0165424SHong Zhang 
620c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",lu->id.ICNTL(24));CHKERRQ(ierr);
621c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",lu->id.ICNTL(25));CHKERRQ(ierr);
622c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",lu->id.ICNTL(26));CHKERRQ(ierr);
623c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",lu->id.ICNTL(27));CHKERRQ(ierr);
624f6c57405SHong Zhang 
625f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
626f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
627f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
628f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
629c0165424SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",lu->id.CNTL(5));CHKERRQ(ierr);
630f6c57405SHong Zhang 
631f6c57405SHong Zhang   /* infomation local to each processor */
632f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
6337adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
6347adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
635f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
6367adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
6377adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
638f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
6397adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
6407adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
641f6c57405SHong Zhang 
642f6c57405SHong 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);}
6437adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
6447adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
645f6c57405SHong Zhang 
646f6c57405SHong 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);}
6477adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
6487adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
649f6c57405SHong Zhang 
650f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
6517adad957SLisandro Dalcin   ierr = PetscSynchronizedPrintf(((PetscObject)A)->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
6527adad957SLisandro Dalcin   ierr = PetscSynchronizedFlush(((PetscObject)A)->comm);
653f6c57405SHong Zhang 
654f6c57405SHong Zhang   if (!lu->myid){ /* information from the host */
655f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
656f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
657f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
658f6c57405SHong Zhang 
659f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
660f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
661f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
662f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
663f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
664f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
665f6c57405SHong 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);
666f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
667f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
668f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
669f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
670f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
671f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
672f6c57405SHong 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);
673f6c57405SHong 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);
674f6c57405SHong 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);
675f6c57405SHong 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);
676f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
677f6c57405SHong 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);
678f6c57405SHong 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);
679f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
680f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
681f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
682f6c57405SHong Zhang   }
683f6c57405SHong Zhang   PetscFunctionReturn(0);
684f6c57405SHong Zhang }
685f6c57405SHong Zhang 
686f6c57405SHong Zhang #undef __FUNCT__
687f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS"
688b24902e0SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
689b24902e0SBarry Smith {
690f6c57405SHong Zhang   PetscErrorCode    ierr;
691f6c57405SHong Zhang   PetscTruth        iascii;
692f6c57405SHong Zhang   PetscViewerFormat format;
693f6c57405SHong Zhang 
694f6c57405SHong Zhang   PetscFunctionBegin;
695f6c57405SHong Zhang     ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
696f6c57405SHong Zhang   if (iascii) {
697f6c57405SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
698f6c57405SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO){
699f6c57405SHong Zhang       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
700f6c57405SHong Zhang     }
701f6c57405SHong Zhang   }
702f6c57405SHong Zhang   PetscFunctionReturn(0);
703f6c57405SHong Zhang }
704f6c57405SHong Zhang 
70535bd34faSBarry Smith #undef __FUNCT__
70635bd34faSBarry Smith #define __FUNCT__ "MatGetInfo_MUMPS"
70735bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
70835bd34faSBarry Smith {
70935bd34faSBarry Smith     Mat_MUMPS  *lu =(Mat_MUMPS*)A->spptr;
71035bd34faSBarry Smith 
71135bd34faSBarry Smith   PetscFunctionBegin;
71235bd34faSBarry Smith   info->block_size        = 1.0;
71335bd34faSBarry Smith   info->nz_allocated      = lu->id.INFOG(20);
71435bd34faSBarry Smith   info->nz_used           = lu->id.INFOG(20);
71535bd34faSBarry Smith   info->nz_unneeded       = 0.0;
71635bd34faSBarry Smith   info->assemblies        = 0.0;
71735bd34faSBarry Smith   info->mallocs           = 0.0;
71835bd34faSBarry Smith   info->memory            = 0.0;
71935bd34faSBarry Smith   info->fill_ratio_given  = 0;
72035bd34faSBarry Smith   info->fill_ratio_needed = 0;
72135bd34faSBarry Smith   info->factor_mallocs    = 0;
72235bd34faSBarry Smith   PetscFunctionReturn(0);
72335bd34faSBarry Smith }
72435bd34faSBarry Smith 
72524b6179bSKris Buschelman /*MC
72641c8de11SBarry Smith   MAT_SOLVER_MUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
72724b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
72824b6179bSKris Buschelman 
72941c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
73024b6179bSKris Buschelman 
73124b6179bSKris Buschelman   Options Database Keys:
73241c8de11SBarry Smith + -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
73324b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
73424b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
73524b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
73624b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
73724b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
73894b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
73924b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
74024b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
74124b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
74224b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
74324b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
74424b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
74524b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
74624b6179bSKris Buschelman 
74724b6179bSKris Buschelman   Level: beginner
74824b6179bSKris Buschelman 
74941c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
75041c8de11SBarry Smith 
75124b6179bSKris Buschelman M*/
75224b6179bSKris Buschelman 
7532877fffaSHong Zhang EXTERN_C_BEGIN
75435bd34faSBarry Smith #undef __FUNCT__
75535bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_mumps"
75635bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_mumps(Mat A,const MatSolverPackage *type)
75735bd34faSBarry Smith {
75835bd34faSBarry Smith   PetscFunctionBegin;
75935bd34faSBarry Smith   *type = MAT_SOLVER_MUMPS;
76035bd34faSBarry Smith   PetscFunctionReturn(0);
76135bd34faSBarry Smith }
76235bd34faSBarry Smith EXTERN_C_END
76335bd34faSBarry Smith 
76435bd34faSBarry Smith EXTERN_C_BEGIN
7652877fffaSHong Zhang /*
7662877fffaSHong Zhang     The seq and mpi versions of this function are the same
7672877fffaSHong Zhang */
7682877fffaSHong Zhang #undef __FUNCT__
7692877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_seqaij_mumps"
7702877fffaSHong Zhang PetscErrorCode MatGetFactor_seqaij_mumps(Mat A,MatFactorType ftype,Mat *F)
7712877fffaSHong Zhang {
7722877fffaSHong Zhang   Mat            B;
7732877fffaSHong Zhang   PetscErrorCode ierr;
7742877fffaSHong Zhang   Mat_MUMPS      *mumps;
7752877fffaSHong Zhang 
7762877fffaSHong Zhang   PetscFunctionBegin;
7772877fffaSHong Zhang   if (ftype != MAT_FACTOR_LU) {
7782877fffaSHong Zhang     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix");
7792877fffaSHong Zhang   }
7802877fffaSHong Zhang   /* Create the factorization matrix */
7812877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
7822877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
7832877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
7842877fffaSHong Zhang   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
7852877fffaSHong Zhang 
7862877fffaSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
7872877fffaSHong Zhang   B->ops->view             = MatView_MUMPS;
78835bd34faSBarry Smith   B->ops->getinfo          = MatGetInfo_MUMPS;
78935bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
79097969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
7912877fffaSHong Zhang   B->factor                = MAT_FACTOR_LU;
7922877fffaSHong Zhang 
7932877fffaSHong Zhang   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
7942877fffaSHong Zhang   mumps->CleanUpMUMPS              = PETSC_FALSE;
7952877fffaSHong Zhang   mumps->isAIJ                     = PETSC_TRUE;
7962877fffaSHong Zhang   mumps->scat_rhs                  = PETSC_NULL;
7972877fffaSHong Zhang   mumps->scat_sol                  = PETSC_NULL;
7982877fffaSHong Zhang   mumps->nSolve                    = 0;
7992877fffaSHong Zhang   mumps->MatDestroy                = B->ops->destroy;
8002877fffaSHong Zhang   B->ops->destroy                  = MatDestroy_MUMPS;
8012877fffaSHong Zhang   B->spptr                         = (void*)mumps;
8022877fffaSHong Zhang 
8032877fffaSHong Zhang   *F = B;
8042877fffaSHong Zhang   PetscFunctionReturn(0);
8052877fffaSHong Zhang }
8062877fffaSHong Zhang EXTERN_C_END
8072877fffaSHong Zhang 
8082877fffaSHong Zhang EXTERN_C_BEGIN
8092877fffaSHong Zhang #undef __FUNCT__
8102877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_mpiaij_mumps"
8112877fffaSHong Zhang PetscErrorCode MatGetFactor_mpiaij_mumps(Mat A,MatFactorType ftype,Mat *F)
8122877fffaSHong Zhang {
8132877fffaSHong Zhang   Mat            B;
8142877fffaSHong Zhang   PetscErrorCode ierr;
8152877fffaSHong Zhang   Mat_MUMPS      *mumps;
8162877fffaSHong Zhang 
8172877fffaSHong Zhang   PetscFunctionBegin;
8182877fffaSHong Zhang   if (ftype != MAT_FACTOR_LU) {
8192877fffaSHong Zhang     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with MUMPS Cholesky, use SBAIJ matrix");
8202877fffaSHong Zhang   }
8212877fffaSHong Zhang   /* Create the factorization matrix */
8222877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
8232877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
8242877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
8252877fffaSHong Zhang   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
8262877fffaSHong Zhang   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
8272877fffaSHong Zhang 
8282877fffaSHong Zhang   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
8292877fffaSHong Zhang   B->ops->view             = MatView_MUMPS;
83035bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
83197969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
8322877fffaSHong Zhang   B->factor                = MAT_FACTOR_LU;
8332877fffaSHong Zhang 
8342877fffaSHong Zhang   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
8352877fffaSHong Zhang   mumps->CleanUpMUMPS              = PETSC_FALSE;
8362877fffaSHong Zhang   mumps->isAIJ                     = PETSC_TRUE;
8372877fffaSHong Zhang   mumps->scat_rhs                  = PETSC_NULL;
8382877fffaSHong Zhang   mumps->scat_sol                  = PETSC_NULL;
8392877fffaSHong Zhang   mumps->nSolve                    = 0;
840f3c0ef26SHong Zhang   mumps->MatDestroy                = B->ops->destroy;
841f3c0ef26SHong Zhang   B->ops->destroy                  = MatDestroy_MUMPS;
8422877fffaSHong Zhang   B->spptr                         = (void*)mumps;
8432877fffaSHong Zhang 
8442877fffaSHong Zhang   *F = B;
8452877fffaSHong Zhang   PetscFunctionReturn(0);
8462877fffaSHong Zhang }
8472877fffaSHong Zhang EXTERN_C_END
8482877fffaSHong Zhang 
8492877fffaSHong Zhang EXTERN_C_BEGIN
8502877fffaSHong Zhang #undef __FUNCT__
8512877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_seqsbaij_mumps"
8522877fffaSHong Zhang PetscErrorCode MatGetFactor_seqsbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
8532877fffaSHong Zhang {
8542877fffaSHong Zhang   Mat            B;
8552877fffaSHong Zhang   PetscErrorCode ierr;
8562877fffaSHong Zhang   Mat_MUMPS      *mumps;
8572877fffaSHong Zhang 
8582877fffaSHong Zhang   PetscFunctionBegin;
8592877fffaSHong Zhang   if (ftype != MAT_FACTOR_CHOLESKY) {
8602877fffaSHong Zhang     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
8612877fffaSHong Zhang   }
8622877fffaSHong Zhang   /* Create the factorization matrix */
8632877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
8642877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
8652877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
8662877fffaSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
8672877fffaSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
8682877fffaSHong Zhang 
8692877fffaSHong Zhang   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
8702877fffaSHong Zhang   B->ops->view                   = MatView_MUMPS;
87135bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
87297969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
8732877fffaSHong Zhang   B->factor                      = MAT_FACTOR_CHOLESKY;
8742877fffaSHong Zhang 
8752877fffaSHong Zhang   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
8762877fffaSHong Zhang   mumps->CleanUpMUMPS              = PETSC_FALSE;
8772877fffaSHong Zhang   mumps->isAIJ                     = PETSC_TRUE;
8782877fffaSHong Zhang   mumps->scat_rhs                  = PETSC_NULL;
8792877fffaSHong Zhang   mumps->scat_sol                  = PETSC_NULL;
8802877fffaSHong Zhang   mumps->nSolve                    = 0;
8812877fffaSHong Zhang   mumps->MatDestroy                = B->ops->destroy;
8822877fffaSHong Zhang   B->ops->destroy                  = MatDestroy_MUMPS;
8832877fffaSHong Zhang   B->spptr                         = (void*)mumps;
884f3c0ef26SHong Zhang 
8852877fffaSHong Zhang   *F = B;
8862877fffaSHong Zhang   PetscFunctionReturn(0);
8872877fffaSHong Zhang }
8882877fffaSHong Zhang EXTERN_C_END
8892877fffaSHong Zhang 
8902877fffaSHong Zhang EXTERN_C_BEGIN
8912877fffaSHong Zhang #undef __FUNCT__
8922877fffaSHong Zhang #define __FUNCT__ "MatGetFactor_mpisbaij_mumps"
8932877fffaSHong Zhang PetscErrorCode MatGetFactor_mpisbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
8942877fffaSHong Zhang {
8952877fffaSHong Zhang   Mat            B;
8962877fffaSHong Zhang   PetscErrorCode ierr;
8972877fffaSHong Zhang   Mat_MUMPS      *mumps;
8982877fffaSHong Zhang 
8992877fffaSHong Zhang   PetscFunctionBegin;
9002877fffaSHong Zhang   if (ftype != MAT_FACTOR_CHOLESKY) {
9012877fffaSHong Zhang     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
9022877fffaSHong Zhang   }
9032877fffaSHong Zhang   /* Create the factorization matrix */
9042877fffaSHong Zhang   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
9052877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
9062877fffaSHong Zhang   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
9072877fffaSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
9082877fffaSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
9092877fffaSHong Zhang 
9102877fffaSHong Zhang   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJMUMPS;
9112877fffaSHong Zhang   B->ops->view                   = MatView_MUMPS;
91235bd34faSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_mumps",MatFactorGetSolverPackage_mumps);CHKERRQ(ierr);
91397969023SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMumpsSetIcntl_C","MatMumpsSetIcntl",MatMumpsSetIcntl);CHKERRQ(ierr);
9142877fffaSHong Zhang   B->factor                      = MAT_FACTOR_CHOLESKY;
9152877fffaSHong Zhang 
9162877fffaSHong Zhang   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
9172877fffaSHong Zhang   mumps->CleanUpMUMPS              = PETSC_FALSE;
9182877fffaSHong Zhang   mumps->isAIJ                     = PETSC_TRUE;
9192877fffaSHong Zhang   mumps->scat_rhs                  = PETSC_NULL;
9202877fffaSHong Zhang   mumps->scat_sol                  = PETSC_NULL;
9212877fffaSHong Zhang   mumps->nSolve                    = 0;
922f3c0ef26SHong Zhang   mumps->MatDestroy                = B->ops->destroy;
923f3c0ef26SHong Zhang   B->ops->destroy                  = MatDestroy_MUMPS;
9242877fffaSHong Zhang   B->spptr                         = (void*)mumps;
925f3c0ef26SHong Zhang 
9262877fffaSHong Zhang   *F = B;
9272877fffaSHong Zhang   PetscFunctionReturn(0);
9282877fffaSHong Zhang }
9292877fffaSHong Zhang EXTERN_C_END
93097969023SHong Zhang 
93161288eaaSHong Zhang /* -------------------------------------------------------------------------------------------*/
93261288eaaSHong Zhang /*@
93361288eaaSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
93461288eaaSHong Zhang 
93561288eaaSHong Zhang    Collective on Mat
93661288eaaSHong Zhang 
93761288eaaSHong Zhang    Input Parameters:
93861288eaaSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
93961288eaaSHong Zhang .  idx - index of MUMPS parameter array ICNTL()
94061288eaaSHong Zhang -  icntl - value of MUMPS ICNTL(imumps)
94161288eaaSHong Zhang 
94261288eaaSHong Zhang   Options Database:
94361288eaaSHong Zhang .   -mat_mumps_icntl_<idx> <icntl>
94461288eaaSHong Zhang 
94561288eaaSHong Zhang    Level: beginner
94661288eaaSHong Zhang 
94761288eaaSHong Zhang    References: MUMPS Users' Guide
94861288eaaSHong Zhang 
94961288eaaSHong Zhang .seealso: MatGetFactor()
95061288eaaSHong Zhang @*/
95197969023SHong Zhang #undef __FUNCT__
95297969023SHong Zhang #define __FUNCT__ "MatMumpsSetIcntl"
95361288eaaSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMumpsSetIcntl(Mat F,PetscInt idx,PetscInt icntl)
95497969023SHong Zhang {
95597969023SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)(F)->spptr;
95697969023SHong Zhang 
95797969023SHong Zhang   PetscFunctionBegin;
95861288eaaSHong Zhang   lu->id.ICNTL(idx) = icntl;
95997969023SHong Zhang   PetscFunctionReturn(0);
96097969023SHong Zhang }
96197969023SHong Zhang 
962