xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 028e57e81cd5976a9a2ed6be4d5a7e7af0db655e)
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 #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   MatType        basetype;
37c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
38f0c56d0fSKris Buschelman   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
39c338a77dSKris Buschelman   int (*MatView)(Mat,PetscViewer);
40c338a77dSKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
41c338a77dSKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
42c338a77dSKris Buschelman   int (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
43c338a77dSKris Buschelman   int (*MatDestroy)(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 
98*028e57e8SHong 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);
1643924e44cSKris Buschelman   ierr = PetscStrfree(mumps->basetype);CHKERRQ(ierr);
165f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
166c338a77dSKris Buschelman   *newmat = B;
167c338a77dSKris Buschelman   PetscFunctionReturn(0);
168c338a77dSKris Buschelman }
169c338a77dSKris Buschelman EXTERN_C_END
170c338a77dSKris Buschelman 
171397b6df1SKris Buschelman #undef __FUNCT__
1723924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
1733924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) {
174f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
175c338a77dSKris Buschelman   int       ierr,size=lu->size;
176397b6df1SKris Buschelman 
177397b6df1SKris Buschelman   PetscFunctionBegin;
178397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
179397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
180397b6df1SKris Buschelman     lu->id.job=JOB_END;
181397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
182397b6df1SKris Buschelman     zmumps_c(&lu->id);
183397b6df1SKris Buschelman #else
184397b6df1SKris Buschelman     dmumps_c(&lu->id);
185397b6df1SKris Buschelman #endif
186c338a77dSKris Buschelman     if (lu->irn) {
187c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
188c338a77dSKris Buschelman     }
189c338a77dSKris Buschelman     if (lu->jcn) {
190c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
191c338a77dSKris Buschelman     }
192c338a77dSKris Buschelman     if (size>1 && lu->val) {
193c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
194c338a77dSKris Buschelman     }
195397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
196397b6df1SKris Buschelman   }
197c338a77dSKris Buschelman   ierr = MatConvert_MUMPS_Base(A,lu->basetype,&A);CHKERRQ(ierr);
198c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
199397b6df1SKris Buschelman   PetscFunctionReturn(0);
200397b6df1SKris Buschelman }
201397b6df1SKris Buschelman 
202397b6df1SKris Buschelman #undef __FUNCT__
203c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
204f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
205f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
206397b6df1SKris Buschelman   int       ierr;
207397b6df1SKris Buschelman 
208397b6df1SKris Buschelman   PetscFunctionBegin;
209c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
210c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
211c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
212c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
213c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
214c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
215c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
216c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
217c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
218c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
219c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
220c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
221c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
222c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
223c338a77dSKris 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);
224c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
225c338a77dSKris 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);
226c338a77dSKris Buschelman 
227c338a77dSKris Buschelman   }
228c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):     %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
229c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):     %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
230c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (efficiency control):     %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
231c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):     %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
232c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):       %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
233c338a77dSKris Buschelman 
234c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
235c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
236c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
237397b6df1SKris Buschelman   PetscFunctionReturn(0);
238397b6df1SKris Buschelman }
239397b6df1SKris Buschelman 
240397b6df1SKris Buschelman #undef __FUNCT__
241f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
242f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
243397b6df1SKris Buschelman   int               ierr;
244397b6df1SKris Buschelman   PetscTruth        isascii;
245397b6df1SKris Buschelman   PetscViewerFormat format;
246f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
247397b6df1SKris Buschelman 
248397b6df1SKris Buschelman   PetscFunctionBegin;
249397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
250397b6df1SKris Buschelman 
251397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
252397b6df1SKris Buschelman   if (isascii) {
253397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
254397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
255397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
256397b6df1SKris Buschelman     }
257397b6df1SKris Buschelman   }
258397b6df1SKris Buschelman   PetscFunctionReturn(0);
259397b6df1SKris Buschelman }
260397b6df1SKris Buschelman 
261397b6df1SKris Buschelman #undef __FUNCT__
262f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
263f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
264f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
265d54de34fSKris Buschelman   PetscScalar *array;
266397b6df1SKris Buschelman   Vec         x_seq;
267397b6df1SKris Buschelman   IS          iden;
268397b6df1SKris Buschelman   VecScatter  scat;
269397b6df1SKris Buschelman   int         ierr;
270397b6df1SKris Buschelman 
271397b6df1SKris Buschelman   PetscFunctionBegin;
272397b6df1SKris Buschelman   if (lu->size > 1){
273397b6df1SKris Buschelman     if (!lu->myid){
274397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
275397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
276397b6df1SKris Buschelman     } else {
277397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
278397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
279397b6df1SKris Buschelman     }
280397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
281397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
282397b6df1SKris Buschelman 
283397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
284397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
285397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
286397b6df1SKris Buschelman   } else {  /* size == 1 */
287397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
288397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
289397b6df1SKris Buschelman   }
290397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
291397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
292397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
293397b6df1SKris Buschelman #else
294397b6df1SKris Buschelman     lu->id.rhs = array;
295397b6df1SKris Buschelman #endif
296397b6df1SKris Buschelman   }
297397b6df1SKris Buschelman 
298397b6df1SKris Buschelman   /* solve phase */
299397b6df1SKris Buschelman   lu->id.job=3;
300397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
301397b6df1SKris Buschelman   zmumps_c(&lu->id);
302397b6df1SKris Buschelman #else
303397b6df1SKris Buschelman   dmumps_c(&lu->id);
304397b6df1SKris Buschelman #endif
305397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
306397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
307397b6df1SKris Buschelman   }
308397b6df1SKris Buschelman 
309397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
310397b6df1SKris Buschelman   if (lu->size > 1) {
311397b6df1SKris Buschelman     if (!lu->myid){
312397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
313397b6df1SKris Buschelman     }
314397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
315397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
316397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
317397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
318397b6df1SKris Buschelman   } else {
319397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
320397b6df1SKris Buschelman   }
321397b6df1SKris Buschelman 
322397b6df1SKris Buschelman   PetscFunctionReturn(0);
323397b6df1SKris Buschelman }
324397b6df1SKris Buschelman 
325397b6df1SKris Buschelman #undef __FUNCT__
326f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
327f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
328f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
329f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
330397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
331397b6df1SKris Buschelman   PetscTruth valOnly,flg;
332397b6df1SKris Buschelman 
333397b6df1SKris Buschelman   PetscFunctionBegin;
334397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
335f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
336397b6df1SKris Buschelman 
337397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
338397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
339397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
34075747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
341397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
342397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
343397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
344397b6df1SKris Buschelman 
345397b6df1SKris Buschelman     /* Set mumps options */
346397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
347397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
348397b6df1SKris Buschelman     lu->id.sym=lu->sym;
349397b6df1SKris Buschelman     if (lu->sym == 2){
350397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
351397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
352397b6df1SKris Buschelman     }
353397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
354397b6df1SKris Buschelman   zmumps_c(&lu->id);
355397b6df1SKris Buschelman #else
356397b6df1SKris Buschelman   dmumps_c(&lu->id);
357397b6df1SKris Buschelman #endif
358397b6df1SKris Buschelman 
359397b6df1SKris Buschelman     if (lu->size == 1){
360397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
361397b6df1SKris Buschelman     } else {
362397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
363397b6df1SKris Buschelman     }
364397b6df1SKris Buschelman 
365397b6df1SKris Buschelman     icntl=-1;
366397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
367397b6df1SKris Buschelman     if (flg && icntl > 0) {
368397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
369397b6df1SKris Buschelman     } else { /* no output */
370397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
371397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
372397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
373397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
374397b6df1SKris Buschelman     }
375397b6df1SKris 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);
376397b6df1SKris Buschelman     icntl=-1;
377397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
378397b6df1SKris Buschelman     if (flg) {
379397b6df1SKris Buschelman       if (icntl== 1){
380397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
381397b6df1SKris Buschelman       } else {
382397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
383397b6df1SKris Buschelman       }
384397b6df1SKris Buschelman     }
385397b6df1SKris 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);
386397b6df1SKris 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);
387397b6df1SKris 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);
388397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
389397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
390397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): efficiency control","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
391397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
392397b6df1SKris Buschelman 
393397b6df1SKris Buschelman     /*
394397b6df1SKris 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);
395397b6df1SKris Buschelman     if (flg){
396397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
397397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
398397b6df1SKris Buschelman       } else {
399397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
400397b6df1SKris Buschelman       }
401397b6df1SKris Buschelman     }
402397b6df1SKris Buschelman     */
403397b6df1SKris Buschelman 
404397b6df1SKris 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);
405397b6df1SKris 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);
406397b6df1SKris 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);
407397b6df1SKris Buschelman     PetscOptionsEnd();
408397b6df1SKris Buschelman   }
409397b6df1SKris Buschelman 
410397b6df1SKris Buschelman   /* define matrix A */
411397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
412397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
413397b6df1SKris Buschelman     if (!lu->myid) {
414c36ead0aSKris Buschelman       if (lua->isAIJ){
415397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
416397b6df1SKris Buschelman         nz               = aa->nz;
417397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
418397b6df1SKris Buschelman       } else {
419397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
420397b6df1SKris Buschelman         nz                  =  aa->s_nz;
421397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
422397b6df1SKris Buschelman       }
423397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
424397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
425397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
426397b6df1SKris Buschelman         nz = 0;
427397b6df1SKris Buschelman         for (i=0; i<M; i++){
428397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
429397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
430397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
431397b6df1SKris Buschelman           }
432397b6df1SKris Buschelman         }
433397b6df1SKris Buschelman       }
434397b6df1SKris Buschelman     }
435397b6df1SKris Buschelman     break;
436397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
437397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
438397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
439397b6df1SKris Buschelman     } else {
440397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
441397b6df1SKris Buschelman     }
442397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
443397b6df1SKris Buschelman     break;
444397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
445397b6df1SKris Buschelman   }
446397b6df1SKris Buschelman 
447397b6df1SKris Buschelman   /* analysis phase */
448397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
449397b6df1SKris Buschelman      lu->id.n = M;
450397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
451397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
452397b6df1SKris Buschelman       if (!lu->myid) {
453397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
454397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
455397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
456397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
457397b6df1SKris Buschelman #else
458397b6df1SKris Buschelman           lu->id.a = lu->val;
459397b6df1SKris Buschelman #endif
460397b6df1SKris Buschelman         }
461397b6df1SKris Buschelman       }
462397b6df1SKris Buschelman       break;
463397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
464397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
465397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
466397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
467397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
468397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
469397b6df1SKris Buschelman #else
470397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
471397b6df1SKris Buschelman #endif
472397b6df1SKris Buschelman       }
473397b6df1SKris Buschelman       break;
474397b6df1SKris Buschelman     }
475397b6df1SKris Buschelman     lu->id.job=1;
476397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
477397b6df1SKris Buschelman   zmumps_c(&lu->id);
478397b6df1SKris Buschelman #else
479397b6df1SKris Buschelman   dmumps_c(&lu->id);
480397b6df1SKris Buschelman #endif
481397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
482397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
483397b6df1SKris Buschelman     }
484397b6df1SKris Buschelman   }
485397b6df1SKris Buschelman 
486397b6df1SKris Buschelman   /* numerical factorization phase */
487397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
488397b6df1SKris Buschelman     if (lu->myid == 0) {
489397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
490397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
491397b6df1SKris Buschelman #else
492397b6df1SKris Buschelman       lu->id.a = lu->val;
493397b6df1SKris Buschelman #endif
494397b6df1SKris Buschelman     }
495397b6df1SKris Buschelman   } else {
496397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
497397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
498397b6df1SKris Buschelman #else
499397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
500397b6df1SKris Buschelman #endif
501397b6df1SKris Buschelman   }
502397b6df1SKris Buschelman   lu->id.job=2;
503397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
504397b6df1SKris Buschelman   zmumps_c(&lu->id);
505397b6df1SKris Buschelman #else
506397b6df1SKris Buschelman   dmumps_c(&lu->id);
507397b6df1SKris Buschelman #endif
508397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
509397b6df1SKris Buschelman     SETERRQ1(1,"1, Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d\n",lu->id.INFOG(1));
510397b6df1SKris Buschelman   }
511397b6df1SKris Buschelman 
512397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
513397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
514397b6df1SKris Buschelman   }
515397b6df1SKris Buschelman 
516397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
517397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
518ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
519397b6df1SKris Buschelman   PetscFunctionReturn(0);
520397b6df1SKris Buschelman }
521397b6df1SKris Buschelman 
522397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
523397b6df1SKris Buschelman #undef __FUNCT__
524f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
525f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
526397b6df1SKris Buschelman   Mat       B;
527f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
528397b6df1SKris Buschelman   int       ierr;
529397b6df1SKris Buschelman 
530397b6df1SKris Buschelman   PetscFunctionBegin;
531397b6df1SKris Buschelman 
532397b6df1SKris Buschelman   /* Create the factorization matrix */
533397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
534397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
535397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
536397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
537397b6df1SKris Buschelman 
538f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
539397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
540f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
541397b6df1SKris Buschelman   lu->sym                 = 0;
542397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
543397b6df1SKris Buschelman 
544397b6df1SKris Buschelman   *F = B;
545397b6df1SKris Buschelman   PetscFunctionReturn(0);
546397b6df1SKris Buschelman }
547397b6df1SKris Buschelman 
548397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
549397b6df1SKris Buschelman #undef __FUNCT__
550f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
551f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
552397b6df1SKris Buschelman   Mat       B;
553f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
554397b6df1SKris Buschelman   int       ierr;
555397b6df1SKris Buschelman 
556397b6df1SKris Buschelman   PetscFunctionBegin;
557397b6df1SKris Buschelman 
558397b6df1SKris Buschelman   /* Create the factorization matrix */
559397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
560397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
561397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
562397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
563397b6df1SKris Buschelman 
564f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
565397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
566f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
567397b6df1SKris Buschelman   lu->sym                       = 2;
568397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
569397b6df1SKris Buschelman 
570397b6df1SKris Buschelman   *F = B;
571397b6df1SKris Buschelman   PetscFunctionReturn(0);
572397b6df1SKris Buschelman }
573397b6df1SKris Buschelman 
574397b6df1SKris Buschelman #undef __FUNCT__
575f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
576f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
577c338a77dSKris Buschelman   int       ierr;
578f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
579c338a77dSKris Buschelman 
580397b6df1SKris Buschelman   PetscFunctionBegin;
581c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
582f0c56d0fSKris Buschelman 
583c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
584c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
585f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
586397b6df1SKris Buschelman   PetscFunctionReturn(0);
587397b6df1SKris Buschelman }
588397b6df1SKris Buschelman 
589c338a77dSKris Buschelman EXTERN_C_BEGIN
590c338a77dSKris Buschelman #undef __FUNCT__
591f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
592f0c56d0fSKris Buschelman int MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,Mat *newmat) {
593c338a77dSKris Buschelman   int       ierr,size;
594c338a77dSKris Buschelman   MPI_Comm  comm;
595c338a77dSKris Buschelman   Mat       B=*newmat;
596f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
597397b6df1SKris Buschelman 
598397b6df1SKris Buschelman   PetscFunctionBegin;
599c338a77dSKris Buschelman   if (B != A) {
600c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
601397b6df1SKris Buschelman   }
602397b6df1SKris Buschelman 
603c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
604f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
605c338a77dSKris Buschelman 
606f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
607c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
608c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
609c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
610c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
611c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
612c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
613f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
614c338a77dSKris Buschelman 
6154b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
616f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_AIJMUMPS;
617f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
618f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
619f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
6203924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
621c338a77dSKris Buschelman 
622c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
623c338a77dSKris Buschelman   if (size == 1) {
6243924e44cSKris Buschelman     ierr = PetscStrallocpy(MATSEQAIJ,&(mumps->basetype));CHKERRQ(ierr);
625c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
626f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
627c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
628c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
629c338a77dSKris Buschelman   } else {
6303924e44cSKris Buschelman     ierr = PetscStrallocpy(MATMPIAIJ,&(mumps->basetype));CHKERRQ(ierr);
631c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
632f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
633c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
634c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
635c338a77dSKris Buschelman   }
636c338a77dSKris Buschelman 
637f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
638c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
639c338a77dSKris Buschelman   *newmat = B;
640397b6df1SKris Buschelman   PetscFunctionReturn(0);
641397b6df1SKris Buschelman }
642c338a77dSKris Buschelman EXTERN_C_END
643397b6df1SKris Buschelman 
644f0c56d0fSKris Buschelman #undef __FUNCT__
645f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS"
646f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
647f0c56d0fSKris Buschelman   int ierr;
648f0c56d0fSKris Buschelman   PetscFunctionBegin;
649f0c56d0fSKris Buschelman   ierr = (*A->ops->duplicate)(A,op,M);CHKERRQ(ierr);
650f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr);
651f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
652f0c56d0fSKris Buschelman }
653f0c56d0fSKris Buschelman 
65424b6179bSKris Buschelman /*MC
655fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
65624b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
65724b6179bSKris Buschelman 
65824b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
65924b6179bSKris Buschelman   on how to declare the existence of external packages),
66024b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
66124b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
66224b6179bSKris Buschelman   This matrix type is only supported for double precision real.
66324b6179bSKris Buschelman 
66424b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
66524b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
66624b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
66724b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
66828b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
66928b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
67028b08bd3SKris Buschelman   without data copy.
67124b6179bSKris Buschelman 
67224b6179bSKris Buschelman   Options Database Keys:
6730bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
67424b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
67524b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
67624b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
67724b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
67824b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
67924b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
68024b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view
68124b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
68224b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
68324b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
68424b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
68524b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
68624b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
68724b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
68824b6179bSKris Buschelman 
68924b6179bSKris Buschelman   Level: beginner
69024b6179bSKris Buschelman 
69124b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
69224b6179bSKris Buschelman M*/
69324b6179bSKris Buschelman 
694397b6df1SKris Buschelman EXTERN_C_BEGIN
695397b6df1SKris Buschelman #undef __FUNCT__
696f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
697f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
698397b6df1SKris Buschelman   int           ierr,size;
699397b6df1SKris Buschelman   MPI_Comm      comm;
700397b6df1SKris Buschelman 
701397b6df1SKris Buschelman   PetscFunctionBegin;
7025441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
7035441df8eSKris Buschelman   /*   and AIJMUMPS types */
7045441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
705397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
706397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
707397b6df1SKris Buschelman   if (size == 1) {
708397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
709397b6df1SKris Buschelman   } else {
710397b6df1SKris Buschelman     ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
711397b6df1SKris Buschelman   }
712f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
713397b6df1SKris Buschelman   PetscFunctionReturn(0);
714397b6df1SKris Buschelman }
715397b6df1SKris Buschelman EXTERN_C_END
716397b6df1SKris Buschelman 
717f579278aSKris Buschelman #undef __FUNCT__
718f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
719f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
720f579278aSKris Buschelman   int       ierr;
721f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
722f579278aSKris Buschelman 
723f579278aSKris Buschelman   PetscFunctionBegin;
724f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
725f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
726f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
727f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
728f579278aSKris Buschelman   PetscFunctionReturn(0);
729f579278aSKris Buschelman }
730f579278aSKris Buschelman 
731f579278aSKris Buschelman EXTERN_C_BEGIN
732f579278aSKris Buschelman #undef __FUNCT__
733f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
734f0c56d0fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,Mat *newmat) {
735f579278aSKris Buschelman   int       ierr,size;
736f579278aSKris Buschelman   MPI_Comm  comm;
737f579278aSKris Buschelman   Mat       B=*newmat;
738f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
739f579278aSKris Buschelman 
740f579278aSKris Buschelman   PetscFunctionBegin;
741f579278aSKris Buschelman   if (B != A) {
742f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
743f579278aSKris Buschelman   }
744f579278aSKris Buschelman 
745f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
746f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
747f579278aSKris Buschelman 
748f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
749f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
750f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
751f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
752f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
753f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
754f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
755f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
756f579278aSKris Buschelman 
757f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
758f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_SBAIJMUMPS;
759f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
760f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
761f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
7623924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
763f579278aSKris Buschelman 
764f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
765f579278aSKris Buschelman   if (size == 1) {
7663924e44cSKris Buschelman     ierr = PetscStrallocpy(MATSEQSBAIJ,&(mumps->basetype));CHKERRQ(ierr);
767f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
768f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
769f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
770f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
771f579278aSKris Buschelman   } else {
7723924e44cSKris Buschelman     ierr = PetscStrallocpy(MATMPISBAIJ,&(mumps->basetype));CHKERRQ(ierr);
773f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
774f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
775f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
776f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
777f579278aSKris Buschelman   }
778f579278aSKris Buschelman 
779f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
780f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
781f579278aSKris Buschelman   *newmat = B;
782f579278aSKris Buschelman   PetscFunctionReturn(0);
783f579278aSKris Buschelman }
784f579278aSKris Buschelman EXTERN_C_END
785f579278aSKris Buschelman 
786f0c56d0fSKris Buschelman #undef __FUNCT__
787f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS"
788f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
789f0c56d0fSKris Buschelman   int ierr;
790f0c56d0fSKris Buschelman   PetscFunctionBegin;
791f0c56d0fSKris Buschelman   ierr = (*A->ops->duplicate)(A,op,M);CHKERRQ(ierr);
792f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr);
793f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
794f0c56d0fSKris Buschelman }
795f0c56d0fSKris Buschelman 
79624b6179bSKris Buschelman /*MC
797fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
79824b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
79924b6179bSKris Buschelman 
80024b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
80124b6179bSKris Buschelman   on how to declare the existence of external packages),
80224b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
80324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
80424b6179bSKris Buschelman   This matrix type is only supported for double precision real.
80524b6179bSKris Buschelman 
80624b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
80724b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
80824b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
80924b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
81028b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
81128b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
81228b08bd3SKris Buschelman   without data copy.
81324b6179bSKris Buschelman 
81424b6179bSKris Buschelman   Options Database Keys:
8150bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
81624b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
81724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
81824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
81924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
82024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
82124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
82224b6179bSKris Buschelman . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -sles_view
82324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
82424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
82524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
82624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
82724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
82824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
82924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
83024b6179bSKris Buschelman 
83124b6179bSKris Buschelman   Level: beginner
83224b6179bSKris Buschelman 
83324b6179bSKris Buschelman .seealso: MATAIJMUMPS
83424b6179bSKris Buschelman M*/
83524b6179bSKris Buschelman 
836397b6df1SKris Buschelman EXTERN_C_BEGIN
837397b6df1SKris Buschelman #undef __FUNCT__
838f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
839f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
840397b6df1SKris Buschelman   int ierr,size;
841397b6df1SKris Buschelman 
842397b6df1SKris Buschelman   PetscFunctionBegin;
8435441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
8445441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
8455441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
8465441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
847397b6df1SKris Buschelman   if (size == 1) {
848397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
849397b6df1SKris Buschelman   } else {
850397b6df1SKris Buschelman     ierr = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
851397b6df1SKris Buschelman   }
852f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
853397b6df1SKris Buschelman   PetscFunctionReturn(0);
854397b6df1SKris Buschelman }
855397b6df1SKris Buschelman EXTERN_C_END
856