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