xref: /petsc/src/mat/impls/aij/mpi/mumps/mumps.c (revision c0be336448177ce2e241891bb6a25ca04d02feb8)
11c2a3de1SBarry Smith 
2397b6df1SKris Buschelman /*
3c2b5dc30SHong Zhang     Provides an interface to the MUMPS sparse solver
4397b6df1SKris Buschelman */
551d5961aSHong Zhang 
6c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
87ee00b23SStefano Zampini #include <../src/mat/impls/sell/mpi/mpisell.h>
9397b6df1SKris Buschelman 
10397b6df1SKris Buschelman EXTERN_C_BEGIN
11397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
122907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
132907cef9SHong Zhang #include <cmumps_c.h>
142907cef9SHong Zhang #else
15c6db04a5SJed Brown #include <zmumps_c.h>
162907cef9SHong Zhang #endif
172907cef9SHong Zhang #else
182907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
192907cef9SHong Zhang #include <smumps_c.h>
20397b6df1SKris Buschelman #else
21c6db04a5SJed Brown #include <dmumps_c.h>
22397b6df1SKris Buschelman #endif
232907cef9SHong Zhang #endif
24397b6df1SKris Buschelman EXTERN_C_END
25397b6df1SKris Buschelman #define JOB_INIT -1
263d472b54SHong Zhang #define JOB_FACTSYMBOLIC 1
273d472b54SHong Zhang #define JOB_FACTNUMERIC 2
283d472b54SHong Zhang #define JOB_SOLVE 3
29397b6df1SKris Buschelman #define JOB_END -2
303d472b54SHong Zhang 
312907cef9SHong Zhang /* calls to MUMPS */
322907cef9SHong Zhang #if defined(PETSC_USE_COMPLEX)
332907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
342907cef9SHong Zhang #define PetscMUMPS_c cmumps_c
352907cef9SHong Zhang #else
362907cef9SHong Zhang #define PetscMUMPS_c zmumps_c
372907cef9SHong Zhang #endif
382907cef9SHong Zhang #else
392907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
402907cef9SHong Zhang #define PetscMUMPS_c smumps_c
412907cef9SHong Zhang #else
422907cef9SHong Zhang #define PetscMUMPS_c dmumps_c
432907cef9SHong Zhang #endif
442907cef9SHong Zhang #endif
452907cef9SHong Zhang 
46940cd9d6SSatish Balay /* declare MumpsScalar */
47940cd9d6SSatish Balay #if defined(PETSC_USE_COMPLEX)
48940cd9d6SSatish Balay #if defined(PETSC_USE_REAL_SINGLE)
49940cd9d6SSatish Balay #define MumpsScalar mumps_complex
50940cd9d6SSatish Balay #else
51940cd9d6SSatish Balay #define MumpsScalar mumps_double_complex
52940cd9d6SSatish Balay #endif
53940cd9d6SSatish Balay #else
54940cd9d6SSatish Balay #define MumpsScalar PetscScalar
55940cd9d6SSatish Balay #endif
563d472b54SHong Zhang 
57397b6df1SKris Buschelman /* macros s.t. indices match MUMPS documentation */
58397b6df1SKris Buschelman #define ICNTL(I) icntl[(I)-1]
59397b6df1SKris Buschelman #define CNTL(I) cntl[(I)-1]
60397b6df1SKris Buschelman #define INFOG(I) infog[(I)-1]
61a7aca84bSHong Zhang #define INFO(I) info[(I)-1]
62397b6df1SKris Buschelman #define RINFOG(I) rinfog[(I)-1]
63adc1d99fSHong Zhang #define RINFO(I) rinfo[(I)-1]
64397b6df1SKris Buschelman 
65397b6df1SKris Buschelman typedef struct {
66397b6df1SKris Buschelman #if defined(PETSC_USE_COMPLEX)
672907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
682907cef9SHong Zhang   CMUMPS_STRUC_C id;
692907cef9SHong Zhang #else
70397b6df1SKris Buschelman   ZMUMPS_STRUC_C id;
712907cef9SHong Zhang #endif
722907cef9SHong Zhang #else
732907cef9SHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
742907cef9SHong Zhang   SMUMPS_STRUC_C id;
75397b6df1SKris Buschelman #else
76397b6df1SKris Buschelman   DMUMPS_STRUC_C id;
77397b6df1SKris Buschelman #endif
782907cef9SHong Zhang #endif
792907cef9SHong Zhang 
80397b6df1SKris Buschelman   MatStructure matstruc;
81c1490034SHong Zhang   PetscMPIInt  myid,size;
82a5e57a09SHong Zhang   PetscInt     *irn,*jcn,nz,sym;
83397b6df1SKris Buschelman   PetscScalar  *val;
84397b6df1SKris Buschelman   MPI_Comm     comm_mumps;
85a5e57a09SHong Zhang   PetscInt     ICNTL9_pre;           /* check if ICNTL(9) is changed from previous MatSolve */
86801fbe65SHong Zhang   VecScatter   scat_rhs, scat_sol;   /* used by MatSolve() */
87801fbe65SHong Zhang   Vec          b_seq,x_seq;
88b34f08ffSHong Zhang   PetscInt     ninfo,*info;          /* display INFO */
89b5fa320bSStefano Zampini   PetscInt     sizeredrhs;
9059ac8732SStefano Zampini   PetscScalar  *schur_sol;
9159ac8732SStefano Zampini   PetscInt     schur_sizesol;
922205254eSKarl Rupp 
93bccb9932SShri Abhyankar   PetscErrorCode (*ConvertToTriples)(Mat, int, MatReuse, int*, int**, int**, PetscScalar**);
94f0c56d0fSKris Buschelman } Mat_MUMPS;
95f0c56d0fSKris Buschelman 
9609573ac7SBarry Smith extern PetscErrorCode MatDuplicate_MUMPS(Mat,MatDuplicateOption,Mat*);
97b24902e0SBarry Smith 
9859ac8732SStefano Zampini static PetscErrorCode MatMumpsResetSchur_Private(Mat_MUMPS* mumps)
99b5fa320bSStefano Zampini {
100b5fa320bSStefano Zampini   PetscErrorCode ierr;
101b5fa320bSStefano Zampini 
102b5fa320bSStefano Zampini   PetscFunctionBegin;
10359ac8732SStefano Zampini   ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
10459ac8732SStefano Zampini   ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
10559ac8732SStefano Zampini   ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
10659ac8732SStefano Zampini   mumps->id.size_schur = 0;
107b3cb21ddSStefano Zampini   mumps->id.schur_lld  = 0;
10859ac8732SStefano Zampini   mumps->id.ICNTL(19)  = 0;
10959ac8732SStefano Zampini   PetscFunctionReturn(0);
11059ac8732SStefano Zampini }
11159ac8732SStefano Zampini 
112b3cb21ddSStefano Zampini /* solve with rhs in mumps->id.redrhs and return in the same location */
113b3cb21ddSStefano Zampini static PetscErrorCode MatMumpsSolveSchur_Private(Mat F)
11459ac8732SStefano Zampini {
115b3cb21ddSStefano Zampini   Mat_MUMPS            *mumps=(Mat_MUMPS*)F->data;
116b3cb21ddSStefano Zampini   Mat                  S,B,X;
117b3cb21ddSStefano Zampini   MatFactorSchurStatus schurstatus;
118b3cb21ddSStefano Zampini   PetscInt             sizesol;
11959ac8732SStefano Zampini   PetscErrorCode       ierr;
12059ac8732SStefano Zampini 
12159ac8732SStefano Zampini   PetscFunctionBegin;
122b3cb21ddSStefano Zampini   ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
123b3cb21ddSStefano Zampini   ierr = MatFactorGetSchurComplement(F,&S,&schurstatus);CHKERRQ(ierr);
124b3cb21ddSStefano Zampini   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&B);CHKERRQ(ierr);
125b3cb21ddSStefano Zampini   switch (schurstatus) {
126b3cb21ddSStefano Zampini   case MAT_FACTOR_SCHUR_FACTORED:
127b3cb21ddSStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,(PetscScalar*)mumps->id.redrhs,&X);CHKERRQ(ierr);
128b3cb21ddSStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
129b3cb21ddSStefano Zampini       ierr = MatMatSolveTranspose(S,B,X);CHKERRQ(ierr);
13059ac8732SStefano Zampini     } else {
131b3cb21ddSStefano Zampini       ierr = MatMatSolve(S,B,X);CHKERRQ(ierr);
13259ac8732SStefano Zampini     }
133b3cb21ddSStefano Zampini     break;
134b3cb21ddSStefano Zampini   case MAT_FACTOR_SCHUR_INVERTED:
135b3cb21ddSStefano Zampini     sizesol = mumps->id.nrhs*mumps->id.size_schur;
13659ac8732SStefano Zampini     if (!mumps->schur_sol || sizesol > mumps->schur_sizesol) {
13759ac8732SStefano Zampini       ierr = PetscFree(mumps->schur_sol);CHKERRQ(ierr);
13859ac8732SStefano Zampini       ierr = PetscMalloc1(sizesol,&mumps->schur_sol);CHKERRQ(ierr);
13959ac8732SStefano Zampini       mumps->schur_sizesol = sizesol;
140b5fa320bSStefano Zampini     }
141b3cb21ddSStefano Zampini     ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.nrhs,mumps->schur_sol,&X);CHKERRQ(ierr);
14259ac8732SStefano Zampini     if (!mumps->id.ICNTL(9)) { /* transpose solve */
143b3cb21ddSStefano Zampini       ierr = MatTransposeMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
144b5fa320bSStefano Zampini     } else {
145b3cb21ddSStefano Zampini       ierr = MatMatMult(S,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr);
146b5fa320bSStefano Zampini     }
147b3cb21ddSStefano Zampini     ierr = MatCopy(X,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
148b3cb21ddSStefano Zampini     break;
149b3cb21ddSStefano Zampini   default:
150b3cb21ddSStefano Zampini     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
151b3cb21ddSStefano Zampini     break;
15259ac8732SStefano Zampini   }
153b3cb21ddSStefano Zampini   ierr = MatFactorRestoreSchurComplement(F,&S,schurstatus);CHKERRQ(ierr);
154b3cb21ddSStefano Zampini   ierr = MatDestroy(&B);CHKERRQ(ierr);
155b3cb21ddSStefano Zampini   ierr = MatDestroy(&X);CHKERRQ(ierr);
156b5fa320bSStefano Zampini   PetscFunctionReturn(0);
157b5fa320bSStefano Zampini }
158b5fa320bSStefano Zampini 
159b3cb21ddSStefano Zampini static PetscErrorCode MatMumpsHandleSchur_Private(Mat F, PetscBool expansion)
160b5fa320bSStefano Zampini {
161b3cb21ddSStefano Zampini   Mat_MUMPS     *mumps=(Mat_MUMPS*)F->data;
162b5fa320bSStefano Zampini   PetscErrorCode ierr;
163b5fa320bSStefano Zampini 
164b5fa320bSStefano Zampini   PetscFunctionBegin;
165b5fa320bSStefano Zampini   if (!mumps->id.ICNTL(19)) { /* do nothing when Schur complement has not been computed */
166b5fa320bSStefano Zampini     PetscFunctionReturn(0);
167b5fa320bSStefano Zampini   }
168b8f61ee1SStefano Zampini   if (!expansion) { /* prepare for the condensation step */
169b5fa320bSStefano Zampini     PetscInt sizeredrhs = mumps->id.nrhs*mumps->id.size_schur;
170b5fa320bSStefano Zampini     /* allocate MUMPS internal array to store reduced right-hand sides */
171b5fa320bSStefano Zampini     if (!mumps->id.redrhs || sizeredrhs > mumps->sizeredrhs) {
172b5fa320bSStefano Zampini       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
173b5fa320bSStefano Zampini       mumps->id.lredrhs = mumps->id.size_schur;
174b5fa320bSStefano Zampini       ierr = PetscMalloc1(mumps->id.nrhs*mumps->id.lredrhs,&mumps->id.redrhs);CHKERRQ(ierr);
175b5fa320bSStefano Zampini       mumps->sizeredrhs = mumps->id.nrhs*mumps->id.lredrhs;
176b5fa320bSStefano Zampini     }
177b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 1; /* condensation phase */
178b5fa320bSStefano Zampini   } else { /* prepare for the expansion step */
179b8f61ee1SStefano Zampini     /* solve Schur complement (this has to be done by the MUMPS user, so basically us) */
180b3cb21ddSStefano Zampini     ierr = MatMumpsSolveSchur_Private(F);CHKERRQ(ierr);
181b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = 2; /* expansion phase */
182b5fa320bSStefano Zampini     PetscMUMPS_c(&mumps->id);
183b5fa320bSStefano Zampini     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
184b5fa320bSStefano Zampini     /* restore defaults */
185b5fa320bSStefano Zampini     mumps->id.ICNTL(26) = -1;
186d3d598ffSStefano Zampini     /* free MUMPS internal array for redrhs if we have solved for multiple rhs in order to save memory space */
187d3d598ffSStefano Zampini     if (mumps->id.nrhs > 1) {
188d3d598ffSStefano Zampini       ierr = PetscFree(mumps->id.redrhs);CHKERRQ(ierr);
189d3d598ffSStefano Zampini       mumps->id.lredrhs = 0;
190d3d598ffSStefano Zampini       mumps->sizeredrhs = 0;
191d3d598ffSStefano Zampini     }
192b5fa320bSStefano Zampini   }
193b5fa320bSStefano Zampini   PetscFunctionReturn(0);
194b5fa320bSStefano Zampini }
195b5fa320bSStefano Zampini 
196397b6df1SKris Buschelman /*
197d341cd04SHong Zhang   MatConvertToTriples_A_B - convert Petsc matrix to triples: row[nz], col[nz], val[nz]
198d341cd04SHong Zhang 
199397b6df1SKris Buschelman   input:
20067877ebaSShri Abhyankar     A       - matrix in aij,baij or sbaij (bs=1) format
201397b6df1SKris Buschelman     shift   - 0: C style output triple; 1: Fortran style output triple.
202bccb9932SShri Abhyankar     reuse   - MAT_INITIAL_MATRIX: spaces are allocated and values are set for the triple
203bccb9932SShri Abhyankar               MAT_REUSE_MATRIX:   only the values in v array are updated
204397b6df1SKris Buschelman   output:
205397b6df1SKris Buschelman     nnz     - dim of r, c, and v (number of local nonzero entries of A)
206397b6df1SKris Buschelman     r, c, v - row and col index, matrix values (matrix triples)
207eb9baa12SBarry Smith 
208eb9baa12SBarry Smith   The returned values r, c, and sometimes v are obtained in a single PetscMalloc(). Then in MatDestroy_MUMPS() it is
2097ee00b23SStefano Zampini   freed with PetscFree(mumps->irn);  This is not ideal code, the fact that v is ONLY sometimes part of mumps->irn means
210eb9baa12SBarry Smith   that the PetscMalloc() cannot easily be replaced with a PetscMalloc3().
211eb9baa12SBarry Smith 
212397b6df1SKris Buschelman  */
21316ebf90aSShri Abhyankar 
214bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
215b24902e0SBarry Smith {
216185f6596SHong Zhang   const PetscInt *ai,*aj,*ajj,M=A->rmap->n;
21767877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
218dfbe8321SBarry Smith   PetscErrorCode ierr;
219c1490034SHong Zhang   PetscInt       *row,*col;
22016ebf90aSShri Abhyankar   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
221397b6df1SKris Buschelman 
222397b6df1SKris Buschelman   PetscFunctionBegin;
22316ebf90aSShri Abhyankar   *v=aa->a;
224bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
2252205254eSKarl Rupp     nz   = aa->nz;
2262205254eSKarl Rupp     ai   = aa->i;
2272205254eSKarl Rupp     aj   = aa->j;
22816ebf90aSShri Abhyankar     *nnz = nz;
229785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
230185f6596SHong Zhang     col  = row + nz;
231185f6596SHong Zhang 
23216ebf90aSShri Abhyankar     nz = 0;
23316ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
23416ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
23567877ebaSShri Abhyankar       ajj = aj + ai[i];
23667877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
23767877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
23816ebf90aSShri Abhyankar       }
23916ebf90aSShri Abhyankar     }
24016ebf90aSShri Abhyankar     *r = row; *c = col;
24116ebf90aSShri Abhyankar   }
24216ebf90aSShri Abhyankar   PetscFunctionReturn(0);
24316ebf90aSShri Abhyankar }
244397b6df1SKris Buschelman 
2457ee00b23SStefano Zampini PetscErrorCode MatConvertToTriples_seqsell_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
2467ee00b23SStefano Zampini {
2477ee00b23SStefano Zampini   Mat_SeqSELL *a=(Mat_SeqSELL*)A->data;
2487ee00b23SStefano Zampini   PetscInt    *ptr;
2497ee00b23SStefano Zampini 
2507ee00b23SStefano Zampini   PetscFunctionBegin;
2517ee00b23SStefano Zampini   *v = a->val;
2527ee00b23SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
2537ee00b23SStefano Zampini     PetscInt       nz,i,j,row;
2547ee00b23SStefano Zampini     PetscErrorCode ierr;
2557ee00b23SStefano Zampini 
2567ee00b23SStefano Zampini     nz   = a->sliidx[a->totalslices];
2577ee00b23SStefano Zampini     *nnz = nz;
2587ee00b23SStefano Zampini     ierr = PetscMalloc1(2*nz, &ptr);CHKERRQ(ierr);
2597ee00b23SStefano Zampini     *r   = ptr;
2607ee00b23SStefano Zampini     *c   = ptr + nz;
2617ee00b23SStefano Zampini 
2627ee00b23SStefano Zampini     for (i=0; i<a->totalslices; i++) {
2637ee00b23SStefano Zampini       for (j=a->sliidx[i],row=0; j<a->sliidx[i+1]; j++,row=((row+1)&0x07)) {
2647ee00b23SStefano Zampini         *ptr++ = 8*i + row + shift;
2657ee00b23SStefano Zampini       }
2667ee00b23SStefano Zampini     }
2677ee00b23SStefano Zampini     for (i=0;i<nz;i++) *ptr++ = a->colidx[i] + shift;
2687ee00b23SStefano Zampini   }
2697ee00b23SStefano Zampini   PetscFunctionReturn(0);
2707ee00b23SStefano Zampini }
2717ee00b23SStefano Zampini 
272bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqbaij_seqaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
27367877ebaSShri Abhyankar {
27467877ebaSShri Abhyankar   Mat_SeqBAIJ    *aa=(Mat_SeqBAIJ*)A->data;
27533d57670SJed Brown   const PetscInt *ai,*aj,*ajj,bs2 = aa->bs2;
27633d57670SJed Brown   PetscInt       bs,M,nz,idx=0,rnz,i,j,k,m;
27767877ebaSShri Abhyankar   PetscErrorCode ierr;
27867877ebaSShri Abhyankar   PetscInt       *row,*col;
27967877ebaSShri Abhyankar 
28067877ebaSShri Abhyankar   PetscFunctionBegin;
28133d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
28233d57670SJed Brown   M = A->rmap->N/bs;
283cf3759fdSShri Abhyankar   *v = aa->a;
284bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
285cf3759fdSShri Abhyankar     ai   = aa->i; aj = aa->j;
28667877ebaSShri Abhyankar     nz   = bs2*aa->nz;
28767877ebaSShri Abhyankar     *nnz = nz;
288785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
289185f6596SHong Zhang     col  = row + nz;
290185f6596SHong Zhang 
29167877ebaSShri Abhyankar     for (i=0; i<M; i++) {
29267877ebaSShri Abhyankar       ajj = aj + ai[i];
29367877ebaSShri Abhyankar       rnz = ai[i+1] - ai[i];
29467877ebaSShri Abhyankar       for (k=0; k<rnz; k++) {
29567877ebaSShri Abhyankar         for (j=0; j<bs; j++) {
29667877ebaSShri Abhyankar           for (m=0; m<bs; m++) {
29767877ebaSShri Abhyankar             row[idx]   = i*bs + m + shift;
298cf3759fdSShri Abhyankar             col[idx++] = bs*(ajj[k]) + j + shift;
29967877ebaSShri Abhyankar           }
30067877ebaSShri Abhyankar         }
30167877ebaSShri Abhyankar       }
30267877ebaSShri Abhyankar     }
303cf3759fdSShri Abhyankar     *r = row; *c = col;
30467877ebaSShri Abhyankar   }
30567877ebaSShri Abhyankar   PetscFunctionReturn(0);
30667877ebaSShri Abhyankar }
30767877ebaSShri Abhyankar 
308bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqsbaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
30916ebf90aSShri Abhyankar {
31067877ebaSShri Abhyankar   const PetscInt *ai, *aj,*ajj,M=A->rmap->n;
31167877ebaSShri Abhyankar   PetscInt       nz,rnz,i,j;
31216ebf90aSShri Abhyankar   PetscErrorCode ierr;
31316ebf90aSShri Abhyankar   PetscInt       *row,*col;
31416ebf90aSShri Abhyankar   Mat_SeqSBAIJ   *aa=(Mat_SeqSBAIJ*)A->data;
31516ebf90aSShri Abhyankar 
31616ebf90aSShri Abhyankar   PetscFunctionBegin;
317882afa5aSHong Zhang   *v = aa->a;
318bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
3192205254eSKarl Rupp     nz   = aa->nz;
3202205254eSKarl Rupp     ai   = aa->i;
3212205254eSKarl Rupp     aj   = aa->j;
3222205254eSKarl Rupp     *v   = aa->a;
32316ebf90aSShri Abhyankar     *nnz = nz;
324785e854fSJed Brown     ierr = PetscMalloc1(2*nz, &row);CHKERRQ(ierr);
325185f6596SHong Zhang     col  = row + nz;
326185f6596SHong Zhang 
32716ebf90aSShri Abhyankar     nz = 0;
32816ebf90aSShri Abhyankar     for (i=0; i<M; i++) {
32916ebf90aSShri Abhyankar       rnz = ai[i+1] - ai[i];
33067877ebaSShri Abhyankar       ajj = aj + ai[i];
33167877ebaSShri Abhyankar       for (j=0; j<rnz; j++) {
33267877ebaSShri Abhyankar         row[nz] = i+shift; col[nz++] = ajj[j] + shift;
33316ebf90aSShri Abhyankar       }
33416ebf90aSShri Abhyankar     }
33516ebf90aSShri Abhyankar     *r = row; *c = col;
33616ebf90aSShri Abhyankar   }
33716ebf90aSShri Abhyankar   PetscFunctionReturn(0);
33816ebf90aSShri Abhyankar }
33916ebf90aSShri Abhyankar 
340bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_seqaij_seqsbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
34116ebf90aSShri Abhyankar {
34267877ebaSShri Abhyankar   const PetscInt    *ai,*aj,*ajj,*adiag,M=A->rmap->n;
34367877ebaSShri Abhyankar   PetscInt          nz,rnz,i,j;
34467877ebaSShri Abhyankar   const PetscScalar *av,*v1;
34516ebf90aSShri Abhyankar   PetscScalar       *val;
34616ebf90aSShri Abhyankar   PetscErrorCode    ierr;
34716ebf90aSShri Abhyankar   PetscInt          *row,*col;
348829b1710SHong Zhang   Mat_SeqAIJ        *aa=(Mat_SeqAIJ*)A->data;
34929b521d4Sstefano_zampini   PetscBool         missing;
35016ebf90aSShri Abhyankar 
35116ebf90aSShri Abhyankar   PetscFunctionBegin;
35216ebf90aSShri Abhyankar   ai    = aa->i; aj = aa->j; av = aa->a;
35316ebf90aSShri Abhyankar   adiag = aa->diag;
35429b521d4Sstefano_zampini   ierr  = MatMissingDiagonal_SeqAIJ(A,&missing,&i);CHKERRQ(ierr);
355bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
3567ee00b23SStefano Zampini     /* count nz in the upper triangular part of A */
357829b1710SHong Zhang     nz = 0;
35829b521d4Sstefano_zampini     if (missing) {
35929b521d4Sstefano_zampini       for (i=0; i<M; i++) {
36029b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
36129b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
36229b521d4Sstefano_zampini             if (aj[j] < i) continue;
36329b521d4Sstefano_zampini             nz++;
36429b521d4Sstefano_zampini           }
36529b521d4Sstefano_zampini         } else {
36629b521d4Sstefano_zampini           nz += ai[i+1] - adiag[i];
36729b521d4Sstefano_zampini         }
36829b521d4Sstefano_zampini       }
36929b521d4Sstefano_zampini     } else {
370829b1710SHong Zhang       for (i=0; i<M; i++) nz += ai[i+1] - adiag[i];
37129b521d4Sstefano_zampini     }
37216ebf90aSShri Abhyankar     *nnz = nz;
373829b1710SHong Zhang 
374185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
375185f6596SHong Zhang     col  = row + nz;
376185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
377185f6596SHong Zhang 
37816ebf90aSShri Abhyankar     nz = 0;
37929b521d4Sstefano_zampini     if (missing) {
38029b521d4Sstefano_zampini       for (i=0; i<M; i++) {
38129b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
38229b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
38329b521d4Sstefano_zampini             if (aj[j] < i) continue;
38429b521d4Sstefano_zampini             row[nz] = i+shift;
38529b521d4Sstefano_zampini             col[nz] = aj[j]+shift;
38629b521d4Sstefano_zampini             val[nz] = av[j];
38729b521d4Sstefano_zampini             nz++;
38829b521d4Sstefano_zampini           }
38929b521d4Sstefano_zampini         } else {
39029b521d4Sstefano_zampini           rnz = ai[i+1] - adiag[i];
39129b521d4Sstefano_zampini           ajj = aj + adiag[i];
39229b521d4Sstefano_zampini           v1  = av + adiag[i];
39329b521d4Sstefano_zampini           for (j=0; j<rnz; j++) {
39429b521d4Sstefano_zampini             row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
39529b521d4Sstefano_zampini           }
39629b521d4Sstefano_zampini         }
39729b521d4Sstefano_zampini       }
39829b521d4Sstefano_zampini     } else {
39916ebf90aSShri Abhyankar       for (i=0; i<M; i++) {
40016ebf90aSShri Abhyankar         rnz = ai[i+1] - adiag[i];
40167877ebaSShri Abhyankar         ajj = aj + adiag[i];
402cf3759fdSShri Abhyankar         v1  = av + adiag[i];
40367877ebaSShri Abhyankar         for (j=0; j<rnz; j++) {
40467877ebaSShri Abhyankar           row[nz] = i+shift; col[nz] = ajj[j] + shift; val[nz++] = v1[j];
40516ebf90aSShri Abhyankar         }
40616ebf90aSShri Abhyankar       }
40729b521d4Sstefano_zampini     }
40816ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
409397b6df1SKris Buschelman   } else {
41016ebf90aSShri Abhyankar     nz = 0; val = *v;
41129b521d4Sstefano_zampini     if (missing) {
41216ebf90aSShri Abhyankar       for (i=0; i <M; i++) {
41329b521d4Sstefano_zampini         if (PetscUnlikely(adiag[i] >= ai[i+1])) {
41429b521d4Sstefano_zampini           for (j=ai[i];j<ai[i+1];j++) {
41529b521d4Sstefano_zampini             if (aj[j] < i) continue;
41629b521d4Sstefano_zampini             val[nz++] = av[j];
41729b521d4Sstefano_zampini           }
41829b521d4Sstefano_zampini         } else {
41916ebf90aSShri Abhyankar           rnz = ai[i+1] - adiag[i];
42067877ebaSShri Abhyankar           v1  = av + adiag[i];
42167877ebaSShri Abhyankar           for (j=0; j<rnz; j++) {
42267877ebaSShri Abhyankar             val[nz++] = v1[j];
42316ebf90aSShri Abhyankar           }
42416ebf90aSShri Abhyankar         }
42516ebf90aSShri Abhyankar       }
42629b521d4Sstefano_zampini     } else {
42716ebf90aSShri Abhyankar       for (i=0; i <M; i++) {
42816ebf90aSShri Abhyankar         rnz = ai[i+1] - adiag[i];
42916ebf90aSShri Abhyankar         v1  = av + adiag[i];
43016ebf90aSShri Abhyankar         for (j=0; j<rnz; j++) {
43116ebf90aSShri Abhyankar           val[nz++] = v1[j];
43216ebf90aSShri Abhyankar         }
43316ebf90aSShri Abhyankar       }
43416ebf90aSShri Abhyankar     }
43529b521d4Sstefano_zampini   }
43616ebf90aSShri Abhyankar   PetscFunctionReturn(0);
43716ebf90aSShri Abhyankar }
43816ebf90aSShri Abhyankar 
439bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
44016ebf90aSShri Abhyankar {
44116ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
44216ebf90aSShri Abhyankar   PetscErrorCode    ierr;
44316ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
44416ebf90aSShri Abhyankar   PetscInt          *row,*col;
44516ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
44616ebf90aSShri Abhyankar   PetscScalar       *val;
447397b6df1SKris Buschelman   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ*)A->data;
448397b6df1SKris Buschelman   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ*)(mat->A)->data;
449397b6df1SKris Buschelman   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ*)(mat->B)->data;
45016ebf90aSShri Abhyankar 
45116ebf90aSShri Abhyankar   PetscFunctionBegin;
452d0f46423SBarry Smith   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
453397b6df1SKris Buschelman   av=aa->a; bv=bb->a;
454397b6df1SKris Buschelman 
4552205254eSKarl Rupp   garray = mat->garray;
4562205254eSKarl Rupp 
457bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
45816ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
45916ebf90aSShri Abhyankar     *nnz = nz;
460185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
461185f6596SHong Zhang     col  = row + nz;
462185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
463185f6596SHong Zhang 
464397b6df1SKris Buschelman     *r = row; *c = col; *v = val;
465397b6df1SKris Buschelman   } else {
466397b6df1SKris Buschelman     row = *r; col = *c; val = *v;
467397b6df1SKris Buschelman   }
468397b6df1SKris Buschelman 
469028e57e8SHong Zhang   jj = 0; irow = rstart;
470397b6df1SKris Buschelman   for (i=0; i<m; i++) {
471397b6df1SKris Buschelman     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
472397b6df1SKris Buschelman     countA = ai[i+1] - ai[i];
473397b6df1SKris Buschelman     countB = bi[i+1] - bi[i];
474397b6df1SKris Buschelman     bjj    = bj + bi[i];
47516ebf90aSShri Abhyankar     v1     = av + ai[i];
47616ebf90aSShri Abhyankar     v2     = bv + bi[i];
477397b6df1SKris Buschelman 
478397b6df1SKris Buschelman     /* A-part */
479397b6df1SKris Buschelman     for (j=0; j<countA; j++) {
480bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
481397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
482397b6df1SKris Buschelman       }
48316ebf90aSShri Abhyankar       val[jj++] = v1[j];
484397b6df1SKris Buschelman     }
48516ebf90aSShri Abhyankar 
48616ebf90aSShri Abhyankar     /* B-part */
48716ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
488bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
489397b6df1SKris Buschelman         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
490397b6df1SKris Buschelman       }
49116ebf90aSShri Abhyankar       val[jj++] = v2[j];
49216ebf90aSShri Abhyankar     }
49316ebf90aSShri Abhyankar     irow++;
49416ebf90aSShri Abhyankar   }
49516ebf90aSShri Abhyankar   PetscFunctionReturn(0);
49616ebf90aSShri Abhyankar }
49716ebf90aSShri Abhyankar 
498bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
49916ebf90aSShri Abhyankar {
50016ebf90aSShri Abhyankar   const PetscInt    *ai, *aj, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
50116ebf90aSShri Abhyankar   PetscErrorCode    ierr;
50216ebf90aSShri Abhyankar   PetscInt          rstart,nz,i,j,jj,irow,countA,countB;
50316ebf90aSShri Abhyankar   PetscInt          *row,*col;
50416ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
50516ebf90aSShri Abhyankar   PetscScalar       *val;
50616ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat = (Mat_MPIAIJ*)A->data;
50716ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ*)(mat->A)->data;
50816ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ*)(mat->B)->data;
50916ebf90aSShri Abhyankar 
51016ebf90aSShri Abhyankar   PetscFunctionBegin;
51116ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; bi=bb->i; bj=bb->j; rstart= A->rmap->rstart;
51216ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
51316ebf90aSShri Abhyankar 
5142205254eSKarl Rupp   garray = mat->garray;
5152205254eSKarl Rupp 
516bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
51716ebf90aSShri Abhyankar     nz   = aa->nz + bb->nz;
51816ebf90aSShri Abhyankar     *nnz = nz;
519185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
520185f6596SHong Zhang     col  = row + nz;
521185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
522185f6596SHong Zhang 
52316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
52416ebf90aSShri Abhyankar   } else {
52516ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
52616ebf90aSShri Abhyankar   }
52716ebf90aSShri Abhyankar 
52816ebf90aSShri Abhyankar   jj = 0; irow = rstart;
52916ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
53016ebf90aSShri Abhyankar     ajj    = aj + ai[i];                 /* ptr to the beginning of this row */
53116ebf90aSShri Abhyankar     countA = ai[i+1] - ai[i];
53216ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
53316ebf90aSShri Abhyankar     bjj    = bj + bi[i];
53416ebf90aSShri Abhyankar     v1     = av + ai[i];
53516ebf90aSShri Abhyankar     v2     = bv + bi[i];
53616ebf90aSShri Abhyankar 
53716ebf90aSShri Abhyankar     /* A-part */
53816ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
539bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
54016ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
54116ebf90aSShri Abhyankar       }
54216ebf90aSShri Abhyankar       val[jj++] = v1[j];
54316ebf90aSShri Abhyankar     }
54416ebf90aSShri Abhyankar 
54516ebf90aSShri Abhyankar     /* B-part */
54616ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
547bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
54816ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
54916ebf90aSShri Abhyankar       }
55016ebf90aSShri Abhyankar       val[jj++] = v2[j];
55116ebf90aSShri Abhyankar     }
55216ebf90aSShri Abhyankar     irow++;
55316ebf90aSShri Abhyankar   }
55416ebf90aSShri Abhyankar   PetscFunctionReturn(0);
55516ebf90aSShri Abhyankar }
55616ebf90aSShri Abhyankar 
557bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpibaij_mpiaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
55867877ebaSShri Abhyankar {
55967877ebaSShri Abhyankar   Mat_MPIBAIJ       *mat    = (Mat_MPIBAIJ*)A->data;
56067877ebaSShri Abhyankar   Mat_SeqBAIJ       *aa     = (Mat_SeqBAIJ*)(mat->A)->data;
56167877ebaSShri Abhyankar   Mat_SeqBAIJ       *bb     = (Mat_SeqBAIJ*)(mat->B)->data;
56267877ebaSShri Abhyankar   const PetscInt    *ai     = aa->i, *bi = bb->i, *aj = aa->j, *bj = bb->j,*ajj, *bjj;
563d985c460SShri Abhyankar   const PetscInt    *garray = mat->garray,mbs=mat->mbs,rstart=A->rmap->rstart;
56433d57670SJed Brown   const PetscInt    bs2=mat->bs2;
56567877ebaSShri Abhyankar   PetscErrorCode    ierr;
56633d57670SJed Brown   PetscInt          bs,nz,i,j,k,n,jj,irow,countA,countB,idx;
56767877ebaSShri Abhyankar   PetscInt          *row,*col;
56867877ebaSShri Abhyankar   const PetscScalar *av=aa->a, *bv=bb->a,*v1,*v2;
56967877ebaSShri Abhyankar   PetscScalar       *val;
57067877ebaSShri Abhyankar 
57167877ebaSShri Abhyankar   PetscFunctionBegin;
57233d57670SJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
573bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
57467877ebaSShri Abhyankar     nz   = bs2*(aa->nz + bb->nz);
57567877ebaSShri Abhyankar     *nnz = nz;
576185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
577185f6596SHong Zhang     col  = row + nz;
578185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
579185f6596SHong Zhang 
58067877ebaSShri Abhyankar     *r = row; *c = col; *v = val;
58167877ebaSShri Abhyankar   } else {
58267877ebaSShri Abhyankar     row = *r; col = *c; val = *v;
58367877ebaSShri Abhyankar   }
58467877ebaSShri Abhyankar 
585d985c460SShri Abhyankar   jj = 0; irow = rstart;
58667877ebaSShri Abhyankar   for (i=0; i<mbs; i++) {
58767877ebaSShri Abhyankar     countA = ai[i+1] - ai[i];
58867877ebaSShri Abhyankar     countB = bi[i+1] - bi[i];
58967877ebaSShri Abhyankar     ajj    = aj + ai[i];
59067877ebaSShri Abhyankar     bjj    = bj + bi[i];
59167877ebaSShri Abhyankar     v1     = av + bs2*ai[i];
59267877ebaSShri Abhyankar     v2     = bv + bs2*bi[i];
59367877ebaSShri Abhyankar 
59467877ebaSShri Abhyankar     idx = 0;
59567877ebaSShri Abhyankar     /* A-part */
59667877ebaSShri Abhyankar     for (k=0; k<countA; k++) {
59767877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
59867877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
599bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
600d985c460SShri Abhyankar             row[jj] = irow + n + shift;
601d985c460SShri Abhyankar             col[jj] = rstart + bs*ajj[k] + j + shift;
60267877ebaSShri Abhyankar           }
60367877ebaSShri Abhyankar           val[jj++] = v1[idx++];
60467877ebaSShri Abhyankar         }
60567877ebaSShri Abhyankar       }
60667877ebaSShri Abhyankar     }
60767877ebaSShri Abhyankar 
60867877ebaSShri Abhyankar     idx = 0;
60967877ebaSShri Abhyankar     /* B-part */
61067877ebaSShri Abhyankar     for (k=0; k<countB; k++) {
61167877ebaSShri Abhyankar       for (j=0; j<bs; j++) {
61267877ebaSShri Abhyankar         for (n=0; n<bs; n++) {
613bccb9932SShri Abhyankar           if (reuse == MAT_INITIAL_MATRIX) {
614d985c460SShri Abhyankar             row[jj] = irow + n + shift;
615d985c460SShri Abhyankar             col[jj] = bs*garray[bjj[k]] + j + shift;
61667877ebaSShri Abhyankar           }
617d985c460SShri Abhyankar           val[jj++] = v2[idx++];
61867877ebaSShri Abhyankar         }
61967877ebaSShri Abhyankar       }
62067877ebaSShri Abhyankar     }
621d985c460SShri Abhyankar     irow += bs;
62267877ebaSShri Abhyankar   }
62367877ebaSShri Abhyankar   PetscFunctionReturn(0);
62467877ebaSShri Abhyankar }
62567877ebaSShri Abhyankar 
626bccb9932SShri Abhyankar PetscErrorCode MatConvertToTriples_mpiaij_mpisbaij(Mat A,int shift,MatReuse reuse,int *nnz,int **r, int **c, PetscScalar **v)
62716ebf90aSShri Abhyankar {
62816ebf90aSShri Abhyankar   const PetscInt    *ai, *aj,*adiag, *bi, *bj,*garray,m=A->rmap->n,*ajj,*bjj;
62916ebf90aSShri Abhyankar   PetscErrorCode    ierr;
630e0bace9bSHong Zhang   PetscInt          rstart,nz,nza,nzb,i,j,jj,irow,countA,countB;
63116ebf90aSShri Abhyankar   PetscInt          *row,*col;
63216ebf90aSShri Abhyankar   const PetscScalar *av, *bv,*v1,*v2;
63316ebf90aSShri Abhyankar   PetscScalar       *val;
63416ebf90aSShri Abhyankar   Mat_MPIAIJ        *mat =  (Mat_MPIAIJ*)A->data;
63516ebf90aSShri Abhyankar   Mat_SeqAIJ        *aa  =(Mat_SeqAIJ*)(mat->A)->data;
63616ebf90aSShri Abhyankar   Mat_SeqAIJ        *bb  =(Mat_SeqAIJ*)(mat->B)->data;
63716ebf90aSShri Abhyankar 
63816ebf90aSShri Abhyankar   PetscFunctionBegin;
63916ebf90aSShri Abhyankar   ai=aa->i; aj=aa->j; adiag=aa->diag;
64016ebf90aSShri Abhyankar   bi=bb->i; bj=bb->j; garray = mat->garray;
64116ebf90aSShri Abhyankar   av=aa->a; bv=bb->a;
6422205254eSKarl Rupp 
64316ebf90aSShri Abhyankar   rstart = A->rmap->rstart;
64416ebf90aSShri Abhyankar 
645bccb9932SShri Abhyankar   if (reuse == MAT_INITIAL_MATRIX) {
646e0bace9bSHong Zhang     nza = 0;    /* num of upper triangular entries in mat->A, including diagonals */
647e0bace9bSHong Zhang     nzb = 0;    /* num of upper triangular entries in mat->B */
64816ebf90aSShri Abhyankar     for (i=0; i<m; i++) {
649e0bace9bSHong Zhang       nza   += (ai[i+1] - adiag[i]);
65016ebf90aSShri Abhyankar       countB = bi[i+1] - bi[i];
65116ebf90aSShri Abhyankar       bjj    = bj + bi[i];
652e0bace9bSHong Zhang       for (j=0; j<countB; j++) {
653e0bace9bSHong Zhang         if (garray[bjj[j]] > rstart) nzb++;
654e0bace9bSHong Zhang       }
655e0bace9bSHong Zhang     }
65616ebf90aSShri Abhyankar 
657e0bace9bSHong Zhang     nz   = nza + nzb; /* total nz of upper triangular part of mat */
65816ebf90aSShri Abhyankar     *nnz = nz;
659185f6596SHong Zhang     ierr = PetscMalloc((2*nz*sizeof(PetscInt)+nz*sizeof(PetscScalar)), &row);CHKERRQ(ierr);
660185f6596SHong Zhang     col  = row + nz;
661185f6596SHong Zhang     val  = (PetscScalar*)(col + nz);
662185f6596SHong Zhang 
66316ebf90aSShri Abhyankar     *r = row; *c = col; *v = val;
66416ebf90aSShri Abhyankar   } else {
66516ebf90aSShri Abhyankar     row = *r; col = *c; val = *v;
66616ebf90aSShri Abhyankar   }
66716ebf90aSShri Abhyankar 
66816ebf90aSShri Abhyankar   jj = 0; irow = rstart;
66916ebf90aSShri Abhyankar   for (i=0; i<m; i++) {
67016ebf90aSShri Abhyankar     ajj    = aj + adiag[i];                 /* ptr to the beginning of the diagonal of this row */
67116ebf90aSShri Abhyankar     v1     = av + adiag[i];
67216ebf90aSShri Abhyankar     countA = ai[i+1] - adiag[i];
67316ebf90aSShri Abhyankar     countB = bi[i+1] - bi[i];
67416ebf90aSShri Abhyankar     bjj    = bj + bi[i];
67516ebf90aSShri Abhyankar     v2     = bv + bi[i];
67616ebf90aSShri Abhyankar 
67716ebf90aSShri Abhyankar     /* A-part */
67816ebf90aSShri Abhyankar     for (j=0; j<countA; j++) {
679bccb9932SShri Abhyankar       if (reuse == MAT_INITIAL_MATRIX) {
68016ebf90aSShri Abhyankar         row[jj] = irow + shift; col[jj] = rstart + ajj[j] + shift;
68116ebf90aSShri Abhyankar       }
68216ebf90aSShri Abhyankar       val[jj++] = v1[j];
68316ebf90aSShri Abhyankar     }
68416ebf90aSShri Abhyankar 
68516ebf90aSShri Abhyankar     /* B-part */
68616ebf90aSShri Abhyankar     for (j=0; j < countB; j++) {
68716ebf90aSShri Abhyankar       if (garray[bjj[j]] > rstart) {
688bccb9932SShri Abhyankar         if (reuse == MAT_INITIAL_MATRIX) {
68916ebf90aSShri Abhyankar           row[jj] = irow + shift; col[jj] = garray[bjj[j]] + shift;
69016ebf90aSShri Abhyankar         }
69116ebf90aSShri Abhyankar         val[jj++] = v2[j];
69216ebf90aSShri Abhyankar       }
693397b6df1SKris Buschelman     }
694397b6df1SKris Buschelman     irow++;
695397b6df1SKris Buschelman   }
696397b6df1SKris Buschelman   PetscFunctionReturn(0);
697397b6df1SKris Buschelman }
698397b6df1SKris Buschelman 
699dfbe8321SBarry Smith PetscErrorCode MatDestroy_MUMPS(Mat A)
700dfbe8321SBarry Smith {
701e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
702dfbe8321SBarry Smith   PetscErrorCode ierr;
703b24902e0SBarry Smith 
704397b6df1SKris Buschelman   PetscFunctionBegin;
705a5e57a09SHong Zhang   ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
706a5e57a09SHong Zhang   ierr = VecScatterDestroy(&mumps->scat_rhs);CHKERRQ(ierr);
707a5e57a09SHong Zhang   ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
708801fbe65SHong Zhang   ierr = VecDestroy(&mumps->b_seq);CHKERRQ(ierr);
709a5e57a09SHong Zhang   ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
710a5e57a09SHong Zhang   ierr = PetscFree(mumps->id.perm_in);CHKERRQ(ierr);
711a5e57a09SHong Zhang   ierr = PetscFree(mumps->irn);CHKERRQ(ierr);
712b34f08ffSHong Zhang   ierr = PetscFree(mumps->info);CHKERRQ(ierr);
71359ac8732SStefano Zampini   ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
714a5e57a09SHong Zhang   mumps->id.job = JOB_END;
715a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
7166f3cc6f9SBarry Smith   ierr = MPI_Comm_free(&mumps->comm_mumps);CHKERRQ(ierr);
717e69c285eSBarry Smith   ierr = PetscFree(A->data);CHKERRQ(ierr);
718bf0cc555SLisandro Dalcin 
71997969023SHong Zhang   /* clear composed functions */
7203ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
7215a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr);
7225a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr);
723bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetIcntl_C",NULL);CHKERRQ(ierr);
724bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetIcntl_C",NULL);CHKERRQ(ierr);
725bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsSetCntl_C",NULL);CHKERRQ(ierr);
726bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetCntl_C",NULL);CHKERRQ(ierr);
727ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfo_C",NULL);CHKERRQ(ierr);
728ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInfog_C",NULL);CHKERRQ(ierr);
729ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfo_C",NULL);CHKERRQ(ierr);
730ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetRinfog_C",NULL);CHKERRQ(ierr);
73189a9c03aSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMumpsGetInverse_C",NULL);CHKERRQ(ierr);
732397b6df1SKris Buschelman   PetscFunctionReturn(0);
733397b6df1SKris Buschelman }
734397b6df1SKris Buschelman 
735b24902e0SBarry Smith PetscErrorCode MatSolve_MUMPS(Mat A,Vec b,Vec x)
736b24902e0SBarry Smith {
737e69c285eSBarry Smith   Mat_MUMPS        *mumps=(Mat_MUMPS*)A->data;
738d54de34fSKris Buschelman   PetscScalar      *array;
73967877ebaSShri Abhyankar   Vec              b_seq;
740329ec9b3SHong Zhang   IS               is_iden,is_petsc;
741dfbe8321SBarry Smith   PetscErrorCode   ierr;
742329ec9b3SHong Zhang   PetscInt         i;
743cc86f929SStefano Zampini   PetscBool        second_solve = PETSC_FALSE;
744883f2eb9SBarry Smith   static PetscBool cite1 = PETSC_FALSE,cite2 = PETSC_FALSE;
745397b6df1SKris Buschelman 
746397b6df1SKris Buschelman   PetscFunctionBegin;
747883f2eb9SBarry Smith   ierr = PetscCitationsRegister("@article{MUMPS01,\n  author = {P.~R. Amestoy and I.~S. Duff and J.-Y. L'Excellent and J. Koster},\n  title = {A fully asynchronous multifrontal solver using distributed dynamic scheduling},\n  journal = {SIAM Journal on Matrix Analysis and Applications},\n  volume = {23},\n  number = {1},\n  pages = {15--41},\n  year = {2001}\n}\n",&cite1);CHKERRQ(ierr);
748883f2eb9SBarry Smith   ierr = PetscCitationsRegister("@article{MUMPS02,\n  author = {P.~R. Amestoy and A. Guermouche and J.-Y. L'Excellent and S. Pralet},\n  title = {Hybrid scheduling for the parallel solution of linear systems},\n  journal = {Parallel Computing},\n  volume = {32},\n  number = {2},\n  pages = {136--156},\n  year = {2006}\n}\n",&cite2);CHKERRQ(ierr);
7492aca8efcSHong Zhang 
750603e8f96SBarry Smith   if (A->factorerrortype) {
7512aca8efcSHong Zhang     ierr = PetscInfo2(A,"MatSolve is called with singular matrix factor, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
7522aca8efcSHong Zhang     ierr = VecSetInf(x);CHKERRQ(ierr);
7532aca8efcSHong Zhang     PetscFunctionReturn(0);
7542aca8efcSHong Zhang   }
7552aca8efcSHong Zhang 
756be818407SHong Zhang   mumps->id.ICNTL(20)= 0; /* dense RHS */
757a5e57a09SHong Zhang   mumps->id.nrhs     = 1;
758a5e57a09SHong Zhang   b_seq          = mumps->b_seq;
759a5e57a09SHong Zhang   if (mumps->size > 1) {
760329ec9b3SHong Zhang     /* MUMPS only supports centralized rhs. Scatter b into a seqential rhs vector */
761a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
762a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_rhs,b,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
763a5e57a09SHong Zhang     if (!mumps->myid) {ierr = VecGetArray(b_seq,&array);CHKERRQ(ierr);}
764397b6df1SKris Buschelman   } else {  /* size == 1 */
765397b6df1SKris Buschelman     ierr = VecCopy(b,x);CHKERRQ(ierr);
766397b6df1SKris Buschelman     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
767397b6df1SKris Buschelman   }
768a5e57a09SHong Zhang   if (!mumps->myid) { /* define rhs on the host */
769a5e57a09SHong Zhang     mumps->id.nrhs = 1;
770940cd9d6SSatish Balay     mumps->id.rhs = (MumpsScalar*)array;
771397b6df1SKris Buschelman   }
772397b6df1SKris Buschelman 
773cc86f929SStefano Zampini   /*
774cc86f929SStefano Zampini      handle condensation step of Schur complement (if any)
775cc86f929SStefano Zampini      We set by default ICNTL(26) == -1 when Schur indices have been provided by the user.
776cc86f929SStefano Zampini      According to MUMPS (5.0.0) manual, any value should be harmful during the factorization phase
777cc86f929SStefano Zampini      Unless the user provides a valid value for ICNTL(26), MatSolve and MatMatSolve routines solve the full system.
778cc86f929SStefano Zampini      This requires an extra call to PetscMUMPS_c and the computation of the factors for S
779cc86f929SStefano Zampini   */
780cc86f929SStefano Zampini   if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
781241dbb5eSStefano Zampini     if (mumps->size > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
782cc86f929SStefano Zampini     second_solve = PETSC_TRUE;
783b3cb21ddSStefano Zampini     ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr);
784cc86f929SStefano Zampini   }
785397b6df1SKris Buschelman   /* solve phase */
786329ec9b3SHong Zhang   /*-------------*/
787a5e57a09SHong Zhang   mumps->id.job = JOB_SOLVE;
788a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
789a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
790397b6df1SKris Buschelman 
791b5fa320bSStefano Zampini   /* handle expansion step of Schur complement (if any) */
792cc86f929SStefano Zampini   if (second_solve) {
793b3cb21ddSStefano Zampini     ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr);
794cc86f929SStefano Zampini   }
795b5fa320bSStefano Zampini 
796a5e57a09SHong Zhang   if (mumps->size > 1) { /* convert mumps distributed solution to petsc mpi x */
797a5e57a09SHong Zhang     if (mumps->scat_sol && mumps->ICNTL9_pre != mumps->id.ICNTL(9)) {
798a5e57a09SHong Zhang       /* when id.ICNTL(9) changes, the contents of lsol_loc may change (not its size, lsol_loc), recreates scat_sol */
799a5e57a09SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
800397b6df1SKris Buschelman     }
801a5e57a09SHong Zhang     if (!mumps->scat_sol) { /* create scatter scat_sol */
802a5e57a09SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,mumps->id.lsol_loc,0,1,&is_iden);CHKERRQ(ierr); /* from */
803a5e57a09SHong Zhang       for (i=0; i<mumps->id.lsol_loc; i++) {
804a5e57a09SHong Zhang         mumps->id.isol_loc[i] -= 1; /* change Fortran style to C style */
805a5e57a09SHong Zhang       }
806a5e57a09SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,mumps->id.lsol_loc,mumps->id.isol_loc,PETSC_COPY_VALUES,&is_petsc);CHKERRQ(ierr);  /* to */
807a5e57a09SHong Zhang       ierr = VecScatterCreate(mumps->x_seq,is_iden,x,is_petsc,&mumps->scat_sol);CHKERRQ(ierr);
8086bf464f9SBarry Smith       ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
8096bf464f9SBarry Smith       ierr = ISDestroy(&is_petsc);CHKERRQ(ierr);
8102205254eSKarl Rupp 
811a5e57a09SHong Zhang       mumps->ICNTL9_pre = mumps->id.ICNTL(9); /* save current value of id.ICNTL(9) */
812397b6df1SKris Buschelman     }
813a5e57a09SHong Zhang 
814a5e57a09SHong Zhang     ierr = VecScatterBegin(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
815a5e57a09SHong Zhang     ierr = VecScatterEnd(mumps->scat_sol,mumps->x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
816329ec9b3SHong Zhang   }
817397b6df1SKris Buschelman   PetscFunctionReturn(0);
818397b6df1SKris Buschelman }
819397b6df1SKris Buschelman 
82051d5961aSHong Zhang PetscErrorCode MatSolveTranspose_MUMPS(Mat A,Vec b,Vec x)
82151d5961aSHong Zhang {
822e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
82351d5961aSHong Zhang   PetscErrorCode ierr;
82451d5961aSHong Zhang 
82551d5961aSHong Zhang   PetscFunctionBegin;
826a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 0;
8270ad0caddSJed Brown   ierr = MatSolve_MUMPS(A,b,x);CHKERRQ(ierr);
828a5e57a09SHong Zhang   mumps->id.ICNTL(9) = 1;
82951d5961aSHong Zhang   PetscFunctionReturn(0);
83051d5961aSHong Zhang }
83151d5961aSHong Zhang 
832e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_MUMPS(Mat A,Mat B,Mat X)
833e0b74bf9SHong Zhang {
834bda8bf91SBarry Smith   PetscErrorCode ierr;
835b8491c3eSStefano Zampini   Mat            Bt = NULL;
836b8491c3eSStefano Zampini   PetscBool      flg, flgT;
837e69c285eSBarry Smith   Mat_MUMPS      *mumps=(Mat_MUMPS*)A->data;
838334c5f61SHong Zhang   PetscInt       i,nrhs,M;
8392cd7d884SHong Zhang   PetscScalar    *array,*bray;
840be818407SHong Zhang   PetscInt       lsol_loc,nlsol_loc,*isol_loc,*idx,*iidx,*idxx,*isol_loc_save;
841be818407SHong Zhang   MumpsScalar    *sol_loc,*sol_loc_save;
842be818407SHong Zhang   IS             is_to,is_from;
843be818407SHong Zhang   PetscInt       k,proc,j,m;
844be818407SHong Zhang   const PetscInt *rstart;
845be818407SHong Zhang   Vec            v_mpi,b_seq,x_seq;
846be818407SHong Zhang   VecScatter     scat_rhs,scat_sol;
847be818407SHong Zhang   PetscScalar    *aa;
848be818407SHong Zhang   PetscInt       spnr,*ia,*ja;
849be818407SHong Zhang   Mat_MPIAIJ     *b;
850bda8bf91SBarry Smith 
851e0b74bf9SHong Zhang   PetscFunctionBegin;
852be818407SHong Zhang   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
853be818407SHong Zhang   if (!flg) SETERRQ(PetscObjectComm((PetscObject)X),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
854be818407SHong Zhang 
8550298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
856be818407SHong Zhang   if (flg) { /* dense B */
857*c0be3364SHong Zhang     if (B->rmap->n != X->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Matrix B and X must have same row distribution");
858be818407SHong Zhang     mumps->id.ICNTL(20)= 0; /* dense RHS */
859be818407SHong Zhang   } else {
860be818407SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&flgT);CHKERRQ(ierr);
861be818407SHong Zhang     if (flgT) { /* input B is transpose of actural RHS matrix, because mumps requires sparse compressed COLUMN storage! */
862be818407SHong Zhang       ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
863be818407SHong Zhang     } else SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
864be818407SHong Zhang     mumps->id.ICNTL(20)= 1; /* sparse RHS */
865b8491c3eSStefano Zampini   }
86687b22cf4SHong Zhang 
8679481e6e9SHong Zhang   ierr = MatGetSize(B,&M,&nrhs);CHKERRQ(ierr);
8689481e6e9SHong Zhang   mumps->id.nrhs = nrhs;
8699481e6e9SHong Zhang   mumps->id.lrhs = M;
8702b691707SHong Zhang   mumps->id.rhs  = NULL;
8719481e6e9SHong Zhang 
8722cd7d884SHong Zhang   if (mumps->size == 1) {
873b8491c3eSStefano Zampini     PetscScalar *aa;
874b8491c3eSStefano Zampini     PetscInt    spnr,*ia,*ja;
875e94cce23SStefano Zampini     PetscBool   second_solve = PETSC_FALSE;
876b8491c3eSStefano Zampini 
8772cd7d884SHong Zhang     ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
878b8491c3eSStefano Zampini     mumps->id.rhs = (MumpsScalar*)array;
8792b691707SHong Zhang 
8802b691707SHong Zhang     if (!Bt) { /* dense B */
8812b691707SHong Zhang       /* copy B to X */
882b8491c3eSStefano Zampini       ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
8836444a565SStefano Zampini       ierr = PetscMemcpy(array,bray,M*nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
8842cd7d884SHong Zhang       ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
8852b691707SHong Zhang     } else { /* sparse B */
886b8491c3eSStefano Zampini       ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr);
887be818407SHong Zhang       ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);CHKERRQ(ierr);
888*c0be3364SHong Zhang       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
8892b691707SHong Zhang       /* mumps requires ia and ja start at 1! */
890b8491c3eSStefano Zampini       mumps->id.irhs_ptr    = ia;
891b8491c3eSStefano Zampini       mumps->id.irhs_sparse = ja;
892b8491c3eSStefano Zampini       mumps->id.nz_rhs      = ia[spnr] - 1;
893b8491c3eSStefano Zampini       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
894b8491c3eSStefano Zampini     }
895e94cce23SStefano Zampini     /* handle condensation step of Schur complement (if any) */
896e94cce23SStefano Zampini     if (mumps->id.ICNTL(26) < 0 || mumps->id.ICNTL(26) > 2) {
897e94cce23SStefano Zampini       second_solve = PETSC_TRUE;
898b3cb21ddSStefano Zampini       ierr = MatMumpsHandleSchur_Private(A,PETSC_FALSE);CHKERRQ(ierr);
899e94cce23SStefano Zampini     }
9002cd7d884SHong Zhang     /* solve phase */
9012cd7d884SHong Zhang     /*-------------*/
9022cd7d884SHong Zhang     mumps->id.job = JOB_SOLVE;
9032cd7d884SHong Zhang     PetscMUMPS_c(&mumps->id);
9042cd7d884SHong Zhang     if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
905b5fa320bSStefano Zampini 
906b5fa320bSStefano Zampini     /* handle expansion step of Schur complement (if any) */
907e94cce23SStefano Zampini     if (second_solve) {
908b3cb21ddSStefano Zampini       ierr = MatMumpsHandleSchur_Private(A,PETSC_TRUE);CHKERRQ(ierr);
909e94cce23SStefano Zampini     }
9102b691707SHong Zhang     if (Bt) { /* sparse B */
911b8491c3eSStefano Zampini       ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr);
912be818407SHong Zhang       ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);CHKERRQ(ierr);
913*c0be3364SHong Zhang       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
914b8491c3eSStefano Zampini     }
9152cd7d884SHong Zhang     ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
916be818407SHong Zhang     PetscFunctionReturn(0);
917be818407SHong Zhang   }
918801fbe65SHong Zhang 
919be818407SHong Zhang   /*--------- parallel case: MUMPS requires rhs B to be centralized on the host! --------*/
92038be02acSStefano Zampini   if (mumps->size > 1 && mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Parallel Schur complements not yet supported from PETSc\n");
921241dbb5eSStefano Zampini 
922be818407SHong Zhang   /* create x_seq to hold mumps local solution */
92371aed81dSHong Zhang   isol_loc_save = mumps->id.isol_loc; /* save it for MatSovle() */
92471aed81dSHong Zhang   sol_loc_save  = mumps->id.sol_loc;
925801fbe65SHong Zhang 
92671aed81dSHong Zhang   lsol_loc  = mumps->id.INFO(23);
92771aed81dSHong Zhang   nlsol_loc = nrhs*lsol_loc;     /* length of sol_loc */
92871aed81dSHong Zhang   ierr = PetscMalloc2(nlsol_loc,&sol_loc,nlsol_loc,&isol_loc);CHKERRQ(ierr);
929940cd9d6SSatish Balay   mumps->id.sol_loc  = (MumpsScalar*)sol_loc;
930801fbe65SHong Zhang   mumps->id.isol_loc = isol_loc;
931801fbe65SHong Zhang 
9321070efccSSatish Balay   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nlsol_loc,(PetscScalar*)sol_loc,&x_seq);CHKERRQ(ierr);
9332cd7d884SHong Zhang 
934334c5f61SHong Zhang   /* scatter v_mpi to b_seq because MUMPS only supports centralized rhs */
93574f0fcc7SHong Zhang   /* idx: maps from k-th index of v_mpi to (i,j)-th global entry of B;
9362b691707SHong Zhang      iidx: inverse of idx, will be used by scattering mumps x_seq -> petsc X */
937be818407SHong Zhang   ierr = PetscMalloc1(nrhs*M,&idx);CHKERRQ(ierr);
9382cd7d884SHong Zhang 
9392b691707SHong Zhang   if (!Bt) { /* dense B */
940be818407SHong Zhang     /* wrap dense rhs matrix B into a vector v_mpi */
9412b691707SHong Zhang     ierr = MatGetLocalSize(B,&m,NULL);CHKERRQ(ierr);
9422b691707SHong Zhang     ierr = MatDenseGetArray(B,&bray);CHKERRQ(ierr);
9432b691707SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)B),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
9442b691707SHong Zhang     ierr = MatDenseRestoreArray(B,&bray);CHKERRQ(ierr);
9452b691707SHong Zhang 
946be818407SHong Zhang     /* scatter v_mpi to b_seq in proc[0]. MUMPS requires rhs to be centralized on the host! */
947801fbe65SHong Zhang     if (!mumps->myid) {
948be818407SHong Zhang       ierr = MatGetOwnershipRanges(B,&rstart);CHKERRQ(ierr);
949be818407SHong Zhang       k = 0;
950be818407SHong Zhang       for (proc=0; proc<mumps->size; proc++){
951be818407SHong Zhang         for (j=0; j<nrhs; j++){
952be818407SHong Zhang           for (i=rstart[proc]; i<rstart[proc+1]; i++){
953be818407SHong Zhang             idx[k++]      = j*M + i;
954be818407SHong Zhang           }
955be818407SHong Zhang         }
956be818407SHong Zhang       }
957be818407SHong Zhang 
958334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,nrhs*M,&b_seq);CHKERRQ(ierr);
959801fbe65SHong Zhang       ierr = ISCreateGeneral(PETSC_COMM_SELF,nrhs*M,idx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
960801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,nrhs*M,0,1,&is_from);CHKERRQ(ierr);
961801fbe65SHong Zhang     } else {
962334c5f61SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&b_seq);CHKERRQ(ierr);
963801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_to);CHKERRQ(ierr);
964801fbe65SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_from);CHKERRQ(ierr);
965801fbe65SHong Zhang     }
966334c5f61SHong Zhang     ierr = VecScatterCreate(v_mpi,is_from,b_seq,is_to,&scat_rhs);CHKERRQ(ierr);
967334c5f61SHong Zhang     ierr = VecScatterBegin(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
968801fbe65SHong Zhang     ierr = ISDestroy(&is_to);CHKERRQ(ierr);
969801fbe65SHong Zhang     ierr = ISDestroy(&is_from);CHKERRQ(ierr);
970334c5f61SHong Zhang     ierr = VecScatterEnd(scat_rhs,v_mpi,b_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
971801fbe65SHong Zhang 
972801fbe65SHong Zhang     if (!mumps->myid) { /* define rhs on the host */
973334c5f61SHong Zhang       ierr = VecGetArray(b_seq,&bray);CHKERRQ(ierr);
974940cd9d6SSatish Balay       mumps->id.rhs = (MumpsScalar*)bray;
975334c5f61SHong Zhang       ierr = VecRestoreArray(b_seq,&bray);CHKERRQ(ierr);
976801fbe65SHong Zhang     }
977801fbe65SHong Zhang 
9782b691707SHong Zhang   } else { /* sparse B */
9792b691707SHong Zhang     b = (Mat_MPIAIJ*)Bt->data;
9802b691707SHong Zhang 
981be818407SHong Zhang     /* wrap dense X into a vector v_mpi */
9822b691707SHong Zhang     ierr = MatGetLocalSize(X,&m,NULL);CHKERRQ(ierr);
9832b691707SHong Zhang     ierr = MatDenseGetArray(X,&bray);CHKERRQ(ierr);
9842b691707SHong Zhang     ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),1,nrhs*m,nrhs*M,(const PetscScalar*)bray,&v_mpi);CHKERRQ(ierr);
9852b691707SHong Zhang     ierr = MatDenseRestoreArray(X,&bray);CHKERRQ(ierr);
9862b691707SHong Zhang 
9872b691707SHong Zhang     if (!mumps->myid) {
9882b691707SHong Zhang       ierr = MatSeqAIJGetArray(b->A,&aa);CHKERRQ(ierr);
989be818407SHong Zhang       ierr = MatGetRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);CHKERRQ(ierr);
990*c0be3364SHong Zhang       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
9912b691707SHong Zhang       /* mumps requires ia and ja start at 1! */
9922b691707SHong Zhang       mumps->id.irhs_ptr    = ia;
9932b691707SHong Zhang       mumps->id.irhs_sparse = ja;
9942b691707SHong Zhang       mumps->id.nz_rhs      = ia[spnr] - 1;
9952b691707SHong Zhang       mumps->id.rhs_sparse  = (MumpsScalar*)aa;
9962b691707SHong Zhang     } else {
9972b691707SHong Zhang       mumps->id.irhs_ptr    = NULL;
9982b691707SHong Zhang       mumps->id.irhs_sparse = NULL;
9992b691707SHong Zhang       mumps->id.nz_rhs      = 0;
10002b691707SHong Zhang       mumps->id.rhs_sparse  = NULL;
10012b691707SHong Zhang     }
10022b691707SHong Zhang   }
10032b691707SHong Zhang 
1004801fbe65SHong Zhang   /* solve phase */
1005801fbe65SHong Zhang   /*-------------*/
1006801fbe65SHong Zhang   mumps->id.job = JOB_SOLVE;
1007801fbe65SHong Zhang   PetscMUMPS_c(&mumps->id);
1008801fbe65SHong Zhang   if (mumps->id.INFOG(1) < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
1009801fbe65SHong Zhang 
1010334c5f61SHong Zhang   /* scatter mumps distributed solution to petsc vector v_mpi, which shares local arrays with solution matrix X */
101174f0fcc7SHong Zhang   ierr = MatDenseGetArray(X,&array);CHKERRQ(ierr);
101274f0fcc7SHong Zhang   ierr = VecPlaceArray(v_mpi,array);CHKERRQ(ierr);
1013801fbe65SHong Zhang 
1014334c5f61SHong Zhang   /* create scatter scat_sol */
1015be818407SHong Zhang   ierr = MatGetOwnershipRanges(X,&rstart);CHKERRQ(ierr);
1016be818407SHong Zhang   /* iidx: inverse of idx computed above, used for scattering mumps x_seq to petsc X */
1017be818407SHong Zhang   iidx = idx;
1018be818407SHong Zhang   k    = 0;
1019be818407SHong Zhang   for (proc=0; proc<mumps->size; proc++){
1020be818407SHong Zhang     for (j=0; j<nrhs; j++){
1021be818407SHong Zhang       for (i=rstart[proc]; i<rstart[proc+1]; i++) iidx[j*M + i] = k++;
1022be818407SHong Zhang     }
1023be818407SHong Zhang   }
1024be818407SHong Zhang 
102571aed81dSHong Zhang   ierr = PetscMalloc1(nlsol_loc,&idxx);CHKERRQ(ierr);
102671aed81dSHong Zhang   ierr = ISCreateStride(PETSC_COMM_SELF,nlsol_loc,0,1,&is_from);CHKERRQ(ierr);
102771aed81dSHong Zhang   for (i=0; i<lsol_loc; i++) {
1028334c5f61SHong Zhang     isol_loc[i] -= 1; /* change Fortran style to C style */
1029334c5f61SHong Zhang     idxx[i] = iidx[isol_loc[i]];
1030801fbe65SHong Zhang     for (j=1; j<nrhs; j++){
1031334c5f61SHong Zhang       idxx[j*lsol_loc+i] = iidx[isol_loc[i]+j*M];
1032801fbe65SHong Zhang     }
1033801fbe65SHong Zhang   }
1034be818407SHong Zhang   ierr = ISCreateGeneral(PETSC_COMM_SELF,nlsol_loc,idxx,PETSC_COPY_VALUES,&is_to);CHKERRQ(ierr);
1035334c5f61SHong Zhang   ierr = VecScatterCreate(x_seq,is_from,v_mpi,is_to,&scat_sol);CHKERRQ(ierr);
1036334c5f61SHong Zhang   ierr = VecScatterBegin(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1037801fbe65SHong Zhang   ierr = ISDestroy(&is_from);CHKERRQ(ierr);
1038801fbe65SHong Zhang   ierr = ISDestroy(&is_to);CHKERRQ(ierr);
1039334c5f61SHong Zhang   ierr = VecScatterEnd(scat_sol,x_seq,v_mpi,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1040801fbe65SHong Zhang   ierr = MatDenseRestoreArray(X,&array);CHKERRQ(ierr);
104171aed81dSHong Zhang 
104271aed81dSHong Zhang   /* free spaces */
104371aed81dSHong Zhang   mumps->id.sol_loc  = sol_loc_save;
104471aed81dSHong Zhang   mumps->id.isol_loc = isol_loc_save;
104571aed81dSHong Zhang 
104671aed81dSHong Zhang   ierr = PetscFree2(sol_loc,isol_loc);CHKERRQ(ierr);
1047be818407SHong Zhang   ierr = PetscFree(idx);CHKERRQ(ierr);
1048801fbe65SHong Zhang   ierr = PetscFree(idxx);CHKERRQ(ierr);
104971aed81dSHong Zhang   ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
105074f0fcc7SHong Zhang   ierr = VecDestroy(&v_mpi);CHKERRQ(ierr);
10512b691707SHong Zhang   if (Bt) {
10522b691707SHong Zhang     if (!mumps->myid) {
10532b691707SHong Zhang       ierr = MatSeqAIJRestoreArray(b->A,&aa);CHKERRQ(ierr);
1054be818407SHong Zhang       ierr = MatRestoreRowIJ(b->A,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&flg);CHKERRQ(ierr);
1055*c0be3364SHong Zhang       if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot restore IJ structure");
10562b691707SHong Zhang     }
10572b691707SHong Zhang   } else {
1058334c5f61SHong Zhang     ierr = VecDestroy(&b_seq);CHKERRQ(ierr);
1059334c5f61SHong Zhang     ierr = VecScatterDestroy(&scat_rhs);CHKERRQ(ierr);
10602b691707SHong Zhang   }
1061334c5f61SHong Zhang   ierr = VecScatterDestroy(&scat_sol);CHKERRQ(ierr);
1062e0b74bf9SHong Zhang   PetscFunctionReturn(0);
1063e0b74bf9SHong Zhang }
1064e0b74bf9SHong Zhang 
1065ace3df97SHong Zhang #if !defined(PETSC_USE_COMPLEX)
1066a58c3f20SHong Zhang /*
1067a58c3f20SHong Zhang   input:
1068a58c3f20SHong Zhang    F:        numeric factor
1069a58c3f20SHong Zhang   output:
1070a58c3f20SHong Zhang    nneg:     total number of negative pivots
107119d49a3bSHong Zhang    nzero:    total number of zero pivots
107219d49a3bSHong Zhang    npos:     (global dimension of F) - nneg - nzero
1073a58c3f20SHong Zhang */
1074dfbe8321SBarry Smith PetscErrorCode MatGetInertia_SBAIJMUMPS(Mat F,int *nneg,int *nzero,int *npos)
1075a58c3f20SHong Zhang {
1076e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
1077dfbe8321SBarry Smith   PetscErrorCode ierr;
1078c1490034SHong Zhang   PetscMPIInt    size;
1079a58c3f20SHong Zhang 
1080a58c3f20SHong Zhang   PetscFunctionBegin;
1081ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&size);CHKERRQ(ierr);
1082bcb30aebSHong 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 */
1083a5e57a09SHong Zhang   if (size > 1 && mumps->id.ICNTL(13) != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_mumps_icntl_13 must be set as 1 for correct global matrix inertia\n",mumps->id.INFOG(13));
1084ed85ac9fSHong Zhang 
1085710ac8efSHong Zhang   if (nneg) *nneg = mumps->id.INFOG(12);
1086ed85ac9fSHong Zhang   if (nzero || npos) {
1087ed85ac9fSHong Zhang     if (mumps->id.ICNTL(24) != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"-mat_mumps_icntl_24 must be set as 1 for null pivot row detection");
1088710ac8efSHong Zhang     if (nzero) *nzero = mumps->id.INFOG(28);
1089710ac8efSHong Zhang     if (npos) *npos   = F->rmap->N - (mumps->id.INFOG(12) + mumps->id.INFOG(28));
1090a58c3f20SHong Zhang   }
1091a58c3f20SHong Zhang   PetscFunctionReturn(0);
1092a58c3f20SHong Zhang }
109319d49a3bSHong Zhang #endif
1094a58c3f20SHong Zhang 
10950481f469SBarry Smith PetscErrorCode MatFactorNumeric_MUMPS(Mat F,Mat A,const MatFactorInfo *info)
1096af281ebdSHong Zhang {
1097e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)(F)->data;
10986849ba73SBarry Smith   PetscErrorCode ierr;
1099ace3abfcSBarry Smith   PetscBool      isMPIAIJ;
1100397b6df1SKris Buschelman 
1101397b6df1SKris Buschelman   PetscFunctionBegin;
11026baea169SHong Zhang   if (mumps->id.INFOG(1) < 0) {
11032aca8efcSHong Zhang     if (mumps->id.INFOG(1) == -6) {
11042aca8efcSHong Zhang       ierr = PetscInfo2(A,"MatFactorNumeric is called with singular matrix structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
11056baea169SHong Zhang     }
11066baea169SHong Zhang     ierr = PetscInfo2(A,"MatFactorNumeric is called after analysis phase fails, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
11072aca8efcSHong Zhang     PetscFunctionReturn(0);
11082aca8efcSHong Zhang   }
11096baea169SHong Zhang 
1110a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_REUSE_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1111397b6df1SKris Buschelman 
1112397b6df1SKris Buschelman   /* numerical factorization phase */
1113329ec9b3SHong Zhang   /*-------------------------------*/
1114a5e57a09SHong Zhang   mumps->id.job = JOB_FACTNUMERIC;
11154e34a73bSHong Zhang   if (!mumps->id.ICNTL(18)) { /* A is centralized */
1116a5e57a09SHong Zhang     if (!mumps->myid) {
1117940cd9d6SSatish Balay       mumps->id.a = (MumpsScalar*)mumps->val;
1118397b6df1SKris Buschelman     }
1119397b6df1SKris Buschelman   } else {
1120940cd9d6SSatish Balay     mumps->id.a_loc = (MumpsScalar*)mumps->val;
1121397b6df1SKris Buschelman   }
1122a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
1123a5e57a09SHong Zhang   if (mumps->id.INFOG(1) < 0) {
1124c0d63f2fSHong Zhang     if (A->erroriffailure) {
1125c0d63f2fSHong Zhang       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
1126151787a6SHong Zhang     } else {
1127c0d63f2fSHong Zhang       if (mumps->id.INFOG(1) == -10) { /* numerically singular matrix */
11282aca8efcSHong Zhang         ierr = PetscInfo2(F,"matrix is numerically singular, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1129603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
1130c0d63f2fSHong Zhang       } else if (mumps->id.INFOG(1) == -13) {
1131c0d63f2fSHong Zhang         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, cannot allocate required memory %d megabytes\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1132603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
1133c0d63f2fSHong Zhang       } else if (mumps->id.INFOG(1) == -8 || mumps->id.INFOG(1) == -9 || (-16 < mumps->id.INFOG(1) && mumps->id.INFOG(1) < -10) ) {
1134c0d63f2fSHong Zhang         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d, problem with workarray \n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1135603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
11362aca8efcSHong Zhang       } else {
1137c0d63f2fSHong Zhang         ierr = PetscInfo2(F,"MUMPS in numerical factorization phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1138603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OTHER;
1139151787a6SHong Zhang       }
11402aca8efcSHong Zhang     }
1141397b6df1SKris Buschelman   }
1142a5e57a09SHong Zhang   if (!mumps->myid && mumps->id.ICNTL(16) > 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"  mumps->id.ICNTL(16):=%d\n",mumps->id.INFOG(16));
1143397b6df1SKris Buschelman 
1144b3cb21ddSStefano Zampini   F->assembled    = PETSC_TRUE;
1145a5e57a09SHong Zhang   mumps->matstruc = SAME_NONZERO_PATTERN;
1146b3cb21ddSStefano Zampini   if (F->schur) { /* reset Schur status to unfactored */
1147b3cb21ddSStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
1148b3cb21ddSStefano Zampini       mumps->id.ICNTL(19) = 2;
1149b3cb21ddSStefano Zampini       ierr = MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);CHKERRQ(ierr);
1150b3cb21ddSStefano Zampini     }
1151b3cb21ddSStefano Zampini     ierr = MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
1152b3cb21ddSStefano Zampini   }
115367877ebaSShri Abhyankar 
1154066565c5SStefano Zampini   /* just to be sure that ICNTL(19) value returned by a call from MatMumpsGetIcntl is always consistent */
1155066565c5SStefano Zampini   if (!mumps->sym && mumps->id.ICNTL(19) && mumps->id.ICNTL(19) != 1) mumps->id.ICNTL(19) = 3;
1156066565c5SStefano Zampini 
1157a5e57a09SHong Zhang   if (mumps->size > 1) {
115867877ebaSShri Abhyankar     PetscInt    lsol_loc;
115967877ebaSShri Abhyankar     PetscScalar *sol_loc;
11602205254eSKarl Rupp 
1161c2093ab7SHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
1162c2093ab7SHong Zhang 
1163c2093ab7SHong Zhang     /* distributed solution; Create x_seq=sol_loc for repeated use */
1164c2093ab7SHong Zhang     if (mumps->x_seq) {
1165c2093ab7SHong Zhang       ierr = VecScatterDestroy(&mumps->scat_sol);CHKERRQ(ierr);
1166c2093ab7SHong Zhang       ierr = PetscFree2(mumps->id.sol_loc,mumps->id.isol_loc);CHKERRQ(ierr);
1167c2093ab7SHong Zhang       ierr = VecDestroy(&mumps->x_seq);CHKERRQ(ierr);
1168c2093ab7SHong Zhang     }
1169a5e57a09SHong Zhang     lsol_loc = mumps->id.INFO(23); /* length of sol_loc */
1170dcca6d9dSJed Brown     ierr = PetscMalloc2(lsol_loc,&sol_loc,lsol_loc,&mumps->id.isol_loc);CHKERRQ(ierr);
1171a5e57a09SHong Zhang     mumps->id.lsol_loc = lsol_loc;
1172940cd9d6SSatish Balay     mumps->id.sol_loc = (MumpsScalar*)sol_loc;
1173a5e57a09SHong Zhang     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,lsol_loc,sol_loc,&mumps->x_seq);CHKERRQ(ierr);
117467877ebaSShri Abhyankar   }
1175397b6df1SKris Buschelman   PetscFunctionReturn(0);
1176397b6df1SKris Buschelman }
1177397b6df1SKris Buschelman 
11789a2535b5SHong Zhang /* Sets MUMPS options from the options database */
11799a2535b5SHong Zhang PetscErrorCode PetscSetMUMPSFromOptions(Mat F, Mat A)
1180dcd589f8SShri Abhyankar {
1181e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1182dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1183b34f08ffSHong Zhang   PetscInt       icntl,info[40],i,ninfo=40;
1184ace3abfcSBarry Smith   PetscBool      flg;
1185dcd589f8SShri Abhyankar 
1186dcd589f8SShri Abhyankar   PetscFunctionBegin;
1187ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MUMPS Options","Mat");CHKERRQ(ierr);
11889a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_1","ICNTL(1): output stream for error messages","None",mumps->id.ICNTL(1),&icntl,&flg);CHKERRQ(ierr);
11899a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(1) = icntl;
11909a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_2","ICNTL(2): output stream for diagnostic printing, statistics, and warning","None",mumps->id.ICNTL(2),&icntl,&flg);CHKERRQ(ierr);
11919a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(2) = icntl;
11929a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_3","ICNTL(3): output stream for global information, collected on the host","None",mumps->id.ICNTL(3),&icntl,&flg);CHKERRQ(ierr);
11939a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(3) = icntl;
1194dcd589f8SShri Abhyankar 
11959a2535b5SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_4","ICNTL(4): level of printing (0 to 4)","None",mumps->id.ICNTL(4),&icntl,&flg);CHKERRQ(ierr);
11969a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(4) = icntl;
11979a2535b5SHong Zhang   if (mumps->id.ICNTL(4) || PetscLogPrintInfo) mumps->id.ICNTL(3) = 6; /* resume MUMPS default id.ICNTL(3) = 6 */
11989a2535b5SHong Zhang 
1199d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_6","ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)","None",mumps->id.ICNTL(6),&icntl,&flg);CHKERRQ(ierr);
12009a2535b5SHong Zhang   if (flg) mumps->id.ICNTL(6) = icntl;
12019a2535b5SHong Zhang 
1202d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_7","ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis","None",mumps->id.ICNTL(7),&icntl,&flg);CHKERRQ(ierr);
1203dcd589f8SShri Abhyankar   if (flg) {
12042205254eSKarl Rupp     if (icntl== 1 && mumps->size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"pivot order be set by the user in PERM_IN -- not supported by the PETSc/MUMPS interface\n");
12052205254eSKarl Rupp     else mumps->id.ICNTL(7) = icntl;
1206dcd589f8SShri Abhyankar   }
1207e0b74bf9SHong Zhang 
12080298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_8","ICNTL(8): scaling strategy (-2 to 8 or 77)","None",mumps->id.ICNTL(8),&mumps->id.ICNTL(8),NULL);CHKERRQ(ierr);
1209d341cd04SHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_9","ICNTL(9): computes the solution using A or A^T","None",mumps->id.ICNTL(9),&mumps->id.ICNTL(9),NULL);CHKERRQ(ierr); handled by MatSolveTranspose_MUMPS() */
12100298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_10","ICNTL(10): max num of refinements","None",mumps->id.ICNTL(10),&mumps->id.ICNTL(10),NULL);CHKERRQ(ierr);
1211d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_11","ICNTL(11): statistics related to an error analysis (via -ksp_view)","None",mumps->id.ICNTL(11),&mumps->id.ICNTL(11),NULL);CHKERRQ(ierr);
1212d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_12","ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)","None",mumps->id.ICNTL(12),&mumps->id.ICNTL(12),NULL);CHKERRQ(ierr);
1213d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_13","ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting","None",mumps->id.ICNTL(13),&mumps->id.ICNTL(13),NULL);CHKERRQ(ierr);
1214d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_14","ICNTL(14): percentage increase in the estimated working space","None",mumps->id.ICNTL(14),&mumps->id.ICNTL(14),NULL);CHKERRQ(ierr);
1215d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_19","ICNTL(19): computes the Schur complement","None",mumps->id.ICNTL(19),&mumps->id.ICNTL(19),NULL);CHKERRQ(ierr);
121659ac8732SStefano Zampini   if (mumps->id.ICNTL(19) <= 0 || mumps->id.ICNTL(19) > 3) { /* reset any schur data (if any) */
1217b3cb21ddSStefano Zampini     ierr = MatDestroy(&F->schur);CHKERRQ(ierr);
121859ac8732SStefano Zampini     ierr = MatMumpsResetSchur_Private(mumps);CHKERRQ(ierr);
121959ac8732SStefano Zampini   }
12204e34a73bSHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_20","ICNTL(20): the format (dense or sparse) of the right-hand sides","None",mumps->id.ICNTL(20),&mumps->id.ICNTL(20),NULL);CHKERRQ(ierr); -- sparse rhs is not supported in PETSc API */
1221d341cd04SHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_21","ICNTL(21): the distribution (centralized or distributed) of the solution vectors","None",mumps->id.ICNTL(21),&mumps->id.ICNTL(21),NULL);CHKERRQ(ierr); we only use distributed solution vector */
12229a2535b5SHong Zhang 
1223d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_22","ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)","None",mumps->id.ICNTL(22),&mumps->id.ICNTL(22),NULL);CHKERRQ(ierr);
12240298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_23","ICNTL(23): max size of the working memory (MB) that can allocate per processor","None",mumps->id.ICNTL(23),&mumps->id.ICNTL(23),NULL);CHKERRQ(ierr);
12250298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_24","ICNTL(24): detection of null pivot rows (0 or 1)","None",mumps->id.ICNTL(24),&mumps->id.ICNTL(24),NULL);CHKERRQ(ierr);
12269a2535b5SHong Zhang   if (mumps->id.ICNTL(24)) {
12279a2535b5SHong Zhang     mumps->id.ICNTL(13) = 1; /* turn-off ScaLAPACK to help with the correct detection of null pivots */
1228d7ebd59bSHong Zhang   }
1229d7ebd59bSHong Zhang 
1230b4ed93dbSHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_25","ICNTL(25): computes a solution of a deficient matrix and a null space basis","None",mumps->id.ICNTL(25),&mumps->id.ICNTL(25),NULL);CHKERRQ(ierr);
1231d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_26","ICNTL(26): drives the solution phase if a Schur complement matrix","None",mumps->id.ICNTL(26),&mumps->id.ICNTL(26),NULL);CHKERRQ(ierr);
12322cd7d884SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_27","ICNTL(27): the blocking size for multiple right-hand sides","None",mumps->id.ICNTL(27),&mumps->id.ICNTL(27),NULL);CHKERRQ(ierr);
12330298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_28","ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering","None",mumps->id.ICNTL(28),&mumps->id.ICNTL(28),NULL);CHKERRQ(ierr);
1234d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_29","ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis","None",mumps->id.ICNTL(29),&mumps->id.ICNTL(29),NULL);CHKERRQ(ierr);
123589a9c03aSHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_30","ICNTL(30): compute user-specified set of entries in inv(A)","None",mumps->id.ICNTL(30),&mumps->id.ICNTL(30),NULL);CHKERRQ(ierr); */ /* call MatMumpsGetInverse() directly */
1236d341cd04SHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_31","ICNTL(31): indicates which factors may be discarded during factorization","None",mumps->id.ICNTL(31),&mumps->id.ICNTL(31),NULL);CHKERRQ(ierr);
12374e34a73bSHong Zhang   /* ierr = PetscOptionsInt("-mat_mumps_icntl_32","ICNTL(32): performs the forward elemination of the right-hand sides during factorization","None",mumps->id.ICNTL(32),&mumps->id.ICNTL(32),NULL);CHKERRQ(ierr);  -- not supported by PETSc API */
12380298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_mumps_icntl_33","ICNTL(33): compute determinant","None",mumps->id.ICNTL(33),&mumps->id.ICNTL(33),NULL);CHKERRQ(ierr);
1239b4ed93dbSHong Zhang   ierr = PetscOptionsInt("-mat_mumps_icntl_35","ICNTL(35): activates Block Lock Rank (BLR) based factorization","None",mumps->id.ICNTL(35),&mumps->id.ICNTL(35),NULL);CHKERRQ(ierr);
1240dcd589f8SShri Abhyankar 
12410298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_1","CNTL(1): relative pivoting threshold","None",mumps->id.CNTL(1),&mumps->id.CNTL(1),NULL);CHKERRQ(ierr);
12420298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_2","CNTL(2): stopping criterion of refinement","None",mumps->id.CNTL(2),&mumps->id.CNTL(2),NULL);CHKERRQ(ierr);
12430298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_3","CNTL(3): absolute pivoting threshold","None",mumps->id.CNTL(3),&mumps->id.CNTL(3),NULL);CHKERRQ(ierr);
12440298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_4","CNTL(4): value for static pivoting","None",mumps->id.CNTL(4),&mumps->id.CNTL(4),NULL);CHKERRQ(ierr);
12450298fd71SBarry Smith   ierr = PetscOptionsReal("-mat_mumps_cntl_5","CNTL(5): fixation for null pivots","None",mumps->id.CNTL(5),&mumps->id.CNTL(5),NULL);CHKERRQ(ierr);
1246b4ed93dbSHong Zhang   ierr = PetscOptionsReal("-mat_mumps_cntl_7","CNTL(7): dropping parameter used during BLR","None",mumps->id.CNTL(7),&mumps->id.CNTL(7),NULL);CHKERRQ(ierr);
1247e5bb22a1SHong Zhang 
12482a808120SBarry Smith   ierr = PetscOptionsString("-mat_mumps_ooc_tmpdir", "out of core directory", "None", mumps->id.ooc_tmpdir, mumps->id.ooc_tmpdir, 256, NULL);CHKERRQ(ierr);
1249b34f08ffSHong Zhang 
125016d797efSHong Zhang   ierr = PetscOptionsIntArray("-mat_mumps_view_info","request INFO local to each processor","",info,&ninfo,NULL);CHKERRQ(ierr);
1251b34f08ffSHong Zhang   if (ninfo) {
1252b34f08ffSHong Zhang     if (ninfo > 40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"number of INFO %d must <= 40\n",ninfo);
1253b34f08ffSHong Zhang     ierr = PetscMalloc1(ninfo,&mumps->info);CHKERRQ(ierr);
1254b34f08ffSHong Zhang     mumps->ninfo = ninfo;
1255b34f08ffSHong Zhang     for (i=0; i<ninfo; i++) {
12566c4ed002SBarry Smith       if (info[i] < 0 || info[i]>40) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"index of INFO %d must between 1 and 40\n",ninfo);
12572a808120SBarry Smith       else  mumps->info[i] = info[i];
1258b34f08ffSHong Zhang     }
1259b34f08ffSHong Zhang   }
1260b34f08ffSHong Zhang 
12612a808120SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1262dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1263dcd589f8SShri Abhyankar }
1264dcd589f8SShri Abhyankar 
1265f697e70eSHong Zhang PetscErrorCode PetscInitializeMUMPS(Mat A,Mat_MUMPS *mumps)
1266dcd589f8SShri Abhyankar {
1267dcd589f8SShri Abhyankar   PetscErrorCode ierr;
1268dcd589f8SShri Abhyankar 
1269dcd589f8SShri Abhyankar   PetscFunctionBegin;
12702a808120SBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)A), &mumps->myid);CHKERRQ(ierr);
1271ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&mumps->size);CHKERRQ(ierr);
1272ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(mumps->comm_mumps));CHKERRQ(ierr);
12732205254eSKarl Rupp 
1274f697e70eSHong Zhang   mumps->id.comm_fortran = MPI_Comm_c2f(mumps->comm_mumps);
1275f697e70eSHong Zhang 
1276f697e70eSHong Zhang   mumps->id.job = JOB_INIT;
1277f697e70eSHong Zhang   mumps->id.par = 1;  /* host participates factorizaton and solve */
1278f697e70eSHong Zhang   mumps->id.sym = mumps->sym;
12792907cef9SHong Zhang   PetscMUMPS_c(&mumps->id);
1280f697e70eSHong Zhang 
12810298fd71SBarry Smith   mumps->scat_rhs     = NULL;
12820298fd71SBarry Smith   mumps->scat_sol     = NULL;
12839a2535b5SHong Zhang 
128470544d5fSHong Zhang   /* set PETSc-MUMPS default options - override MUMPS default */
12859a2535b5SHong Zhang   mumps->id.ICNTL(3) = 0;
12869a2535b5SHong Zhang   mumps->id.ICNTL(4) = 0;
12879a2535b5SHong Zhang   if (mumps->size == 1) {
12889a2535b5SHong Zhang     mumps->id.ICNTL(18) = 0;   /* centralized assembled matrix input */
12899a2535b5SHong Zhang   } else {
12909a2535b5SHong Zhang     mumps->id.ICNTL(18) = 3;   /* distributed assembled matrix input */
12914e34a73bSHong Zhang     mumps->id.ICNTL(20) = 0;   /* rhs is in dense format */
129270544d5fSHong Zhang     mumps->id.ICNTL(21) = 1;   /* distributed solution */
12939a2535b5SHong Zhang   }
12946444a565SStefano Zampini 
12956444a565SStefano Zampini   /* schur */
12966444a565SStefano Zampini   mumps->id.size_schur      = 0;
12976444a565SStefano Zampini   mumps->id.listvar_schur   = NULL;
12986444a565SStefano Zampini   mumps->id.schur           = NULL;
1299b5fa320bSStefano Zampini   mumps->sizeredrhs         = 0;
130059ac8732SStefano Zampini   mumps->schur_sol          = NULL;
130159ac8732SStefano Zampini   mumps->schur_sizesol      = 0;
1302dcd589f8SShri Abhyankar   PetscFunctionReturn(0);
1303dcd589f8SShri Abhyankar }
1304dcd589f8SShri Abhyankar 
13059a625307SHong Zhang PetscErrorCode MatFactorSymbolic_MUMPS_ReportIfError(Mat F,Mat A,const MatFactorInfo *info,Mat_MUMPS *mumps)
13065cd7cf9dSHong Zhang {
13075cd7cf9dSHong Zhang   PetscErrorCode ierr;
13085cd7cf9dSHong Zhang 
13095cd7cf9dSHong Zhang   PetscFunctionBegin;
13105cd7cf9dSHong Zhang   if (mumps->id.INFOG(1) < 0) {
13115cd7cf9dSHong Zhang     if (A->erroriffailure) {
13125cd7cf9dSHong Zhang       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in analysis phase: INFOG(1)=%d\n",mumps->id.INFOG(1));
13135cd7cf9dSHong Zhang     } else {
13145cd7cf9dSHong Zhang       if (mumps->id.INFOG(1) == -6) {
13155cd7cf9dSHong Zhang         ierr = PetscInfo2(F,"matrix is singular in structure, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1316603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_STRUCT_ZEROPIVOT;
13175cd7cf9dSHong Zhang       } else if (mumps->id.INFOG(1) == -5 || mumps->id.INFOG(1) == -7) {
13185cd7cf9dSHong Zhang         ierr = PetscInfo2(F,"problem of workspace, INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1319603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
13205cd7cf9dSHong Zhang       } else {
13215cd7cf9dSHong Zhang         ierr = PetscInfo2(F,"Error reported by MUMPS in analysis phase: INFOG(1)=%d, INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));CHKERRQ(ierr);
1322603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OTHER;
13235cd7cf9dSHong Zhang       }
13245cd7cf9dSHong Zhang     }
13255cd7cf9dSHong Zhang   }
13265cd7cf9dSHong Zhang   PetscFunctionReturn(0);
13275cd7cf9dSHong Zhang }
13285cd7cf9dSHong Zhang 
1329a5e57a09SHong Zhang /* Note Petsc r(=c) permutation is used when mumps->id.ICNTL(7)==1 with centralized assembled matrix input; otherwise r and c are ignored */
13300481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_AIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1331b24902e0SBarry Smith {
1332e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1333dcd589f8SShri Abhyankar   PetscErrorCode ierr;
133467877ebaSShri Abhyankar   Vec            b;
133567877ebaSShri Abhyankar   IS             is_iden;
133667877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1337397b6df1SKris Buschelman 
1338397b6df1SKris Buschelman   PetscFunctionBegin;
1339a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1340dcd589f8SShri Abhyankar 
13419a2535b5SHong Zhang   /* Set MUMPS options from the options database */
13429a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1343dcd589f8SShri Abhyankar 
1344a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1345dcd589f8SShri Abhyankar 
134667877ebaSShri Abhyankar   /* analysis phase */
134767877ebaSShri Abhyankar   /*----------------*/
1348a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1349a5e57a09SHong Zhang   mumps->id.n   = M;
1350a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
135167877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1352a5e57a09SHong Zhang     if (!mumps->myid) {
1353a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1354a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1355940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
135667877ebaSShri Abhyankar       }
1357a5e57a09SHong Zhang       if (mumps->id.ICNTL(7) == 1) { /* use user-provide matrix ordering - assuming r = c ordering */
13585248a706SHong Zhang         /*
13595248a706SHong Zhang         PetscBool      flag;
13605248a706SHong Zhang         ierr = ISEqual(r,c,&flag);CHKERRQ(ierr);
13615248a706SHong Zhang         if (!flag) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"row_perm != col_perm");
13625248a706SHong Zhang         ierr = ISView(r,PETSC_VIEWER_STDOUT_SELF);
13635248a706SHong Zhang          */
1364a5e57a09SHong Zhang         if (!mumps->myid) {
1365e0b74bf9SHong Zhang           const PetscInt *idx;
1366e0b74bf9SHong Zhang           PetscInt       i,*perm_in;
13672205254eSKarl Rupp 
1368785e854fSJed Brown           ierr = PetscMalloc1(M,&perm_in);CHKERRQ(ierr);
1369e0b74bf9SHong Zhang           ierr = ISGetIndices(r,&idx);CHKERRQ(ierr);
13702205254eSKarl Rupp 
1371a5e57a09SHong Zhang           mumps->id.perm_in = perm_in;
1372e0b74bf9SHong Zhang           for (i=0; i<M; i++) perm_in[i] = idx[i]+1; /* perm_in[]: start from 1, not 0! */
1373e0b74bf9SHong Zhang           ierr = ISRestoreIndices(r,&idx);CHKERRQ(ierr);
1374e0b74bf9SHong Zhang         }
1375e0b74bf9SHong Zhang       }
137667877ebaSShri Abhyankar     }
137767877ebaSShri Abhyankar     break;
137867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1379a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1380a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1381a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1382940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
138367877ebaSShri Abhyankar     }
138467877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1385a5e57a09SHong Zhang     if (!mumps->myid) {
13862cd7d884SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->rmap->N,&mumps->b_seq);CHKERRQ(ierr);
13872cd7d884SHong Zhang       ierr = ISCreateStride(PETSC_COMM_SELF,A->rmap->N,0,1,&is_iden);CHKERRQ(ierr);
138867877ebaSShri Abhyankar     } else {
1389a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
139067877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
139167877ebaSShri Abhyankar     }
13922a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1393a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
13946bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
13956bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
139667877ebaSShri Abhyankar     break;
139767877ebaSShri Abhyankar   }
1398a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
13995cd7cf9dSHong Zhang   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
140067877ebaSShri Abhyankar 
1401719d5645SBarry Smith   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1402dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
140351d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
14044e34a73bSHong Zhang   F->ops->matsolve        = MatMatSolve_MUMPS;
1405b24902e0SBarry Smith   PetscFunctionReturn(0);
1406b24902e0SBarry Smith }
1407b24902e0SBarry Smith 
1408450b117fSShri Abhyankar /* Note the Petsc r and c permutations are ignored */
1409450b117fSShri Abhyankar PetscErrorCode MatLUFactorSymbolic_BAIJMUMPS(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1410450b117fSShri Abhyankar {
1411e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1412dcd589f8SShri Abhyankar   PetscErrorCode ierr;
141367877ebaSShri Abhyankar   Vec            b;
141467877ebaSShri Abhyankar   IS             is_iden;
141567877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1416450b117fSShri Abhyankar 
1417450b117fSShri Abhyankar   PetscFunctionBegin;
1418a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1419dcd589f8SShri Abhyankar 
14209a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14219a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1422dcd589f8SShri Abhyankar 
1423a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
142467877ebaSShri Abhyankar 
142567877ebaSShri Abhyankar   /* analysis phase */
142667877ebaSShri Abhyankar   /*----------------*/
1427a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1428a5e57a09SHong Zhang   mumps->id.n   = M;
1429a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
143067877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1431a5e57a09SHong Zhang     if (!mumps->myid) {
1432a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1433a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1434940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
143567877ebaSShri Abhyankar       }
143667877ebaSShri Abhyankar     }
143767877ebaSShri Abhyankar     break;
143867877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1439a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1440a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1441a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1442940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
144367877ebaSShri Abhyankar     }
144467877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1445a5e57a09SHong Zhang     if (!mumps->myid) {
1446a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
144767877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
144867877ebaSShri Abhyankar     } else {
1449a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
145067877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
145167877ebaSShri Abhyankar     }
14522a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1453a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
14546bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
14556bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
145667877ebaSShri Abhyankar     break;
145767877ebaSShri Abhyankar   }
1458a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
14595cd7cf9dSHong Zhang   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
146067877ebaSShri Abhyankar 
1461450b117fSShri Abhyankar   F->ops->lufactornumeric = MatFactorNumeric_MUMPS;
1462dcd589f8SShri Abhyankar   F->ops->solve           = MatSolve_MUMPS;
146351d5961aSHong Zhang   F->ops->solvetranspose  = MatSolveTranspose_MUMPS;
1464450b117fSShri Abhyankar   PetscFunctionReturn(0);
1465450b117fSShri Abhyankar }
1466b24902e0SBarry Smith 
1467141f4205SHong Zhang /* Note the Petsc r permutation and factor info are ignored */
146867877ebaSShri Abhyankar PetscErrorCode MatCholeskyFactorSymbolic_MUMPS(Mat F,Mat A,IS r,const MatFactorInfo *info)
1469b24902e0SBarry Smith {
1470e69c285eSBarry Smith   Mat_MUMPS      *mumps = (Mat_MUMPS*)F->data;
1471dcd589f8SShri Abhyankar   PetscErrorCode ierr;
147267877ebaSShri Abhyankar   Vec            b;
147367877ebaSShri Abhyankar   IS             is_iden;
147467877ebaSShri Abhyankar   const PetscInt M = A->rmap->N;
1475397b6df1SKris Buschelman 
1476397b6df1SKris Buschelman   PetscFunctionBegin;
1477a5e57a09SHong Zhang   mumps->matstruc = DIFFERENT_NONZERO_PATTERN;
1478dcd589f8SShri Abhyankar 
14799a2535b5SHong Zhang   /* Set MUMPS options from the options database */
14809a2535b5SHong Zhang   ierr = PetscSetMUMPSFromOptions(F,A);CHKERRQ(ierr);
1481dcd589f8SShri Abhyankar 
1482a5e57a09SHong Zhang   ierr = (*mumps->ConvertToTriples)(A, 1, MAT_INITIAL_MATRIX, &mumps->nz, &mumps->irn, &mumps->jcn, &mumps->val);CHKERRQ(ierr);
1483dcd589f8SShri Abhyankar 
148467877ebaSShri Abhyankar   /* analysis phase */
148567877ebaSShri Abhyankar   /*----------------*/
1486a5e57a09SHong Zhang   mumps->id.job = JOB_FACTSYMBOLIC;
1487a5e57a09SHong Zhang   mumps->id.n   = M;
1488a5e57a09SHong Zhang   switch (mumps->id.ICNTL(18)) {
148967877ebaSShri Abhyankar   case 0:  /* centralized assembled matrix input */
1490a5e57a09SHong Zhang     if (!mumps->myid) {
1491a5e57a09SHong Zhang       mumps->id.nz =mumps->nz; mumps->id.irn=mumps->irn; mumps->id.jcn=mumps->jcn;
1492a5e57a09SHong Zhang       if (mumps->id.ICNTL(6)>1) {
1493940cd9d6SSatish Balay         mumps->id.a = (MumpsScalar*)mumps->val;
149467877ebaSShri Abhyankar       }
149567877ebaSShri Abhyankar     }
149667877ebaSShri Abhyankar     break;
149767877ebaSShri Abhyankar   case 3:  /* distributed assembled matrix input (size>1) */
1498a5e57a09SHong Zhang     mumps->id.nz_loc = mumps->nz;
1499a5e57a09SHong Zhang     mumps->id.irn_loc=mumps->irn; mumps->id.jcn_loc=mumps->jcn;
1500a5e57a09SHong Zhang     if (mumps->id.ICNTL(6)>1) {
1501940cd9d6SSatish Balay       mumps->id.a_loc = (MumpsScalar*)mumps->val;
150267877ebaSShri Abhyankar     }
150367877ebaSShri Abhyankar     /* MUMPS only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
1504a5e57a09SHong Zhang     if (!mumps->myid) {
1505a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&mumps->b_seq);CHKERRQ(ierr);
150667877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
150767877ebaSShri Abhyankar     } else {
1508a5e57a09SHong Zhang       ierr = VecCreateSeq(PETSC_COMM_SELF,0,&mumps->b_seq);CHKERRQ(ierr);
150967877ebaSShri Abhyankar       ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&is_iden);CHKERRQ(ierr);
151067877ebaSShri Abhyankar     }
15112a7a6963SBarry Smith     ierr = MatCreateVecs(A,NULL,&b);CHKERRQ(ierr);
1512a5e57a09SHong Zhang     ierr = VecScatterCreate(b,is_iden,mumps->b_seq,is_iden,&mumps->scat_rhs);CHKERRQ(ierr);
15136bf464f9SBarry Smith     ierr = ISDestroy(&is_iden);CHKERRQ(ierr);
15146bf464f9SBarry Smith     ierr = VecDestroy(&b);CHKERRQ(ierr);
151567877ebaSShri Abhyankar     break;
151667877ebaSShri Abhyankar   }
1517a5e57a09SHong Zhang   PetscMUMPS_c(&mumps->id);
15185cd7cf9dSHong Zhang   ierr = MatFactorSymbolic_MUMPS_ReportIfError(F,A,info,mumps);CHKERRQ(ierr);
15195cd7cf9dSHong Zhang 
15202792810eSHong Zhang   F->ops->choleskyfactornumeric = MatFactorNumeric_MUMPS;
1521dcd589f8SShri Abhyankar   F->ops->solve                 = MatSolve_MUMPS;
152251d5961aSHong Zhang   F->ops->solvetranspose        = MatSolve_MUMPS;
15234e34a73bSHong Zhang   F->ops->matsolve              = MatMatSolve_MUMPS;
15244e34a73bSHong Zhang #if defined(PETSC_USE_COMPLEX)
15250298fd71SBarry Smith   F->ops->getinertia = NULL;
15264e34a73bSHong Zhang #else
15274e34a73bSHong Zhang   F->ops->getinertia = MatGetInertia_SBAIJMUMPS;
1528db4efbfdSBarry Smith #endif
1529b24902e0SBarry Smith   PetscFunctionReturn(0);
1530b24902e0SBarry Smith }
1531b24902e0SBarry Smith 
153264e6c443SBarry Smith PetscErrorCode MatView_MUMPS(Mat A,PetscViewer viewer)
153374ed9c26SBarry Smith {
1534f6c57405SHong Zhang   PetscErrorCode    ierr;
153564e6c443SBarry Smith   PetscBool         iascii;
153664e6c443SBarry Smith   PetscViewerFormat format;
1537e69c285eSBarry Smith   Mat_MUMPS         *mumps=(Mat_MUMPS*)A->data;
1538f6c57405SHong Zhang 
1539f6c57405SHong Zhang   PetscFunctionBegin;
154064e6c443SBarry Smith   /* check if matrix is mumps type */
154164e6c443SBarry Smith   if (A->ops->solve != MatSolve_MUMPS) PetscFunctionReturn(0);
154264e6c443SBarry Smith 
1543251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
154464e6c443SBarry Smith   if (iascii) {
154564e6c443SBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
154664e6c443SBarry Smith     if (format == PETSC_VIEWER_ASCII_INFO) {
154764e6c443SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"MUMPS run parameters:\n");CHKERRQ(ierr);
1548a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  SYM (matrix type):                   %d \n",mumps->id.sym);CHKERRQ(ierr);
1549a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  PAR (host participation):            %d \n",mumps->id.par);CHKERRQ(ierr);
1550a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(1) (output for error):         %d \n",mumps->id.ICNTL(1));CHKERRQ(ierr);
1551a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(2) (output of diagnostic msg): %d \n",mumps->id.ICNTL(2));CHKERRQ(ierr);
1552a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(3) (output for global info):   %d \n",mumps->id.ICNTL(3));CHKERRQ(ierr);
1553a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(4) (level of printing):        %d \n",mumps->id.ICNTL(4));CHKERRQ(ierr);
1554a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(5) (input mat struct):         %d \n",mumps->id.ICNTL(5));CHKERRQ(ierr);
1555a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(6) (matrix prescaling):        %d \n",mumps->id.ICNTL(6));CHKERRQ(ierr);
1556d4ed72bbSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(7) (sequential matrix ordering):%d \n",mumps->id.ICNTL(7));CHKERRQ(ierr);
1557d4ed72bbSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(8) (scaling strategy):        %d \n",mumps->id.ICNTL(8));CHKERRQ(ierr);
1558a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(10) (max num of refinements):  %d \n",mumps->id.ICNTL(10));CHKERRQ(ierr);
1559a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(11) (error analysis):          %d \n",mumps->id.ICNTL(11));CHKERRQ(ierr);
1560a5e57a09SHong Zhang       if (mumps->id.ICNTL(11)>0) {
1561a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(4) (inf norm of input mat):        %g\n",mumps->id.RINFOG(4));CHKERRQ(ierr);
1562a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(5) (inf norm of solution):         %g\n",mumps->id.RINFOG(5));CHKERRQ(ierr);
1563a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(6) (inf norm of residual):         %g\n",mumps->id.RINFOG(6));CHKERRQ(ierr);
1564a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(7),RINFOG(8) (backward error est): %g, %g\n",mumps->id.RINFOG(7),mumps->id.RINFOG(8));CHKERRQ(ierr);
1565a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(9) (error estimate):               %g \n",mumps->id.RINFOG(9));CHKERRQ(ierr);
1566a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"    RINFOG(10),RINFOG(11)(condition numbers): %g, %g\n",mumps->id.RINFOG(10),mumps->id.RINFOG(11));CHKERRQ(ierr);
1567f6c57405SHong Zhang       }
1568a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(12) (efficiency control):                         %d \n",mumps->id.ICNTL(12));CHKERRQ(ierr);
1569a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(13) (efficiency control):                         %d \n",mumps->id.ICNTL(13));CHKERRQ(ierr);
1570a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(14) (percentage of estimated workspace increase): %d \n",mumps->id.ICNTL(14));CHKERRQ(ierr);
1571f6c57405SHong Zhang       /* ICNTL(15-17) not used */
1572a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(18) (input mat struct):                           %d \n",mumps->id.ICNTL(18));CHKERRQ(ierr);
1573d4ed72bbSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(19) (Schur complement info):                       %d \n",mumps->id.ICNTL(19));CHKERRQ(ierr);
1574a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(20) (rhs sparse pattern):                         %d \n",mumps->id.ICNTL(20));CHKERRQ(ierr);
1575ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(21) (solution struct):                            %d \n",mumps->id.ICNTL(21));CHKERRQ(ierr);
1576a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(22) (in-core/out-of-core facility):               %d \n",mumps->id.ICNTL(22));CHKERRQ(ierr);
1577a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(23) (max size of memory can be allocated locally):%d \n",mumps->id.ICNTL(23));CHKERRQ(ierr);
1578c0165424SHong Zhang 
1579a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(24) (detection of null pivot rows):               %d \n",mumps->id.ICNTL(24));CHKERRQ(ierr);
1580a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(25) (computation of a null space basis):          %d \n",mumps->id.ICNTL(25));CHKERRQ(ierr);
1581a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(26) (Schur options for rhs or solution):          %d \n",mumps->id.ICNTL(26));CHKERRQ(ierr);
1582a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(27) (experimental parameter):                     %d \n",mumps->id.ICNTL(27));CHKERRQ(ierr);
1583a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(28) (use parallel or sequential ordering):        %d \n",mumps->id.ICNTL(28));CHKERRQ(ierr);
1584a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(29) (parallel ordering):                          %d \n",mumps->id.ICNTL(29));CHKERRQ(ierr);
158542179a6aSHong Zhang 
1586a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(30) (user-specified set of entries in inv(A)):    %d \n",mumps->id.ICNTL(30));CHKERRQ(ierr);
1587a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(31) (factors is discarded in the solve phase):    %d \n",mumps->id.ICNTL(31));CHKERRQ(ierr);
1588a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(33) (compute determinant):                        %d \n",mumps->id.ICNTL(33));CHKERRQ(ierr);
15896e32de5dSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  ICNTL(35) (activate BLR based factorization):           %d \n",mumps->id.ICNTL(35));CHKERRQ(ierr);
1590f6c57405SHong Zhang 
1591a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(1) (relative pivoting threshold):      %g \n",mumps->id.CNTL(1));CHKERRQ(ierr);
1592a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(2) (stopping criterion of refinement): %g \n",mumps->id.CNTL(2));CHKERRQ(ierr);
1593ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(3) (absolute pivoting threshold):      %g \n",mumps->id.CNTL(3));CHKERRQ(ierr);
1594ca3dc52bSPierre Jolivet       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(4) (value of static pivoting):         %g \n",mumps->id.CNTL(4));CHKERRQ(ierr);
1595a5e57a09SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(5) (fixation for null pivots):         %g \n",mumps->id.CNTL(5));CHKERRQ(ierr);
15966e32de5dSHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  CNTL(7) (dropping parameter for BLR):       %g \n",mumps->id.CNTL(7));CHKERRQ(ierr);
1597f6c57405SHong Zhang 
1598f6c57405SHong Zhang       /* infomation local to each processor */
159934ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(1) (local estimated flops for the elimination after analysis): \n");CHKERRQ(ierr);
16001575c14dSBarry Smith       ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
1601a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %g \n",mumps->myid,mumps->id.RINFO(1));CHKERRQ(ierr);
16022a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
160334ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(2) (local estimated flops for the assembly after factorization): \n");CHKERRQ(ierr);
1604a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(2));CHKERRQ(ierr);
16052a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
160634ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  RINFO(3) (local estimated flops for the elimination after factorization): \n");CHKERRQ(ierr);
1607a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d]  %g \n",mumps->myid,mumps->id.RINFO(3));CHKERRQ(ierr);
16082a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1609f6c57405SHong Zhang 
161034ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(15) (estimated size of (in MB) MUMPS internal data for running numerical factorization): \n");CHKERRQ(ierr);
1611a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] %d \n",mumps->myid,mumps->id.INFO(15));CHKERRQ(ierr);
16122a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1613f6c57405SHong Zhang 
161434ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(16) (size of (in MB) MUMPS internal data used during numerical factorization): \n");CHKERRQ(ierr);
1615a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(16));CHKERRQ(ierr);
16162a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1617f6c57405SHong Zhang 
161834ed7027SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  INFO(23) (num of pivots eliminated on this processor after factorization): \n");CHKERRQ(ierr);
1619a5e57a09SHong Zhang       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(23));CHKERRQ(ierr);
16202a808120SBarry Smith       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1621b34f08ffSHong Zhang 
1622b34f08ffSHong Zhang       if (mumps->ninfo && mumps->ninfo <= 40){
1623b34f08ffSHong Zhang         PetscInt i;
1624b34f08ffSHong Zhang         for (i=0; i<mumps->ninfo; i++){
1625b34f08ffSHong Zhang           ierr = PetscViewerASCIIPrintf(viewer, "  INFO(%d): \n",mumps->info[i]);CHKERRQ(ierr);
1626b34f08ffSHong Zhang           ierr = PetscViewerASCIISynchronizedPrintf(viewer,"    [%d] %d \n",mumps->myid,mumps->id.INFO(mumps->info[i]));CHKERRQ(ierr);
16272a808120SBarry Smith           ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1628b34f08ffSHong Zhang         }
1629b34f08ffSHong Zhang       }
1630b34f08ffSHong Zhang 
1631b34f08ffSHong Zhang 
16321575c14dSBarry Smith       ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
1633f6c57405SHong Zhang 
1634a5e57a09SHong Zhang       if (!mumps->myid) { /* information from the host */
1635a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(1) (global estimated flops for the elimination after analysis): %g \n",mumps->id.RINFOG(1));CHKERRQ(ierr);
1636a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(2) (global estimated flops for the assembly after factorization): %g \n",mumps->id.RINFOG(2));CHKERRQ(ierr);
1637a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  RINFOG(3) (global estimated flops for the elimination after factorization): %g \n",mumps->id.RINFOG(3));CHKERRQ(ierr);
1638a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant): (%g,%g)*(2^%d)\n",mumps->id.RINFOG(12),mumps->id.RINFOG(13),mumps->id.INFOG(34));CHKERRQ(ierr);
1639f6c57405SHong Zhang 
1640a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(3) (estimated real workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(3));CHKERRQ(ierr);
1641a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(4) (estimated integer workspace for factors on all processors after analysis): %d \n",mumps->id.INFOG(4));CHKERRQ(ierr);
1642a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(5) (estimated maximum front size in the complete tree): %d \n",mumps->id.INFOG(5));CHKERRQ(ierr);
1643a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(6) (number of nodes in the complete tree): %d \n",mumps->id.INFOG(6));CHKERRQ(ierr);
1644a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(7) (ordering option effectively use after analysis): %d \n",mumps->id.INFOG(7));CHKERRQ(ierr);
1645a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(8) (structural symmetry in percent of the permuted matrix after analysis): %d \n",mumps->id.INFOG(8));CHKERRQ(ierr);
1646a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(9) (total real/complex workspace to store the matrix factors after factorization): %d \n",mumps->id.INFOG(9));CHKERRQ(ierr);
1647a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(10) (total integer space store the matrix factors after factorization): %d \n",mumps->id.INFOG(10));CHKERRQ(ierr);
1648a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(11) (order of largest frontal matrix after factorization): %d \n",mumps->id.INFOG(11));CHKERRQ(ierr);
1649a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(12) (number of off-diagonal pivots): %d \n",mumps->id.INFOG(12));CHKERRQ(ierr);
1650a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(13) (number of delayed pivots after factorization): %d \n",mumps->id.INFOG(13));CHKERRQ(ierr);
1651a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(14) (number of memory compress after factorization): %d \n",mumps->id.INFOG(14));CHKERRQ(ierr);
1652a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(15) (number of steps of iterative refinement after solution): %d \n",mumps->id.INFOG(15));CHKERRQ(ierr);
1653a5e57a09SHong 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",mumps->id.INFOG(16));CHKERRQ(ierr);
1654a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(17) (estimated size of all MUMPS internal data for factorization after analysis: sum over all processors): %d \n",mumps->id.INFOG(17));CHKERRQ(ierr);
1655a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(18) (size of all MUMPS internal data allocated during factorization: value on the most memory consuming processor): %d \n",mumps->id.INFOG(18));CHKERRQ(ierr);
1656a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(19) (size of all MUMPS internal data allocated during factorization: sum over all processors): %d \n",mumps->id.INFOG(19));CHKERRQ(ierr);
1657a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(20) (estimated number of entries in the factors): %d \n",mumps->id.INFOG(20));CHKERRQ(ierr);
1658a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(21) (size in MB of memory effectively used during factorization - value on the most memory consuming processor): %d \n",mumps->id.INFOG(21));CHKERRQ(ierr);
1659a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(22) (size in MB of memory effectively used during factorization - sum over all processors): %d \n",mumps->id.INFOG(22));CHKERRQ(ierr);
1660a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(23) (after analysis: value of ICNTL(6) effectively used): %d \n",mumps->id.INFOG(23));CHKERRQ(ierr);
1661a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(24) (after analysis: value of ICNTL(12) effectively used): %d \n",mumps->id.INFOG(24));CHKERRQ(ierr);
1662a5e57a09SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(25) (after factorization: number of pivots modified by static pivoting): %d \n",mumps->id.INFOG(25));CHKERRQ(ierr);
166340d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(28) (after factorization: number of null pivots encountered): %d\n",mumps->id.INFOG(28));CHKERRQ(ierr);
166440d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(29) (after factorization: effective number of entries in the factors (sum over all processors)): %d\n",mumps->id.INFOG(29));CHKERRQ(ierr);
166540d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(30, 31) (after solution: size in Mbytes of memory used during solution phase): %d, %d\n",mumps->id.INFOG(30),mumps->id.INFOG(31));CHKERRQ(ierr);
166640d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(32) (after analysis: type of analysis done): %d\n",mumps->id.INFOG(32));CHKERRQ(ierr);
166740d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(33) (value used for ICNTL(8)): %d\n",mumps->id.INFOG(33));CHKERRQ(ierr);
166840d435e3SHong Zhang         ierr = PetscViewerASCIIPrintf(viewer,"  INFOG(34) (exponent of the determinant if determinant is requested): %d\n",mumps->id.INFOG(34));CHKERRQ(ierr);
1669f6c57405SHong Zhang       }
1670f6c57405SHong Zhang     }
1671cb828f0fSHong Zhang   }
1672f6c57405SHong Zhang   PetscFunctionReturn(0);
1673f6c57405SHong Zhang }
1674f6c57405SHong Zhang 
167535bd34faSBarry Smith PetscErrorCode MatGetInfo_MUMPS(Mat A,MatInfoType flag,MatInfo *info)
167635bd34faSBarry Smith {
1677e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)A->data;
167835bd34faSBarry Smith 
167935bd34faSBarry Smith   PetscFunctionBegin;
168035bd34faSBarry Smith   info->block_size        = 1.0;
1681cb828f0fSHong Zhang   info->nz_allocated      = mumps->id.INFOG(20);
1682cb828f0fSHong Zhang   info->nz_used           = mumps->id.INFOG(20);
168335bd34faSBarry Smith   info->nz_unneeded       = 0.0;
168435bd34faSBarry Smith   info->assemblies        = 0.0;
168535bd34faSBarry Smith   info->mallocs           = 0.0;
168635bd34faSBarry Smith   info->memory            = 0.0;
168735bd34faSBarry Smith   info->fill_ratio_given  = 0;
168835bd34faSBarry Smith   info->fill_ratio_needed = 0;
168935bd34faSBarry Smith   info->factor_mallocs    = 0;
169035bd34faSBarry Smith   PetscFunctionReturn(0);
169135bd34faSBarry Smith }
169235bd34faSBarry Smith 
16935ccb76cbSHong Zhang /* -------------------------------------------------------------------------------------------*/
16948e7ba810SStefano Zampini PetscErrorCode MatFactorSetSchurIS_MUMPS(Mat F, IS is)
16956444a565SStefano Zampini {
1696e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
16978e7ba810SStefano Zampini   const PetscInt *idxs;
16988e7ba810SStefano Zampini   PetscInt       size,i;
16996444a565SStefano Zampini   PetscErrorCode ierr;
17006444a565SStefano Zampini 
17016444a565SStefano Zampini   PetscFunctionBegin;
1702b3cb21ddSStefano Zampini   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
1703241dbb5eSStefano Zampini   if (mumps->size > 1) {
1704241dbb5eSStefano Zampini     PetscBool ls,gs;
1705241dbb5eSStefano Zampini 
17064c644ebcSSatish Balay     ls   = mumps->myid ? (size ? PETSC_FALSE : PETSC_TRUE) : PETSC_TRUE;
1707241dbb5eSStefano Zampini     ierr = MPI_Allreduce(&ls,&gs,1,MPIU_BOOL,MPI_LAND,mumps->comm_mumps);CHKERRQ(ierr);
1708241dbb5eSStefano Zampini     if (!gs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MUMPS distributed parallel Schur complements not yet supported from PETSc\n");
1709241dbb5eSStefano Zampini   }
17106444a565SStefano Zampini   if (mumps->id.size_schur != size) {
17116444a565SStefano Zampini     ierr = PetscFree2(mumps->id.listvar_schur,mumps->id.schur);CHKERRQ(ierr);
17126444a565SStefano Zampini     mumps->id.size_schur = size;
17136444a565SStefano Zampini     mumps->id.schur_lld  = size;
17146444a565SStefano Zampini     ierr = PetscMalloc2(size,&mumps->id.listvar_schur,size*size,&mumps->id.schur);CHKERRQ(ierr);
17156444a565SStefano Zampini   }
1716b3cb21ddSStefano Zampini 
1717b3cb21ddSStefano Zampini   /* Schur complement matrix */
1718b3cb21ddSStefano Zampini   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mumps->id.size_schur,mumps->id.size_schur,(PetscScalar*)mumps->id.schur,&F->schur);CHKERRQ(ierr);
1719b3cb21ddSStefano Zampini   if (mumps->sym == 1) {
1720b3cb21ddSStefano Zampini     ierr = MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
1721b3cb21ddSStefano Zampini   }
1722b3cb21ddSStefano Zampini 
1723b3cb21ddSStefano Zampini   /* MUMPS expects Fortran style indices */
17248e7ba810SStefano Zampini   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
17256444a565SStefano Zampini   ierr = PetscMemcpy(mumps->id.listvar_schur,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
17268e7ba810SStefano Zampini   for (i=0;i<size;i++) mumps->id.listvar_schur[i]++;
17278e7ba810SStefano Zampini   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
1728241dbb5eSStefano Zampini   if (mumps->size > 1) {
1729241dbb5eSStefano Zampini     mumps->id.ICNTL(19) = 1; /* MUMPS returns Schur centralized on the host */
1730241dbb5eSStefano Zampini   } else {
17316444a565SStefano Zampini     if (F->factortype == MAT_FACTOR_LU) {
173259ac8732SStefano Zampini       mumps->id.ICNTL(19) = 3; /* MUMPS returns full matrix */
17336444a565SStefano Zampini     } else {
173459ac8732SStefano Zampini       mumps->id.ICNTL(19) = 2; /* MUMPS returns lower triangular part */
17356444a565SStefano Zampini     }
1736241dbb5eSStefano Zampini   }
173759ac8732SStefano Zampini   /* set a special value of ICNTL (not handled my MUMPS) to be used in the solve phase by PETSc */
1738b5fa320bSStefano Zampini   mumps->id.ICNTL(26) = -1;
17396444a565SStefano Zampini   PetscFunctionReturn(0);
17406444a565SStefano Zampini }
174159ac8732SStefano Zampini 
17426444a565SStefano Zampini /* -------------------------------------------------------------------------------------------*/
17435a05ddb0SStefano Zampini PetscErrorCode MatFactorCreateSchurComplement_MUMPS(Mat F,Mat* S)
17446444a565SStefano Zampini {
17456444a565SStefano Zampini   Mat            St;
1746e69c285eSBarry Smith   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
17476444a565SStefano Zampini   PetscScalar    *array;
17486444a565SStefano Zampini #if defined(PETSC_USE_COMPLEX)
17498ac429a0SStefano Zampini   PetscScalar    im = PetscSqrtScalar((PetscScalar)-1.0);
17506444a565SStefano Zampini #endif
17516444a565SStefano Zampini   PetscErrorCode ierr;
17526444a565SStefano Zampini 
17536444a565SStefano Zampini   PetscFunctionBegin;
17545a05ddb0SStefano Zampini   if (!mumps->id.ICNTL(19)) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
1755241dbb5eSStefano Zampini   ierr = MatCreate(PETSC_COMM_SELF,&St);CHKERRQ(ierr);
17566444a565SStefano Zampini   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mumps->id.size_schur,mumps->id.size_schur);CHKERRQ(ierr);
17576444a565SStefano Zampini   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
17586444a565SStefano Zampini   ierr = MatSetUp(St);CHKERRQ(ierr);
17596444a565SStefano Zampini   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
176059ac8732SStefano Zampini   if (!mumps->sym) { /* MUMPS always return a full matrix */
17616444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 1) { /* stored by rows */
17626444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17636444a565SStefano Zampini       for (i=0;i<N;i++) {
17646444a565SStefano Zampini         for (j=0;j<N;j++) {
17656444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17666444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17676444a565SStefano Zampini #else
17686444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
17696444a565SStefano Zampini #endif
17706444a565SStefano Zampini           array[j*N+i] = val;
17716444a565SStefano Zampini         }
17726444a565SStefano Zampini       }
17736444a565SStefano Zampini     } else { /* stored by columns */
17746444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
17756444a565SStefano Zampini     }
17766444a565SStefano Zampini   } else { /* either full or lower-triangular (not packed) */
17776444a565SStefano Zampini     if (mumps->id.ICNTL(19) == 2) { /* lower triangular stored by columns */
17786444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17796444a565SStefano Zampini       for (i=0;i<N;i++) {
17806444a565SStefano Zampini         for (j=i;j<N;j++) {
17816444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17826444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17836444a565SStefano Zampini #else
17846444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
17856444a565SStefano Zampini #endif
17866444a565SStefano Zampini           array[i*N+j] = val;
17876444a565SStefano Zampini           array[j*N+i] = val;
17886444a565SStefano Zampini         }
17896444a565SStefano Zampini       }
17906444a565SStefano Zampini     } else if (mumps->id.ICNTL(19) == 3) { /* full matrix */
17916444a565SStefano Zampini       ierr = PetscMemcpy(array,mumps->id.schur,mumps->id.size_schur*mumps->id.size_schur*sizeof(PetscScalar));CHKERRQ(ierr);
17926444a565SStefano Zampini     } else { /* ICNTL(19) == 1 lower triangular stored by rows */
17936444a565SStefano Zampini       PetscInt i,j,N=mumps->id.size_schur;
17946444a565SStefano Zampini       for (i=0;i<N;i++) {
17956444a565SStefano Zampini         for (j=0;j<i+1;j++) {
17966444a565SStefano Zampini #if !defined(PETSC_USE_COMPLEX)
17976444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j];
17986444a565SStefano Zampini #else
17996444a565SStefano Zampini           PetscScalar val = mumps->id.schur[i*N+j].r + im*mumps->id.schur[i*N+j].i;
18006444a565SStefano Zampini #endif
18016444a565SStefano Zampini           array[i*N+j] = val;
18026444a565SStefano Zampini           array[j*N+i] = val;
18036444a565SStefano Zampini         }
18046444a565SStefano Zampini       }
18056444a565SStefano Zampini     }
18066444a565SStefano Zampini   }
18076444a565SStefano Zampini   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
18086444a565SStefano Zampini   *S   = St;
18096444a565SStefano Zampini   PetscFunctionReturn(0);
18106444a565SStefano Zampini }
18116444a565SStefano Zampini 
181259ac8732SStefano Zampini /* -------------------------------------------------------------------------------------------*/
18135ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt ival)
18145ccb76cbSHong Zhang {
1815e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
18165ccb76cbSHong Zhang 
18175ccb76cbSHong Zhang   PetscFunctionBegin;
1818a5e57a09SHong Zhang   mumps->id.ICNTL(icntl) = ival;
18195ccb76cbSHong Zhang   PetscFunctionReturn(0);
18205ccb76cbSHong Zhang }
18215ccb76cbSHong Zhang 
1822bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl_MUMPS(Mat F,PetscInt icntl,PetscInt *ival)
1823bc6112feSHong Zhang {
1824e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1825bc6112feSHong Zhang 
1826bc6112feSHong Zhang   PetscFunctionBegin;
1827bc6112feSHong Zhang   *ival = mumps->id.ICNTL(icntl);
1828bc6112feSHong Zhang   PetscFunctionReturn(0);
1829bc6112feSHong Zhang }
1830bc6112feSHong Zhang 
18315ccb76cbSHong Zhang /*@
18325ccb76cbSHong Zhang   MatMumpsSetIcntl - Set MUMPS parameter ICNTL()
18335ccb76cbSHong Zhang 
18345ccb76cbSHong Zhang    Logically Collective on Mat
18355ccb76cbSHong Zhang 
18365ccb76cbSHong Zhang    Input Parameters:
18375ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
18385ccb76cbSHong Zhang .  icntl - index of MUMPS parameter array ICNTL()
18395ccb76cbSHong Zhang -  ival - value of MUMPS ICNTL(icntl)
18405ccb76cbSHong Zhang 
18415ccb76cbSHong Zhang   Options Database:
18425ccb76cbSHong Zhang .   -mat_mumps_icntl_<icntl> <ival>
18435ccb76cbSHong Zhang 
18445ccb76cbSHong Zhang    Level: beginner
18455ccb76cbSHong Zhang 
184696a0c994SBarry Smith    References:
184796a0c994SBarry Smith .     MUMPS Users' Guide
18485ccb76cbSHong Zhang 
18499fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
18505ccb76cbSHong Zhang  @*/
18515ccb76cbSHong Zhang PetscErrorCode MatMumpsSetIcntl(Mat F,PetscInt icntl,PetscInt ival)
18525ccb76cbSHong Zhang {
18535ccb76cbSHong Zhang   PetscErrorCode ierr;
18545ccb76cbSHong Zhang 
18555ccb76cbSHong Zhang   PetscFunctionBegin;
18562989dfd4SHong Zhang   PetscValidType(F,1);
18572989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
18585ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
18595ccb76cbSHong Zhang   PetscValidLogicalCollectiveInt(F,ival,3);
18605ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetIcntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
18615ccb76cbSHong Zhang   PetscFunctionReturn(0);
18625ccb76cbSHong Zhang }
18635ccb76cbSHong Zhang 
1864a21f80fcSHong Zhang /*@
1865a21f80fcSHong Zhang   MatMumpsGetIcntl - Get MUMPS parameter ICNTL()
1866a21f80fcSHong Zhang 
1867a21f80fcSHong Zhang    Logically Collective on Mat
1868a21f80fcSHong Zhang 
1869a21f80fcSHong Zhang    Input Parameters:
1870a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1871a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array ICNTL()
1872a21f80fcSHong Zhang 
1873a21f80fcSHong Zhang   Output Parameter:
1874a21f80fcSHong Zhang .  ival - value of MUMPS ICNTL(icntl)
1875a21f80fcSHong Zhang 
1876a21f80fcSHong Zhang    Level: beginner
1877a21f80fcSHong Zhang 
187896a0c994SBarry Smith    References:
187996a0c994SBarry Smith .     MUMPS Users' Guide
1880a21f80fcSHong Zhang 
18819fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1882a21f80fcSHong Zhang @*/
1883bc6112feSHong Zhang PetscErrorCode MatMumpsGetIcntl(Mat F,PetscInt icntl,PetscInt *ival)
1884bc6112feSHong Zhang {
1885bc6112feSHong Zhang   PetscErrorCode ierr;
1886bc6112feSHong Zhang 
1887bc6112feSHong Zhang   PetscFunctionBegin;
18882989dfd4SHong Zhang   PetscValidType(F,1);
18892989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1890bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1891bc6112feSHong Zhang   PetscValidIntPointer(ival,3);
18922989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetIcntl_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
1893bc6112feSHong Zhang   PetscFunctionReturn(0);
1894bc6112feSHong Zhang }
1895bc6112feSHong Zhang 
18968928b65cSHong Zhang /* -------------------------------------------------------------------------------------------*/
18978928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal val)
18988928b65cSHong Zhang {
1899e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
19008928b65cSHong Zhang 
19018928b65cSHong Zhang   PetscFunctionBegin;
19028928b65cSHong Zhang   mumps->id.CNTL(icntl) = val;
19038928b65cSHong Zhang   PetscFunctionReturn(0);
19048928b65cSHong Zhang }
19058928b65cSHong Zhang 
1906bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl_MUMPS(Mat F,PetscInt icntl,PetscReal *val)
1907bc6112feSHong Zhang {
1908e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1909bc6112feSHong Zhang 
1910bc6112feSHong Zhang   PetscFunctionBegin;
1911bc6112feSHong Zhang   *val = mumps->id.CNTL(icntl);
1912bc6112feSHong Zhang   PetscFunctionReturn(0);
1913bc6112feSHong Zhang }
1914bc6112feSHong Zhang 
19158928b65cSHong Zhang /*@
19168928b65cSHong Zhang   MatMumpsSetCntl - Set MUMPS parameter CNTL()
19178928b65cSHong Zhang 
19188928b65cSHong Zhang    Logically Collective on Mat
19198928b65cSHong Zhang 
19208928b65cSHong Zhang    Input Parameters:
19218928b65cSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
19228928b65cSHong Zhang .  icntl - index of MUMPS parameter array CNTL()
19238928b65cSHong Zhang -  val - value of MUMPS CNTL(icntl)
19248928b65cSHong Zhang 
19258928b65cSHong Zhang   Options Database:
19268928b65cSHong Zhang .   -mat_mumps_cntl_<icntl> <val>
19278928b65cSHong Zhang 
19288928b65cSHong Zhang    Level: beginner
19298928b65cSHong Zhang 
193096a0c994SBarry Smith    References:
193196a0c994SBarry Smith .     MUMPS Users' Guide
19328928b65cSHong Zhang 
19339fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
19348928b65cSHong Zhang @*/
19358928b65cSHong Zhang PetscErrorCode MatMumpsSetCntl(Mat F,PetscInt icntl,PetscReal val)
19368928b65cSHong Zhang {
19378928b65cSHong Zhang   PetscErrorCode ierr;
19388928b65cSHong Zhang 
19398928b65cSHong Zhang   PetscFunctionBegin;
19402989dfd4SHong Zhang   PetscValidType(F,1);
19412989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
19428928b65cSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1943bc6112feSHong Zhang   PetscValidLogicalCollectiveReal(F,val,3);
19448928b65cSHong Zhang   ierr = PetscTryMethod(F,"MatMumpsSetCntl_C",(Mat,PetscInt,PetscReal),(F,icntl,val));CHKERRQ(ierr);
19458928b65cSHong Zhang   PetscFunctionReturn(0);
19468928b65cSHong Zhang }
19478928b65cSHong Zhang 
1948a21f80fcSHong Zhang /*@
1949a21f80fcSHong Zhang   MatMumpsGetCntl - Get MUMPS parameter CNTL()
1950a21f80fcSHong Zhang 
1951a21f80fcSHong Zhang    Logically Collective on Mat
1952a21f80fcSHong Zhang 
1953a21f80fcSHong Zhang    Input Parameters:
1954a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
1955a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array CNTL()
1956a21f80fcSHong Zhang 
1957a21f80fcSHong Zhang   Output Parameter:
1958a21f80fcSHong Zhang .  val - value of MUMPS CNTL(icntl)
1959a21f80fcSHong Zhang 
1960a21f80fcSHong Zhang    Level: beginner
1961a21f80fcSHong Zhang 
196296a0c994SBarry Smith    References:
196396a0c994SBarry Smith .      MUMPS Users' Guide
1964a21f80fcSHong Zhang 
19659fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
1966a21f80fcSHong Zhang @*/
1967bc6112feSHong Zhang PetscErrorCode MatMumpsGetCntl(Mat F,PetscInt icntl,PetscReal *val)
1968bc6112feSHong Zhang {
1969bc6112feSHong Zhang   PetscErrorCode ierr;
1970bc6112feSHong Zhang 
1971bc6112feSHong Zhang   PetscFunctionBegin;
19722989dfd4SHong Zhang   PetscValidType(F,1);
19732989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
1974bc6112feSHong Zhang   PetscValidLogicalCollectiveInt(F,icntl,2);
1975bc6112feSHong Zhang   PetscValidRealPointer(val,3);
19762989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetCntl_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
1977bc6112feSHong Zhang   PetscFunctionReturn(0);
1978bc6112feSHong Zhang }
1979bc6112feSHong Zhang 
1980ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo_MUMPS(Mat F,PetscInt icntl,PetscInt *info)
1981bc6112feSHong Zhang {
1982e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1983bc6112feSHong Zhang 
1984bc6112feSHong Zhang   PetscFunctionBegin;
1985bc6112feSHong Zhang   *info = mumps->id.INFO(icntl);
1986bc6112feSHong Zhang   PetscFunctionReturn(0);
1987bc6112feSHong Zhang }
1988bc6112feSHong Zhang 
1989ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog_MUMPS(Mat F,PetscInt icntl,PetscInt *infog)
1990bc6112feSHong Zhang {
1991e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
1992bc6112feSHong Zhang 
1993bc6112feSHong Zhang   PetscFunctionBegin;
1994bc6112feSHong Zhang   *infog = mumps->id.INFOG(icntl);
1995bc6112feSHong Zhang   PetscFunctionReturn(0);
1996bc6112feSHong Zhang }
1997bc6112feSHong Zhang 
1998ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfo)
1999bc6112feSHong Zhang {
2000e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2001bc6112feSHong Zhang 
2002bc6112feSHong Zhang   PetscFunctionBegin;
2003bc6112feSHong Zhang   *rinfo = mumps->id.RINFO(icntl);
2004bc6112feSHong Zhang   PetscFunctionReturn(0);
2005bc6112feSHong Zhang }
2006bc6112feSHong Zhang 
2007ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog_MUMPS(Mat F,PetscInt icntl,PetscReal *rinfog)
2008bc6112feSHong Zhang {
2009e69c285eSBarry Smith   Mat_MUMPS *mumps =(Mat_MUMPS*)F->data;
2010bc6112feSHong Zhang 
2011bc6112feSHong Zhang   PetscFunctionBegin;
2012bc6112feSHong Zhang   *rinfog = mumps->id.RINFOG(icntl);
2013bc6112feSHong Zhang   PetscFunctionReturn(0);
2014bc6112feSHong Zhang }
2015bc6112feSHong Zhang 
201689a9c03aSHong Zhang PetscErrorCode MatMumpsGetInverse_MUMPS(Mat F,Mat spRHS)
2017bb599dfdSHong Zhang {
2018bb599dfdSHong Zhang   PetscErrorCode ierr;
2019bb599dfdSHong Zhang   Mat            Bt = NULL;
2020bb599dfdSHong Zhang   PetscBool      flgT;
2021bb599dfdSHong Zhang   Mat_MUMPS      *mumps =(Mat_MUMPS*)F->data;
2022bb599dfdSHong Zhang   PetscBool      done;
2023bb599dfdSHong Zhang   PetscScalar    *aa;
2024bb599dfdSHong Zhang   PetscInt       spnr,*ia,*ja;
2025bb599dfdSHong Zhang 
2026bb599dfdSHong Zhang   PetscFunctionBegin;
2027e3f2db6aSHong Zhang   if (!mumps->myid) {
2028e3f2db6aSHong Zhang     PetscValidIntPointer(spRHS,2);
2029bb599dfdSHong Zhang     ierr = PetscObjectTypeCompare((PetscObject)spRHS,MATTRANSPOSEMAT,&flgT);CHKERRQ(ierr);
2030bb599dfdSHong Zhang     if (flgT) {
2031bb599dfdSHong Zhang       ierr = MatTransposeGetMat(spRHS,&Bt);CHKERRQ(ierr);
2032bb599dfdSHong Zhang     } else {
2033bb599dfdSHong Zhang       SETERRQ(PetscObjectComm((PetscObject)spRHS),PETSC_ERR_ARG_WRONG,"Matrix spRHS must be type MATTRANSPOSEMAT matrix");
2034bb599dfdSHong Zhang     }
2035e3f2db6aSHong Zhang   }
2036bb599dfdSHong Zhang 
2037bb599dfdSHong Zhang   ierr = MatMumpsSetIcntl(F,30,1);CHKERRQ(ierr);
2038bb599dfdSHong Zhang 
2039e3f2db6aSHong Zhang   if (!mumps->myid) {
2040bb599dfdSHong Zhang     ierr = MatSeqAIJGetArray(Bt,&aa);CHKERRQ(ierr);
2041bb599dfdSHong Zhang     ierr = MatGetRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
2042bb599dfdSHong Zhang     if (!done) SETERRQ(PetscObjectComm((PetscObject)Bt),PETSC_ERR_ARG_WRONG,"Cannot get IJ structure");
2043bb599dfdSHong Zhang 
2044bb599dfdSHong Zhang     mumps->id.irhs_ptr    = ia;
2045bb599dfdSHong Zhang     mumps->id.irhs_sparse = ja;
2046bb599dfdSHong Zhang     mumps->id.nz_rhs      = ia[spnr] - 1;
2047bb599dfdSHong Zhang     mumps->id.rhs_sparse  = (MumpsScalar*)aa;
2048e3f2db6aSHong Zhang   } else {
2049e3f2db6aSHong Zhang     mumps->id.irhs_ptr    = NULL;
2050e3f2db6aSHong Zhang     mumps->id.irhs_sparse = NULL;
2051e3f2db6aSHong Zhang     mumps->id.nz_rhs      = 0;
2052e3f2db6aSHong Zhang     mumps->id.rhs_sparse  = NULL;
2053e3f2db6aSHong Zhang   }
2054bb599dfdSHong Zhang   mumps->id.ICNTL(20)   = 1; /* rhs is sparse */
2055e3f2db6aSHong Zhang   mumps->id.ICNTL(21)   = 0; /* solution is in assembled centralized format */
2056bb599dfdSHong Zhang 
2057bb599dfdSHong Zhang   /* solve phase */
2058bb599dfdSHong Zhang   /*-------------*/
2059bb599dfdSHong Zhang   mumps->id.job = JOB_SOLVE;
2060bb599dfdSHong Zhang   PetscMUMPS_c(&mumps->id);
2061e3f2db6aSHong Zhang   if (mumps->id.INFOG(1) < 0)
2062e3f2db6aSHong Zhang     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MUMPS in solve phase: INFOG(1)=%d INFO(2)=%d\n",mumps->id.INFOG(1),mumps->id.INFO(2));
206314267174SHong Zhang 
2064e3f2db6aSHong Zhang   if (!mumps->myid) {
206514267174SHong Zhang     ierr = MatSeqAIJRestoreArray(Bt,&aa);CHKERRQ(ierr);
206614267174SHong Zhang     ierr = MatRestoreRowIJ(Bt,1,PETSC_FALSE,PETSC_FALSE,&spnr,(const PetscInt**)&ia,(const PetscInt**)&ja,&done);CHKERRQ(ierr);
2067e3f2db6aSHong Zhang   }
2068bb599dfdSHong Zhang   PetscFunctionReturn(0);
2069bb599dfdSHong Zhang }
2070bb599dfdSHong Zhang 
2071bb599dfdSHong Zhang /*@
207289a9c03aSHong Zhang   MatMumpsGetInverse - Get user-specified set of entries in inverse of A
2073bb599dfdSHong Zhang 
2074bb599dfdSHong Zhang    Logically Collective on Mat
2075bb599dfdSHong Zhang 
2076bb599dfdSHong Zhang    Input Parameters:
2077bb599dfdSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2078e3f2db6aSHong Zhang -  spRHS - sequential sparse matrix in MATTRANSPOSEMAT format holding specified indices in processor[0]
2079bb599dfdSHong Zhang 
2080bb599dfdSHong Zhang   Output Parameter:
2081e3f2db6aSHong Zhang . spRHS - requested entries of inverse of A
2082bb599dfdSHong Zhang 
2083bb599dfdSHong Zhang    Level: beginner
2084bb599dfdSHong Zhang 
2085bb599dfdSHong Zhang    References:
2086bb599dfdSHong Zhang .      MUMPS Users' Guide
2087bb599dfdSHong Zhang 
2088bb599dfdSHong Zhang .seealso: MatGetFactor(), MatCreateTranspose()
2089bb599dfdSHong Zhang @*/
209089a9c03aSHong Zhang PetscErrorCode MatMumpsGetInverse(Mat F,Mat spRHS)
2091bb599dfdSHong Zhang {
2092bb599dfdSHong Zhang   PetscErrorCode ierr;
2093bb599dfdSHong Zhang 
2094bb599dfdSHong Zhang   PetscFunctionBegin;
2095bb599dfdSHong Zhang   PetscValidType(F,1);
2096bb599dfdSHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
209789a9c03aSHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetInverse_C",(Mat,Mat),(F,spRHS));CHKERRQ(ierr);
2098bb599dfdSHong Zhang   PetscFunctionReturn(0);
2099bb599dfdSHong Zhang }
2100bb599dfdSHong Zhang 
2101a21f80fcSHong Zhang /*@
2102a21f80fcSHong Zhang   MatMumpsGetInfo - Get MUMPS parameter INFO()
2103a21f80fcSHong Zhang 
2104a21f80fcSHong Zhang    Logically Collective on Mat
2105a21f80fcSHong Zhang 
2106a21f80fcSHong Zhang    Input Parameters:
2107a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2108a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFO()
2109a21f80fcSHong Zhang 
2110a21f80fcSHong Zhang   Output Parameter:
2111a21f80fcSHong Zhang .  ival - value of MUMPS INFO(icntl)
2112a21f80fcSHong Zhang 
2113a21f80fcSHong Zhang    Level: beginner
2114a21f80fcSHong Zhang 
211596a0c994SBarry Smith    References:
211696a0c994SBarry Smith .      MUMPS Users' Guide
2117a21f80fcSHong Zhang 
21189fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2119a21f80fcSHong Zhang @*/
2120ca810319SHong Zhang PetscErrorCode MatMumpsGetInfo(Mat F,PetscInt icntl,PetscInt *ival)
2121bc6112feSHong Zhang {
2122bc6112feSHong Zhang   PetscErrorCode ierr;
2123bc6112feSHong Zhang 
2124bc6112feSHong Zhang   PetscFunctionBegin;
21252989dfd4SHong Zhang   PetscValidType(F,1);
21262989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2127ca810319SHong Zhang   PetscValidIntPointer(ival,3);
21282989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetInfo_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2129bc6112feSHong Zhang   PetscFunctionReturn(0);
2130bc6112feSHong Zhang }
2131bc6112feSHong Zhang 
2132a21f80fcSHong Zhang /*@
2133a21f80fcSHong Zhang   MatMumpsGetInfog - Get MUMPS parameter INFOG()
2134a21f80fcSHong Zhang 
2135a21f80fcSHong Zhang    Logically Collective on Mat
2136a21f80fcSHong Zhang 
2137a21f80fcSHong Zhang    Input Parameters:
2138a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2139a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array INFOG()
2140a21f80fcSHong Zhang 
2141a21f80fcSHong Zhang   Output Parameter:
2142a21f80fcSHong Zhang .  ival - value of MUMPS INFOG(icntl)
2143a21f80fcSHong Zhang 
2144a21f80fcSHong Zhang    Level: beginner
2145a21f80fcSHong Zhang 
214696a0c994SBarry Smith    References:
214796a0c994SBarry Smith .      MUMPS Users' Guide
2148a21f80fcSHong Zhang 
21499fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2150a21f80fcSHong Zhang @*/
2151ca810319SHong Zhang PetscErrorCode MatMumpsGetInfog(Mat F,PetscInt icntl,PetscInt *ival)
2152bc6112feSHong Zhang {
2153bc6112feSHong Zhang   PetscErrorCode ierr;
2154bc6112feSHong Zhang 
2155bc6112feSHong Zhang   PetscFunctionBegin;
21562989dfd4SHong Zhang   PetscValidType(F,1);
21572989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2158ca810319SHong Zhang   PetscValidIntPointer(ival,3);
21592989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetInfog_C",(Mat,PetscInt,PetscInt*),(F,icntl,ival));CHKERRQ(ierr);
2160bc6112feSHong Zhang   PetscFunctionReturn(0);
2161bc6112feSHong Zhang }
2162bc6112feSHong Zhang 
2163a21f80fcSHong Zhang /*@
2164a21f80fcSHong Zhang   MatMumpsGetRinfo - Get MUMPS parameter RINFO()
2165a21f80fcSHong Zhang 
2166a21f80fcSHong Zhang    Logically Collective on Mat
2167a21f80fcSHong Zhang 
2168a21f80fcSHong Zhang    Input Parameters:
2169a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2170a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFO()
2171a21f80fcSHong Zhang 
2172a21f80fcSHong Zhang   Output Parameter:
2173a21f80fcSHong Zhang .  val - value of MUMPS RINFO(icntl)
2174a21f80fcSHong Zhang 
2175a21f80fcSHong Zhang    Level: beginner
2176a21f80fcSHong Zhang 
217796a0c994SBarry Smith    References:
217896a0c994SBarry Smith .       MUMPS Users' Guide
2179a21f80fcSHong Zhang 
21809fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2181a21f80fcSHong Zhang @*/
2182ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfo(Mat F,PetscInt icntl,PetscReal *val)
2183bc6112feSHong Zhang {
2184bc6112feSHong Zhang   PetscErrorCode ierr;
2185bc6112feSHong Zhang 
2186bc6112feSHong Zhang   PetscFunctionBegin;
21872989dfd4SHong Zhang   PetscValidType(F,1);
21882989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2189bc6112feSHong Zhang   PetscValidRealPointer(val,3);
21902989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetRinfo_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2191bc6112feSHong Zhang   PetscFunctionReturn(0);
2192bc6112feSHong Zhang }
2193bc6112feSHong Zhang 
2194a21f80fcSHong Zhang /*@
2195a21f80fcSHong Zhang   MatMumpsGetRinfog - Get MUMPS parameter RINFOG()
2196a21f80fcSHong Zhang 
2197a21f80fcSHong Zhang    Logically Collective on Mat
2198a21f80fcSHong Zhang 
2199a21f80fcSHong Zhang    Input Parameters:
2200a21f80fcSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-MUMPS interface
2201a21f80fcSHong Zhang -  icntl - index of MUMPS parameter array RINFOG()
2202a21f80fcSHong Zhang 
2203a21f80fcSHong Zhang   Output Parameter:
2204a21f80fcSHong Zhang .  val - value of MUMPS RINFOG(icntl)
2205a21f80fcSHong Zhang 
2206a21f80fcSHong Zhang    Level: beginner
2207a21f80fcSHong Zhang 
220896a0c994SBarry Smith    References:
220996a0c994SBarry Smith .      MUMPS Users' Guide
2210a21f80fcSHong Zhang 
22119fc87aa7SBarry Smith .seealso: MatGetFactor(), MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog()
2212a21f80fcSHong Zhang @*/
2213ca810319SHong Zhang PetscErrorCode MatMumpsGetRinfog(Mat F,PetscInt icntl,PetscReal *val)
2214bc6112feSHong Zhang {
2215bc6112feSHong Zhang   PetscErrorCode ierr;
2216bc6112feSHong Zhang 
2217bc6112feSHong Zhang   PetscFunctionBegin;
22182989dfd4SHong Zhang   PetscValidType(F,1);
22192989dfd4SHong Zhang   if (!F->factortype) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ARG_WRONGSTATE,"Only for factored matrix");
2220bc6112feSHong Zhang   PetscValidRealPointer(val,3);
22212989dfd4SHong Zhang   ierr = PetscUseMethod(F,"MatMumpsGetRinfog_C",(Mat,PetscInt,PetscReal*),(F,icntl,val));CHKERRQ(ierr);
2222bc6112feSHong Zhang   PetscFunctionReturn(0);
2223bc6112feSHong Zhang }
2224bc6112feSHong Zhang 
222524b6179bSKris Buschelman /*MC
22262692d6eeSBarry Smith   MATSOLVERMUMPS -  A matrix type providing direct solvers (LU and Cholesky) for
222724b6179bSKris Buschelman   distributed and sequential matrices via the external package MUMPS.
222824b6179bSKris Buschelman 
222941c8de11SBarry Smith   Works with MATAIJ and MATSBAIJ matrices
223024b6179bSKris Buschelman 
2231c2b89b5dSBarry Smith   Use ./configure --download-mumps --download-scalapack --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with MUMPS
2232c2b89b5dSBarry Smith 
22333ca39a21SBarry Smith   Use -pc_type cholesky or lu -pc_factor_mat_solver_type mumps to use this direct solver
2234c2b89b5dSBarry Smith 
223524b6179bSKris Buschelman   Options Database Keys:
22364422a9fcSPatrick Sanan +  -mat_mumps_icntl_1 - ICNTL(1): output stream for error messages
22374422a9fcSPatrick Sanan .  -mat_mumps_icntl_2 - ICNTL(2): output stream for diagnostic printing, statistics, and warning
22384422a9fcSPatrick Sanan .  -mat_mumps_icntl_3 -  ICNTL(3): output stream for global information, collected on the host
22394422a9fcSPatrick Sanan .  -mat_mumps_icntl_4 -  ICNTL(4): level of printing (0 to 4)
22404422a9fcSPatrick Sanan .  -mat_mumps_icntl_6 - ICNTL(6): permutes to a zero-free diagonal and/or scale the matrix (0 to 7)
22414422a9fcSPatrick Sanan .  -mat_mumps_icntl_7 - ICNTL(7): computes a symmetric permutation in sequential analysis (0 to 7). 3=Scotch, 4=PORD, 5=Metis
22424422a9fcSPatrick Sanan .  -mat_mumps_icntl_8  - ICNTL(8): scaling strategy (-2 to 8 or 77)
22434422a9fcSPatrick Sanan .  -mat_mumps_icntl_10  - ICNTL(10): max num of refinements
22444422a9fcSPatrick Sanan .  -mat_mumps_icntl_11  - ICNTL(11): statistics related to an error analysis (via -ksp_view)
22454422a9fcSPatrick Sanan .  -mat_mumps_icntl_12  - ICNTL(12): an ordering strategy for symmetric matrices (0 to 3)
22464422a9fcSPatrick Sanan .  -mat_mumps_icntl_13  - ICNTL(13): parallelism of the root node (enable ScaLAPACK) and its splitting
22474422a9fcSPatrick Sanan .  -mat_mumps_icntl_14  - ICNTL(14): percentage increase in the estimated working space
22484422a9fcSPatrick Sanan .  -mat_mumps_icntl_19  - ICNTL(19): computes the Schur complement
22494422a9fcSPatrick Sanan .  -mat_mumps_icntl_22  - ICNTL(22): in-core/out-of-core factorization and solve (0 or 1)
22504422a9fcSPatrick Sanan .  -mat_mumps_icntl_23  - ICNTL(23): max size of the working memory (MB) that can allocate per processor
22514422a9fcSPatrick Sanan .  -mat_mumps_icntl_24  - ICNTL(24): detection of null pivot rows (0 or 1)
22524422a9fcSPatrick Sanan .  -mat_mumps_icntl_25  - ICNTL(25): compute a solution of a deficient matrix and a null space basis
22534422a9fcSPatrick Sanan .  -mat_mumps_icntl_26  - ICNTL(26): drives the solution phase if a Schur complement matrix
22544422a9fcSPatrick Sanan .  -mat_mumps_icntl_28  - ICNTL(28): use 1 for sequential analysis and ictnl(7) ordering, or 2 for parallel analysis and ictnl(29) ordering
22554422a9fcSPatrick Sanan .  -mat_mumps_icntl_29 - ICNTL(29): parallel ordering 1 = ptscotch, 2 = parmetis
22564422a9fcSPatrick Sanan .  -mat_mumps_icntl_30 - ICNTL(30): compute user-specified set of entries in inv(A)
22574422a9fcSPatrick Sanan .  -mat_mumps_icntl_31 - ICNTL(31): indicates which factors may be discarded during factorization
22584422a9fcSPatrick Sanan .  -mat_mumps_icntl_33 - ICNTL(33): compute determinant
22594422a9fcSPatrick Sanan .  -mat_mumps_cntl_1  - CNTL(1): relative pivoting threshold
22604422a9fcSPatrick Sanan .  -mat_mumps_cntl_2  -  CNTL(2): stopping criterion of refinement
22614422a9fcSPatrick Sanan .  -mat_mumps_cntl_3 - CNTL(3): absolute pivoting threshold
22624422a9fcSPatrick Sanan .  -mat_mumps_cntl_4 - CNTL(4): value for static pivoting
22634422a9fcSPatrick Sanan -  -mat_mumps_cntl_5 - CNTL(5): fixation for null pivots
226424b6179bSKris Buschelman 
226524b6179bSKris Buschelman   Level: beginner
226624b6179bSKris Buschelman 
226795452b02SPatrick Sanan     Notes:
226895452b02SPatrick Sanan     When a MUMPS factorization fails inside a KSP solve, for example with a KSP_DIVERGED_PCSETUP_FAILED, one can find the MUMPS information about the failure by calling
22699fc87aa7SBarry Smith $          KSPGetPC(ksp,&pc);
22709fc87aa7SBarry Smith $          PCFactorGetMatrix(pc,&mat);
22719fc87aa7SBarry Smith $          MatMumpsGetInfo(mat,....);
22729fc87aa7SBarry Smith $          MatMumpsGetInfog(mat,....); etc.
22739fc87aa7SBarry Smith            Or you can run with -ksp_error_if_not_converged and the program will be stopped and the information printed in the error message.
22749fc87aa7SBarry Smith 
22753ca39a21SBarry Smith .seealso: PCFactorSetMatSolverType(), MatSolverType, MatMumpsSetICntl(), MatMumpsGetIcntl(), MatMumpsSetCntl(), MatMumpsGetCntl(), MatMumpsGetInfo(), MatMumpsGetInfog(), MatMumpsGetRinfo(), MatMumpsGetRinfog(), KSPGetPC(), PCGetFactor(), PCFactorGetMatrix()
227641c8de11SBarry Smith 
227724b6179bSKris Buschelman M*/
227824b6179bSKris Buschelman 
2279ea799195SBarry Smith static PetscErrorCode MatFactorGetSolverType_mumps(Mat A,MatSolverType *type)
228035bd34faSBarry Smith {
228135bd34faSBarry Smith   PetscFunctionBegin;
22822692d6eeSBarry Smith   *type = MATSOLVERMUMPS;
228335bd34faSBarry Smith   PetscFunctionReturn(0);
228435bd34faSBarry Smith }
228535bd34faSBarry Smith 
2286bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI AIJ matrices */
2287cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_aij_mumps(Mat A,MatFactorType ftype,Mat *F)
22882877fffaSHong Zhang {
22892877fffaSHong Zhang   Mat            B;
22902877fffaSHong Zhang   PetscErrorCode ierr;
22912877fffaSHong Zhang   Mat_MUMPS      *mumps;
2292ace3abfcSBarry Smith   PetscBool      isSeqAIJ;
22932877fffaSHong Zhang 
22942877fffaSHong Zhang   PetscFunctionBegin;
22952877fffaSHong Zhang   /* Create the factorization matrix */
2296251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
2297ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
22982877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2299e69c285eSBarry Smith   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2300e69c285eSBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
23012877fffaSHong Zhang 
2302b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
23032205254eSKarl Rupp 
23042877fffaSHong Zhang   B->ops->view        = MatView_MUMPS;
230535bd34faSBarry Smith   B->ops->getinfo     = MatGetInfo_MUMPS;
23062205254eSKarl Rupp 
23073ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
23085a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
23095a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2310bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2311bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2312bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2313bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2314ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2315ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2316ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2317ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
231889a9c03aSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);CHKERRQ(ierr);
23196444a565SStefano Zampini 
2320450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2321450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
2322d5f3da31SBarry Smith     B->factortype            = MAT_FACTOR_LU;
2323bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqaij;
2324bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij;
2325746480a1SHong Zhang     mumps->sym = 0;
2326dcd589f8SShri Abhyankar   } else {
232767877ebaSShri Abhyankar     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2328450b117fSShri Abhyankar     B->factortype                  = MAT_FACTOR_CHOLESKY;
2329bccb9932SShri Abhyankar     if (isSeqAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqaij_seqsbaij;
2330bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpiaij_mpisbaij;
233159ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
233259ac8732SStefano Zampini     mumps->sym = 2;
233359ac8732SStefano Zampini #else
23346fdc2a6dSBarry Smith     if (A->spd_set && A->spd) mumps->sym = 1;
23356fdc2a6dSBarry Smith     else                      mumps->sym = 2;
233659ac8732SStefano Zampini #endif
2337450b117fSShri Abhyankar   }
23382877fffaSHong Zhang 
233900c67f3bSHong Zhang   /* set solvertype */
234000c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
234100c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
234200c67f3bSHong Zhang 
23432877fffaSHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
2344e69c285eSBarry Smith   B->data         = (void*)mumps;
23452205254eSKarl Rupp 
2346f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2347746480a1SHong Zhang 
23482877fffaSHong Zhang   *F = B;
23492877fffaSHong Zhang   PetscFunctionReturn(0);
23502877fffaSHong Zhang }
23512877fffaSHong Zhang 
2352bccb9932SShri Abhyankar /* MatGetFactor for Seq and MPI SBAIJ matrices */
2353cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_sbaij_mumps(Mat A,MatFactorType ftype,Mat *F)
23542877fffaSHong Zhang {
23552877fffaSHong Zhang   Mat            B;
23562877fffaSHong Zhang   PetscErrorCode ierr;
23572877fffaSHong Zhang   Mat_MUMPS      *mumps;
2358ace3abfcSBarry Smith   PetscBool      isSeqSBAIJ;
23592877fffaSHong Zhang 
23602877fffaSHong Zhang   PetscFunctionBegin;
2361ce94432eSBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with MUMPS LU, use AIJ matrix");
2362ce94432eSBarry Smith   if (A->rmap->bs > 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with block size > 1 with MUMPS Cholesky, use AIJ matrix instead");
2363251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
23642877fffaSHong Zhang   /* Create the factorization matrix */
2365ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
23662877fffaSHong Zhang   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2367e69c285eSBarry Smith   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2368e69c285eSBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
2369e69c285eSBarry Smith 
2370b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2371bccb9932SShri Abhyankar   if (isSeqSBAIJ) {
237216ebf90aSShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_seqsbaij_seqsbaij;
2373dcd589f8SShri Abhyankar   } else {
2374bccb9932SShri Abhyankar     mumps->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij;
2375bccb9932SShri Abhyankar   }
2376bccb9932SShri Abhyankar 
2377e69c285eSBarry Smith   B->ops->getinfo                = MatGetInfo_External;
237867877ebaSShri Abhyankar   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_MUMPS;
2379bccb9932SShri Abhyankar   B->ops->view                   = MatView_MUMPS;
23802205254eSKarl Rupp 
23813ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
23825a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
23835a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2384b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2385b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2386b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2387b13644aeSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2388ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2389ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2390ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2391ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
239289a9c03aSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);CHKERRQ(ierr);
23932205254eSKarl Rupp 
2394f4762488SHong Zhang   B->factortype = MAT_FACTOR_CHOLESKY;
239559ac8732SStefano Zampini #if defined(PETSC_USE_COMPLEX)
239659ac8732SStefano Zampini   mumps->sym = 2;
239759ac8732SStefano Zampini #else
23986fdc2a6dSBarry Smith   if (A->spd_set && A->spd) mumps->sym = 1;
23996fdc2a6dSBarry Smith   else                      mumps->sym = 2;
240059ac8732SStefano Zampini #endif
2401a214ac2aSShri Abhyankar 
240200c67f3bSHong Zhang   /* set solvertype */
240300c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
240400c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
240500c67f3bSHong Zhang 
2406f3c0ef26SHong Zhang   B->ops->destroy = MatDestroy_MUMPS;
2407e69c285eSBarry Smith   B->data         = (void*)mumps;
24082205254eSKarl Rupp 
2409f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2410746480a1SHong Zhang 
24112877fffaSHong Zhang   *F = B;
24122877fffaSHong Zhang   PetscFunctionReturn(0);
24132877fffaSHong Zhang }
241497969023SHong Zhang 
2415cc2e6a90SBarry Smith static PetscErrorCode MatGetFactor_baij_mumps(Mat A,MatFactorType ftype,Mat *F)
241667877ebaSShri Abhyankar {
241767877ebaSShri Abhyankar   Mat            B;
241867877ebaSShri Abhyankar   PetscErrorCode ierr;
241967877ebaSShri Abhyankar   Mat_MUMPS      *mumps;
2420ace3abfcSBarry Smith   PetscBool      isSeqBAIJ;
242167877ebaSShri Abhyankar 
242267877ebaSShri Abhyankar   PetscFunctionBegin;
242367877ebaSShri Abhyankar   /* Create the factorization matrix */
2424251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
2425ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
242667877ebaSShri Abhyankar   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
2427e69c285eSBarry Smith   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
2428e69c285eSBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
2429450b117fSShri Abhyankar 
2430b00a9115SJed Brown   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
2431450b117fSShri Abhyankar   if (ftype == MAT_FACTOR_LU) {
2432450b117fSShri Abhyankar     B->ops->lufactorsymbolic = MatLUFactorSymbolic_BAIJMUMPS;
2433450b117fSShri Abhyankar     B->factortype            = MAT_FACTOR_LU;
2434bccb9932SShri Abhyankar     if (isSeqBAIJ) mumps->ConvertToTriples = MatConvertToTriples_seqbaij_seqaij;
2435bccb9932SShri Abhyankar     else mumps->ConvertToTriples = MatConvertToTriples_mpibaij_mpiaij;
2436746480a1SHong Zhang     mumps->sym = 0;
2437f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc BAIJ matrices with MUMPS Cholesky, use SBAIJ or AIJ matrix instead\n");
2438bccb9932SShri Abhyankar 
2439e69c285eSBarry Smith   B->ops->getinfo     = MatGetInfo_External;
2440450b117fSShri Abhyankar   B->ops->view        = MatView_MUMPS;
24412205254eSKarl Rupp 
24423ca39a21SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
24435a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
24445a05ddb0SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
2445bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
2446bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
2447bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
2448bc6112feSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
2449ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
2450ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
2451ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
2452ca810319SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
245389a9c03aSHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInverse_C",MatMumpsGetInverse_MUMPS);CHKERRQ(ierr);
2454450b117fSShri Abhyankar 
245500c67f3bSHong Zhang   /* set solvertype */
245600c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
245700c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
245800c67f3bSHong Zhang 
24597ee00b23SStefano Zampini   B->ops->destroy = MatDestroy_MUMPS;
24607ee00b23SStefano Zampini   B->data         = (void*)mumps;
24617ee00b23SStefano Zampini 
24627ee00b23SStefano Zampini   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
24637ee00b23SStefano Zampini 
24647ee00b23SStefano Zampini   *F = B;
24657ee00b23SStefano Zampini   PetscFunctionReturn(0);
24667ee00b23SStefano Zampini }
24677ee00b23SStefano Zampini 
24687ee00b23SStefano Zampini /* MatGetFactor for Seq and MPI SELL matrices */
24697ee00b23SStefano Zampini static PetscErrorCode MatGetFactor_sell_mumps(Mat A,MatFactorType ftype,Mat *F)
24707ee00b23SStefano Zampini {
24717ee00b23SStefano Zampini   Mat            B;
24727ee00b23SStefano Zampini   PetscErrorCode ierr;
24737ee00b23SStefano Zampini   Mat_MUMPS      *mumps;
24747ee00b23SStefano Zampini   PetscBool      isSeqSELL;
24757ee00b23SStefano Zampini 
24767ee00b23SStefano Zampini   PetscFunctionBegin;
24777ee00b23SStefano Zampini   /* Create the factorization matrix */
24787ee00b23SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSELL,&isSeqSELL);CHKERRQ(ierr);
24797ee00b23SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
24807ee00b23SStefano Zampini   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
24817ee00b23SStefano Zampini   ierr = PetscStrallocpy("mumps",&((PetscObject)B)->type_name);CHKERRQ(ierr);
24827ee00b23SStefano Zampini   ierr = MatSetUp(B);CHKERRQ(ierr);
24837ee00b23SStefano Zampini 
24847ee00b23SStefano Zampini   ierr = PetscNewLog(B,&mumps);CHKERRQ(ierr);
24857ee00b23SStefano Zampini 
24867ee00b23SStefano Zampini   B->ops->view        = MatView_MUMPS;
24877ee00b23SStefano Zampini   B->ops->getinfo     = MatGetInfo_MUMPS;
24887ee00b23SStefano Zampini 
24897ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mumps);CHKERRQ(ierr);
24907ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MUMPS);CHKERRQ(ierr);
24917ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MUMPS);CHKERRQ(ierr);
24927ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetIcntl_C",MatMumpsSetIcntl_MUMPS);CHKERRQ(ierr);
24937ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetIcntl_C",MatMumpsGetIcntl_MUMPS);CHKERRQ(ierr);
24947ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsSetCntl_C",MatMumpsSetCntl_MUMPS);CHKERRQ(ierr);
24957ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetCntl_C",MatMumpsGetCntl_MUMPS);CHKERRQ(ierr);
24967ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfo_C",MatMumpsGetInfo_MUMPS);CHKERRQ(ierr);
24977ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetInfog_C",MatMumpsGetInfog_MUMPS);CHKERRQ(ierr);
24987ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfo_C",MatMumpsGetRinfo_MUMPS);CHKERRQ(ierr);
24997ee00b23SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMumpsGetRinfog_C",MatMumpsGetRinfog_MUMPS);CHKERRQ(ierr);
25007ee00b23SStefano Zampini 
25017ee00b23SStefano Zampini   if (ftype == MAT_FACTOR_LU) {
25027ee00b23SStefano Zampini     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMUMPS;
25037ee00b23SStefano Zampini     B->factortype            = MAT_FACTOR_LU;
25047ee00b23SStefano Zampini     if (isSeqSELL) mumps->ConvertToTriples = MatConvertToTriples_seqsell_seqaij;
25057ee00b23SStefano Zampini     else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
25067ee00b23SStefano Zampini     mumps->sym = 0;
25077ee00b23SStefano Zampini   } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
25087ee00b23SStefano Zampini 
25097ee00b23SStefano Zampini   /* set solvertype */
25107ee00b23SStefano Zampini   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
25117ee00b23SStefano Zampini   ierr = PetscStrallocpy(MATSOLVERMUMPS,&B->solvertype);CHKERRQ(ierr);
25127ee00b23SStefano Zampini 
2513450b117fSShri Abhyankar   B->ops->destroy = MatDestroy_MUMPS;
2514e69c285eSBarry Smith   B->data         = (void*)mumps;
25152205254eSKarl Rupp 
2516f697e70eSHong Zhang   ierr = PetscInitializeMUMPS(A,mumps);CHKERRQ(ierr);
2517746480a1SHong Zhang 
2518450b117fSShri Abhyankar   *F = B;
2519450b117fSShri Abhyankar   PetscFunctionReturn(0);
2520450b117fSShri Abhyankar }
252142c9c57cSBarry Smith 
25223ca39a21SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MUMPS(void)
252342c9c57cSBarry Smith {
252442c9c57cSBarry Smith   PetscErrorCode ierr;
252542c9c57cSBarry Smith 
252642c9c57cSBarry Smith   PetscFunctionBegin;
25273ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
25283ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
25293ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
25303ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPIBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
25313ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATMPISBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
25323ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mumps);CHKERRQ(ierr);
25333ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mumps);CHKERRQ(ierr);
25343ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_baij_mumps);CHKERRQ(ierr);
25353ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_baij_mumps);CHKERRQ(ierr);
25363ca39a21SBarry Smith   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mumps);CHKERRQ(ierr);
25377ee00b23SStefano Zampini   ierr = MatSolverTypeRegister(MATSOLVERMUMPS,MATSEQSELL,MAT_FACTOR_LU,MatGetFactor_sell_mumps);CHKERRQ(ierr);
253842c9c57cSBarry Smith   PetscFunctionReturn(0);
253942c9c57cSBarry Smith }
254042c9c57cSBarry Smith 
2541