xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision f6c5740593df9f17f0e8e1ad8f01ca5cd635062c)
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;
426849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
436849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
446849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
456849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
466849ba73SBarry Smith   PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
476849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
486849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
496849ba73SBarry Smith   PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
50f0c56d0fSKris Buschelman } Mat_MUMPS;
51f0c56d0fSKris Buschelman 
52dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
53892f6c3fSKris Buschelman EXTERN_C_BEGIN
5475179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*);
55892f6c3fSKris Buschelman EXTERN_C_END
56397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
57397b6df1SKris Buschelman /*
58397b6df1SKris Buschelman   input:
5975747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
60397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
61397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
62397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
63397b6df1SKris Buschelman   output:
64397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
65397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
66397b6df1SKris Buschelman  */
67dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
68c1490034SHong Zhang   PetscInt       *ai, *aj, *bi, *bj, rstart,nz, *garray;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
702a4c71feSBarry Smith   PetscInt       i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol;
71c1490034SHong Zhang   PetscInt       *row,*col;
72397b6df1SKris Buschelman   PetscScalar    *av, *bv,*val;
73f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
74397b6df1SKris Buschelman 
75397b6df1SKris Buschelman   PetscFunctionBegin;
76397b6df1SKris Buschelman   if (mumps->isAIJ){
77397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
78397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
79397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
80397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
812a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
82397b6df1SKris Buschelman     garray = mat->garray;
83397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
84397b6df1SKris Buschelman 
85397b6df1SKris Buschelman   } else {
86397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
87397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
88397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
890c0e133fSBarry Smith     if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs);
906c6c5352SBarry Smith     nz = aa->nz + bb->nz;
912a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
92397b6df1SKris Buschelman     garray = mat->garray;
93397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
94397b6df1SKris Buschelman   }
95397b6df1SKris Buschelman 
96397b6df1SKris Buschelman   if (!valOnly){
977c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr);
987c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr);
99397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
100397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
101397b6df1SKris Buschelman   } else {
102397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
103397b6df1SKris Buschelman   }
104397b6df1SKris Buschelman   *nnz = nz;
105397b6df1SKris Buschelman 
106028e57e8SHong Zhang   jj = 0; irow = rstart;
107397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
108397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
109397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
110397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
111397b6df1SKris Buschelman     bjj = bj + bi[i];
112397b6df1SKris Buschelman 
113397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
114397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11575747be1SHong Zhang     j=-1;
11675747be1SHong Zhang     do {
11775747be1SHong Zhang       j++;
11875747be1SHong Zhang       if (j == countB) break;
119397b6df1SKris Buschelman       jcol = garray[bjj[j]];
12075747be1SHong Zhang     } while (jcol < colA_start);
12175747be1SHong Zhang     jB = j;
122397b6df1SKris Buschelman 
123397b6df1SKris Buschelman     /* B-part, smaller col index */
124397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
125397b6df1SKris Buschelman     for (j=0; j<jB; j++){
126397b6df1SKris Buschelman       jcol = garray[bjj[j]];
127397b6df1SKris Buschelman       if (!valOnly){
128397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12975747be1SHong Zhang 
130397b6df1SKris Buschelman       }
131397b6df1SKris Buschelman       val[jj++] = *bv++;
132397b6df1SKris Buschelman     }
133397b6df1SKris Buschelman     /* A-part */
134397b6df1SKris Buschelman     for (j=0; j<countA; j++){
135397b6df1SKris Buschelman       if (!valOnly){
136397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
137397b6df1SKris Buschelman       }
138397b6df1SKris Buschelman       val[jj++] = *av++;
139397b6df1SKris Buschelman     }
140397b6df1SKris Buschelman     /* B-part, larger col index */
141397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
142397b6df1SKris Buschelman       if (!valOnly){
143397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
144397b6df1SKris Buschelman       }
145397b6df1SKris Buschelman       val[jj++] = *bv++;
146397b6df1SKris Buschelman     }
147397b6df1SKris Buschelman     irow++;
148397b6df1SKris Buschelman   }
149397b6df1SKris Buschelman 
150397b6df1SKris Buschelman   PetscFunctionReturn(0);
151397b6df1SKris Buschelman }
152397b6df1SKris Buschelman 
153c338a77dSKris Buschelman EXTERN_C_BEGIN
154c338a77dSKris Buschelman #undef __FUNCT__
155c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
15617667f90SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat)
157521d7252SBarry Smith {
158dfbe8321SBarry Smith   PetscErrorCode ierr;
159c338a77dSKris Buschelman   Mat            B=*newmat;
160f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
161901853e0SKris Buschelman   void           (*f)(void);
162c338a77dSKris Buschelman 
163c338a77dSKris Buschelman   PetscFunctionBegin;
164ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
165c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
166c338a77dSKris Buschelman   }
167f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
168f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
169f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
170f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
171f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
172f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
173901853e0SKris Buschelman 
17417667f90SBarry Smith   /* put back original composed preallocation function */
175f68b968cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
176901853e0SKris Buschelman   if (f) {
177f68b968cSBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr);
178901853e0SKris Buschelman   }
179901853e0SKris Buschelman 
180901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
181901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
182901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
183901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
1842895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
1852895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
186901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
187901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
188901853e0SKris Buschelman 
189901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
190c338a77dSKris Buschelman   *newmat = B;
191c338a77dSKris Buschelman   PetscFunctionReturn(0);
192c338a77dSKris Buschelman }
193c338a77dSKris Buschelman EXTERN_C_END
194c338a77dSKris Buschelman 
195397b6df1SKris Buschelman #undef __FUNCT__
1963924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
197dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
198dfbe8321SBarry Smith {
199f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
200dfbe8321SBarry Smith   PetscErrorCode ierr;
201c1490034SHong Zhang   PetscMPIInt    size=lu->size;
2026849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
203397b6df1SKris Buschelman   PetscFunctionBegin;
204397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
205397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
206329ec9b3SHong Zhang     if (size > 1){
207329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
208329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
209329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
210329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
211329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
212329ec9b3SHong Zhang       ierr = PetscFree(lu->val);CHKERRQ(ierr);
213329ec9b3SHong Zhang     }
214397b6df1SKris Buschelman     lu->id.job=JOB_END;
215397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
216397b6df1SKris Buschelman     zmumps_c(&lu->id);
217397b6df1SKris Buschelman #else
218397b6df1SKris Buschelman     dmumps_c(&lu->id);
219397b6df1SKris Buschelman #endif
220c338a77dSKris Buschelman     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
221c338a77dSKris Buschelman     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
222397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
223397b6df1SKris Buschelman   }
224a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
225a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
226c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
227397b6df1SKris Buschelman   PetscFunctionReturn(0);
228397b6df1SKris Buschelman }
229397b6df1SKris Buschelman 
230397b6df1SKris Buschelman #undef __FUNCT__
231a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
232dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
233dfbe8321SBarry Smith {
2346849ba73SBarry Smith   PetscErrorCode ierr;
235329ec9b3SHong Zhang   PetscMPIInt    size;
236a39386dcSKris Buschelman 
237a39386dcSKris Buschelman   PetscFunctionBegin;
238a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
239a39386dcSKris Buschelman   if (size==1) {
240ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
241a39386dcSKris Buschelman   } else {
242ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
243a39386dcSKris Buschelman   }
244a39386dcSKris Buschelman   PetscFunctionReturn(0);
245a39386dcSKris Buschelman }
246a39386dcSKris Buschelman 
247a39386dcSKris Buschelman #undef __FUNCT__
248a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
249dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
250dfbe8321SBarry Smith {
2516849ba73SBarry Smith   PetscErrorCode ierr;
252329ec9b3SHong Zhang   PetscMPIInt    size;
253a39386dcSKris Buschelman 
254a39386dcSKris Buschelman   PetscFunctionBegin;
255a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
256a39386dcSKris Buschelman   if (size==1) {
257ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
258a39386dcSKris Buschelman   } else {
259ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
260a39386dcSKris Buschelman   }
261a39386dcSKris Buschelman   PetscFunctionReturn(0);
262a39386dcSKris Buschelman }
263a39386dcSKris Buschelman 
264a39386dcSKris Buschelman #undef __FUNCT__
265*f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
266*f6c57405SHong Zhang PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) {
267f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
268d54de34fSKris Buschelman   PetscScalar    *array;
269397b6df1SKris Buschelman   Vec            x_seq;
270329ec9b3SHong Zhang   IS             is_iden,is_petsc;
271329ec9b3SHong Zhang   VecScatter     scat_rhs=lu->scat_rhs,scat_sol=lu->scat_sol;
272dfbe8321SBarry Smith   PetscErrorCode ierr;
273329ec9b3SHong Zhang   PetscInt       i;
274397b6df1SKris Buschelman 
275397b6df1SKris Buschelman   PetscFunctionBegin;
276329ec9b3SHong Zhang   lu->id.nrhs = 1;
277329ec9b3SHong Zhang   x_seq = lu->b_seq;
278397b6df1SKris Buschelman   if (lu->size > 1){
279329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
280329ec9b3SHong Zhang     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
281329ec9b3SHong Zhang     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
282397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
283397b6df1SKris Buschelman   } else {  /* size == 1 */
284397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
285397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
286397b6df1SKris Buschelman   }
287397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
288397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
289397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
290397b6df1SKris Buschelman #else
291397b6df1SKris Buschelman     lu->id.rhs = array;
292397b6df1SKris Buschelman #endif
293397b6df1SKris Buschelman   }
294329ec9b3SHong Zhang   if (lu->size == 1){
295329ec9b3SHong Zhang     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
296329ec9b3SHong Zhang   } else if (!lu->myid){
297329ec9b3SHong Zhang     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
298329ec9b3SHong Zhang   }
299329ec9b3SHong Zhang 
300329ec9b3SHong Zhang   if (lu->size > 1){
301329ec9b3SHong Zhang     /* distributed solution */
302329ec9b3SHong Zhang     lu->id.ICNTL(21) = 1;
303329ec9b3SHong Zhang     if (!lu->nSolve){
304329ec9b3SHong Zhang       /* Create x_seq=sol_loc for repeated use */
305329ec9b3SHong Zhang       PetscInt    lsol_loc;
306329ec9b3SHong Zhang       PetscScalar *sol_loc;
307329ec9b3SHong Zhang       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
308329ec9b3SHong Zhang       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
309329ec9b3SHong Zhang       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
310329ec9b3SHong Zhang       lu->id.lsol_loc = lsol_loc;
311ea5584a4SSatish Balay       lu->id.sol_loc  = (F_DOUBLE *)sol_loc;
312329ec9b3SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
313329ec9b3SHong Zhang     }
314329ec9b3SHong Zhang   }
315397b6df1SKris Buschelman 
316397b6df1SKris Buschelman   /* solve phase */
317329ec9b3SHong Zhang   /*-------------*/
318397b6df1SKris Buschelman   lu->id.job = 3;
319397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
320397b6df1SKris Buschelman   zmumps_c(&lu->id);
321397b6df1SKris Buschelman #else
322397b6df1SKris Buschelman   dmumps_c(&lu->id);
323397b6df1SKris Buschelman #endif
324397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
32579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
326397b6df1SKris Buschelman   }
327397b6df1SKris Buschelman 
328329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
329329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
330329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
331329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
332329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
333397b6df1SKris Buschelman       }
334329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
335329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
336329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
337329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
338397b6df1SKris Buschelman     }
339329ec9b3SHong Zhang     ierr = VecScatterBegin(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
340329ec9b3SHong Zhang     ierr = VecScatterEnd(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
341329ec9b3SHong Zhang   }
342329ec9b3SHong Zhang   lu->nSolve++;
343397b6df1SKris Buschelman   PetscFunctionReturn(0);
344397b6df1SKris Buschelman }
345397b6df1SKris Buschelman 
346a58c3f20SHong Zhang /*
347a58c3f20SHong Zhang   input:
348a58c3f20SHong Zhang    F:        numeric factor
349a58c3f20SHong Zhang   output:
350a58c3f20SHong Zhang    nneg:     total number of negative pivots
351a58c3f20SHong Zhang    nzero:    0
352a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
353a58c3f20SHong Zhang */
354a58c3f20SHong Zhang 
355a58c3f20SHong Zhang #undef __FUNCT__
356a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
357dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
358a58c3f20SHong Zhang {
359a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
360dfbe8321SBarry Smith   PetscErrorCode ierr;
361c1490034SHong Zhang   PetscMPIInt    size;
362a58c3f20SHong Zhang 
363a58c3f20SHong Zhang   PetscFunctionBegin;
364bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
365bcb30aebSHong 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 */
366bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
36779a5c55eSBarry 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));
368bcb30aebSHong Zhang   }
369a58c3f20SHong Zhang   if (nneg){
370a58c3f20SHong Zhang     if (!lu->myid){
371a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
372a58c3f20SHong Zhang     }
373bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
374a58c3f20SHong Zhang   }
375a58c3f20SHong Zhang   if (nzero) *nzero = 0;
3762a4c71feSBarry Smith   if (npos)  *npos  = F->rmap.N - (*nneg);
377a58c3f20SHong Zhang   PetscFunctionReturn(0);
378a58c3f20SHong Zhang }
379a58c3f20SHong Zhang 
380397b6df1SKris Buschelman #undef __FUNCT__
381*f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
382*f6c57405SHong Zhang PetscErrorCode MatFactorNumeric_MUMPS(Mat A,MatFactorInfo *info,Mat *F)
383af281ebdSHong Zhang {
384f0c56d0fSKris Buschelman   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
385f0c56d0fSKris Buschelman   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
3866849ba73SBarry Smith   PetscErrorCode ierr;
3872a4c71feSBarry Smith   PetscInt       rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
388397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
389e09efc27SHong Zhang   Mat            F_diag;
390397b6df1SKris Buschelman 
391397b6df1SKris Buschelman   PetscFunctionBegin;
392397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
393*f6c57405SHong Zhang     (*F)->ops->solve    = MatSolve_MUMPS;
394397b6df1SKris Buschelman 
395397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
396397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
397397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
39875747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
399397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
400397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
401a0e2756fSSatish Balay     ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr);
402397b6df1SKris Buschelman 
403397b6df1SKris Buschelman     /* Set mumps options */
404397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
405397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
406397b6df1SKris Buschelman     lu->id.sym=lu->sym;
407397b6df1SKris Buschelman     if (lu->sym == 2){
408397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
409397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
410397b6df1SKris Buschelman     }
411397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
412397b6df1SKris Buschelman     zmumps_c(&lu->id);
413397b6df1SKris Buschelman #else
414397b6df1SKris Buschelman     dmumps_c(&lu->id);
415397b6df1SKris Buschelman #endif
416397b6df1SKris Buschelman 
417397b6df1SKris Buschelman     if (lu->size == 1){
418397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
419397b6df1SKris Buschelman     } else {
420397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
421397b6df1SKris Buschelman     }
422397b6df1SKris Buschelman 
423397b6df1SKris Buschelman     icntl=-1;
42421f4b680SHong Zhang     lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
425397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
42619facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
427397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
428397b6df1SKris Buschelman     } else { /* no output */
429397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
430397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
431397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
432397b6df1SKris Buschelman     }
433397b6df1SKris 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);
434397b6df1SKris Buschelman     icntl=-1;
435397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
436397b6df1SKris Buschelman     if (flg) {
437397b6df1SKris Buschelman       if (icntl== 1){
438397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
439397b6df1SKris Buschelman       } else {
440397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
441397b6df1SKris Buschelman       }
442397b6df1SKris Buschelman     }
443397b6df1SKris 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);
444397b6df1SKris 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);
44594b7f48cSBarry 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);
446397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
447397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
448adc1d99fSHong 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);
449397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
450397b6df1SKris Buschelman 
451397b6df1SKris 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);
452397b6df1SKris 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);
453397b6df1SKris 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);
45425f9c88cSHong 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);
455397b6df1SKris Buschelman     PetscOptionsEnd();
456397b6df1SKris Buschelman   }
457397b6df1SKris Buschelman 
458397b6df1SKris Buschelman   /* define matrix A */
459397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
460397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
461397b6df1SKris Buschelman     if (!lu->myid) {
462c36ead0aSKris Buschelman       if (lua->isAIJ){
463397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
464397b6df1SKris Buschelman         nz               = aa->nz;
465397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
466397b6df1SKris Buschelman       } else {
467397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
4686c6c5352SBarry Smith         nz                  =  aa->nz;
469397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
470397b6df1SKris Buschelman       }
471397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
4727c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
4737c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
474397b6df1SKris Buschelman         nz = 0;
475397b6df1SKris Buschelman         for (i=0; i<M; i++){
476397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
477397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
478397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
479397b6df1SKris Buschelman           }
480397b6df1SKris Buschelman         }
481397b6df1SKris Buschelman       }
482397b6df1SKris Buschelman     }
483397b6df1SKris Buschelman     break;
484397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
485397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
486397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
487397b6df1SKris Buschelman     } else {
488397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
489397b6df1SKris Buschelman     }
490397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
491397b6df1SKris Buschelman     break;
492397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
493397b6df1SKris Buschelman   }
494397b6df1SKris Buschelman 
495397b6df1SKris Buschelman   /* analysis phase */
496329ec9b3SHong Zhang   /*----------------*/
497397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
498329ec9b3SHong Zhang     lu->id.job = 1;
499329ec9b3SHong Zhang 
500397b6df1SKris Buschelman     lu->id.n = M;
501397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
502397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
503397b6df1SKris Buschelman       if (!lu->myid) {
504397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
505397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
506397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
507397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
508397b6df1SKris Buschelman #else
509397b6df1SKris Buschelman           lu->id.a = lu->val;
510397b6df1SKris Buschelman #endif
511397b6df1SKris Buschelman         }
512397b6df1SKris Buschelman       }
513397b6df1SKris Buschelman       break;
514397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
515397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
516397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
517397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
518397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
519397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
520397b6df1SKris Buschelman #else
521397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
522397b6df1SKris Buschelman #endif
523397b6df1SKris Buschelman       }
524329ec9b3SHong Zhang       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
525329ec9b3SHong Zhang       IS  is_iden;
526329ec9b3SHong Zhang       Vec b;
527329ec9b3SHong Zhang       if (!lu->myid){
528329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr);
529329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr);
530329ec9b3SHong Zhang       } else {
531329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
532329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
533329ec9b3SHong Zhang       }
534329ec9b3SHong Zhang       ierr = VecCreate(A->comm,&b);CHKERRQ(ierr);
535329ec9b3SHong Zhang       ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
536329ec9b3SHong Zhang       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
537329ec9b3SHong Zhang 
538329ec9b3SHong Zhang       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
539329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
540329ec9b3SHong Zhang       ierr = VecDestroy(b);CHKERRQ(ierr);
541397b6df1SKris Buschelman       break;
542397b6df1SKris Buschelman     }
543397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
544397b6df1SKris Buschelman     zmumps_c(&lu->id);
545397b6df1SKris Buschelman #else
546397b6df1SKris Buschelman     dmumps_c(&lu->id);
547397b6df1SKris Buschelman #endif
548397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
54979a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
550397b6df1SKris Buschelman     }
551397b6df1SKris Buschelman   }
552397b6df1SKris Buschelman 
553397b6df1SKris Buschelman   /* numerical factorization phase */
554329ec9b3SHong Zhang   /*-------------------------------*/
555329ec9b3SHong Zhang   lu->id.job = 2;
556958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
557a7aca84bSHong Zhang     if (!lu->myid) {
558397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
559397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
560397b6df1SKris Buschelman #else
561397b6df1SKris Buschelman       lu->id.a = lu->val;
562397b6df1SKris Buschelman #endif
563397b6df1SKris Buschelman     }
564397b6df1SKris Buschelman   } else {
565397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
566397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
567397b6df1SKris Buschelman #else
568397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
569397b6df1SKris Buschelman #endif
570397b6df1SKris Buschelman   }
571397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
572397b6df1SKris Buschelman   zmumps_c(&lu->id);
573397b6df1SKris Buschelman #else
574397b6df1SKris Buschelman   dmumps_c(&lu->id);
575397b6df1SKris Buschelman #endif
576397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
57719facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
57819facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
57919facb7aSBarry Smith     } else {
58079a5c55eSBarry 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));
581397b6df1SKris Buschelman     }
58219facb7aSBarry Smith   }
583397b6df1SKris Buschelman 
58419facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
58579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
586397b6df1SKris Buschelman   }
587397b6df1SKris Buschelman 
5888ada1bb4SHong Zhang   if (lu->size > 1){
589e09efc27SHong Zhang     if ((*F)->factor == FACTOR_LU){
590e09efc27SHong Zhang       F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
591e09efc27SHong Zhang     } else {
592e09efc27SHong Zhang       F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
593e09efc27SHong Zhang     }
594e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
595329ec9b3SHong Zhang     if (lu->nSolve){
596329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
597329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
598329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
599329ec9b3SHong Zhang     }
6008ada1bb4SHong Zhang   }
601397b6df1SKris Buschelman   (*F)->assembled   = PETSC_TRUE;
602397b6df1SKris Buschelman   lu->matstruc      = SAME_NONZERO_PATTERN;
603ace87b0dSHong Zhang   lu->CleanUpMUMPS  = PETSC_TRUE;
604329ec9b3SHong Zhang   lu->nSolve        = 0;
605397b6df1SKris Buschelman   PetscFunctionReturn(0);
606397b6df1SKris Buschelman }
607397b6df1SKris Buschelman 
608397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
609397b6df1SKris Buschelman #undef __FUNCT__
610f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
611dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
612397b6df1SKris Buschelman   Mat            B;
613f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
614dfbe8321SBarry Smith   PetscErrorCode ierr;
615397b6df1SKris Buschelman 
616397b6df1SKris Buschelman   PetscFunctionBegin;
617397b6df1SKris Buschelman   /* Create the factorization matrix */
618f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
6192a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
620be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
621397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
622397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
623397b6df1SKris Buschelman 
624*f6c57405SHong Zhang   B->ops->lufactornumeric = MatFactorNumeric_MUMPS;
625397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
626f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
627397b6df1SKris Buschelman   lu->sym                 = 0;
628397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
629397b6df1SKris Buschelman 
630397b6df1SKris Buschelman   *F = B;
631397b6df1SKris Buschelman   PetscFunctionReturn(0);
632397b6df1SKris Buschelman }
633397b6df1SKris Buschelman 
634397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
635397b6df1SKris Buschelman #undef __FUNCT__
636f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
637dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
638397b6df1SKris Buschelman   Mat            B;
639f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
640dfbe8321SBarry Smith   PetscErrorCode ierr;
641397b6df1SKris Buschelman 
642397b6df1SKris Buschelman   PetscFunctionBegin;
643397b6df1SKris Buschelman   /* Create the factorization matrix */
644f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
6452a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
646be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
647efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
648efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
649397b6df1SKris Buschelman 
650*f6c57405SHong Zhang   B->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
651a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
652397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
653f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
654397b6df1SKris Buschelman   lu->sym                       = 2;
655397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
656397b6df1SKris Buschelman 
657397b6df1SKris Buschelman   *F = B;
658397b6df1SKris Buschelman   PetscFunctionReturn(0);
659397b6df1SKris Buschelman }
660397b6df1SKris Buschelman 
661397b6df1SKris Buschelman #undef __FUNCT__
662*f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS"
663*f6c57405SHong Zhang PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
664*f6c57405SHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
665*f6c57405SHong Zhang   PetscErrorCode ierr;
666*f6c57405SHong Zhang 
667*f6c57405SHong Zhang   PetscFunctionBegin;
668*f6c57405SHong Zhang   /* check if matrix is mumps type */
669*f6c57405SHong Zhang   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
670*f6c57405SHong Zhang 
671*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
672*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
673*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
674*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
675*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
676*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
677*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
678*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
679*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
680*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
681*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
682*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
683*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
684*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
685*f6c57405SHong Zhang   if (!lu->myid && lu->id.ICNTL(11)>0) {
686*f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
687*f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
688*f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
689*f6c57405SHong 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);
690*f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
691*f6c57405SHong 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);
692*f6c57405SHong Zhang 
693*f6c57405SHong Zhang   }
694*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
695*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
696*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
697*f6c57405SHong Zhang   /* ICNTL(15-17) not used */
698*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
699*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
700*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
701*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
702*f6c57405SHong Zhang 
703*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
704*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
705*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
706*f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
707*f6c57405SHong Zhang 
708*f6c57405SHong Zhang   /* infomation local to each processor */
709*f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
710*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
711*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
712*f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
713*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
714*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
715*f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
716*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
717*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
718*f6c57405SHong Zhang   /*
719*f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);}
720*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr);
721*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
722*f6c57405SHong Zhang   */
723*f6c57405SHong Zhang 
724*f6c57405SHong 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);}
725*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
726*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
727*f6c57405SHong Zhang 
728*f6c57405SHong 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);}
729*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
730*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
731*f6c57405SHong Zhang 
732*f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
733*f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
734*f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
735*f6c57405SHong Zhang 
736*f6c57405SHong Zhang   if (!lu->myid){ /* information from the host */
737*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
738*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
739*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
740*f6c57405SHong Zhang 
741*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
742*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
743*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
744*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
745*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
746*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
747*f6c57405SHong 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);
748*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
749*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
750*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
751*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
752*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
753*f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
754*f6c57405SHong 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);
755*f6c57405SHong 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);
756*f6c57405SHong 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);
757*f6c57405SHong 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);
758*f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
759*f6c57405SHong 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);
760*f6c57405SHong 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);
761*f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
762*f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
763*f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
764*f6c57405SHong Zhang   }
765*f6c57405SHong Zhang 
766*f6c57405SHong Zhang   PetscFunctionReturn(0);
767*f6c57405SHong Zhang }
768*f6c57405SHong Zhang 
769*f6c57405SHong Zhang #undef __FUNCT__
770*f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS"
771*f6c57405SHong Zhang PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
772*f6c57405SHong Zhang   PetscErrorCode    ierr;
773*f6c57405SHong Zhang   PetscTruth        iascii;
774*f6c57405SHong Zhang   PetscViewerFormat format;
775*f6c57405SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
776*f6c57405SHong Zhang 
777*f6c57405SHong Zhang   PetscFunctionBegin;
778*f6c57405SHong Zhang   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
779*f6c57405SHong Zhang 
780*f6c57405SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
781*f6c57405SHong Zhang   if (iascii) {
782*f6c57405SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
783*f6c57405SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO){
784*f6c57405SHong Zhang       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
785*f6c57405SHong Zhang     }
786*f6c57405SHong Zhang   }
787*f6c57405SHong Zhang   PetscFunctionReturn(0);
788*f6c57405SHong Zhang }
789*f6c57405SHong Zhang 
790*f6c57405SHong Zhang #undef __FUNCT__
791f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
792dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
793dfbe8321SBarry Smith   PetscErrorCode ierr;
794f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
795c338a77dSKris Buschelman 
796397b6df1SKris Buschelman   PetscFunctionBegin;
797c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
798f0c56d0fSKris Buschelman 
799c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
800c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
801f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
802397b6df1SKris Buschelman   PetscFunctionReturn(0);
803397b6df1SKris Buschelman }
804397b6df1SKris Buschelman 
805c338a77dSKris Buschelman EXTERN_C_BEGIN
806c338a77dSKris Buschelman #undef __FUNCT__
807f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
80874a14a02SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
809dfbe8321SBarry Smith {
810dfbe8321SBarry Smith   PetscErrorCode ierr;
811521d7252SBarry Smith   PetscMPIInt    size;
812c338a77dSKris Buschelman   MPI_Comm       comm;
813c338a77dSKris Buschelman   Mat            B=*newmat;
814f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps;
815397b6df1SKris Buschelman 
816397b6df1SKris Buschelman   PetscFunctionBegin;
817c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
818f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
819c338a77dSKris Buschelman 
8205c088420SHong Zhang   if (reuse == MAT_INITIAL_MATRIX) {
8215c088420SHong Zhang     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
8225c088420SHong Zhang     /* A may have special container that is not duplicated,
8235c088420SHong Zhang        e.g., A is obtainted from MatMatMult(,&A). Save B->ops instead */
8245c088420SHong Zhang     mumps->MatDuplicate              = B->ops->duplicate;
8255c088420SHong Zhang     mumps->MatView                   = B->ops->view;
8265c088420SHong Zhang     mumps->MatAssemblyEnd            = B->ops->assemblyend;
8275c088420SHong Zhang     mumps->MatLUFactorSymbolic       = B->ops->lufactorsymbolic;
8285c088420SHong Zhang     mumps->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic;
8295c088420SHong Zhang     mumps->MatDestroy                = B->ops->destroy;
8305c088420SHong Zhang   } else {
831f0c56d0fSKris Buschelman     mumps->MatDuplicate              = A->ops->duplicate;
832c338a77dSKris Buschelman     mumps->MatView                   = A->ops->view;
833c338a77dSKris Buschelman     mumps->MatAssemblyEnd            = A->ops->assemblyend;
834c338a77dSKris Buschelman     mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
835c338a77dSKris Buschelman     mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
836c338a77dSKris Buschelman     mumps->MatDestroy                = A->ops->destroy;
8375c088420SHong Zhang   }
838a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
839c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
840f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
841c338a77dSKris Buschelman 
8424b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
843422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
844c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
845f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
846f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
8473924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
848c338a77dSKris Buschelman 
849c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
850c338a77dSKris Buschelman   if (size == 1) {
851c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
852f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
853c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
854c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
855c338a77dSKris Buschelman   } else {
856c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
857f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
858c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
859c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
860c338a77dSKris Buschelman   }
861c338a77dSKris Buschelman 
862ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr);
863c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
864c338a77dSKris Buschelman   *newmat = B;
865397b6df1SKris Buschelman   PetscFunctionReturn(0);
866397b6df1SKris Buschelman }
867c338a77dSKris Buschelman EXTERN_C_END
868397b6df1SKris Buschelman 
86924b6179bSKris Buschelman /*MC
870fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
87124b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
87224b6179bSKris Buschelman 
87324b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
87424b6179bSKris Buschelman   on how to declare the existence of external packages),
87524b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
87624b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
87724b6179bSKris Buschelman 
87824b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
87924b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
88024b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
88124b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
88228b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
88328b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
88428b08bd3SKris Buschelman   without data copy.
88524b6179bSKris Buschelman 
88624b6179bSKris Buschelman   Options Database Keys:
8870bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
88824b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
88924b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
89024b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
89124b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
89224b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
89324b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
89494b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
89524b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
89624b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
89724b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
89824b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
89924b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
90024b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
90124b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
90224b6179bSKris Buschelman 
90324b6179bSKris Buschelman   Level: beginner
90424b6179bSKris Buschelman 
90524b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
90624b6179bSKris Buschelman M*/
90724b6179bSKris Buschelman 
908397b6df1SKris Buschelman EXTERN_C_BEGIN
909397b6df1SKris Buschelman #undef __FUNCT__
910f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
911be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A)
912dfbe8321SBarry Smith {
913dfbe8321SBarry Smith   PetscErrorCode ierr;
914c1490034SHong Zhang   PetscMPIInt    size;
915397b6df1SKris Buschelman 
916397b6df1SKris Buschelman   PetscFunctionBegin;
91717667f90SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
918397b6df1SKris Buschelman   if (size == 1) {
919397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
920397b6df1SKris Buschelman   } else {
921397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
92200ff2a26SHong Zhang     /*
92300ff2a26SHong Zhang     Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
924ceb03754SKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
92500ff2a26SHong Zhang     */
926397b6df1SKris Buschelman   }
927ceb03754SKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
928397b6df1SKris Buschelman   PetscFunctionReturn(0);
929397b6df1SKris Buschelman }
930397b6df1SKris Buschelman EXTERN_C_END
931397b6df1SKris Buschelman 
932f579278aSKris Buschelman #undef __FUNCT__
933f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
934dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
935dfbe8321SBarry Smith {
936dfbe8321SBarry Smith   PetscErrorCode ierr;
937f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
938f579278aSKris Buschelman 
939f579278aSKris Buschelman   PetscFunctionBegin;
940f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
941f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
942f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
943f579278aSKris Buschelman   PetscFunctionReturn(0);
944f579278aSKris Buschelman }
945f579278aSKris Buschelman 
946f579278aSKris Buschelman EXTERN_C_BEGIN
947f579278aSKris Buschelman #undef __FUNCT__
9489c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
949c1490034SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
9509c097c71SKris Buschelman {
9519c097c71SKris Buschelman   Mat       A;
95202217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
953dfbe8321SBarry Smith   PetscErrorCode ierr;
9549c097c71SKris Buschelman 
9559c097c71SKris Buschelman   PetscFunctionBegin;
9569c097c71SKris Buschelman   /*
9579c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
9589c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
9599c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
9609c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
9619c097c71SKris Buschelman     in the preallocation routine.
9629c097c71SKris Buschelman   */
9639c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
9649c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
965ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
9669c097c71SKris Buschelman   PetscFunctionReturn(0);
9679c097c71SKris Buschelman }
9689c097c71SKris Buschelman EXTERN_C_END
9699c097c71SKris Buschelman 
9709c097c71SKris Buschelman EXTERN_C_BEGIN
9719c097c71SKris Buschelman #undef __FUNCT__
972f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
97375179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
974dfbe8321SBarry Smith {
975dfbe8321SBarry Smith   PetscErrorCode ierr;
976521d7252SBarry Smith   PetscMPIInt    size;
977f579278aSKris Buschelman   MPI_Comm       comm;
978f579278aSKris Buschelman   Mat            B=*newmat;
979422e82a1SHong Zhang   Mat_MUMPS      *mumps;
9809c097c71SKris Buschelman   void           (*f)(void);
981f579278aSKris Buschelman 
982f579278aSKris Buschelman   PetscFunctionBegin;
983ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
984f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
985f579278aSKris Buschelman   }
986f579278aSKris Buschelman 
987f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
988f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
989f579278aSKris Buschelman 
990f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
991f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
992f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
993f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
994f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
995f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
996a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
997f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
998f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
999f579278aSKris Buschelman 
1000f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
1001422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
1002c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
1003f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
1004f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
10053924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
1006f579278aSKris Buschelman 
1007f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1008f579278aSKris Buschelman   if (size == 1) {
1009f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
1010f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1011f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
1012f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1013f579278aSKris Buschelman   } else {
10149c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
1015f68b968cSBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
1016901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
10176849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
10189c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
10199c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
10209c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
10219c097c71SKris Buschelman     }
1022f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
1023f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1024f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
1025f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1026f579278aSKris Buschelman   }
1027f579278aSKris Buschelman 
1028ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr);
1029f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
1030f579278aSKris Buschelman   *newmat = B;
1031f579278aSKris Buschelman   PetscFunctionReturn(0);
1032f579278aSKris Buschelman }
1033f579278aSKris Buschelman EXTERN_C_END
1034f579278aSKris Buschelman 
1035f0c56d0fSKris Buschelman #undef __FUNCT__
1036422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
1037dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1038dfbe8321SBarry Smith   PetscErrorCode ierr;
10398e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
10408f340917SKris Buschelman 
1041f0c56d0fSKris Buschelman   PetscFunctionBegin;
10428f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
10438e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
1044f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
1045f0c56d0fSKris Buschelman }
1046f0c56d0fSKris Buschelman 
104724b6179bSKris Buschelman /*MC
1048fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
104924b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
105024b6179bSKris Buschelman 
105124b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
105224b6179bSKris Buschelman   on how to declare the existence of external packages),
105324b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
105424b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
105524b6179bSKris Buschelman 
105624b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
105724b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
105824b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
105924b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
106028b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
106128b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
106228b08bd3SKris Buschelman   without data copy.
106324b6179bSKris Buschelman 
106424b6179bSKris Buschelman   Options Database Keys:
10650bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
106624b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
106724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
106824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
106924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
107024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
107124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
107294b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
107324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
107424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
107524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
107624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
107724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
107824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
107924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
108024b6179bSKris Buschelman 
108124b6179bSKris Buschelman   Level: beginner
108224b6179bSKris Buschelman 
108324b6179bSKris Buschelman .seealso: MATAIJMUMPS
108424b6179bSKris Buschelman M*/
108524b6179bSKris Buschelman 
1086397b6df1SKris Buschelman EXTERN_C_BEGIN
1087397b6df1SKris Buschelman #undef __FUNCT__
1088f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1089be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A)
1090dfbe8321SBarry Smith {
10916849ba73SBarry Smith   PetscErrorCode ierr;
10928ada1bb4SHong Zhang   PetscMPIInt    size;
1093397b6df1SKris Buschelman 
1094397b6df1SKris Buschelman   PetscFunctionBegin;
10955441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1096397b6df1SKris Buschelman   if (size == 1) {
1097397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1098397b6df1SKris Buschelman   } else {
1099397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1100397b6df1SKris Buschelman   }
1101ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1102397b6df1SKris Buschelman   PetscFunctionReturn(0);
1103397b6df1SKris Buschelman }
1104397b6df1SKris Buschelman EXTERN_C_END
1105