xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision be5d1d56a128fdbca06f8d9818f1d611ccde2ba2)
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);
469c097c71SKris Buschelman   int (*MatPreallocate)(Mat,int,int,int*,int,int*);
47f0c56d0fSKris Buschelman } Mat_MUMPS;
48f0c56d0fSKris Buschelman 
49f0c56d0fSKris Buschelman EXTERN int MatDuplicate_AIJMUMPS(Mat,MatDuplicateOption,Mat*);
50f0c56d0fSKris Buschelman EXTERN int MatDuplicate_SBAIJMUMPS(Mat,MatDuplicateOption,Mat*);
51892f6c3fSKris Buschelman EXTERN_C_BEGIN
52892f6c3fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,Mat*);
53892f6c3fSKris Buschelman EXTERN_C_END
54397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
55397b6df1SKris Buschelman /*
56397b6df1SKris Buschelman   input:
5775747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
58397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
59397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
60397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
61397b6df1SKris Buschelman   output:
62397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
63397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
64397b6df1SKris Buschelman  */
65f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
66397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
67397b6df1SKris Buschelman   int         ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
68d54de34fSKris Buschelman   int         *row,*col;
69397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
70f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
71397b6df1SKris Buschelman 
72397b6df1SKris Buschelman   PetscFunctionBegin;
73397b6df1SKris Buschelman   if (mumps->isAIJ){
74397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
75397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
76397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
77397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
78397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
79397b6df1SKris Buschelman     garray = mat->garray;
80397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
81397b6df1SKris Buschelman 
82397b6df1SKris Buschelman   } else {
83397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
84397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
85397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
86847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs);
876c6c5352SBarry Smith     nz = aa->nz + bb->nz;
88397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
89397b6df1SKris Buschelman     garray = mat->garray;
90397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
91397b6df1SKris Buschelman   }
92397b6df1SKris Buschelman 
93397b6df1SKris Buschelman   if (!valOnly){
94397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
95397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
96397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
97397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
98397b6df1SKris Buschelman   } else {
99397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
100397b6df1SKris Buschelman   }
101397b6df1SKris Buschelman   *nnz = nz;
102397b6df1SKris Buschelman 
103028e57e8SHong Zhang   jj = 0; irow = rstart;
104397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
105397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
106397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
107397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
108397b6df1SKris Buschelman     bjj = bj + bi[i];
109397b6df1SKris Buschelman 
110397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
111397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11275747be1SHong Zhang     j=-1;
11375747be1SHong Zhang     do {
11475747be1SHong Zhang       j++;
11575747be1SHong Zhang       if (j == countB) break;
116397b6df1SKris Buschelman       jcol = garray[bjj[j]];
11775747be1SHong Zhang     } while (jcol < colA_start);
11875747be1SHong Zhang     jB = j;
119397b6df1SKris Buschelman 
120397b6df1SKris Buschelman     /* B-part, smaller col index */
121397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
122397b6df1SKris Buschelman     for (j=0; j<jB; j++){
123397b6df1SKris Buschelman       jcol = garray[bjj[j]];
124397b6df1SKris Buschelman       if (!valOnly){
125397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12675747be1SHong Zhang 
127397b6df1SKris Buschelman       }
128397b6df1SKris Buschelman       val[jj++] = *bv++;
129397b6df1SKris Buschelman     }
130397b6df1SKris Buschelman     /* A-part */
131397b6df1SKris Buschelman     for (j=0; j<countA; j++){
132397b6df1SKris Buschelman       if (!valOnly){
133397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
134397b6df1SKris Buschelman       }
135397b6df1SKris Buschelman       val[jj++] = *av++;
136397b6df1SKris Buschelman     }
137397b6df1SKris Buschelman     /* B-part, larger col index */
138397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
139397b6df1SKris Buschelman       if (!valOnly){
140397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
141397b6df1SKris Buschelman       }
142397b6df1SKris Buschelman       val[jj++] = *bv++;
143397b6df1SKris Buschelman     }
144397b6df1SKris Buschelman     irow++;
145397b6df1SKris Buschelman   }
146397b6df1SKris Buschelman 
147397b6df1SKris Buschelman   PetscFunctionReturn(0);
148397b6df1SKris Buschelman }
149397b6df1SKris Buschelman 
150c338a77dSKris Buschelman EXTERN_C_BEGIN
151c338a77dSKris Buschelman #undef __FUNCT__
152c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
1538e9aea5cSBarry Smith int MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) {
154c338a77dSKris Buschelman   int       ierr;
155c338a77dSKris Buschelman   Mat       B=*newmat;
156f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
157c338a77dSKris Buschelman 
158c338a77dSKris Buschelman   PetscFunctionBegin;
159c338a77dSKris Buschelman   if (B != A) {
160c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
161c338a77dSKris Buschelman   }
162f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
163f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
164f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
165f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
166f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
167f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
1683924e44cSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
169f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
170c338a77dSKris Buschelman   *newmat = B;
171c338a77dSKris Buschelman   PetscFunctionReturn(0);
172c338a77dSKris Buschelman }
173c338a77dSKris Buschelman EXTERN_C_END
174c338a77dSKris Buschelman 
175397b6df1SKris Buschelman #undef __FUNCT__
1763924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
1773924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) {
178f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
179c338a77dSKris Buschelman   int       ierr,size=lu->size;
180a39386dcSKris Buschelman   int       (*specialdestroy)(Mat);
181397b6df1SKris Buschelman   PetscFunctionBegin;
182397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
183397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
184397b6df1SKris Buschelman     lu->id.job=JOB_END;
185397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
186397b6df1SKris Buschelman     zmumps_c(&lu->id);
187397b6df1SKris Buschelman #else
188397b6df1SKris Buschelman     dmumps_c(&lu->id);
189397b6df1SKris Buschelman #endif
190c338a77dSKris Buschelman     if (lu->irn) {
191c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
192c338a77dSKris Buschelman     }
193c338a77dSKris Buschelman     if (lu->jcn) {
194c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
195c338a77dSKris Buschelman     }
196c338a77dSKris Buschelman     if (size>1 && lu->val) {
197c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
198c338a77dSKris Buschelman     }
199397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
200397b6df1SKris Buschelman   }
201a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
202a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
203c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
204397b6df1SKris Buschelman   PetscFunctionReturn(0);
205397b6df1SKris Buschelman }
206397b6df1SKris Buschelman 
207397b6df1SKris Buschelman #undef __FUNCT__
208a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
209a39386dcSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) {
210a39386dcSKris Buschelman   int ierr, size;
211a39386dcSKris Buschelman 
212a39386dcSKris Buschelman   PetscFunctionBegin;
213a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
214a39386dcSKris Buschelman   if (size==1) {
215a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr);
216a39386dcSKris Buschelman   } else {
217a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr);
218a39386dcSKris Buschelman   }
219a39386dcSKris Buschelman   PetscFunctionReturn(0);
220a39386dcSKris Buschelman }
221a39386dcSKris Buschelman 
222a39386dcSKris Buschelman #undef __FUNCT__
223a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
224a39386dcSKris Buschelman int MatDestroy_SBAIJMUMPS(Mat A) {
225a39386dcSKris Buschelman   int ierr, size;
226a39386dcSKris Buschelman 
227a39386dcSKris Buschelman   PetscFunctionBegin;
228a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
229a39386dcSKris Buschelman   if (size==1) {
230a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr);
231a39386dcSKris Buschelman   } else {
232a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr);
233a39386dcSKris Buschelman   }
234a39386dcSKris Buschelman   PetscFunctionReturn(0);
235a39386dcSKris Buschelman }
236a39386dcSKris Buschelman 
237a39386dcSKris Buschelman #undef __FUNCT__
238c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
239f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
240f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
241397b6df1SKris Buschelman   int       ierr;
242397b6df1SKris Buschelman 
243397b6df1SKris Buschelman   PetscFunctionBegin;
244c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
245c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
246c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
247c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
248c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
249c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
250c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
251c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
252c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
253c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
254c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
255c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
256c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
257c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
258c338a77dSKris 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);
259c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
260c338a77dSKris 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);
261c338a77dSKris Buschelman 
262c338a77dSKris Buschelman   }
263c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
264c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
265adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
266c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
267c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
268c338a77dSKris Buschelman 
269c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
270c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
271c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
27257f0c58bSHong Zhang 
27357f0c58bSHong Zhang   /* infomation local to each processor */
27457f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
27557f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
27657f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27757f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
27857f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
27957f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
28057f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
28157f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
28257f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
283adc1d99fSHong Zhang 
284adc1d99fSHong Zhang   if (lu->myid == 0){ /* information from the host */
285adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
286adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
287adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
288adc1d99fSHong Zhang 
289adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
290adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
291adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
292adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
293adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
294adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
295adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
296adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
297adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
298adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
299adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
300adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
301adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
302adc1d99fSHong 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);
303adc1d99fSHong 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);
304adc1d99fSHong 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);
305adc1d99fSHong 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);
306adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
307adc1d99fSHong Zhang   }
308adc1d99fSHong Zhang 
309397b6df1SKris Buschelman   PetscFunctionReturn(0);
310397b6df1SKris Buschelman }
311397b6df1SKris Buschelman 
312397b6df1SKris Buschelman #undef __FUNCT__
313f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
314f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
315397b6df1SKris Buschelman   int               ierr;
316397b6df1SKris Buschelman   PetscTruth        isascii;
317397b6df1SKris Buschelman   PetscViewerFormat format;
318f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
319397b6df1SKris Buschelman 
320397b6df1SKris Buschelman   PetscFunctionBegin;
321397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
322397b6df1SKris Buschelman 
323397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
324397b6df1SKris Buschelman   if (isascii) {
325397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
326397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
327397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
328397b6df1SKris Buschelman     }
329397b6df1SKris Buschelman   }
330397b6df1SKris Buschelman   PetscFunctionReturn(0);
331397b6df1SKris Buschelman }
332397b6df1SKris Buschelman 
333397b6df1SKris Buschelman #undef __FUNCT__
334f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
335f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
336f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
337d54de34fSKris Buschelman   PetscScalar *array;
338397b6df1SKris Buschelman   Vec         x_seq;
339397b6df1SKris Buschelman   IS          iden;
340397b6df1SKris Buschelman   VecScatter  scat;
341397b6df1SKris Buschelman   int         ierr;
342397b6df1SKris Buschelman 
343397b6df1SKris Buschelman   PetscFunctionBegin;
344397b6df1SKris Buschelman   if (lu->size > 1){
345397b6df1SKris Buschelman     if (!lu->myid){
346397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
347397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
348397b6df1SKris Buschelman     } else {
349397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
350397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
351397b6df1SKris Buschelman     }
352397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
353397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
354397b6df1SKris Buschelman 
355397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
356397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
357397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
358397b6df1SKris Buschelman   } else {  /* size == 1 */
359397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
360397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
361397b6df1SKris Buschelman   }
362397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
363397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
364397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
365397b6df1SKris Buschelman #else
366397b6df1SKris Buschelman     lu->id.rhs = array;
367397b6df1SKris Buschelman #endif
368397b6df1SKris Buschelman   }
369397b6df1SKris Buschelman 
370397b6df1SKris Buschelman   /* solve phase */
371397b6df1SKris Buschelman   lu->id.job=3;
372397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
373397b6df1SKris Buschelman   zmumps_c(&lu->id);
374397b6df1SKris Buschelman #else
375397b6df1SKris Buschelman   dmumps_c(&lu->id);
376397b6df1SKris Buschelman #endif
377397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
378397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
379397b6df1SKris Buschelman   }
380397b6df1SKris Buschelman 
381397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
382397b6df1SKris Buschelman   if (lu->size > 1) {
383397b6df1SKris Buschelman     if (!lu->myid){
384397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
385397b6df1SKris Buschelman     }
386397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
387397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
388397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
389397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
390397b6df1SKris Buschelman   } else {
391397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
392397b6df1SKris Buschelman   }
393397b6df1SKris Buschelman 
394397b6df1SKris Buschelman   PetscFunctionReturn(0);
395397b6df1SKris Buschelman }
396397b6df1SKris Buschelman 
397a58c3f20SHong Zhang /*
398a58c3f20SHong Zhang   input:
399a58c3f20SHong Zhang    F:        numeric factor
400a58c3f20SHong Zhang   output:
401a58c3f20SHong Zhang    nneg:     total number of negative pivots
402a58c3f20SHong Zhang    nzero:    0
403a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
404a58c3f20SHong Zhang */
405a58c3f20SHong Zhang 
406a58c3f20SHong Zhang #undef __FUNCT__
407a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
408a58c3f20SHong Zhang int MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
409a58c3f20SHong Zhang {
410a58c3f20SHong Zhang   Mat_MUMPS  *lu =(Mat_MUMPS*)F->spptr;
411bcb30aebSHong Zhang   int        ierr,neg,zero,pos,size;
412a58c3f20SHong Zhang 
413a58c3f20SHong Zhang   PetscFunctionBegin;
414bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
415bcb30aebSHong 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 */
416bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
417bcb30aebSHong 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));
418bcb30aebSHong Zhang   }
419a58c3f20SHong Zhang   if (nneg){
420a58c3f20SHong Zhang     if (!lu->myid){
421a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
422a58c3f20SHong Zhang     }
423bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
424a58c3f20SHong Zhang   }
425a58c3f20SHong Zhang   if (nzero) *nzero = 0;
426a58c3f20SHong Zhang   if (npos)  *npos  = F->M - (*nneg);
427a58c3f20SHong Zhang   PetscFunctionReturn(0);
428a58c3f20SHong Zhang }
429a58c3f20SHong Zhang 
430397b6df1SKris Buschelman #undef __FUNCT__
431f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
432f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
433f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
434f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
435397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
436397b6df1SKris Buschelman   PetscTruth valOnly,flg;
437397b6df1SKris Buschelman 
438397b6df1SKris Buschelman   PetscFunctionBegin;
439397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
440f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
441397b6df1SKris Buschelman 
442397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
443397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
444397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
44575747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
446397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
447397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
448397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
449397b6df1SKris Buschelman 
450397b6df1SKris Buschelman     /* Set mumps options */
451397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
452397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
453397b6df1SKris Buschelman     lu->id.sym=lu->sym;
454397b6df1SKris Buschelman     if (lu->sym == 2){
455397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
456397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
457397b6df1SKris Buschelman     }
458397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
459397b6df1SKris Buschelman   zmumps_c(&lu->id);
460397b6df1SKris Buschelman #else
461397b6df1SKris Buschelman   dmumps_c(&lu->id);
462397b6df1SKris Buschelman #endif
463397b6df1SKris Buschelman 
464397b6df1SKris Buschelman     if (lu->size == 1){
465397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
466397b6df1SKris Buschelman     } else {
467397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
468397b6df1SKris Buschelman     }
469397b6df1SKris Buschelman 
470397b6df1SKris Buschelman     icntl=-1;
471397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
472397b6df1SKris Buschelman     if (flg && icntl > 0) {
473397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
474397b6df1SKris Buschelman     } else { /* no output */
475397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
476397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
477397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
478397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
479397b6df1SKris Buschelman     }
480397b6df1SKris 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);
481397b6df1SKris Buschelman     icntl=-1;
482397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
483397b6df1SKris Buschelman     if (flg) {
484397b6df1SKris Buschelman       if (icntl== 1){
485397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
486397b6df1SKris Buschelman       } else {
487397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
488397b6df1SKris Buschelman       }
489397b6df1SKris Buschelman     }
490397b6df1SKris 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);
491397b6df1SKris 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);
49294b7f48cSBarry 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);
493397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
494397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
495adc1d99fSHong 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);
496397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
497397b6df1SKris Buschelman 
498397b6df1SKris Buschelman     /*
499397b6df1SKris 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);
500397b6df1SKris Buschelman     if (flg){
501397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
502397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
503397b6df1SKris Buschelman       } else {
504397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
505397b6df1SKris Buschelman       }
506397b6df1SKris Buschelman     }
507397b6df1SKris Buschelman     */
508397b6df1SKris Buschelman 
509397b6df1SKris 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);
510397b6df1SKris 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);
511397b6df1SKris 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);
512397b6df1SKris Buschelman     PetscOptionsEnd();
513397b6df1SKris Buschelman   }
514397b6df1SKris Buschelman 
515397b6df1SKris Buschelman   /* define matrix A */
516397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
517397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
518397b6df1SKris Buschelman     if (!lu->myid) {
519c36ead0aSKris Buschelman       if (lua->isAIJ){
520397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
521397b6df1SKris Buschelman         nz               = aa->nz;
522397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
523397b6df1SKris Buschelman       } else {
524397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5256c6c5352SBarry Smith         nz                  =  aa->nz;
526397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
527397b6df1SKris Buschelman       }
528397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
529397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
530397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
531397b6df1SKris Buschelman         nz = 0;
532397b6df1SKris Buschelman         for (i=0; i<M; i++){
533397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
534397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
535397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
536397b6df1SKris Buschelman           }
537397b6df1SKris Buschelman         }
538397b6df1SKris Buschelman       }
539397b6df1SKris Buschelman     }
540397b6df1SKris Buschelman     break;
541397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
542397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
543397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
544397b6df1SKris Buschelman     } else {
545397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
546397b6df1SKris Buschelman     }
547397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
548397b6df1SKris Buschelman     break;
549397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
550397b6df1SKris Buschelman   }
551397b6df1SKris Buschelman 
552397b6df1SKris Buschelman   /* analysis phase */
553397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
554397b6df1SKris Buschelman      lu->id.n = M;
555397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
556397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
557397b6df1SKris Buschelman       if (!lu->myid) {
558397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
559397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
560397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
561397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
562397b6df1SKris Buschelman #else
563397b6df1SKris Buschelman           lu->id.a = lu->val;
564397b6df1SKris Buschelman #endif
565397b6df1SKris Buschelman         }
566397b6df1SKris Buschelman       }
567397b6df1SKris Buschelman       break;
568397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
569397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
570397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
571397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
572397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
573397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
574397b6df1SKris Buschelman #else
575397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
576397b6df1SKris Buschelman #endif
577397b6df1SKris Buschelman       }
578397b6df1SKris Buschelman       break;
579397b6df1SKris Buschelman     }
580397b6df1SKris Buschelman     lu->id.job=1;
581397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
582397b6df1SKris Buschelman   zmumps_c(&lu->id);
583397b6df1SKris Buschelman #else
584397b6df1SKris Buschelman   dmumps_c(&lu->id);
585397b6df1SKris Buschelman #endif
586397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
587397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
588397b6df1SKris Buschelman     }
589397b6df1SKris Buschelman   }
590397b6df1SKris Buschelman 
591397b6df1SKris Buschelman   /* numerical factorization phase */
592397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
593a7aca84bSHong Zhang     if (!lu->myid) {
594397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
595397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
596397b6df1SKris Buschelman #else
597397b6df1SKris Buschelman       lu->id.a = lu->val;
598397b6df1SKris Buschelman #endif
599397b6df1SKris Buschelman     }
600397b6df1SKris Buschelman   } else {
601397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
602397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
603397b6df1SKris Buschelman #else
604397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
605397b6df1SKris Buschelman #endif
606397b6df1SKris Buschelman   }
607397b6df1SKris Buschelman   lu->id.job=2;
608397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
609397b6df1SKris Buschelman   zmumps_c(&lu->id);
610397b6df1SKris Buschelman #else
611397b6df1SKris Buschelman   dmumps_c(&lu->id);
612397b6df1SKris Buschelman #endif
613397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
614a7aca84bSHong 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));
615397b6df1SKris Buschelman   }
616397b6df1SKris Buschelman 
617397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
618397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
619397b6df1SKris Buschelman   }
620397b6df1SKris Buschelman 
621397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
622397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
623ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
624397b6df1SKris Buschelman   PetscFunctionReturn(0);
625397b6df1SKris Buschelman }
626397b6df1SKris Buschelman 
627397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
628397b6df1SKris Buschelman #undef __FUNCT__
629f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
630f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
631397b6df1SKris Buschelman   Mat       B;
632f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
633397b6df1SKris Buschelman   int       ierr;
634397b6df1SKris Buschelman 
635397b6df1SKris Buschelman   PetscFunctionBegin;
636397b6df1SKris Buschelman 
637397b6df1SKris Buschelman   /* Create the factorization matrix */
638397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
639*be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
640397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
641397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
642397b6df1SKris Buschelman 
643f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
644397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
645f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
646397b6df1SKris Buschelman   lu->sym                 = 0;
647397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
648397b6df1SKris Buschelman 
649397b6df1SKris Buschelman   *F = B;
650397b6df1SKris Buschelman   PetscFunctionReturn(0);
651397b6df1SKris Buschelman }
652397b6df1SKris Buschelman 
653397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
654397b6df1SKris Buschelman #undef __FUNCT__
655f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
656f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
657397b6df1SKris Buschelman   Mat       B;
658f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
659397b6df1SKris Buschelman   int       ierr;
660397b6df1SKris Buschelman 
661397b6df1SKris Buschelman   PetscFunctionBegin;
662397b6df1SKris Buschelman 
663397b6df1SKris Buschelman   /* Create the factorization matrix */
664397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
665*be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
666397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
667397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
668397b6df1SKris Buschelman 
669f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
670a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
671397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
672f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
673397b6df1SKris Buschelman   lu->sym                       = 2;
674397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
675397b6df1SKris Buschelman 
676397b6df1SKris Buschelman   *F = B;
677397b6df1SKris Buschelman   PetscFunctionReturn(0);
678397b6df1SKris Buschelman }
679397b6df1SKris Buschelman 
680397b6df1SKris Buschelman #undef __FUNCT__
681f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
682f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
683c338a77dSKris Buschelman   int       ierr;
684f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
685c338a77dSKris Buschelman 
686397b6df1SKris Buschelman   PetscFunctionBegin;
687c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
688f0c56d0fSKris Buschelman 
689c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
690c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
691f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
692397b6df1SKris Buschelman   PetscFunctionReturn(0);
693397b6df1SKris Buschelman }
694397b6df1SKris Buschelman 
695c338a77dSKris Buschelman EXTERN_C_BEGIN
696c338a77dSKris Buschelman #undef __FUNCT__
697f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
6988e9aea5cSBarry Smith int MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
699c338a77dSKris Buschelman   int       ierr,size;
700c338a77dSKris Buschelman   MPI_Comm  comm;
701c338a77dSKris Buschelman   Mat       B=*newmat;
702f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
703397b6df1SKris Buschelman 
704397b6df1SKris Buschelman   PetscFunctionBegin;
705c338a77dSKris Buschelman   if (B != A) {
706c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
707397b6df1SKris Buschelman   }
708397b6df1SKris Buschelman 
709c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
710f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
711c338a77dSKris Buschelman 
712f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
713c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
714c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
715c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
716c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
717c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
718a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
719c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
720f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
721c338a77dSKris Buschelman 
7224b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
723f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_AIJMUMPS;
724f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
725f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
726f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
7273924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
728c338a77dSKris Buschelman 
729c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
730c338a77dSKris Buschelman   if (size == 1) {
731c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
732f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
733c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
734c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
735c338a77dSKris Buschelman   } else {
736c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
737f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
738c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
739c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
740c338a77dSKris Buschelman   }
741c338a77dSKris Buschelman 
742f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
743c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
744c338a77dSKris Buschelman   *newmat = B;
745397b6df1SKris Buschelman   PetscFunctionReturn(0);
746397b6df1SKris Buschelman }
747c338a77dSKris Buschelman EXTERN_C_END
748397b6df1SKris Buschelman 
749f0c56d0fSKris Buschelman #undef __FUNCT__
750f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS"
751f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
752f0c56d0fSKris Buschelman   int       ierr;
7538f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
7548f340917SKris Buschelman 
755f0c56d0fSKris Buschelman   PetscFunctionBegin;
7568f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
757f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr);
758a39386dcSKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
759f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
760f0c56d0fSKris Buschelman }
761f0c56d0fSKris Buschelman 
76224b6179bSKris Buschelman /*MC
763fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
76424b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
76524b6179bSKris Buschelman 
76624b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
76724b6179bSKris Buschelman   on how to declare the existence of external packages),
76824b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
76924b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
77024b6179bSKris Buschelman   This matrix type is only supported for double precision real.
77124b6179bSKris Buschelman 
77224b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
77324b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
77424b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
77524b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
77628b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
77728b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
77828b08bd3SKris Buschelman   without data copy.
77924b6179bSKris Buschelman 
78024b6179bSKris Buschelman   Options Database Keys:
7810bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
78224b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
78324b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
78424b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
78524b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
78624b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
78724b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
78894b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
78924b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
79024b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
79124b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
79224b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
79324b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
79424b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
79524b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
79624b6179bSKris Buschelman 
79724b6179bSKris Buschelman   Level: beginner
79824b6179bSKris Buschelman 
79924b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
80024b6179bSKris Buschelman M*/
80124b6179bSKris Buschelman 
802397b6df1SKris Buschelman EXTERN_C_BEGIN
803397b6df1SKris Buschelman #undef __FUNCT__
804f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
805f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
806397b6df1SKris Buschelman   int      ierr,size;
807e2d9671bSKris Buschelman   Mat      A_diag;
808397b6df1SKris Buschelman   MPI_Comm comm;
809397b6df1SKris Buschelman 
810397b6df1SKris Buschelman   PetscFunctionBegin;
8115441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
8125441df8eSKris Buschelman   /*   and AIJMUMPS types */
8135441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
814397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
815397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
816397b6df1SKris Buschelman   if (size == 1) {
817397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
818397b6df1SKris Buschelman   } else {
819397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
820e2d9671bSKris Buschelman     A_diag = ((Mat_MPIAIJ *)A->data)->A;
821e2d9671bSKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr);
822397b6df1SKris Buschelman   }
823f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
824397b6df1SKris Buschelman   PetscFunctionReturn(0);
825397b6df1SKris Buschelman }
826397b6df1SKris Buschelman EXTERN_C_END
827397b6df1SKris Buschelman 
828f579278aSKris Buschelman #undef __FUNCT__
829f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
830f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
831f579278aSKris Buschelman   int       ierr;
832f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
833f579278aSKris Buschelman 
834f579278aSKris Buschelman   PetscFunctionBegin;
835f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
836f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
837f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
838f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
839f579278aSKris Buschelman   PetscFunctionReturn(0);
840f579278aSKris Buschelman }
841f579278aSKris Buschelman 
842f579278aSKris Buschelman EXTERN_C_BEGIN
843f579278aSKris Buschelman #undef __FUNCT__
8449c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
8459c097c71SKris Buschelman int MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat  B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz)
8469c097c71SKris Buschelman {
8479c097c71SKris Buschelman   Mat       A;
84802217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
8499c097c71SKris Buschelman   int       ierr;
8509c097c71SKris Buschelman 
8519c097c71SKris Buschelman   PetscFunctionBegin;
8529c097c71SKris Buschelman   /*
8539c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
8549c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
8559c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
8569c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
8579c097c71SKris Buschelman     in the preallocation routine.
8589c097c71SKris Buschelman   */
8599c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
8609c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
8619c097c71SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
8629c097c71SKris Buschelman   PetscFunctionReturn(0);
8639c097c71SKris Buschelman }
8649c097c71SKris Buschelman EXTERN_C_END
8659c097c71SKris Buschelman 
8669c097c71SKris Buschelman EXTERN_C_BEGIN
8679c097c71SKris Buschelman #undef __FUNCT__
868f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
8698e9aea5cSBarry Smith int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
870f579278aSKris Buschelman   int       ierr,size;
871f579278aSKris Buschelman   MPI_Comm  comm;
872f579278aSKris Buschelman   Mat       B=*newmat;
87302217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
8749c097c71SKris Buschelman   void      (*f)(void);
875f579278aSKris Buschelman 
876f579278aSKris Buschelman   PetscFunctionBegin;
877f579278aSKris Buschelman   if (B != A) {
878f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
879f579278aSKris Buschelman   }
880f579278aSKris Buschelman 
881f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
882f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
883f579278aSKris Buschelman 
884f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
885f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
886f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
887f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
888f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
889f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
890a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
891f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
892f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
893f579278aSKris Buschelman 
894f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
895f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_SBAIJMUMPS;
896f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
897f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
898f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
8993924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
900f579278aSKris Buschelman 
901f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
902f579278aSKris Buschelman   if (size == 1) {
903f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
904f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
905f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
906f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
907f579278aSKris Buschelman   } else {
9089c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
9099c097c71SKris Buschelman     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
9109c097c71SKris Buschelman     if (f) {
9119c097c71SKris Buschelman       mumps->MatPreallocate = (int (*)(Mat,int,int,int*,int,int*))f;
9129c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
9139c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
9149c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
9159c097c71SKris Buschelman     }
916f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
917f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
918f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
919f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
920f579278aSKris Buschelman   }
921f579278aSKris Buschelman 
922f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
923f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
924f579278aSKris Buschelman   *newmat = B;
925f579278aSKris Buschelman   PetscFunctionReturn(0);
926f579278aSKris Buschelman }
927f579278aSKris Buschelman EXTERN_C_END
928f579278aSKris Buschelman 
929f0c56d0fSKris Buschelman #undef __FUNCT__
930f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS"
931f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
932f0c56d0fSKris Buschelman   int       ierr;
9338f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
9348f340917SKris Buschelman 
935f0c56d0fSKris Buschelman   PetscFunctionBegin;
9368f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
937f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr);
9383f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
939f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
940f0c56d0fSKris Buschelman }
941f0c56d0fSKris Buschelman 
94224b6179bSKris Buschelman /*MC
943fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
94424b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
94524b6179bSKris Buschelman 
94624b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
94724b6179bSKris Buschelman   on how to declare the existence of external packages),
94824b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
94924b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
95024b6179bSKris Buschelman   This matrix type is only supported for double precision real.
95124b6179bSKris Buschelman 
95224b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
95324b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
95424b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
95524b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
95628b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
95728b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
95828b08bd3SKris Buschelman   without data copy.
95924b6179bSKris Buschelman 
96024b6179bSKris Buschelman   Options Database Keys:
9610bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
96224b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
96324b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
96424b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
96524b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
96624b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
96724b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
96894b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
96924b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
97024b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
97124b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
97224b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
97324b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
97424b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
97524b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
97624b6179bSKris Buschelman 
97724b6179bSKris Buschelman   Level: beginner
97824b6179bSKris Buschelman 
97924b6179bSKris Buschelman .seealso: MATAIJMUMPS
98024b6179bSKris Buschelman M*/
98124b6179bSKris Buschelman 
982397b6df1SKris Buschelman EXTERN_C_BEGIN
983397b6df1SKris Buschelman #undef __FUNCT__
984f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
985f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
986397b6df1SKris Buschelman   int ierr,size;
987397b6df1SKris Buschelman 
988397b6df1SKris Buschelman   PetscFunctionBegin;
9895441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
9905441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
9915441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
9925441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
993397b6df1SKris Buschelman   if (size == 1) {
994397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
995397b6df1SKris Buschelman   } else {
996397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
997397b6df1SKris Buschelman   }
998f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
999397b6df1SKris Buschelman   PetscFunctionReturn(0);
1000397b6df1SKris Buschelman }
1001397b6df1SKris Buschelman EXTERN_C_END
1002