xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision d1390fe99e56d011c149b7273cb450901cf54230)
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
51dfbe8321SBarry Smith PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,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;
86847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->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"
153dfbe8321SBarry Smith PetscErrorCode MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) {
154dfbe8321SBarry Smith   PetscErrorCode ierr;
155c338a77dSKris Buschelman   Mat       B=*newmat;
156f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
157901853e0SKris Buschelman   void      (*f)(void);
158c338a77dSKris Buschelman 
159c338a77dSKris Buschelman   PetscFunctionBegin;
160c338a77dSKris Buschelman   if (B != A) {
161c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
162c338a77dSKris Buschelman   }
163f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
164f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
165f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
166f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
167f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
168f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
169901853e0SKris Buschelman 
170901853e0SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
171901853e0SKris Buschelman   if (f) {
172*d1390fe9SSatish Balay     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(FCNVOID)mumps->MatPreallocate);CHKERRQ(ierr);
173901853e0SKris Buschelman   }
174f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
175901853e0SKris Buschelman 
176901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
177901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
178901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
179901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
180901853e0SKris Buschelman   ierr = PetscObjectComposeFunction(PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
181901853e0SKris Buschelman   ierr = PetscObjectComposeFunction(PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
182901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
183901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
184901853e0SKris Buschelman 
185901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
186c338a77dSKris Buschelman   *newmat = B;
187c338a77dSKris Buschelman   PetscFunctionReturn(0);
188c338a77dSKris Buschelman }
189c338a77dSKris Buschelman EXTERN_C_END
190c338a77dSKris Buschelman 
191397b6df1SKris Buschelman #undef __FUNCT__
1923924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
193dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
194dfbe8321SBarry Smith {
195f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
196dfbe8321SBarry Smith   PetscErrorCode ierr;
197dfbe8321SBarry Smith   int       size=lu->size;
1986849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
199397b6df1SKris Buschelman   PetscFunctionBegin;
200397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
201397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
202397b6df1SKris Buschelman     lu->id.job=JOB_END;
203397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
204397b6df1SKris Buschelman     zmumps_c(&lu->id);
205397b6df1SKris Buschelman #else
206397b6df1SKris Buschelman     dmumps_c(&lu->id);
207397b6df1SKris Buschelman #endif
208c338a77dSKris Buschelman     if (lu->irn) {
209c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
210c338a77dSKris Buschelman     }
211c338a77dSKris Buschelman     if (lu->jcn) {
212c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
213c338a77dSKris Buschelman     }
214c338a77dSKris Buschelman     if (size>1 && lu->val) {
215c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
216c338a77dSKris Buschelman     }
217397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
218397b6df1SKris Buschelman   }
219a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
220a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
221c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
222397b6df1SKris Buschelman   PetscFunctionReturn(0);
223397b6df1SKris Buschelman }
224397b6df1SKris Buschelman 
225397b6df1SKris Buschelman #undef __FUNCT__
226a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
227dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
228dfbe8321SBarry Smith {
2296849ba73SBarry Smith   PetscErrorCode ierr;
2306849ba73SBarry Smith   int  size;
231a39386dcSKris Buschelman 
232a39386dcSKris Buschelman   PetscFunctionBegin;
233a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
234a39386dcSKris Buschelman   if (size==1) {
235a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr);
236a39386dcSKris Buschelman   } else {
237a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr);
238a39386dcSKris Buschelman   }
239a39386dcSKris Buschelman   PetscFunctionReturn(0);
240a39386dcSKris Buschelman }
241a39386dcSKris Buschelman 
242a39386dcSKris Buschelman #undef __FUNCT__
243a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
244dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
245dfbe8321SBarry Smith {
2466849ba73SBarry Smith   PetscErrorCode ierr;
2476849ba73SBarry Smith   int  size;
248a39386dcSKris Buschelman 
249a39386dcSKris Buschelman   PetscFunctionBegin;
250a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
251a39386dcSKris Buschelman   if (size==1) {
252a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr);
253a39386dcSKris Buschelman   } else {
254a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr);
255a39386dcSKris Buschelman   }
256a39386dcSKris Buschelman   PetscFunctionReturn(0);
257a39386dcSKris Buschelman }
258a39386dcSKris Buschelman 
259a39386dcSKris Buschelman #undef __FUNCT__
260c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
261dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
262f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
263dfbe8321SBarry Smith   PetscErrorCode ierr;
264397b6df1SKris Buschelman 
265397b6df1SKris Buschelman   PetscFunctionBegin;
266c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
267c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
268c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
269c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
270c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
271c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
272c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
273c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
274c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
275c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
276958c9bccSBarry Smith   if (!lu->myid && lu->id.ICNTL(11)>0) {
277c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
278c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
279c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
280c338a77dSKris 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);
281c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
282c338a77dSKris 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);
283c338a77dSKris Buschelman 
284c338a77dSKris Buschelman   }
285c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
286c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
287adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
288c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
289c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
290c338a77dSKris Buschelman 
291c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
292c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
293c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
29457f0c58bSHong Zhang 
29557f0c58bSHong Zhang   /* infomation local to each processor */
296958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
29757f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
29857f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
299958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
30057f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
30157f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
302958c9bccSBarry Smith   if (!lu->myid) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
30357f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
30457f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
305adc1d99fSHong Zhang 
306958c9bccSBarry Smith   if (!lu->myid){ /* information from the host */
307adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
308adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
309adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
310adc1d99fSHong Zhang 
311adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
312adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
313adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
314adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
315adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
316adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
317adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
318adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
319adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
320adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
321adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
322adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
323adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
324adc1d99fSHong 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);
325adc1d99fSHong 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);
326adc1d99fSHong 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);
327adc1d99fSHong 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);
328adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
329adc1d99fSHong Zhang   }
330adc1d99fSHong Zhang 
331397b6df1SKris Buschelman   PetscFunctionReturn(0);
332397b6df1SKris Buschelman }
333397b6df1SKris Buschelman 
334397b6df1SKris Buschelman #undef __FUNCT__
335f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
336dfbe8321SBarry Smith PetscErrorCode MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
337dfbe8321SBarry Smith   PetscErrorCode ierr;
33832077d6dSBarry Smith   PetscTruth        iascii;
339397b6df1SKris Buschelman   PetscViewerFormat format;
340f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
341397b6df1SKris Buschelman 
342397b6df1SKris Buschelman   PetscFunctionBegin;
343397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
344397b6df1SKris Buschelman 
34532077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
34632077d6dSBarry Smith   if (iascii) {
347397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
348397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
349397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
350397b6df1SKris Buschelman     }
351397b6df1SKris Buschelman   }
352397b6df1SKris Buschelman   PetscFunctionReturn(0);
353397b6df1SKris Buschelman }
354397b6df1SKris Buschelman 
355397b6df1SKris Buschelman #undef __FUNCT__
356f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
357dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
358f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
359d54de34fSKris Buschelman   PetscScalar *array;
360397b6df1SKris Buschelman   Vec         x_seq;
361397b6df1SKris Buschelman   IS          iden;
362397b6df1SKris Buschelman   VecScatter  scat;
363dfbe8321SBarry Smith   PetscErrorCode ierr;
364397b6df1SKris Buschelman 
365397b6df1SKris Buschelman   PetscFunctionBegin;
366397b6df1SKris Buschelman   if (lu->size > 1){
367397b6df1SKris Buschelman     if (!lu->myid){
368397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
369397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
370397b6df1SKris Buschelman     } else {
371397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
372397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
373397b6df1SKris Buschelman     }
374397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
375397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
376397b6df1SKris Buschelman 
377397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
378397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
379397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
380397b6df1SKris Buschelman   } else {  /* size == 1 */
381397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
382397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
383397b6df1SKris Buschelman   }
384397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
385397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
386397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
387397b6df1SKris Buschelman #else
388397b6df1SKris Buschelman     lu->id.rhs = array;
389397b6df1SKris Buschelman #endif
390397b6df1SKris Buschelman   }
391397b6df1SKris Buschelman 
392397b6df1SKris Buschelman   /* solve phase */
393397b6df1SKris Buschelman   lu->id.job=3;
394397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
395397b6df1SKris Buschelman   zmumps_c(&lu->id);
396397b6df1SKris Buschelman #else
397397b6df1SKris Buschelman   dmumps_c(&lu->id);
398397b6df1SKris Buschelman #endif
399397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
400397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
401397b6df1SKris Buschelman   }
402397b6df1SKris Buschelman 
403397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
404397b6df1SKris Buschelman   if (lu->size > 1) {
405397b6df1SKris Buschelman     if (!lu->myid){
406397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
407397b6df1SKris Buschelman     }
408397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
409397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
410397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
411397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
412397b6df1SKris Buschelman   } else {
413397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
414397b6df1SKris Buschelman   }
415397b6df1SKris Buschelman 
416397b6df1SKris Buschelman   PetscFunctionReturn(0);
417397b6df1SKris Buschelman }
418397b6df1SKris Buschelman 
419a58c3f20SHong Zhang /*
420a58c3f20SHong Zhang   input:
421a58c3f20SHong Zhang    F:        numeric factor
422a58c3f20SHong Zhang   output:
423a58c3f20SHong Zhang    nneg:     total number of negative pivots
424a58c3f20SHong Zhang    nzero:    0
425a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
426a58c3f20SHong Zhang */
427a58c3f20SHong Zhang 
428a58c3f20SHong Zhang #undef __FUNCT__
429a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
430dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
431a58c3f20SHong Zhang {
432a58c3f20SHong Zhang   Mat_MUMPS  *lu =(Mat_MUMPS*)F->spptr;
433dfbe8321SBarry Smith   PetscErrorCode ierr;
434dfbe8321SBarry Smith   int        size;
435a58c3f20SHong Zhang 
436a58c3f20SHong Zhang   PetscFunctionBegin;
437bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
438bcb30aebSHong 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 */
439bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
440bcb30aebSHong Zhang     SETERRQ1(1,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
441bcb30aebSHong Zhang   }
442a58c3f20SHong Zhang   if (nneg){
443a58c3f20SHong Zhang     if (!lu->myid){
444a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
445a58c3f20SHong Zhang     }
446bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
447a58c3f20SHong Zhang   }
448a58c3f20SHong Zhang   if (nzero) *nzero = 0;
449a58c3f20SHong Zhang   if (npos)  *npos  = F->M - (*nneg);
450a58c3f20SHong Zhang   PetscFunctionReturn(0);
451a58c3f20SHong Zhang }
452a58c3f20SHong Zhang 
453397b6df1SKris Buschelman #undef __FUNCT__
454f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
455dfbe8321SBarry Smith PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
456f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
457f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
4586849ba73SBarry Smith   PetscErrorCode ierr;
4596849ba73SBarry Smith   int        rnz,nnz,nz,i,M=A->M,*ai,*aj,icntl;
460397b6df1SKris Buschelman   PetscTruth valOnly,flg;
461397b6df1SKris Buschelman 
462397b6df1SKris Buschelman   PetscFunctionBegin;
463397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
464f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
465397b6df1SKris Buschelman 
466397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
467397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
468397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
46975747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
470397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
471397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
472397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
473397b6df1SKris Buschelman 
474397b6df1SKris Buschelman     /* Set mumps options */
475397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
476397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
477397b6df1SKris Buschelman     lu->id.sym=lu->sym;
478397b6df1SKris Buschelman     if (lu->sym == 2){
479397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
480397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
481397b6df1SKris Buschelman     }
482397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
483397b6df1SKris Buschelman   zmumps_c(&lu->id);
484397b6df1SKris Buschelman #else
485397b6df1SKris Buschelman   dmumps_c(&lu->id);
486397b6df1SKris Buschelman #endif
487397b6df1SKris Buschelman 
488397b6df1SKris Buschelman     if (lu->size == 1){
489397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
490397b6df1SKris Buschelman     } else {
491397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
492397b6df1SKris Buschelman     }
493397b6df1SKris Buschelman 
494397b6df1SKris Buschelman     icntl=-1;
495397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
496397b6df1SKris Buschelman     if (flg && icntl > 0) {
497397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
498397b6df1SKris Buschelman     } else { /* no output */
499397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
500397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
501397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
502397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
503397b6df1SKris Buschelman     }
504397b6df1SKris 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);
505397b6df1SKris Buschelman     icntl=-1;
506397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
507397b6df1SKris Buschelman     if (flg) {
508397b6df1SKris Buschelman       if (icntl== 1){
509397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
510397b6df1SKris Buschelman       } else {
511397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
512397b6df1SKris Buschelman       }
513397b6df1SKris Buschelman     }
514397b6df1SKris 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);
515397b6df1SKris 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);
51694b7f48cSBarry 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);
517397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
518397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
519adc1d99fSHong 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);
520397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
521397b6df1SKris Buschelman 
522397b6df1SKris Buschelman     /*
523397b6df1SKris 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);
524397b6df1SKris Buschelman     if (flg){
525397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
526397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
527397b6df1SKris Buschelman       } else {
528397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
529397b6df1SKris Buschelman       }
530397b6df1SKris Buschelman     }
531397b6df1SKris Buschelman     */
532397b6df1SKris Buschelman 
533397b6df1SKris 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);
534397b6df1SKris 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);
535397b6df1SKris 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);
536397b6df1SKris Buschelman     PetscOptionsEnd();
537397b6df1SKris Buschelman   }
538397b6df1SKris Buschelman 
539397b6df1SKris Buschelman   /* define matrix A */
540397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
541397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
542397b6df1SKris Buschelman     if (!lu->myid) {
543c36ead0aSKris Buschelman       if (lua->isAIJ){
544397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
545397b6df1SKris Buschelman         nz               = aa->nz;
546397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
547397b6df1SKris Buschelman       } else {
548397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5496c6c5352SBarry Smith         nz                  =  aa->nz;
550397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
551397b6df1SKris Buschelman       }
552397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
553397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
554397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
555397b6df1SKris Buschelman         nz = 0;
556397b6df1SKris Buschelman         for (i=0; i<M; i++){
557397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
558397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
559397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
560397b6df1SKris Buschelman           }
561397b6df1SKris Buschelman         }
562397b6df1SKris Buschelman       }
563397b6df1SKris Buschelman     }
564397b6df1SKris Buschelman     break;
565397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
566397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
567397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
568397b6df1SKris Buschelman     } else {
569397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
570397b6df1SKris Buschelman     }
571397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
572397b6df1SKris Buschelman     break;
573397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
574397b6df1SKris Buschelman   }
575397b6df1SKris Buschelman 
576397b6df1SKris Buschelman   /* analysis phase */
577397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
578397b6df1SKris Buschelman      lu->id.n = M;
579397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
580397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
581397b6df1SKris Buschelman       if (!lu->myid) {
582397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
583397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
584397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
585397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
586397b6df1SKris Buschelman #else
587397b6df1SKris Buschelman           lu->id.a = lu->val;
588397b6df1SKris Buschelman #endif
589397b6df1SKris Buschelman         }
590397b6df1SKris Buschelman       }
591397b6df1SKris Buschelman       break;
592397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
593397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
594397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
595397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
596397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
597397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
598397b6df1SKris Buschelman #else
599397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
600397b6df1SKris Buschelman #endif
601397b6df1SKris Buschelman       }
602397b6df1SKris Buschelman       break;
603397b6df1SKris Buschelman     }
604397b6df1SKris Buschelman     lu->id.job=1;
605397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
606397b6df1SKris Buschelman   zmumps_c(&lu->id);
607397b6df1SKris Buschelman #else
608397b6df1SKris Buschelman   dmumps_c(&lu->id);
609397b6df1SKris Buschelman #endif
610397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
611397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
612397b6df1SKris Buschelman     }
613397b6df1SKris Buschelman   }
614397b6df1SKris Buschelman 
615397b6df1SKris Buschelman   /* numerical factorization phase */
616958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
617a7aca84bSHong Zhang     if (!lu->myid) {
618397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
619397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
620397b6df1SKris Buschelman #else
621397b6df1SKris Buschelman       lu->id.a = lu->val;
622397b6df1SKris Buschelman #endif
623397b6df1SKris Buschelman     }
624397b6df1SKris Buschelman   } else {
625397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
626397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
627397b6df1SKris Buschelman #else
628397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
629397b6df1SKris Buschelman #endif
630397b6df1SKris Buschelman   }
631397b6df1SKris Buschelman   lu->id.job=2;
632397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
633397b6df1SKris Buschelman   zmumps_c(&lu->id);
634397b6df1SKris Buschelman #else
635397b6df1SKris Buschelman   dmumps_c(&lu->id);
636397b6df1SKris Buschelman #endif
637397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
638a7aca84bSHong Zhang     SETERRQ2(1,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
639397b6df1SKris Buschelman   }
640397b6df1SKris Buschelman 
641397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
642397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
643397b6df1SKris Buschelman   }
644397b6df1SKris Buschelman 
645397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
646397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
647ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
648397b6df1SKris Buschelman   PetscFunctionReturn(0);
649397b6df1SKris Buschelman }
650397b6df1SKris Buschelman 
651397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
652397b6df1SKris Buschelman #undef __FUNCT__
653f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
654dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
655397b6df1SKris Buschelman   Mat       B;
656f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
657dfbe8321SBarry Smith   PetscErrorCode ierr;
658397b6df1SKris Buschelman 
659397b6df1SKris Buschelman   PetscFunctionBegin;
660397b6df1SKris Buschelman 
661397b6df1SKris Buschelman   /* Create the factorization matrix */
662397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
663be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
664397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
665397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
666397b6df1SKris Buschelman 
667f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
668397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
669f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
670397b6df1SKris Buschelman   lu->sym                 = 0;
671397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
672397b6df1SKris Buschelman 
673397b6df1SKris Buschelman   *F = B;
674397b6df1SKris Buschelman   PetscFunctionReturn(0);
675397b6df1SKris Buschelman }
676397b6df1SKris Buschelman 
677397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
678397b6df1SKris Buschelman #undef __FUNCT__
679f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
680dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
681397b6df1SKris Buschelman   Mat       B;
682f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
683dfbe8321SBarry Smith   PetscErrorCode ierr;
684397b6df1SKris Buschelman 
685397b6df1SKris Buschelman   PetscFunctionBegin;
686397b6df1SKris Buschelman 
687397b6df1SKris Buschelman   /* Create the factorization matrix */
688397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
689be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
690efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
691efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
692397b6df1SKris Buschelman 
693f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
694a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
695397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
696f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
697397b6df1SKris Buschelman   lu->sym                       = 2;
698397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
699397b6df1SKris Buschelman 
700397b6df1SKris Buschelman   *F = B;
701397b6df1SKris Buschelman   PetscFunctionReturn(0);
702397b6df1SKris Buschelman }
703397b6df1SKris Buschelman 
704397b6df1SKris Buschelman #undef __FUNCT__
705f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
706dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
707dfbe8321SBarry Smith   PetscErrorCode ierr;
708f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
709c338a77dSKris Buschelman 
710397b6df1SKris Buschelman   PetscFunctionBegin;
711c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
712f0c56d0fSKris Buschelman 
713c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
714c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
715f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
716397b6df1SKris Buschelman   PetscFunctionReturn(0);
717397b6df1SKris Buschelman }
718397b6df1SKris Buschelman 
719c338a77dSKris Buschelman EXTERN_C_BEGIN
720c338a77dSKris Buschelman #undef __FUNCT__
721f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
722dfbe8321SBarry Smith PetscErrorCode MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat)\
723dfbe8321SBarry Smith {
724dfbe8321SBarry Smith   PetscErrorCode ierr;
725dfbe8321SBarry Smith   int       size;
726c338a77dSKris Buschelman   MPI_Comm  comm;
727c338a77dSKris Buschelman   Mat       B=*newmat;
728f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
729397b6df1SKris Buschelman 
730397b6df1SKris Buschelman   PetscFunctionBegin;
731c338a77dSKris Buschelman   if (B != A) {
732c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
733397b6df1SKris Buschelman   }
734397b6df1SKris Buschelman 
735c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
736f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
737c338a77dSKris Buschelman 
738f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
739c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
740c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
741c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
742c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
743c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
744a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
745c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
746f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
747c338a77dSKris Buschelman 
7484b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
749422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
750f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
751f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
752f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
7533924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
754c338a77dSKris Buschelman 
755c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
756c338a77dSKris Buschelman   if (size == 1) {
757c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
758f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
759c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
760c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
761c338a77dSKris Buschelman   } else {
762c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
763f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
764c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
765c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
766c338a77dSKris Buschelman   }
767c338a77dSKris Buschelman 
768f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
769c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
770c338a77dSKris Buschelman   *newmat = B;
771397b6df1SKris Buschelman   PetscFunctionReturn(0);
772397b6df1SKris Buschelman }
773c338a77dSKris Buschelman EXTERN_C_END
774397b6df1SKris Buschelman 
77524b6179bSKris Buschelman /*MC
776fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
77724b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
77824b6179bSKris Buschelman 
77924b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
78024b6179bSKris Buschelman   on how to declare the existence of external packages),
78124b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
78224b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
78324b6179bSKris Buschelman   This matrix type is only supported for double precision real.
78424b6179bSKris Buschelman 
78524b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
78624b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
78724b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
78824b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
78928b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
79028b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
79128b08bd3SKris Buschelman   without data copy.
79224b6179bSKris Buschelman 
79324b6179bSKris Buschelman   Options Database Keys:
7940bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
79524b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
79624b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
79724b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
79824b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
79924b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
80024b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
80194b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
80224b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
80324b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
80424b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
80524b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
80624b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
80724b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
80824b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
80924b6179bSKris Buschelman 
81024b6179bSKris Buschelman   Level: beginner
81124b6179bSKris Buschelman 
81224b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
81324b6179bSKris Buschelman M*/
81424b6179bSKris Buschelman 
815397b6df1SKris Buschelman EXTERN_C_BEGIN
816397b6df1SKris Buschelman #undef __FUNCT__
817f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
818dfbe8321SBarry Smith PetscErrorCode MatCreate_AIJMUMPS(Mat A)
819dfbe8321SBarry Smith {
820dfbe8321SBarry Smith   PetscErrorCode ierr;
821dfbe8321SBarry Smith   int      size;
822e2d9671bSKris Buschelman   Mat      A_diag;
823397b6df1SKris Buschelman   MPI_Comm comm;
824397b6df1SKris Buschelman 
825397b6df1SKris Buschelman   PetscFunctionBegin;
8265441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
8275441df8eSKris Buschelman   /*   and AIJMUMPS types */
8285441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
829397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
830397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
831397b6df1SKris Buschelman   if (size == 1) {
832397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
833397b6df1SKris Buschelman   } else {
834397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
835e2d9671bSKris Buschelman     A_diag = ((Mat_MPIAIJ *)A->data)->A;
836e2d9671bSKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr);
837397b6df1SKris Buschelman   }
838f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
839397b6df1SKris Buschelman   PetscFunctionReturn(0);
840397b6df1SKris Buschelman }
841397b6df1SKris Buschelman EXTERN_C_END
842397b6df1SKris Buschelman 
843f579278aSKris Buschelman #undef __FUNCT__
844f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
845dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
846dfbe8321SBarry Smith {
847dfbe8321SBarry Smith   PetscErrorCode ierr;
848f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
849f579278aSKris Buschelman 
850f579278aSKris Buschelman   PetscFunctionBegin;
851f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
852f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
853f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
854f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
855f579278aSKris Buschelman   PetscFunctionReturn(0);
856f579278aSKris Buschelman }
857f579278aSKris Buschelman 
858f579278aSKris Buschelman EXTERN_C_BEGIN
859f579278aSKris Buschelman #undef __FUNCT__
8609c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
861dfbe8321SBarry Smith PetscErrorCode MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat  B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz)
8629c097c71SKris Buschelman {
8639c097c71SKris Buschelman   Mat       A;
86402217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
865dfbe8321SBarry Smith   PetscErrorCode ierr;
8669c097c71SKris Buschelman 
8679c097c71SKris Buschelman   PetscFunctionBegin;
8689c097c71SKris Buschelman   /*
8699c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
8709c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
8719c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
8729c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
8739c097c71SKris Buschelman     in the preallocation routine.
8749c097c71SKris Buschelman   */
8759c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
8769c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
8779c097c71SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
8789c097c71SKris Buschelman   PetscFunctionReturn(0);
8799c097c71SKris Buschelman }
8809c097c71SKris Buschelman EXTERN_C_END
8819c097c71SKris Buschelman 
8829c097c71SKris Buschelman EXTERN_C_BEGIN
8839c097c71SKris Buschelman #undef __FUNCT__
884f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
885dfbe8321SBarry Smith PetscErrorCode MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat)
886dfbe8321SBarry Smith {
887dfbe8321SBarry Smith   PetscErrorCode ierr;
888dfbe8321SBarry Smith   int       size;
889f579278aSKris Buschelman   MPI_Comm  comm;
890f579278aSKris Buschelman   Mat       B=*newmat;
891422e82a1SHong Zhang   Mat_MUMPS *mumps;
8929c097c71SKris Buschelman   void      (*f)(void);
893f579278aSKris Buschelman 
894f579278aSKris Buschelman   PetscFunctionBegin;
895f579278aSKris Buschelman   if (B != A) {
896f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
897f579278aSKris Buschelman   }
898f579278aSKris Buschelman 
899f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
900f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
901f579278aSKris Buschelman 
902f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
903f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
904f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
905f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
906f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
907f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
908a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
909f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
910f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
911f579278aSKris Buschelman 
912f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
913422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
914f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
915f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
916f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
9173924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
918f579278aSKris Buschelman 
919f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
920f579278aSKris Buschelman   if (size == 1) {
921f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
922f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
923f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
924f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
925f579278aSKris Buschelman   } else {
9269c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
9279c097c71SKris Buschelman     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
928901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
9296849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
9309c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
9319c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
9329c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
9339c097c71SKris Buschelman     }
934f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
935f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
936f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
937f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
938f579278aSKris Buschelman   }
939f579278aSKris Buschelman 
940f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
941f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
942f579278aSKris Buschelman   *newmat = B;
943f579278aSKris Buschelman   PetscFunctionReturn(0);
944f579278aSKris Buschelman }
945f579278aSKris Buschelman EXTERN_C_END
946f579278aSKris Buschelman 
947f0c56d0fSKris Buschelman #undef __FUNCT__
948422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
949dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
950dfbe8321SBarry Smith   PetscErrorCode ierr;
9518e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
9528f340917SKris Buschelman 
953f0c56d0fSKris Buschelman   PetscFunctionBegin;
9548f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
9558e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
956f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
957f0c56d0fSKris Buschelman }
958f0c56d0fSKris Buschelman 
95924b6179bSKris Buschelman /*MC
960fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
96124b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
96224b6179bSKris Buschelman 
96324b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
96424b6179bSKris Buschelman   on how to declare the existence of external packages),
96524b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
96624b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
96724b6179bSKris Buschelman   This matrix type is only supported for double precision real.
96824b6179bSKris Buschelman 
96924b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
97024b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
97124b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
97224b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
97328b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
97428b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
97528b08bd3SKris Buschelman   without data copy.
97624b6179bSKris Buschelman 
97724b6179bSKris Buschelman   Options Database Keys:
9780bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
97924b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
98024b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
98124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
98224b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
98324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
98424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
98594b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
98624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
98724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
98824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
98924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
99024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
99124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
99224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
99324b6179bSKris Buschelman 
99424b6179bSKris Buschelman   Level: beginner
99524b6179bSKris Buschelman 
99624b6179bSKris Buschelman .seealso: MATAIJMUMPS
99724b6179bSKris Buschelman M*/
99824b6179bSKris Buschelman 
999397b6df1SKris Buschelman EXTERN_C_BEGIN
1000397b6df1SKris Buschelman #undef __FUNCT__
1001f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1002dfbe8321SBarry Smith PetscErrorCode MatCreate_SBAIJMUMPS(Mat A)
1003dfbe8321SBarry Smith {
10046849ba73SBarry Smith   PetscErrorCode ierr;
10056849ba73SBarry Smith   int size;
1006397b6df1SKris Buschelman 
1007397b6df1SKris Buschelman   PetscFunctionBegin;
10085441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
10095441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
10105441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
10115441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1012397b6df1SKris Buschelman   if (size == 1) {
1013397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1014397b6df1SKris Buschelman   } else {
1015397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1016397b6df1SKris Buschelman   }
1017f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
1018397b6df1SKris Buschelman   PetscFunctionReturn(0);
1019397b6df1SKris Buschelman }
1020397b6df1SKris Buschelman EXTERN_C_END
1021