xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision fafad7475eccdf22e42c12aa21085e705729f888)
1397b6df1SKris Buschelman /*$Id: mumps.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2397b6df1SKris Buschelman /*
3397b6df1SKris Buschelman     Provides an interface to the MUMPS_4.2_beta sparse solver
4397b6df1SKris Buschelman */
5397b6df1SKris Buschelman 
6397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
7397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h"
8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h"
9397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
10397b6df1SKris Buschelman 
11397b6df1SKris Buschelman EXTERN_C_BEGIN
12397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
13397b6df1SKris Buschelman #include "zmumps_c.h"
14397b6df1SKris Buschelman #else
15397b6df1SKris Buschelman #include "dmumps_c.h"
16397b6df1SKris Buschelman #endif
17397b6df1SKris Buschelman EXTERN_C_END
18397b6df1SKris Buschelman #define JOB_INIT -1
19397b6df1SKris Buschelman #define JOB_END -2
20397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
21397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
22397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
23397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
24397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
25397b6df1SKris Buschelman 
26397b6df1SKris Buschelman typedef struct {
27397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
28397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
29397b6df1SKris Buschelman #else
30397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
31397b6df1SKris Buschelman #endif
32397b6df1SKris Buschelman   MatStructure   matstruc;
33397b6df1SKris Buschelman   int            myid,size,*irn,*jcn,sym;
34397b6df1SKris Buschelman   PetscScalar    *val;
35397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
36397b6df1SKris Buschelman 
37c338a77dSKris Buschelman   MatType        basetype;
38c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
39f0c56d0fSKris Buschelman   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
40c338a77dSKris Buschelman   int (*MatView)(Mat,PetscViewer);
41c338a77dSKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
42c338a77dSKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
43c338a77dSKris Buschelman   int (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
44c338a77dSKris Buschelman   int (*MatDestroy)(Mat);
45f0c56d0fSKris Buschelman } Mat_MUMPS;
46f0c56d0fSKris Buschelman 
47f0c56d0fSKris Buschelman EXTERN int MatDuplicate_AIJMUMPS(Mat,MatDuplicateOption,Mat*);
48f0c56d0fSKris Buschelman EXTERN int MatDuplicate_SBAIJMUMPS(Mat,MatDuplicateOption,Mat*);
490e3434eeSKris Buschelman 
50397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
51397b6df1SKris Buschelman /*
52397b6df1SKris Buschelman   input:
53397b6df1SKris Buschelman     A       - matrix in mpiaij format
54397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
55397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
56397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
57397b6df1SKris Buschelman   output:
58397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
59397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
60397b6df1SKris Buschelman  */
61f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
62397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
63397b6df1SKris Buschelman   int         ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
64d54de34fSKris Buschelman   int         *row,*col;
65397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
66f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
67397b6df1SKris Buschelman 
68397b6df1SKris Buschelman   PetscFunctionBegin;
69397b6df1SKris Buschelman 
70397b6df1SKris Buschelman   if (mumps->isAIJ){
71397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
72397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
73397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
74397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
75397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
76397b6df1SKris Buschelman     garray = mat->garray;
77397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
78397b6df1SKris Buschelman 
79397b6df1SKris Buschelman   } else {
80397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
81397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
82397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
83847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs);
84397b6df1SKris Buschelman     nz = aa->s_nz + bb->nz;
85397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
86397b6df1SKris Buschelman     garray = mat->garray;
87397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
88397b6df1SKris Buschelman   }
89397b6df1SKris Buschelman 
90397b6df1SKris Buschelman   if (!valOnly){
91397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
92397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
93397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
94397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
95397b6df1SKris Buschelman   } else {
96397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
97397b6df1SKris Buschelman   }
98397b6df1SKris Buschelman   *nnz = nz;
99397b6df1SKris Buschelman 
100397b6df1SKris Buschelman   jj = 0; jB = 0; irow = rstart;
101397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
102397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
103397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
104397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
105397b6df1SKris Buschelman     bjj = bj + bi[i];
106397b6df1SKris Buschelman 
107397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
108397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
109397b6df1SKris Buschelman     for (j=0; j<countB; j++){
110397b6df1SKris Buschelman       jcol = garray[bjj[j]];
111397b6df1SKris Buschelman       if (jcol > colA_start) { jB = j; break; }
112397b6df1SKris Buschelman       if (j==countB-1) jB = countB;
113397b6df1SKris Buschelman     }
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;
121397b6df1SKris Buschelman       }
122397b6df1SKris Buschelman       val[jj++] = *bv++;
123397b6df1SKris Buschelman     }
124397b6df1SKris Buschelman     /* A-part */
125397b6df1SKris Buschelman     for (j=0; j<countA; j++){
126397b6df1SKris Buschelman       if (!valOnly){
127397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
128397b6df1SKris Buschelman       }
129397b6df1SKris Buschelman       val[jj++] = *av++;
130397b6df1SKris Buschelman     }
131397b6df1SKris Buschelman     /* B-part, larger col index */
132397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
133397b6df1SKris Buschelman       if (!valOnly){
134397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
135397b6df1SKris Buschelman       }
136397b6df1SKris Buschelman       val[jj++] = *bv++;
137397b6df1SKris Buschelman     }
138397b6df1SKris Buschelman     irow++;
139397b6df1SKris Buschelman   }
140397b6df1SKris Buschelman 
141397b6df1SKris Buschelman   PetscFunctionReturn(0);
142397b6df1SKris Buschelman }
143397b6df1SKris Buschelman 
144c338a77dSKris Buschelman EXTERN_C_BEGIN
145c338a77dSKris Buschelman #undef __FUNCT__
146c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
147c338a77dSKris Buschelman int MatConvert_MUMPS_Base(Mat A,MatType type,Mat *newmat) {
148c338a77dSKris Buschelman   int       ierr;
149c338a77dSKris Buschelman   Mat       B=*newmat;
150f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
151c338a77dSKris Buschelman 
152c338a77dSKris Buschelman   PetscFunctionBegin;
153c338a77dSKris Buschelman   if (B != A) {
154c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
155c338a77dSKris Buschelman   }
156f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
157f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
158f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
159f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
160f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
161f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
1623924e44cSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
1633924e44cSKris Buschelman   ierr = PetscStrfree(mumps->basetype);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;
175397b6df1SKris Buschelman 
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   }
196c338a77dSKris Buschelman   ierr = MatConvert_MUMPS_Base(A,lu->basetype,&A);CHKERRQ(ierr);
197c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
198397b6df1SKris Buschelman   PetscFunctionReturn(0);
199397b6df1SKris Buschelman }
200397b6df1SKris Buschelman 
201397b6df1SKris Buschelman #undef __FUNCT__
202c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
203f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
204f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
205397b6df1SKris Buschelman   int       ierr;
206397b6df1SKris Buschelman 
207397b6df1SKris Buschelman   PetscFunctionBegin;
208c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
209c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
210c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
211c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
212c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
213c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
214c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
215c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
216c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
217c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
218c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
219c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
220c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
221c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
222c338a77dSKris 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);
223c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
224c338a77dSKris 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);
225c338a77dSKris Buschelman 
226c338a77dSKris Buschelman   }
227c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):     %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
228c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):     %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
229c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (efficiency control):     %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
230c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):     %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
231c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):       %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
232c338a77dSKris Buschelman 
233c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
234c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
235c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
236397b6df1SKris Buschelman   PetscFunctionReturn(0);
237397b6df1SKris Buschelman }
238397b6df1SKris Buschelman 
239397b6df1SKris Buschelman #undef __FUNCT__
240f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
241f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
242397b6df1SKris Buschelman   int               ierr;
243397b6df1SKris Buschelman   PetscTruth        isascii;
244397b6df1SKris Buschelman   PetscViewerFormat format;
245f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
246397b6df1SKris Buschelman 
247397b6df1SKris Buschelman   PetscFunctionBegin;
248397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
249397b6df1SKris Buschelman 
250397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
251397b6df1SKris Buschelman   if (isascii) {
252397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
253397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
254397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
255397b6df1SKris Buschelman     }
256397b6df1SKris Buschelman   }
257397b6df1SKris Buschelman   PetscFunctionReturn(0);
258397b6df1SKris Buschelman }
259397b6df1SKris Buschelman 
260397b6df1SKris Buschelman #undef __FUNCT__
261f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
262f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
263f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
264d54de34fSKris Buschelman   PetscScalar *array;
265397b6df1SKris Buschelman   Vec         x_seq;
266397b6df1SKris Buschelman   IS          iden;
267397b6df1SKris Buschelman   VecScatter  scat;
268397b6df1SKris Buschelman   int         ierr;
269397b6df1SKris Buschelman 
270397b6df1SKris Buschelman   PetscFunctionBegin;
271397b6df1SKris Buschelman   if (lu->size > 1){
272397b6df1SKris Buschelman     if (!lu->myid){
273397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
274397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
275397b6df1SKris Buschelman     } else {
276397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
277397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
278397b6df1SKris Buschelman     }
279397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
280397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
281397b6df1SKris Buschelman 
282397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
283397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
284397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
285397b6df1SKris Buschelman   } else {  /* size == 1 */
286397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
287397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
288397b6df1SKris Buschelman   }
289397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
290397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
291397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
292397b6df1SKris Buschelman #else
293397b6df1SKris Buschelman     lu->id.rhs = array;
294397b6df1SKris Buschelman #endif
295397b6df1SKris Buschelman   }
296397b6df1SKris Buschelman 
297397b6df1SKris Buschelman   /* solve phase */
298397b6df1SKris Buschelman   lu->id.job=3;
299397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
300397b6df1SKris Buschelman   zmumps_c(&lu->id);
301397b6df1SKris Buschelman #else
302397b6df1SKris Buschelman   dmumps_c(&lu->id);
303397b6df1SKris Buschelman #endif
304397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
305397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
306397b6df1SKris Buschelman   }
307397b6df1SKris Buschelman 
308397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
309397b6df1SKris Buschelman   if (lu->size > 1) {
310397b6df1SKris Buschelman     if (!lu->myid){
311397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
312397b6df1SKris Buschelman     }
313397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
314397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
315397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
316397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
317397b6df1SKris Buschelman   } else {
318397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
319397b6df1SKris Buschelman   }
320397b6df1SKris Buschelman 
321397b6df1SKris Buschelman   PetscFunctionReturn(0);
322397b6df1SKris Buschelman }
323397b6df1SKris Buschelman 
324397b6df1SKris Buschelman #undef __FUNCT__
325f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
326f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
327f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
328f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
329397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
330397b6df1SKris Buschelman   PetscTruth valOnly,flg;
331397b6df1SKris Buschelman 
332397b6df1SKris Buschelman   PetscFunctionBegin;
333397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
334f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
335397b6df1SKris Buschelman 
336397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
337397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
338397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
339397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
340397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
341397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
342397b6df1SKris Buschelman 
343397b6df1SKris Buschelman     /* Set mumps options */
344397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
345397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
346397b6df1SKris Buschelman     lu->id.sym=lu->sym;
347397b6df1SKris Buschelman     if (lu->sym == 2){
348397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
349397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
350397b6df1SKris Buschelman     }
351397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
352397b6df1SKris Buschelman   zmumps_c(&lu->id);
353397b6df1SKris Buschelman #else
354397b6df1SKris Buschelman   dmumps_c(&lu->id);
355397b6df1SKris Buschelman #endif
356397b6df1SKris Buschelman 
357397b6df1SKris Buschelman     if (lu->size == 1){
358397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
359397b6df1SKris Buschelman     } else {
360397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
361397b6df1SKris Buschelman     }
362397b6df1SKris Buschelman 
363397b6df1SKris Buschelman     icntl=-1;
364397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
365397b6df1SKris Buschelman     if (flg && icntl > 0) {
366397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
367397b6df1SKris Buschelman     } else { /* no output */
368397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
369397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
370397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
371397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
372397b6df1SKris Buschelman     }
373397b6df1SKris 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);
374397b6df1SKris Buschelman     icntl=-1;
375397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
376397b6df1SKris Buschelman     if (flg) {
377397b6df1SKris Buschelman       if (icntl== 1){
378397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
379397b6df1SKris Buschelman       } else {
380397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
381397b6df1SKris Buschelman       }
382397b6df1SKris Buschelman     }
383397b6df1SKris 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);
384397b6df1SKris 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);
385397b6df1SKris 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);
386397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
387397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
388397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): efficiency control","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
389397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
390397b6df1SKris Buschelman 
391397b6df1SKris Buschelman     /*
392397b6df1SKris 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);
393397b6df1SKris Buschelman     if (flg){
394397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
395397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
396397b6df1SKris Buschelman       } else {
397397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
398397b6df1SKris Buschelman       }
399397b6df1SKris Buschelman     }
400397b6df1SKris Buschelman     */
401397b6df1SKris Buschelman 
402397b6df1SKris 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);
403397b6df1SKris 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);
404397b6df1SKris 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);
405397b6df1SKris Buschelman     PetscOptionsEnd();
406397b6df1SKris Buschelman   }
407397b6df1SKris Buschelman 
408397b6df1SKris Buschelman   /* define matrix A */
409397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
410397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
411397b6df1SKris Buschelman     if (!lu->myid) {
412c36ead0aSKris Buschelman       if (lua->isAIJ){
413397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
414397b6df1SKris Buschelman         nz               = aa->nz;
415397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
416397b6df1SKris Buschelman       } else {
417397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
418397b6df1SKris Buschelman         nz                  =  aa->s_nz;
419397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
420397b6df1SKris Buschelman       }
421397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
422397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
423397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
424397b6df1SKris Buschelman         nz = 0;
425397b6df1SKris Buschelman         for (i=0; i<M; i++){
426397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
427397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
428397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
429397b6df1SKris Buschelman           }
430397b6df1SKris Buschelman         }
431397b6df1SKris Buschelman       }
432397b6df1SKris Buschelman     }
433397b6df1SKris Buschelman     break;
434397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
435397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
436397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
437397b6df1SKris Buschelman     } else {
438397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
439397b6df1SKris Buschelman     }
440397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
441397b6df1SKris Buschelman     break;
442397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
443397b6df1SKris Buschelman   }
444397b6df1SKris Buschelman 
445397b6df1SKris Buschelman   /* analysis phase */
446397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
447397b6df1SKris Buschelman      lu->id.n = M;
448397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
449397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
450397b6df1SKris Buschelman       if (!lu->myid) {
451397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
452397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
453397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
454397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
455397b6df1SKris Buschelman #else
456397b6df1SKris Buschelman           lu->id.a = lu->val;
457397b6df1SKris Buschelman #endif
458397b6df1SKris Buschelman         }
459397b6df1SKris Buschelman       }
460397b6df1SKris Buschelman       break;
461397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
462397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
463397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
464397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
465397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
466397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
467397b6df1SKris Buschelman #else
468397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
469397b6df1SKris Buschelman #endif
470397b6df1SKris Buschelman       }
471397b6df1SKris Buschelman       break;
472397b6df1SKris Buschelman     }
473397b6df1SKris Buschelman     lu->id.job=1;
474397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
475397b6df1SKris Buschelman   zmumps_c(&lu->id);
476397b6df1SKris Buschelman #else
477397b6df1SKris Buschelman   dmumps_c(&lu->id);
478397b6df1SKris Buschelman #endif
479397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
480397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
481397b6df1SKris Buschelman     }
482397b6df1SKris Buschelman   }
483397b6df1SKris Buschelman 
484397b6df1SKris Buschelman   /* numerical factorization phase */
485397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
486397b6df1SKris Buschelman     if (lu->myid == 0) {
487397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
488397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
489397b6df1SKris Buschelman #else
490397b6df1SKris Buschelman       lu->id.a = lu->val;
491397b6df1SKris Buschelman #endif
492397b6df1SKris Buschelman     }
493397b6df1SKris Buschelman   } else {
494397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
495397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
496397b6df1SKris Buschelman #else
497397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
498397b6df1SKris Buschelman #endif
499397b6df1SKris Buschelman   }
500397b6df1SKris Buschelman   lu->id.job=2;
501397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
502397b6df1SKris Buschelman   zmumps_c(&lu->id);
503397b6df1SKris Buschelman #else
504397b6df1SKris Buschelman   dmumps_c(&lu->id);
505397b6df1SKris Buschelman #endif
506397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
507397b6df1SKris Buschelman     SETERRQ1(1,"1, Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d\n",lu->id.INFOG(1));
508397b6df1SKris Buschelman   }
509397b6df1SKris Buschelman 
510397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
511397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
512397b6df1SKris Buschelman   }
513397b6df1SKris Buschelman 
514397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
515397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
516ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
517397b6df1SKris Buschelman   PetscFunctionReturn(0);
518397b6df1SKris Buschelman }
519397b6df1SKris Buschelman 
520397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
521397b6df1SKris Buschelman #undef __FUNCT__
522f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
523f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
524397b6df1SKris Buschelman   Mat       B;
525f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
526397b6df1SKris Buschelman   int       ierr;
527397b6df1SKris Buschelman 
528397b6df1SKris Buschelman   PetscFunctionBegin;
529397b6df1SKris Buschelman 
530397b6df1SKris Buschelman   /* Create the factorization matrix */
531397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
532397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
533397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
534397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
535397b6df1SKris Buschelman 
536f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
537397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
538f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
539397b6df1SKris Buschelman   lu->sym                 = 0;
540397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
541397b6df1SKris Buschelman 
542397b6df1SKris Buschelman   *F = B;
543397b6df1SKris Buschelman   PetscFunctionReturn(0);
544397b6df1SKris Buschelman }
545397b6df1SKris Buschelman 
546397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
547397b6df1SKris Buschelman #undef __FUNCT__
548f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
549f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
550397b6df1SKris Buschelman   Mat       B;
551f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
552397b6df1SKris Buschelman   int       ierr;
553397b6df1SKris Buschelman 
554397b6df1SKris Buschelman   PetscFunctionBegin;
555397b6df1SKris Buschelman 
556397b6df1SKris Buschelman   /* Create the factorization matrix */
557397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
558397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
559397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
560397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
561397b6df1SKris Buschelman 
562f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
563397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
564f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
565397b6df1SKris Buschelman   lu->sym                       = 2;
566397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
567397b6df1SKris Buschelman 
568397b6df1SKris Buschelman   *F = B;
569397b6df1SKris Buschelman   PetscFunctionReturn(0);
570397b6df1SKris Buschelman }
571397b6df1SKris Buschelman 
572397b6df1SKris Buschelman #undef __FUNCT__
573f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
574f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
575c338a77dSKris Buschelman   int       ierr;
576f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
577c338a77dSKris Buschelman 
578397b6df1SKris Buschelman   PetscFunctionBegin;
579c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
580f0c56d0fSKris Buschelman 
581c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
582c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
583f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
584397b6df1SKris Buschelman   PetscFunctionReturn(0);
585397b6df1SKris Buschelman }
586397b6df1SKris Buschelman 
587c338a77dSKris Buschelman EXTERN_C_BEGIN
588c338a77dSKris Buschelman #undef __FUNCT__
589f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
590f0c56d0fSKris Buschelman int MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,Mat *newmat) {
591c338a77dSKris Buschelman   int       ierr,size;
592c338a77dSKris Buschelman   MPI_Comm  comm;
593c338a77dSKris Buschelman   Mat       B=*newmat;
594f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
595397b6df1SKris Buschelman 
596397b6df1SKris Buschelman   PetscFunctionBegin;
597c338a77dSKris Buschelman   if (B != A) {
598c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
599397b6df1SKris Buschelman   }
600397b6df1SKris Buschelman 
601c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
602f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
603c338a77dSKris Buschelman 
604f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
605c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
606c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
607c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
608c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
609c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
610c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
611f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
612c338a77dSKris Buschelman 
6134b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
614f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_AIJMUMPS;
615f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
616f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
617f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
6183924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
619c338a77dSKris Buschelman 
620c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
621c338a77dSKris Buschelman   if (size == 1) {
6223924e44cSKris Buschelman     ierr = PetscStrallocpy(MATSEQAIJ,&(mumps->basetype));CHKERRQ(ierr);
623c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
624f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
625c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
626c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
627c338a77dSKris Buschelman   } else {
6283924e44cSKris Buschelman     ierr = PetscStrallocpy(MATMPIAIJ,&(mumps->basetype));CHKERRQ(ierr);
629c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
630f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
631c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
632c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
633c338a77dSKris Buschelman   }
634c338a77dSKris Buschelman 
635f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
636c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
637c338a77dSKris Buschelman   *newmat = B;
638397b6df1SKris Buschelman   PetscFunctionReturn(0);
639397b6df1SKris Buschelman }
640c338a77dSKris Buschelman EXTERN_C_END
641397b6df1SKris Buschelman 
642f0c56d0fSKris Buschelman #undef __FUNCT__
643f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS"
644f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
645f0c56d0fSKris Buschelman   int ierr;
646f0c56d0fSKris Buschelman   PetscFunctionBegin;
647f0c56d0fSKris Buschelman   ierr = (*A->ops->duplicate)(A,op,M);CHKERRQ(ierr);
648f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr);
649f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
650f0c56d0fSKris Buschelman }
651f0c56d0fSKris Buschelman 
65224b6179bSKris Buschelman /*MC
653*fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
65424b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
65524b6179bSKris Buschelman 
65624b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
65724b6179bSKris Buschelman   on how to declare the existence of external packages),
65824b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
65924b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
66024b6179bSKris Buschelman   This matrix type is only supported for double precision real.
66124b6179bSKris Buschelman 
66224b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
66324b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
66424b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
66524b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
66628b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
66728b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
66828b08bd3SKris Buschelman   without data copy.
66924b6179bSKris Buschelman 
67024b6179bSKris Buschelman   Options Database Keys:
6710bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
67224b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
67324b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
67424b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
67524b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
67624b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
67724b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
67824b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view
67924b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
68024b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
68124b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
68224b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
68324b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
68424b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
68524b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
68624b6179bSKris Buschelman 
68724b6179bSKris Buschelman   Level: beginner
68824b6179bSKris Buschelman 
68924b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
69024b6179bSKris Buschelman M*/
69124b6179bSKris Buschelman 
692397b6df1SKris Buschelman EXTERN_C_BEGIN
693397b6df1SKris Buschelman #undef __FUNCT__
694f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
695f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
696397b6df1SKris Buschelman   int           ierr,size;
697397b6df1SKris Buschelman   MPI_Comm      comm;
698397b6df1SKris Buschelman 
699397b6df1SKris Buschelman   PetscFunctionBegin;
7005441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
7015441df8eSKris Buschelman   /*   and AIJMUMPS types */
7025441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
703397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
704397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
705397b6df1SKris Buschelman   if (size == 1) {
706397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
707397b6df1SKris Buschelman   } else {
708397b6df1SKris Buschelman     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
709397b6df1SKris Buschelman   }
710f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
711397b6df1SKris Buschelman   PetscFunctionReturn(0);
712397b6df1SKris Buschelman }
713397b6df1SKris Buschelman EXTERN_C_END
714397b6df1SKris Buschelman 
715f579278aSKris Buschelman #undef __FUNCT__
716f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
717f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
718f579278aSKris Buschelman   int       ierr;
719f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
720f579278aSKris Buschelman 
721f579278aSKris Buschelman   PetscFunctionBegin;
722f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
723f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
724f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
725f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
726f579278aSKris Buschelman   PetscFunctionReturn(0);
727f579278aSKris Buschelman }
728f579278aSKris Buschelman 
729f579278aSKris Buschelman EXTERN_C_BEGIN
730f579278aSKris Buschelman #undef __FUNCT__
731f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
732f0c56d0fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,Mat *newmat) {
733f579278aSKris Buschelman   int       ierr,size;
734f579278aSKris Buschelman   MPI_Comm  comm;
735f579278aSKris Buschelman   Mat       B=*newmat;
736f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
737f579278aSKris Buschelman 
738f579278aSKris Buschelman   PetscFunctionBegin;
739f579278aSKris Buschelman   if (B != A) {
740f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
741f579278aSKris Buschelman   }
742f579278aSKris Buschelman 
743f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
744f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
745f579278aSKris Buschelman 
746f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
747f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
748f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
749f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
750f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
751f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
752f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
753f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
754f579278aSKris Buschelman 
755f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
756f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_SBAIJMUMPS;
757f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
758f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
759f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
7603924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
761f579278aSKris Buschelman 
762f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
763f579278aSKris Buschelman   if (size == 1) {
7643924e44cSKris Buschelman     ierr = PetscStrallocpy(MATSEQSBAIJ,&(mumps->basetype));CHKERRQ(ierr);
765f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
766f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
767f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
768f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
769f579278aSKris Buschelman   } else {
7703924e44cSKris Buschelman     ierr = PetscStrallocpy(MATMPISBAIJ,&(mumps->basetype));CHKERRQ(ierr);
771f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
772f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
773f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
774f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
775f579278aSKris Buschelman   }
776f579278aSKris Buschelman 
777f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
778f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
779f579278aSKris Buschelman   *newmat = B;
780f579278aSKris Buschelman   PetscFunctionReturn(0);
781f579278aSKris Buschelman }
782f579278aSKris Buschelman EXTERN_C_END
783f579278aSKris Buschelman 
784f0c56d0fSKris Buschelman #undef __FUNCT__
785f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS"
786f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
787f0c56d0fSKris Buschelman   int ierr;
788f0c56d0fSKris Buschelman   PetscFunctionBegin;
789f0c56d0fSKris Buschelman   ierr = (*A->ops->duplicate)(A,op,M);CHKERRQ(ierr);
790f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr);
791f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
792f0c56d0fSKris Buschelman }
793f0c56d0fSKris Buschelman 
79424b6179bSKris Buschelman /*MC
795*fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
79624b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
79724b6179bSKris Buschelman 
79824b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
79924b6179bSKris Buschelman   on how to declare the existence of external packages),
80024b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
80124b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
80224b6179bSKris Buschelman   This matrix type is only supported for double precision real.
80324b6179bSKris Buschelman 
80424b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
80524b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
80624b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
80724b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
80828b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
80928b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
81028b08bd3SKris Buschelman   without data copy.
81124b6179bSKris Buschelman 
81224b6179bSKris Buschelman   Options Database Keys:
8130bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
81424b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
81524b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
81624b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
81724b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
81824b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
81924b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
82024b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view
82124b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
82224b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
82324b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
82424b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
82524b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
82624b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
82724b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
82824b6179bSKris Buschelman 
82924b6179bSKris Buschelman   Level: beginner
83024b6179bSKris Buschelman 
83124b6179bSKris Buschelman .seealso: MATAIJMUMPS
83224b6179bSKris Buschelman M*/
83324b6179bSKris Buschelman 
834397b6df1SKris Buschelman EXTERN_C_BEGIN
835397b6df1SKris Buschelman #undef __FUNCT__
836f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
837f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
838397b6df1SKris Buschelman   int ierr,size;
839397b6df1SKris Buschelman 
840397b6df1SKris Buschelman   PetscFunctionBegin;
8415441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
8425441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
8435441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
8445441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
845397b6df1SKris Buschelman   if (size == 1) {
846397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
847397b6df1SKris Buschelman   } else {
848397b6df1SKris Buschelman     ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
849397b6df1SKris Buschelman   }
850f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
851397b6df1SKris Buschelman   PetscFunctionReturn(0);
852397b6df1SKris Buschelman }
853397b6df1SKris Buschelman EXTERN_C_END
854