xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision ceb03754fa7dd360e20fd24445c18760198643f9)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3a58c3f20SHong Zhang     Provides an interface to the MUMPS_4.3.1 sparse solver
4397b6df1SKris Buschelman */
5397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
6397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h"
7397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h"
8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
12397b6df1SKris Buschelman #include "zmumps_c.h"
13397b6df1SKris Buschelman #else
14397b6df1SKris Buschelman #include "dmumps_c.h"
15397b6df1SKris Buschelman #endif
16397b6df1SKris Buschelman EXTERN_C_END
17397b6df1SKris Buschelman #define JOB_INIT -1
18397b6df1SKris Buschelman #define JOB_END -2
19397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
20397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
21397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
22397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
23a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
24397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
25adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
26397b6df1SKris Buschelman 
27397b6df1SKris Buschelman typedef struct {
28397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
29397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
30397b6df1SKris Buschelman #else
31397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
32397b6df1SKris Buschelman #endif
33397b6df1SKris Buschelman   MatStructure   matstruc;
34397b6df1SKris Buschelman   int            myid,size,*irn,*jcn,sym;
35397b6df1SKris Buschelman   PetscScalar    *val;
36397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
37397b6df1SKris Buschelman 
38c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
396849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
406849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
416849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
426849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
436849ba73SBarry Smith   PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
446849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
456849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
466849ba73SBarry Smith   PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
47f0c56d0fSKris Buschelman } Mat_MUMPS;
48f0c56d0fSKris Buschelman 
49dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
50892f6c3fSKris Buschelman EXTERN_C_BEGIN
51*ceb03754SKris Buschelman PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,MatReuse,Mat*);
52892f6c3fSKris Buschelman EXTERN_C_END
53397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
54397b6df1SKris Buschelman /*
55397b6df1SKris Buschelman   input:
5675747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
57397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
58397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
59397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
60397b6df1SKris Buschelman   output:
61397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
62397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
63397b6df1SKris Buschelman  */
64dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
65397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
66dfbe8321SBarry Smith   PetscErrorCode ierr;
67dfbe8321SBarry Smith   int         i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
68d54de34fSKris Buschelman   int         *row,*col;
69397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
70f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
71397b6df1SKris Buschelman 
72397b6df1SKris Buschelman   PetscFunctionBegin;
73397b6df1SKris Buschelman   if (mumps->isAIJ){
74397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
75397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
76397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
77397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
78397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
79397b6df1SKris Buschelman     garray = mat->garray;
80397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
81397b6df1SKris Buschelman 
82397b6df1SKris Buschelman   } else {
83397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
84397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
85397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
86521d7252SBarry Smith     if (A->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->bs);
876c6c5352SBarry Smith     nz = aa->nz + bb->nz;
88397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
89397b6df1SKris Buschelman     garray = mat->garray;
90397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
91397b6df1SKris Buschelman   }
92397b6df1SKris Buschelman 
93397b6df1SKris Buschelman   if (!valOnly){
94397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
95397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
96397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
97397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
98397b6df1SKris Buschelman   } else {
99397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
100397b6df1SKris Buschelman   }
101397b6df1SKris Buschelman   *nnz = nz;
102397b6df1SKris Buschelman 
103028e57e8SHong Zhang   jj = 0; irow = rstart;
104397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
105397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
106397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
107397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
108397b6df1SKris Buschelman     bjj = bj + bi[i];
109397b6df1SKris Buschelman 
110397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
111397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11275747be1SHong Zhang     j=-1;
11375747be1SHong Zhang     do {
11475747be1SHong Zhang       j++;
11575747be1SHong Zhang       if (j == countB) break;
116397b6df1SKris Buschelman       jcol = garray[bjj[j]];
11775747be1SHong Zhang     } while (jcol < colA_start);
11875747be1SHong Zhang     jB = j;
119397b6df1SKris Buschelman 
120397b6df1SKris Buschelman     /* B-part, smaller col index */
121397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
122397b6df1SKris Buschelman     for (j=0; j<jB; j++){
123397b6df1SKris Buschelman       jcol = garray[bjj[j]];
124397b6df1SKris Buschelman       if (!valOnly){
125397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12675747be1SHong Zhang 
127397b6df1SKris Buschelman       }
128397b6df1SKris Buschelman       val[jj++] = *bv++;
129397b6df1SKris Buschelman     }
130397b6df1SKris Buschelman     /* A-part */
131397b6df1SKris Buschelman     for (j=0; j<countA; j++){
132397b6df1SKris Buschelman       if (!valOnly){
133397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
134397b6df1SKris Buschelman       }
135397b6df1SKris Buschelman       val[jj++] = *av++;
136397b6df1SKris Buschelman     }
137397b6df1SKris Buschelman     /* B-part, larger col index */
138397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
139397b6df1SKris Buschelman       if (!valOnly){
140397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
141397b6df1SKris Buschelman       }
142397b6df1SKris Buschelman       val[jj++] = *bv++;
143397b6df1SKris Buschelman     }
144397b6df1SKris Buschelman     irow++;
145397b6df1SKris Buschelman   }
146397b6df1SKris Buschelman 
147397b6df1SKris Buschelman   PetscFunctionReturn(0);
148397b6df1SKris Buschelman }
149397b6df1SKris Buschelman 
150c338a77dSKris Buschelman EXTERN_C_BEGIN
151c338a77dSKris Buschelman #undef __FUNCT__
152c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
153*ceb03754SKris Buschelman PetscErrorCode MatConvert_MUMPS_Base(Mat A,const MatType type,MatReuse reuse,Mat *newmat) \
154521d7252SBarry Smith {
155dfbe8321SBarry Smith   PetscErrorCode ierr;
156c338a77dSKris Buschelman   Mat            B=*newmat;
157f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
158901853e0SKris Buschelman   void           (*f)(void);
159c338a77dSKris Buschelman 
160c338a77dSKris Buschelman   PetscFunctionBegin;
161*ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
162c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
163c338a77dSKris Buschelman   }
164f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
165f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
166f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
167f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
168f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
169f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
170901853e0SKris Buschelman 
171901853e0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
172901853e0SKris Buschelman   if (f) {
173d1390fe9SSatish Balay     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(FCNVOID)mumps->MatPreallocate);CHKERRQ(ierr);
174901853e0SKris Buschelman   }
175f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
176901853e0SKris Buschelman 
177901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
178901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
179901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
180901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
1812895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
1822895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
183901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
184901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
185901853e0SKris Buschelman 
186901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
187c338a77dSKris Buschelman   *newmat = B;
188c338a77dSKris Buschelman   PetscFunctionReturn(0);
189c338a77dSKris Buschelman }
190c338a77dSKris Buschelman EXTERN_C_END
191c338a77dSKris Buschelman 
192397b6df1SKris Buschelman #undef __FUNCT__
1933924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
194dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
195dfbe8321SBarry Smith {
196f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
197dfbe8321SBarry Smith   PetscErrorCode ierr;
198dfbe8321SBarry Smith   int       size=lu->size;
1996849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
200397b6df1SKris Buschelman   PetscFunctionBegin;
201397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
202397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
203397b6df1SKris Buschelman     lu->id.job=JOB_END;
204397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
205397b6df1SKris Buschelman     zmumps_c(&lu->id);
206397b6df1SKris Buschelman #else
207397b6df1SKris Buschelman     dmumps_c(&lu->id);
208397b6df1SKris Buschelman #endif
209c338a77dSKris Buschelman     if (lu->irn) {
210c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
211c338a77dSKris Buschelman     }
212c338a77dSKris Buschelman     if (lu->jcn) {
213c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
214c338a77dSKris Buschelman     }
215c338a77dSKris Buschelman     if (size>1 && lu->val) {
216c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
217c338a77dSKris Buschelman     }
218397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
219397b6df1SKris Buschelman   }
220a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
221a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
222c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
223397b6df1SKris Buschelman   PetscFunctionReturn(0);
224397b6df1SKris Buschelman }
225397b6df1SKris Buschelman 
226397b6df1SKris Buschelman #undef __FUNCT__
227a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
228dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
229dfbe8321SBarry Smith {
2306849ba73SBarry Smith   PetscErrorCode ierr;
2316849ba73SBarry Smith   int  size;
232a39386dcSKris Buschelman 
233a39386dcSKris Buschelman   PetscFunctionBegin;
234a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
235a39386dcSKris Buschelman   if (size==1) {
236*ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
237a39386dcSKris Buschelman   } else {
238*ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
239a39386dcSKris Buschelman   }
240a39386dcSKris Buschelman   PetscFunctionReturn(0);
241a39386dcSKris Buschelman }
242a39386dcSKris Buschelman 
243a39386dcSKris Buschelman #undef __FUNCT__
244a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
245dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
246dfbe8321SBarry Smith {
2476849ba73SBarry Smith   PetscErrorCode ierr;
2486849ba73SBarry Smith   int  size;
249a39386dcSKris Buschelman 
250a39386dcSKris Buschelman   PetscFunctionBegin;
251a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
252a39386dcSKris Buschelman   if (size==1) {
253*ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
254a39386dcSKris Buschelman   } else {
255*ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
256a39386dcSKris Buschelman   }
257a39386dcSKris Buschelman   PetscFunctionReturn(0);
258a39386dcSKris Buschelman }
259a39386dcSKris Buschelman 
260a39386dcSKris Buschelman #undef __FUNCT__
261c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
262dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
263f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
264dfbe8321SBarry Smith   PetscErrorCode ierr;
265397b6df1SKris Buschelman 
266397b6df1SKris Buschelman   PetscFunctionBegin;
267c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
268c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
269c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
270c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
271c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
272c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
273c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
274c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
275c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
276c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
277958c9bccSBarry Smith   if (!lu->myid && lu->id.ICNTL(11)>0) {
278c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
279c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
280c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
281c338a77dSKris Buschelman     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);
282c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
283c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",lu->id.RINFOG(10),lu->id.RINFOG(11));CHKERRQ(ierr);
284c338a77dSKris Buschelman 
285c338a77dSKris Buschelman   }
286c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
287c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
288adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
289c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
290c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
291c338a77dSKris Buschelman 
292c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
293c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
294c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
29557f0c58bSHong Zhang 
29657f0c58bSHong Zhang   /* infomation local to each processor */
297958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
29857f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
29957f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
300958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
30157f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
30257f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
303958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
30457f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
30557f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
306adc1d99fSHong Zhang 
307958c9bccSBarry Smith   if (!lu->myid){ /* information from the host */
308adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
309adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
310adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
311adc1d99fSHong Zhang 
312adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
313adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
314adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
315adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
316adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
317adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
318adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
319adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
320adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
321adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
322adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
323adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
324adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
325adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(16) (estimated size (in million of bytes) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr);
326adc1d99fSHong 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);
327adc1d99fSHong 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);
328adc1d99fSHong 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);
329adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
330adc1d99fSHong Zhang   }
331adc1d99fSHong Zhang 
332397b6df1SKris Buschelman   PetscFunctionReturn(0);
333397b6df1SKris Buschelman }
334397b6df1SKris Buschelman 
335397b6df1SKris Buschelman #undef __FUNCT__
336f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
337dfbe8321SBarry Smith PetscErrorCode MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
338dfbe8321SBarry Smith   PetscErrorCode ierr;
33932077d6dSBarry Smith   PetscTruth        iascii;
340397b6df1SKris Buschelman   PetscViewerFormat format;
341f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
342397b6df1SKris Buschelman 
343397b6df1SKris Buschelman   PetscFunctionBegin;
344397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
345397b6df1SKris Buschelman 
34632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
34732077d6dSBarry Smith   if (iascii) {
348397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
349397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
350397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
351397b6df1SKris Buschelman     }
352397b6df1SKris Buschelman   }
353397b6df1SKris Buschelman   PetscFunctionReturn(0);
354397b6df1SKris Buschelman }
355397b6df1SKris Buschelman 
356397b6df1SKris Buschelman #undef __FUNCT__
357f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
358dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
359f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
360d54de34fSKris Buschelman   PetscScalar *array;
361397b6df1SKris Buschelman   Vec         x_seq;
362397b6df1SKris Buschelman   IS          iden;
363397b6df1SKris Buschelman   VecScatter  scat;
364dfbe8321SBarry Smith   PetscErrorCode ierr;
365397b6df1SKris Buschelman 
366397b6df1SKris Buschelman   PetscFunctionBegin;
367397b6df1SKris Buschelman   if (lu->size > 1){
368397b6df1SKris Buschelman     if (!lu->myid){
369397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
370397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
371397b6df1SKris Buschelman     } else {
372397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
373397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
374397b6df1SKris Buschelman     }
375397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
376397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
377397b6df1SKris Buschelman 
378397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
379397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
380397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
381397b6df1SKris Buschelman   } else {  /* size == 1 */
382397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
383397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
384397b6df1SKris Buschelman   }
385397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
386397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
387397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
388397b6df1SKris Buschelman #else
389397b6df1SKris Buschelman     lu->id.rhs = array;
390397b6df1SKris Buschelman #endif
391397b6df1SKris Buschelman   }
392397b6df1SKris Buschelman 
393397b6df1SKris Buschelman   /* solve phase */
394397b6df1SKris Buschelman   lu->id.job=3;
395397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
396397b6df1SKris Buschelman   zmumps_c(&lu->id);
397397b6df1SKris Buschelman #else
398397b6df1SKris Buschelman   dmumps_c(&lu->id);
399397b6df1SKris Buschelman #endif
400397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
40179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
402397b6df1SKris Buschelman   }
403397b6df1SKris Buschelman 
404397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
405397b6df1SKris Buschelman   if (lu->size > 1) {
406397b6df1SKris Buschelman     if (!lu->myid){
407397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
408397b6df1SKris Buschelman     }
409397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
410397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
411397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
412397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
413397b6df1SKris Buschelman   } else {
414397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
415397b6df1SKris Buschelman   }
416397b6df1SKris Buschelman 
417397b6df1SKris Buschelman   PetscFunctionReturn(0);
418397b6df1SKris Buschelman }
419397b6df1SKris Buschelman 
420a58c3f20SHong Zhang /*
421a58c3f20SHong Zhang   input:
422a58c3f20SHong Zhang    F:        numeric factor
423a58c3f20SHong Zhang   output:
424a58c3f20SHong Zhang    nneg:     total number of negative pivots
425a58c3f20SHong Zhang    nzero:    0
426a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
427a58c3f20SHong Zhang */
428a58c3f20SHong Zhang 
429a58c3f20SHong Zhang #undef __FUNCT__
430a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
431dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
432a58c3f20SHong Zhang {
433a58c3f20SHong Zhang   Mat_MUMPS  *lu =(Mat_MUMPS*)F->spptr;
434dfbe8321SBarry Smith   PetscErrorCode ierr;
435dfbe8321SBarry Smith   int        size;
436a58c3f20SHong Zhang 
437a58c3f20SHong Zhang   PetscFunctionBegin;
438bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
439bcb30aebSHong 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 */
440bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
44179a5c55eSBarry 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));
442bcb30aebSHong Zhang   }
443a58c3f20SHong Zhang   if (nneg){
444a58c3f20SHong Zhang     if (!lu->myid){
445a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
446a58c3f20SHong Zhang     }
447bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
448a58c3f20SHong Zhang   }
449a58c3f20SHong Zhang   if (nzero) *nzero = 0;
450a58c3f20SHong Zhang   if (npos)  *npos  = F->M - (*nneg);
451a58c3f20SHong Zhang   PetscFunctionReturn(0);
452a58c3f20SHong Zhang }
453a58c3f20SHong Zhang 
454397b6df1SKris Buschelman #undef __FUNCT__
45519facb7aSBarry Smith #define __FUNCT__ "MatFactorNumeric_AIJMUMPS"
456af281ebdSHong Zhang PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F)
457af281ebdSHong Zhang {
458f0c56d0fSKris Buschelman   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
459f0c56d0fSKris Buschelman   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
4606849ba73SBarry Smith   PetscErrorCode ierr;
461af281ebdSHong Zhang   PetscInt       rnz,nnz,nz,i,M=A->M,*ai,*aj,icntl;
462397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
463397b6df1SKris Buschelman 
464397b6df1SKris Buschelman   PetscFunctionBegin;
465397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
466f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
467397b6df1SKris Buschelman 
468397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
469397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
470397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
47175747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
472397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
473397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
474397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
475397b6df1SKris Buschelman 
476397b6df1SKris Buschelman     /* Set mumps options */
477397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
478397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
479397b6df1SKris Buschelman     lu->id.sym=lu->sym;
480397b6df1SKris Buschelman     if (lu->sym == 2){
481397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
482397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
483397b6df1SKris Buschelman     }
484397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
485397b6df1SKris Buschelman     zmumps_c(&lu->id);
486397b6df1SKris Buschelman #else
487397b6df1SKris Buschelman     dmumps_c(&lu->id);
488397b6df1SKris Buschelman #endif
489397b6df1SKris Buschelman 
490397b6df1SKris Buschelman     if (lu->size == 1){
491397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
492397b6df1SKris Buschelman     } else {
493397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
494397b6df1SKris Buschelman     }
495397b6df1SKris Buschelman 
496397b6df1SKris Buschelman     icntl=-1;
497397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
49819facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
499397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
500397b6df1SKris Buschelman     } else { /* no output */
501397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
502397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
503397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
504397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
505397b6df1SKris Buschelman     }
506397b6df1SKris 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);
507397b6df1SKris Buschelman     icntl=-1;
508397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
509397b6df1SKris Buschelman     if (flg) {
510397b6df1SKris Buschelman       if (icntl== 1){
511397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
512397b6df1SKris Buschelman       } else {
513397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
514397b6df1SKris Buschelman       }
515397b6df1SKris Buschelman     }
516397b6df1SKris 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);
517397b6df1SKris 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);
51894b7f48cSBarry 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);
519397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
520397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
521adc1d99fSHong 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);
522397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
523397b6df1SKris Buschelman 
524397b6df1SKris Buschelman     /*
525397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_16","ICNTL(16): 1: rank detection; 2: rank detection and nullspace","None",lu->id.ICNTL(16),&icntl,&flg);CHKERRQ(ierr);
526397b6df1SKris Buschelman     if (flg){
527397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
528397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
529397b6df1SKris Buschelman       } else {
530397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
531397b6df1SKris Buschelman       }
532397b6df1SKris Buschelman     }
533397b6df1SKris Buschelman     */
534397b6df1SKris Buschelman 
535397b6df1SKris 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);
536397b6df1SKris 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);
537397b6df1SKris 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);
538397b6df1SKris Buschelman     PetscOptionsEnd();
539397b6df1SKris Buschelman   }
540397b6df1SKris Buschelman 
541397b6df1SKris Buschelman   /* define matrix A */
542397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
543397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
544397b6df1SKris Buschelman     if (!lu->myid) {
545c36ead0aSKris Buschelman       if (lua->isAIJ){
546397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
547397b6df1SKris Buschelman         nz               = aa->nz;
548397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
549397b6df1SKris Buschelman       } else {
550397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5516c6c5352SBarry Smith         nz                  =  aa->nz;
552397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
553397b6df1SKris Buschelman       }
554397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
555397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
556397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
557397b6df1SKris Buschelman         nz = 0;
558397b6df1SKris Buschelman         for (i=0; i<M; i++){
559397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
560397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
561397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
562397b6df1SKris Buschelman           }
563397b6df1SKris Buschelman         }
564397b6df1SKris Buschelman       }
565397b6df1SKris Buschelman     }
566397b6df1SKris Buschelman     break;
567397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
568397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
569397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
570397b6df1SKris Buschelman     } else {
571397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
572397b6df1SKris Buschelman     }
573397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
574397b6df1SKris Buschelman     break;
575397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
576397b6df1SKris Buschelman   }
577397b6df1SKris Buschelman 
578397b6df1SKris Buschelman   /* analysis phase */
579397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
580397b6df1SKris Buschelman      lu->id.n = M;
581397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
582397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
583397b6df1SKris Buschelman       if (!lu->myid) {
584397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
585397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
586397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
587397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
588397b6df1SKris Buschelman #else
589397b6df1SKris Buschelman           lu->id.a = lu->val;
590397b6df1SKris Buschelman #endif
591397b6df1SKris Buschelman         }
592397b6df1SKris Buschelman       }
593397b6df1SKris Buschelman       break;
594397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
595397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
596397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
597397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
598397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
599397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
600397b6df1SKris Buschelman #else
601397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
602397b6df1SKris Buschelman #endif
603397b6df1SKris Buschelman       }
604397b6df1SKris Buschelman       break;
605397b6df1SKris Buschelman     }
606397b6df1SKris Buschelman     lu->id.job=1;
607397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
608397b6df1SKris Buschelman     zmumps_c(&lu->id);
609397b6df1SKris Buschelman #else
610397b6df1SKris Buschelman     dmumps_c(&lu->id);
611397b6df1SKris Buschelman #endif
612397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
61379a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
614397b6df1SKris Buschelman     }
615397b6df1SKris Buschelman   }
616397b6df1SKris Buschelman 
617397b6df1SKris Buschelman   /* numerical factorization phase */
618958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
619a7aca84bSHong Zhang     if (!lu->myid) {
620397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
621397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
622397b6df1SKris Buschelman #else
623397b6df1SKris Buschelman       lu->id.a = lu->val;
624397b6df1SKris Buschelman #endif
625397b6df1SKris Buschelman     }
626397b6df1SKris Buschelman   } else {
627397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
628397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
629397b6df1SKris Buschelman #else
630397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
631397b6df1SKris Buschelman #endif
632397b6df1SKris Buschelman   }
633397b6df1SKris Buschelman   lu->id.job=2;
634397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
635397b6df1SKris Buschelman   zmumps_c(&lu->id);
636397b6df1SKris Buschelman #else
637397b6df1SKris Buschelman   dmumps_c(&lu->id);
638397b6df1SKris Buschelman #endif
639397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
64019facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
64119facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
64219facb7aSBarry Smith     } else {
64379a5c55eSBarry 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));
644397b6df1SKris Buschelman     }
64519facb7aSBarry Smith   }
646397b6df1SKris Buschelman 
64719facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
64879a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
649397b6df1SKris Buschelman   }
650397b6df1SKris Buschelman 
651397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
652397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
653ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
654397b6df1SKris Buschelman   PetscFunctionReturn(0);
655397b6df1SKris Buschelman }
656397b6df1SKris Buschelman 
657397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
658397b6df1SKris Buschelman #undef __FUNCT__
659f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
660dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
661397b6df1SKris Buschelman   Mat            B;
662f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
663dfbe8321SBarry Smith   PetscErrorCode ierr;
664397b6df1SKris Buschelman 
665397b6df1SKris Buschelman   PetscFunctionBegin;
666397b6df1SKris Buschelman 
667397b6df1SKris Buschelman   /* Create the factorization matrix */
668397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
669be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
670397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
671397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
672397b6df1SKris Buschelman 
673f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
674397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
675f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
676397b6df1SKris Buschelman   lu->sym                 = 0;
677397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
678397b6df1SKris Buschelman 
679397b6df1SKris Buschelman   *F = B;
680397b6df1SKris Buschelman   PetscFunctionReturn(0);
681397b6df1SKris Buschelman }
682397b6df1SKris Buschelman 
683397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
684397b6df1SKris Buschelman #undef __FUNCT__
685f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
686dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
687397b6df1SKris Buschelman   Mat       B;
688f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
689dfbe8321SBarry Smith   PetscErrorCode ierr;
690397b6df1SKris Buschelman 
691397b6df1SKris Buschelman   PetscFunctionBegin;
692397b6df1SKris Buschelman   /* Create the factorization matrix */
693397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
694be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
695efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
696efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
697397b6df1SKris Buschelman 
698f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
699a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
700397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
701f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
702397b6df1SKris Buschelman   lu->sym                       = 2;
703397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
704397b6df1SKris Buschelman 
705397b6df1SKris Buschelman   *F = B;
706397b6df1SKris Buschelman   PetscFunctionReturn(0);
707397b6df1SKris Buschelman }
708397b6df1SKris Buschelman 
709397b6df1SKris Buschelman #undef __FUNCT__
710f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
711dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
712dfbe8321SBarry Smith   PetscErrorCode ierr;
713f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
714c338a77dSKris Buschelman 
715397b6df1SKris Buschelman   PetscFunctionBegin;
716c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
717f0c56d0fSKris Buschelman 
718c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
719c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
720f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
721397b6df1SKris Buschelman   PetscFunctionReturn(0);
722397b6df1SKris Buschelman }
723397b6df1SKris Buschelman 
724c338a77dSKris Buschelman EXTERN_C_BEGIN
725c338a77dSKris Buschelman #undef __FUNCT__
726f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
727*ceb03754SKris Buschelman PetscErrorCode MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
728dfbe8321SBarry Smith {
729dfbe8321SBarry Smith   PetscErrorCode ierr;
730521d7252SBarry Smith   PetscMPIInt    size;
731c338a77dSKris Buschelman   MPI_Comm       comm;
732c338a77dSKris Buschelman   Mat            B=*newmat;
733f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps;
734397b6df1SKris Buschelman 
735397b6df1SKris Buschelman   PetscFunctionBegin;
736*ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
737c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
738397b6df1SKris Buschelman   }
739397b6df1SKris Buschelman 
740c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
741f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
742c338a77dSKris Buschelman 
743f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
744c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
745c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
746c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
747c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
748c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
749a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
750c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
751f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
752c338a77dSKris Buschelman 
7534b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
754422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
755f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
756f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
757f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
7583924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
759c338a77dSKris Buschelman 
760c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
761c338a77dSKris Buschelman   if (size == 1) {
762c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
763f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
764c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
765c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
766c338a77dSKris Buschelman   } else {
767c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
768f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
769c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
770c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
771c338a77dSKris Buschelman   }
772c338a77dSKris Buschelman 
77352e6d16bSBarry Smith   PetscLogInfo(0,"MatConvert_AIJ_AIJMUMPS:Using MUMPS for LU factorization and solves.");
774c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
775c338a77dSKris Buschelman   *newmat = B;
776397b6df1SKris Buschelman   PetscFunctionReturn(0);
777397b6df1SKris Buschelman }
778c338a77dSKris Buschelman EXTERN_C_END
779397b6df1SKris Buschelman 
78024b6179bSKris Buschelman /*MC
781fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
78224b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
78324b6179bSKris Buschelman 
78424b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
78524b6179bSKris Buschelman   on how to declare the existence of external packages),
78624b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
78724b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
78824b6179bSKris Buschelman 
78924b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
79024b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
79124b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
79224b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
79328b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
79428b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
79528b08bd3SKris Buschelman   without data copy.
79624b6179bSKris Buschelman 
79724b6179bSKris Buschelman   Options Database Keys:
7980bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
79924b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
80024b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
80124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
80224b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
80324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
80424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
80594b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
80624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
80724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
80824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
80924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
81024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
81124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
81224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
81324b6179bSKris Buschelman 
81424b6179bSKris Buschelman   Level: beginner
81524b6179bSKris Buschelman 
81624b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
81724b6179bSKris Buschelman M*/
81824b6179bSKris Buschelman 
819397b6df1SKris Buschelman EXTERN_C_BEGIN
820397b6df1SKris Buschelman #undef __FUNCT__
821f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
822dfbe8321SBarry Smith PetscErrorCode MatCreate_AIJMUMPS(Mat A)
823dfbe8321SBarry Smith {
824dfbe8321SBarry Smith   PetscErrorCode ierr;
825dfbe8321SBarry Smith   int      size;
826e2d9671bSKris Buschelman   Mat      A_diag;
827397b6df1SKris Buschelman   MPI_Comm comm;
828397b6df1SKris Buschelman 
829397b6df1SKris Buschelman   PetscFunctionBegin;
8305441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
8315441df8eSKris Buschelman   /*   and AIJMUMPS types */
8325441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
833397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
834397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
835397b6df1SKris Buschelman   if (size == 1) {
836397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
837397b6df1SKris Buschelman   } else {
838397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
839e2d9671bSKris Buschelman     A_diag = ((Mat_MPIAIJ *)A->data)->A;
840*ceb03754SKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
841397b6df1SKris Buschelman   }
842*ceb03754SKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
843397b6df1SKris Buschelman   PetscFunctionReturn(0);
844397b6df1SKris Buschelman }
845397b6df1SKris Buschelman EXTERN_C_END
846397b6df1SKris Buschelman 
847f579278aSKris Buschelman #undef __FUNCT__
848f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
849dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
850dfbe8321SBarry Smith {
851dfbe8321SBarry Smith   PetscErrorCode ierr;
852f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
853f579278aSKris Buschelman 
854f579278aSKris Buschelman   PetscFunctionBegin;
855f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
856f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
857f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
858f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
859f579278aSKris Buschelman   PetscFunctionReturn(0);
860f579278aSKris Buschelman }
861f579278aSKris Buschelman 
862f579278aSKris Buschelman EXTERN_C_BEGIN
863f579278aSKris Buschelman #undef __FUNCT__
8649c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
865dfbe8321SBarry Smith PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat  B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz)
8669c097c71SKris Buschelman {
8679c097c71SKris Buschelman   Mat       A;
86802217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
869dfbe8321SBarry Smith   PetscErrorCode ierr;
8709c097c71SKris Buschelman 
8719c097c71SKris Buschelman   PetscFunctionBegin;
8729c097c71SKris Buschelman   /*
8739c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
8749c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
8759c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
8769c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
8779c097c71SKris Buschelman     in the preallocation routine.
8789c097c71SKris Buschelman   */
8799c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
8809c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
881*ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
8829c097c71SKris Buschelman   PetscFunctionReturn(0);
8839c097c71SKris Buschelman }
8849c097c71SKris Buschelman EXTERN_C_END
8859c097c71SKris Buschelman 
8869c097c71SKris Buschelman EXTERN_C_BEGIN
8879c097c71SKris Buschelman #undef __FUNCT__
888f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
889*ceb03754SKris Buschelman PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,MatReuse reuse,Mat *newmat)
890dfbe8321SBarry Smith {
891dfbe8321SBarry Smith   PetscErrorCode ierr;
892521d7252SBarry Smith   PetscMPIInt    size;
893f579278aSKris Buschelman   MPI_Comm       comm;
894f579278aSKris Buschelman   Mat            B=*newmat;
895422e82a1SHong Zhang   Mat_MUMPS      *mumps;
8969c097c71SKris Buschelman   void           (*f)(void);
897f579278aSKris Buschelman 
898f579278aSKris Buschelman   PetscFunctionBegin;
899*ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
900f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
901f579278aSKris Buschelman   }
902f579278aSKris Buschelman 
903f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
904f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
905f579278aSKris Buschelman 
906f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
907f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
908f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
909f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
910f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
911f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
912a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
913f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
914f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
915f579278aSKris Buschelman 
916f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
917422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
918f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
919f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
920f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
9213924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
922f579278aSKris Buschelman 
923f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
924f579278aSKris Buschelman   if (size == 1) {
925f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
926f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
927f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
928f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
929f579278aSKris Buschelman   } else {
9309c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
9319c097c71SKris Buschelman     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
932901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
9336849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
9349c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
9359c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
9369c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
9379c097c71SKris Buschelman     }
938f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
939f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
940f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
941f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
942f579278aSKris Buschelman   }
943f579278aSKris Buschelman 
94452e6d16bSBarry Smith   PetscLogInfo(0,"MatConvert_AIJ_AIJMUMPS:Using MUMPS for Cholesky factorization and solves.");
945f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
946f579278aSKris Buschelman   *newmat = B;
947f579278aSKris Buschelman   PetscFunctionReturn(0);
948f579278aSKris Buschelman }
949f579278aSKris Buschelman EXTERN_C_END
950f579278aSKris Buschelman 
951f0c56d0fSKris Buschelman #undef __FUNCT__
952422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
953dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
954dfbe8321SBarry Smith   PetscErrorCode ierr;
9558e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
9568f340917SKris Buschelman 
957f0c56d0fSKris Buschelman   PetscFunctionBegin;
9588f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
9598e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
960f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
961f0c56d0fSKris Buschelman }
962f0c56d0fSKris Buschelman 
96324b6179bSKris Buschelman /*MC
964fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
96524b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
96624b6179bSKris Buschelman 
96724b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
96824b6179bSKris Buschelman   on how to declare the existence of external packages),
96924b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
97024b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
97124b6179bSKris Buschelman 
97224b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
97324b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
97424b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
97524b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
97628b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
97728b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
97828b08bd3SKris Buschelman   without data copy.
97924b6179bSKris Buschelman 
98024b6179bSKris Buschelman   Options Database Keys:
9810bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
98224b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
98324b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
98424b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
98524b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
98624b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
98724b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
98894b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
98924b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
99024b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
99124b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
99224b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
99324b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
99424b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
99524b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
99624b6179bSKris Buschelman 
99724b6179bSKris Buschelman   Level: beginner
99824b6179bSKris Buschelman 
99924b6179bSKris Buschelman .seealso: MATAIJMUMPS
100024b6179bSKris Buschelman M*/
100124b6179bSKris Buschelman 
1002397b6df1SKris Buschelman EXTERN_C_BEGIN
1003397b6df1SKris Buschelman #undef __FUNCT__
1004f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1005dfbe8321SBarry Smith PetscErrorCode MatCreate_SBAIJMUMPS(Mat A)
1006dfbe8321SBarry Smith {
10076849ba73SBarry Smith   PetscErrorCode ierr;
10086849ba73SBarry Smith   int size;
1009397b6df1SKris Buschelman 
1010397b6df1SKris Buschelman   PetscFunctionBegin;
10115441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
10125441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
10135441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
10145441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1015397b6df1SKris Buschelman   if (size == 1) {
1016397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1017397b6df1SKris Buschelman   } else {
1018397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1019397b6df1SKris Buschelman   }
1020*ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1021397b6df1SKris Buschelman   PetscFunctionReturn(0);
1022397b6df1SKris Buschelman }
1023397b6df1SKris Buschelman EXTERN_C_END
1024