xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 8e393735c4facc1ad349ab70cfcf515eeda33366)
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 
49422e82a1SHong Zhang EXTERN int MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
50892f6c3fSKris Buschelman EXTERN_C_BEGIN
51892f6c3fSKris Buschelman int MatConvert_SBAIJ_SBAIJMUMPS(Mat,const MatType,Mat*);
52892f6c3fSKris Buschelman EXTERN_C_END
53397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
54397b6df1SKris Buschelman /*
55397b6df1SKris Buschelman   input:
5675747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
57397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
58397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
59397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
60397b6df1SKris Buschelman   output:
61397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
62397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
63397b6df1SKris Buschelman  */
64f0c56d0fSKris Buschelman int MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
65397b6df1SKris Buschelman   int         *ai, *aj, *bi, *bj, rstart,nz, *garray;
66397b6df1SKris Buschelman   int         ierr,i,j,jj,jB,irow,m=A->m,*ajj,*bjj,countA,countB,colA_start,jcol;
67d54de34fSKris Buschelman   int         *row,*col;
68397b6df1SKris Buschelman   PetscScalar *av, *bv,*val;
69f0c56d0fSKris Buschelman   Mat_MUMPS   *mumps=(Mat_MUMPS*)A->spptr;
70397b6df1SKris Buschelman 
71397b6df1SKris Buschelman   PetscFunctionBegin;
72397b6df1SKris Buschelman   if (mumps->isAIJ){
73397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
74397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
75397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
76397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
77397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
78397b6df1SKris Buschelman     garray = mat->garray;
79397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
80397b6df1SKris Buschelman 
81397b6df1SKris Buschelman   } else {
82397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
83397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
84397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
85847143adSKris Buschelman     if (mat->bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", mat->bs);
866c6c5352SBarry Smith     nz = aa->nz + bb->nz;
87397b6df1SKris Buschelman     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= mat->rstart;
88397b6df1SKris Buschelman     garray = mat->garray;
89397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
90397b6df1SKris Buschelman   }
91397b6df1SKris Buschelman 
92397b6df1SKris Buschelman   if (!valOnly){
93397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&row);CHKERRQ(ierr);
94397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(int),&col);CHKERRQ(ierr);
95397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
96397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
97397b6df1SKris Buschelman   } else {
98397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
99397b6df1SKris Buschelman   }
100397b6df1SKris Buschelman   *nnz = nz;
101397b6df1SKris Buschelman 
102028e57e8SHong Zhang   jj = 0; irow = rstart;
103397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
104397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
105397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
106397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
107397b6df1SKris Buschelman     bjj = bj + bi[i];
108397b6df1SKris Buschelman 
109397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
110397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11175747be1SHong Zhang     j=-1;
11275747be1SHong Zhang     do {
11375747be1SHong Zhang       j++;
11475747be1SHong Zhang       if (j == countB) break;
115397b6df1SKris Buschelman       jcol = garray[bjj[j]];
11675747be1SHong Zhang     } while (jcol < colA_start);
11775747be1SHong Zhang     jB = j;
118397b6df1SKris Buschelman 
119397b6df1SKris Buschelman     /* B-part, smaller col index */
120397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
121397b6df1SKris Buschelman     for (j=0; j<jB; j++){
122397b6df1SKris Buschelman       jcol = garray[bjj[j]];
123397b6df1SKris Buschelman       if (!valOnly){
124397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12575747be1SHong Zhang 
126397b6df1SKris Buschelman       }
127397b6df1SKris Buschelman       val[jj++] = *bv++;
128397b6df1SKris Buschelman     }
129397b6df1SKris Buschelman     /* A-part */
130397b6df1SKris Buschelman     for (j=0; j<countA; j++){
131397b6df1SKris Buschelman       if (!valOnly){
132397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
133397b6df1SKris Buschelman       }
134397b6df1SKris Buschelman       val[jj++] = *av++;
135397b6df1SKris Buschelman     }
136397b6df1SKris Buschelman     /* B-part, larger col index */
137397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
138397b6df1SKris Buschelman       if (!valOnly){
139397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
140397b6df1SKris Buschelman       }
141397b6df1SKris Buschelman       val[jj++] = *bv++;
142397b6df1SKris Buschelman     }
143397b6df1SKris Buschelman     irow++;
144397b6df1SKris Buschelman   }
145397b6df1SKris Buschelman 
146397b6df1SKris Buschelman   PetscFunctionReturn(0);
147397b6df1SKris Buschelman }
148397b6df1SKris Buschelman 
149c338a77dSKris Buschelman EXTERN_C_BEGIN
150c338a77dSKris Buschelman #undef __FUNCT__
151c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
1528e9aea5cSBarry Smith int MatConvert_MUMPS_Base(Mat A,const MatType type,Mat *newmat) {
153c338a77dSKris Buschelman   int       ierr;
154c338a77dSKris Buschelman   Mat       B=*newmat;
155f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
156c338a77dSKris Buschelman 
157c338a77dSKris Buschelman   PetscFunctionBegin;
158c338a77dSKris Buschelman   if (B != A) {
159c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
160c338a77dSKris Buschelman   }
161f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
162f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
163f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
164f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
165f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
166f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
1673924e44cSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
168f0c56d0fSKris Buschelman   ierr = PetscFree(mumps);CHKERRQ(ierr);
169c338a77dSKris Buschelman   *newmat = B;
170c338a77dSKris Buschelman   PetscFunctionReturn(0);
171c338a77dSKris Buschelman }
172c338a77dSKris Buschelman EXTERN_C_END
173c338a77dSKris Buschelman 
174397b6df1SKris Buschelman #undef __FUNCT__
1753924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
1763924e44cSKris Buschelman int MatDestroy_MUMPS(Mat A) {
177f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
178c338a77dSKris Buschelman   int       ierr,size=lu->size;
179a39386dcSKris Buschelman   int       (*specialdestroy)(Mat);
180397b6df1SKris Buschelman   PetscFunctionBegin;
181397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
182397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
183397b6df1SKris Buschelman     lu->id.job=JOB_END;
184397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
185397b6df1SKris Buschelman     zmumps_c(&lu->id);
186397b6df1SKris Buschelman #else
187397b6df1SKris Buschelman     dmumps_c(&lu->id);
188397b6df1SKris Buschelman #endif
189c338a77dSKris Buschelman     if (lu->irn) {
190c338a77dSKris Buschelman       ierr = PetscFree(lu->irn);CHKERRQ(ierr);
191c338a77dSKris Buschelman     }
192c338a77dSKris Buschelman     if (lu->jcn) {
193c338a77dSKris Buschelman       ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
194c338a77dSKris Buschelman     }
195c338a77dSKris Buschelman     if (size>1 && lu->val) {
196c338a77dSKris Buschelman       ierr = PetscFree(lu->val);CHKERRQ(ierr);
197c338a77dSKris Buschelman     }
198397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
199397b6df1SKris Buschelman   }
200a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
201a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
202c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
203397b6df1SKris Buschelman   PetscFunctionReturn(0);
204397b6df1SKris Buschelman }
205397b6df1SKris Buschelman 
206397b6df1SKris Buschelman #undef __FUNCT__
207a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
208a39386dcSKris Buschelman int MatDestroy_AIJMUMPS(Mat A) {
209a39386dcSKris Buschelman   int ierr, size;
210a39386dcSKris Buschelman 
211a39386dcSKris Buschelman   PetscFunctionBegin;
212a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
213a39386dcSKris Buschelman   if (size==1) {
214a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,&A);CHKERRQ(ierr);
215a39386dcSKris Buschelman   } else {
216a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,&A);CHKERRQ(ierr);
217a39386dcSKris Buschelman   }
218a39386dcSKris Buschelman   PetscFunctionReturn(0);
219a39386dcSKris Buschelman }
220a39386dcSKris Buschelman 
221a39386dcSKris Buschelman #undef __FUNCT__
222a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
223a39386dcSKris Buschelman int MatDestroy_SBAIJMUMPS(Mat A) {
224a39386dcSKris Buschelman   int ierr, size;
225a39386dcSKris Buschelman 
226a39386dcSKris Buschelman   PetscFunctionBegin;
227a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
228a39386dcSKris Buschelman   if (size==1) {
229a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,&A);CHKERRQ(ierr);
230a39386dcSKris Buschelman   } else {
231a39386dcSKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,&A);CHKERRQ(ierr);
232a39386dcSKris Buschelman   }
233a39386dcSKris Buschelman   PetscFunctionReturn(0);
234a39386dcSKris Buschelman }
235a39386dcSKris Buschelman 
236a39386dcSKris Buschelman #undef __FUNCT__
237c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
238f0c56d0fSKris Buschelman int MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
239f0c56d0fSKris Buschelman   Mat_MUMPS *lu=(Mat_MUMPS*)A->spptr;
240397b6df1SKris Buschelman   int       ierr;
241397b6df1SKris Buschelman 
242397b6df1SKris Buschelman   PetscFunctionBegin;
243c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
244c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
245c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
246c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
247c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
248c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
249c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
250c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
251c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
252c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
253c338a77dSKris Buschelman   if (lu->myid == 0 && lu->id.ICNTL(11)>0) {
254c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
255c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
256c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
257c338a77dSKris 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);
258c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
259c338a77dSKris 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);
260c338a77dSKris Buschelman 
261c338a77dSKris Buschelman   }
262c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
263c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
264adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
265c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(15) (efficiency control):                         %d \n",lu->id.ICNTL(15));CHKERRQ(ierr);
266c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
267c338a77dSKris Buschelman 
268c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
269c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
270c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
27157f0c58bSHong Zhang 
27257f0c58bSHong Zhang   /* infomation local to each processor */
27357f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
27457f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
27557f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27657f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
27757f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
27857f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
27957f0c58bSHong Zhang   if (lu->myid == 0) ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
28057f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
28157f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
282adc1d99fSHong Zhang 
283adc1d99fSHong Zhang   if (lu->myid == 0){ /* information from the host */
284adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
285adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
286adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
287adc1d99fSHong Zhang 
288adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
289adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
290adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
291adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
292adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
293adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
294adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real space store the matrix factors after analysis): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
295adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after analysis): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
296adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
297adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
298adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
299adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
300adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
301adc1d99fSHong 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);
302adc1d99fSHong 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);
303adc1d99fSHong 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);
304adc1d99fSHong 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);
305adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
306adc1d99fSHong Zhang   }
307adc1d99fSHong Zhang 
308397b6df1SKris Buschelman   PetscFunctionReturn(0);
309397b6df1SKris Buschelman }
310397b6df1SKris Buschelman 
311397b6df1SKris Buschelman #undef __FUNCT__
312f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_AIJMUMPS"
313f0c56d0fSKris Buschelman int MatView_AIJMUMPS(Mat A,PetscViewer viewer) {
314397b6df1SKris Buschelman   int               ierr;
315397b6df1SKris Buschelman   PetscTruth        isascii;
316397b6df1SKris Buschelman   PetscViewerFormat format;
317f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
318397b6df1SKris Buschelman 
319397b6df1SKris Buschelman   PetscFunctionBegin;
320397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
321397b6df1SKris Buschelman 
322397b6df1SKris Buschelman   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr);
323397b6df1SKris Buschelman   if (isascii) {
324397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
325397b6df1SKris Buschelman     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
326397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
327397b6df1SKris Buschelman     }
328397b6df1SKris Buschelman   }
329397b6df1SKris Buschelman   PetscFunctionReturn(0);
330397b6df1SKris Buschelman }
331397b6df1SKris Buschelman 
332397b6df1SKris Buschelman #undef __FUNCT__
333f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
334f0c56d0fSKris Buschelman int MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
335f0c56d0fSKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS*)A->spptr;
336d54de34fSKris Buschelman   PetscScalar *array;
337397b6df1SKris Buschelman   Vec         x_seq;
338397b6df1SKris Buschelman   IS          iden;
339397b6df1SKris Buschelman   VecScatter  scat;
340397b6df1SKris Buschelman   int         ierr;
341397b6df1SKris Buschelman 
342397b6df1SKris Buschelman   PetscFunctionBegin;
343397b6df1SKris Buschelman   if (lu->size > 1){
344397b6df1SKris Buschelman     if (!lu->myid){
345397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,A->N,&x_seq);CHKERRQ(ierr);
346397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,A->N,0,1,&iden);CHKERRQ(ierr);
347397b6df1SKris Buschelman     } else {
348397b6df1SKris Buschelman       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&x_seq);CHKERRQ(ierr);
349397b6df1SKris Buschelman       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&iden);CHKERRQ(ierr);
350397b6df1SKris Buschelman     }
351397b6df1SKris Buschelman     ierr = VecScatterCreate(b,iden,x_seq,iden,&scat);CHKERRQ(ierr);
352397b6df1SKris Buschelman     ierr = ISDestroy(iden);CHKERRQ(ierr);
353397b6df1SKris Buschelman 
354397b6df1SKris Buschelman     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
355397b6df1SKris Buschelman     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat);CHKERRQ(ierr);
356397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
357397b6df1SKris Buschelman   } else {  /* size == 1 */
358397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
359397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
360397b6df1SKris Buschelman   }
361397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
362397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
363397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
364397b6df1SKris Buschelman #else
365397b6df1SKris Buschelman     lu->id.rhs = array;
366397b6df1SKris Buschelman #endif
367397b6df1SKris Buschelman   }
368397b6df1SKris Buschelman 
369397b6df1SKris Buschelman   /* solve phase */
370397b6df1SKris Buschelman   lu->id.job=3;
371397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
372397b6df1SKris Buschelman   zmumps_c(&lu->id);
373397b6df1SKris Buschelman #else
374397b6df1SKris Buschelman   dmumps_c(&lu->id);
375397b6df1SKris Buschelman #endif
376397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
377397b6df1SKris Buschelman     SETERRQ1(1,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
378397b6df1SKris Buschelman   }
379397b6df1SKris Buschelman 
380397b6df1SKris Buschelman   /* convert mumps solution x_seq to petsc mpi x */
381397b6df1SKris Buschelman   if (lu->size > 1) {
382397b6df1SKris Buschelman     if (!lu->myid){
383397b6df1SKris Buschelman       ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
384397b6df1SKris Buschelman     }
385397b6df1SKris Buschelman     ierr = VecScatterBegin(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
386397b6df1SKris Buschelman     ierr = VecScatterEnd(x_seq,x,INSERT_VALUES,SCATTER_REVERSE,scat);CHKERRQ(ierr);
387397b6df1SKris Buschelman     ierr = VecScatterDestroy(scat);CHKERRQ(ierr);
388397b6df1SKris Buschelman     ierr = VecDestroy(x_seq);CHKERRQ(ierr);
389397b6df1SKris Buschelman   } else {
390397b6df1SKris Buschelman     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
391397b6df1SKris Buschelman   }
392397b6df1SKris Buschelman 
393397b6df1SKris Buschelman   PetscFunctionReturn(0);
394397b6df1SKris Buschelman }
395397b6df1SKris Buschelman 
396a58c3f20SHong Zhang /*
397a58c3f20SHong Zhang   input:
398a58c3f20SHong Zhang    F:        numeric factor
399a58c3f20SHong Zhang   output:
400a58c3f20SHong Zhang    nneg:     total number of negative pivots
401a58c3f20SHong Zhang    nzero:    0
402a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
403a58c3f20SHong Zhang */
404a58c3f20SHong Zhang 
405a58c3f20SHong Zhang #undef __FUNCT__
406a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
407a58c3f20SHong Zhang int MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
408a58c3f20SHong Zhang {
409a58c3f20SHong Zhang   Mat_MUMPS  *lu =(Mat_MUMPS*)F->spptr;
410bcb30aebSHong Zhang   int        ierr,neg,zero,pos,size;
411a58c3f20SHong Zhang 
412a58c3f20SHong Zhang   PetscFunctionBegin;
413bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
414bcb30aebSHong 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 */
415bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
416bcb30aebSHong 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));
417bcb30aebSHong Zhang   }
418a58c3f20SHong Zhang   if (nneg){
419a58c3f20SHong Zhang     if (!lu->myid){
420a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
421a58c3f20SHong Zhang     }
422bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
423a58c3f20SHong Zhang   }
424a58c3f20SHong Zhang   if (nzero) *nzero = 0;
425a58c3f20SHong Zhang   if (npos)  *npos  = F->M - (*nneg);
426a58c3f20SHong Zhang   PetscFunctionReturn(0);
427a58c3f20SHong Zhang }
428a58c3f20SHong Zhang 
429397b6df1SKris Buschelman #undef __FUNCT__
430f0c56d0fSKris Buschelman #define __FUNCT__ "MatFactorNumeric_MPIAIJMUMPS"
431f0c56d0fSKris Buschelman int MatFactorNumeric_AIJMUMPS(Mat A,Mat *F) {
432f0c56d0fSKris Buschelman   Mat_MUMPS  *lu =(Mat_MUMPS*)(*F)->spptr;
433f0c56d0fSKris Buschelman   Mat_MUMPS  *lua=(Mat_MUMPS*)(A)->spptr;
434397b6df1SKris Buschelman   int        rnz,nnz,ierr,nz,i,M=A->M,*ai,*aj,icntl;
435397b6df1SKris Buschelman   PetscTruth valOnly,flg;
436397b6df1SKris Buschelman 
437397b6df1SKris Buschelman   PetscFunctionBegin;
438397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
439f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
440397b6df1SKris Buschelman 
441397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
442397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
443397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
44475747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
445397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
446397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
447397b6df1SKris Buschelman     lu->id.comm_fortran = lu->comm_mumps;
448397b6df1SKris Buschelman 
449397b6df1SKris Buschelman     /* Set mumps options */
450397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
451397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
452397b6df1SKris Buschelman     lu->id.sym=lu->sym;
453397b6df1SKris Buschelman     if (lu->sym == 2){
454397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
455397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
456397b6df1SKris Buschelman     }
457397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
458397b6df1SKris Buschelman   zmumps_c(&lu->id);
459397b6df1SKris Buschelman #else
460397b6df1SKris Buschelman   dmumps_c(&lu->id);
461397b6df1SKris Buschelman #endif
462397b6df1SKris Buschelman 
463397b6df1SKris Buschelman     if (lu->size == 1){
464397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
465397b6df1SKris Buschelman     } else {
466397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
467397b6df1SKris Buschelman     }
468397b6df1SKris Buschelman 
469397b6df1SKris Buschelman     icntl=-1;
470397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
471397b6df1SKris Buschelman     if (flg && icntl > 0) {
472397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
473397b6df1SKris Buschelman     } else { /* no output */
474397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
475397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
476397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
477397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
478397b6df1SKris Buschelman     }
479397b6df1SKris 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);
480397b6df1SKris Buschelman     icntl=-1;
481397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
482397b6df1SKris Buschelman     if (flg) {
483397b6df1SKris Buschelman       if (icntl== 1){
484397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
485397b6df1SKris Buschelman       } else {
486397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
487397b6df1SKris Buschelman       }
488397b6df1SKris Buschelman     }
489397b6df1SKris 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);
490397b6df1SKris 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);
49194b7f48cSBarry 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);
492397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
493397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
494adc1d99fSHong 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);
495397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
496397b6df1SKris Buschelman 
497397b6df1SKris Buschelman     /*
498397b6df1SKris 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);
499397b6df1SKris Buschelman     if (flg){
500397b6df1SKris Buschelman       if (icntl >-1 && icntl <3 ){
501397b6df1SKris Buschelman         if (lu->myid==0) lu->id.ICNTL(16) = icntl;
502397b6df1SKris Buschelman       } else {
503397b6df1SKris Buschelman         SETERRQ1(PETSC_ERR_SUP,"ICNTL(16)=%d -- not supported\n",icntl);
504397b6df1SKris Buschelman       }
505397b6df1SKris Buschelman     }
506397b6df1SKris Buschelman     */
507397b6df1SKris Buschelman 
508397b6df1SKris 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);
509397b6df1SKris 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);
510397b6df1SKris 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);
511397b6df1SKris Buschelman     PetscOptionsEnd();
512397b6df1SKris Buschelman   }
513397b6df1SKris Buschelman 
514397b6df1SKris Buschelman   /* define matrix A */
515397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
516397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
517397b6df1SKris Buschelman     if (!lu->myid) {
518c36ead0aSKris Buschelman       if (lua->isAIJ){
519397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
520397b6df1SKris Buschelman         nz               = aa->nz;
521397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
522397b6df1SKris Buschelman       } else {
523397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5246c6c5352SBarry Smith         nz                  =  aa->nz;
525397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
526397b6df1SKris Buschelman       }
527397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
528397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->irn);CHKERRQ(ierr);
529397b6df1SKris Buschelman         ierr = PetscMalloc(nz*sizeof(int),&lu->jcn);CHKERRQ(ierr);
530397b6df1SKris Buschelman         nz = 0;
531397b6df1SKris Buschelman         for (i=0; i<M; i++){
532397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
533397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
534397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
535397b6df1SKris Buschelman           }
536397b6df1SKris Buschelman         }
537397b6df1SKris Buschelman       }
538397b6df1SKris Buschelman     }
539397b6df1SKris Buschelman     break;
540397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
541397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
542397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
543397b6df1SKris Buschelman     } else {
544397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
545397b6df1SKris Buschelman     }
546397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
547397b6df1SKris Buschelman     break;
548397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
549397b6df1SKris Buschelman   }
550397b6df1SKris Buschelman 
551397b6df1SKris Buschelman   /* analysis phase */
552397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
553397b6df1SKris Buschelman      lu->id.n = M;
554397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
555397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
556397b6df1SKris Buschelman       if (!lu->myid) {
557397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
558397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
559397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
560397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
561397b6df1SKris Buschelman #else
562397b6df1SKris Buschelman           lu->id.a = lu->val;
563397b6df1SKris Buschelman #endif
564397b6df1SKris Buschelman         }
565397b6df1SKris Buschelman       }
566397b6df1SKris Buschelman       break;
567397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
568397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
569397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
570397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
571397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
572397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
573397b6df1SKris Buschelman #else
574397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
575397b6df1SKris Buschelman #endif
576397b6df1SKris Buschelman       }
577397b6df1SKris Buschelman       break;
578397b6df1SKris Buschelman     }
579397b6df1SKris Buschelman     lu->id.job=1;
580397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
581397b6df1SKris Buschelman   zmumps_c(&lu->id);
582397b6df1SKris Buschelman #else
583397b6df1SKris Buschelman   dmumps_c(&lu->id);
584397b6df1SKris Buschelman #endif
585397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
586397b6df1SKris Buschelman       SETERRQ1(1,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
587397b6df1SKris Buschelman     }
588397b6df1SKris Buschelman   }
589397b6df1SKris Buschelman 
590397b6df1SKris Buschelman   /* numerical factorization phase */
591397b6df1SKris Buschelman   if(lu->id.ICNTL(18) == 0) {
592a7aca84bSHong Zhang     if (!lu->myid) {
593397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
594397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
595397b6df1SKris Buschelman #else
596397b6df1SKris Buschelman       lu->id.a = lu->val;
597397b6df1SKris Buschelman #endif
598397b6df1SKris Buschelman     }
599397b6df1SKris Buschelman   } else {
600397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
601397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
602397b6df1SKris Buschelman #else
603397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
604397b6df1SKris Buschelman #endif
605397b6df1SKris Buschelman   }
606397b6df1SKris Buschelman   lu->id.job=2;
607397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
608397b6df1SKris Buschelman   zmumps_c(&lu->id);
609397b6df1SKris Buschelman #else
610397b6df1SKris Buschelman   dmumps_c(&lu->id);
611397b6df1SKris Buschelman #endif
612397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
613a7aca84bSHong 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));
614397b6df1SKris Buschelman   }
615397b6df1SKris Buschelman 
616397b6df1SKris Buschelman   if (lu->myid==0 && lu->id.ICNTL(16) > 0){
617397b6df1SKris Buschelman     SETERRQ1(1,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
618397b6df1SKris Buschelman   }
619397b6df1SKris Buschelman 
620397b6df1SKris Buschelman   (*F)->assembled  = PETSC_TRUE;
621397b6df1SKris Buschelman   lu->matstruc     = SAME_NONZERO_PATTERN;
622ace87b0dSHong Zhang   lu->CleanUpMUMPS = PETSC_TRUE;
623397b6df1SKris Buschelman   PetscFunctionReturn(0);
624397b6df1SKris Buschelman }
625397b6df1SKris Buschelman 
626397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
627397b6df1SKris Buschelman #undef __FUNCT__
628f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
629f0c56d0fSKris Buschelman int MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
630397b6df1SKris Buschelman   Mat       B;
631f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
632397b6df1SKris Buschelman   int       ierr;
633397b6df1SKris Buschelman 
634397b6df1SKris Buschelman   PetscFunctionBegin;
635397b6df1SKris Buschelman 
636397b6df1SKris Buschelman   /* Create the factorization matrix */
637397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
638be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
639397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
640397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
641397b6df1SKris Buschelman 
642f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
643397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
644f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
645397b6df1SKris Buschelman   lu->sym                 = 0;
646397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
647397b6df1SKris Buschelman 
648397b6df1SKris Buschelman   *F = B;
649397b6df1SKris Buschelman   PetscFunctionReturn(0);
650397b6df1SKris Buschelman }
651397b6df1SKris Buschelman 
652397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
653397b6df1SKris Buschelman #undef __FUNCT__
654f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
655f0c56d0fSKris Buschelman int MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
656397b6df1SKris Buschelman   Mat       B;
657f0c56d0fSKris Buschelman   Mat_MUMPS *lu;
658397b6df1SKris Buschelman   int       ierr;
659397b6df1SKris Buschelman 
660397b6df1SKris Buschelman   PetscFunctionBegin;
661397b6df1SKris Buschelman 
662397b6df1SKris Buschelman   /* Create the factorization matrix */
663397b6df1SKris Buschelman   ierr = MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);CHKERRQ(ierr);
664be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
665efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
666efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
667397b6df1SKris Buschelman 
668f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
669a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
670397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
671f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
672397b6df1SKris Buschelman   lu->sym                       = 2;
673397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
674397b6df1SKris Buschelman 
675397b6df1SKris Buschelman   *F = B;
676397b6df1SKris Buschelman   PetscFunctionReturn(0);
677397b6df1SKris Buschelman }
678397b6df1SKris Buschelman 
679397b6df1SKris Buschelman #undef __FUNCT__
680f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
681f0c56d0fSKris Buschelman int MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
682c338a77dSKris Buschelman   int       ierr;
683f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
684c338a77dSKris Buschelman 
685397b6df1SKris Buschelman   PetscFunctionBegin;
686c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
687f0c56d0fSKris Buschelman 
688c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
689c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
690f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
691397b6df1SKris Buschelman   PetscFunctionReturn(0);
692397b6df1SKris Buschelman }
693397b6df1SKris Buschelman 
694c338a77dSKris Buschelman EXTERN_C_BEGIN
695c338a77dSKris Buschelman #undef __FUNCT__
696f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
6978e9aea5cSBarry Smith int MatConvert_AIJ_AIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
698c338a77dSKris Buschelman   int       ierr,size;
699c338a77dSKris Buschelman   MPI_Comm  comm;
700c338a77dSKris Buschelman   Mat       B=*newmat;
701f0c56d0fSKris Buschelman   Mat_MUMPS *mumps;
702397b6df1SKris Buschelman 
703397b6df1SKris Buschelman   PetscFunctionBegin;
704c338a77dSKris Buschelman   if (B != A) {
705c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
706397b6df1SKris Buschelman   }
707397b6df1SKris Buschelman 
708c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
709f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
710c338a77dSKris Buschelman 
711f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
712c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
713c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
714c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
715c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
716c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
717a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
718c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
719f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
720c338a77dSKris Buschelman 
7214b68dd72SKris Buschelman   B->spptr                         = (void *)mumps;
722422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
723f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
724f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
725f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
7263924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
727c338a77dSKris Buschelman 
728c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
729c338a77dSKris Buschelman   if (size == 1) {
730c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
731f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
732c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
733c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
734c338a77dSKris Buschelman   } else {
735c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
736f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
737c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
738c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
739c338a77dSKris Buschelman   }
740c338a77dSKris Buschelman 
741f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for LU factorization and solves.");
742c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
743c338a77dSKris Buschelman   *newmat = B;
744397b6df1SKris Buschelman   PetscFunctionReturn(0);
745397b6df1SKris Buschelman }
746c338a77dSKris Buschelman EXTERN_C_END
747397b6df1SKris Buschelman 
74824b6179bSKris Buschelman /*MC
749fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
75024b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
75124b6179bSKris Buschelman 
75224b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
75324b6179bSKris Buschelman   on how to declare the existence of external packages),
75424b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
75524b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
75624b6179bSKris Buschelman   This matrix type is only supported for double precision real.
75724b6179bSKris Buschelman 
75824b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
75924b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
76024b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
76124b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
76228b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
76328b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
76428b08bd3SKris Buschelman   without data copy.
76524b6179bSKris Buschelman 
76624b6179bSKris Buschelman   Options Database Keys:
7670bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
76824b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
76924b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
77024b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
77124b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
77224b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
77324b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
77494b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
77524b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
77624b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
77724b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
77824b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
77924b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
78024b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
78124b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
78224b6179bSKris Buschelman 
78324b6179bSKris Buschelman   Level: beginner
78424b6179bSKris Buschelman 
78524b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
78624b6179bSKris Buschelman M*/
78724b6179bSKris Buschelman 
788397b6df1SKris Buschelman EXTERN_C_BEGIN
789397b6df1SKris Buschelman #undef __FUNCT__
790f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
791f0c56d0fSKris Buschelman int MatCreate_AIJMUMPS(Mat A) {
792397b6df1SKris Buschelman   int      ierr,size;
793e2d9671bSKris Buschelman   Mat      A_diag;
794397b6df1SKris Buschelman   MPI_Comm comm;
795397b6df1SKris Buschelman 
796397b6df1SKris Buschelman   PetscFunctionBegin;
7975441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
7985441df8eSKris Buschelman   /*   and AIJMUMPS types */
7995441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
800397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
801397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
802397b6df1SKris Buschelman   if (size == 1) {
803397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
804397b6df1SKris Buschelman   } else {
805397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
806e2d9671bSKris Buschelman     A_diag = ((Mat_MPIAIJ *)A->data)->A;
807e2d9671bSKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,&A_diag);CHKERRQ(ierr);
808397b6df1SKris Buschelman   }
809f0c56d0fSKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,&A);CHKERRQ(ierr);
810397b6df1SKris Buschelman   PetscFunctionReturn(0);
811397b6df1SKris Buschelman }
812397b6df1SKris Buschelman EXTERN_C_END
813397b6df1SKris Buschelman 
814f579278aSKris Buschelman #undef __FUNCT__
815f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
816f0c56d0fSKris Buschelman int MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode) {
817f579278aSKris Buschelman   int       ierr;
818f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
819f579278aSKris Buschelman 
820f579278aSKris Buschelman   PetscFunctionBegin;
821f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
822f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
823f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
824f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
825f579278aSKris Buschelman   PetscFunctionReturn(0);
826f579278aSKris Buschelman }
827f579278aSKris Buschelman 
828f579278aSKris Buschelman EXTERN_C_BEGIN
829f579278aSKris Buschelman #undef __FUNCT__
8309c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
8319c097c71SKris Buschelman int MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat  B,int bs,int d_nz,int *d_nnz,int o_nz,int *o_nnz)
8329c097c71SKris Buschelman {
8339c097c71SKris Buschelman   Mat       A;
83402217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
8359c097c71SKris Buschelman   int       ierr;
8369c097c71SKris Buschelman 
8379c097c71SKris Buschelman   PetscFunctionBegin;
8389c097c71SKris Buschelman   /*
8399c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
8409c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
8419c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
8429c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
8439c097c71SKris Buschelman     in the preallocation routine.
8449c097c71SKris Buschelman   */
8459c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
8469c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
8479c097c71SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
8489c097c71SKris Buschelman   PetscFunctionReturn(0);
8499c097c71SKris Buschelman }
8509c097c71SKris Buschelman EXTERN_C_END
8519c097c71SKris Buschelman 
8529c097c71SKris Buschelman EXTERN_C_BEGIN
8539c097c71SKris Buschelman #undef __FUNCT__
854f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
8558e9aea5cSBarry Smith int MatConvert_SBAIJ_SBAIJMUMPS(Mat A,const MatType newtype,Mat *newmat) {
856f579278aSKris Buschelman   int       ierr,size;
857f579278aSKris Buschelman   MPI_Comm  comm;
858f579278aSKris Buschelman   Mat       B=*newmat;
859422e82a1SHong Zhang   Mat_MUMPS *mumps;
8609c097c71SKris Buschelman   void      (*f)(void);
861f579278aSKris Buschelman 
862f579278aSKris Buschelman   PetscFunctionBegin;
863f579278aSKris Buschelman   if (B != A) {
864f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
865f579278aSKris Buschelman   }
866f579278aSKris Buschelman 
867f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
868f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
869f579278aSKris Buschelman 
870f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
871f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
872f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
873f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
874f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
875f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
876a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
877f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
878f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
879f579278aSKris Buschelman 
880f579278aSKris Buschelman   B->spptr                         = (void *)mumps;
881422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
882f0c56d0fSKris Buschelman   B->ops->view                     = MatView_AIJMUMPS;
883f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
884f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
8853924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
886f579278aSKris Buschelman 
887f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
888f579278aSKris Buschelman   if (size == 1) {
889f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
890f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
891f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
892f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
893f579278aSKris Buschelman   } else {
8949c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
8959c097c71SKris Buschelman     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);CHKERRQ(ierr);
8969c097c71SKris Buschelman     if (f) {
8979c097c71SKris Buschelman       mumps->MatPreallocate = (int (*)(Mat,int,int,int*,int,int*))f;
8989c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
8999c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
9009c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
9019c097c71SKris Buschelman     }
902f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
903f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
904f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
905f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
906f579278aSKris Buschelman   }
907f579278aSKris Buschelman 
908f579278aSKris Buschelman   PetscLogInfo(0,"Using MUMPS for Cholesky factorization and solves.");
909f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
910f579278aSKris Buschelman   *newmat = B;
911f579278aSKris Buschelman   PetscFunctionReturn(0);
912f579278aSKris Buschelman }
913f579278aSKris Buschelman EXTERN_C_END
914f579278aSKris Buschelman 
915f0c56d0fSKris Buschelman #undef __FUNCT__
916422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
917422e82a1SHong Zhang int MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
918f0c56d0fSKris Buschelman   int         ierr;
919*8e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
9208f340917SKris Buschelman 
921f0c56d0fSKris Buschelman   PetscFunctionBegin;
9228f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
923*8e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
924f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
925f0c56d0fSKris Buschelman }
926f0c56d0fSKris Buschelman 
92724b6179bSKris Buschelman /*MC
928fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
92924b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
93024b6179bSKris Buschelman 
93124b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
93224b6179bSKris Buschelman   on how to declare the existence of external packages),
93324b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
93424b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
93524b6179bSKris Buschelman   This matrix type is only supported for double precision real.
93624b6179bSKris Buschelman 
93724b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
93824b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
93924b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
94024b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
94128b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
94228b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
94328b08bd3SKris Buschelman   without data copy.
94424b6179bSKris Buschelman 
94524b6179bSKris Buschelman   Options Database Keys:
9460bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
94724b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
94824b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
94924b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
95024b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
95124b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
95224b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
95394b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
95424b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
95524b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
95624b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
95724b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
95824b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
95924b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
96024b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
96124b6179bSKris Buschelman 
96224b6179bSKris Buschelman   Level: beginner
96324b6179bSKris Buschelman 
96424b6179bSKris Buschelman .seealso: MATAIJMUMPS
96524b6179bSKris Buschelman M*/
96624b6179bSKris Buschelman 
967397b6df1SKris Buschelman EXTERN_C_BEGIN
968397b6df1SKris Buschelman #undef __FUNCT__
969f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
970f0c56d0fSKris Buschelman int MatCreate_SBAIJMUMPS(Mat A) {
971397b6df1SKris Buschelman   int ierr,size;
972397b6df1SKris Buschelman 
973397b6df1SKris Buschelman   PetscFunctionBegin;
9745441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
9755441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
9765441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
9775441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
978397b6df1SKris Buschelman   if (size == 1) {
979397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
980397b6df1SKris Buschelman   } else {
981397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
982397b6df1SKris Buschelman   }
983f0c56d0fSKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,&A);CHKERRQ(ierr);
984397b6df1SKris Buschelman   PetscFunctionReturn(0);
985397b6df1SKris Buschelman }
986397b6df1SKris Buschelman EXTERN_C_END
987