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