xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 38f2d2fdb3b6f522a3102c6eb796cebecf3224c0)
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   }
17915d897acSHong Zhang   ierr     = PetscFree(mumps);CHKERRQ(ierr);
18015d897acSHong Zhang   A->spptr = PETSC_NULL;
181901853e0SKris Buschelman 
182901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
183901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
184901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
185901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
1862895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
1872895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
188901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
189901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
190901853e0SKris Buschelman 
191901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
192c338a77dSKris Buschelman   *newmat = B;
193c338a77dSKris Buschelman   PetscFunctionReturn(0);
194c338a77dSKris Buschelman }
195c338a77dSKris Buschelman EXTERN_C_END
196c338a77dSKris Buschelman 
197397b6df1SKris Buschelman #undef __FUNCT__
1983924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
199dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
200dfbe8321SBarry Smith {
201f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
202dfbe8321SBarry Smith   PetscErrorCode ierr;
203c1490034SHong Zhang   PetscMPIInt    size=lu->size;
2046849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
205397b6df1SKris Buschelman   PetscFunctionBegin;
206397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
207397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
208329ec9b3SHong Zhang     if (size > 1){
209329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
210329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
211329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
212329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
213329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
214329ec9b3SHong Zhang       ierr = PetscFree(lu->val);CHKERRQ(ierr);
215329ec9b3SHong Zhang     }
216397b6df1SKris Buschelman     lu->id.job=JOB_END;
217397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
218397b6df1SKris Buschelman     zmumps_c(&lu->id);
219397b6df1SKris Buschelman #else
220397b6df1SKris Buschelman     dmumps_c(&lu->id);
221397b6df1SKris Buschelman #endif
222c338a77dSKris Buschelman     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
223c338a77dSKris Buschelman     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
224397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
225397b6df1SKris Buschelman   }
226a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
227a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
228c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
229397b6df1SKris Buschelman   PetscFunctionReturn(0);
230397b6df1SKris Buschelman }
231397b6df1SKris Buschelman 
232397b6df1SKris Buschelman #undef __FUNCT__
233a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
234dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
235dfbe8321SBarry Smith {
2366849ba73SBarry Smith   PetscErrorCode ierr;
237329ec9b3SHong Zhang   PetscMPIInt    size;
238a39386dcSKris Buschelman 
239a39386dcSKris Buschelman   PetscFunctionBegin;
240a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
241a39386dcSKris Buschelman   if (size==1) {
242ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
243a39386dcSKris Buschelman   } else {
244ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
245a39386dcSKris Buschelman   }
246a39386dcSKris Buschelman   PetscFunctionReturn(0);
247a39386dcSKris Buschelman }
248a39386dcSKris Buschelman 
249a39386dcSKris Buschelman #undef __FUNCT__
250a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
251dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
252dfbe8321SBarry Smith {
2536849ba73SBarry Smith   PetscErrorCode ierr;
254329ec9b3SHong Zhang   PetscMPIInt    size;
255a39386dcSKris Buschelman 
256a39386dcSKris Buschelman   PetscFunctionBegin;
257a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
258a39386dcSKris Buschelman   if (size==1) {
259ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
260a39386dcSKris Buschelman   } else {
261ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
262a39386dcSKris Buschelman   }
263a39386dcSKris Buschelman   PetscFunctionReturn(0);
264a39386dcSKris Buschelman }
265a39386dcSKris Buschelman 
266a39386dcSKris Buschelman #undef __FUNCT__
267f6c57405SHong Zhang #define __FUNCT__ "MatSolve_MUMPS"
268f6c57405SHong Zhang PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x) {
269f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
270d54de34fSKris Buschelman   PetscScalar    *array;
271397b6df1SKris Buschelman   Vec            x_seq;
272329ec9b3SHong Zhang   IS             is_iden,is_petsc;
273ca93a389SHong Zhang   VecScatter     scat_rhs=lu->scat_rhs;
274dfbe8321SBarry Smith   PetscErrorCode ierr;
275329ec9b3SHong Zhang   PetscInt       i;
276397b6df1SKris Buschelman 
277397b6df1SKris Buschelman   PetscFunctionBegin;
278329ec9b3SHong Zhang   lu->id.nrhs = 1;
279329ec9b3SHong Zhang   x_seq = lu->b_seq;
280397b6df1SKris Buschelman   if (lu->size > 1){
281329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
282ca9f406cSSatish Balay     ierr = VecScatterBegin(scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
283ca9f406cSSatish Balay     ierr = VecScatterEnd(scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
284397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
285397b6df1SKris Buschelman   } else {  /* size == 1 */
286397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
287397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
288397b6df1SKris Buschelman   }
289397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
290397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
291397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
292397b6df1SKris Buschelman #else
293397b6df1SKris Buschelman     lu->id.rhs = array;
294397b6df1SKris Buschelman #endif
295397b6df1SKris Buschelman   }
296329ec9b3SHong Zhang   if (lu->size == 1){
297329ec9b3SHong Zhang     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
298329ec9b3SHong Zhang   } else if (!lu->myid){
299329ec9b3SHong Zhang     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
300329ec9b3SHong Zhang   }
301329ec9b3SHong Zhang 
302329ec9b3SHong Zhang   if (lu->size > 1){
303329ec9b3SHong Zhang     /* distributed solution */
304329ec9b3SHong Zhang     lu->id.ICNTL(21) = 1;
305329ec9b3SHong Zhang     if (!lu->nSolve){
306329ec9b3SHong Zhang       /* Create x_seq=sol_loc for repeated use */
307329ec9b3SHong Zhang       PetscInt    lsol_loc;
308329ec9b3SHong Zhang       PetscScalar *sol_loc;
309329ec9b3SHong Zhang       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
310329ec9b3SHong Zhang       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
311329ec9b3SHong Zhang       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
312329ec9b3SHong Zhang       lu->id.lsol_loc = lsol_loc;
3136f8312c5SHong Zhang #if defined(PETSC_USE_COMPLEX)
3146f8312c5SHong Zhang       lu->id.sol_loc  = (ZMUMPS_DOUBLE *)sol_loc;
3156f8312c5SHong Zhang #else
3166f8312c5SHong Zhang       lu->id.sol_loc  = (DMUMPS_DOUBLE *)sol_loc;
3176f8312c5SHong Zhang #endif
318329ec9b3SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
319329ec9b3SHong Zhang     }
320329ec9b3SHong Zhang   }
321397b6df1SKris Buschelman 
322397b6df1SKris Buschelman   /* solve phase */
323329ec9b3SHong Zhang   /*-------------*/
324397b6df1SKris Buschelman   lu->id.job = 3;
325397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
326397b6df1SKris Buschelman   zmumps_c(&lu->id);
327397b6df1SKris Buschelman #else
328397b6df1SKris Buschelman   dmumps_c(&lu->id);
329397b6df1SKris Buschelman #endif
330397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
33179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
332397b6df1SKris Buschelman   }
333397b6df1SKris Buschelman 
334329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
335329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
336329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
337329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
338329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
339397b6df1SKris Buschelman       }
340329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
341329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
342329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
343329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
344397b6df1SKris Buschelman     }
345ca9f406cSSatish Balay     ierr = VecScatterBegin(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
346ca9f406cSSatish Balay     ierr = VecScatterEnd(lu->scat_sol,lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
347329ec9b3SHong Zhang   }
348329ec9b3SHong Zhang   lu->nSolve++;
349397b6df1SKris Buschelman   PetscFunctionReturn(0);
350397b6df1SKris Buschelman }
351397b6df1SKris Buschelman 
352a58c3f20SHong Zhang /*
353a58c3f20SHong Zhang   input:
354a58c3f20SHong Zhang    F:        numeric factor
355a58c3f20SHong Zhang   output:
356a58c3f20SHong Zhang    nneg:     total number of negative pivots
357a58c3f20SHong Zhang    nzero:    0
358a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
359a58c3f20SHong Zhang */
360a58c3f20SHong Zhang 
361a58c3f20SHong Zhang #undef __FUNCT__
362a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
363dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
364a58c3f20SHong Zhang {
365a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
366dfbe8321SBarry Smith   PetscErrorCode ierr;
367c1490034SHong Zhang   PetscMPIInt    size;
368a58c3f20SHong Zhang 
369a58c3f20SHong Zhang   PetscFunctionBegin;
370bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
371bcb30aebSHong 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 */
372bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
37379a5c55eSBarry 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));
374bcb30aebSHong Zhang   }
375a58c3f20SHong Zhang   if (nneg){
376a58c3f20SHong Zhang     if (!lu->myid){
377a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
378a58c3f20SHong Zhang     }
379bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
380a58c3f20SHong Zhang   }
381a58c3f20SHong Zhang   if (nzero) *nzero = 0;
3822a4c71feSBarry Smith   if (npos)  *npos  = F->rmap.N - (*nneg);
383a58c3f20SHong Zhang   PetscFunctionReturn(0);
384a58c3f20SHong Zhang }
385a58c3f20SHong Zhang 
386397b6df1SKris Buschelman #undef __FUNCT__
387f6c57405SHong Zhang #define __FUNCT__ "MatFactorNumeric_MUMPS"
388f6c57405SHong Zhang PetscErrorCode MatFactorNumeric_MUMPS(Mat A,MatFactorInfo *info,Mat *F)
389af281ebdSHong Zhang {
390f0c56d0fSKris Buschelman   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
391f0c56d0fSKris Buschelman   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
3926849ba73SBarry Smith   PetscErrorCode ierr;
3932a4c71feSBarry Smith   PetscInt       rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
394397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
395e09efc27SHong Zhang   Mat            F_diag;
396397b6df1SKris Buschelman 
397397b6df1SKris Buschelman   PetscFunctionBegin;
398397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
399f6c57405SHong Zhang     (*F)->ops->solve    = MatSolve_MUMPS;
400397b6df1SKris Buschelman 
401397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
402397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
403397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
40475747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
405397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
406397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
407a0e2756fSSatish Balay     ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr);
408397b6df1SKris Buschelman 
409397b6df1SKris Buschelman     /* Set mumps options */
410397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
411397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
412397b6df1SKris Buschelman     lu->id.sym=lu->sym;
413397b6df1SKris Buschelman     if (lu->sym == 2){
414397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
415397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
416397b6df1SKris Buschelman     }
417397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
418397b6df1SKris Buschelman     zmumps_c(&lu->id);
419397b6df1SKris Buschelman #else
420397b6df1SKris Buschelman     dmumps_c(&lu->id);
421397b6df1SKris Buschelman #endif
422397b6df1SKris Buschelman 
423397b6df1SKris Buschelman     if (lu->size == 1){
424397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
425397b6df1SKris Buschelman     } else {
426397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
427397b6df1SKris Buschelman     }
428397b6df1SKris Buschelman 
429397b6df1SKris Buschelman     icntl=-1;
43021f4b680SHong Zhang     lu->id.ICNTL(4) = 0;  /* level of printing; overwrite mumps default ICNTL(4)=2 */
431397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
43219facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
433397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
434397b6df1SKris Buschelman     } else { /* no output */
435397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
436397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
437397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
438397b6df1SKris Buschelman     }
439397b6df1SKris 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);
440397b6df1SKris Buschelman     icntl=-1;
441397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
442397b6df1SKris Buschelman     if (flg) {
443397b6df1SKris Buschelman       if (icntl== 1){
444397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
445397b6df1SKris Buschelman       } else {
446397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
447397b6df1SKris Buschelman       }
448397b6df1SKris Buschelman     }
449397b6df1SKris 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);
450397b6df1SKris 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);
45194b7f48cSBarry 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);
452397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
453397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
454adc1d99fSHong 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);
455397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
456397b6df1SKris Buschelman 
457397b6df1SKris 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);
458397b6df1SKris 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);
459397b6df1SKris 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);
46025f9c88cSHong 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);
461397b6df1SKris Buschelman     PetscOptionsEnd();
462397b6df1SKris Buschelman   }
463397b6df1SKris Buschelman 
464397b6df1SKris Buschelman   /* define matrix A */
465397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
466397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
467397b6df1SKris Buschelman     if (!lu->myid) {
468c36ead0aSKris Buschelman       if (lua->isAIJ){
469397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
470397b6df1SKris Buschelman         nz               = aa->nz;
471397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
472397b6df1SKris Buschelman       } else {
473397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
4746c6c5352SBarry Smith         nz                  =  aa->nz;
475397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
476397b6df1SKris Buschelman       }
477397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
4787c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
4797c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
480397b6df1SKris Buschelman         nz = 0;
481397b6df1SKris Buschelman         for (i=0; i<M; i++){
482397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
483397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
484397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
485397b6df1SKris Buschelman           }
486397b6df1SKris Buschelman         }
487397b6df1SKris Buschelman       }
488397b6df1SKris Buschelman     }
489397b6df1SKris Buschelman     break;
490397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
491397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
492397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
493397b6df1SKris Buschelman     } else {
494397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
495397b6df1SKris Buschelman     }
496397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
497397b6df1SKris Buschelman     break;
498397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
499397b6df1SKris Buschelman   }
500397b6df1SKris Buschelman 
501397b6df1SKris Buschelman   /* analysis phase */
502329ec9b3SHong Zhang   /*----------------*/
503397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
504329ec9b3SHong Zhang     lu->id.job = 1;
505329ec9b3SHong Zhang 
506397b6df1SKris Buschelman     lu->id.n = M;
507397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
508397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
509397b6df1SKris Buschelman       if (!lu->myid) {
510397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
511397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
512397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
513397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
514397b6df1SKris Buschelman #else
515397b6df1SKris Buschelman           lu->id.a = lu->val;
516397b6df1SKris Buschelman #endif
517397b6df1SKris Buschelman         }
518397b6df1SKris Buschelman       }
519397b6df1SKris Buschelman       break;
520397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
521397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
522397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
523397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
524397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
525397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
526397b6df1SKris Buschelman #else
527397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
528397b6df1SKris Buschelman #endif
529397b6df1SKris Buschelman       }
530329ec9b3SHong Zhang       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
531329ec9b3SHong Zhang       IS  is_iden;
532329ec9b3SHong Zhang       Vec b;
533329ec9b3SHong Zhang       if (!lu->myid){
534329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr);
535329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr);
536329ec9b3SHong Zhang       } else {
537329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
538329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
539329ec9b3SHong Zhang       }
540329ec9b3SHong Zhang       ierr = VecCreate(A->comm,&b);CHKERRQ(ierr);
541329ec9b3SHong Zhang       ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
542329ec9b3SHong Zhang       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
543329ec9b3SHong Zhang 
544329ec9b3SHong Zhang       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
545329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
546329ec9b3SHong Zhang       ierr = VecDestroy(b);CHKERRQ(ierr);
547397b6df1SKris Buschelman       break;
548397b6df1SKris Buschelman     }
549397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
550397b6df1SKris Buschelman     zmumps_c(&lu->id);
551397b6df1SKris Buschelman #else
552397b6df1SKris Buschelman     dmumps_c(&lu->id);
553397b6df1SKris Buschelman #endif
554397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
55579a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
556397b6df1SKris Buschelman     }
557397b6df1SKris Buschelman   }
558397b6df1SKris Buschelman 
559397b6df1SKris Buschelman   /* numerical factorization phase */
560329ec9b3SHong Zhang   /*-------------------------------*/
561329ec9b3SHong Zhang   lu->id.job = 2;
562958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
563a7aca84bSHong Zhang     if (!lu->myid) {
564397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
565397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
566397b6df1SKris Buschelman #else
567397b6df1SKris Buschelman       lu->id.a = lu->val;
568397b6df1SKris Buschelman #endif
569397b6df1SKris Buschelman     }
570397b6df1SKris Buschelman   } else {
571397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
572397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
573397b6df1SKris Buschelman #else
574397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
575397b6df1SKris Buschelman #endif
576397b6df1SKris Buschelman   }
577397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
578397b6df1SKris Buschelman   zmumps_c(&lu->id);
579397b6df1SKris Buschelman #else
580397b6df1SKris Buschelman   dmumps_c(&lu->id);
581397b6df1SKris Buschelman #endif
582397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
58319facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
58419facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
58519facb7aSBarry Smith     } else {
58679a5c55eSBarry 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));
587397b6df1SKris Buschelman     }
58819facb7aSBarry Smith   }
589397b6df1SKris Buschelman 
59019facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
59179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
592397b6df1SKris Buschelman   }
593397b6df1SKris Buschelman 
5948ada1bb4SHong Zhang   if (lu->size > 1){
595e09efc27SHong Zhang     if ((*F)->factor == FACTOR_LU){
596e09efc27SHong Zhang       F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
597e09efc27SHong Zhang     } else {
598e09efc27SHong Zhang       F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
599e09efc27SHong Zhang     }
600e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
601329ec9b3SHong Zhang     if (lu->nSolve){
602329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
603329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
604329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
605329ec9b3SHong Zhang     }
6068ada1bb4SHong Zhang   }
607397b6df1SKris Buschelman   (*F)->assembled   = PETSC_TRUE;
608397b6df1SKris Buschelman   lu->matstruc      = SAME_NONZERO_PATTERN;
609ace87b0dSHong Zhang   lu->CleanUpMUMPS  = PETSC_TRUE;
610329ec9b3SHong Zhang   lu->nSolve        = 0;
611397b6df1SKris Buschelman   PetscFunctionReturn(0);
612397b6df1SKris Buschelman }
613397b6df1SKris Buschelman 
614397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
615397b6df1SKris Buschelman #undef __FUNCT__
616f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
617dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
618397b6df1SKris Buschelman   Mat            B;
619f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
620dfbe8321SBarry Smith   PetscErrorCode ierr;
621397b6df1SKris Buschelman 
622397b6df1SKris Buschelman   PetscFunctionBegin;
623397b6df1SKris Buschelman   /* Create the factorization matrix */
624f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
6252a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
626be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
627397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
628397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
629397b6df1SKris Buschelman 
630f6c57405SHong Zhang   B->ops->lufactornumeric = MatFactorNumeric_MUMPS;
631397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
632f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
633397b6df1SKris Buschelman   lu->sym                 = 0;
634397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
635397b6df1SKris Buschelman 
636397b6df1SKris Buschelman   *F = B;
637397b6df1SKris Buschelman   PetscFunctionReturn(0);
638397b6df1SKris Buschelman }
639397b6df1SKris Buschelman 
640397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
641397b6df1SKris Buschelman #undef __FUNCT__
642f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
643dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
644397b6df1SKris Buschelman   Mat            B;
645f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
646dfbe8321SBarry Smith   PetscErrorCode ierr;
647397b6df1SKris Buschelman 
648397b6df1SKris Buschelman   PetscFunctionBegin;
649397b6df1SKris Buschelman   /* Create the factorization matrix */
650f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
6512a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
652be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
653efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
654efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
655397b6df1SKris Buschelman 
656f6c57405SHong Zhang   B->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
657a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
658397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
659f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
660397b6df1SKris Buschelman   lu->sym                       = 2;
661397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
662397b6df1SKris Buschelman 
663397b6df1SKris Buschelman   *F = B;
664397b6df1SKris Buschelman   PetscFunctionReturn(0);
665397b6df1SKris Buschelman }
666397b6df1SKris Buschelman 
667397b6df1SKris Buschelman #undef __FUNCT__
668f6c57405SHong Zhang #define __FUNCT__ "MatFactorInfo_MUMPS"
669f6c57405SHong Zhang PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
670f6c57405SHong Zhang   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
671f6c57405SHong Zhang   PetscErrorCode ierr;
672f6c57405SHong Zhang 
673f6c57405SHong Zhang   PetscFunctionBegin;
674f6c57405SHong Zhang   /* check if matrix is mumps type */
675f6c57405SHong Zhang   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
676f6c57405SHong Zhang 
677f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
678f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
679f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
680f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
681f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
682f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
683f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
684f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
685f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
686f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
687f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
688f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
689f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
690f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
691f6c57405SHong Zhang   if (!lu->myid && lu->id.ICNTL(11)>0) {
692f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
693f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
694f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
695f6c57405SHong 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);
696f6c57405SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
697f6c57405SHong 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);
698f6c57405SHong Zhang 
699f6c57405SHong Zhang   }
700f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
701f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
702f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
703f6c57405SHong Zhang   /* ICNTL(15-17) not used */
704f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
705f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
706f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
707f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
708f6c57405SHong Zhang 
709f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
710f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
711f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
712f6c57405SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
713f6c57405SHong Zhang 
714f6c57405SHong Zhang   /* infomation local to each processor */
715f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
716f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
717f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
718f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
719f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
720f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
721f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
722f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
723f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
724f6c57405SHong Zhang   /*
725f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);}
726f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr);
727f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
728f6c57405SHong Zhang   */
729f6c57405SHong Zhang 
730f6c57405SHong 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);}
731f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
732f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
733f6c57405SHong Zhang 
734f6c57405SHong 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);}
735f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
736f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
737f6c57405SHong Zhang 
738f6c57405SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
739f6c57405SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
740f6c57405SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
741f6c57405SHong Zhang 
742f6c57405SHong Zhang   if (!lu->myid){ /* information from the host */
743f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
744f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
745f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
746f6c57405SHong Zhang 
747f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
748f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
749f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
750f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
751f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
752f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
753f6c57405SHong 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);
754f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
755f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
756f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
757f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
758f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
759f6c57405SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
760f6c57405SHong 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);
761f6c57405SHong 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);
762f6c57405SHong 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);
763f6c57405SHong 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);
764f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
765f6c57405SHong 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);
766f6c57405SHong 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);
767f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
768f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
769f6c57405SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
770f6c57405SHong Zhang   }
771f6c57405SHong Zhang 
772f6c57405SHong Zhang   PetscFunctionReturn(0);
773f6c57405SHong Zhang }
774f6c57405SHong Zhang 
775f6c57405SHong Zhang #undef __FUNCT__
776f6c57405SHong Zhang #define __FUNCT__ "MatView_MUMPS"
777f6c57405SHong Zhang PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
778f6c57405SHong Zhang   PetscErrorCode    ierr;
779f6c57405SHong Zhang   PetscTruth        iascii;
780f6c57405SHong Zhang   PetscViewerFormat format;
781f6c57405SHong Zhang   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
782f6c57405SHong Zhang 
783f6c57405SHong Zhang   PetscFunctionBegin;
784f6c57405SHong Zhang   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
785f6c57405SHong Zhang 
786f6c57405SHong Zhang   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
787f6c57405SHong Zhang   if (iascii) {
788f6c57405SHong Zhang     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
789f6c57405SHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO){
790f6c57405SHong Zhang       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
791f6c57405SHong Zhang     }
792f6c57405SHong Zhang   }
793f6c57405SHong Zhang   PetscFunctionReturn(0);
794f6c57405SHong Zhang }
795f6c57405SHong Zhang 
796f6c57405SHong Zhang #undef __FUNCT__
797f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
798dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
799dfbe8321SBarry Smith   PetscErrorCode ierr;
800f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
801c338a77dSKris Buschelman 
802397b6df1SKris Buschelman   PetscFunctionBegin;
803c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
804f0c56d0fSKris Buschelman 
805c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
806c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
807f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
808397b6df1SKris Buschelman   PetscFunctionReturn(0);
809397b6df1SKris Buschelman }
810397b6df1SKris Buschelman 
811c338a77dSKris Buschelman EXTERN_C_BEGIN
812c338a77dSKris Buschelman #undef __FUNCT__
813f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
81474a14a02SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
815dfbe8321SBarry Smith {
816dfbe8321SBarry Smith   PetscErrorCode ierr;
817521d7252SBarry Smith   PetscMPIInt    size;
818c338a77dSKris Buschelman   MPI_Comm       comm;
819c338a77dSKris Buschelman   Mat            B=*newmat;
820f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps;
821397b6df1SKris Buschelman 
822397b6df1SKris Buschelman   PetscFunctionBegin;
823c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
824*38f2d2fdSLisandro Dalcin   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
825c338a77dSKris Buschelman 
8265c088420SHong Zhang   if (reuse == MAT_INITIAL_MATRIX) {
8275c088420SHong Zhang     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
8285c088420SHong Zhang     /* A may have special container that is not duplicated,
8295c088420SHong Zhang        e.g., A is obtainted from MatMatMult(,&A). Save B->ops instead */
8305c088420SHong Zhang     mumps->MatDuplicate              = B->ops->duplicate;
8315c088420SHong Zhang     mumps->MatView                   = B->ops->view;
8325c088420SHong Zhang     mumps->MatAssemblyEnd            = B->ops->assemblyend;
8335c088420SHong Zhang     mumps->MatLUFactorSymbolic       = B->ops->lufactorsymbolic;
8345c088420SHong Zhang     mumps->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic;
8355c088420SHong Zhang     mumps->MatDestroy                = B->ops->destroy;
8365c088420SHong Zhang   } else {
837f0c56d0fSKris Buschelman     mumps->MatDuplicate              = A->ops->duplicate;
838c338a77dSKris Buschelman     mumps->MatView                   = A->ops->view;
839c338a77dSKris Buschelman     mumps->MatAssemblyEnd            = A->ops->assemblyend;
840c338a77dSKris Buschelman     mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
841c338a77dSKris Buschelman     mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
842c338a77dSKris Buschelman     mumps->MatDestroy                = A->ops->destroy;
8435c088420SHong Zhang   }
844a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
845c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
846f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
847c338a77dSKris Buschelman 
8484b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
849422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
850c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
851f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
852f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
8533924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
854c338a77dSKris Buschelman 
855c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
856c338a77dSKris Buschelman   if (size == 1) {
857c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
858f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
859c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
860c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
861c338a77dSKris Buschelman   } else {
862c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
863f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
864c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
865c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
866c338a77dSKris Buschelman   }
867c338a77dSKris Buschelman 
868ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr);
869c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
870c338a77dSKris Buschelman   *newmat = B;
871397b6df1SKris Buschelman   PetscFunctionReturn(0);
872397b6df1SKris Buschelman }
873c338a77dSKris Buschelman EXTERN_C_END
874397b6df1SKris Buschelman 
87524b6179bSKris Buschelman /*MC
876fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
87724b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
87824b6179bSKris Buschelman 
87924b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
88024b6179bSKris Buschelman   on how to declare the existence of external packages),
88124b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
88224b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
88324b6179bSKris Buschelman 
88424b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
88524b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
8863ec795f1SBarry Smith   MatSeqAIJSetPreallocation() is supported, and similarly MatMPIAIJSetPreallocation() is supported
88724b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
8883ec795f1SBarry Smith   the above preallocation routines for simplicity.  One can also call MatConvert() for an inplace
88928b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
89028b08bd3SKris Buschelman   without data copy.
89124b6179bSKris Buschelman 
89224b6179bSKris Buschelman   Options Database Keys:
8930bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
89424b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
89524b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
89624b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
89724b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
89824b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
89924b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
90094b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
90124b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
90224b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
90324b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
90424b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
90524b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
90624b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
90724b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
90824b6179bSKris Buschelman 
90924b6179bSKris Buschelman   Level: beginner
91024b6179bSKris Buschelman 
91124b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
91224b6179bSKris Buschelman M*/
91324b6179bSKris Buschelman 
914397b6df1SKris Buschelman EXTERN_C_BEGIN
915397b6df1SKris Buschelman #undef __FUNCT__
916f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
917be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A)
918dfbe8321SBarry Smith {
919dfbe8321SBarry Smith   PetscErrorCode ierr;
920c1490034SHong Zhang   PetscMPIInt    size;
921397b6df1SKris Buschelman 
922397b6df1SKris Buschelman   PetscFunctionBegin;
92317667f90SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
924397b6df1SKris Buschelman   if (size == 1) {
925397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
926397b6df1SKris Buschelman   } else {
927397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
92800ff2a26SHong Zhang     /*
92900ff2a26SHong Zhang     Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
930ceb03754SKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
93100ff2a26SHong Zhang     */
932397b6df1SKris Buschelman   }
933ceb03754SKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
934397b6df1SKris Buschelman   PetscFunctionReturn(0);
935397b6df1SKris Buschelman }
936397b6df1SKris Buschelman EXTERN_C_END
937397b6df1SKris Buschelman 
938f579278aSKris Buschelman #undef __FUNCT__
939f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
940dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
941dfbe8321SBarry Smith {
942dfbe8321SBarry Smith   PetscErrorCode ierr;
943f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
944f579278aSKris Buschelman 
945f579278aSKris Buschelman   PetscFunctionBegin;
946f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
947f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
948f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
949f579278aSKris Buschelman   PetscFunctionReturn(0);
950f579278aSKris Buschelman }
951f579278aSKris Buschelman 
952f579278aSKris Buschelman EXTERN_C_BEGIN
953f579278aSKris Buschelman #undef __FUNCT__
9549c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
955c1490034SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
9569c097c71SKris Buschelman {
9579c097c71SKris Buschelman   Mat       A;
95802217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
959dfbe8321SBarry Smith   PetscErrorCode ierr;
9609c097c71SKris Buschelman 
9619c097c71SKris Buschelman   PetscFunctionBegin;
9629c097c71SKris Buschelman   /*
9639c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
9649c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
9659c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
9669c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
9679c097c71SKris Buschelman     in the preallocation routine.
9689c097c71SKris Buschelman   */
9699c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
9709c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
971ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
9729c097c71SKris Buschelman   PetscFunctionReturn(0);
9739c097c71SKris Buschelman }
9749c097c71SKris Buschelman EXTERN_C_END
9759c097c71SKris Buschelman 
9769c097c71SKris Buschelman EXTERN_C_BEGIN
9779c097c71SKris Buschelman #undef __FUNCT__
978f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
97975179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
980dfbe8321SBarry Smith {
981dfbe8321SBarry Smith   PetscErrorCode ierr;
982521d7252SBarry Smith   PetscMPIInt    size;
983f579278aSKris Buschelman   MPI_Comm       comm;
984f579278aSKris Buschelman   Mat            B=*newmat;
985422e82a1SHong Zhang   Mat_MUMPS      *mumps;
9869c097c71SKris Buschelman   void           (*f)(void);
987f579278aSKris Buschelman 
988f579278aSKris Buschelman   PetscFunctionBegin;
989ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
990f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
991f579278aSKris Buschelman   }
992f579278aSKris Buschelman 
993f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
994*38f2d2fdSLisandro Dalcin   ierr = PetscNewLog(B,Mat_MUMPS,&mumps);CHKERRQ(ierr);
995f579278aSKris Buschelman 
996f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
997f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
998f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
999f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
1000f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
1001f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
1002a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
1003f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
1004f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
1005f579278aSKris Buschelman 
1006f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
1007422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
1008c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
1009f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
1010f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
10113924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
1012f579278aSKris Buschelman 
1013f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1014f579278aSKris Buschelman   if (size == 1) {
1015f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
1016f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1017f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
1018f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1019f579278aSKris Buschelman   } else {
10209c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
1021f68b968cSBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
1022901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
10236849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
10249c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
10259c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
10269c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
10279c097c71SKris Buschelman     }
1028f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
1029f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1030f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
1031f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1032f579278aSKris Buschelman   }
1033f579278aSKris Buschelman 
1034ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr);
1035f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
1036f579278aSKris Buschelman   *newmat = B;
1037f579278aSKris Buschelman   PetscFunctionReturn(0);
1038f579278aSKris Buschelman }
1039f579278aSKris Buschelman EXTERN_C_END
1040f579278aSKris Buschelman 
1041f0c56d0fSKris Buschelman #undef __FUNCT__
1042422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
1043dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1044dfbe8321SBarry Smith   PetscErrorCode ierr;
10458e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
10468f340917SKris Buschelman 
1047f0c56d0fSKris Buschelman   PetscFunctionBegin;
10488f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
10498e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
1050f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
1051f0c56d0fSKris Buschelman }
1052f0c56d0fSKris Buschelman 
105324b6179bSKris Buschelman /*MC
1054fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
105524b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
105624b6179bSKris Buschelman 
105724b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
105824b6179bSKris Buschelman   on how to declare the existence of external packages),
105924b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
106024b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
106124b6179bSKris Buschelman 
106224b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
106324b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
106424b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
106524b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
106628b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
106728b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
106828b08bd3SKris Buschelman   without data copy.
106924b6179bSKris Buschelman 
107024b6179bSKris Buschelman   Options Database Keys:
10710bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
107224b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
107324b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
107424b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
107524b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
107624b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
107724b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
107894b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
107924b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
108024b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
108124b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
108224b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
108324b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
108424b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
108524b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
108624b6179bSKris Buschelman 
108724b6179bSKris Buschelman   Level: beginner
108824b6179bSKris Buschelman 
108924b6179bSKris Buschelman .seealso: MATAIJMUMPS
109024b6179bSKris Buschelman M*/
109124b6179bSKris Buschelman 
1092397b6df1SKris Buschelman EXTERN_C_BEGIN
1093397b6df1SKris Buschelman #undef __FUNCT__
1094f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1095be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A)
1096dfbe8321SBarry Smith {
10976849ba73SBarry Smith   PetscErrorCode ierr;
10988ada1bb4SHong Zhang   PetscMPIInt    size;
1099397b6df1SKris Buschelman 
1100397b6df1SKris Buschelman   PetscFunctionBegin;
11015441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1102397b6df1SKris Buschelman   if (size == 1) {
1103397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1104397b6df1SKris Buschelman   } else {
1105397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1106397b6df1SKris Buschelman   }
1107ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1108397b6df1SKris Buschelman   PetscFunctionReturn(0);
1109397b6df1SKris Buschelman }
1110397b6df1SKris Buschelman EXTERN_C_END
1111