xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 9c097c7127bb68c85bf5f7a436f8a4b636fd1b64)
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);
46*9c097c71SKris 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*);
51*9c097c71SKris Buschelman EXTERN int MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,Mat*);
52397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
53397b6df1SKris Buschelman /*
54397b6df1SKris Buschelman   input:
5575747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
56397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
57397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
58397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
59397b6df1SKris Buschelman   output:
60397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
61397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
62397b6df1SKris Buschelman  */
63f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
64397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
65397b6df1SKris Buschelman   int         ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
66d54de34fSKris Buschelman   int         *row,*col;
67397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
68f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
69397b6df1SKris Buschelman 
70397b6df1SKris Buschelman   PetscFunctionBegin;
71397b6df1SKris Buschelman   if (mumps->isAIJ){
72397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
73397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
74397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
75397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
76397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
77397b6df1SKris Buschelman     garray = mat->garray;
78397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
79397b6df1SKris Buschelman 
80397b6df1SKris Buschelman   } else {
81397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
82397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
83397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
84847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs);
856c6c5352SBarry Smith     nz = aa->nz + bb->nz;
86397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
87397b6df1SKris Buschelman     garray = mat->garray;
88397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
89397b6df1SKris Buschelman   }
90397b6df1SKris Buschelman 
91397b6df1SKris Buschelman   if (!valOnly){
92397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
93397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
94397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
95397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
96397b6df1SKris Buschelman   } else {
97397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
98397b6df1SKris Buschelman   }
99397b6df1SKris Buschelman   *nnz = nz;
100397b6df1SKris Buschelman 
101028e57e8SHong Zhang   jj = 0; irow = rstart;
102397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
103397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
104397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
105397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
106397b6df1SKris Buschelman     bjj = bj + bi[i];
107397b6df1SKris Buschelman 
108397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
109397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11075747be1SHong Zhang     j=-1;
11175747be1SHong Zhang     do {
11275747be1SHong Zhang       j++;
11375747be1SHong Zhang       if (j == countB) break;
114397b6df1SKris Buschelman       jcol = garray[bjj[j]];
11575747be1SHong Zhang     } while (jcol < colA_start);
11675747be1SHong Zhang     jB = j;
117397b6df1SKris Buschelman 
118397b6df1SKris Buschelman     /* B-part, smaller col index */
119397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
120397b6df1SKris Buschelman     for (j=0; j<jB; j++){
121397b6df1SKris Buschelman       jcol = garray[bjj[j]];
122397b6df1SKris Buschelman       if (!valOnly){
123397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12475747be1SHong Zhang 
125397b6df1SKris Buschelman       }
126397b6df1SKris Buschelman       val[jj++] = *bv++;
127397b6df1SKris Buschelman     }
128397b6df1SKris Buschelman     /* A-part */
129397b6df1SKris Buschelman     for (j=0; j<countA; j++){
130397b6df1SKris Buschelman       if (!valOnly){
131397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
132397b6df1SKris Buschelman       }
133397b6df1SKris Buschelman       val[jj++] = *av++;
134397b6df1SKris Buschelman     }
135397b6df1SKris Buschelman     /* B-part, larger col index */
136397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
137397b6df1SKris Buschelman       if (!valOnly){
138397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
139397b6df1SKris Buschelman       }
140397b6df1SKris Buschelman       val[jj++] = *bv++;
141397b6df1SKris Buschelman     }
142397b6df1SKris Buschelman     irow++;
143397b6df1SKris Buschelman   }
144397b6df1SKris Buschelman 
145397b6df1SKris Buschelman   PetscFunctionReturn(0);
146397b6df1SKris Buschelman }
147397b6df1SKris Buschelman 
148c338a77dSKris Buschelman EXTERN_C_BEGIN
149c338a77dSKris Buschelman #undef __FUNCT__
150c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
1518e9aea5cSBarry Smith int MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) {
152c338a77dSKris Buschelman   int       ierr;
153c338a77dSKris Buschelman   Mat       B=*newmat;
154f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
155c338a77dSKris Buschelman 
156c338a77dSKris Buschelman   PetscFunctionBegin;
157c338a77dSKris Buschelman   if (B != A) {
158c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
159c338a77dSKris Buschelman   }
160f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
161f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
162f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
163f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
164f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
165f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
1663924e44cSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
167f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
168c338a77dSKris Buschelman   *newmat = B;
169c338a77dSKris Buschelman   PetscFunctionReturn(0);
170c338a77dSKris Buschelman }
171c338a77dSKris Buschelman EXTERN_C_END
172c338a77dSKris Buschelman 
173397b6df1SKris Buschelman #undef __FUNCT__
1743924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
1753924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) {
176f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
177c338a77dSKris Buschelman   int       ierr,size=lu->size;
178a39386dcSKris Buschelman   int       (*specialdestroy)(Mat);
179397b6df1SKris Buschelman   PetscFunctionBegin;
180397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
181397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
182397b6df1SKris Buschelman     lu->id.job=JOB_END;
183397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
184397b6df1SKris Buschelman     zmumps_c(&lu->id);
185397b6df1SKris Buschelman #else
186397b6df1SKris Buschelman     dmumps_c(&lu->id);
187397b6df1SKris Buschelman #endif
188c338a77dSKris Buschelman     if (lu->irn) {
189c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
190c338a77dSKris Buschelman     }
191c338a77dSKris Buschelman     if (lu->jcn) {
192c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
193c338a77dSKris Buschelman     }
194c338a77dSKris Buschelman     if (size>1 && lu->val) {
195c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
196c338a77dSKris Buschelman     }
197397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
198397b6df1SKris Buschelman   }
199a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
200a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
201c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
202397b6df1SKris Buschelman   PetscFunctionReturn(0);
203397b6df1SKris Buschelman }
204397b6df1SKris Buschelman 
205397b6df1SKris Buschelman #undef __FUNCT__
206a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
207a39386dcSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) {
208a39386dcSKris Buschelman   int ierr, size;
209a39386dcSKris Buschelman 
210a39386dcSKris Buschelman   PetscFunctionBegin;
211a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
212a39386dcSKris Buschelman   if (size==1) {
213a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr);
214a39386dcSKris Buschelman   } else {
215a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr);
216a39386dcSKris Buschelman   }
217a39386dcSKris Buschelman   PetscFunctionReturn(0);
218a39386dcSKris Buschelman }
219a39386dcSKris Buschelman 
220a39386dcSKris Buschelman #undef __FUNCT__
221a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
222a39386dcSKris Buschelman int MatDestroy_SBAIJMUMPS(Mat A) {
223a39386dcSKris Buschelman   int ierr, size;
224a39386dcSKris Buschelman 
225a39386dcSKris Buschelman   PetscFunctionBegin;
226a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
227a39386dcSKris Buschelman   if (size==1) {
228a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr);
229a39386dcSKris Buschelman   } else {
230a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr);
231a39386dcSKris Buschelman   }
232a39386dcSKris Buschelman   PetscFunctionReturn(0);
233a39386dcSKris Buschelman }
234a39386dcSKris Buschelman 
235a39386dcSKris Buschelman #undef __FUNCT__
236c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
237f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
238f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
239397b6df1SKris Buschelman   int       ierr;
240397b6df1SKris Buschelman 
241397b6df1SKris Buschelman   PetscFunctionBegin;
242c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
243c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
244c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
245c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
246c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
247c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
248c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
249c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
250c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
251c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
252c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
253c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
254c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
255c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
256c338a77dSKris 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);
257c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
258c338a77dSKris 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);
259c338a77dSKris Buschelman 
260c338a77dSKris Buschelman   }
261c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
262c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
263adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
264c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
265c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
266c338a77dSKris Buschelman 
267c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
268c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
269c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
27057f0c58bSHong Zhang 
27157f0c58bSHong Zhang   /* infomation local to each processor */
27257f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
27357f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
27457f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27557f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
27657f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
27757f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27857f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
27957f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
28057f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
281adc1d99fSHong Zhang 
282adc1d99fSHong Zhang   if (lu->myid == 0){ /* information from the host */
283adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
284adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
285adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
286adc1d99fSHong Zhang 
287adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
288adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
289adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
290adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
291adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
292adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
293adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
294adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
295adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
296adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
297adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
298adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
299adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
300adc1d99fSHong 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);
301adc1d99fSHong 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);
302adc1d99fSHong 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);
303adc1d99fSHong 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);
304adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
305adc1d99fSHong Zhang   }
306adc1d99fSHong Zhang 
307397b6df1SKris Buschelman   PetscFunctionReturn(0);
308397b6df1SKris Buschelman }
309397b6df1SKris Buschelman 
310397b6df1SKris Buschelman #undef __FUNCT__
311f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
312f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
313397b6df1SKris Buschelman   int               ierr;
314397b6df1SKris Buschelman   PetscTruth        isascii;
315397b6df1SKris Buschelman   PetscViewerFormat format;
316f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
317397b6df1SKris Buschelman 
318397b6df1SKris Buschelman   PetscFunctionBegin;
319397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
320397b6df1SKris Buschelman 
321397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
322397b6df1SKris Buschelman   if (isascii) {
323397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
324397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
325397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
326397b6df1SKris Buschelman     }
327397b6df1SKris Buschelman   }
328397b6df1SKris Buschelman   PetscFunctionReturn(0);
329397b6df1SKris Buschelman }
330397b6df1SKris Buschelman 
331397b6df1SKris Buschelman #undef __FUNCT__
332f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
333f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
334f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
335d54de34fSKris Buschelman   PetscScalar *array;
336397b6df1SKris Buschelman   Vec         x_seq;
337397b6df1SKris Buschelman   IS          iden;
338397b6df1SKris Buschelman   VecScatter  scat;
339397b6df1SKris Buschelman   int         ierr;
340397b6df1SKris Buschelman 
341397b6df1SKris Buschelman   PetscFunctionBegin;
342397b6df1SKris Buschelman   if (lu->size > 1){
343397b6df1SKris Buschelman     if (!lu->myid){
344397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
345397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
346397b6df1SKris Buschelman     } else {
347397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
348397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
349397b6df1SKris Buschelman     }
350397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
351397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
352397b6df1SKris Buschelman 
353397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
354397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
355397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
356397b6df1SKris Buschelman   } else {  /* size == 1 */
357397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
358397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
359397b6df1SKris Buschelman   }
360397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
361397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
362397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
363397b6df1SKris Buschelman #else
364397b6df1SKris Buschelman     lu->id.rhs = array;
365397b6df1SKris Buschelman #endif
366397b6df1SKris Buschelman   }
367397b6df1SKris Buschelman 
368397b6df1SKris Buschelman   /* solve phase */
369397b6df1SKris Buschelman   lu->id.job=3;
370397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
371397b6df1SKris Buschelman   zmumps_c(&lu->id);
372397b6df1SKris Buschelman #else
373397b6df1SKris Buschelman   dmumps_c(&lu->id);
374397b6df1SKris Buschelman #endif
375397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
376397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
377397b6df1SKris Buschelman   }
378397b6df1SKris Buschelman 
379397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
380397b6df1SKris Buschelman   if (lu->size > 1) {
381397b6df1SKris Buschelman     if (!lu->myid){
382397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
383397b6df1SKris Buschelman     }
384397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
385397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
386397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
387397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
388397b6df1SKris Buschelman   } else {
389397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
390397b6df1SKris Buschelman   }
391397b6df1SKris Buschelman 
392397b6df1SKris Buschelman   PetscFunctionReturn(0);
393397b6df1SKris Buschelman }
394397b6df1SKris Buschelman 
395a58c3f20SHong Zhang /*
396a58c3f20SHong Zhang   input:
397a58c3f20SHong Zhang    F:        numeric factor
398a58c3f20SHong Zhang   output:
399a58c3f20SHong Zhang    nneg:     total number of negative pivots
400a58c3f20SHong Zhang    nzero:    0
401a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
402a58c3f20SHong Zhang */
403a58c3f20SHong Zhang 
404a58c3f20SHong Zhang #undef __FUNCT__
405a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
406a58c3f20SHong Zhang int MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
407a58c3f20SHong Zhang {
408a58c3f20SHong Zhang   Mat_MUMPS  *lu =(Mat_MUMPS*)F->spptr;
409bcb30aebSHong Zhang   int        ierr,neg,zero,pos,size;
410a58c3f20SHong Zhang 
411a58c3f20SHong Zhang   PetscFunctionBegin;
412bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
413bcb30aebSHong 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 */
414bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
415bcb30aebSHong 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));
416bcb30aebSHong Zhang   }
417a58c3f20SHong Zhang   if (nneg){
418a58c3f20SHong Zhang     if (!lu->myid){
419a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
420a58c3f20SHong Zhang     }
421bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
422a58c3f20SHong Zhang   }
423a58c3f20SHong Zhang   if (nzero) *nzero = 0;
424a58c3f20SHong Zhang   if (npos)  *npos  = F->M - (*nneg);
425a58c3f20SHong Zhang   PetscFunctionReturn(0);
426a58c3f20SHong Zhang }
427a58c3f20SHong Zhang 
428397b6df1SKris Buschelman #undef __FUNCT__
429f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
430f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
431f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
432f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
433397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
434397b6df1SKris Buschelman   PetscTruth valOnly,flg;
435397b6df1SKris Buschelman 
436397b6df1SKris Buschelman   PetscFunctionBegin;
437397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
438f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
439397b6df1SKris Buschelman 
440397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
441397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
442397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
44375747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
444397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
445397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
446397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
447397b6df1SKris Buschelman 
448397b6df1SKris Buschelman     /* Set mumps options */
449397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
450397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
451397b6df1SKris Buschelman     lu->id.sym=lu->sym;
452397b6df1SKris Buschelman     if (lu->sym == 2){
453397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
454397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
455397b6df1SKris Buschelman     }
456397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
457397b6df1SKris Buschelman   zmumps_c(&lu->id);
458397b6df1SKris Buschelman #else
459397b6df1SKris Buschelman   dmumps_c(&lu->id);
460397b6df1SKris Buschelman #endif
461397b6df1SKris Buschelman 
462397b6df1SKris Buschelman     if (lu->size == 1){
463397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
464397b6df1SKris Buschelman     } else {
465397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
466397b6df1SKris Buschelman     }
467397b6df1SKris Buschelman 
468397b6df1SKris Buschelman     icntl=-1;
469397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
470397b6df1SKris Buschelman     if (flg && icntl > 0) {
471397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
472397b6df1SKris Buschelman     } else { /* no output */
473397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
474397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
475397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
476397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
477397b6df1SKris Buschelman     }
478397b6df1SKris 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);
479397b6df1SKris Buschelman     icntl=-1;
480397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
481397b6df1SKris Buschelman     if (flg) {
482397b6df1SKris Buschelman       if (icntl== 1){
483397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
484397b6df1SKris Buschelman       } else {
485397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
486397b6df1SKris Buschelman       }
487397b6df1SKris Buschelman     }
488397b6df1SKris 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);
489397b6df1SKris 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);
49094b7f48cSBarry 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);
491397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
492397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
493adc1d99fSHong 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);
494397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
495397b6df1SKris Buschelman 
496397b6df1SKris Buschelman     /*
497397b6df1SKris 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);
498397b6df1SKris Buschelman     if (flg){
499397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
500397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
501397b6df1SKris Buschelman       } else {
502397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
503397b6df1SKris Buschelman       }
504397b6df1SKris Buschelman     }
505397b6df1SKris Buschelman     */
506397b6df1SKris Buschelman 
507397b6df1SKris 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);
508397b6df1SKris 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);
509397b6df1SKris 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);
510397b6df1SKris Buschelman     PetscOptionsEnd();
511397b6df1SKris Buschelman   }
512397b6df1SKris Buschelman 
513397b6df1SKris Buschelman   /* define matrix A */
514397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
515397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
516397b6df1SKris Buschelman     if (!lu->myid) {
517c36ead0aSKris Buschelman       if (lua->isAIJ){
518397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
519397b6df1SKris Buschelman         nz               = aa->nz;
520397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
521397b6df1SKris Buschelman       } else {
522397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5236c6c5352SBarry Smith         nz                  =  aa->nz;
524397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
525397b6df1SKris Buschelman       }
526397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
527397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
528397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
529397b6df1SKris Buschelman         nz = 0;
530397b6df1SKris Buschelman         for (i=0; i<M; i++){
531397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
532397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
533397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
534397b6df1SKris Buschelman           }
535397b6df1SKris Buschelman         }
536397b6df1SKris Buschelman       }
537397b6df1SKris Buschelman     }
538397b6df1SKris Buschelman     break;
539397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
540397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
541397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
542397b6df1SKris Buschelman     } else {
543397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
544397b6df1SKris Buschelman     }
545397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
546397b6df1SKris Buschelman     break;
547397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
548397b6df1SKris Buschelman   }
549397b6df1SKris Buschelman 
550397b6df1SKris Buschelman   /* analysis phase */
551397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
552397b6df1SKris Buschelman      lu->id.n = M;
553397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
554397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
555397b6df1SKris Buschelman       if (!lu->myid) {
556397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
557397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
558397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
559397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
560397b6df1SKris Buschelman #else
561397b6df1SKris Buschelman           lu->id.a = lu->val;
562397b6df1SKris Buschelman #endif
563397b6df1SKris Buschelman         }
564397b6df1SKris Buschelman       }
565397b6df1SKris Buschelman       break;
566397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
567397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
568397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
569397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
570397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
571397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
572397b6df1SKris Buschelman #else
573397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
574397b6df1SKris Buschelman #endif
575397b6df1SKris Buschelman       }
576397b6df1SKris Buschelman       break;
577397b6df1SKris Buschelman     }
578397b6df1SKris Buschelman     lu->id.job=1;
579397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
580397b6df1SKris Buschelman   zmumps_c(&lu->id);
581397b6df1SKris Buschelman #else
582397b6df1SKris Buschelman   dmumps_c(&lu->id);
583397b6df1SKris Buschelman #endif
584397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
585397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
586397b6df1SKris Buschelman     }
587397b6df1SKris Buschelman   }
588397b6df1SKris Buschelman 
589397b6df1SKris Buschelman   /* numerical factorization phase */
590397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
591a7aca84bSHong Zhang     if (!lu->myid) {
592397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
593397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
594397b6df1SKris Buschelman #else
595397b6df1SKris Buschelman       lu->id.a = lu->val;
596397b6df1SKris Buschelman #endif
597397b6df1SKris Buschelman     }
598397b6df1SKris Buschelman   } else {
599397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
600397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
601397b6df1SKris Buschelman #else
602397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
603397b6df1SKris Buschelman #endif
604397b6df1SKris Buschelman   }
605397b6df1SKris Buschelman   lu->id.job=2;
606397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
607397b6df1SKris Buschelman   zmumps_c(&lu->id);
608397b6df1SKris Buschelman #else
609397b6df1SKris Buschelman   dmumps_c(&lu->id);
610397b6df1SKris Buschelman #endif
611397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
612a7aca84bSHong 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));
613397b6df1SKris Buschelman   }
614397b6df1SKris Buschelman 
615397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
616397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
617397b6df1SKris Buschelman   }
618397b6df1SKris Buschelman 
619397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
620397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
621ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
622397b6df1SKris Buschelman   PetscFunctionReturn(0);
623397b6df1SKris Buschelman }
624397b6df1SKris Buschelman 
625397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
626397b6df1SKris Buschelman #undef __FUNCT__
627f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
628f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
629397b6df1SKris Buschelman   Mat       B;
630f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
631397b6df1SKris Buschelman   int       ierr;
632397b6df1SKris Buschelman 
633397b6df1SKris Buschelman   PetscFunctionBegin;
634397b6df1SKris Buschelman 
635397b6df1SKris Buschelman   /* Create the factorization matrix */
636397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
637397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
638397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
639397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
640397b6df1SKris Buschelman 
641f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
642397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
643f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
644397b6df1SKris Buschelman   lu->sym                 = 0;
645397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
646397b6df1SKris Buschelman 
647397b6df1SKris Buschelman   *F = B;
648397b6df1SKris Buschelman   PetscFunctionReturn(0);
649397b6df1SKris Buschelman }
650397b6df1SKris Buschelman 
651397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
652397b6df1SKris Buschelman #undef __FUNCT__
653f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
654f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
655397b6df1SKris Buschelman   Mat       B;
656f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
657397b6df1SKris Buschelman   int       ierr;
658397b6df1SKris Buschelman 
659397b6df1SKris Buschelman   PetscFunctionBegin;
660397b6df1SKris Buschelman 
661397b6df1SKris Buschelman   /* Create the factorization matrix */
662397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
663397b6df1SKris Buschelman   ierr = MatSetType(B,MATAIJMUMPS);CHKERRQ(ierr);
664397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
665397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
666397b6df1SKris Buschelman 
667f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
668a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
669397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
670f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
671397b6df1SKris Buschelman   lu->sym                       = 2;
672397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
673397b6df1SKris Buschelman 
674397b6df1SKris Buschelman   *F = B;
675397b6df1SKris Buschelman   PetscFunctionReturn(0);
676397b6df1SKris Buschelman }
677397b6df1SKris Buschelman 
678397b6df1SKris Buschelman #undef __FUNCT__
679f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
680f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
681c338a77dSKris Buschelman   int       ierr;
682f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
683c338a77dSKris Buschelman 
684397b6df1SKris Buschelman   PetscFunctionBegin;
685c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
686f0c56d0fSKris Buschelman 
687c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
688c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
689f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
690397b6df1SKris Buschelman   PetscFunctionReturn(0);
691397b6df1SKris Buschelman }
692397b6df1SKris Buschelman 
693c338a77dSKris Buschelman EXTERN_C_BEGIN
694c338a77dSKris Buschelman #undef __FUNCT__
695f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
6968e9aea5cSBarry Smith int MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
697c338a77dSKris Buschelman   int       ierr,size;
698c338a77dSKris Buschelman   MPI_Comm  comm;
699c338a77dSKris Buschelman   Mat       B=*newmat;
700f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
701397b6df1SKris Buschelman 
702397b6df1SKris Buschelman   PetscFunctionBegin;
703c338a77dSKris Buschelman   if (B != A) {
704c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
705397b6df1SKris Buschelman   }
706397b6df1SKris Buschelman 
707c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
708f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
709c338a77dSKris Buschelman 
710f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
711c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
712c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
713c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
714c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
715c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
716a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
717c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
718f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
719c338a77dSKris Buschelman 
7204b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
721f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_AIJMUMPS;
722f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
723f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
724f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
7253924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
726c338a77dSKris Buschelman 
727c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
728c338a77dSKris Buschelman   if (size == 1) {
729c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
730f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
731c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
732c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
733c338a77dSKris Buschelman   } else {
734c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
735f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
736c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
737c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
738c338a77dSKris Buschelman   }
739c338a77dSKris Buschelman 
740f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
741c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
742c338a77dSKris Buschelman   *newmat = B;
743397b6df1SKris Buschelman   PetscFunctionReturn(0);
744397b6df1SKris Buschelman }
745c338a77dSKris Buschelman EXTERN_C_END
746397b6df1SKris Buschelman 
747f0c56d0fSKris Buschelman #undef __FUNCT__
748f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_AIJMUMPS"
749f0c56d0fSKris Buschelman int MatDuplicate_AIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
750f0c56d0fSKris Buschelman   int       ierr;
7518f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
7528f340917SKris Buschelman 
753f0c56d0fSKris Buschelman   PetscFunctionBegin;
7548f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
755f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(*M,MATAIJMUMPS,M);CHKERRQ(ierr);
756a39386dcSKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
757f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
758f0c56d0fSKris Buschelman }
759f0c56d0fSKris Buschelman 
76024b6179bSKris Buschelman /*MC
761fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
76224b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
76324b6179bSKris Buschelman 
76424b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
76524b6179bSKris Buschelman   on how to declare the existence of external packages),
76624b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
76724b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
76824b6179bSKris Buschelman   This matrix type is only supported for double precision real.
76924b6179bSKris Buschelman 
77024b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
77124b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
77224b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
77324b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
77428b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
77528b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
77628b08bd3SKris Buschelman   without data copy.
77724b6179bSKris Buschelman 
77824b6179bSKris Buschelman   Options Database Keys:
7790bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
78024b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
78124b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
78224b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
78324b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
78424b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
78524b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
78694b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
78724b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
78824b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
78924b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
79024b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
79124b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
79224b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
79324b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
79424b6179bSKris Buschelman 
79524b6179bSKris Buschelman   Level: beginner
79624b6179bSKris Buschelman 
79724b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
79824b6179bSKris Buschelman M*/
79924b6179bSKris Buschelman 
800397b6df1SKris Buschelman EXTERN_C_BEGIN
801397b6df1SKris Buschelman #undef __FUNCT__
802f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
803f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
804397b6df1SKris Buschelman   int      ierr,size;
805e2d9671bSKris Buschelman   Mat      A_diag;
806397b6df1SKris Buschelman   MPI_Comm comm;
807397b6df1SKris Buschelman 
808397b6df1SKris Buschelman   PetscFunctionBegin;
8095441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
8105441df8eSKris Buschelman   /*   and AIJMUMPS types */
8115441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
812397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
813397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
814397b6df1SKris Buschelman   if (size == 1) {
815397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
816397b6df1SKris Buschelman   } else {
817397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
818e2d9671bSKris Buschelman     A_diag = ((Mat_MPIAIJ *)A->data)->A;
819e2d9671bSKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr);
820397b6df1SKris Buschelman   }
821f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
822397b6df1SKris Buschelman   PetscFunctionReturn(0);
823397b6df1SKris Buschelman }
824397b6df1SKris Buschelman EXTERN_C_END
825397b6df1SKris Buschelman 
826f579278aSKris Buschelman #undef __FUNCT__
827f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
828f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
829f579278aSKris Buschelman   int       ierr;
830f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
831f579278aSKris Buschelman 
832f579278aSKris Buschelman   PetscFunctionBegin;
833f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
834f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
835f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
836f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
837f579278aSKris Buschelman   PetscFunctionReturn(0);
838f579278aSKris Buschelman }
839f579278aSKris Buschelman 
840f579278aSKris Buschelman EXTERN_C_BEGIN
841f579278aSKris Buschelman #undef __FUNCT__
842*9c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
843*9c097c71SKris Buschelman int MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat  B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz)
844*9c097c71SKris Buschelman {
845*9c097c71SKris Buschelman   Mat       A;
846*9c097c71SKris Buschelman   Mat_MUMPS *mumps;
847*9c097c71SKris Buschelman   int       ierr;
848*9c097c71SKris Buschelman 
849*9c097c71SKris Buschelman   PetscFunctionBegin;
850*9c097c71SKris Buschelman   /*
851*9c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
852*9c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
853*9c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
854*9c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
855*9c097c71SKris Buschelman     in the preallocation routine.
856*9c097c71SKris Buschelman   */
857*9c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
858*9c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
859*9c097c71SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
860*9c097c71SKris Buschelman   PetscFunctionReturn(0);
861*9c097c71SKris Buschelman }
862*9c097c71SKris Buschelman EXTERN_C_END
863*9c097c71SKris Buschelman 
864*9c097c71SKris Buschelman EXTERN_C_BEGIN
865*9c097c71SKris Buschelman #undef __FUNCT__
866f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
8678e9aea5cSBarry Smith int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
868f579278aSKris Buschelman   int       ierr,size;
869f579278aSKris Buschelman   MPI_Comm  comm;
870f579278aSKris Buschelman   Mat       B=*newmat;
871f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
872*9c097c71SKris Buschelman   void      (*f)(void);
873f579278aSKris Buschelman 
874f579278aSKris Buschelman   PetscFunctionBegin;
875f579278aSKris Buschelman   if (B != A) {
876f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
877f579278aSKris Buschelman   }
878f579278aSKris Buschelman 
879f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
880f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
881f579278aSKris Buschelman 
882f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
883f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
884f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
885f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
886f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
887f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
888a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
889f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
890f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
891f579278aSKris Buschelman 
892f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
893f0c56d0fSKris Buschelman   B->ops->duplicate                = MatDuplicate_SBAIJMUMPS;
894f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
895f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
896f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
8973924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
898f579278aSKris Buschelman 
899f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
900f579278aSKris Buschelman   if (size == 1) {
901f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
902f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
903f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
904f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
905f579278aSKris Buschelman   } else {
906*9c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
907*9c097c71SKris Buschelman     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
908*9c097c71SKris Buschelman     if (f) {
909*9c097c71SKris Buschelman       mumps->MatPreallocate = (int (*)(Mat,int,int,int*,int,int*))f;
910*9c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
911*9c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
912*9c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
913*9c097c71SKris Buschelman     }
914f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
915f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
916f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
917f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
918f579278aSKris Buschelman   }
919f579278aSKris Buschelman 
920f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
921f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
922f579278aSKris Buschelman   *newmat = B;
923f579278aSKris Buschelman   PetscFunctionReturn(0);
924f579278aSKris Buschelman }
925f579278aSKris Buschelman EXTERN_C_END
926f579278aSKris Buschelman 
927f0c56d0fSKris Buschelman #undef __FUNCT__
928f0c56d0fSKris Buschelman #define __FUNCT__ "MatDuplicate_SBAIJMUMPS"
929f0c56d0fSKris Buschelman int MatDuplicate_SBAIJMUMPS(Mat A, MatDuplicateOption op, Mat *M) {
930f0c56d0fSKris Buschelman   int       ierr;
9318f340917SKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS *)A->spptr;
9328f340917SKris Buschelman 
933f0c56d0fSKris Buschelman   PetscFunctionBegin;
9348f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
935f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(*M,MATSBAIJMUMPS,M);CHKERRQ(ierr);
9363f953163SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
937f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
938f0c56d0fSKris Buschelman }
939f0c56d0fSKris Buschelman 
94024b6179bSKris Buschelman /*MC
941fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
94224b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
94324b6179bSKris Buschelman 
94424b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
94524b6179bSKris Buschelman   on how to declare the existence of external packages),
94624b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
94724b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
94824b6179bSKris Buschelman   This matrix type is only supported for double precision real.
94924b6179bSKris Buschelman 
95024b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
95124b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
95224b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
95324b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
95428b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
95528b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
95628b08bd3SKris Buschelman   without data copy.
95724b6179bSKris Buschelman 
95824b6179bSKris Buschelman   Options Database Keys:
9590bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
96024b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
96124b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
96224b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
96324b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
96424b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
96524b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
96694b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
96724b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
96824b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
96924b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
97024b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
97124b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
97224b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
97324b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
97424b6179bSKris Buschelman 
97524b6179bSKris Buschelman   Level: beginner
97624b6179bSKris Buschelman 
97724b6179bSKris Buschelman .seealso: MATAIJMUMPS
97824b6179bSKris Buschelman M*/
97924b6179bSKris Buschelman 
980397b6df1SKris Buschelman EXTERN_C_BEGIN
981397b6df1SKris Buschelman #undef __FUNCT__
982f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
983f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
984397b6df1SKris Buschelman   int ierr,size;
985397b6df1SKris Buschelman 
986397b6df1SKris Buschelman   PetscFunctionBegin;
9875441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
9885441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
9895441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
9905441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
991397b6df1SKris Buschelman   if (size == 1) {
992397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
993397b6df1SKris Buschelman   } else {
994397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
995397b6df1SKris Buschelman   }
996f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
997397b6df1SKris Buschelman   PetscFunctionReturn(0);
998397b6df1SKris Buschelman }
999397b6df1SKris Buschelman EXTERN_C_END
1000