xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision e2d9671b86a8697a4c2c5fff9ed2002d1fb085ae)
1397b6df1SKris Buschelman /*$Id: mumps.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2397b6df1SKris Buschelman /*
3a58c3f20SHong Zhang     Provides an interface to the MUMPS_4.3.1 sparse solver
4397b6df1SKris Buschelman */
5397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
6397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h"
7397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h"
8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
12397b6df1SKris Buschelman #include "zmumps_c.h"
13397b6df1SKris Buschelman #else
14397b6df1SKris Buschelman #include "dmumps_c.h"
15397b6df1SKris Buschelman #endif
16397b6df1SKris Buschelman EXTERN_C_END
17397b6df1SKris Buschelman #define JOB_INIT -1
18397b6df1SKris Buschelman #define JOB_END -2
19397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
20397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
21397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
22397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
23a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
24397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
25adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
26397b6df1SKris Buschelman 
27397b6df1SKris Buschelman typedef struct {
28397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
29397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
30397b6df1SKris Buschelman #else
31397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
32397b6df1SKris Buschelman #endif
33397b6df1SKris Buschelman   MatStructure   matstruc;
34397b6df1SKris Buschelman   int            myid,size,*irn,*jcn,sym;
35397b6df1SKris Buschelman   PetscScalar    *val;
36397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
37397b6df1SKris Buschelman 
38c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
39f0c56d0fSKris Buschelman   int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
40c338a77dSKris Buschelman   int (*MatView)(Mat,PetscViewer);
41c338a77dSKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
42c338a77dSKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
43c338a77dSKris Buschelman   int (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
44c338a77dSKris Buschelman   int (*MatDestroy)(Mat);
45a39386dcSKris Buschelman   int (*specialdestroy)(Mat);
46f0c56d0fSKris Buschelman } Mat_MUMPS;
47f0c56d0fSKris Buschelman 
48f0c56d0fSKris Buschelman EXTERN int MatDuplicate_AIJMUMPS(Mat,MatDuplicateOption,Mat*);
49f0c56d0fSKris Buschelman EXTERN int MatDuplicate_SBAIJMUMPS(Mat,MatDuplicateOption,Mat*);
500e3434eeSKris Buschelman 
51397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
52397b6df1SKris Buschelman /*
53397b6df1SKris Buschelman   input:
5475747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
55397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
56397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
57397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
58397b6df1SKris Buschelman   output:
59397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
60397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
61397b6df1SKris Buschelman  */
62f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
63397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
64397b6df1SKris Buschelman   int         ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
65d54de34fSKris Buschelman   int         *row,*col;
66397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
67f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
68397b6df1SKris Buschelman 
69397b6df1SKris Buschelman   PetscFunctionBegin;
70397b6df1SKris Buschelman   if (mumps->isAIJ){
71397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
72397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
73397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
74397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
75397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
76397b6df1SKris Buschelman     garray = mat->garray;
77397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
78397b6df1SKris Buschelman 
79397b6df1SKris Buschelman   } else {
80397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
81397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
82397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
83847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs);
846c6c5352SBarry Smith     nz = aa->nz + bb->nz;
85397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
86397b6df1SKris Buschelman     garray = mat->garray;
87397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
88397b6df1SKris Buschelman   }
89397b6df1SKris Buschelman 
90397b6df1SKris Buschelman   if (!valOnly){
91397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
92397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
93397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
94397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
95397b6df1SKris Buschelman   } else {
96397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
97397b6df1SKris Buschelman   }
98397b6df1SKris Buschelman   *nnz = nz;
99397b6df1SKris Buschelman 
100028e57e8SHong Zhang   jj = 0; irow = rstart;
101397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
102397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
103397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
104397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
105397b6df1SKris Buschelman     bjj = bj + bi[i];
106397b6df1SKris Buschelman 
107397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
108397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
10975747be1SHong Zhang     j=-1;
11075747be1SHong Zhang     do {
11175747be1SHong Zhang       j++;
11275747be1SHong Zhang       if (j == countB) break;
113397b6df1SKris Buschelman       jcol = garray[bjj[j]];
11475747be1SHong Zhang     } while (jcol < colA_start);
11575747be1SHong Zhang     jB = j;
116397b6df1SKris Buschelman 
117397b6df1SKris Buschelman     /* B-part, smaller col index */
118397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
119397b6df1SKris Buschelman     for (j=0; j<jB; j++){
120397b6df1SKris Buschelman       jcol = garray[bjj[j]];
121397b6df1SKris Buschelman       if (!valOnly){
122397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12375747be1SHong Zhang 
124397b6df1SKris Buschelman       }
125397b6df1SKris Buschelman       val[jj++] = *bv++;
126397b6df1SKris Buschelman     }
127397b6df1SKris Buschelman     /* A-part */
128397b6df1SKris Buschelman     for (j=0; j<countA; j++){
129397b6df1SKris Buschelman       if (!valOnly){
130397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
131397b6df1SKris Buschelman       }
132397b6df1SKris Buschelman       val[jj++] = *av++;
133397b6df1SKris Buschelman     }
134397b6df1SKris Buschelman     /* B-part, larger col index */
135397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
136397b6df1SKris Buschelman       if (!valOnly){
137397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
138397b6df1SKris Buschelman       }
139397b6df1SKris Buschelman       val[jj++] = *bv++;
140397b6df1SKris Buschelman     }
141397b6df1SKris Buschelman     irow++;
142397b6df1SKris Buschelman   }
143397b6df1SKris Buschelman 
144397b6df1SKris Buschelman   PetscFunctionReturn(0);
145397b6df1SKris Buschelman }
146397b6df1SKris Buschelman 
147c338a77dSKris Buschelman EXTERN_C_BEGIN
148c338a77dSKris Buschelman #undef __FUNCT__
149c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
1508e9aea5cSBarry Smith int MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) {
151c338a77dSKris Buschelman   int       ierr;
152c338a77dSKris Buschelman   Mat       B=*newmat;
153f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
154c338a77dSKris Buschelman 
155c338a77dSKris Buschelman   PetscFunctionBegin;
156c338a77dSKris Buschelman   if (B != A) {
157c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
158c338a77dSKris Buschelman   }
159f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
160f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
161f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
162f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
163f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
164f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
1653924e44cSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
166f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
167c338a77dSKris Buschelman   *newmat = B;
168c338a77dSKris Buschelman   PetscFunctionReturn(0);
169c338a77dSKris Buschelman }
170c338a77dSKris Buschelman EXTERN_C_END
171c338a77dSKris Buschelman 
172397b6df1SKris Buschelman #undef __FUNCT__
1733924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
1743924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) {
175f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
176c338a77dSKris Buschelman   int       ierr,size=lu->size;
177a39386dcSKris Buschelman   int       (*specialdestroy)(Mat);
178397b6df1SKris Buschelman   PetscFunctionBegin;
179397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
180397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
181397b6df1SKris Buschelman     lu->id.job=JOB_END;
182397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
183397b6df1SKris Buschelman     zmumps_c(&lu->id);
184397b6df1SKris Buschelman #else
185397b6df1SKris Buschelman     dmumps_c(&lu->id);
186397b6df1SKris Buschelman #endif
187c338a77dSKris Buschelman     if (lu->irn) {
188c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
189c338a77dSKris Buschelman     }
190c338a77dSKris Buschelman     if (lu->jcn) {
191c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
192c338a77dSKris Buschelman     }
193c338a77dSKris Buschelman     if (size>1 && lu->val) {
194c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
195c338a77dSKris Buschelman     }
196397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
197397b6df1SKris Buschelman   }
198a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
199a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
200c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
201397b6df1SKris Buschelman   PetscFunctionReturn(0);
202397b6df1SKris Buschelman }
203397b6df1SKris Buschelman 
204397b6df1SKris Buschelman #undef __FUNCT__
205a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
206a39386dcSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) {
207a39386dcSKris Buschelman   int ierr, size;
208a39386dcSKris Buschelman 
209a39386dcSKris Buschelman   PetscFunctionBegin;
210a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
211a39386dcSKris Buschelman   if (size==1) {
212a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr);
213a39386dcSKris Buschelman   } else {
214a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr);
215a39386dcSKris Buschelman   }
216a39386dcSKris Buschelman   PetscFunctionReturn(0);
217a39386dcSKris Buschelman }
218a39386dcSKris Buschelman 
219a39386dcSKris Buschelman #undef __FUNCT__
220a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
221a39386dcSKris Buschelman int MatDestroy_SBAIJMUMPS(Mat A) {
222a39386dcSKris Buschelman   int ierr, size;
223a39386dcSKris Buschelman 
224a39386dcSKris Buschelman   PetscFunctionBegin;
225a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
226a39386dcSKris Buschelman   if (size==1) {
227a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr);
228a39386dcSKris Buschelman   } else {
229a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr);
230a39386dcSKris Buschelman   }
231a39386dcSKris Buschelman   PetscFunctionReturn(0);
232a39386dcSKris Buschelman }
233a39386dcSKris Buschelman 
234a39386dcSKris Buschelman #undef __FUNCT__
235c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
236f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
237f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
238397b6df1SKris Buschelman   int       ierr;
239397b6df1SKris Buschelman 
240397b6df1SKris Buschelman   PetscFunctionBegin;
241c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
242c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
243c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
244c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
245c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
246c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
247c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
248c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
249c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
250c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
251c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
252c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
253c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
254c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
255c338a77dSKris 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);
256c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
257c338a77dSKris 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);
258c338a77dSKris Buschelman 
259c338a77dSKris Buschelman   }
260c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
261c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
262adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
263c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
264c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
265c338a77dSKris Buschelman 
266c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
267c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
268c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
26957f0c58bSHong Zhang 
27057f0c58bSHong Zhang   /* infomation local to each processor */
27157f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
27257f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
27357f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27457f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
27557f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
27657f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27757f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
27857f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
27957f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
280adc1d99fSHong Zhang 
281adc1d99fSHong Zhang   if (lu->myid == 0){ /* information from the host */
282adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
283adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
284adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
285adc1d99fSHong Zhang 
286adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
287adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
288adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
289adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
290adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
291adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
292adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
293adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
294adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
295adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
296adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
297adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
298adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
299adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(16) (estimated size (in million of bytes) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr);
300adc1d99fSHong 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);
301adc1d99fSHong 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);
302adc1d99fSHong 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);
303adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
304adc1d99fSHong Zhang   }
305adc1d99fSHong Zhang 
306397b6df1SKris Buschelman   PetscFunctionReturn(0);
307397b6df1SKris Buschelman }
308397b6df1SKris Buschelman 
309397b6df1SKris Buschelman #undef __FUNCT__
310f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
311f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
312397b6df1SKris Buschelman   int               ierr;
313397b6df1SKris Buschelman   PetscTruth        isascii;
314397b6df1SKris Buschelman   PetscViewerFormat format;
315f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
316397b6df1SKris Buschelman 
317397b6df1SKris Buschelman   PetscFunctionBegin;
318397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
319397b6df1SKris Buschelman 
320397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
321397b6df1SKris Buschelman   if (isascii) {
322397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
323397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
324397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
325397b6df1SKris Buschelman     }
326397b6df1SKris Buschelman   }
327397b6df1SKris Buschelman   PetscFunctionReturn(0);
328397b6df1SKris Buschelman }
329397b6df1SKris Buschelman 
330397b6df1SKris Buschelman #undef __FUNCT__
331f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
332f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
333f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
334d54de34fSKris Buschelman   PetscScalar *array;
335397b6df1SKris Buschelman   Vec         x_seq;
336397b6df1SKris Buschelman   IS          iden;
337397b6df1SKris Buschelman   VecScatter  scat;
338397b6df1SKris Buschelman   int         ierr;
339397b6df1SKris Buschelman 
340397b6df1SKris Buschelman   PetscFunctionBegin;
341397b6df1SKris Buschelman   if (lu->size > 1){
342397b6df1SKris Buschelman     if (!lu->myid){
343397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
344397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
345397b6df1SKris Buschelman     } else {
346397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
347397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
348397b6df1SKris Buschelman     }
349397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
350397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
351397b6df1SKris Buschelman 
352397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
353397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
354397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
355397b6df1SKris Buschelman   } else {  /* size == 1 */
356397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
357397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
358397b6df1SKris Buschelman   }
359397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
360397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
361397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
362397b6df1SKris Buschelman #else
363397b6df1SKris Buschelman     lu->id.rhs = array;
364397b6df1SKris Buschelman #endif
365397b6df1SKris Buschelman   }
366397b6df1SKris Buschelman 
367397b6df1SKris Buschelman   /* solve phase */
368397b6df1SKris Buschelman   lu->id.job=3;
369397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
370397b6df1SKris Buschelman   zmumps_c(&lu->id);
371397b6df1SKris Buschelman #else
372397b6df1SKris Buschelman   dmumps_c(&lu->id);
373397b6df1SKris Buschelman #endif
374397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
375397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
376397b6df1SKris Buschelman   }
377397b6df1SKris Buschelman 
378397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
379397b6df1SKris Buschelman   if (lu->size > 1) {
380397b6df1SKris Buschelman     if (!lu->myid){
381397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
382397b6df1SKris Buschelman     }
383397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
384397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
385397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
386397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
387397b6df1SKris Buschelman   } else {
388397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
389397b6df1SKris Buschelman   }
390397b6df1SKris Buschelman 
391397b6df1SKris Buschelman   PetscFunctionReturn(0);
392397b6df1SKris Buschelman }
393397b6df1SKris Buschelman 
394a58c3f20SHong Zhang /*
395a58c3f20SHong Zhang   input:
396a58c3f20SHong Zhang    F:        numeric factor
397a58c3f20SHong Zhang   output:
398a58c3f20SHong Zhang    nneg:     total number of negative pivots
399a58c3f20SHong Zhang    nzero:    0
400a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
401a58c3f20SHong Zhang */
402a58c3f20SHong Zhang 
403a58c3f20SHong Zhang #undef __FUNCT__
404a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
405a58c3f20SHong Zhang int MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
406a58c3f20SHong Zhang {
407a58c3f20SHong Zhang   Mat_MUMPS  *lu =(Mat_MUMPS*)F->spptr;
408bcb30aebSHong Zhang   int        ierr,neg,zero,pos,size;
409a58c3f20SHong Zhang 
410a58c3f20SHong Zhang   PetscFunctionBegin;
411bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
412bcb30aebSHong 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 */
413bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
414bcb30aebSHong Zhang     SETERRQ1(1,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
415bcb30aebSHong Zhang   }
416a58c3f20SHong Zhang   if (nneg){
417a58c3f20SHong Zhang     if (!lu->myid){
418a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
419a58c3f20SHong Zhang     }
420bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
421a58c3f20SHong Zhang   }
422a58c3f20SHong Zhang   if (nzero) *nzero = 0;
423a58c3f20SHong Zhang   if (npos)  *npos  = F->M - (*nneg);
424a58c3f20SHong Zhang   PetscFunctionReturn(0);
425a58c3f20SHong Zhang }
426a58c3f20SHong Zhang 
427397b6df1SKris Buschelman #undef __FUNCT__
428f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
429f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
430f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
431f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
432397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
433397b6df1SKris Buschelman   PetscTruth valOnly,flg;
434397b6df1SKris Buschelman 
435397b6df1SKris Buschelman   PetscFunctionBegin;
436397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
437f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
438397b6df1SKris Buschelman 
439397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
440397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
441397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
44275747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
443397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
444397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
445397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
446397b6df1SKris Buschelman 
447397b6df1SKris Buschelman     /* Set mumps options */
448397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
449397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
450397b6df1SKris Buschelman     lu->id.sym=lu->sym;
451397b6df1SKris Buschelman     if (lu->sym == 2){
452397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
453397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
454397b6df1SKris Buschelman     }
455397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
456397b6df1SKris Buschelman   zmumps_c(&lu->id);
457397b6df1SKris Buschelman #else
458397b6df1SKris Buschelman   dmumps_c(&lu->id);
459397b6df1SKris Buschelman #endif
460397b6df1SKris Buschelman 
461397b6df1SKris Buschelman     if (lu->size == 1){
462397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
463397b6df1SKris Buschelman     } else {
464397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
465397b6df1SKris Buschelman     }
466397b6df1SKris Buschelman 
467397b6df1SKris Buschelman     icntl=-1;
468397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
469397b6df1SKris Buschelman     if (flg && icntl > 0) {
470397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
471397b6df1SKris Buschelman     } else { /* no output */
472397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
473397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
474397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
475397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
476397b6df1SKris Buschelman     }
477397b6df1SKris 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);
478397b6df1SKris Buschelman     icntl=-1;
479397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
480397b6df1SKris Buschelman     if (flg) {
481397b6df1SKris Buschelman       if (icntl== 1){
482397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
483397b6df1SKris Buschelman       } else {
484397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
485397b6df1SKris Buschelman       }
486397b6df1SKris Buschelman     }
487397b6df1SKris 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);
488397b6df1SKris 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);
48994b7f48cSBarry 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);
490397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
491397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
492adc1d99fSHong 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);
493397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
494397b6df1SKris Buschelman 
495397b6df1SKris Buschelman     /*
496397b6df1SKris 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);
497397b6df1SKris Buschelman     if (flg){
498397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
499397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
500397b6df1SKris Buschelman       } else {
501397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
502397b6df1SKris Buschelman       }
503397b6df1SKris Buschelman     }
504397b6df1SKris Buschelman     */
505397b6df1SKris Buschelman 
506397b6df1SKris 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);
507397b6df1SKris 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);
508397b6df1SKris 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);
509397b6df1SKris Buschelman     PetscOptionsEnd();
510397b6df1SKris Buschelman   }
511397b6df1SKris Buschelman 
512397b6df1SKris Buschelman   /* define matrix A */
513397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
514397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
515397b6df1SKris Buschelman     if (!lu->myid) {
516c36ead0aSKris Buschelman       if (lua->isAIJ){
517397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
518397b6df1SKris Buschelman         nz               = aa->nz;
519397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
520397b6df1SKris Buschelman       } else {
521397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5226c6c5352SBarry Smith         nz                  =  aa->nz;
523397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
524397b6df1SKris Buschelman       }
525397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
526397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
527397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
528397b6df1SKris Buschelman         nz = 0;
529397b6df1SKris Buschelman         for (i=0; i<M; i++){
530397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
531397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
532397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
533397b6df1SKris Buschelman           }
534397b6df1SKris Buschelman         }
535397b6df1SKris Buschelman       }
536397b6df1SKris Buschelman     }
537397b6df1SKris Buschelman     break;
538397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
539397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
540397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
541397b6df1SKris Buschelman     } else {
542397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
543397b6df1SKris Buschelman     }
544397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
545397b6df1SKris Buschelman     break;
546397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
547397b6df1SKris Buschelman   }
548397b6df1SKris Buschelman 
549397b6df1SKris Buschelman   /* analysis phase */
550397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
551397b6df1SKris Buschelman      lu->id.n = M;
552397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
553397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
554397b6df1SKris Buschelman       if (!lu->myid) {
555397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
556397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
557397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
558397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
559397b6df1SKris Buschelman #else
560397b6df1SKris Buschelman           lu->id.a = lu->val;
561397b6df1SKris Buschelman #endif
562397b6df1SKris Buschelman         }
563397b6df1SKris Buschelman       }
564397b6df1SKris Buschelman       break;
565397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
566397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
567397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
568397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
569397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
570397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
571397b6df1SKris Buschelman #else
572397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
573397b6df1SKris Buschelman #endif
574397b6df1SKris Buschelman       }
575397b6df1SKris Buschelman       break;
576397b6df1SKris Buschelman     }
577397b6df1SKris Buschelman     lu->id.job=1;
578397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
579397b6df1SKris Buschelman   zmumps_c(&lu->id);
580397b6df1SKris Buschelman #else
581397b6df1SKris Buschelman   dmumps_c(&lu->id);
582397b6df1SKris Buschelman #endif
583397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
584397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
585397b6df1SKris Buschelman     }
586397b6df1SKris Buschelman   }
587397b6df1SKris Buschelman 
588397b6df1SKris Buschelman   /* numerical factorization phase */
589397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
590a7aca84bSHong Zhang     if (!lu->myid) {
591397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
592397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
593397b6df1SKris Buschelman #else
594397b6df1SKris Buschelman       lu->id.a = lu->val;
595397b6df1SKris Buschelman #endif
596397b6df1SKris Buschelman     }
597397b6df1SKris Buschelman   } else {
598397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
599397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
600397b6df1SKris Buschelman #else
601397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
602397b6df1SKris Buschelman #endif
603397b6df1SKris Buschelman   }
604397b6df1SKris Buschelman   lu->id.job=2;
605397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
606397b6df1SKris Buschelman   zmumps_c(&lu->id);
607397b6df1SKris Buschelman #else
608397b6df1SKris Buschelman   dmumps_c(&lu->id);
609397b6df1SKris Buschelman #endif
610397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
611a7aca84bSHong Zhang     SETERRQ2(1,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
612397b6df1SKris Buschelman   }
613397b6df1SKris Buschelman 
614397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
615397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
616397b6df1SKris Buschelman   }
617397b6df1SKris Buschelman 
618397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
619397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
620ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
621397b6df1SKris Buschelman   PetscFunctionReturn(0);
622397b6df1SKris Buschelman }
623397b6df1SKris Buschelman 
624397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
625397b6df1SKris Buschelman #undef __FUNCT__
626f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
627f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
628397b6df1SKris Buschelman   Mat       B;
629f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
630397b6df1SKris Buschelman   int       ierr;
631397b6df1SKris Buschelman 
632397b6df1SKris Buschelman   PetscFunctionBegin;
633397b6df1SKris Buschelman 
634397b6df1SKris Buschelman   /* Create the factorization matrix */
635397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
636397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
637397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
638397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
639397b6df1SKris Buschelman 
640f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
641397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
642f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
643397b6df1SKris Buschelman   lu->sym                 = 0;
644397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
645397b6df1SKris Buschelman 
646397b6df1SKris Buschelman   *F = B;
647397b6df1SKris Buschelman   PetscFunctionReturn(0);
648397b6df1SKris Buschelman }
649397b6df1SKris Buschelman 
650397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
651397b6df1SKris Buschelman #undef __FUNCT__
652f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
653f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
654397b6df1SKris Buschelman   Mat       B;
655f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
656397b6df1SKris Buschelman   int       ierr;
657397b6df1SKris Buschelman 
658397b6df1SKris Buschelman   PetscFunctionBegin;
659397b6df1SKris Buschelman 
660397b6df1SKris Buschelman   /* Create the factorization matrix */
661397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
662397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
663397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
664397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
665397b6df1SKris Buschelman 
666f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
667a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
668397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
669f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
670397b6df1SKris Buschelman   lu->sym                       = 2;
671397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
672397b6df1SKris Buschelman 
673397b6df1SKris Buschelman   *F = B;
674397b6df1SKris Buschelman   PetscFunctionReturn(0);
675397b6df1SKris Buschelman }
676397b6df1SKris Buschelman 
677397b6df1SKris Buschelman #undef __FUNCT__
678f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
679f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
680c338a77dSKris Buschelman   int       ierr;
681f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
682c338a77dSKris Buschelman 
683397b6df1SKris Buschelman   PetscFunctionBegin;
684c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
685f0c56d0fSKris Buschelman 
686c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
687c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
688f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
689397b6df1SKris Buschelman   PetscFunctionReturn(0);
690397b6df1SKris Buschelman }
691397b6df1SKris Buschelman 
692c338a77dSKris Buschelman EXTERN_C_BEGIN
693c338a77dSKris Buschelman #undef __FUNCT__
694f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
6958e9aea5cSBarry Smith int MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
696c338a77dSKris Buschelman   int       ierr,size;
697c338a77dSKris Buschelman   MPI_Comm  comm;
698c338a77dSKris Buschelman   Mat       B=*newmat;
699f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
700397b6df1SKris Buschelman 
701397b6df1SKris Buschelman   PetscFunctionBegin;
702c338a77dSKris Buschelman   if (B != A) {
703c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
704397b6df1SKris Buschelman   }
705397b6df1SKris Buschelman 
706c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
707f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
708c338a77dSKris Buschelman 
709f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
710c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
711c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
712c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
713c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
714c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
715a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
716c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
717f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
718c338a77dSKris Buschelman 
7194b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
720f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_AIJMUMPS;
721f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
722f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
723f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
7243924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
725c338a77dSKris Buschelman 
726c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
727c338a77dSKris Buschelman   if (size == 1) {
728c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
729f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
730c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
731c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
732c338a77dSKris Buschelman   } else {
733c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
734f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
735c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
736c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
737c338a77dSKris Buschelman   }
738c338a77dSKris Buschelman 
739f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
740c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
741c338a77dSKris Buschelman   *newmat = B;
742397b6df1SKris Buschelman   PetscFunctionReturn(0);
743397b6df1SKris Buschelman }
744c338a77dSKris Buschelman EXTERN_C_END
745397b6df1SKris Buschelman 
746f0c56d0fSKris Buschelman #undef __FUNCT__
747f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS"
748f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
749f0c56d0fSKris Buschelman   int       ierr;
7508f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
7518f340917SKris Buschelman 
752f0c56d0fSKris Buschelman   PetscFunctionBegin;
7538f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
754f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr);
755a39386dcSKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
756f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
757f0c56d0fSKris Buschelman }
758f0c56d0fSKris Buschelman 
75924b6179bSKris Buschelman /*MC
760fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
76124b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
76224b6179bSKris Buschelman 
76324b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
76424b6179bSKris Buschelman   on how to declare the existence of external packages),
76524b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
76624b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
76724b6179bSKris Buschelman   This matrix type is only supported for double precision real.
76824b6179bSKris Buschelman 
76924b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
77024b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
77124b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
77224b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
77328b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
77428b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
77528b08bd3SKris Buschelman   without data copy.
77624b6179bSKris Buschelman 
77724b6179bSKris Buschelman   Options Database Keys:
7780bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
77924b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
78024b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
78124b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
78224b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
78324b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
78424b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
78594b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
78624b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
78724b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
78824b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
78924b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
79024b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
79124b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
79224b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
79324b6179bSKris Buschelman 
79424b6179bSKris Buschelman   Level: beginner
79524b6179bSKris Buschelman 
79624b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
79724b6179bSKris Buschelman M*/
79824b6179bSKris Buschelman 
799397b6df1SKris Buschelman EXTERN_C_BEGIN
800397b6df1SKris Buschelman #undef __FUNCT__
801f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
802f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
803397b6df1SKris Buschelman   int      ierr,size;
804*e2d9671bSKris Buschelman   Mat      A_diag;
805397b6df1SKris Buschelman   MPI_Comm comm;
806397b6df1SKris Buschelman 
807397b6df1SKris Buschelman   PetscFunctionBegin;
8085441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
8095441df8eSKris Buschelman   /*   and AIJMUMPS types */
8105441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
811397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
812397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
813397b6df1SKris Buschelman   if (size == 1) {
814397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
815397b6df1SKris Buschelman   } else {
816397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
817*e2d9671bSKris Buschelman     A_diag = ((Mat_MPIAIJ *)A->data)->A;
818*e2d9671bSKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr);
819397b6df1SKris Buschelman   }
820f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
821397b6df1SKris Buschelman   PetscFunctionReturn(0);
822397b6df1SKris Buschelman }
823397b6df1SKris Buschelman EXTERN_C_END
824397b6df1SKris Buschelman 
825f579278aSKris Buschelman #undef __FUNCT__
826f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
827f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
828f579278aSKris Buschelman   int       ierr;
829f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
830f579278aSKris Buschelman 
831f579278aSKris Buschelman   PetscFunctionBegin;
832f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
833f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
834f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
835f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
836f579278aSKris Buschelman   PetscFunctionReturn(0);
837f579278aSKris Buschelman }
838f579278aSKris Buschelman 
839f579278aSKris Buschelman EXTERN_C_BEGIN
840f579278aSKris Buschelman #undef __FUNCT__
841f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
8428e9aea5cSBarry Smith int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
843f579278aSKris Buschelman   int       ierr,size;
844f579278aSKris Buschelman   MPI_Comm  comm;
845f579278aSKris Buschelman   Mat       B=*newmat;
846f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
847f579278aSKris Buschelman 
848f579278aSKris Buschelman   PetscFunctionBegin;
849f579278aSKris Buschelman   if (B != A) {
850f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
851f579278aSKris Buschelman   }
852f579278aSKris Buschelman 
853f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
854f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
855f579278aSKris Buschelman 
856f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
857f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
858f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
859f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
860f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
861f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
862a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
863f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
864f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
865f579278aSKris Buschelman 
866f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
867f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_SBAIJMUMPS;
868f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
869f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
870f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
8713924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
872f579278aSKris Buschelman 
873f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
874f579278aSKris Buschelman   if (size == 1) {
875f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
876f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
877f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
878f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
879f579278aSKris Buschelman   } else {
880f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
881f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
882f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
883f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
884f579278aSKris Buschelman   }
885f579278aSKris Buschelman 
886f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
887f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
888f579278aSKris Buschelman   *newmat = B;
889f579278aSKris Buschelman   PetscFunctionReturn(0);
890f579278aSKris Buschelman }
891f579278aSKris Buschelman EXTERN_C_END
892f579278aSKris Buschelman 
893f0c56d0fSKris Buschelman #undef __FUNCT__
894f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS"
895f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
896f0c56d0fSKris Buschelman   int       ierr;
8978f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
8988f340917SKris Buschelman 
899f0c56d0fSKris Buschelman   PetscFunctionBegin;
9008f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
901f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr);
9023f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
903f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
904f0c56d0fSKris Buschelman }
905f0c56d0fSKris Buschelman 
90624b6179bSKris Buschelman /*MC
907fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
90824b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
90924b6179bSKris Buschelman 
91024b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
91124b6179bSKris Buschelman   on how to declare the existence of external packages),
91224b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
91324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
91424b6179bSKris Buschelman   This matrix type is only supported for double precision real.
91524b6179bSKris Buschelman 
91624b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
91724b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
91824b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
91924b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
92028b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
92128b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
92228b08bd3SKris Buschelman   without data copy.
92324b6179bSKris Buschelman 
92424b6179bSKris Buschelman   Options Database Keys:
9250bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
92624b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
92724b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
92824b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
92924b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
93024b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
93124b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
93294b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
93324b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
93424b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
93524b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
93624b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
93724b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
93824b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
93924b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
94024b6179bSKris Buschelman 
94124b6179bSKris Buschelman   Level: beginner
94224b6179bSKris Buschelman 
94324b6179bSKris Buschelman .seealso: MATAIJMUMPS
94424b6179bSKris Buschelman M*/
94524b6179bSKris Buschelman 
946397b6df1SKris Buschelman EXTERN_C_BEGIN
947397b6df1SKris Buschelman #undef __FUNCT__
948f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
949f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
950397b6df1SKris Buschelman   int ierr,size;
951*e2d9671bSKris Buschelman   Mat A_diag;
952397b6df1SKris Buschelman 
953397b6df1SKris Buschelman   PetscFunctionBegin;
9545441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
9555441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
9565441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
9575441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
958397b6df1SKris Buschelman   if (size == 1) {
959397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
960397b6df1SKris Buschelman   } else {
961397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
962*e2d9671bSKris Buschelman     A_diag = ((Mat_MPISBAIJ *)A->data)->A;
963*e2d9671bSKris Buschelman     ierr   = MatConvert_SBAIJ_SBAIJMUMPS(A_diag,MATSBAIJMUMPS,&A_diag);CHKERRQ(ierr);
964397b6df1SKris Buschelman   }
965f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
966397b6df1SKris Buschelman   PetscFunctionReturn(0);
967397b6df1SKris Buschelman }
968397b6df1SKris Buschelman EXTERN_C_END
969