xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 25f9c88cff50e825c0ea89c3ecd71095ccaf5bff)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
3397b6df1SKris Buschelman /*
4c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
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]
24a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
25397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
26adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
27397b6df1SKris Buschelman 
28397b6df1SKris Buschelman typedef struct {
29397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
30397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
31397b6df1SKris Buschelman #else
32397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
33397b6df1SKris Buschelman #endif
34397b6df1SKris Buschelman   MatStructure   matstruc;
35c1490034SHong Zhang   PetscMPIInt    myid,size;
36329ec9b3SHong Zhang   PetscInt       *irn,*jcn,sym,nSolve;
37397b6df1SKris Buschelman   PetscScalar    *val;
38397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
39329ec9b3SHong Zhang   VecScatter     scat_rhs, scat_sol;
40c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
41329ec9b3SHong Zhang   Vec            b_seq,x_seq;
426849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
436849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
446849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
456849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
466849ba73SBarry Smith   PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
476849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
486849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
496849ba73SBarry Smith   PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
50f0c56d0fSKris Buschelman } Mat_MUMPS;
51f0c56d0fSKris Buschelman 
52dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
53892f6c3fSKris Buschelman EXTERN_C_BEGIN
5475179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*);
55892f6c3fSKris Buschelman EXTERN_C_END
56397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
57397b6df1SKris Buschelman /*
58397b6df1SKris Buschelman   input:
5975747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
60397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
61397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
62397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
63397b6df1SKris Buschelman   output:
64397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
65397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
66397b6df1SKris Buschelman  */
67dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
68c1490034SHong Zhang   PetscInt       *ai, *aj, *bi, *bj, rstart,nz, *garray;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
702a4c71feSBarry Smith   PetscInt       i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol;
71c1490034SHong Zhang   PetscInt       *row,*col;
72397b6df1SKris Buschelman   PetscScalar    *av, *bv,*val;
73f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
74397b6df1SKris Buschelman 
75397b6df1SKris Buschelman   PetscFunctionBegin;
76397b6df1SKris Buschelman   if (mumps->isAIJ){
77397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
78397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
79397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
80397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
812a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
82397b6df1SKris Buschelman     garray = mat->garray;
83397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
84397b6df1SKris Buschelman 
85397b6df1SKris Buschelman   } else {
86397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
87397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
88397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
890c0e133fSBarry Smith     if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs);
906c6c5352SBarry Smith     nz = aa->nz + bb->nz;
912a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
92397b6df1SKris Buschelman     garray = mat->garray;
93397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
94397b6df1SKris Buschelman   }
95397b6df1SKris Buschelman 
96397b6df1SKris Buschelman   if (!valOnly){
977c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr);
987c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr);
99397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
100397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
101397b6df1SKris Buschelman   } else {
102397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
103397b6df1SKris Buschelman   }
104397b6df1SKris Buschelman   *nnz = nz;
105397b6df1SKris Buschelman 
106028e57e8SHong Zhang   jj = 0; irow = rstart;
107397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
108397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
109397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
110397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
111397b6df1SKris Buschelman     bjj = bj + bi[i];
112397b6df1SKris Buschelman 
113397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
114397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11575747be1SHong Zhang     j=-1;
11675747be1SHong Zhang     do {
11775747be1SHong Zhang       j++;
11875747be1SHong Zhang       if (j == countB) break;
119397b6df1SKris Buschelman       jcol = garray[bjj[j]];
12075747be1SHong Zhang     } while (jcol < colA_start);
12175747be1SHong Zhang     jB = j;
122397b6df1SKris Buschelman 
123397b6df1SKris Buschelman     /* B-part, smaller col index */
124397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
125397b6df1SKris Buschelman     for (j=0; j<jB; j++){
126397b6df1SKris Buschelman       jcol = garray[bjj[j]];
127397b6df1SKris Buschelman       if (!valOnly){
128397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12975747be1SHong Zhang 
130397b6df1SKris Buschelman       }
131397b6df1SKris Buschelman       val[jj++] = *bv++;
132397b6df1SKris Buschelman     }
133397b6df1SKris Buschelman     /* A-part */
134397b6df1SKris Buschelman     for (j=0; j<countA; j++){
135397b6df1SKris Buschelman       if (!valOnly){
136397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
137397b6df1SKris Buschelman       }
138397b6df1SKris Buschelman       val[jj++] = *av++;
139397b6df1SKris Buschelman     }
140397b6df1SKris Buschelman     /* B-part, larger col index */
141397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
142397b6df1SKris Buschelman       if (!valOnly){
143397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
144397b6df1SKris Buschelman       }
145397b6df1SKris Buschelman       val[jj++] = *bv++;
146397b6df1SKris Buschelman     }
147397b6df1SKris Buschelman     irow++;
148397b6df1SKris Buschelman   }
149397b6df1SKris Buschelman 
150397b6df1SKris Buschelman   PetscFunctionReturn(0);
151397b6df1SKris Buschelman }
152397b6df1SKris Buschelman 
153c338a77dSKris Buschelman EXTERN_C_BEGIN
154c338a77dSKris Buschelman #undef __FUNCT__
155c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
15617667f90SBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat)
157521d7252SBarry Smith {
158dfbe8321SBarry Smith   PetscErrorCode ierr;
159c338a77dSKris Buschelman   Mat            B=*newmat;
160f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
161901853e0SKris Buschelman   void           (*f)(void);
162c338a77dSKris Buschelman 
163c338a77dSKris Buschelman   PetscFunctionBegin;
164ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
165c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
166c338a77dSKris Buschelman   }
167f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
168f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
169f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
170f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
171f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
172f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
173901853e0SKris Buschelman 
17417667f90SBarry Smith   /* put back original composed preallocation function */
175f68b968cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
176901853e0SKris Buschelman   if (f) {
177f68b968cSBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr);
178901853e0SKris Buschelman   }
179901853e0SKris Buschelman 
180901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
181901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
182901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
183901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
1842895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
1852895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
186901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
187901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
188901853e0SKris Buschelman 
189901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
190c338a77dSKris Buschelman   *newmat = B;
191c338a77dSKris Buschelman   PetscFunctionReturn(0);
192c338a77dSKris Buschelman }
193c338a77dSKris Buschelman EXTERN_C_END
194c338a77dSKris Buschelman 
195397b6df1SKris Buschelman #undef __FUNCT__
1963924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
197dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
198dfbe8321SBarry Smith {
199f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
200dfbe8321SBarry Smith   PetscErrorCode ierr;
201c1490034SHong Zhang   PetscMPIInt    size=lu->size;
2026849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
203397b6df1SKris Buschelman   PetscFunctionBegin;
204397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
205397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
206329ec9b3SHong Zhang     if (size > 1){
207329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
208329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
209329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
210329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
211329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
212329ec9b3SHong Zhang       ierr = PetscFree(lu->val);CHKERRQ(ierr);
213329ec9b3SHong Zhang     }
214397b6df1SKris Buschelman     lu->id.job=JOB_END;
215397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
216397b6df1SKris Buschelman     zmumps_c(&lu->id);
217397b6df1SKris Buschelman #else
218397b6df1SKris Buschelman     dmumps_c(&lu->id);
219397b6df1SKris Buschelman #endif
220c338a77dSKris Buschelman     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
221c338a77dSKris Buschelman     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
222397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
223397b6df1SKris Buschelman   }
224a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
225a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
226c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
227397b6df1SKris Buschelman   PetscFunctionReturn(0);
228397b6df1SKris Buschelman }
229397b6df1SKris Buschelman 
230397b6df1SKris Buschelman #undef __FUNCT__
231a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
232dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
233dfbe8321SBarry Smith {
2346849ba73SBarry Smith   PetscErrorCode ierr;
235329ec9b3SHong Zhang   PetscMPIInt    size;
236a39386dcSKris Buschelman 
237a39386dcSKris Buschelman   PetscFunctionBegin;
238a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
239a39386dcSKris Buschelman   if (size==1) {
240ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
241a39386dcSKris Buschelman   } else {
242ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
243a39386dcSKris Buschelman   }
244a39386dcSKris Buschelman   PetscFunctionReturn(0);
245a39386dcSKris Buschelman }
246a39386dcSKris Buschelman 
247a39386dcSKris Buschelman #undef __FUNCT__
248a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
249dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
250dfbe8321SBarry Smith {
2516849ba73SBarry Smith   PetscErrorCode ierr;
252329ec9b3SHong Zhang   PetscMPIInt    size;
253a39386dcSKris Buschelman 
254a39386dcSKris Buschelman   PetscFunctionBegin;
255a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
256a39386dcSKris Buschelman   if (size==1) {
257ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
258a39386dcSKris Buschelman   } else {
259ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
260a39386dcSKris Buschelman   }
261a39386dcSKris Buschelman   PetscFunctionReturn(0);
262a39386dcSKris Buschelman }
263a39386dcSKris Buschelman 
264a39386dcSKris Buschelman #undef __FUNCT__
265c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
266dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
267f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
268dfbe8321SBarry Smith   PetscErrorCode ierr;
269397b6df1SKris Buschelman 
270397b6df1SKris Buschelman   PetscFunctionBegin;
271c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
272c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
273c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
27468533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
27568533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
27668533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
277c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
278c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
279c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
280c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
28168533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
282c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
283c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
284c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
285958c9bccSBarry Smith   if (!lu->myid && lu->id.ICNTL(11)>0) {
286c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
287c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
288c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
289c338a77dSKris 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);
290c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
291c338a77dSKris 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);
292c338a77dSKris Buschelman 
293c338a77dSKris Buschelman   }
294c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
295c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
296adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
29768533f03SHong Zhang   /* ICNTL(15-17) not used */
298c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
29968533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
30068533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
30168533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
302c338a77dSKris Buschelman 
303c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
304c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
305c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
30668533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
30757f0c58bSHong Zhang 
30857f0c58bSHong Zhang   /* infomation local to each processor */
30968533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
31057f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
31157f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
31268533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
31357f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
31457f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
31568533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
31657f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
31757f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
31868533f03SHong Zhang   /*
31968533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);}
32068533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr);
32168533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
32268533f03SHong Zhang   */
32368533f03SHong Zhang 
32468533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);}
32568533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
32668533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
32768533f03SHong Zhang 
32868533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);}
32968533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
33068533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
33168533f03SHong Zhang 
33268533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
33368533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
33468533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
335adc1d99fSHong Zhang 
336958c9bccSBarry Smith   if (!lu->myid){ /* information from the host */
337adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
338adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
339adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
340adc1d99fSHong Zhang 
341adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
342adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
343adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
344adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
345adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
346adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
34768533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
34868533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
34968533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
350adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
351adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
352adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
353adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
35468533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr);
355adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",lu->id.INFOG(17));CHKERRQ(ierr);
356adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",lu->id.INFOG(18));CHKERRQ(ierr);
357adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",lu->id.INFOG(19));CHKERRQ(ierr);
358adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
35968533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr);
36068533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr);
36168533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
36268533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
36368533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
364adc1d99fSHong Zhang   }
365adc1d99fSHong Zhang 
366397b6df1SKris Buschelman   PetscFunctionReturn(0);
367397b6df1SKris Buschelman }
368397b6df1SKris Buschelman 
369397b6df1SKris Buschelman #undef __FUNCT__
370c1490034SHong Zhang #define __FUNCT__ "MatView_MUMPS"
371c1490034SHong Zhang PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
372dfbe8321SBarry Smith   PetscErrorCode    ierr;
37332077d6dSBarry Smith   PetscTruth        iascii;
374397b6df1SKris Buschelman   PetscViewerFormat format;
375f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
376397b6df1SKris Buschelman 
377397b6df1SKris Buschelman   PetscFunctionBegin;
378397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
379397b6df1SKris Buschelman 
38032077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
38132077d6dSBarry Smith   if (iascii) {
382397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3832f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO){
384397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
385397b6df1SKris Buschelman     }
386397b6df1SKris Buschelman   }
387397b6df1SKris Buschelman   PetscFunctionReturn(0);
388397b6df1SKris Buschelman }
389397b6df1SKris Buschelman 
390397b6df1SKris Buschelman #undef __FUNCT__
391f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
392dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
393f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
394d54de34fSKris Buschelman   PetscScalar    *array;
395397b6df1SKris Buschelman   Vec            x_seq;
396329ec9b3SHong Zhang   IS             is_iden,is_petsc;
397329ec9b3SHong Zhang   VecScatter     scat_rhs=lu->scat_rhs,scat_sol=lu->scat_sol;
398dfbe8321SBarry Smith   PetscErrorCode ierr;
399329ec9b3SHong Zhang   PetscInt       i;
400397b6df1SKris Buschelman 
401397b6df1SKris Buschelman   PetscFunctionBegin;
402329ec9b3SHong Zhang   lu->id.nrhs = 1;
403329ec9b3SHong Zhang   x_seq = lu->b_seq;
404397b6df1SKris Buschelman   if (lu->size > 1){
405329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
406329ec9b3SHong Zhang     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
407329ec9b3SHong Zhang     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
408397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
409397b6df1SKris Buschelman   } else {  /* size == 1 */
410397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
411397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
412397b6df1SKris Buschelman   }
413397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
414397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
415397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
416397b6df1SKris Buschelman #else
417397b6df1SKris Buschelman     lu->id.rhs = array;
418397b6df1SKris Buschelman #endif
419397b6df1SKris Buschelman   }
420329ec9b3SHong Zhang   if (lu->size == 1){
421329ec9b3SHong Zhang     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
422329ec9b3SHong Zhang   } else if (!lu->myid){
423329ec9b3SHong Zhang     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
424329ec9b3SHong Zhang   }
425329ec9b3SHong Zhang 
426329ec9b3SHong Zhang   if (lu->size > 1){
427329ec9b3SHong Zhang     /* distributed solution */
428329ec9b3SHong Zhang     lu->id.ICNTL(21) = 1;
429329ec9b3SHong Zhang     if (!lu->nSolve){
430329ec9b3SHong Zhang       /* Create x_seq=sol_loc for repeated use */
431329ec9b3SHong Zhang       PetscInt    lsol_loc;
432329ec9b3SHong Zhang       PetscScalar *sol_loc;
433329ec9b3SHong Zhang       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
434329ec9b3SHong Zhang       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
435329ec9b3SHong Zhang       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
436329ec9b3SHong Zhang       lu->id.lsol_loc = lsol_loc;
437ea5584a4SSatish Balay       lu->id.sol_loc  = (F_DOUBLE *)sol_loc;
438329ec9b3SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
439329ec9b3SHong Zhang     }
440329ec9b3SHong Zhang   }
441397b6df1SKris Buschelman 
442397b6df1SKris Buschelman   /* solve phase */
443329ec9b3SHong Zhang   /*-------------*/
444397b6df1SKris Buschelman   lu->id.job = 3;
445397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
446397b6df1SKris Buschelman   zmumps_c(&lu->id);
447397b6df1SKris Buschelman #else
448397b6df1SKris Buschelman   dmumps_c(&lu->id);
449397b6df1SKris Buschelman #endif
450397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
45179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
452397b6df1SKris Buschelman   }
453397b6df1SKris Buschelman 
454329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
455329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
456329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
457329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
458329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
459397b6df1SKris Buschelman       }
460329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
461329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
462329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
463329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
464397b6df1SKris Buschelman     }
465329ec9b3SHong Zhang     ierr = VecScatterBegin(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
466329ec9b3SHong Zhang     ierr = VecScatterEnd(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
467329ec9b3SHong Zhang   }
468329ec9b3SHong Zhang   lu->nSolve++;
469397b6df1SKris Buschelman   PetscFunctionReturn(0);
470397b6df1SKris Buschelman }
471397b6df1SKris Buschelman 
472a58c3f20SHong Zhang /*
473a58c3f20SHong Zhang   input:
474a58c3f20SHong Zhang    F:        numeric factor
475a58c3f20SHong Zhang   output:
476a58c3f20SHong Zhang    nneg:     total number of negative pivots
477a58c3f20SHong Zhang    nzero:    0
478a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
479a58c3f20SHong Zhang */
480a58c3f20SHong Zhang 
481a58c3f20SHong Zhang #undef __FUNCT__
482a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
483dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
484a58c3f20SHong Zhang {
485a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
486dfbe8321SBarry Smith   PetscErrorCode ierr;
487c1490034SHong Zhang   PetscMPIInt    size;
488a58c3f20SHong Zhang 
489a58c3f20SHong Zhang   PetscFunctionBegin;
490bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
491bcb30aebSHong Zhang   /* MUMPS 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK */
492bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
49379a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
494bcb30aebSHong Zhang   }
495a58c3f20SHong Zhang   if (nneg){
496a58c3f20SHong Zhang     if (!lu->myid){
497a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
498a58c3f20SHong Zhang     }
499bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
500a58c3f20SHong Zhang   }
501a58c3f20SHong Zhang   if (nzero) *nzero = 0;
5022a4c71feSBarry Smith   if (npos)  *npos  = F->rmap.N - (*nneg);
503a58c3f20SHong Zhang   PetscFunctionReturn(0);
504a58c3f20SHong Zhang }
505a58c3f20SHong Zhang 
506397b6df1SKris Buschelman #undef __FUNCT__
50719facb7aSBarry Smith #define __FUNCT__ "MatFactorNumeric_AIJMUMPS"
508af281ebdSHong Zhang PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F)
509af281ebdSHong Zhang {
510f0c56d0fSKris Buschelman   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
511f0c56d0fSKris Buschelman   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
5126849ba73SBarry Smith   PetscErrorCode ierr;
5132a4c71feSBarry Smith   PetscInt       rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
514397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
515e09efc27SHong Zhang   Mat            F_diag;
516397b6df1SKris Buschelman 
517397b6df1SKris Buschelman   PetscFunctionBegin;
518397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
519f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
520397b6df1SKris Buschelman 
521397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
522397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
523397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
52475747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
525397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
526397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
527a0e2756fSSatish Balay     ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr);
528397b6df1SKris Buschelman 
529397b6df1SKris Buschelman     /* Set mumps options */
530397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
531397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
532397b6df1SKris Buschelman     lu->id.sym=lu->sym;
533397b6df1SKris Buschelman     if (lu->sym == 2){
534397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
535397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
536397b6df1SKris Buschelman     }
537397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
538397b6df1SKris Buschelman     zmumps_c(&lu->id);
539397b6df1SKris Buschelman #else
540397b6df1SKris Buschelman     dmumps_c(&lu->id);
541397b6df1SKris Buschelman #endif
542397b6df1SKris Buschelman 
543397b6df1SKris Buschelman     if (lu->size == 1){
544397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
545397b6df1SKris Buschelman     } else {
546397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
547397b6df1SKris Buschelman     }
548397b6df1SKris Buschelman 
549397b6df1SKris Buschelman     icntl=-1;
550397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
55119facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
552397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
553397b6df1SKris Buschelman     } else { /* no output */
554397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
555397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
556397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
557397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
558397b6df1SKris Buschelman     }
559397b6df1SKris 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);
560397b6df1SKris Buschelman     icntl=-1;
561397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
562397b6df1SKris Buschelman     if (flg) {
563397b6df1SKris Buschelman       if (icntl== 1){
564397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
565397b6df1SKris Buschelman       } else {
566397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
567397b6df1SKris Buschelman       }
568397b6df1SKris Buschelman     }
569397b6df1SKris 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);
570397b6df1SKris 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);
57194b7f48cSBarry Smith     ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): error analysis, a positive value returns statistics (by -ksp_view)","None",lu->id.ICNTL(11),&lu->id.ICNTL(11),PETSC_NULL);CHKERRQ(ierr);
572397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
573397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
574adc1d99fSHong Zhang     ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage of estimated workspace increase","None",lu->id.ICNTL(14),&lu->id.ICNTL(14),PETSC_NULL);CHKERRQ(ierr);
575397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
576397b6df1SKris Buschelman 
577397b6df1SKris 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);
578397b6df1SKris 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);
579397b6df1SKris 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);
580*25f9c88cSHong Zhang     ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",lu->id.CNTL(4),&lu->id.CNTL(4),PETSC_NULL);CHKERRQ(ierr);
581397b6df1SKris Buschelman     PetscOptionsEnd();
582397b6df1SKris Buschelman   }
583397b6df1SKris Buschelman 
584397b6df1SKris Buschelman   /* define matrix A */
585397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
586397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
587397b6df1SKris Buschelman     if (!lu->myid) {
588c36ead0aSKris Buschelman       if (lua->isAIJ){
589397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
590397b6df1SKris Buschelman         nz               = aa->nz;
591397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
592397b6df1SKris Buschelman       } else {
593397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5946c6c5352SBarry Smith         nz                  =  aa->nz;
595397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
596397b6df1SKris Buschelman       }
597397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
5987c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
5997c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
600397b6df1SKris Buschelman         nz = 0;
601397b6df1SKris Buschelman         for (i=0; i<M; i++){
602397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
603397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
604397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
605397b6df1SKris Buschelman           }
606397b6df1SKris Buschelman         }
607397b6df1SKris Buschelman       }
608397b6df1SKris Buschelman     }
609397b6df1SKris Buschelman     break;
610397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
611397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
612397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
613397b6df1SKris Buschelman     } else {
614397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
615397b6df1SKris Buschelman     }
616397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
617397b6df1SKris Buschelman     break;
618397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
619397b6df1SKris Buschelman   }
620397b6df1SKris Buschelman 
621397b6df1SKris Buschelman   /* analysis phase */
622329ec9b3SHong Zhang   /*----------------*/
623397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
624329ec9b3SHong Zhang     lu->id.job = 1;
625329ec9b3SHong Zhang 
626397b6df1SKris Buschelman     lu->id.n = M;
627397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
628397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
629397b6df1SKris Buschelman       if (!lu->myid) {
630397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
631397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
632397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
633397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
634397b6df1SKris Buschelman #else
635397b6df1SKris Buschelman           lu->id.a = lu->val;
636397b6df1SKris Buschelman #endif
637397b6df1SKris Buschelman         }
638397b6df1SKris Buschelman       }
639397b6df1SKris Buschelman       break;
640397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
641397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
642397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
643397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
644397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
645397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
646397b6df1SKris Buschelman #else
647397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
648397b6df1SKris Buschelman #endif
649397b6df1SKris Buschelman       }
650329ec9b3SHong Zhang       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
651329ec9b3SHong Zhang       IS  is_iden;
652329ec9b3SHong Zhang       Vec b;
653329ec9b3SHong Zhang       if (!lu->myid){
654329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr);
655329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr);
656329ec9b3SHong Zhang       } else {
657329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
658329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
659329ec9b3SHong Zhang       }
660329ec9b3SHong Zhang       ierr = VecCreate(A->comm,&b);CHKERRQ(ierr);
661329ec9b3SHong Zhang       ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
662329ec9b3SHong Zhang       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
663329ec9b3SHong Zhang 
664329ec9b3SHong Zhang       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
665329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
666329ec9b3SHong Zhang       ierr = VecDestroy(b);CHKERRQ(ierr);
667397b6df1SKris Buschelman       break;
668397b6df1SKris Buschelman     }
669397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
670397b6df1SKris Buschelman     zmumps_c(&lu->id);
671397b6df1SKris Buschelman #else
672397b6df1SKris Buschelman     dmumps_c(&lu->id);
673397b6df1SKris Buschelman #endif
674397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
67579a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
676397b6df1SKris Buschelman     }
677397b6df1SKris Buschelman   }
678397b6df1SKris Buschelman 
679397b6df1SKris Buschelman   /* numerical factorization phase */
680329ec9b3SHong Zhang   /*-------------------------------*/
681329ec9b3SHong Zhang   lu->id.job = 2;
682958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
683a7aca84bSHong Zhang     if (!lu->myid) {
684397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
685397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
686397b6df1SKris Buschelman #else
687397b6df1SKris Buschelman       lu->id.a = lu->val;
688397b6df1SKris Buschelman #endif
689397b6df1SKris Buschelman     }
690397b6df1SKris Buschelman   } else {
691397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
692397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
693397b6df1SKris Buschelman #else
694397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
695397b6df1SKris Buschelman #endif
696397b6df1SKris Buschelman   }
697397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
698397b6df1SKris Buschelman   zmumps_c(&lu->id);
699397b6df1SKris Buschelman #else
700397b6df1SKris Buschelman   dmumps_c(&lu->id);
701397b6df1SKris Buschelman #endif
702397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
70319facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
70419facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
70519facb7aSBarry Smith     } else {
70679a5c55eSBarry Smith       SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
707397b6df1SKris Buschelman     }
70819facb7aSBarry Smith   }
709397b6df1SKris Buschelman 
71019facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
71179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
712397b6df1SKris Buschelman   }
713397b6df1SKris Buschelman 
7148ada1bb4SHong Zhang   if (lu->size > 1){
715e09efc27SHong Zhang     if ((*F)->factor == FACTOR_LU){
716e09efc27SHong Zhang       F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
717e09efc27SHong Zhang     } else {
718e09efc27SHong Zhang       F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
719e09efc27SHong Zhang     }
720e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
721329ec9b3SHong Zhang     if (lu->nSolve){
722329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
723329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
724329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
725329ec9b3SHong Zhang     }
7268ada1bb4SHong Zhang   }
727397b6df1SKris Buschelman   (*F)->assembled   = PETSC_TRUE;
728397b6df1SKris Buschelman   lu->matstruc      = SAME_NONZERO_PATTERN;
729ace87b0dSHong Zhang   lu->CleanUpMUMPS  = PETSC_TRUE;
730329ec9b3SHong Zhang   lu->nSolve        = 0;
731397b6df1SKris Buschelman   PetscFunctionReturn(0);
732397b6df1SKris Buschelman }
733397b6df1SKris Buschelman 
734397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
735397b6df1SKris Buschelman #undef __FUNCT__
736f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
737dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
738397b6df1SKris Buschelman   Mat            B;
739f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
740dfbe8321SBarry Smith   PetscErrorCode ierr;
741397b6df1SKris Buschelman 
742397b6df1SKris Buschelman   PetscFunctionBegin;
743397b6df1SKris Buschelman   /* Create the factorization matrix */
744f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
7452a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
746be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
747397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
748397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
749397b6df1SKris Buschelman 
750f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
751397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
752f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
753397b6df1SKris Buschelman   lu->sym                 = 0;
754397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
755397b6df1SKris Buschelman 
756397b6df1SKris Buschelman   *F = B;
757397b6df1SKris Buschelman   PetscFunctionReturn(0);
758397b6df1SKris Buschelman }
759397b6df1SKris Buschelman 
760397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
761397b6df1SKris Buschelman #undef __FUNCT__
762f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
763dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
764397b6df1SKris Buschelman   Mat            B;
765f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
766dfbe8321SBarry Smith   PetscErrorCode ierr;
767397b6df1SKris Buschelman 
768397b6df1SKris Buschelman   PetscFunctionBegin;
769397b6df1SKris Buschelman   /* Create the factorization matrix */
770f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
7712a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
772be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
773efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
774efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
775397b6df1SKris Buschelman 
776f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
777a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
778397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
779f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
780397b6df1SKris Buschelman   lu->sym                       = 2;
781397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
782397b6df1SKris Buschelman 
783397b6df1SKris Buschelman   *F = B;
784397b6df1SKris Buschelman   PetscFunctionReturn(0);
785397b6df1SKris Buschelman }
786397b6df1SKris Buschelman 
787397b6df1SKris Buschelman #undef __FUNCT__
788f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
789dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
790dfbe8321SBarry Smith   PetscErrorCode ierr;
791f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
792c338a77dSKris Buschelman 
793397b6df1SKris Buschelman   PetscFunctionBegin;
794c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
795f0c56d0fSKris Buschelman 
796c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
797c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
798f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
799397b6df1SKris Buschelman   PetscFunctionReturn(0);
800397b6df1SKris Buschelman }
801397b6df1SKris Buschelman 
802c338a77dSKris Buschelman EXTERN_C_BEGIN
803c338a77dSKris Buschelman #undef __FUNCT__
804f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
80574a14a02SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
806dfbe8321SBarry Smith {
807dfbe8321SBarry Smith   PetscErrorCode ierr;
808521d7252SBarry Smith   PetscMPIInt    size;
809c338a77dSKris Buschelman   MPI_Comm       comm;
810c338a77dSKris Buschelman   Mat            B=*newmat;
811f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps;
812397b6df1SKris Buschelman 
813397b6df1SKris Buschelman   PetscFunctionBegin;
814ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
815c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
816397b6df1SKris Buschelman   }
817397b6df1SKris Buschelman 
818c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
819f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
820c338a77dSKris Buschelman 
821f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
822c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
823c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
824c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
825c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
826c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
827a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
828c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
829f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
830c338a77dSKris Buschelman 
8314b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
832422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
833c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
834f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
835f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
8363924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
837c338a77dSKris Buschelman 
838c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
839c338a77dSKris Buschelman   if (size == 1) {
840c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
841f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
842c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
843c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
844c338a77dSKris Buschelman   } else {
845c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
846f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
847c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
848c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
849c338a77dSKris Buschelman   }
850c338a77dSKris Buschelman 
851ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr);
852c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
853c338a77dSKris Buschelman   *newmat = B;
854397b6df1SKris Buschelman   PetscFunctionReturn(0);
855397b6df1SKris Buschelman }
856c338a77dSKris Buschelman EXTERN_C_END
857397b6df1SKris Buschelman 
85824b6179bSKris Buschelman /*MC
859fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
86024b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
86124b6179bSKris Buschelman 
86224b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
86324b6179bSKris Buschelman   on how to declare the existence of external packages),
86424b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
86524b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
86624b6179bSKris Buschelman 
86724b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
86824b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
86924b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
87024b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
87128b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
87228b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
87328b08bd3SKris Buschelman   without data copy.
87424b6179bSKris Buschelman 
87524b6179bSKris Buschelman   Options Database Keys:
8760bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
87724b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
87824b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
87924b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
88024b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
88124b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
88224b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
88394b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
88424b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
88524b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
88624b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
88724b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
88824b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
88924b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
89024b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
89124b6179bSKris Buschelman 
89224b6179bSKris Buschelman   Level: beginner
89324b6179bSKris Buschelman 
89424b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
89524b6179bSKris Buschelman M*/
89624b6179bSKris Buschelman 
897397b6df1SKris Buschelman EXTERN_C_BEGIN
898397b6df1SKris Buschelman #undef __FUNCT__
899f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
900be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A)
901dfbe8321SBarry Smith {
902dfbe8321SBarry Smith   PetscErrorCode ierr;
903c1490034SHong Zhang   PetscMPIInt    size;
904397b6df1SKris Buschelman 
905397b6df1SKris Buschelman   PetscFunctionBegin;
90617667f90SBarry Smith   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
907397b6df1SKris Buschelman   if (size == 1) {
908397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
909397b6df1SKris Buschelman   } else {
910397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
91100ff2a26SHong Zhang     /*
91200ff2a26SHong Zhang     Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
913ceb03754SKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
91400ff2a26SHong Zhang     */
915397b6df1SKris Buschelman   }
916ceb03754SKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
917397b6df1SKris Buschelman   PetscFunctionReturn(0);
918397b6df1SKris Buschelman }
919397b6df1SKris Buschelman EXTERN_C_END
920397b6df1SKris Buschelman 
921f579278aSKris Buschelman #undef __FUNCT__
922f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
923dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
924dfbe8321SBarry Smith {
925dfbe8321SBarry Smith   PetscErrorCode ierr;
926f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
927f579278aSKris Buschelman 
928f579278aSKris Buschelman   PetscFunctionBegin;
929f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
930f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
931f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
932f579278aSKris Buschelman   PetscFunctionReturn(0);
933f579278aSKris Buschelman }
934f579278aSKris Buschelman 
935f579278aSKris Buschelman EXTERN_C_BEGIN
936f579278aSKris Buschelman #undef __FUNCT__
9379c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
938c1490034SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
9399c097c71SKris Buschelman {
9409c097c71SKris Buschelman   Mat       A;
94102217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
942dfbe8321SBarry Smith   PetscErrorCode ierr;
9439c097c71SKris Buschelman 
9449c097c71SKris Buschelman   PetscFunctionBegin;
9459c097c71SKris Buschelman   /*
9469c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
9479c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
9489c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
9499c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
9509c097c71SKris Buschelman     in the preallocation routine.
9519c097c71SKris Buschelman   */
9529c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
9539c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
954ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
9559c097c71SKris Buschelman   PetscFunctionReturn(0);
9569c097c71SKris Buschelman }
9579c097c71SKris Buschelman EXTERN_C_END
9589c097c71SKris Buschelman 
9599c097c71SKris Buschelman EXTERN_C_BEGIN
9609c097c71SKris Buschelman #undef __FUNCT__
961f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
96275179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
963dfbe8321SBarry Smith {
964dfbe8321SBarry Smith   PetscErrorCode ierr;
965521d7252SBarry Smith   PetscMPIInt    size;
966f579278aSKris Buschelman   MPI_Comm       comm;
967f579278aSKris Buschelman   Mat            B=*newmat;
968422e82a1SHong Zhang   Mat_MUMPS      *mumps;
9699c097c71SKris Buschelman   void           (*f)(void);
970f579278aSKris Buschelman 
971f579278aSKris Buschelman   PetscFunctionBegin;
972ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
973f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
974f579278aSKris Buschelman   }
975f579278aSKris Buschelman 
976f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
977f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
978f579278aSKris Buschelman 
979f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
980f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
981f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
982f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
983f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
984f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
985a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
986f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
987f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
988f579278aSKris Buschelman 
989f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
990422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
991c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
992f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
993f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
9943924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
995f579278aSKris Buschelman 
996f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
997f579278aSKris Buschelman   if (size == 1) {
998f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
999f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1000f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
1001f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1002f579278aSKris Buschelman   } else {
10039c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
1004f68b968cSBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
1005901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
10066849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
10079c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
10089c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
10099c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
10109c097c71SKris Buschelman     }
1011f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
1012f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1013f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
1014f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1015f579278aSKris Buschelman   }
1016f579278aSKris Buschelman 
1017ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr);
1018f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
1019f579278aSKris Buschelman   *newmat = B;
1020f579278aSKris Buschelman   PetscFunctionReturn(0);
1021f579278aSKris Buschelman }
1022f579278aSKris Buschelman EXTERN_C_END
1023f579278aSKris Buschelman 
1024f0c56d0fSKris Buschelman #undef __FUNCT__
1025422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
1026dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1027dfbe8321SBarry Smith   PetscErrorCode ierr;
10288e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
10298f340917SKris Buschelman 
1030f0c56d0fSKris Buschelman   PetscFunctionBegin;
10318f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
10328e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
1033f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
1034f0c56d0fSKris Buschelman }
1035f0c56d0fSKris Buschelman 
103624b6179bSKris Buschelman /*MC
1037fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
103824b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
103924b6179bSKris Buschelman 
104024b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
104124b6179bSKris Buschelman   on how to declare the existence of external packages),
104224b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
104324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
104424b6179bSKris Buschelman 
104524b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
104624b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
104724b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
104824b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
104928b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
105028b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
105128b08bd3SKris Buschelman   without data copy.
105224b6179bSKris Buschelman 
105324b6179bSKris Buschelman   Options Database Keys:
10540bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
105524b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
105624b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
105724b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
105824b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
105924b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
106024b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
106194b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
106224b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
106324b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
106424b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
106524b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
106624b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
106724b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
106824b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
106924b6179bSKris Buschelman 
107024b6179bSKris Buschelman   Level: beginner
107124b6179bSKris Buschelman 
107224b6179bSKris Buschelman .seealso: MATAIJMUMPS
107324b6179bSKris Buschelman M*/
107424b6179bSKris Buschelman 
1075397b6df1SKris Buschelman EXTERN_C_BEGIN
1076397b6df1SKris Buschelman #undef __FUNCT__
1077f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1078be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A)
1079dfbe8321SBarry Smith {
10806849ba73SBarry Smith   PetscErrorCode ierr;
10818ada1bb4SHong Zhang   PetscMPIInt    size;
1082397b6df1SKris Buschelman 
1083397b6df1SKris Buschelman   PetscFunctionBegin;
10845441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1085397b6df1SKris Buschelman   if (size == 1) {
1086397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1087397b6df1SKris Buschelman   } else {
1088397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1089397b6df1SKris Buschelman   }
1090ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1091397b6df1SKris Buschelman   PetscFunctionReturn(0);
1092397b6df1SKris Buschelman }
1093397b6df1SKris Buschelman EXTERN_C_END
1094