xref: /petsc/src/mat/impls/aij/seq/mkl_pardiso/mkl_pardiso.c (revision 1873fc9fedfc3bd8fc2f82ccd29ab8b44a80f1db)
1 #include <../src/mat/impls/aij/seq/aij.h>        /*I "petscmat.h" I*/
2 #include <../src/mat/impls/sbaij/seq/sbaij.h>
3 #include <../src/mat/impls/dense/seq/dense.h>
4 
5 #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
6 #define MKL_ILP64
7 #endif
8 #include <mkl_pardiso.h>
9 
10 PETSC_EXTERN void PetscSetMKL_PARDISOThreads(int);
11 
12 /*
13  *  Possible mkl_pardiso phases that controls the execution of the solver.
14  *  For more information check mkl_pardiso manual.
15  */
16 #define JOB_ANALYSIS 11
17 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION 12
18 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13
19 #define JOB_NUMERICAL_FACTORIZATION 22
20 #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 23
21 #define JOB_SOLVE_ITERATIVE_REFINEMENT 33
22 #define JOB_SOLVE_FORWARD_SUBSTITUTION 331
23 #define JOB_SOLVE_DIAGONAL_SUBSTITUTION 332
24 #define JOB_SOLVE_BACKWARD_SUBSTITUTION 333
25 #define JOB_RELEASE_OF_LU_MEMORY 0
26 #define JOB_RELEASE_OF_ALL_MEMORY -1
27 
28 #define IPARM_SIZE 64
29 
30 #if defined(PETSC_USE_64BIT_INDICES)
31  #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
32   #define INT_TYPE long long int
33   #define MKL_PARDISO pardiso
34   #define MKL_PARDISO_INIT pardisoinit
35  #else
36   /* this is the case where the MKL BLAS/LAPACK are 32 bit integers but the 64 bit integer version of
37      of Pardiso code is used; hence the need for the 64 below*/
38   #define INT_TYPE long long int
39   #define MKL_PARDISO pardiso_64
40   #define MKL_PARDISO_INIT pardiso_64init
41 void pardiso_64init(void *pt, INT_TYPE *mtype, INT_TYPE iparm [])
42 {
43   int iparm_copy[IPARM_SIZE], mtype_copy, i;
44 
45   mtype_copy = *mtype;
46   pardisoinit(pt, &mtype_copy, iparm_copy);
47   for (i=0; i<IPARM_SIZE; i++) iparm[i] = iparm_copy[i];
48 }
49  #endif
50 #else
51  #define INT_TYPE int
52  #define MKL_PARDISO pardiso
53  #define MKL_PARDISO_INIT pardisoinit
54 #endif
55 
56 
57 /*
58  *  Internal data structure.
59  *  For more information check mkl_pardiso manual.
60  */
61 typedef struct {
62 
63   /* Configuration vector*/
64   INT_TYPE     iparm[IPARM_SIZE];
65 
66   /*
67    * Internal mkl_pardiso memory location.
68    * After the first call to mkl_pardiso do not modify pt, as that could cause a serious memory leak.
69    */
70   void         *pt[IPARM_SIZE];
71 
72   /* Basic mkl_pardiso info*/
73   INT_TYPE     phase, maxfct, mnum, mtype, n, nrhs, msglvl, err;
74 
75   /* Matrix structure*/
76   void         *a;
77   INT_TYPE     *ia, *ja;
78 
79   /* Number of non-zero elements*/
80   INT_TYPE     nz;
81 
82   /* Row permutaton vector*/
83   INT_TYPE     *perm;
84 
85   /* Define if matrix preserves sparse structure.*/
86   MatStructure matstruc;
87 
88   PetscBool    needsym;
89   PetscBool    freeaij;
90 
91   /* Schur complement */
92   PetscScalar  *schur;
93   PetscInt     schur_size;
94   PetscInt     *schur_idxs;
95   PetscScalar  *schur_work;
96   PetscBLASInt schur_work_size;
97   PetscBool    solve_interior;
98 
99   /* True if mkl_pardiso function have been used.*/
100   PetscBool CleanUp;
101 
102   /* Conversion to a format suitable for MKL */
103   PetscErrorCode (*Convert)(Mat, PetscBool, MatReuse, PetscBool*, INT_TYPE*, INT_TYPE**, INT_TYPE**, PetscScalar**);
104 } Mat_MKL_PARDISO;
105 
106 PetscErrorCode MatMKLPardiso_Convert_seqsbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
107 {
108   Mat_SeqSBAIJ   *aa = (Mat_SeqSBAIJ*)A->data;
109   PetscInt       bs  = A->rmap->bs,i;
110   PetscErrorCode ierr;
111 
112   PetscFunctionBegin;
113   if (!sym) {
114     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
115   }
116   *v      = aa->a;
117   if (bs == 1) { /* already in the correct format */
118     /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
119     *r    = (INT_TYPE*)aa->i;
120     *c    = (INT_TYPE*)aa->j;
121     *nnz  = (INT_TYPE)aa->nz;
122     *free = PETSC_FALSE;
123   } else if (reuse == MAT_INITIAL_MATRIX) {
124     PetscInt m = A->rmap->n,nz = aa->nz;
125     PetscInt *row,*col;
126     ierr = PetscMalloc2(m+1,&row,nz,&col);CHKERRQ(ierr);
127     for (i=0; i<m+1; i++) {
128       row[i] = aa->i[i]+1;
129     }
130     for (i=0; i<nz; i++) {
131       col[i] = aa->j[i]+1;
132     }
133     *r    = (INT_TYPE*)row;
134     *c    = (INT_TYPE*)col;
135     *nnz  = (INT_TYPE)nz;
136     *free = PETSC_TRUE;
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 PetscErrorCode MatMKLPardiso_Convert_seqbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
142 {
143   Mat_SeqBAIJ    *aa = (Mat_SeqBAIJ*)A->data;
144   PetscInt       bs  = A->rmap->bs,i;
145   PetscErrorCode ierr;
146 
147   PetscFunctionBegin;
148   if (!sym) {
149     *v      = aa->a;
150     if (bs == 1) { /* already in the correct format */
151       /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
152       *r    = (INT_TYPE*)aa->i;
153       *c    = (INT_TYPE*)aa->j;
154       *nnz  = (INT_TYPE)aa->nz;
155       *free = PETSC_FALSE;
156       PetscFunctionReturn(0);
157     } else if (reuse == MAT_INITIAL_MATRIX) {
158       PetscInt m = A->rmap->n,nz = aa->nz;
159       PetscInt *row,*col;
160       ierr = PetscMalloc2(m+1,&row,nz,&col);CHKERRQ(ierr);
161       for (i=0; i<m+1; i++) {
162         row[i] = aa->i[i]+1;
163       }
164       for (i=0; i<nz; i++) {
165         col[i] = aa->j[i]+1;
166       }
167       *r    = (INT_TYPE*)row;
168       *c    = (INT_TYPE*)col;
169       *nnz  = (INT_TYPE)nz;
170     }
171     *free = PETSC_TRUE;
172   } else {
173     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
174   }
175   PetscFunctionReturn(0);
176 }
177 
178 PetscErrorCode MatMKLPardiso_Convert_seqaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
179 {
180   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)A->data;
181   PetscScalar    *aav;
182   PetscErrorCode ierr;
183 
184   PetscFunctionBegin;
185   ierr = MatSeqAIJGetArrayRead(A,(const PetscScalar**)&aav);CHKERRQ(ierr);
186   if (!sym) { /* already in the correct format */
187     *v    = aav;
188     *r    = (INT_TYPE*)aa->i;
189     *c    = (INT_TYPE*)aa->j;
190     *nnz  = (INT_TYPE)aa->nz;
191     *free = PETSC_FALSE;
192   } else if (reuse == MAT_INITIAL_MATRIX) { /* need to get the triangular part */
193     PetscScalar *vals,*vv;
194     PetscInt    *row,*col,*jj;
195     PetscInt    m = A->rmap->n,nz,i;
196 
197     nz = 0;
198     for (i=0; i<m; i++) nz += aa->i[i+1] - aa->diag[i];
199     ierr = PetscMalloc2(m+1,&row,nz,&col);CHKERRQ(ierr);
200     ierr = PetscMalloc1(nz,&vals);CHKERRQ(ierr);
201     jj = col;
202     vv = vals;
203 
204     row[0] = 0;
205     for (i=0; i<m; i++) {
206       PetscInt    *aj = aa->j + aa->diag[i];
207       PetscScalar *av = aav + aa->diag[i];
208       PetscInt    rl  = aa->i[i+1] - aa->diag[i],j;
209 
210       for (j=0; j<rl; j++) {
211         *jj = *aj; jj++; aj++;
212         *vv = *av; vv++; av++;
213       }
214       row[i+1] = row[i] + rl;
215     }
216     *v    = vals;
217     *r    = (INT_TYPE*)row;
218     *c    = (INT_TYPE*)col;
219     *nnz  = (INT_TYPE)nz;
220     *free = PETSC_TRUE;
221   } else {
222     PetscScalar *vv;
223     PetscInt    m = A->rmap->n,i;
224 
225     vv = *v;
226     for (i=0; i<m; i++) {
227       PetscScalar *av = aav + aa->diag[i];
228       PetscInt    rl  = aa->i[i+1] - aa->diag[i],j;
229       for (j=0; j<rl; j++) {
230         *vv = *av; vv++; av++;
231       }
232     }
233     *free = PETSC_TRUE;
234   }
235   ierr = MatSeqAIJRestoreArrayRead(A,(const PetscScalar**)&aav);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 
239 
240 static PetscErrorCode MatMKLPardisoSolveSchur_Private(Mat F, PetscScalar *B, PetscScalar *X)
241 {
242   Mat_MKL_PARDISO      *mpardiso = (Mat_MKL_PARDISO*)F->data;
243   Mat                  S,Xmat,Bmat;
244   MatFactorSchurStatus schurstatus;
245   PetscErrorCode       ierr;
246 
247   PetscFunctionBegin;
248   ierr = MatFactorFactorizeSchurComplement(F);CHKERRQ(ierr);
249   ierr = MatFactorGetSchurComplement(F,&S,&schurstatus);CHKERRQ(ierr);
250   if (X == B && schurstatus == MAT_FACTOR_SCHUR_INVERTED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"X and B cannot point to the same address");
251   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mpardiso->schur_size,mpardiso->nrhs,B,&Bmat);CHKERRQ(ierr);
252   ierr = MatCreateSeqDense(PETSC_COMM_SELF,mpardiso->schur_size,mpardiso->nrhs,X,&Xmat);CHKERRQ(ierr);
253   ierr = MatSetType(Bmat,((PetscObject)S)->type_name);CHKERRQ(ierr);
254   ierr = MatSetType(Xmat,((PetscObject)S)->type_name);CHKERRQ(ierr);
255 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
256   ierr = MatPinToCPU(Xmat,S->pinnedtocpu);CHKERRQ(ierr);
257   ierr = MatPinToCPU(Bmat,S->pinnedtocpu);CHKERRQ(ierr);
258 #endif
259 
260 #if defined(PETSC_USE_COMPLEX)
261   if (mpardiso->iparm[12-1] == 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Hermitian solve not implemented yet");
262 #endif
263 
264   switch (schurstatus) {
265   case MAT_FACTOR_SCHUR_FACTORED:
266     if (!mpardiso->iparm[12-1]) {
267       ierr = MatMatSolve(S,Bmat,Xmat);CHKERRQ(ierr);
268     } else { /* transpose solve */
269       ierr = MatMatSolveTranspose(S,Bmat,Xmat);CHKERRQ(ierr);
270     }
271     break;
272   case MAT_FACTOR_SCHUR_INVERTED:
273     if (!mpardiso->iparm[12-1]) {
274       ierr = MatMatMult(S,Bmat,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Xmat);CHKERRQ(ierr);
275     } else { /* transpose solve */
276       ierr = MatTransposeMatMult(S,Bmat,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Xmat);CHKERRQ(ierr);
277     }
278     break;
279   default:
280     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
281     break;
282   }
283   ierr = MatFactorRestoreSchurComplement(F,&S,schurstatus);CHKERRQ(ierr);
284   ierr = MatDestroy(&Bmat);CHKERRQ(ierr);
285   ierr = MatDestroy(&Xmat);CHKERRQ(ierr);
286   PetscFunctionReturn(0);
287 }
288 
289 PetscErrorCode MatFactorSetSchurIS_MKL_PARDISO(Mat F, IS is)
290 {
291   Mat_MKL_PARDISO   *mpardiso = (Mat_MKL_PARDISO*)F->data;
292   const PetscScalar *arr;
293   const PetscInt    *idxs;
294   PetscInt          size,i;
295   PetscMPIInt       csize;
296   PetscBool         sorted;
297   PetscErrorCode    ierr;
298 
299   PetscFunctionBegin;
300   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&csize);CHKERRQ(ierr);
301   if (csize > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MKL_PARDISO parallel Schur complements not yet supported from PETSc");
302   ierr = ISSorted(is,&sorted);CHKERRQ(ierr);
303   if (!sorted) {
304     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IS for MKL_PARDISO Schur complements needs to be sorted");
305   }
306   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
307   ierr = PetscFree(mpardiso->schur_work);CHKERRQ(ierr);
308   ierr = PetscBLASIntCast(PetscMax(mpardiso->n,2*size),&mpardiso->schur_work_size);CHKERRQ(ierr);
309   ierr = PetscMalloc1(mpardiso->schur_work_size,&mpardiso->schur_work);CHKERRQ(ierr);
310   ierr = MatDestroy(&F->schur);CHKERRQ(ierr);
311   ierr = MatCreateSeqDense(PETSC_COMM_SELF,size,size,NULL,&F->schur);CHKERRQ(ierr);
312   ierr = MatDenseGetArrayRead(F->schur,&arr);CHKERRQ(ierr);
313   mpardiso->schur      = (PetscScalar*)arr;
314   mpardiso->schur_size = size;
315   ierr = MatDenseRestoreArrayRead(F->schur,&arr);CHKERRQ(ierr);
316   if (mpardiso->mtype == 2) {
317     ierr = MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr);
318   }
319 
320   ierr = PetscFree(mpardiso->schur_idxs);CHKERRQ(ierr);
321   ierr = PetscMalloc1(size,&mpardiso->schur_idxs);CHKERRQ(ierr);
322   ierr = PetscArrayzero(mpardiso->perm,mpardiso->n);CHKERRQ(ierr);
323   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
324   ierr = PetscArraycpy(mpardiso->schur_idxs,idxs,size);CHKERRQ(ierr);
325   for (i=0;i<size;i++) mpardiso->perm[idxs[i]] = 1;
326   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
327   if (size) { /* turn on Schur switch if the set of indices is not empty */
328     mpardiso->iparm[36-1] = 2;
329   }
330   PetscFunctionReturn(0);
331 }
332 
333 PetscErrorCode MatDestroy_MKL_PARDISO(Mat A)
334 {
335   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
336   PetscErrorCode  ierr;
337 
338   PetscFunctionBegin;
339   if (mat_mkl_pardiso->CleanUp) {
340     mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;
341 
342     MKL_PARDISO (mat_mkl_pardiso->pt,
343       &mat_mkl_pardiso->maxfct,
344       &mat_mkl_pardiso->mnum,
345       &mat_mkl_pardiso->mtype,
346       &mat_mkl_pardiso->phase,
347       &mat_mkl_pardiso->n,
348       NULL,
349       NULL,
350       NULL,
351       NULL,
352       &mat_mkl_pardiso->nrhs,
353       mat_mkl_pardiso->iparm,
354       &mat_mkl_pardiso->msglvl,
355       NULL,
356       NULL,
357       &mat_mkl_pardiso->err);
358   }
359   ierr = PetscFree(mat_mkl_pardiso->perm);CHKERRQ(ierr);
360   ierr = PetscFree(mat_mkl_pardiso->schur_work);CHKERRQ(ierr);
361   ierr = PetscFree(mat_mkl_pardiso->schur_idxs);CHKERRQ(ierr);
362   if (mat_mkl_pardiso->freeaij) {
363     ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr);
364     if (mat_mkl_pardiso->iparm[34] == 1) {
365       ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr);
366     }
367   }
368   ierr = PetscFree(A->data);CHKERRQ(ierr);
369 
370   /* clear composed functions */
371   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);CHKERRQ(ierr);
372   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr);
373   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);CHKERRQ(ierr);
374   PetscFunctionReturn(0);
375 }
376 
377 static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce)
378 {
379   PetscFunctionBegin;
380   if (reduce) { /* data given for the whole matrix */
381     PetscInt i,m=0,p=0;
382     for (i=0;i<mpardiso->nrhs;i++) {
383       PetscInt j;
384       for (j=0;j<mpardiso->schur_size;j++) {
385         schur[p+j] = whole[m+mpardiso->schur_idxs[j]];
386       }
387       m += mpardiso->n;
388       p += mpardiso->schur_size;
389     }
390   } else { /* from Schur to whole */
391     PetscInt i,m=0,p=0;
392     for (i=0;i<mpardiso->nrhs;i++) {
393       PetscInt j;
394       for (j=0;j<mpardiso->schur_size;j++) {
395         whole[m+mpardiso->schur_idxs[j]] = schur[p+j];
396       }
397       m += mpardiso->n;
398       p += mpardiso->schur_size;
399     }
400   }
401   PetscFunctionReturn(0);
402 }
403 
404 PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x)
405 {
406   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
407   PetscErrorCode    ierr;
408   PetscScalar       *xarray;
409   const PetscScalar *barray;
410 
411   PetscFunctionBegin;
412   mat_mkl_pardiso->nrhs = 1;
413   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
414   ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
415 
416   if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
417   else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
418 
419   if (barray == xarray) { /* if the two vectors share the same memory */
420     PetscScalar *work;
421     if (!mat_mkl_pardiso->schur_work) {
422       ierr = PetscMalloc1(mat_mkl_pardiso->n,&work);CHKERRQ(ierr);
423     } else {
424       work = mat_mkl_pardiso->schur_work;
425     }
426     mat_mkl_pardiso->iparm[6-1] = 1;
427     MKL_PARDISO (mat_mkl_pardiso->pt,
428       &mat_mkl_pardiso->maxfct,
429       &mat_mkl_pardiso->mnum,
430       &mat_mkl_pardiso->mtype,
431       &mat_mkl_pardiso->phase,
432       &mat_mkl_pardiso->n,
433       mat_mkl_pardiso->a,
434       mat_mkl_pardiso->ia,
435       mat_mkl_pardiso->ja,
436       NULL,
437       &mat_mkl_pardiso->nrhs,
438       mat_mkl_pardiso->iparm,
439       &mat_mkl_pardiso->msglvl,
440       (void*)xarray,
441       (void*)work,
442       &mat_mkl_pardiso->err);
443     if (!mat_mkl_pardiso->schur_work) {
444       ierr = PetscFree(work);CHKERRQ(ierr);
445     }
446   } else {
447     mat_mkl_pardiso->iparm[6-1] = 0;
448     MKL_PARDISO (mat_mkl_pardiso->pt,
449       &mat_mkl_pardiso->maxfct,
450       &mat_mkl_pardiso->mnum,
451       &mat_mkl_pardiso->mtype,
452       &mat_mkl_pardiso->phase,
453       &mat_mkl_pardiso->n,
454       mat_mkl_pardiso->a,
455       mat_mkl_pardiso->ia,
456       mat_mkl_pardiso->ja,
457       mat_mkl_pardiso->perm,
458       &mat_mkl_pardiso->nrhs,
459       mat_mkl_pardiso->iparm,
460       &mat_mkl_pardiso->msglvl,
461       (void*)barray,
462       (void*)xarray,
463       &mat_mkl_pardiso->err);
464   }
465   ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
466 
467   if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
468 
469   if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
470     PetscInt shift = mat_mkl_pardiso->schur_size;
471 
472     /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
473     if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;
474 
475     if (!mat_mkl_pardiso->solve_interior) {
476       /* solve Schur complement */
477       ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr);
478       ierr = MatMKLPardisoSolveSchur_Private(A,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr);
479       ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr);
480     } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substituted to xarray[schur] will be neglected */
481       PetscInt i;
482       for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
483         xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.;
484       }
485     }
486 
487     /* expansion phase */
488     mat_mkl_pardiso->iparm[6-1] = 1;
489     mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
490     MKL_PARDISO (mat_mkl_pardiso->pt,
491       &mat_mkl_pardiso->maxfct,
492       &mat_mkl_pardiso->mnum,
493       &mat_mkl_pardiso->mtype,
494       &mat_mkl_pardiso->phase,
495       &mat_mkl_pardiso->n,
496       mat_mkl_pardiso->a,
497       mat_mkl_pardiso->ia,
498       mat_mkl_pardiso->ja,
499       mat_mkl_pardiso->perm,
500       &mat_mkl_pardiso->nrhs,
501       mat_mkl_pardiso->iparm,
502       &mat_mkl_pardiso->msglvl,
503       (void*)xarray,
504       (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
505       &mat_mkl_pardiso->err);
506 
507     if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
508     mat_mkl_pardiso->iparm[6-1] = 0;
509   }
510   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
511   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
512   PetscFunctionReturn(0);
513 }
514 
515 PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x)
516 {
517   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
518   PetscInt        oiparm12;
519   PetscErrorCode  ierr;
520 
521   PetscFunctionBegin;
522   oiparm12 = mat_mkl_pardiso->iparm[12 - 1];
523   mat_mkl_pardiso->iparm[12 - 1] = 2;
524   ierr = MatSolve_MKL_PARDISO(A,b,x);CHKERRQ(ierr);
525   mat_mkl_pardiso->iparm[12 - 1] = oiparm12;
526   PetscFunctionReturn(0);
527 }
528 
529 PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X)
530 {
531   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->data;
532   PetscErrorCode    ierr;
533   const PetscScalar *barray;
534   PetscScalar       *xarray;
535   PetscBool         flg;
536 
537   PetscFunctionBegin;
538   ierr = PetscObjectBaseTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr);
539   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix");
540   if (X != B) {
541     ierr = PetscObjectBaseTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr);
542     if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix");
543   }
544 
545   ierr = MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);CHKERRQ(ierr);
546 
547   if (mat_mkl_pardiso->nrhs > 0) {
548     ierr = MatDenseGetArrayRead(B,&barray);CHKERRQ(ierr);
549     ierr = MatDenseGetArray(X,&xarray);CHKERRQ(ierr);
550 
551     if (barray == xarray) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"B and X cannot share the same memory location");
552     if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
553     else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
554 
555     MKL_PARDISO (mat_mkl_pardiso->pt,
556       &mat_mkl_pardiso->maxfct,
557       &mat_mkl_pardiso->mnum,
558       &mat_mkl_pardiso->mtype,
559       &mat_mkl_pardiso->phase,
560       &mat_mkl_pardiso->n,
561       mat_mkl_pardiso->a,
562       mat_mkl_pardiso->ia,
563       mat_mkl_pardiso->ja,
564       mat_mkl_pardiso->perm,
565       &mat_mkl_pardiso->nrhs,
566       mat_mkl_pardiso->iparm,
567       &mat_mkl_pardiso->msglvl,
568       (void*)barray,
569       (void*)xarray,
570       &mat_mkl_pardiso->err);
571     if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
572 
573     ierr = MatDenseRestoreArrayRead(B,&barray);CHKERRQ(ierr);
574     if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
575       PetscScalar *o_schur_work = NULL;
576       PetscInt    shift = mat_mkl_pardiso->schur_size*mat_mkl_pardiso->nrhs,scale;
577       PetscInt    mem = mat_mkl_pardiso->n*mat_mkl_pardiso->nrhs;
578 
579       /* allocate extra memory if it is needed */
580       scale = 1;
581       if (A->schur_status == MAT_FACTOR_SCHUR_INVERTED) scale = 2;
582 
583       mem *= scale;
584       if (mem > mat_mkl_pardiso->schur_work_size) {
585         o_schur_work = mat_mkl_pardiso->schur_work;
586         ierr = PetscMalloc1(mem,&mat_mkl_pardiso->schur_work);CHKERRQ(ierr);
587       }
588 
589       /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
590       if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;
591 
592       /* solve Schur complement */
593       if (!mat_mkl_pardiso->solve_interior) {
594         ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr);
595         ierr = MatMKLPardisoSolveSchur_Private(A,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr);
596         ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr);
597       } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substituted to xarray[schur,n] will be neglected */
598         PetscInt i,n,m=0;
599         for (n=0;n<mat_mkl_pardiso->nrhs;n++) {
600           for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
601             xarray[mat_mkl_pardiso->schur_idxs[i]+m] = 0.;
602           }
603           m += mat_mkl_pardiso->n;
604         }
605       }
606 
607       /* expansion phase */
608       mat_mkl_pardiso->iparm[6-1] = 1;
609       mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
610       MKL_PARDISO (mat_mkl_pardiso->pt,
611         &mat_mkl_pardiso->maxfct,
612         &mat_mkl_pardiso->mnum,
613         &mat_mkl_pardiso->mtype,
614         &mat_mkl_pardiso->phase,
615         &mat_mkl_pardiso->n,
616         mat_mkl_pardiso->a,
617         mat_mkl_pardiso->ia,
618         mat_mkl_pardiso->ja,
619         mat_mkl_pardiso->perm,
620         &mat_mkl_pardiso->nrhs,
621         mat_mkl_pardiso->iparm,
622         &mat_mkl_pardiso->msglvl,
623         (void*)xarray,
624         (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
625         &mat_mkl_pardiso->err);
626       if (o_schur_work) { /* restore original schur_work (minimal size) */
627         ierr = PetscFree(mat_mkl_pardiso->schur_work);CHKERRQ(ierr);
628         mat_mkl_pardiso->schur_work = o_schur_work;
629       }
630       if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
631       mat_mkl_pardiso->iparm[6-1] = 0;
632     }
633   }
634   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
635   PetscFunctionReturn(0);
636 }
637 
638 PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info)
639 {
640   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->data;
641   PetscErrorCode  ierr;
642 
643   PetscFunctionBegin;
644   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
645   ierr = (*mat_mkl_pardiso->Convert)(A,mat_mkl_pardiso->needsym,MAT_REUSE_MATRIX,&mat_mkl_pardiso->freeaij,&mat_mkl_pardiso->nz,&mat_mkl_pardiso->ia,&mat_mkl_pardiso->ja,(PetscScalar**)&mat_mkl_pardiso->a);CHKERRQ(ierr);
646 
647   mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION;
648   MKL_PARDISO (mat_mkl_pardiso->pt,
649     &mat_mkl_pardiso->maxfct,
650     &mat_mkl_pardiso->mnum,
651     &mat_mkl_pardiso->mtype,
652     &mat_mkl_pardiso->phase,
653     &mat_mkl_pardiso->n,
654     mat_mkl_pardiso->a,
655     mat_mkl_pardiso->ia,
656     mat_mkl_pardiso->ja,
657     mat_mkl_pardiso->perm,
658     &mat_mkl_pardiso->nrhs,
659     mat_mkl_pardiso->iparm,
660     &mat_mkl_pardiso->msglvl,
661     NULL,
662     (void*)mat_mkl_pardiso->schur,
663     &mat_mkl_pardiso->err);
664   if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
665 
666   /* report flops */
667   if (mat_mkl_pardiso->iparm[18] > 0) {
668     ierr = PetscLogFlops(PetscPowRealInt(10.,6)*mat_mkl_pardiso->iparm[18]);CHKERRQ(ierr);
669   }
670 
671   if (F->schur) { /* schur output from pardiso is in row major format */
672 #if defined(PETSC_HAVE_CUDA)
673     F->schur->offloadmask = PETSC_OFFLOAD_CPU;
674 #endif
675     ierr = MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);CHKERRQ(ierr);
676     ierr = MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);CHKERRQ(ierr);
677   }
678   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
679   mat_mkl_pardiso->CleanUp  = PETSC_TRUE;
680   PetscFunctionReturn(0);
681 }
682 
683 PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A)
684 {
685   Mat_MKL_PARDISO     *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
686   PetscErrorCode      ierr;
687   PetscInt            icntl,threads=1;
688   PetscBool           flg;
689 
690   PetscFunctionBegin;
691   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");CHKERRQ(ierr);
692 
693   ierr = PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use within PARDISO","None",threads,&threads,&flg);CHKERRQ(ierr);
694   if (flg) PetscSetMKL_PARDISOThreads((int)threads);
695 
696   ierr = PetscOptionsInt("-mat_mkl_pardiso_66","Maximum number of factors with identical sparsity structure that must be kept in memory at the same time","None",mat_mkl_pardiso->maxfct,&icntl,&flg);CHKERRQ(ierr);
697   if (flg) mat_mkl_pardiso->maxfct = icntl;
698 
699   ierr = PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);CHKERRQ(ierr);
700   if (flg) mat_mkl_pardiso->mnum = icntl;
701 
702   ierr = PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);CHKERRQ(ierr);
703   if (flg) mat_mkl_pardiso->msglvl = icntl;
704 
705   ierr = PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);CHKERRQ(ierr);
706   if (flg) {
707     void *pt[IPARM_SIZE];
708     mat_mkl_pardiso->mtype = icntl;
709     MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
710 #if defined(PETSC_USE_REAL_SINGLE)
711     mat_mkl_pardiso->iparm[27] = 1;
712 #else
713     mat_mkl_pardiso->iparm[27] = 0;
714 #endif
715     mat_mkl_pardiso->iparm[34] = 1; /* use 0-based indexing */
716   }
717   ierr = PetscOptionsInt("-mat_mkl_pardiso_1","Use default values (if 0)","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);CHKERRQ(ierr);
718 
719   if (flg && icntl != 0) {
720     ierr = PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);CHKERRQ(ierr);
721     if (flg) mat_mkl_pardiso->iparm[1] = icntl;
722 
723     ierr = PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);CHKERRQ(ierr);
724     if (flg) mat_mkl_pardiso->iparm[3] = icntl;
725 
726     ierr = PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);CHKERRQ(ierr);
727     if (flg) mat_mkl_pardiso->iparm[4] = icntl;
728 
729     ierr = PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);CHKERRQ(ierr);
730     if (flg) mat_mkl_pardiso->iparm[5] = icntl;
731 
732     ierr = PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);CHKERRQ(ierr);
733     if (flg) mat_mkl_pardiso->iparm[7] = icntl;
734 
735     ierr = PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);CHKERRQ(ierr);
736     if (flg) mat_mkl_pardiso->iparm[9] = icntl;
737 
738     ierr = PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);CHKERRQ(ierr);
739     if (flg) mat_mkl_pardiso->iparm[10] = icntl;
740 
741     ierr = PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);CHKERRQ(ierr);
742     if (flg) mat_mkl_pardiso->iparm[11] = icntl;
743 
744     ierr = PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);CHKERRQ(ierr);
745     if (flg) mat_mkl_pardiso->iparm[12] = icntl;
746 
747     ierr = PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);CHKERRQ(ierr);
748     if (flg) mat_mkl_pardiso->iparm[17] = icntl;
749 
750     ierr = PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations (0 to disable)","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);CHKERRQ(ierr);
751     if (flg) mat_mkl_pardiso->iparm[18] = icntl;
752 
753     ierr = PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);CHKERRQ(ierr);
754     if (flg) mat_mkl_pardiso->iparm[20] = icntl;
755 
756     ierr = PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);CHKERRQ(ierr);
757     if (flg) mat_mkl_pardiso->iparm[23] = icntl;
758 
759     ierr = PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);CHKERRQ(ierr);
760     if (flg) mat_mkl_pardiso->iparm[24] = icntl;
761 
762     ierr = PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);CHKERRQ(ierr);
763     if (flg) mat_mkl_pardiso->iparm[26] = icntl;
764 
765     ierr = PetscOptionsInt("-mat_mkl_pardiso_31","Partial solve and computing selected components of the solution vectors","None",mat_mkl_pardiso->iparm[30],&icntl,&flg);CHKERRQ(ierr);
766     if (flg) mat_mkl_pardiso->iparm[30] = icntl;
767 
768     ierr = PetscOptionsInt("-mat_mkl_pardiso_34","Optimal number of threads for conditional numerical reproducibility (CNR) mode","None",mat_mkl_pardiso->iparm[33],&icntl,&flg);CHKERRQ(ierr);
769     if (flg) mat_mkl_pardiso->iparm[33] = icntl;
770 
771     ierr = PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);CHKERRQ(ierr);
772     if (flg) mat_mkl_pardiso->iparm[59] = icntl;
773   }
774   PetscOptionsEnd();
775   PetscFunctionReturn(0);
776 }
777 
778 PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso)
779 {
780   PetscErrorCode ierr;
781   PetscInt       i,bs;
782   PetscBool      match;
783 
784   PetscFunctionBegin;
785   for (i=0; i<IPARM_SIZE; i++) mat_mkl_pardiso->iparm[i] = 0;
786   for (i=0; i<IPARM_SIZE; i++) mat_mkl_pardiso->pt[i] = 0;
787   /* Default options for both sym and unsym */
788   mat_mkl_pardiso->iparm[ 0] =  1; /* Solver default parameters overriden with provided by iparm */
789   mat_mkl_pardiso->iparm[ 1] =  2; /* Metis reordering */
790   mat_mkl_pardiso->iparm[ 5] =  0; /* Write solution into x */
791   mat_mkl_pardiso->iparm[ 7] =  0; /* Max number of iterative refinement steps */
792   mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
793   mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
794 #if 0
795   mat_mkl_pardiso->iparm[23] =  1; /* Parallel factorization control*/
796 #endif
797   ierr = PetscObjectTypeCompareAny((PetscObject)A,&match,MATSEQBAIJ,MATSEQSBAIJ,"");CHKERRQ(ierr);
798   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
799   if (!match || bs == 1) {
800     mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */
801     mat_mkl_pardiso->n         = A->rmap->N;
802   } else {
803     mat_mkl_pardiso->iparm[34] = 0; /* Cluster Sparse Solver use Fortran-style indexing for ia and ja arrays */
804     mat_mkl_pardiso->iparm[36] = bs;
805     mat_mkl_pardiso->n         = A->rmap->N/bs;
806   }
807   mat_mkl_pardiso->iparm[39] =  0; /* Input: matrix/rhs/solution stored on master */
808 
809   mat_mkl_pardiso->CleanUp   = PETSC_FALSE;
810   mat_mkl_pardiso->maxfct    = 1; /* Maximum number of numerical factorizations. */
811   mat_mkl_pardiso->mnum      = 1; /* Which factorization to use. */
812   mat_mkl_pardiso->msglvl    = 0; /* 0: do not print 1: Print statistical information in file */
813   mat_mkl_pardiso->phase     = -1;
814   mat_mkl_pardiso->err       = 0;
815 
816   mat_mkl_pardiso->nrhs      = 1;
817   mat_mkl_pardiso->err       = 0;
818   mat_mkl_pardiso->phase     = -1;
819 
820   if (ftype == MAT_FACTOR_LU) {
821     mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */
822     mat_mkl_pardiso->iparm[10] =  1; /* Use nonsymmetric permutation and scaling MPS */
823     mat_mkl_pardiso->iparm[12] =  1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
824   } else {
825     mat_mkl_pardiso->iparm[ 9] = 8; /* Perturb the pivot elements with 1E-8 */
826     mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */
827     mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
828 #if defined(PETSC_USE_DEBUG)
829     mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */
830 #endif
831   }
832   ierr = PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);CHKERRQ(ierr);
833   for (i=0; i<A->rmap->N; i++) {
834     mat_mkl_pardiso->perm[i] = 0;
835   }
836   mat_mkl_pardiso->schur_size = 0;
837   PetscFunctionReturn(0);
838 }
839 
840 PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info)
841 {
842   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
843   PetscErrorCode  ierr;
844 
845   PetscFunctionBegin;
846   mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
847   ierr = PetscSetMKL_PARDISOFromOptions(F,A);CHKERRQ(ierr);
848 
849   /* throw away any previously computed structure */
850   if (mat_mkl_pardiso->freeaij) {
851     ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr);
852     if (mat_mkl_pardiso->iparm[34] == 1) {
853       ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr);
854     }
855   }
856   ierr = (*mat_mkl_pardiso->Convert)(A,mat_mkl_pardiso->needsym,MAT_INITIAL_MATRIX,&mat_mkl_pardiso->freeaij,&mat_mkl_pardiso->nz,&mat_mkl_pardiso->ia,&mat_mkl_pardiso->ja,(PetscScalar**)&mat_mkl_pardiso->a);CHKERRQ(ierr);
857   if (mat_mkl_pardiso->iparm[34] == 1) mat_mkl_pardiso->n = A->rmap->N;
858   else mat_mkl_pardiso->n = A->rmap->N/A->rmap->bs;
859 
860   mat_mkl_pardiso->phase = JOB_ANALYSIS;
861 
862   /* reset flops counting if requested */
863   if (mat_mkl_pardiso->iparm[18]) mat_mkl_pardiso->iparm[18] = -1;
864 
865   MKL_PARDISO (mat_mkl_pardiso->pt,
866     &mat_mkl_pardiso->maxfct,
867     &mat_mkl_pardiso->mnum,
868     &mat_mkl_pardiso->mtype,
869     &mat_mkl_pardiso->phase,
870     &mat_mkl_pardiso->n,
871     mat_mkl_pardiso->a,
872     mat_mkl_pardiso->ia,
873     mat_mkl_pardiso->ja,
874     mat_mkl_pardiso->perm,
875     &mat_mkl_pardiso->nrhs,
876     mat_mkl_pardiso->iparm,
877     &mat_mkl_pardiso->msglvl,
878     NULL,
879     NULL,
880     &mat_mkl_pardiso->err);
881   if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
882 
883   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
884 
885   if (F->factortype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO;
886   else F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO;
887 
888   F->ops->solve           = MatSolve_MKL_PARDISO;
889   F->ops->solvetranspose  = MatSolveTranspose_MKL_PARDISO;
890   F->ops->matsolve        = MatMatSolve_MKL_PARDISO;
891   PetscFunctionReturn(0);
892 }
893 
894 PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
895 {
896   PetscErrorCode ierr;
897 
898   PetscFunctionBegin;
899   ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr);
900   PetscFunctionReturn(0);
901 }
902 
903 #if !defined(PETSC_USE_COMPLEX)
904 PetscErrorCode MatGetInertia_MKL_PARDISO(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
905 {
906   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)F->data;
907 
908   PetscFunctionBegin;
909   if (nneg) *nneg = mat_mkl_pardiso->iparm[22];
910   if (npos) *npos = mat_mkl_pardiso->iparm[21];
911   if (nzero) *nzero = F->rmap->N - (mat_mkl_pardiso->iparm[22] + mat_mkl_pardiso->iparm[21]);
912   PetscFunctionReturn(0);
913 }
914 #endif
915 
916 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info)
917 {
918   PetscErrorCode ierr;
919 
920   PetscFunctionBegin;
921   ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr);
922 #if defined(PETSC_USE_COMPLEX)
923   F->ops->getinertia = NULL;
924 #else
925   F->ops->getinertia = MatGetInertia_MKL_PARDISO;
926 #endif
927   PetscFunctionReturn(0);
928 }
929 
930 PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer)
931 {
932   PetscErrorCode    ierr;
933   PetscBool         iascii;
934   PetscViewerFormat format;
935   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
936   PetscInt          i;
937 
938   PetscFunctionBegin;
939   if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(0);
940 
941   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
942   if (iascii) {
943     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
944     if (format == PETSC_VIEWER_ASCII_INFO) {
945       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");CHKERRQ(ierr);
946       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase:             %d \n",mat_mkl_pardiso->phase);CHKERRQ(ierr);
947       for (i=1; i<=64; i++) {
948         ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]:     %d \n",i, mat_mkl_pardiso->iparm[i - 1]);CHKERRQ(ierr);
949       }
950       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct:     %d \n", mat_mkl_pardiso->maxfct);CHKERRQ(ierr);
951       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum:     %d \n", mat_mkl_pardiso->mnum);CHKERRQ(ierr);
952       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype:     %d \n", mat_mkl_pardiso->mtype);CHKERRQ(ierr);
953       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n:     %d \n", mat_mkl_pardiso->n);CHKERRQ(ierr);
954       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs:     %d \n", mat_mkl_pardiso->nrhs);CHKERRQ(ierr);
955       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl:     %d \n", mat_mkl_pardiso->msglvl);CHKERRQ(ierr);
956     }
957   }
958   PetscFunctionReturn(0);
959 }
960 
961 
962 PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info)
963 {
964   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)A->data;
965 
966   PetscFunctionBegin;
967   info->block_size        = 1.0;
968   info->nz_used           = mat_mkl_pardiso->iparm[17];
969   info->nz_allocated      = mat_mkl_pardiso->iparm[17];
970   info->nz_unneeded       = 0.0;
971   info->assemblies        = 0.0;
972   info->mallocs           = 0.0;
973   info->memory            = 0.0;
974   info->fill_ratio_given  = 0;
975   info->fill_ratio_needed = 0;
976   info->factor_mallocs    = 0;
977   PetscFunctionReturn(0);
978 }
979 
980 PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival)
981 {
982   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
983 
984   PetscFunctionBegin;
985   if (icntl <= 64) {
986     mat_mkl_pardiso->iparm[icntl - 1] = ival;
987   } else {
988     if (icntl == 65) PetscSetMKL_PARDISOThreads(ival);
989     else if (icntl == 66) mat_mkl_pardiso->maxfct = ival;
990     else if (icntl == 67) mat_mkl_pardiso->mnum = ival;
991     else if (icntl == 68) mat_mkl_pardiso->msglvl = ival;
992     else if (icntl == 69) {
993       void *pt[IPARM_SIZE];
994       mat_mkl_pardiso->mtype = ival;
995       MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
996 #if defined(PETSC_USE_REAL_SINGLE)
997       mat_mkl_pardiso->iparm[27] = 1;
998 #else
999       mat_mkl_pardiso->iparm[27] = 0;
1000 #endif
1001       mat_mkl_pardiso->iparm[34] = 1;
1002     } else if (icntl==70) mat_mkl_pardiso->solve_interior = (PetscBool)!!ival;
1003   }
1004   PetscFunctionReturn(0);
1005 }
1006 
1007 /*@
1008   MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters
1009 
1010    Logically Collective on Mat
1011 
1012    Input Parameters:
1013 +  F - the factored matrix obtained by calling MatGetFactor()
1014 .  icntl - index of Mkl_Pardiso parameter
1015 -  ival - value of Mkl_Pardiso parameter
1016 
1017   Options Database:
1018 .   -mat_mkl_pardiso_<icntl> <ival>
1019 
1020    Level: beginner
1021 
1022    References:
1023 .      Mkl_Pardiso Users' Guide
1024 
1025 .seealso: MatGetFactor()
1026 @*/
1027 PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival)
1028 {
1029   PetscErrorCode ierr;
1030 
1031   PetscFunctionBegin;
1032   ierr = PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 /*MC
1037   MATSOLVERMKL_PARDISO -  A matrix type providing direct solvers (LU) for
1038   sequential matrices via the external package MKL_PARDISO.
1039 
1040   Works with MATSEQAIJ matrices
1041 
1042   Use -pc_type lu -pc_factor_mat_solver_type mkl_pardiso to use this direct solver
1043 
1044   Options Database Keys:
1045 + -mat_mkl_pardiso_65 - Number of threads to use within MKL_PARDISO
1046 . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
1047 . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase
1048 . -mat_mkl_pardiso_68 - Message level information
1049 . -mat_mkl_pardiso_69 - Defines the matrix type. IMPORTANT: When you set this flag, iparm parameters are going to be set to the default ones for the matrix type
1050 . -mat_mkl_pardiso_1  - Use default values
1051 . -mat_mkl_pardiso_2  - Fill-in reducing ordering for the input matrix
1052 . -mat_mkl_pardiso_4  - Preconditioned CGS/CG
1053 . -mat_mkl_pardiso_5  - User permutation
1054 . -mat_mkl_pardiso_6  - Write solution on x
1055 . -mat_mkl_pardiso_8  - Iterative refinement step
1056 . -mat_mkl_pardiso_10 - Pivoting perturbation
1057 . -mat_mkl_pardiso_11 - Scaling vectors
1058 . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A
1059 . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching
1060 . -mat_mkl_pardiso_18 - Numbers of non-zero elements
1061 . -mat_mkl_pardiso_19 - Report number of floating point operations
1062 . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices
1063 . -mat_mkl_pardiso_24 - Parallel factorization control
1064 . -mat_mkl_pardiso_25 - Parallel forward/backward solve control
1065 . -mat_mkl_pardiso_27 - Matrix checker
1066 . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors
1067 . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
1068 - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode
1069 
1070   Level: beginner
1071 
1072   For more information please check  mkl_pardiso manual
1073 
1074 .seealso: PCFactorSetMatSolverType(), MatSolverType
1075 
1076 M*/
1077 static PetscErrorCode MatFactorGetSolverType_mkl_pardiso(Mat A, MatSolverType *type)
1078 {
1079   PetscFunctionBegin;
1080   *type = MATSOLVERMKL_PARDISO;
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F)
1085 {
1086   Mat             B;
1087   PetscErrorCode  ierr;
1088   Mat_MKL_PARDISO *mat_mkl_pardiso;
1089   PetscBool       isSeqAIJ,isSeqBAIJ,isSeqSBAIJ;
1090 
1091   PetscFunctionBegin;
1092   ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1093   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1094   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
1095   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
1096   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1097   ierr = PetscStrallocpy("mkl_pardiso",&((PetscObject)B)->type_name);CHKERRQ(ierr);
1098   ierr = MatSetUp(B);CHKERRQ(ierr);
1099 
1100   ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr);
1101   B->data = mat_mkl_pardiso;
1102 
1103   ierr = MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr);
1104   if (ftype == MAT_FACTOR_LU) {
1105     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO;
1106     B->factortype            = MAT_FACTOR_LU;
1107     mat_mkl_pardiso->needsym = PETSC_FALSE;
1108     if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1109     else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1110     else if (isSeqSBAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead");
1111     else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU with %s format",((PetscObject)A)->type_name);
1112 #if defined(PETSC_USE_COMPLEX)
1113     mat_mkl_pardiso->mtype = 13;
1114 #else
1115     mat_mkl_pardiso->mtype = 11;
1116 #endif
1117   } else {
1118 #if defined(PETSC_USE_COMPLEX)
1119     SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with complex scalars! Use MAT_FACTOR_LU instead",((PetscObject)A)->type_name);
1120 #endif
1121     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO;
1122     B->factortype                  = MAT_FACTOR_CHOLESKY;
1123     if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1124     else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1125     else if (isSeqSBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij;
1126     else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with %s format",((PetscObject)A)->type_name);
1127 
1128     mat_mkl_pardiso->needsym = PETSC_TRUE;
1129     if (A->spd_set && A->spd) mat_mkl_pardiso->mtype = 2;
1130     else                      mat_mkl_pardiso->mtype = -2;
1131   }
1132   B->ops->destroy = MatDestroy_MKL_PARDISO;
1133   B->ops->view    = MatView_MKL_PARDISO;
1134   B->ops->getinfo = MatGetInfo_MKL_PARDISO;
1135   B->factortype   = ftype;
1136   B->assembled    = PETSC_TRUE;
1137 
1138   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
1139   ierr = PetscStrallocpy(MATSOLVERMKL_PARDISO,&B->solvertype);CHKERRQ(ierr);
1140 
1141   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mkl_pardiso);CHKERRQ(ierr);
1142   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MKL_PARDISO);CHKERRQ(ierr);
1143   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr);
1144 
1145   *F = B;
1146   PetscFunctionReturn(0);
1147 }
1148 
1149 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MKL_Pardiso(void)
1150 {
1151   PetscErrorCode ierr;
1152 
1153   PetscFunctionBegin;
1154   ierr = MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1155   ierr = MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1156   ierr = MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1157   ierr = MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1158   PetscFunctionReturn(0);
1159 }
1160 
1161