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