xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision a39386dc5f1e6b8d7f614f7d32f7995628ebc4be)
1397b6df1SKris Buschelman /*$Id: mumps.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2397b6df1SKris Buschelman /*
33a0aa19bSHong Zhang     Provides an interface to the MUMPS_4.3 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]
23397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
24397b6df1SKris Buschelman 
25397b6df1SKris Buschelman typedef struct {
26397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
27397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
28397b6df1SKris Buschelman #else
29397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
30397b6df1SKris Buschelman #endif
31397b6df1SKris Buschelman   MatStructure   matstruc;
32397b6df1SKris Buschelman   int            myid,size,*irn,*jcn,sym;
33397b6df1SKris Buschelman   PetscScalar    *val;
34397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
35397b6df1SKris Buschelman 
36c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
37f0c56d0fSKris Buschelman   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
38c338a77dSKris Buschelman   int (*MatView)(Mat,PetscViewer);
39c338a77dSKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
40c338a77dSKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
41c338a77dSKris Buschelman   int (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
42c338a77dSKris Buschelman   int (*MatDestroy)(Mat);
43*a39386dcSKris Buschelman   int (*specialdestroy)(Mat);
44f0c56d0fSKris Buschelman } Mat_MUMPS;
45f0c56d0fSKris Buschelman 
46f0c56d0fSKris Buschelman EXTERN int MatDuplicate_AIJMUMPS(Mat,MatDuplicateOption,Mat*);
47f0c56d0fSKris Buschelman EXTERN int MatDuplicate_SBAIJMUMPS(Mat,MatDuplicateOption,Mat*);
480e3434eeSKris Buschelman 
49397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
50397b6df1SKris Buschelman /*
51397b6df1SKris Buschelman   input:
5275747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
53397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
54397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
55397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
56397b6df1SKris Buschelman   output:
57397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
58397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
59397b6df1SKris Buschelman  */
60f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
61397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
62397b6df1SKris Buschelman   int         ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
63d54de34fSKris Buschelman   int         *row,*col;
64397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
65f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
66397b6df1SKris Buschelman 
67397b6df1SKris Buschelman   PetscFunctionBegin;
68397b6df1SKris Buschelman   if (mumps->isAIJ){
69397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
70397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
71397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
72397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
73397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
74397b6df1SKris Buschelman     garray = mat->garray;
75397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
76397b6df1SKris Buschelman 
77397b6df1SKris Buschelman   } else {
78397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
79397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
80397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
81847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs);
82397b6df1SKris Buschelman     nz = aa->s_nz + bb->nz;
83397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
84397b6df1SKris Buschelman     garray = mat->garray;
85397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
86397b6df1SKris Buschelman   }
87397b6df1SKris Buschelman 
88397b6df1SKris Buschelman   if (!valOnly){
89397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
90397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
91397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
92397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
93397b6df1SKris Buschelman   } else {
94397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
95397b6df1SKris Buschelman   }
96397b6df1SKris Buschelman   *nnz = nz;
97397b6df1SKris Buschelman 
98028e57e8SHong Zhang   jj = 0; irow = rstart;
99397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
100397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
101397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
102397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
103397b6df1SKris Buschelman     bjj = bj + bi[i];
104397b6df1SKris Buschelman 
105397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
106397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
10775747be1SHong Zhang     j=-1;
10875747be1SHong Zhang     do {
10975747be1SHong Zhang       j++;
11075747be1SHong Zhang       if (j == countB) break;
111397b6df1SKris Buschelman       jcol = garray[bjj[j]];
11275747be1SHong Zhang     } while (jcol < colA_start);
11375747be1SHong Zhang     jB = j;
114397b6df1SKris Buschelman 
115397b6df1SKris Buschelman     /* B-part, smaller col index */
116397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
117397b6df1SKris Buschelman     for (j=0; j<jB; j++){
118397b6df1SKris Buschelman       jcol = garray[bjj[j]];
119397b6df1SKris Buschelman       if (!valOnly){
120397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12175747be1SHong Zhang 
122397b6df1SKris Buschelman       }
123397b6df1SKris Buschelman       val[jj++] = *bv++;
124397b6df1SKris Buschelman     }
125397b6df1SKris Buschelman     /* A-part */
126397b6df1SKris Buschelman     for (j=0; j<countA; j++){
127397b6df1SKris Buschelman       if (!valOnly){
128397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
129397b6df1SKris Buschelman       }
130397b6df1SKris Buschelman       val[jj++] = *av++;
131397b6df1SKris Buschelman     }
132397b6df1SKris Buschelman     /* B-part, larger col index */
133397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
134397b6df1SKris Buschelman       if (!valOnly){
135397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
136397b6df1SKris Buschelman       }
137397b6df1SKris Buschelman       val[jj++] = *bv++;
138397b6df1SKris Buschelman     }
139397b6df1SKris Buschelman     irow++;
140397b6df1SKris Buschelman   }
141397b6df1SKris Buschelman 
142397b6df1SKris Buschelman   PetscFunctionReturn(0);
143397b6df1SKris Buschelman }
144397b6df1SKris Buschelman 
145c338a77dSKris Buschelman EXTERN_C_BEGIN
146c338a77dSKris Buschelman #undef __FUNCT__
147c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
148c338a77dSKris Buschelman int MatConvert_MUMPS_Base(Mat A,MatType type,Mat *newmat) {
149c338a77dSKris Buschelman   int       ierr;
150c338a77dSKris Buschelman   Mat       B=*newmat;
151f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
152c338a77dSKris Buschelman 
153c338a77dSKris Buschelman   PetscFunctionBegin;
154c338a77dSKris Buschelman   if (B != A) {
155c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
156c338a77dSKris Buschelman   }
157f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
158f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
159f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
160f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
161f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
162f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
1633924e44cSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
164f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
165c338a77dSKris Buschelman   *newmat = B;
166c338a77dSKris Buschelman   PetscFunctionReturn(0);
167c338a77dSKris Buschelman }
168c338a77dSKris Buschelman EXTERN_C_END
169c338a77dSKris Buschelman 
170397b6df1SKris Buschelman #undef __FUNCT__
1713924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
1723924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) {
173f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
174c338a77dSKris Buschelman   int       ierr,size=lu->size;
175*a39386dcSKris Buschelman   int       (*specialdestroy)(Mat);
176397b6df1SKris Buschelman   PetscFunctionBegin;
177397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
178397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
179397b6df1SKris Buschelman     lu->id.job=JOB_END;
180397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
181397b6df1SKris Buschelman     zmumps_c(&lu->id);
182397b6df1SKris Buschelman #else
183397b6df1SKris Buschelman     dmumps_c(&lu->id);
184397b6df1SKris Buschelman #endif
185c338a77dSKris Buschelman     if (lu->irn) {
186c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
187c338a77dSKris Buschelman     }
188c338a77dSKris Buschelman     if (lu->jcn) {
189c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
190c338a77dSKris Buschelman     }
191c338a77dSKris Buschelman     if (size>1 && lu->val) {
192c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
193c338a77dSKris Buschelman     }
194397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
195397b6df1SKris Buschelman   }
196*a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
197*a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
198c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
199397b6df1SKris Buschelman   PetscFunctionReturn(0);
200397b6df1SKris Buschelman }
201397b6df1SKris Buschelman 
202397b6df1SKris Buschelman #undef __FUNCT__
203*a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
204*a39386dcSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) {
205*a39386dcSKris Buschelman   int ierr, size;
206*a39386dcSKris Buschelman 
207*a39386dcSKris Buschelman   PetscFunctionBegin;
208*a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
209*a39386dcSKris Buschelman   if (size==1) {
210*a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr);
211*a39386dcSKris Buschelman   } else {
212*a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr);
213*a39386dcSKris Buschelman   }
214*a39386dcSKris Buschelman   PetscFunctionReturn(0);
215*a39386dcSKris Buschelman }
216*a39386dcSKris Buschelman 
217*a39386dcSKris Buschelman #undef __FUNCT__
218*a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
219*a39386dcSKris Buschelman int MatDestroy_SBAIJMUMPS(Mat A) {
220*a39386dcSKris Buschelman   int ierr, size;
221*a39386dcSKris Buschelman 
222*a39386dcSKris Buschelman   PetscFunctionBegin;
223*a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
224*a39386dcSKris Buschelman   if (size==1) {
225*a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr);
226*a39386dcSKris Buschelman   } else {
227*a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr);
228*a39386dcSKris Buschelman   }
229*a39386dcSKris Buschelman   PetscFunctionReturn(0);
230*a39386dcSKris Buschelman }
231*a39386dcSKris Buschelman 
232*a39386dcSKris Buschelman #undef __FUNCT__
233c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
234f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
235f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
236397b6df1SKris Buschelman   int       ierr;
237397b6df1SKris Buschelman 
238397b6df1SKris Buschelman   PetscFunctionBegin;
239c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
240c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
241c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
242c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
243c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
244c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
245c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
246c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
247c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
248c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
249c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
250c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
251c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
252c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
253c338a77dSKris 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);
254c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
255c338a77dSKris 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);
256c338a77dSKris Buschelman 
257c338a77dSKris Buschelman   }
258c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):     %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
259c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):     %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
260c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (efficiency control):     %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
261c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):     %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
262c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):       %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
263c338a77dSKris Buschelman 
264c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
265c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
266c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
267397b6df1SKris Buschelman   PetscFunctionReturn(0);
268397b6df1SKris Buschelman }
269397b6df1SKris Buschelman 
270397b6df1SKris Buschelman #undef __FUNCT__
271f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
272f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
273397b6df1SKris Buschelman   int               ierr;
274397b6df1SKris Buschelman   PetscTruth        isascii;
275397b6df1SKris Buschelman   PetscViewerFormat format;
276f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
277397b6df1SKris Buschelman 
278397b6df1SKris Buschelman   PetscFunctionBegin;
279397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
280397b6df1SKris Buschelman 
281397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
282397b6df1SKris Buschelman   if (isascii) {
283397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
284397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
285397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
286397b6df1SKris Buschelman     }
287397b6df1SKris Buschelman   }
288397b6df1SKris Buschelman   PetscFunctionReturn(0);
289397b6df1SKris Buschelman }
290397b6df1SKris Buschelman 
291397b6df1SKris Buschelman #undef __FUNCT__
292f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
293f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
294f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
295d54de34fSKris Buschelman   PetscScalar *array;
296397b6df1SKris Buschelman   Vec         x_seq;
297397b6df1SKris Buschelman   IS          iden;
298397b6df1SKris Buschelman   VecScatter  scat;
299397b6df1SKris Buschelman   int         ierr;
300397b6df1SKris Buschelman 
301397b6df1SKris Buschelman   PetscFunctionBegin;
302397b6df1SKris Buschelman   if (lu->size > 1){
303397b6df1SKris Buschelman     if (!lu->myid){
304397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
305397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
306397b6df1SKris Buschelman     } else {
307397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
308397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
309397b6df1SKris Buschelman     }
310397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
311397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
312397b6df1SKris Buschelman 
313397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
314397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
315397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
316397b6df1SKris Buschelman   } else {  /* size == 1 */
317397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
318397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
319397b6df1SKris Buschelman   }
320397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
321397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
322397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
323397b6df1SKris Buschelman #else
324397b6df1SKris Buschelman     lu->id.rhs = array;
325397b6df1SKris Buschelman #endif
326397b6df1SKris Buschelman   }
327397b6df1SKris Buschelman 
328397b6df1SKris Buschelman   /* solve phase */
329397b6df1SKris Buschelman   lu->id.job=3;
330397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
331397b6df1SKris Buschelman   zmumps_c(&lu->id);
332397b6df1SKris Buschelman #else
333397b6df1SKris Buschelman   dmumps_c(&lu->id);
334397b6df1SKris Buschelman #endif
335397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
336397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
337397b6df1SKris Buschelman   }
338397b6df1SKris Buschelman 
339397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
340397b6df1SKris Buschelman   if (lu->size > 1) {
341397b6df1SKris Buschelman     if (!lu->myid){
342397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
343397b6df1SKris Buschelman     }
344397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
345397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
346397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
347397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
348397b6df1SKris Buschelman   } else {
349397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
350397b6df1SKris Buschelman   }
351397b6df1SKris Buschelman 
352397b6df1SKris Buschelman   PetscFunctionReturn(0);
353397b6df1SKris Buschelman }
354397b6df1SKris Buschelman 
355397b6df1SKris Buschelman #undef __FUNCT__
356f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
357f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
358f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
359f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
360397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
361397b6df1SKris Buschelman   PetscTruth valOnly,flg;
362397b6df1SKris Buschelman 
363397b6df1SKris Buschelman   PetscFunctionBegin;
364397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
365f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
366397b6df1SKris Buschelman 
367397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
368397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
369397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
37075747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
371397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
372397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
373397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
374397b6df1SKris Buschelman 
375397b6df1SKris Buschelman     /* Set mumps options */
376397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
377397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
378397b6df1SKris Buschelman     lu->id.sym=lu->sym;
379397b6df1SKris Buschelman     if (lu->sym == 2){
380397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
381397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
382397b6df1SKris Buschelman     }
383397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
384397b6df1SKris Buschelman   zmumps_c(&lu->id);
385397b6df1SKris Buschelman #else
386397b6df1SKris Buschelman   dmumps_c(&lu->id);
387397b6df1SKris Buschelman #endif
388397b6df1SKris Buschelman 
389397b6df1SKris Buschelman     if (lu->size == 1){
390397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
391397b6df1SKris Buschelman     } else {
392397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
393397b6df1SKris Buschelman     }
394397b6df1SKris Buschelman 
395397b6df1SKris Buschelman     icntl=-1;
396397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
397397b6df1SKris Buschelman     if (flg && icntl > 0) {
398397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
399397b6df1SKris Buschelman     } else { /* no output */
400397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
401397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
402397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
403397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
404397b6df1SKris Buschelman     }
405397b6df1SKris 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);
406397b6df1SKris Buschelman     icntl=-1;
407397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
408397b6df1SKris Buschelman     if (flg) {
409397b6df1SKris Buschelman       if (icntl== 1){
410397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
411397b6df1SKris Buschelman       } else {
412397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
413397b6df1SKris Buschelman       }
414397b6df1SKris Buschelman     }
415397b6df1SKris 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);
416397b6df1SKris 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);
417397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -sles_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr);
418397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
419397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
420397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): efficiency control","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
421397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
422397b6df1SKris Buschelman 
423397b6df1SKris Buschelman     /*
424397b6df1SKris 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);
425397b6df1SKris Buschelman     if (flg){
426397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
427397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
428397b6df1SKris Buschelman       } else {
429397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
430397b6df1SKris Buschelman       }
431397b6df1SKris Buschelman     }
432397b6df1SKris Buschelman     */
433397b6df1SKris Buschelman 
434397b6df1SKris 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);
435397b6df1SKris 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);
436397b6df1SKris 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);
437397b6df1SKris Buschelman     PetscOptionsEnd();
438397b6df1SKris Buschelman   }
439397b6df1SKris Buschelman 
440397b6df1SKris Buschelman   /* define matrix A */
441397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
442397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
443397b6df1SKris Buschelman     if (!lu->myid) {
444c36ead0aSKris Buschelman       if (lua->isAIJ){
445397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
446397b6df1SKris Buschelman         nz               = aa->nz;
447397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
448397b6df1SKris Buschelman       } else {
449397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
450397b6df1SKris Buschelman         nz                  =  aa->s_nz;
451397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
452397b6df1SKris Buschelman       }
453397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
454397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
455397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
456397b6df1SKris Buschelman         nz = 0;
457397b6df1SKris Buschelman         for (i=0; i<M; i++){
458397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
459397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
460397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
461397b6df1SKris Buschelman           }
462397b6df1SKris Buschelman         }
463397b6df1SKris Buschelman       }
464397b6df1SKris Buschelman     }
465397b6df1SKris Buschelman     break;
466397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
467397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
468397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
469397b6df1SKris Buschelman     } else {
470397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
471397b6df1SKris Buschelman     }
472397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
473397b6df1SKris Buschelman     break;
474397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
475397b6df1SKris Buschelman   }
476397b6df1SKris Buschelman 
477397b6df1SKris Buschelman   /* analysis phase */
478397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
479397b6df1SKris Buschelman      lu->id.n = M;
480397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
481397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
482397b6df1SKris Buschelman       if (!lu->myid) {
483397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
484397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
485397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
486397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
487397b6df1SKris Buschelman #else
488397b6df1SKris Buschelman           lu->id.a = lu->val;
489397b6df1SKris Buschelman #endif
490397b6df1SKris Buschelman         }
491397b6df1SKris Buschelman       }
492397b6df1SKris Buschelman       break;
493397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
494397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
495397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
496397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
497397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
498397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
499397b6df1SKris Buschelman #else
500397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
501397b6df1SKris Buschelman #endif
502397b6df1SKris Buschelman       }
503397b6df1SKris Buschelman       break;
504397b6df1SKris Buschelman     }
505397b6df1SKris Buschelman     lu->id.job=1;
506397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
507397b6df1SKris Buschelman   zmumps_c(&lu->id);
508397b6df1SKris Buschelman #else
509397b6df1SKris Buschelman   dmumps_c(&lu->id);
510397b6df1SKris Buschelman #endif
511397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
512397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
513397b6df1SKris Buschelman     }
514397b6df1SKris Buschelman   }
515397b6df1SKris Buschelman 
516397b6df1SKris Buschelman   /* numerical factorization phase */
517397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
518397b6df1SKris Buschelman     if (lu->myid == 0) {
519397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
520397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
521397b6df1SKris Buschelman #else
522397b6df1SKris Buschelman       lu->id.a = lu->val;
523397b6df1SKris Buschelman #endif
524397b6df1SKris Buschelman     }
525397b6df1SKris Buschelman   } else {
526397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
527397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
528397b6df1SKris Buschelman #else
529397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
530397b6df1SKris Buschelman #endif
531397b6df1SKris Buschelman   }
532397b6df1SKris Buschelman   lu->id.job=2;
533397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
534397b6df1SKris Buschelman   zmumps_c(&lu->id);
535397b6df1SKris Buschelman #else
536397b6df1SKris Buschelman   dmumps_c(&lu->id);
537397b6df1SKris Buschelman #endif
538397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
539397b6df1SKris Buschelman     SETERRQ1(1,"1, Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d\n",lu->id.INFOG(1));
540397b6df1SKris Buschelman   }
541397b6df1SKris Buschelman 
542397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
543397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
544397b6df1SKris Buschelman   }
545397b6df1SKris Buschelman 
546397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
547397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
548ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
549397b6df1SKris Buschelman   PetscFunctionReturn(0);
550397b6df1SKris Buschelman }
551397b6df1SKris Buschelman 
552397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
553397b6df1SKris Buschelman #undef __FUNCT__
554f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
555f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
556397b6df1SKris Buschelman   Mat       B;
557f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
558397b6df1SKris Buschelman   int       ierr;
559397b6df1SKris Buschelman 
560397b6df1SKris Buschelman   PetscFunctionBegin;
561397b6df1SKris Buschelman 
562397b6df1SKris Buschelman   /* Create the factorization matrix */
563397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
564397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
565397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
566397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
567397b6df1SKris Buschelman 
568f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
569397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
570f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
571397b6df1SKris Buschelman   lu->sym                 = 0;
572397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
573397b6df1SKris Buschelman 
574397b6df1SKris Buschelman   *F = B;
575397b6df1SKris Buschelman   PetscFunctionReturn(0);
576397b6df1SKris Buschelman }
577397b6df1SKris Buschelman 
578397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
579397b6df1SKris Buschelman #undef __FUNCT__
580f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
581f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
582397b6df1SKris Buschelman   Mat       B;
583f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
584397b6df1SKris Buschelman   int       ierr;
585397b6df1SKris Buschelman 
586397b6df1SKris Buschelman   PetscFunctionBegin;
587397b6df1SKris Buschelman 
588397b6df1SKris Buschelman   /* Create the factorization matrix */
589397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
590397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
591397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
592397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
593397b6df1SKris Buschelman 
594f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
595397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
596f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
597397b6df1SKris Buschelman   lu->sym                       = 2;
598397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
599397b6df1SKris Buschelman 
600397b6df1SKris Buschelman   *F = B;
601397b6df1SKris Buschelman   PetscFunctionReturn(0);
602397b6df1SKris Buschelman }
603397b6df1SKris Buschelman 
604397b6df1SKris Buschelman #undef __FUNCT__
605f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
606f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
607c338a77dSKris Buschelman   int       ierr;
608f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
609c338a77dSKris Buschelman 
610397b6df1SKris Buschelman   PetscFunctionBegin;
611c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
612f0c56d0fSKris Buschelman 
613c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
614c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
615f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
616397b6df1SKris Buschelman   PetscFunctionReturn(0);
617397b6df1SKris Buschelman }
618397b6df1SKris Buschelman 
619c338a77dSKris Buschelman EXTERN_C_BEGIN
620c338a77dSKris Buschelman #undef __FUNCT__
621f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
622f0c56d0fSKris Buschelman int MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,Mat *newmat) {
623c338a77dSKris Buschelman   int       ierr,size;
624c338a77dSKris Buschelman   MPI_Comm  comm;
625c338a77dSKris Buschelman   Mat       B=*newmat;
626f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
627397b6df1SKris Buschelman 
628397b6df1SKris Buschelman   PetscFunctionBegin;
629c338a77dSKris Buschelman   if (B != A) {
630c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
631397b6df1SKris Buschelman   }
632397b6df1SKris Buschelman 
633c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
634f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
635c338a77dSKris Buschelman 
636f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
637c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
638c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
639c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
640c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
641c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
642*a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
643c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
644f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
645c338a77dSKris Buschelman 
6464b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
647f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_AIJMUMPS;
648f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
649f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
650f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
6513924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
652c338a77dSKris Buschelman 
653c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
654c338a77dSKris Buschelman   if (size == 1) {
655c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
656f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
657c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
658c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
659c338a77dSKris Buschelman   } else {
660c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
661f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
662c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
663c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
664c338a77dSKris Buschelman   }
665c338a77dSKris Buschelman 
666f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
667c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
668c338a77dSKris Buschelman   *newmat = B;
669397b6df1SKris Buschelman   PetscFunctionReturn(0);
670397b6df1SKris Buschelman }
671c338a77dSKris Buschelman EXTERN_C_END
672397b6df1SKris Buschelman 
673f0c56d0fSKris Buschelman #undef __FUNCT__
674f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS"
675f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
676f0c56d0fSKris Buschelman   int       ierr;
6778f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
6788f340917SKris Buschelman 
679f0c56d0fSKris Buschelman   PetscFunctionBegin;
6808f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
681f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr);
682*a39386dcSKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
683f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
684f0c56d0fSKris Buschelman }
685f0c56d0fSKris Buschelman 
68624b6179bSKris Buschelman /*MC
687fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
68824b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
68924b6179bSKris Buschelman 
69024b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
69124b6179bSKris Buschelman   on how to declare the existence of external packages),
69224b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
69324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
69424b6179bSKris Buschelman   This matrix type is only supported for double precision real.
69524b6179bSKris Buschelman 
69624b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
69724b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
69824b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
69924b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
70028b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
70128b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
70228b08bd3SKris Buschelman   without data copy.
70324b6179bSKris Buschelman 
70424b6179bSKris Buschelman   Options Database Keys:
7050bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
70624b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
70724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
70824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
70924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
71024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
71124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
71224b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view
71324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
71424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
71524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
71624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
71724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
71824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
71924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
72024b6179bSKris Buschelman 
72124b6179bSKris Buschelman   Level: beginner
72224b6179bSKris Buschelman 
72324b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
72424b6179bSKris Buschelman M*/
72524b6179bSKris Buschelman 
726397b6df1SKris Buschelman EXTERN_C_BEGIN
727397b6df1SKris Buschelman #undef __FUNCT__
728f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
729f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
730397b6df1SKris Buschelman   int           ierr,size;
731397b6df1SKris Buschelman   MPI_Comm      comm;
732397b6df1SKris Buschelman 
733397b6df1SKris Buschelman   PetscFunctionBegin;
7345441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
7355441df8eSKris Buschelman   /*   and AIJMUMPS types */
7365441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
737397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
738397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
739397b6df1SKris Buschelman   if (size == 1) {
740397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
741397b6df1SKris Buschelman   } else {
742397b6df1SKris Buschelman     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
743397b6df1SKris Buschelman   }
744f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
745397b6df1SKris Buschelman   PetscFunctionReturn(0);
746397b6df1SKris Buschelman }
747397b6df1SKris Buschelman EXTERN_C_END
748397b6df1SKris Buschelman 
749f579278aSKris Buschelman #undef __FUNCT__
750f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
751f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
752f579278aSKris Buschelman   int       ierr;
753f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
754f579278aSKris Buschelman 
755f579278aSKris Buschelman   PetscFunctionBegin;
756f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
757f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
758f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
759f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
760f579278aSKris Buschelman   PetscFunctionReturn(0);
761f579278aSKris Buschelman }
762f579278aSKris Buschelman 
763f579278aSKris Buschelman EXTERN_C_BEGIN
764f579278aSKris Buschelman #undef __FUNCT__
765f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
766f0c56d0fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,Mat *newmat) {
767f579278aSKris Buschelman   int       ierr,size;
768f579278aSKris Buschelman   MPI_Comm  comm;
769f579278aSKris Buschelman   Mat       B=*newmat;
770f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
771f579278aSKris Buschelman 
772f579278aSKris Buschelman   PetscFunctionBegin;
773f579278aSKris Buschelman   if (B != A) {
774f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
775f579278aSKris Buschelman   }
776f579278aSKris Buschelman 
777f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
778f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
779f579278aSKris Buschelman 
780f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
781f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
782f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
783f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
784f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
785f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
786*a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
787f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
788f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
789f579278aSKris Buschelman 
790f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
791f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_SBAIJMUMPS;
792f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
793f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
794f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
7953924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
796f579278aSKris Buschelman 
797f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
798f579278aSKris Buschelman   if (size == 1) {
799f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
800f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
801f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
802f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
803f579278aSKris Buschelman   } else {
804f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
805f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
806f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
807f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
808f579278aSKris Buschelman   }
809f579278aSKris Buschelman 
810f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
811f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
812f579278aSKris Buschelman   *newmat = B;
813f579278aSKris Buschelman   PetscFunctionReturn(0);
814f579278aSKris Buschelman }
815f579278aSKris Buschelman EXTERN_C_END
816f579278aSKris Buschelman 
817f0c56d0fSKris Buschelman #undef __FUNCT__
818f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS"
819f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
820f0c56d0fSKris Buschelman   int       ierr;
8218f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
8228f340917SKris Buschelman 
823f0c56d0fSKris Buschelman   PetscFunctionBegin;
8248f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
825f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr);
8263f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
827f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
828f0c56d0fSKris Buschelman }
829f0c56d0fSKris Buschelman 
83024b6179bSKris Buschelman /*MC
831fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
83224b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
83324b6179bSKris Buschelman 
83424b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
83524b6179bSKris Buschelman   on how to declare the existence of external packages),
83624b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
83724b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
83824b6179bSKris Buschelman   This matrix type is only supported for double precision real.
83924b6179bSKris Buschelman 
84024b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
84124b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
84224b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
84324b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
84428b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
84528b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
84628b08bd3SKris Buschelman   without data copy.
84724b6179bSKris Buschelman 
84824b6179bSKris Buschelman   Options Database Keys:
8490bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
85024b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
85124b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
85224b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
85324b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
85424b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
85524b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
85624b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view
85724b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
85824b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
85924b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
86024b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
86124b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
86224b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
86324b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
86424b6179bSKris Buschelman 
86524b6179bSKris Buschelman   Level: beginner
86624b6179bSKris Buschelman 
86724b6179bSKris Buschelman .seealso: MATAIJMUMPS
86824b6179bSKris Buschelman M*/
86924b6179bSKris Buschelman 
870397b6df1SKris Buschelman EXTERN_C_BEGIN
871397b6df1SKris Buschelman #undef __FUNCT__
872f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
873f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
874397b6df1SKris Buschelman   int ierr,size;
875397b6df1SKris Buschelman 
876397b6df1SKris Buschelman   PetscFunctionBegin;
8775441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
8785441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
8795441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
8805441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
881397b6df1SKris Buschelman   if (size == 1) {
882397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
883397b6df1SKris Buschelman   } else {
884397b6df1SKris Buschelman     ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
885397b6df1SKris Buschelman   }
886f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
887397b6df1SKris Buschelman   PetscFunctionReturn(0);
888397b6df1SKris Buschelman }
889397b6df1SKris Buschelman EXTERN_C_END
890