xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision 74a14a02d99cfad1d024587fb3d7711056ccd482)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
3397b6df1SKris Buschelman /*
4c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
5397b6df1SKris Buschelman */
6397b6df1SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
7397b6df1SKris Buschelman #include "src/mat/impls/aij/mpi/mpiaij.h"
8397b6df1SKris Buschelman #include "src/mat/impls/sbaij/seq/sbaij.h"
9397b6df1SKris Buschelman #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
10397b6df1SKris Buschelman 
11397b6df1SKris Buschelman EXTERN_C_BEGIN
12397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
13397b6df1SKris Buschelman #include "zmumps_c.h"
14397b6df1SKris Buschelman #else
15397b6df1SKris Buschelman #include "dmumps_c.h"
16397b6df1SKris Buschelman #endif
17397b6df1SKris Buschelman EXTERN_C_END
18397b6df1SKris Buschelman #define JOB_INIT -1
19397b6df1SKris Buschelman #define JOB_END -2
20397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
21397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
22397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
23397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
24a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
25397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
26adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
27397b6df1SKris Buschelman 
28397b6df1SKris Buschelman typedef struct {
29397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
30397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
31397b6df1SKris Buschelman #else
32397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
33397b6df1SKris Buschelman #endif
34397b6df1SKris Buschelman   MatStructure   matstruc;
35c1490034SHong Zhang   PetscMPIInt    myid,size;
36329ec9b3SHong Zhang   PetscInt       *irn,*jcn,sym,nSolve;
37397b6df1SKris Buschelman   PetscScalar    *val;
38397b6df1SKris Buschelman   MPI_Comm       comm_mumps;
39329ec9b3SHong Zhang   VecScatter     scat_rhs, scat_sol;
40c338a77dSKris Buschelman   PetscTruth     isAIJ,CleanUpMUMPS;
41329ec9b3SHong Zhang   Vec            b_seq,x_seq;
426849ba73SBarry Smith   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
436849ba73SBarry Smith   PetscErrorCode (*MatView)(Mat,PetscViewer);
446849ba73SBarry Smith   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
456849ba73SBarry Smith   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
466849ba73SBarry Smith   PetscErrorCode (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
476849ba73SBarry Smith   PetscErrorCode (*MatDestroy)(Mat);
486849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
496849ba73SBarry Smith   PetscErrorCode (*MatPreallocate)(Mat,int,int,int*,int,int*);
50f0c56d0fSKris Buschelman } Mat_MUMPS;
51f0c56d0fSKris Buschelman 
52dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
53892f6c3fSKris Buschelman EXTERN_C_BEGIN
5475179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat,MatType,MatReuse,Mat*);
55892f6c3fSKris Buschelman EXTERN_C_END
56397b6df1SKris Buschelman /* convert Petsc mpiaij matrix to triples: row[nz], col[nz], val[nz] */
57397b6df1SKris Buschelman /*
58397b6df1SKris Buschelman   input:
5975747be1SHong Zhang     A       - matrix in mpiaij or mpisbaij (bs=1) format
60397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
61397b6df1SKris Buschelman     valOnly - FALSE: spaces are allocated and values are set for the triple
62397b6df1SKris Buschelman               TRUE:  only the values in v array are updated
63397b6df1SKris Buschelman   output:
64397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
65397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
66397b6df1SKris Buschelman  */
67dfbe8321SBarry Smith PetscErrorCode MatConvertToTriples(Mat A,int shift,PetscTruth valOnly,int *nnz,int **r, int **c, PetscScalar **v) {
68c1490034SHong Zhang   PetscInt       *ai, *aj, *bi, *bj, rstart,nz, *garray;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
702a4c71feSBarry Smith   PetscInt       i,j,jj,jB,irow,m=A->rmap.n,*ajj,*bjj,countA,countB,colA_start,jcol;
71c1490034SHong Zhang   PetscInt       *row,*col;
72397b6df1SKris Buschelman   PetscScalar    *av, *bv,*val;
73f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
74397b6df1SKris Buschelman 
75397b6df1SKris Buschelman   PetscFunctionBegin;
76397b6df1SKris Buschelman   if (mumps->isAIJ){
77397b6df1SKris Buschelman     Mat_MPIAIJ    *mat =  (Mat_MPIAIJ*)A->data;
78397b6df1SKris Buschelman     Mat_SeqAIJ    *aa=(Mat_SeqAIJ*)(mat->A)->data;
79397b6df1SKris Buschelman     Mat_SeqAIJ    *bb=(Mat_SeqAIJ*)(mat->B)->data;
80397b6df1SKris Buschelman     nz = aa->nz + bb->nz;
812a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
82397b6df1SKris Buschelman     garray = mat->garray;
83397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
84397b6df1SKris Buschelman 
85397b6df1SKris Buschelman   } else {
86397b6df1SKris Buschelman     Mat_MPISBAIJ  *mat =  (Mat_MPISBAIJ*)A->data;
87397b6df1SKris Buschelman     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
88397b6df1SKris Buschelman     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
890c0e133fSBarry Smith     if (A->rmap.bs > 1) SETERRQ1(PETSC_ERR_SUP," bs=%d is not supported yet\n", A->rmap.bs);
906c6c5352SBarry Smith     nz = aa->nz + bb->nz;
912a4c71feSBarry Smith     ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap.rstart;
92397b6df1SKris Buschelman     garray = mat->garray;
93397b6df1SKris Buschelman     av=aa->a; bv=bb->a;
94397b6df1SKris Buschelman   }
95397b6df1SKris Buschelman 
96397b6df1SKris Buschelman   if (!valOnly){
977c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt) ,&row);CHKERRQ(ierr);
987c307921SBarry Smith     ierr = PetscMalloc(nz*sizeof(PetscInt),&col);CHKERRQ(ierr);
99397b6df1SKris Buschelman     ierr = PetscMalloc(nz*sizeof(PetscScalar),&val);CHKERRQ(ierr);
100397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
101397b6df1SKris Buschelman   } else {
102397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
103397b6df1SKris Buschelman   }
104397b6df1SKris Buschelman   *nnz = nz;
105397b6df1SKris Buschelman 
106028e57e8SHong Zhang   jj = 0; irow = rstart;
107397b6df1SKris Buschelman   for ( i=0; i<m; i++ ) {
108397b6df1SKris Buschelman     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
109397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
110397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
111397b6df1SKris Buschelman     bjj = bj + bi[i];
112397b6df1SKris Buschelman 
113397b6df1SKris Buschelman     /* get jB, the starting local col index for the 2nd B-part */
114397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
11575747be1SHong Zhang     j=-1;
11675747be1SHong Zhang     do {
11775747be1SHong Zhang       j++;
11875747be1SHong Zhang       if (j == countB) break;
119397b6df1SKris Buschelman       jcol = garray[bjj[j]];
12075747be1SHong Zhang     } while (jcol < colA_start);
12175747be1SHong Zhang     jB = j;
122397b6df1SKris Buschelman 
123397b6df1SKris Buschelman     /* B-part, smaller col index */
124397b6df1SKris Buschelman     colA_start = rstart + ajj[0]; /* the smallest col index for A */
125397b6df1SKris Buschelman     for (j=0; j<jB; j++){
126397b6df1SKris Buschelman       jcol = garray[bjj[j]];
127397b6df1SKris Buschelman       if (!valOnly){
128397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = jcol + shift;
12975747be1SHong Zhang 
130397b6df1SKris Buschelman       }
131397b6df1SKris Buschelman       val[jj++] = *bv++;
132397b6df1SKris Buschelman     }
133397b6df1SKris Buschelman     /* A-part */
134397b6df1SKris Buschelman     for (j=0; j<countA; j++){
135397b6df1SKris Buschelman       if (!valOnly){
136397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
137397b6df1SKris Buschelman       }
138397b6df1SKris Buschelman       val[jj++] = *av++;
139397b6df1SKris Buschelman     }
140397b6df1SKris Buschelman     /* B-part, larger col index */
141397b6df1SKris Buschelman     for (j=jB; j<countB; j++){
142397b6df1SKris Buschelman       if (!valOnly){
143397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
144397b6df1SKris Buschelman       }
145397b6df1SKris Buschelman       val[jj++] = *bv++;
146397b6df1SKris Buschelman     }
147397b6df1SKris Buschelman     irow++;
148397b6df1SKris Buschelman   }
149397b6df1SKris Buschelman 
150397b6df1SKris Buschelman   PetscFunctionReturn(0);
151397b6df1SKris Buschelman }
152397b6df1SKris Buschelman 
153c338a77dSKris Buschelman EXTERN_C_BEGIN
154c338a77dSKris Buschelman #undef __FUNCT__
155c338a77dSKris Buschelman #define __FUNCT__ "MatConvert_MUMPS_Base"
15675179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_MUMPS_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat) \
157521d7252SBarry Smith {
158dfbe8321SBarry Smith   PetscErrorCode ierr;
159c338a77dSKris Buschelman   Mat            B=*newmat;
160f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->spptr;
161901853e0SKris Buschelman   void           (*f)(void);
162c338a77dSKris Buschelman 
163c338a77dSKris Buschelman   PetscFunctionBegin;
164ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
165c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
166c338a77dSKris Buschelman   }
167f0c56d0fSKris Buschelman   B->ops->duplicate              = mumps->MatDuplicate;
168f0c56d0fSKris Buschelman   B->ops->view                   = mumps->MatView;
169f0c56d0fSKris Buschelman   B->ops->assemblyend            = mumps->MatAssemblyEnd;
170f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic       = mumps->MatLUFactorSymbolic;
171f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic = mumps->MatCholeskyFactorSymbolic;
172f0c56d0fSKris Buschelman   B->ops->destroy                = mumps->MatDestroy;
173901853e0SKris Buschelman 
174f68b968cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
175901853e0SKris Buschelman   if (f) {
176f68b968cSBarry Smith     ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)mumps->MatPreallocate);CHKERRQ(ierr);
177901853e0SKris Buschelman   }
178901853e0SKris Buschelman 
179901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
180901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_seqaij_C","",PETSC_NULL);CHKERRQ(ierr);
181901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_aijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
182901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_aijmumps_mpiaij_C","",PETSC_NULL);CHKERRQ(ierr);
1832895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
1842895b8caSSatish Balay   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr);
185901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C","",PETSC_NULL);CHKERRQ(ierr);
186901853e0SKris Buschelman   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C","",PETSC_NULL);CHKERRQ(ierr);
187901853e0SKris Buschelman 
188901853e0SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
189c338a77dSKris Buschelman   *newmat = B;
190c338a77dSKris Buschelman   PetscFunctionReturn(0);
191c338a77dSKris Buschelman }
192c338a77dSKris Buschelman EXTERN_C_END
193c338a77dSKris Buschelman 
194397b6df1SKris Buschelman #undef __FUNCT__
1953924e44cSKris Buschelman #define __FUNCT__ "MatDestroy_MUMPS"
196dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
197dfbe8321SBarry Smith {
198f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
199dfbe8321SBarry Smith   PetscErrorCode ierr;
200c1490034SHong Zhang   PetscMPIInt    size=lu->size;
2016849ba73SBarry Smith   PetscErrorCode (*specialdestroy)(Mat);
202397b6df1SKris Buschelman   PetscFunctionBegin;
203397b6df1SKris Buschelman   if (lu->CleanUpMUMPS) {
204397b6df1SKris Buschelman     /* Terminate instance, deallocate memories */
205329ec9b3SHong Zhang     if (size > 1){
206329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
207329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
208329ec9b3SHong Zhang       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
209329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
210329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
211329ec9b3SHong Zhang       ierr = PetscFree(lu->val);CHKERRQ(ierr);
212329ec9b3SHong Zhang     }
213397b6df1SKris Buschelman     lu->id.job=JOB_END;
214397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
215397b6df1SKris Buschelman     zmumps_c(&lu->id);
216397b6df1SKris Buschelman #else
217397b6df1SKris Buschelman     dmumps_c(&lu->id);
218397b6df1SKris Buschelman #endif
219c338a77dSKris Buschelman     ierr = PetscFree(lu->irn);CHKERRQ(ierr);
220c338a77dSKris Buschelman     ierr = PetscFree(lu->jcn);CHKERRQ(ierr);
221397b6df1SKris Buschelman     ierr = MPI_Comm_free(&(lu->comm_mumps));CHKERRQ(ierr);
222397b6df1SKris Buschelman   }
223a39386dcSKris Buschelman   specialdestroy = lu->specialdestroy;
224a39386dcSKris Buschelman   ierr = (*specialdestroy)(A);CHKERRQ(ierr);
225c338a77dSKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
226397b6df1SKris Buschelman   PetscFunctionReturn(0);
227397b6df1SKris Buschelman }
228397b6df1SKris Buschelman 
229397b6df1SKris Buschelman #undef __FUNCT__
230a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_AIJMUMPS"
231dfbe8321SBarry Smith PetscErrorCode MatDestroy_AIJMUMPS(Mat A)
232dfbe8321SBarry Smith {
2336849ba73SBarry Smith   PetscErrorCode ierr;
234329ec9b3SHong Zhang   PetscMPIInt    size;
235a39386dcSKris Buschelman 
236a39386dcSKris Buschelman   PetscFunctionBegin;
237a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
238a39386dcSKris Buschelman   if (size==1) {
239ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
240a39386dcSKris Buschelman   } else {
241ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPIAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
242a39386dcSKris Buschelman   }
243a39386dcSKris Buschelman   PetscFunctionReturn(0);
244a39386dcSKris Buschelman }
245a39386dcSKris Buschelman 
246a39386dcSKris Buschelman #undef __FUNCT__
247a39386dcSKris Buschelman #define __FUNCT__ "MatDestroy_SBAIJMUMPS"
248dfbe8321SBarry Smith PetscErrorCode MatDestroy_SBAIJMUMPS(Mat A)
249dfbe8321SBarry Smith {
2506849ba73SBarry Smith   PetscErrorCode ierr;
251329ec9b3SHong Zhang   PetscMPIInt    size;
252a39386dcSKris Buschelman 
253a39386dcSKris Buschelman   PetscFunctionBegin;
254a39386dcSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);
255a39386dcSKris Buschelman   if (size==1) {
256ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATSEQSBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
257a39386dcSKris Buschelman   } else {
258ceb03754SKris Buschelman     ierr = MatConvert_MUMPS_Base(A,MATMPISBAIJ,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
259a39386dcSKris Buschelman   }
260a39386dcSKris Buschelman   PetscFunctionReturn(0);
261a39386dcSKris Buschelman }
262a39386dcSKris Buschelman 
263a39386dcSKris Buschelman #undef __FUNCT__
264c338a77dSKris Buschelman #define __FUNCT__ "MatFactorInfo_MUMPS"
265dfbe8321SBarry Smith PetscErrorCode MatFactorInfo_MUMPS(Mat A,PetscViewer viewer) {
266f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
267dfbe8321SBarry Smith   PetscErrorCode ierr;
268397b6df1SKris Buschelman 
269397b6df1SKris Buschelman   PetscFunctionBegin;
270c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
271c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                  %d \n",lu->id.sym);CHKERRQ(ierr);
272c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):           %d \n",lu->id.par);CHKERRQ(ierr);
27368533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):        %d \n",lu->id.ICNTL(1));CHKERRQ(ierr);
27468533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg):%d \n",lu->id.ICNTL(2));CHKERRQ(ierr);
27568533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):  %d \n",lu->id.ICNTL(3));CHKERRQ(ierr);
276c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):       %d \n",lu->id.ICNTL(4));CHKERRQ(ierr);
277c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):        %d \n",lu->id.ICNTL(5));CHKERRQ(ierr);
278c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):       %d \n",lu->id.ICNTL(6));CHKERRQ(ierr);
279c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (matrix ordering):         %d \n",lu->id.ICNTL(7));CHKERRQ(ierr);
28068533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scalling strategy):       %d \n",lu->id.ICNTL(8));CHKERRQ(ierr);
281c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(9) (A/A^T x=b is solved):     %d \n",lu->id.ICNTL(9));CHKERRQ(ierr);
282c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements): %d \n",lu->id.ICNTL(10));CHKERRQ(ierr);
283c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):         %d \n",lu->id.ICNTL(11));CHKERRQ(ierr);
284958c9bccSBarry Smith   if (!lu->myid && lu->id.ICNTL(11)>0) {
285c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(4) (inf norm of input mat):        %g\n",lu->id.RINFOG(4));CHKERRQ(ierr);
286c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(5) (inf norm of solution):         %g\n",lu->id.RINFOG(5));CHKERRQ(ierr);
287c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(6) (inf norm of residual):         %g\n",lu->id.RINFOG(6));CHKERRQ(ierr);
288c338a77dSKris 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);
289c338a77dSKris Buschelman     ierr = PetscPrintf(PETSC_COMM_SELF,"        RINFOG(9) (error estimate):               %g \n",lu->id.RINFOG(9));CHKERRQ(ierr);
290c338a77dSKris 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);
291c338a77dSKris Buschelman 
292c338a77dSKris Buschelman   }
293c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",lu->id.ICNTL(12));CHKERRQ(ierr);
294c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",lu->id.ICNTL(13));CHKERRQ(ierr);
295adc1d99fSHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",lu->id.ICNTL(14));CHKERRQ(ierr);
29668533f03SHong Zhang   /* ICNTL(15-17) not used */
297c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",lu->id.ICNTL(18));CHKERRQ(ierr);
29868533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Shur complement info):                       %d \n",lu->id.ICNTL(19));CHKERRQ(ierr);
29968533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",lu->id.ICNTL(20));CHKERRQ(ierr);
30068533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",lu->id.ICNTL(21));CHKERRQ(ierr);
301c338a77dSKris Buschelman 
302c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",lu->id.CNTL(1));CHKERRQ(ierr);
303c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",lu->id.CNTL(2));CHKERRQ(ierr);
304c338a77dSKris Buschelman   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",lu->id.CNTL(3));CHKERRQ(ierr);
30568533f03SHong Zhang   ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",lu->id.CNTL(4));CHKERRQ(ierr);
30657f0c58bSHong Zhang 
30757f0c58bSHong Zhang   /* infomation local to each processor */
30868533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);}
30957f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %g \n",lu->myid,lu->id.RINFO(1));CHKERRQ(ierr);
31057f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
31168533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);}
31257f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(2));CHKERRQ(ierr);
31357f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
31468533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);}
31557f0c58bSHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d]  %g \n",lu->myid,lu->id.RINFO(3));CHKERRQ(ierr);
31657f0c58bSHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
31768533f03SHong Zhang   /*
31868533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(2) (info about error or warning ): \n");CHKERRQ(ierr);}
31968533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(2));CHKERRQ(ierr);
32068533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
32168533f03SHong Zhang   */
32268533f03SHong Zhang 
32368533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);}
32468533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(15));CHKERRQ(ierr);
32568533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
32668533f03SHong Zhang 
32768533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);}
32868533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(16));CHKERRQ(ierr);
32968533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
33068533f03SHong Zhang 
33168533f03SHong Zhang   if (!lu->myid) {ierr = PetscPrintf(PETSC_COMM_SELF, "      INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);}
33268533f03SHong Zhang   ierr = PetscSynchronizedPrintf(A->comm,"             [%d] %d \n",lu->myid,lu->id.INFO(23));CHKERRQ(ierr);
33368533f03SHong Zhang   ierr = PetscSynchronizedFlush(A->comm);
334adc1d99fSHong Zhang 
335958c9bccSBarry Smith   if (!lu->myid){ /* information from the host */
336adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",lu->id.RINFOG(1));CHKERRQ(ierr);
337adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",lu->id.RINFOG(2));CHKERRQ(ierr);
338adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",lu->id.RINFOG(3));CHKERRQ(ierr);
339adc1d99fSHong Zhang 
340adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(3));CHKERRQ(ierr);
341adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",lu->id.INFOG(4));CHKERRQ(ierr);
342adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",lu->id.INFOG(5));CHKERRQ(ierr);
343adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",lu->id.INFOG(6));CHKERRQ(ierr);
344adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively uese after analysis): %d \n",lu->id.INFOG(7));CHKERRQ(ierr);
345adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",lu->id.INFOG(8));CHKERRQ(ierr);
34668533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",lu->id.INFOG(9));CHKERRQ(ierr);
34768533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",lu->id.INFOG(10));CHKERRQ(ierr);
34868533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",lu->id.INFOG(11));CHKERRQ(ierr);
349adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",lu->id.INFOG(12));CHKERRQ(ierr);
350adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",lu->id.INFOG(13));CHKERRQ(ierr);
351adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",lu->id.INFOG(14));CHKERRQ(ierr);
352adc1d99fSHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",lu->id.INFOG(15));CHKERRQ(ierr);
35368533f03SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(16) (estimated size (in MB) of all MUMPS internal data for factorization after analysis: value on the most memory consuming processor): %d \n",lu->id.INFOG(16));CHKERRQ(ierr);
354adc1d99fSHong 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);
355adc1d99fSHong 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);
356adc1d99fSHong 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);
357adc1d99fSHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",lu->id.INFOG(20));CHKERRQ(ierr);
35868533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",lu->id.INFOG(21));CHKERRQ(ierr);
35968533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",lu->id.INFOG(22));CHKERRQ(ierr);
36068533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",lu->id.INFOG(23));CHKERRQ(ierr);
36168533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",lu->id.INFOG(24));CHKERRQ(ierr);
36268533f03SHong Zhang      ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",lu->id.INFOG(25));CHKERRQ(ierr);
363adc1d99fSHong Zhang   }
364adc1d99fSHong Zhang 
365397b6df1SKris Buschelman   PetscFunctionReturn(0);
366397b6df1SKris Buschelman }
367397b6df1SKris Buschelman 
368397b6df1SKris Buschelman #undef __FUNCT__
369c1490034SHong Zhang #define __FUNCT__ "MatView_MUMPS"
370c1490034SHong Zhang PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer) {
371dfbe8321SBarry Smith   PetscErrorCode    ierr;
37232077d6dSBarry Smith   PetscTruth        iascii;
373397b6df1SKris Buschelman   PetscViewerFormat format;
374f0c56d0fSKris Buschelman   Mat_MUMPS         *mumps=(Mat_MUMPS*)(A->spptr);
375397b6df1SKris Buschelman 
376397b6df1SKris Buschelman   PetscFunctionBegin;
377397b6df1SKris Buschelman   ierr = (*mumps->MatView)(A,viewer);CHKERRQ(ierr);
378397b6df1SKris Buschelman 
37932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
38032077d6dSBarry Smith   if (iascii) {
381397b6df1SKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
3822f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO){
383397b6df1SKris Buschelman       ierr = MatFactorInfo_MUMPS(A,viewer);CHKERRQ(ierr);
384397b6df1SKris Buschelman     }
385397b6df1SKris Buschelman   }
386397b6df1SKris Buschelman   PetscFunctionReturn(0);
387397b6df1SKris Buschelman }
388397b6df1SKris Buschelman 
389397b6df1SKris Buschelman #undef __FUNCT__
390f0c56d0fSKris Buschelman #define __FUNCT__ "MatSolve_AIJMUMPS"
391dfbe8321SBarry Smith PetscErrorCode MatSolve_AIJMUMPS(Mat A,Vec b,Vec x) {
392f0c56d0fSKris Buschelman   Mat_MUMPS      *lu=(Mat_MUMPS*)A->spptr;
393d54de34fSKris Buschelman   PetscScalar    *array;
394397b6df1SKris Buschelman   Vec            x_seq;
395329ec9b3SHong Zhang   IS             is_iden,is_petsc;
396329ec9b3SHong Zhang   VecScatter     scat_rhs=lu->scat_rhs,scat_sol=lu->scat_sol;
397dfbe8321SBarry Smith   PetscErrorCode ierr;
398329ec9b3SHong Zhang   PetscInt       i;
399397b6df1SKris Buschelman 
400397b6df1SKris Buschelman   PetscFunctionBegin;
401329ec9b3SHong Zhang   lu->id.nrhs = 1;
402329ec9b3SHong Zhang   x_seq = lu->b_seq;
403397b6df1SKris Buschelman   if (lu->size > 1){
404329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
405329ec9b3SHong Zhang     ierr = VecScatterBegin(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
406329ec9b3SHong Zhang     ierr = VecScatterEnd(b,x_seq,INSERT_VALUES,SCATTER_FORWARD,scat_rhs);CHKERRQ(ierr);
407397b6df1SKris Buschelman     if (!lu->myid) {ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);}
408397b6df1SKris Buschelman   } else {  /* size == 1 */
409397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
410397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
411397b6df1SKris Buschelman   }
412397b6df1SKris Buschelman   if (!lu->myid) { /* define rhs on the host */
413397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
414397b6df1SKris Buschelman     lu->id.rhs = (mumps_double_complex*)array;
415397b6df1SKris Buschelman #else
416397b6df1SKris Buschelman     lu->id.rhs = array;
417397b6df1SKris Buschelman #endif
418397b6df1SKris Buschelman   }
419329ec9b3SHong Zhang   if (lu->size == 1){
420329ec9b3SHong Zhang     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
421329ec9b3SHong Zhang   } else if (!lu->myid){
422329ec9b3SHong Zhang     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
423329ec9b3SHong Zhang   }
424329ec9b3SHong Zhang 
425329ec9b3SHong Zhang   if (lu->size > 1){
426329ec9b3SHong Zhang     /* distributed solution */
427329ec9b3SHong Zhang     lu->id.ICNTL(21) = 1;
428329ec9b3SHong Zhang     if (!lu->nSolve){
429329ec9b3SHong Zhang       /* Create x_seq=sol_loc for repeated use */
430329ec9b3SHong Zhang       PetscInt    lsol_loc;
431329ec9b3SHong Zhang       PetscScalar *sol_loc;
432329ec9b3SHong Zhang       lsol_loc = lu->id.INFO(23); /* length of sol_loc */
433329ec9b3SHong Zhang       ierr = PetscMalloc((1+lsol_loc)*(sizeof(PetscScalar)+sizeof(PetscInt)),&sol_loc);CHKERRQ(ierr);
434329ec9b3SHong Zhang       lu->id.isol_loc = (PetscInt *)(sol_loc + lsol_loc);
435329ec9b3SHong Zhang       lu->id.lsol_loc = lsol_loc;
436ea5584a4SSatish Balay       lu->id.sol_loc  = (F_DOUBLE *)sol_loc;
437329ec9b3SHong Zhang       ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,lsol_loc,sol_loc,&lu->x_seq);CHKERRQ(ierr);
438329ec9b3SHong Zhang     }
439329ec9b3SHong Zhang   }
440397b6df1SKris Buschelman 
441397b6df1SKris Buschelman   /* solve phase */
442329ec9b3SHong Zhang   /*-------------*/
443397b6df1SKris Buschelman   lu->id.job = 3;
444397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
445397b6df1SKris Buschelman   zmumps_c(&lu->id);
446397b6df1SKris Buschelman #else
447397b6df1SKris Buschelman   dmumps_c(&lu->id);
448397b6df1SKris Buschelman #endif
449397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
45079a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",lu->id.INFOG(1));
451397b6df1SKris Buschelman   }
452397b6df1SKris Buschelman 
453329ec9b3SHong Zhang   if (lu->size > 1) { /* convert mumps distributed solution to petsc mpi x */
454329ec9b3SHong Zhang     if (!lu->nSolve){ /* create scatter scat_sol */
455329ec9b3SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,lu->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
456329ec9b3SHong Zhang       for (i=0; i<lu->id.lsol_loc; i++){
457329ec9b3SHong Zhang         lu->id.isol_loc[i] -= 1; /* change Fortran style to C style */
458397b6df1SKris Buschelman       }
459329ec9b3SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,lu->id.lsol_loc,lu->id.isol_loc,&is_petsc);CHKERRQ(ierr);  /* to */
460329ec9b3SHong Zhang       ierr = VecScatterCreate(lu->x_seq,is_iden,x,is_petsc,&lu->scat_sol);CHKERRQ(ierr);
461329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
462329ec9b3SHong Zhang       ierr = ISDestroy(is_petsc);CHKERRQ(ierr);
463397b6df1SKris Buschelman     }
464329ec9b3SHong Zhang     ierr = VecScatterBegin(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
465329ec9b3SHong Zhang     ierr = VecScatterEnd(lu->x_seq,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat_sol);CHKERRQ(ierr);
466329ec9b3SHong Zhang   }
467329ec9b3SHong Zhang   lu->nSolve++;
468397b6df1SKris Buschelman   PetscFunctionReturn(0);
469397b6df1SKris Buschelman }
470397b6df1SKris Buschelman 
471a58c3f20SHong Zhang /*
472a58c3f20SHong Zhang   input:
473a58c3f20SHong Zhang    F:        numeric factor
474a58c3f20SHong Zhang   output:
475a58c3f20SHong Zhang    nneg:     total number of negative pivots
476a58c3f20SHong Zhang    nzero:    0
477a58c3f20SHong Zhang    npos:     (global dimension of F) - nneg
478a58c3f20SHong Zhang */
479a58c3f20SHong Zhang 
480a58c3f20SHong Zhang #undef __FUNCT__
481a58c3f20SHong Zhang #define __FUNCT__ "MatGetInertia_SBAIJMUMPS"
482dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
483a58c3f20SHong Zhang {
484a58c3f20SHong Zhang   Mat_MUMPS      *lu =(Mat_MUMPS*)F->spptr;
485dfbe8321SBarry Smith   PetscErrorCode ierr;
486c1490034SHong Zhang   PetscMPIInt    size;
487a58c3f20SHong Zhang 
488a58c3f20SHong Zhang   PetscFunctionBegin;
489bcb30aebSHong Zhang   ierr = MPI_Comm_size(F->comm,&size);CHKERRQ(ierr);
490bcb30aebSHong 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 */
491bcb30aebSHong Zhang   if (size > 1 && lu->id.ICNTL(13) != 1){
49279a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13));
493bcb30aebSHong Zhang   }
494a58c3f20SHong Zhang   if (nneg){
495a58c3f20SHong Zhang     if (!lu->myid){
496a58c3f20SHong Zhang       *nneg = lu->id.INFOG(12);
497a58c3f20SHong Zhang     }
498bcb30aebSHong Zhang     ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_mumps);CHKERRQ(ierr);
499a58c3f20SHong Zhang   }
500a58c3f20SHong Zhang   if (nzero) *nzero = 0;
5012a4c71feSBarry Smith   if (npos)  *npos  = F->rmap.N - (*nneg);
502a58c3f20SHong Zhang   PetscFunctionReturn(0);
503a58c3f20SHong Zhang }
504a58c3f20SHong Zhang 
505397b6df1SKris Buschelman #undef __FUNCT__
50619facb7aSBarry Smith #define __FUNCT__ "MatFactorNumeric_AIJMUMPS"
507af281ebdSHong Zhang PetscErrorCode MatFactorNumeric_AIJMUMPS(Mat A,MatFactorInfo *info,Mat *F)
508af281ebdSHong Zhang {
509f0c56d0fSKris Buschelman   Mat_MUMPS      *lu =(Mat_MUMPS*)(*F)->spptr;
510f0c56d0fSKris Buschelman   Mat_MUMPS      *lua=(Mat_MUMPS*)(A)->spptr;
5116849ba73SBarry Smith   PetscErrorCode ierr;
5122a4c71feSBarry Smith   PetscInt       rnz,nnz,nz=0,i,M=A->rmap.N,*ai,*aj,icntl;
513397b6df1SKris Buschelman   PetscTruth     valOnly,flg;
514e09efc27SHong Zhang   Mat            F_diag;
515397b6df1SKris Buschelman 
516397b6df1SKris Buschelman   PetscFunctionBegin;
517397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
518f0c56d0fSKris Buschelman     (*F)->ops->solve    = MatSolve_AIJMUMPS;
519397b6df1SKris Buschelman 
520397b6df1SKris Buschelman     /* Initialize a MUMPS instance */
521397b6df1SKris Buschelman     ierr = MPI_Comm_rank(A->comm, &lu->myid);
522397b6df1SKris Buschelman     ierr = MPI_Comm_size(A->comm,&lu->size);CHKERRQ(ierr);
52375747be1SHong Zhang     lua->myid = lu->myid; lua->size = lu->size;
524397b6df1SKris Buschelman     lu->id.job = JOB_INIT;
525397b6df1SKris Buschelman     ierr = MPI_Comm_dup(A->comm,&(lu->comm_mumps));CHKERRQ(ierr);
526a0e2756fSSatish Balay     ierr = MPICCommToFortranComm(lu->comm_mumps,&(lu->id.comm_fortran));CHKERRQ(ierr);
527397b6df1SKris Buschelman 
528397b6df1SKris Buschelman     /* Set mumps options */
529397b6df1SKris Buschelman     ierr = PetscOptionsBegin(A->comm,A->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
530397b6df1SKris Buschelman     lu->id.par=1;  /* host participates factorizaton and solve */
531397b6df1SKris Buschelman     lu->id.sym=lu->sym;
532397b6df1SKris Buschelman     if (lu->sym == 2){
533397b6df1SKris Buschelman       ierr = PetscOptionsInt("-mat_mumps_sym","SYM: (1,2)","None",lu->id.sym,&icntl,&flg);CHKERRQ(ierr);
534397b6df1SKris Buschelman       if (flg && icntl == 1) lu->id.sym=icntl;  /* matrix is spd */
535397b6df1SKris Buschelman     }
536397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
537397b6df1SKris Buschelman     zmumps_c(&lu->id);
538397b6df1SKris Buschelman #else
539397b6df1SKris Buschelman     dmumps_c(&lu->id);
540397b6df1SKris Buschelman #endif
541397b6df1SKris Buschelman 
542397b6df1SKris Buschelman     if (lu->size == 1){
543397b6df1SKris Buschelman       lu->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
544397b6df1SKris Buschelman     } else {
545397b6df1SKris Buschelman       lu->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
546397b6df1SKris Buschelman     }
547397b6df1SKris Buschelman 
548397b6df1SKris Buschelman     icntl=-1;
549397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",lu->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
55019facb7aSBarry Smith     if ((flg && icntl > 0) || PetscLogPrintInfo) {
551397b6df1SKris Buschelman       lu->id.ICNTL(4)=icntl; /* and use mumps default icntl(i), i=1,2,3 */
552397b6df1SKris Buschelman     } else { /* no output */
553397b6df1SKris Buschelman       lu->id.ICNTL(1) = 0;  /* error message, default= 6 */
554397b6df1SKris Buschelman       lu->id.ICNTL(2) = -1; /* output stream for diagnostic printing, statistics, and warning. default=0 */
555397b6df1SKris Buschelman       lu->id.ICNTL(3) = -1; /* output stream for global information, default=6 */
556397b6df1SKris Buschelman       lu->id.ICNTL(4) = 0;  /* level of printing, 0,1,2,3,4, default=2 */
557397b6df1SKris Buschelman     }
558397b6df1SKris 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);
559397b6df1SKris Buschelman     icntl=-1;
560397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): matrix ordering (0 to 7)","None",lu->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
561397b6df1SKris Buschelman     if (flg) {
562397b6df1SKris Buschelman       if (icntl== 1){
563397b6df1SKris Buschelman         SETERRQ(PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
564397b6df1SKris Buschelman       } else {
565397b6df1SKris Buschelman         lu->id.ICNTL(7) = icntl;
566397b6df1SKris Buschelman       }
567397b6df1SKris Buschelman     }
568397b6df1SKris 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);
569397b6df1SKris 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);
57094b7f48cSBarry 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);
571397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): efficiency control","None",lu->id.ICNTL(12),&lu->id.ICNTL(12),PETSC_NULL);CHKERRQ(ierr);
572397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): efficiency control","None",lu->id.ICNTL(13),&lu->id.ICNTL(13),PETSC_NULL);CHKERRQ(ierr);
573adc1d99fSHong 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);
574397b6df1SKris Buschelman     ierr = PetscOptionsInt("-mat_mumps_icntl_15","ICNTL(15): efficiency control","None",lu->id.ICNTL(15),&lu->id.ICNTL(15),PETSC_NULL);CHKERRQ(ierr);
575397b6df1SKris Buschelman 
576397b6df1SKris 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);
577397b6df1SKris 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);
578397b6df1SKris 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);
579397b6df1SKris Buschelman     PetscOptionsEnd();
580397b6df1SKris Buschelman   }
581397b6df1SKris Buschelman 
582397b6df1SKris Buschelman   /* define matrix A */
583397b6df1SKris Buschelman   switch (lu->id.ICNTL(18)){
584397b6df1SKris Buschelman   case 0:  /* centralized assembled matrix input (size=1) */
585397b6df1SKris Buschelman     if (!lu->myid) {
586c36ead0aSKris Buschelman       if (lua->isAIJ){
587397b6df1SKris Buschelman         Mat_SeqAIJ   *aa = (Mat_SeqAIJ*)A->data;
588397b6df1SKris Buschelman         nz               = aa->nz;
589397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
590397b6df1SKris Buschelman       } else {
591397b6df1SKris Buschelman         Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
5926c6c5352SBarry Smith         nz                  =  aa->nz;
593397b6df1SKris Buschelman         ai = aa->i; aj = aa->j; lu->val = aa->a;
594397b6df1SKris Buschelman       }
595397b6df1SKris Buschelman       if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization, get irn and jcn */
5967c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->irn);CHKERRQ(ierr);
5977c307921SBarry Smith         ierr = PetscMalloc(nz*sizeof(PetscInt),&lu->jcn);CHKERRQ(ierr);
598397b6df1SKris Buschelman         nz = 0;
599397b6df1SKris Buschelman         for (i=0; i<M; i++){
600397b6df1SKris Buschelman           rnz = ai[i+1] - ai[i];
601397b6df1SKris Buschelman           while (rnz--) {  /* Fortran row/col index! */
602397b6df1SKris Buschelman             lu->irn[nz] = i+1; lu->jcn[nz] = (*aj)+1; aj++; nz++;
603397b6df1SKris Buschelman           }
604397b6df1SKris Buschelman         }
605397b6df1SKris Buschelman       }
606397b6df1SKris Buschelman     }
607397b6df1SKris Buschelman     break;
608397b6df1SKris Buschelman   case 3:  /* distributed assembled matrix input (size>1) */
609397b6df1SKris Buschelman     if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
610397b6df1SKris Buschelman       valOnly = PETSC_FALSE;
611397b6df1SKris Buschelman     } else {
612397b6df1SKris Buschelman       valOnly = PETSC_TRUE; /* only update mat values, not row and col index */
613397b6df1SKris Buschelman     }
614397b6df1SKris Buschelman     ierr = MatConvertToTriples(A,1,valOnly, &nnz, &lu->irn, &lu->jcn, &lu->val);CHKERRQ(ierr);
615397b6df1SKris Buschelman     break;
616397b6df1SKris Buschelman   default: SETERRQ(PETSC_ERR_SUP,"Matrix input format is not supported by MUMPS.");
617397b6df1SKris Buschelman   }
618397b6df1SKris Buschelman 
619397b6df1SKris Buschelman   /* analysis phase */
620329ec9b3SHong Zhang   /*----------------*/
621397b6df1SKris Buschelman   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
622329ec9b3SHong Zhang     lu->id.job = 1;
623329ec9b3SHong Zhang 
624397b6df1SKris Buschelman     lu->id.n = M;
625397b6df1SKris Buschelman     switch (lu->id.ICNTL(18)){
626397b6df1SKris Buschelman     case 0:  /* centralized assembled matrix input */
627397b6df1SKris Buschelman       if (!lu->myid) {
628397b6df1SKris Buschelman         lu->id.nz =nz; lu->id.irn=lu->irn; lu->id.jcn=lu->jcn;
629397b6df1SKris Buschelman         if (lu->id.ICNTL(6)>1){
630397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
631397b6df1SKris Buschelman           lu->id.a = (mumps_double_complex*)lu->val;
632397b6df1SKris Buschelman #else
633397b6df1SKris Buschelman           lu->id.a = lu->val;
634397b6df1SKris Buschelman #endif
635397b6df1SKris Buschelman         }
636397b6df1SKris Buschelman       }
637397b6df1SKris Buschelman       break;
638397b6df1SKris Buschelman     case 3:  /* distributed assembled matrix input (size>1) */
639397b6df1SKris Buschelman       lu->id.nz_loc = nnz;
640397b6df1SKris Buschelman       lu->id.irn_loc=lu->irn; lu->id.jcn_loc=lu->jcn;
641397b6df1SKris Buschelman       if (lu->id.ICNTL(6)>1) {
642397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
643397b6df1SKris Buschelman         lu->id.a_loc = (mumps_double_complex*)lu->val;
644397b6df1SKris Buschelman #else
645397b6df1SKris Buschelman         lu->id.a_loc = lu->val;
646397b6df1SKris Buschelman #endif
647397b6df1SKris Buschelman       }
648329ec9b3SHong Zhang       /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
649329ec9b3SHong Zhang       IS  is_iden;
650329ec9b3SHong Zhang       Vec b;
651329ec9b3SHong Zhang       if (!lu->myid){
652329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap.N,&lu->b_seq);CHKERRQ(ierr);
653329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap.N,0,1,&is_iden);CHKERRQ(ierr);
654329ec9b3SHong Zhang       } else {
655329ec9b3SHong Zhang         ierr = VecCreateSeq(PETSC_COMM_SELF,0,&lu->b_seq);CHKERRQ(ierr);
656329ec9b3SHong Zhang         ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
657329ec9b3SHong Zhang       }
658329ec9b3SHong Zhang       ierr = VecCreate(A->comm,&b);CHKERRQ(ierr);
659329ec9b3SHong Zhang       ierr = VecSetSizes(b,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr);
660329ec9b3SHong Zhang       ierr = VecSetFromOptions(b);CHKERRQ(ierr);
661329ec9b3SHong Zhang 
662329ec9b3SHong Zhang       ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
663329ec9b3SHong Zhang       ierr = ISDestroy(is_iden);CHKERRQ(ierr);
664329ec9b3SHong Zhang       ierr = VecDestroy(b);CHKERRQ(ierr);
665397b6df1SKris Buschelman       break;
666397b6df1SKris Buschelman     }
667397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
668397b6df1SKris Buschelman     zmumps_c(&lu->id);
669397b6df1SKris Buschelman #else
670397b6df1SKris Buschelman     dmumps_c(&lu->id);
671397b6df1SKris Buschelman #endif
672397b6df1SKris Buschelman     if (lu->id.INFOG(1) < 0) {
67379a5c55eSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",lu->id.INFOG(1));
674397b6df1SKris Buschelman     }
675397b6df1SKris Buschelman   }
676397b6df1SKris Buschelman 
677397b6df1SKris Buschelman   /* numerical factorization phase */
678329ec9b3SHong Zhang   /*-------------------------------*/
679329ec9b3SHong Zhang   lu->id.job = 2;
680958c9bccSBarry Smith   if(!lu->id.ICNTL(18)) {
681a7aca84bSHong Zhang     if (!lu->myid) {
682397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
683397b6df1SKris Buschelman       lu->id.a = (mumps_double_complex*)lu->val;
684397b6df1SKris Buschelman #else
685397b6df1SKris Buschelman       lu->id.a = lu->val;
686397b6df1SKris Buschelman #endif
687397b6df1SKris Buschelman     }
688397b6df1SKris Buschelman   } else {
689397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
690397b6df1SKris Buschelman     lu->id.a_loc = (mumps_double_complex*)lu->val;
691397b6df1SKris Buschelman #else
692397b6df1SKris Buschelman     lu->id.a_loc = lu->val;
693397b6df1SKris Buschelman #endif
694397b6df1SKris Buschelman   }
695397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
696397b6df1SKris Buschelman   zmumps_c(&lu->id);
697397b6df1SKris Buschelman #else
698397b6df1SKris Buschelman   dmumps_c(&lu->id);
699397b6df1SKris Buschelman #endif
700397b6df1SKris Buschelman   if (lu->id.INFOG(1) < 0) {
70119facb7aSBarry Smith     if (lu->id.INFO(1) == -13) {
70219facb7aSBarry Smith       SETERRQ1(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: Cannot allocate required memory %d megabytes\n",lu->id.INFO(2));
70319facb7aSBarry Smith     } else {
70479a5c55eSBarry Smith       SETERRQ2(PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFO(1)=%d, INFO(2)=%d\n",lu->id.INFO(1),lu->id.INFO(2));
705397b6df1SKris Buschelman     }
70619facb7aSBarry Smith   }
707397b6df1SKris Buschelman 
70819facb7aSBarry Smith   if (!lu->myid && lu->id.ICNTL(16) > 0){
70979a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_LIB,"  lu->id.ICNTL(16):=%d\n",lu->id.INFOG(16));
710397b6df1SKris Buschelman   }
711397b6df1SKris Buschelman 
7128ada1bb4SHong Zhang   if (lu->size > 1){
713e09efc27SHong Zhang     if ((*F)->factor == FACTOR_LU){
714e09efc27SHong Zhang       F_diag = ((Mat_MPIAIJ *)(*F)->data)->A;
715e09efc27SHong Zhang     } else {
716e09efc27SHong Zhang       F_diag = ((Mat_MPISBAIJ *)(*F)->data)->A;
717e09efc27SHong Zhang     }
718e09efc27SHong Zhang     F_diag->assembled = PETSC_TRUE;
719329ec9b3SHong Zhang     if (lu->nSolve){
720329ec9b3SHong Zhang       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
721329ec9b3SHong Zhang       ierr = PetscFree(lu->id.sol_loc);CHKERRQ(ierr);
722329ec9b3SHong Zhang       ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);
723329ec9b3SHong Zhang     }
7248ada1bb4SHong Zhang   }
725397b6df1SKris Buschelman   (*F)->assembled   = PETSC_TRUE;
726397b6df1SKris Buschelman   lu->matstruc      = SAME_NONZERO_PATTERN;
727ace87b0dSHong Zhang   lu->CleanUpMUMPS  = PETSC_TRUE;
728329ec9b3SHong Zhang   lu->nSolve        = 0;
729397b6df1SKris Buschelman   PetscFunctionReturn(0);
730397b6df1SKris Buschelman }
731397b6df1SKris Buschelman 
732397b6df1SKris Buschelman /* Note the Petsc r and c permutations are ignored */
733397b6df1SKris Buschelman #undef __FUNCT__
734f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_AIJMUMPS"
735dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
736397b6df1SKris Buschelman   Mat            B;
737f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
738dfbe8321SBarry Smith   PetscErrorCode ierr;
739397b6df1SKris Buschelman 
740397b6df1SKris Buschelman   PetscFunctionBegin;
741397b6df1SKris Buschelman   /* Create the factorization matrix */
742f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
7432a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
744be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
745397b6df1SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
746397b6df1SKris Buschelman   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
747397b6df1SKris Buschelman 
748f0c56d0fSKris Buschelman   B->ops->lufactornumeric = MatFactorNumeric_AIJMUMPS;
749397b6df1SKris Buschelman   B->factor               = FACTOR_LU;
750f0c56d0fSKris Buschelman   lu                      = (Mat_MUMPS*)B->spptr;
751397b6df1SKris Buschelman   lu->sym                 = 0;
752397b6df1SKris Buschelman   lu->matstruc            = DIFFERENT_NONZERO_PATTERN;
753397b6df1SKris Buschelman 
754397b6df1SKris Buschelman   *F = B;
755397b6df1SKris Buschelman   PetscFunctionReturn(0);
756397b6df1SKris Buschelman }
757397b6df1SKris Buschelman 
758397b6df1SKris Buschelman /* Note the Petsc r permutation is ignored */
759397b6df1SKris Buschelman #undef __FUNCT__
760f0c56d0fSKris Buschelman #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJMUMPS"
761dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SBAIJMUMPS(Mat A,IS r,MatFactorInfo *info,Mat *F) {
762397b6df1SKris Buschelman   Mat            B;
763f0c56d0fSKris Buschelman   Mat_MUMPS      *lu;
764dfbe8321SBarry Smith   PetscErrorCode ierr;
765397b6df1SKris Buschelman 
766397b6df1SKris Buschelman   PetscFunctionBegin;
767397b6df1SKris Buschelman   /* Create the factorization matrix */
768f69a0ea3SMatthew Knepley   ierr = MatCreate(A->comm,&B);CHKERRQ(ierr);
7692a4c71feSBarry Smith   ierr = MatSetSizes(B,A->rmap.n,A->cmap.n,A->rmap.N,A->cmap.N);CHKERRQ(ierr);
770be5d1d56SKris Buschelman   ierr = MatSetType(B,A->type_name);CHKERRQ(ierr);
771efc670deSHong Zhang   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
772efc670deSHong Zhang   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
773397b6df1SKris Buschelman 
774f0c56d0fSKris Buschelman   B->ops->choleskyfactornumeric = MatFactorNumeric_AIJMUMPS;
775a58c3f20SHong Zhang   B->ops->getinertia            = MatGetInertia_SBAIJMUMPS;
776397b6df1SKris Buschelman   B->factor                     = FACTOR_CHOLESKY;
777f0c56d0fSKris Buschelman   lu                            = (Mat_MUMPS*)B->spptr;
778397b6df1SKris Buschelman   lu->sym                       = 2;
779397b6df1SKris Buschelman   lu->matstruc                  = DIFFERENT_NONZERO_PATTERN;
780397b6df1SKris Buschelman 
781397b6df1SKris Buschelman   *F = B;
782397b6df1SKris Buschelman   PetscFunctionReturn(0);
783397b6df1SKris Buschelman }
784397b6df1SKris Buschelman 
785397b6df1SKris Buschelman #undef __FUNCT__
786f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_AIJMUMPS"
787dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_AIJMUMPS(Mat A,MatAssemblyType mode) {
788dfbe8321SBarry Smith   PetscErrorCode ierr;
789f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
790c338a77dSKris Buschelman 
791397b6df1SKris Buschelman   PetscFunctionBegin;
792c338a77dSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
793f0c56d0fSKris Buschelman 
794c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
795c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
796f0c56d0fSKris Buschelman   A->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
797397b6df1SKris Buschelman   PetscFunctionReturn(0);
798397b6df1SKris Buschelman }
799397b6df1SKris Buschelman 
800c338a77dSKris Buschelman EXTERN_C_BEGIN
801c338a77dSKris Buschelman #undef __FUNCT__
802f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_AIJ_AIJMUMPS"
803*74a14a02SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_AIJ_AIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
804dfbe8321SBarry Smith {
805dfbe8321SBarry Smith   PetscErrorCode ierr;
806521d7252SBarry Smith   PetscMPIInt    size;
807c338a77dSKris Buschelman   MPI_Comm       comm;
808c338a77dSKris Buschelman   Mat            B=*newmat;
809f0c56d0fSKris Buschelman   Mat_MUMPS      *mumps;
810397b6df1SKris Buschelman 
811397b6df1SKris Buschelman   PetscFunctionBegin;
812ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
813c338a77dSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
814397b6df1SKris Buschelman   }
815397b6df1SKris Buschelman 
816c338a77dSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
817f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
818c338a77dSKris Buschelman 
819f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
820c338a77dSKris Buschelman   mumps->MatView                   = A->ops->view;
821c338a77dSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
822c338a77dSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
823c338a77dSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
824c338a77dSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
825a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_AIJMUMPS;
826c338a77dSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
827f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_TRUE;
828c338a77dSKris Buschelman 
8294b68dd72SKris Buschelman   B->spptr                         = (void*)mumps;
830422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
831c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
832f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_AIJMUMPS;
833f0c56d0fSKris Buschelman   B->ops->lufactorsymbolic         = MatLUFactorSymbolic_AIJMUMPS;
8343924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
835c338a77dSKris Buschelman 
836c338a77dSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
837c338a77dSKris Buschelman   if (size == 1) {
838c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_aijmumps_C",
839f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
840c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_seqaij_C",
841c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
842c338a77dSKris Buschelman   } else {
843c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpiaij_aijmumps_C",
844f0c56d0fSKris Buschelman                                              "MatConvert_AIJ_AIJMUMPS",MatConvert_AIJ_AIJMUMPS);CHKERRQ(ierr);
845c338a77dSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_aijmumps_mpiaij_C",
846c338a77dSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
847c338a77dSKris Buschelman   }
848c338a77dSKris Buschelman 
849ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for LU factorization and solves.\n");CHKERRQ(ierr);
850c338a77dSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
851c338a77dSKris Buschelman   *newmat = B;
852397b6df1SKris Buschelman   PetscFunctionReturn(0);
853397b6df1SKris Buschelman }
854c338a77dSKris Buschelman EXTERN_C_END
855397b6df1SKris Buschelman 
85624b6179bSKris Buschelman /*MC
857fafad747SKris Buschelman   MATAIJMUMPS - MATAIJMUMPS = "aijmumps" - A matrix type providing direct solvers (LU) for distributed
85824b6179bSKris Buschelman   and sequential matrices via the external package MUMPS.
85924b6179bSKris Buschelman 
86024b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
86124b6179bSKris Buschelman   on how to declare the existence of external packages),
86224b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
86324b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJMUMPS).
86424b6179bSKris Buschelman 
86524b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQAIJ.
86624b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPIAIJ.  Hence for single process communicators,
86724b6179bSKris Buschelman   MatSeqAIJSetPreallocation is supported, and similarly MatMPIAIJSetPreallocation is supported
86824b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
86928b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
87028b08bd3SKris Buschelman   conversion to or from the MATSEQAIJ or MATMPIAIJ type (depending on the communicator size)
87128b08bd3SKris Buschelman   without data copy.
87224b6179bSKris Buschelman 
87324b6179bSKris Buschelman   Options Database Keys:
8740bad9183SKris Buschelman + -mat_type aijmumps - sets the matrix type to "aijmumps" during a call to MatSetFromOptions()
87524b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
87624b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,1,2,3,4> - print level
87724b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
87824b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
87924b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
88024b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
88194b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
88224b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
88324b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
88424b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
88524b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
88624b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
88724b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
88824b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
88924b6179bSKris Buschelman 
89024b6179bSKris Buschelman   Level: beginner
89124b6179bSKris Buschelman 
89224b6179bSKris Buschelman .seealso: MATSBAIJMUMPS
89324b6179bSKris Buschelman M*/
89424b6179bSKris Buschelman 
895397b6df1SKris Buschelman EXTERN_C_BEGIN
896397b6df1SKris Buschelman #undef __FUNCT__
897f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_AIJMUMPS"
898be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJMUMPS(Mat A)
899dfbe8321SBarry Smith {
900dfbe8321SBarry Smith   PetscErrorCode ierr;
901c1490034SHong Zhang   PetscMPIInt    size;
902397b6df1SKris Buschelman   MPI_Comm       comm;
903397b6df1SKris Buschelman 
904397b6df1SKris Buschelman   PetscFunctionBegin;
9055441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqAIJ or MPIAIJ */
9065441df8eSKris Buschelman   /*   and AIJMUMPS types */
9075441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATAIJMUMPS);CHKERRQ(ierr);
908397b6df1SKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
909397b6df1SKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
910397b6df1SKris Buschelman   if (size == 1) {
911397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr);
912397b6df1SKris Buschelman   } else {
913397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
91400ff2a26SHong Zhang     /*
91500ff2a26SHong Zhang     Mat A_diag = ((Mat_MPIAIJ *)A->data)->A;
916ceb03754SKris Buschelman     ierr   = MatConvert_AIJ_AIJMUMPS(A_diag,MATAIJMUMPS,MAT_REUSE_MATRIX,&A_diag);CHKERRQ(ierr);
91700ff2a26SHong Zhang     */
918397b6df1SKris Buschelman   }
919ceb03754SKris Buschelman   ierr = MatConvert_AIJ_AIJMUMPS(A,MATAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
920397b6df1SKris Buschelman   PetscFunctionReturn(0);
921397b6df1SKris Buschelman }
922397b6df1SKris Buschelman EXTERN_C_END
923397b6df1SKris Buschelman 
924f579278aSKris Buschelman #undef __FUNCT__
925f0c56d0fSKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SBAIJMUMPS"
926dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_SBAIJMUMPS(Mat A,MatAssemblyType mode)
927dfbe8321SBarry Smith {
928dfbe8321SBarry Smith   PetscErrorCode ierr;
929f0c56d0fSKris Buschelman   Mat_MUMPS *mumps=(Mat_MUMPS*)A->spptr;
930f579278aSKris Buschelman 
931f579278aSKris Buschelman   PetscFunctionBegin;
932f579278aSKris Buschelman   ierr = (*mumps->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
933f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
934f0c56d0fSKris Buschelman   A->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
935f579278aSKris Buschelman   PetscFunctionReturn(0);
936f579278aSKris Buschelman }
937f579278aSKris Buschelman 
938f579278aSKris Buschelman EXTERN_C_BEGIN
939f579278aSKris Buschelman #undef __FUNCT__
9409c097c71SKris Buschelman #define __FUNCT__ "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS"
941c1490034SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation_MPISBAIJMUMPS(Mat B,PetscInt bs,PetscInt d_nz,PetscInt *d_nnz,PetscInt o_nz,PetscInt *o_nnz)
9429c097c71SKris Buschelman {
9439c097c71SKris Buschelman   Mat       A;
94402217bfdSHong Zhang   Mat_MUMPS *mumps=(Mat_MUMPS*)B->spptr;
945dfbe8321SBarry Smith   PetscErrorCode ierr;
9469c097c71SKris Buschelman 
9479c097c71SKris Buschelman   PetscFunctionBegin;
9489c097c71SKris Buschelman   /*
9499c097c71SKris Buschelman     After performing the MPISBAIJ Preallocation, we need to convert the local diagonal block matrix
9509c097c71SKris Buschelman     into MUMPS type so that the block jacobi preconditioner (for example) can use MUMPS.  I would
9519c097c71SKris Buschelman     like this to be done in the MatCreate routine, but the creation of this inner matrix requires
9529c097c71SKris Buschelman     block size info so that PETSc can determine the local size properly.  The block size info is set
9539c097c71SKris Buschelman     in the preallocation routine.
9549c097c71SKris Buschelman   */
9559c097c71SKris Buschelman   ierr = (*mumps->MatPreallocate)(B,bs,d_nz,d_nnz,o_nz,o_nnz);
9569c097c71SKris Buschelman   A    = ((Mat_MPISBAIJ *)B->data)->A;
957ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
9589c097c71SKris Buschelman   PetscFunctionReturn(0);
9599c097c71SKris Buschelman }
9609c097c71SKris Buschelman EXTERN_C_END
9619c097c71SKris Buschelman 
9629c097c71SKris Buschelman EXTERN_C_BEGIN
9639c097c71SKris Buschelman #undef __FUNCT__
964f0c56d0fSKris Buschelman #define __FUNCT__ "MatConvert_SBAIJ_SBAIJMUMPS"
96575179d2cSHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SBAIJ_SBAIJMUMPS(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
966dfbe8321SBarry Smith {
967dfbe8321SBarry Smith   PetscErrorCode ierr;
968521d7252SBarry Smith   PetscMPIInt    size;
969f579278aSKris Buschelman   MPI_Comm       comm;
970f579278aSKris Buschelman   Mat            B=*newmat;
971422e82a1SHong Zhang   Mat_MUMPS      *mumps;
9729c097c71SKris Buschelman   void           (*f)(void);
973f579278aSKris Buschelman 
974f579278aSKris Buschelman   PetscFunctionBegin;
975ceb03754SKris Buschelman   if (reuse == MAT_INITIAL_MATRIX) {
976f579278aSKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
977f579278aSKris Buschelman   }
978f579278aSKris Buschelman 
979f579278aSKris Buschelman   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
980f0c56d0fSKris Buschelman   ierr = PetscNew(Mat_MUMPS,&mumps);CHKERRQ(ierr);
981f579278aSKris Buschelman 
982f0c56d0fSKris Buschelman   mumps->MatDuplicate              = A->ops->duplicate;
983f579278aSKris Buschelman   mumps->MatView                   = A->ops->view;
984f579278aSKris Buschelman   mumps->MatAssemblyEnd            = A->ops->assemblyend;
985f579278aSKris Buschelman   mumps->MatLUFactorSymbolic       = A->ops->lufactorsymbolic;
986f579278aSKris Buschelman   mumps->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
987f579278aSKris Buschelman   mumps->MatDestroy                = A->ops->destroy;
988a39386dcSKris Buschelman   mumps->specialdestroy            = MatDestroy_SBAIJMUMPS;
989f579278aSKris Buschelman   mumps->CleanUpMUMPS              = PETSC_FALSE;
990f579278aSKris Buschelman   mumps->isAIJ                     = PETSC_FALSE;
991f579278aSKris Buschelman 
992f579278aSKris Buschelman   B->spptr                         = (void*)mumps;
993422e82a1SHong Zhang   B->ops->duplicate                = MatDuplicate_MUMPS;
994c1490034SHong Zhang   B->ops->view                     = MatView_MUMPS;
995f0c56d0fSKris Buschelman   B->ops->assemblyend              = MatAssemblyEnd_SBAIJMUMPS;
996f0c56d0fSKris Buschelman   B->ops->choleskyfactorsymbolic   = MatCholeskyFactorSymbolic_SBAIJMUMPS;
9973924e44cSKris Buschelman   B->ops->destroy                  = MatDestroy_MUMPS;
998f579278aSKris Buschelman 
999f579278aSKris Buschelman   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1000f579278aSKris Buschelman   if (size == 1) {
1001f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqsbaij_sbaijmumps_C",
1002f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1003f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_seqsbaij_C",
1004f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1005f579278aSKris Buschelman   } else {
10069c097c71SKris Buschelman   /* I really don't like needing to know the tag: MatMPISBAIJSetPreallocation_C */
1007f68b968cSBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",(PetscVoidStarFunction)&f);CHKERRQ(ierr);
1008901853e0SKris Buschelman     if (f) { /* This case should always be true when this routine is called */
10096849ba73SBarry Smith       mumps->MatPreallocate = (PetscErrorCode (*)(Mat,int,int,int*,int,int*))f;
10109c097c71SKris Buschelman       ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPISBAIJSetPreallocation_C",
10119c097c71SKris Buschelman                                                "MatMPISBAIJSetPreallocation_MPISBAIJMUMPS",
10129c097c71SKris Buschelman                                                MatMPISBAIJSetPreallocation_MPISBAIJMUMPS);CHKERRQ(ierr);
10139c097c71SKris Buschelman     }
1014f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpisbaij_sbaijmumps_C",
1015f0c56d0fSKris Buschelman                                              "MatConvert_SBAIJ_SBAIJMUMPS",MatConvert_SBAIJ_SBAIJMUMPS);CHKERRQ(ierr);
1016f0c56d0fSKris Buschelman     ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_sbaijmumps_mpisbaij_C",
1017f579278aSKris Buschelman                                              "MatConvert_MUMPS_Base",MatConvert_MUMPS_Base);CHKERRQ(ierr);
1018f579278aSKris Buschelman   }
1019f579278aSKris Buschelman 
1020ae15b995SBarry Smith   ierr = PetscInfo(0,"Using MUMPS for Cholesky factorization and solves.\n");CHKERRQ(ierr);
1021f579278aSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,newtype);CHKERRQ(ierr);
1022f579278aSKris Buschelman   *newmat = B;
1023f579278aSKris Buschelman   PetscFunctionReturn(0);
1024f579278aSKris Buschelman }
1025f579278aSKris Buschelman EXTERN_C_END
1026f579278aSKris Buschelman 
1027f0c56d0fSKris Buschelman #undef __FUNCT__
1028422e82a1SHong Zhang #define __FUNCT__ "MatDuplicate_MUMPS"
1029dfbe8321SBarry Smith PetscErrorCode MatDuplicate_MUMPS(Mat A, MatDuplicateOption op, Mat *M) {
1030dfbe8321SBarry Smith   PetscErrorCode ierr;
10318e393735SKris Buschelman   Mat_MUMPS   *lu=(Mat_MUMPS *)A->spptr;
10328f340917SKris Buschelman 
1033f0c56d0fSKris Buschelman   PetscFunctionBegin;
10348f340917SKris Buschelman   ierr = (*lu->MatDuplicate)(A,op,M);CHKERRQ(ierr);
10358e393735SKris Buschelman   ierr = PetscMemcpy((*M)->spptr,lu,sizeof(Mat_MUMPS));CHKERRQ(ierr);
1036f0c56d0fSKris Buschelman   PetscFunctionReturn(0);
1037f0c56d0fSKris Buschelman }
1038f0c56d0fSKris Buschelman 
103924b6179bSKris Buschelman /*MC
1040fafad747SKris Buschelman   MATSBAIJMUMPS - MATSBAIJMUMPS = "sbaijmumps" - A symmetric matrix type providing direct solvers (Cholesky) for
104124b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
104224b6179bSKris Buschelman 
104324b6179bSKris Buschelman   If MUMPS is installed (see the manual for instructions
104424b6179bSKris Buschelman   on how to declare the existence of external packages),
104524b6179bSKris Buschelman   a matrix type can be constructed which invokes MUMPS solvers.
104624b6179bSKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATSBAIJMUMPS).
104724b6179bSKris Buschelman 
104824b6179bSKris Buschelman   If created with a single process communicator, this matrix type inherits from MATSEQSBAIJ.
104924b6179bSKris Buschelman   Otherwise, this matrix type inherits from MATMPISBAIJ.  Hence for single process communicators,
105024b6179bSKris Buschelman   MatSeqSBAIJSetPreallocation is supported, and similarly MatMPISBAIJSetPreallocation is supported
105124b6179bSKris Buschelman   for communicators controlling multiple processes.  It is recommended that you call both of
105228b08bd3SKris Buschelman   the above preallocation routines for simplicity.  One can also call MatConvert for an inplace
105328b08bd3SKris Buschelman   conversion to or from the MATSEQSBAIJ or MATMPISBAIJ type (depending on the communicator size)
105428b08bd3SKris Buschelman   without data copy.
105524b6179bSKris Buschelman 
105624b6179bSKris Buschelman   Options Database Keys:
10570bad9183SKris Buschelman + -mat_type sbaijmumps - sets the matrix type to "sbaijmumps" during a call to MatSetFromOptions()
105824b6179bSKris Buschelman . -mat_mumps_sym <0,1,2> - 0 the matrix is unsymmetric, 1 symmetric positive definite, 2 symmetric
105924b6179bSKris Buschelman . -mat_mumps_icntl_4 <0,...,4> - print level
106024b6179bSKris Buschelman . -mat_mumps_icntl_6 <0,...,7> - matrix prescaling options (see MUMPS User's Guide)
106124b6179bSKris Buschelman . -mat_mumps_icntl_7 <0,...,7> - matrix orderings (see MUMPS User's Guide)
106224b6179bSKris Buschelman . -mat_mumps_icntl_9 <1,2> - A or A^T x=b to be solved: 1 denotes A, 2 denotes A^T
106324b6179bSKris Buschelman . -mat_mumps_icntl_10 <n> - maximum number of iterative refinements
106494b7f48cSBarry Smith . -mat_mumps_icntl_11 <n> - error analysis, a positive value returns statistics during -ksp_view
106524b6179bSKris Buschelman . -mat_mumps_icntl_12 <n> - efficiency control (see MUMPS User's Guide)
106624b6179bSKris Buschelman . -mat_mumps_icntl_13 <n> - efficiency control (see MUMPS User's Guide)
106724b6179bSKris Buschelman . -mat_mumps_icntl_14 <n> - efficiency control (see MUMPS User's Guide)
106824b6179bSKris Buschelman . -mat_mumps_icntl_15 <n> - efficiency control (see MUMPS User's Guide)
106924b6179bSKris Buschelman . -mat_mumps_cntl_1 <delta> - relative pivoting threshold
107024b6179bSKris Buschelman . -mat_mumps_cntl_2 <tol> - stopping criterion for refinement
107124b6179bSKris Buschelman - -mat_mumps_cntl_3 <adelta> - absolute pivoting threshold
107224b6179bSKris Buschelman 
107324b6179bSKris Buschelman   Level: beginner
107424b6179bSKris Buschelman 
107524b6179bSKris Buschelman .seealso: MATAIJMUMPS
107624b6179bSKris Buschelman M*/
107724b6179bSKris Buschelman 
1078397b6df1SKris Buschelman EXTERN_C_BEGIN
1079397b6df1SKris Buschelman #undef __FUNCT__
1080f0c56d0fSKris Buschelman #define __FUNCT__ "MatCreate_SBAIJMUMPS"
1081be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SBAIJMUMPS(Mat A)
1082dfbe8321SBarry Smith {
10836849ba73SBarry Smith   PetscErrorCode ierr;
10848ada1bb4SHong Zhang   PetscMPIInt    size;
1085397b6df1SKris Buschelman 
1086397b6df1SKris Buschelman   PetscFunctionBegin;
10875441df8eSKris Buschelman   /* Change type name before calling MatSetType to force proper construction of SeqSBAIJ or MPISBAIJ */
10885441df8eSKris Buschelman   /*   and SBAIJMUMPS types */
10895441df8eSKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)A,MATSBAIJMUMPS);CHKERRQ(ierr);
10905441df8eSKris Buschelman   ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr);CHKERRQ(ierr);
1091397b6df1SKris Buschelman   if (size == 1) {
1092397b6df1SKris Buschelman     ierr = MatSetType(A,MATSEQSBAIJ);CHKERRQ(ierr);
1093397b6df1SKris Buschelman   } else {
1094397b6df1SKris Buschelman     ierr   = MatSetType(A,MATMPISBAIJ);CHKERRQ(ierr);
1095397b6df1SKris Buschelman   }
1096ceb03754SKris Buschelman   ierr = MatConvert_SBAIJ_SBAIJMUMPS(A,MATSBAIJMUMPS,MAT_REUSE_MATRIX,&A);CHKERRQ(ierr);
1097397b6df1SKris Buschelman   PetscFunctionReturn(0);
1098397b6df1SKris Buschelman }
1099397b6df1SKris Buschelman EXTERN_C_END
1100