xref: /petsc/src/mat/impls/aij/mpi/mkl_cpardiso/mkl_cpardiso.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 #include <petscsys.h>
2 #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I  "petscmat.h"  I*/
3 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
4 
5 #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
6   #define MKL_ILP64
7 #endif
8 #include <mkl.h>
9 #include <mkl_cluster_sparse_solver.h>
10 
11 /*
12  *  Possible mkl_cpardiso phases that controls the execution of the solver.
13  *  For more information check mkl_cpardiso manual.
14  */
15 #define JOB_ANALYSIS                                                    11
16 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION                            12
17 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13
18 #define JOB_NUMERICAL_FACTORIZATION                                     22
19 #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT          23
20 #define JOB_SOLVE_ITERATIVE_REFINEMENT                                  33
21 #define JOB_SOLVE_FORWARD_SUBSTITUTION                                  331
22 #define JOB_SOLVE_DIAGONAL_SUBSTITUTION                                 332
23 #define JOB_SOLVE_BACKWARD_SUBSTITUTION                                 333
24 #define JOB_RELEASE_OF_LU_MEMORY                                        0
25 #define JOB_RELEASE_OF_ALL_MEMORY                                       -1
26 
27 #define IPARM_SIZE 64
28 #define INT_TYPE   MKL_INT
29 
30 static const char *Err_MSG_CPardiso(int errNo)
31 {
32   switch (errNo) {
33   case -1:
34     return "input inconsistent";
35     break;
36   case -2:
37     return "not enough memory";
38     break;
39   case -3:
40     return "reordering problem";
41     break;
42   case -4:
43     return "zero pivot, numerical factorization or iterative refinement problem";
44     break;
45   case -5:
46     return "unclassified (internal) error";
47     break;
48   case -6:
49     return "preordering failed (matrix types 11, 13 only)";
50     break;
51   case -7:
52     return "diagonal matrix problem";
53     break;
54   case -8:
55     return "32-bit integer overflow problem";
56     break;
57   case -9:
58     return "not enough memory for OOC";
59     break;
60   case -10:
61     return "problems with opening OOC temporary files";
62     break;
63   case -11:
64     return "read/write problems with the OOC data file";
65     break;
66   default:
67     return "unknown error";
68   }
69 }
70 
71 /*
72  *  Internal data structure.
73  *  For more information check mkl_cpardiso manual.
74  */
75 
76 typedef struct {
77   /* Configuration vector */
78   INT_TYPE iparm[IPARM_SIZE];
79 
80   /*
81    * Internal mkl_cpardiso memory location.
82    * After the first call to mkl_cpardiso do not modify pt, as that could cause a serious memory leak.
83    */
84   void *pt[IPARM_SIZE];
85 
86   MPI_Fint comm_mkl_cpardiso;
87 
88   /* Basic mkl_cpardiso info*/
89   INT_TYPE phase, maxfct, mnum, mtype, n, nrhs, msglvl, err;
90 
91   /* Matrix structure */
92   PetscScalar *a;
93 
94   INT_TYPE *ia, *ja;
95 
96   /* Number of non-zero elements */
97   INT_TYPE nz;
98 
99   /* Row permutaton vector*/
100   INT_TYPE *perm;
101 
102   /* Define is matrix preserve sparse structure. */
103   MatStructure matstruc;
104 
105   PetscErrorCode (*ConvertToTriples)(Mat, MatReuse, PetscInt *, PetscInt **, PetscInt **, PetscScalar **);
106 
107   /* True if mkl_cpardiso function have been used. */
108   PetscBool CleanUp;
109 } Mat_MKL_CPARDISO;
110 
111 /*
112  * Copy the elements of matrix A.
113  * Input:
114  *   - Mat A: MATSEQAIJ matrix
115  *   - int shift: matrix index.
116  *     - 0 for c representation
117  *     - 1 for fortran representation
118  *   - MatReuse reuse:
119  *     - MAT_INITIAL_MATRIX: Create a new aij representation
120  *     - MAT_REUSE_MATRIX: Reuse all aij representation and just change values
121  * Output:
122  *   - int *nnz: Number of nonzero-elements.
123  *   - int **r pointer to i index
124  *   - int **c pointer to j elements
125  *   - MATRIXTYPE **v: Non-zero elements
126  */
127 static PetscErrorCode MatCopy_seqaij_seqaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v)
128 {
129   Mat_SeqAIJ *aa = (Mat_SeqAIJ *)A->data;
130 
131   PetscFunctionBegin;
132   *v = aa->a;
133   if (reuse == MAT_INITIAL_MATRIX) {
134     *r   = (INT_TYPE *)aa->i;
135     *c   = (INT_TYPE *)aa->j;
136     *nnz = aa->nz;
137   }
138   PetscFunctionReturn(PETSC_SUCCESS);
139 }
140 
141 static PetscErrorCode MatConvertToTriples_mpiaij_mpiaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v)
142 {
143   const PetscInt    *ai, *aj, *bi, *bj, *garray, m = A->rmap->n, *ajj, *bjj;
144   PetscInt           rstart, nz, i, j, countA, countB;
145   PetscInt          *row, *col;
146   const PetscScalar *av, *bv;
147   PetscScalar       *val;
148   Mat_MPIAIJ        *mat = (Mat_MPIAIJ *)A->data;
149   Mat_SeqAIJ        *aa  = (Mat_SeqAIJ *)(mat->A)->data;
150   Mat_SeqAIJ        *bb  = (Mat_SeqAIJ *)(mat->B)->data;
151   PetscInt           colA_start, jB, jcol;
152 
153   PetscFunctionBegin;
154   ai     = aa->i;
155   aj     = aa->j;
156   bi     = bb->i;
157   bj     = bb->j;
158   rstart = A->rmap->rstart;
159   av     = aa->a;
160   bv     = bb->a;
161 
162   garray = mat->garray;
163 
164   if (reuse == MAT_INITIAL_MATRIX) {
165     nz   = aa->nz + bb->nz;
166     *nnz = nz;
167     PetscCall(PetscMalloc3(m + 1, &row, nz, &col, nz, &val));
168     *r = row;
169     *c = col;
170     *v = val;
171   } else {
172     row = *r;
173     col = *c;
174     val = *v;
175   }
176 
177   nz = 0;
178   for (i = 0; i < m; i++) {
179     row[i] = nz;
180     countA = ai[i + 1] - ai[i];
181     countB = bi[i + 1] - bi[i];
182     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
183     bjj    = bj + bi[i];
184 
185     /* B part, smaller col index */
186     colA_start = rstart + ajj[0]; /* the smallest global col index of A */
187     jB         = 0;
188     for (j = 0; j < countB; j++) {
189       jcol = garray[bjj[j]];
190       if (jcol > colA_start) break;
191       col[nz]   = jcol;
192       val[nz++] = *bv++;
193     }
194     jB = j;
195 
196     /* A part */
197     for (j = 0; j < countA; j++) {
198       col[nz]   = rstart + ajj[j];
199       val[nz++] = *av++;
200     }
201 
202     /* B part, larger col index */
203     for (j = jB; j < countB; j++) {
204       col[nz]   = garray[bjj[j]];
205       val[nz++] = *bv++;
206     }
207   }
208   row[m] = nz;
209 
210   PetscFunctionReturn(PETSC_SUCCESS);
211 }
212 
213 static PetscErrorCode MatConvertToTriples_mpibaij_mpibaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v)
214 {
215   const PetscInt    *ai, *aj, *bi, *bj, *garray, bs = A->rmap->bs, bs2 = bs * bs, m = A->rmap->n / bs, *ajj, *bjj;
216   PetscInt           rstart, nz, i, j, countA, countB;
217   PetscInt          *row, *col;
218   const PetscScalar *av, *bv;
219   PetscScalar       *val;
220   Mat_MPIBAIJ       *mat = (Mat_MPIBAIJ *)A->data;
221   Mat_SeqBAIJ       *aa  = (Mat_SeqBAIJ *)(mat->A)->data;
222   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ *)(mat->B)->data;
223   PetscInt           colA_start, jB, jcol;
224 
225   PetscFunctionBegin;
226   ai     = aa->i;
227   aj     = aa->j;
228   bi     = bb->i;
229   bj     = bb->j;
230   rstart = A->rmap->rstart / bs;
231   av     = aa->a;
232   bv     = bb->a;
233 
234   garray = mat->garray;
235 
236   if (reuse == MAT_INITIAL_MATRIX) {
237     nz   = aa->nz + bb->nz;
238     *nnz = nz;
239     PetscCall(PetscMalloc3(m + 1, &row, nz, &col, nz * bs2, &val));
240     *r = row;
241     *c = col;
242     *v = val;
243   } else {
244     row = *r;
245     col = *c;
246     val = *v;
247   }
248 
249   nz = 0;
250   for (i = 0; i < m; i++) {
251     row[i] = nz + 1;
252     countA = ai[i + 1] - ai[i];
253     countB = bi[i + 1] - bi[i];
254     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
255     bjj    = bj + bi[i];
256 
257     /* B part, smaller col index */
258     colA_start = rstart + (countA > 0 ? ajj[0] : 0); /* the smallest global col index of A */
259     jB         = 0;
260     for (j = 0; j < countB; j++) {
261       jcol = garray[bjj[j]];
262       if (jcol > colA_start) break;
263       col[nz++] = jcol + 1;
264     }
265     jB = j;
266     PetscCall(PetscArraycpy(val, bv, jB * bs2));
267     val += jB * bs2;
268     bv += jB * bs2;
269 
270     /* A part */
271     for (j = 0; j < countA; j++) col[nz++] = rstart + ajj[j] + 1;
272     PetscCall(PetscArraycpy(val, av, countA * bs2));
273     val += countA * bs2;
274     av += countA * bs2;
275 
276     /* B part, larger col index */
277     for (j = jB; j < countB; j++) col[nz++] = garray[bjj[j]] + 1;
278     PetscCall(PetscArraycpy(val, bv, (countB - jB) * bs2));
279     val += (countB - jB) * bs2;
280     bv += (countB - jB) * bs2;
281   }
282   row[m] = nz + 1;
283 
284   PetscFunctionReturn(PETSC_SUCCESS);
285 }
286 
287 static PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v)
288 {
289   const PetscInt    *ai, *aj, *bi, *bj, *garray, bs = A->rmap->bs, bs2 = bs * bs, m = A->rmap->n / bs, *ajj, *bjj;
290   PetscInt           rstart, nz, i, j, countA, countB;
291   PetscInt          *row, *col;
292   const PetscScalar *av, *bv;
293   PetscScalar       *val;
294   Mat_MPISBAIJ      *mat = (Mat_MPISBAIJ *)A->data;
295   Mat_SeqSBAIJ      *aa  = (Mat_SeqSBAIJ *)(mat->A)->data;
296   Mat_SeqBAIJ       *bb  = (Mat_SeqBAIJ *)(mat->B)->data;
297 
298   PetscFunctionBegin;
299   ai     = aa->i;
300   aj     = aa->j;
301   bi     = bb->i;
302   bj     = bb->j;
303   rstart = A->rmap->rstart / bs;
304   av     = aa->a;
305   bv     = bb->a;
306 
307   garray = mat->garray;
308 
309   if (reuse == MAT_INITIAL_MATRIX) {
310     nz   = aa->nz + bb->nz;
311     *nnz = nz;
312     PetscCall(PetscMalloc3(m + 1, &row, nz, &col, nz * bs2, &val));
313     *r = row;
314     *c = col;
315     *v = val;
316   } else {
317     row = *r;
318     col = *c;
319     val = *v;
320   }
321 
322   nz = 0;
323   for (i = 0; i < m; i++) {
324     row[i] = nz + 1;
325     countA = ai[i + 1] - ai[i];
326     countB = bi[i + 1] - bi[i];
327     ajj    = aj + ai[i]; /* ptr to the beginning of this row */
328     bjj    = bj + bi[i];
329 
330     /* A part */
331     for (j = 0; j < countA; j++) col[nz++] = rstart + ajj[j] + 1;
332     PetscCall(PetscArraycpy(val, av, countA * bs2));
333     val += countA * bs2;
334     av += countA * bs2;
335 
336     /* B part, larger col index */
337     for (j = 0; j < countB; j++) col[nz++] = garray[bjj[j]] + 1;
338     PetscCall(PetscArraycpy(val, bv, countB * bs2));
339     val += countB * bs2;
340     bv += countB * bs2;
341   }
342   row[m] = nz + 1;
343 
344   PetscFunctionReturn(PETSC_SUCCESS);
345 }
346 
347 /*
348  * Free memory for Mat_MKL_CPARDISO structure and pointers to objects.
349  */
350 static PetscErrorCode MatDestroy_MKL_CPARDISO(Mat A)
351 {
352   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data;
353   MPI_Comm          comm;
354 
355   PetscFunctionBegin;
356   /* Terminate instance, deallocate memories */
357   if (mat_mkl_cpardiso->CleanUp) {
358     mat_mkl_cpardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;
359 
360     cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, NULL, NULL, NULL, mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs,
361                           mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err);
362   }
363 
364   if (mat_mkl_cpardiso->ConvertToTriples != MatCopy_seqaij_seqaij_MKL_CPARDISO) PetscCall(PetscFree3(mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, mat_mkl_cpardiso->a));
365   comm = MPI_Comm_f2c(mat_mkl_cpardiso->comm_mkl_cpardiso);
366   PetscCallMPI(MPI_Comm_free(&comm));
367   PetscCall(PetscFree(A->data));
368 
369   /* clear composed functions */
370   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
371   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMkl_CPardisoSetCntl_C", NULL));
372   PetscFunctionReturn(PETSC_SUCCESS);
373 }
374 
375 /*
376  * Computes Ax = b
377  */
378 static PetscErrorCode MatSolve_MKL_CPARDISO(Mat A, Vec b, Vec x)
379 {
380   Mat_MKL_CPARDISO  *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)(A)->data;
381   PetscScalar       *xarray;
382   const PetscScalar *barray;
383 
384   PetscFunctionBegin;
385   mat_mkl_cpardiso->nrhs = 1;
386   PetscCall(VecGetArray(x, &xarray));
387   PetscCall(VecGetArrayRead(b, &barray));
388 
389   /* solve phase */
390   mat_mkl_cpardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
391   cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja,
392                         mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err);
393 
394   PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL Cluster PARDISO: err=%d, msg = \"%s\". Please check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err));
395 
396   PetscCall(VecRestoreArray(x, &xarray));
397   PetscCall(VecRestoreArrayRead(b, &barray));
398   mat_mkl_cpardiso->CleanUp = PETSC_TRUE;
399   PetscFunctionReturn(PETSC_SUCCESS);
400 }
401 
402 static PetscErrorCode MatSolveTranspose_MKL_CPARDISO(Mat A, Vec b, Vec x)
403 {
404   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data;
405 
406   PetscFunctionBegin;
407 #if defined(PETSC_USE_COMPLEX)
408   mat_mkl_cpardiso->iparm[12 - 1] = 1;
409 #else
410   mat_mkl_cpardiso->iparm[12 - 1] = 2;
411 #endif
412   PetscCall(MatSolve_MKL_CPARDISO(A, b, x));
413   mat_mkl_cpardiso->iparm[12 - 1] = 0;
414   PetscFunctionReturn(PETSC_SUCCESS);
415 }
416 
417 static PetscErrorCode MatMatSolve_MKL_CPARDISO(Mat A, Mat B, Mat X)
418 {
419   Mat_MKL_CPARDISO  *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)(A)->data;
420   PetscScalar       *xarray;
421   const PetscScalar *barray;
422 
423   PetscFunctionBegin;
424   PetscCall(MatGetSize(B, NULL, (PetscInt *)&mat_mkl_cpardiso->nrhs));
425 
426   if (mat_mkl_cpardiso->nrhs > 0) {
427     PetscCall(MatDenseGetArrayRead(B, &barray));
428     PetscCall(MatDenseGetArray(X, &xarray));
429 
430     PetscCheck(barray != xarray, PETSC_COMM_SELF, PETSC_ERR_SUP, "B and X cannot share the same memory location");
431 
432     /* solve phase */
433     mat_mkl_cpardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
434     cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja,
435                           mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err);
436     PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL Cluster PARDISO: err=%d, msg = \"%s\". Please check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err));
437     PetscCall(MatDenseRestoreArrayRead(B, &barray));
438     PetscCall(MatDenseRestoreArray(X, &xarray));
439   }
440   mat_mkl_cpardiso->CleanUp = PETSC_TRUE;
441   PetscFunctionReturn(PETSC_SUCCESS);
442 }
443 
444 /*
445  * LU Decomposition
446  */
447 static PetscErrorCode MatFactorNumeric_MKL_CPARDISO(Mat F, Mat A, const MatFactorInfo *info)
448 {
449   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)(F)->data;
450 
451   PetscFunctionBegin;
452   mat_mkl_cpardiso->matstruc = SAME_NONZERO_PATTERN;
453   PetscCall((*mat_mkl_cpardiso->ConvertToTriples)(A, MAT_REUSE_MATRIX, &mat_mkl_cpardiso->nz, &mat_mkl_cpardiso->ia, &mat_mkl_cpardiso->ja, &mat_mkl_cpardiso->a));
454 
455   mat_mkl_cpardiso->phase = JOB_NUMERICAL_FACTORIZATION;
456   cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja,
457                         mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, &mat_mkl_cpardiso->err);
458   PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL Cluster PARDISO: err=%d, msg = \"%s\". Please check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err));
459 
460   mat_mkl_cpardiso->matstruc = SAME_NONZERO_PATTERN;
461   mat_mkl_cpardiso->CleanUp  = PETSC_TRUE;
462   PetscFunctionReturn(PETSC_SUCCESS);
463 }
464 
465 /* Sets mkl_cpardiso options from the options database */
466 static PetscErrorCode MatSetFromOptions_MKL_CPARDISO(Mat F, Mat A)
467 {
468   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data;
469   PetscInt          icntl, threads;
470   PetscBool         flg;
471 
472   PetscFunctionBegin;
473   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "MKL Cluster PARDISO Options", "Mat");
474   PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_65", "Suggested number of threads to use within MKL Cluster PARDISO", "None", threads, &threads, &flg));
475   if (flg) mkl_set_num_threads((int)threads);
476 
477   PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_66", "Maximum number of factors with identical sparsity structure that must be kept in memory at the same time", "None", mat_mkl_cpardiso->maxfct, &icntl, &flg));
478   if (flg) mat_mkl_cpardiso->maxfct = icntl;
479 
480   PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_67", "Indicates the actual matrix for the solution phase", "None", mat_mkl_cpardiso->mnum, &icntl, &flg));
481   if (flg) mat_mkl_cpardiso->mnum = icntl;
482 
483   PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_68", "Message level information", "None", mat_mkl_cpardiso->msglvl, &icntl, &flg));
484   if (flg) mat_mkl_cpardiso->msglvl = icntl;
485 
486   PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_69", "Defines the matrix type", "None", mat_mkl_cpardiso->mtype, &icntl, &flg));
487   if (flg) mat_mkl_cpardiso->mtype = icntl;
488   PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_1", "Use default values", "None", mat_mkl_cpardiso->iparm[0], &icntl, &flg));
489 
490   if (flg && icntl != 0) {
491     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_2", "Fill-in reducing ordering for the input matrix", "None", mat_mkl_cpardiso->iparm[1], &icntl, &flg));
492     if (flg) mat_mkl_cpardiso->iparm[1] = icntl;
493 
494     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_4", "Preconditioned CGS/CG", "None", mat_mkl_cpardiso->iparm[3], &icntl, &flg));
495     if (flg) mat_mkl_cpardiso->iparm[3] = icntl;
496 
497     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_5", "User permutation", "None", mat_mkl_cpardiso->iparm[4], &icntl, &flg));
498     if (flg) mat_mkl_cpardiso->iparm[4] = icntl;
499 
500     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_6", "Write solution on x", "None", mat_mkl_cpardiso->iparm[5], &icntl, &flg));
501     if (flg) mat_mkl_cpardiso->iparm[5] = icntl;
502 
503     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_8", "Iterative refinement step", "None", mat_mkl_cpardiso->iparm[7], &icntl, &flg));
504     if (flg) mat_mkl_cpardiso->iparm[7] = icntl;
505 
506     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_10", "Pivoting perturbation", "None", mat_mkl_cpardiso->iparm[9], &icntl, &flg));
507     if (flg) mat_mkl_cpardiso->iparm[9] = icntl;
508 
509     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_11", "Scaling vectors", "None", mat_mkl_cpardiso->iparm[10], &icntl, &flg));
510     if (flg) mat_mkl_cpardiso->iparm[10] = icntl;
511 
512     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_12", "Solve with transposed or conjugate transposed matrix A", "None", mat_mkl_cpardiso->iparm[11], &icntl, &flg));
513     if (flg) mat_mkl_cpardiso->iparm[11] = icntl;
514 
515     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_13", "Improved accuracy using (non-) symmetric weighted matching", "None", mat_mkl_cpardiso->iparm[12], &icntl, &flg));
516     if (flg) mat_mkl_cpardiso->iparm[12] = icntl;
517 
518     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_18", "Numbers of non-zero elements", "None", mat_mkl_cpardiso->iparm[17], &icntl, &flg));
519     if (flg) mat_mkl_cpardiso->iparm[17] = icntl;
520 
521     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_19", "Report number of floating point operations", "None", mat_mkl_cpardiso->iparm[18], &icntl, &flg));
522     if (flg) mat_mkl_cpardiso->iparm[18] = icntl;
523 
524     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_21", "Pivoting for symmetric indefinite matrices", "None", mat_mkl_cpardiso->iparm[20], &icntl, &flg));
525     if (flg) mat_mkl_cpardiso->iparm[20] = icntl;
526 
527     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_24", "Parallel factorization control", "None", mat_mkl_cpardiso->iparm[23], &icntl, &flg));
528     if (flg) mat_mkl_cpardiso->iparm[23] = icntl;
529 
530     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_25", "Parallel forward/backward solve control", "None", mat_mkl_cpardiso->iparm[24], &icntl, &flg));
531     if (flg) mat_mkl_cpardiso->iparm[24] = icntl;
532 
533     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_27", "Matrix checker", "None", mat_mkl_cpardiso->iparm[26], &icntl, &flg));
534     if (flg) mat_mkl_cpardiso->iparm[26] = icntl;
535 
536     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_31", "Partial solve and computing selected components of the solution vectors", "None", mat_mkl_cpardiso->iparm[30], &icntl, &flg));
537     if (flg) mat_mkl_cpardiso->iparm[30] = icntl;
538 
539     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_34", "Optimal number of threads for conditional numerical reproducibility (CNR) mode", "None", mat_mkl_cpardiso->iparm[33], &icntl, &flg));
540     if (flg) mat_mkl_cpardiso->iparm[33] = icntl;
541 
542     PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_60", "Intel MKL Cluster PARDISO mode", "None", mat_mkl_cpardiso->iparm[59], &icntl, &flg));
543     if (flg) mat_mkl_cpardiso->iparm[59] = icntl;
544   }
545 
546   PetscOptionsEnd();
547   PetscFunctionReturn(PETSC_SUCCESS);
548 }
549 
550 static PetscErrorCode PetscInitialize_MKL_CPARDISO(Mat A, Mat_MKL_CPARDISO *mat_mkl_cpardiso)
551 {
552   PetscInt    bs;
553   PetscBool   match;
554   PetscMPIInt size;
555   MPI_Comm    comm;
556 
557   PetscFunctionBegin;
558 
559   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)A), &comm));
560   PetscCallMPI(MPI_Comm_size(comm, &size));
561   mat_mkl_cpardiso->comm_mkl_cpardiso = MPI_Comm_c2f(comm);
562 
563   mat_mkl_cpardiso->CleanUp = PETSC_FALSE;
564   mat_mkl_cpardiso->maxfct  = 1;
565   mat_mkl_cpardiso->mnum    = 1;
566   mat_mkl_cpardiso->n       = A->rmap->N;
567   if (mat_mkl_cpardiso->iparm[36]) mat_mkl_cpardiso->n /= mat_mkl_cpardiso->iparm[36];
568   mat_mkl_cpardiso->msglvl = 0;
569   mat_mkl_cpardiso->nrhs   = 1;
570   mat_mkl_cpardiso->err    = 0;
571   mat_mkl_cpardiso->phase  = -1;
572 #if defined(PETSC_USE_COMPLEX)
573   mat_mkl_cpardiso->mtype = 13;
574 #else
575   mat_mkl_cpardiso->mtype         = 11;
576 #endif
577 
578 #if defined(PETSC_USE_REAL_SINGLE)
579   mat_mkl_cpardiso->iparm[27] = 1;
580 #else
581   mat_mkl_cpardiso->iparm[27]     = 0;
582 #endif
583 
584   mat_mkl_cpardiso->iparm[0]  = 1;  /* Solver default parameters overridden with provided by iparm */
585   mat_mkl_cpardiso->iparm[1]  = 2;  /* Use METIS for fill-in reordering */
586   mat_mkl_cpardiso->iparm[5]  = 0;  /* Write solution into x */
587   mat_mkl_cpardiso->iparm[7]  = 2;  /* Max number of iterative refinement steps */
588   mat_mkl_cpardiso->iparm[9]  = 13; /* Perturb the pivot elements with 1E-13 */
589   mat_mkl_cpardiso->iparm[10] = 1;  /* Use nonsymmetric permutation and scaling MPS */
590   mat_mkl_cpardiso->iparm[12] = 1;  /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
591   mat_mkl_cpardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
592   mat_mkl_cpardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
593   mat_mkl_cpardiso->iparm[26] = 1;  /* Check input data for correctness */
594 
595   mat_mkl_cpardiso->iparm[39] = 0;
596   if (size > 1) {
597     mat_mkl_cpardiso->iparm[39] = 2;
598     mat_mkl_cpardiso->iparm[40] = A->rmap->rstart;
599     mat_mkl_cpardiso->iparm[41] = A->rmap->rend - 1;
600   }
601   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &match, MATMPIBAIJ, MATMPISBAIJ, ""));
602   if (match) {
603     PetscCall(MatGetBlockSize(A, &bs));
604     mat_mkl_cpardiso->iparm[36] = bs;
605     mat_mkl_cpardiso->iparm[40] /= bs;
606     mat_mkl_cpardiso->iparm[41] /= bs;
607     mat_mkl_cpardiso->iparm[40]++;
608     mat_mkl_cpardiso->iparm[41]++;
609     mat_mkl_cpardiso->iparm[34] = 0; /* Fortran style */
610   } else {
611     mat_mkl_cpardiso->iparm[34] = 1; /* C style */
612   }
613 
614   mat_mkl_cpardiso->perm = 0;
615   PetscFunctionReturn(PETSC_SUCCESS);
616 }
617 
618 /*
619  * Symbolic decomposition. Mkl_Pardiso analysis phase.
620  */
621 static PetscErrorCode MatLUFactorSymbolic_AIJMKL_CPARDISO(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
622 {
623   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data;
624 
625   PetscFunctionBegin;
626   mat_mkl_cpardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
627 
628   /* Set MKL_CPARDISO options from the options database */
629   PetscCall(MatSetFromOptions_MKL_CPARDISO(F, A));
630   PetscCall((*mat_mkl_cpardiso->ConvertToTriples)(A, MAT_INITIAL_MATRIX, &mat_mkl_cpardiso->nz, &mat_mkl_cpardiso->ia, &mat_mkl_cpardiso->ja, &mat_mkl_cpardiso->a));
631 
632   mat_mkl_cpardiso->n = A->rmap->N;
633   if (mat_mkl_cpardiso->iparm[36]) mat_mkl_cpardiso->n /= mat_mkl_cpardiso->iparm[36];
634 
635   /* analysis phase */
636   mat_mkl_cpardiso->phase = JOB_ANALYSIS;
637 
638   cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja,
639                         mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err);
640 
641   PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL Cluster PARDISO: err=%d, msg = \"%s\".Check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err));
642 
643   mat_mkl_cpardiso->CleanUp = PETSC_TRUE;
644   F->ops->lufactornumeric   = MatFactorNumeric_MKL_CPARDISO;
645   F->ops->solve             = MatSolve_MKL_CPARDISO;
646   F->ops->solvetranspose    = MatSolveTranspose_MKL_CPARDISO;
647   F->ops->matsolve          = MatMatSolve_MKL_CPARDISO;
648   PetscFunctionReturn(PETSC_SUCCESS);
649 }
650 
651 static PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_CPARDISO(Mat F, Mat A, IS perm, const MatFactorInfo *info)
652 {
653   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data;
654 
655   PetscFunctionBegin;
656   mat_mkl_cpardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
657 
658   /* Set MKL_CPARDISO options from the options database */
659   PetscCall(MatSetFromOptions_MKL_CPARDISO(F, A));
660   PetscCall((*mat_mkl_cpardiso->ConvertToTriples)(A, MAT_INITIAL_MATRIX, &mat_mkl_cpardiso->nz, &mat_mkl_cpardiso->ia, &mat_mkl_cpardiso->ja, &mat_mkl_cpardiso->a));
661 
662   mat_mkl_cpardiso->n = A->rmap->N;
663   if (mat_mkl_cpardiso->iparm[36]) mat_mkl_cpardiso->n /= mat_mkl_cpardiso->iparm[36];
664   PetscCheck(!PetscDefined(USE_COMPLEX), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for PARDISO CHOLESKY with complex scalars! Use MAT_FACTOR_LU instead");
665   if (A->spd == PETSC_BOOL3_TRUE) mat_mkl_cpardiso->mtype = 2;
666   else mat_mkl_cpardiso->mtype = -2;
667 
668   /* analysis phase */
669   mat_mkl_cpardiso->phase = JOB_ANALYSIS;
670 
671   cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja,
672                         mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err);
673 
674   PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL Cluster PARDISO: err=%d, msg = \"%s\".Check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err));
675 
676   mat_mkl_cpardiso->CleanUp     = PETSC_TRUE;
677   F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_CPARDISO;
678   F->ops->solve                 = MatSolve_MKL_CPARDISO;
679   F->ops->solvetranspose        = MatSolveTranspose_MKL_CPARDISO;
680   F->ops->matsolve              = MatMatSolve_MKL_CPARDISO;
681   PetscFunctionReturn(PETSC_SUCCESS);
682 }
683 
684 static PetscErrorCode MatView_MKL_CPARDISO(Mat A, PetscViewer viewer)
685 {
686   PetscBool         iascii;
687   PetscViewerFormat format;
688   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data;
689   PetscInt          i;
690 
691   PetscFunctionBegin;
692   /* check if matrix is mkl_cpardiso type */
693   if (A->ops->solve != MatSolve_MKL_CPARDISO) PetscFunctionReturn(PETSC_SUCCESS);
694 
695   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
696   if (iascii) {
697     PetscCall(PetscViewerGetFormat(viewer, &format));
698     if (format == PETSC_VIEWER_ASCII_INFO) {
699       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO run parameters:\n"));
700       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO phase:             %d \n", mat_mkl_cpardiso->phase));
701       for (i = 1; i <= 64; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO iparm[%d]:     %d \n", i, mat_mkl_cpardiso->iparm[i - 1]));
702       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO maxfct:     %d \n", mat_mkl_cpardiso->maxfct));
703       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO mnum:     %d \n", mat_mkl_cpardiso->mnum));
704       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO mtype:     %d \n", mat_mkl_cpardiso->mtype));
705       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO n:     %d \n", mat_mkl_cpardiso->n));
706       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO nrhs:     %d \n", mat_mkl_cpardiso->nrhs));
707       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL Cluster PARDISO msglvl:     %d \n", mat_mkl_cpardiso->msglvl));
708     }
709   }
710   PetscFunctionReturn(PETSC_SUCCESS);
711 }
712 
713 static PetscErrorCode MatGetInfo_MKL_CPARDISO(Mat A, MatInfoType flag, MatInfo *info)
714 {
715   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data;
716 
717   PetscFunctionBegin;
718   info->block_size        = 1.0;
719   info->nz_allocated      = mat_mkl_cpardiso->nz + 0.0;
720   info->nz_unneeded       = 0.0;
721   info->assemblies        = 0.0;
722   info->mallocs           = 0.0;
723   info->memory            = 0.0;
724   info->fill_ratio_given  = 0;
725   info->fill_ratio_needed = 0;
726   info->factor_mallocs    = 0;
727   PetscFunctionReturn(PETSC_SUCCESS);
728 }
729 
730 static PetscErrorCode MatMkl_CPardisoSetCntl_MKL_CPARDISO(Mat F, PetscInt icntl, PetscInt ival)
731 {
732   Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data;
733 
734   PetscFunctionBegin;
735   if (icntl <= 64) {
736     mat_mkl_cpardiso->iparm[icntl - 1] = ival;
737   } else {
738     if (icntl == 65) mkl_set_num_threads((int)ival);
739     else if (icntl == 66) mat_mkl_cpardiso->maxfct = ival;
740     else if (icntl == 67) mat_mkl_cpardiso->mnum = ival;
741     else if (icntl == 68) mat_mkl_cpardiso->msglvl = ival;
742     else if (icntl == 69) mat_mkl_cpardiso->mtype = ival;
743   }
744   PetscFunctionReturn(PETSC_SUCCESS);
745 }
746 
747 /*@
748   MatMkl_CPardisoSetCntl - Set MKL Cluster PARDISO parameters
749   <https://www.intel.com/content/www/us/en/docs/onemkl/developer-reference-c/2023-2/onemkl-pardiso-parallel-direct-sparse-solver-iface.html>
750 
751   Logically Collective
752 
753   Input Parameters:
754 + F     - the factored matrix obtained by calling `MatGetFactor()`
755 . icntl - index of MKL Cluster PARDISO parameter
756 - ival  - value of MKL Cluster PARDISO parameter
757 
758   Options Database Key:
759 . -mat_mkl_cpardiso_<icntl> <ival> - set the option numbered icntl to ival
760 
761   Level: intermediate
762 
763   Note:
764   This routine cannot be used if you are solving the linear system with `TS`, `SNES`, or `KSP`, only if you directly call `MatGetFactor()` so use the options
765   database approach when working with `TS`, `SNES`, or `KSP`. See `MATSOLVERMKL_CPARDISO` for the options
766 
767 .seealso: [](ch_matrices), `Mat`, `MatGetFactor()`, `MATMPIAIJ`, `MATSOLVERMKL_CPARDISO`
768 @*/
769 PetscErrorCode MatMkl_CPardisoSetCntl(Mat F, PetscInt icntl, PetscInt ival)
770 {
771   PetscFunctionBegin;
772   PetscTryMethod(F, "MatMkl_CPardisoSetCntl_C", (Mat, PetscInt, PetscInt), (F, icntl, ival));
773   PetscFunctionReturn(PETSC_SUCCESS);
774 }
775 
776 /*MC
777   MATSOLVERMKL_CPARDISO -  A matrix type providing direct solvers (LU) for parallel matrices via the external package MKL Cluster PARDISO
778   <https://www.intel.com/content/www/us/en/docs/onemkl/developer-reference-c/2023-2/onemkl-pardiso-parallel-direct-sparse-solver-iface.html>
779 
780   Works with `MATMPIAIJ` matrices
781 
782   Use `-pc_type lu` `-pc_factor_mat_solver_type mkl_cpardiso` to use this direct solver
783 
784   Options Database Keys:
785 + -mat_mkl_cpardiso_65 - Suggested number of threads to use within MKL Cluster PARDISO
786 . -mat_mkl_cpardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
787 . -mat_mkl_cpardiso_67 - Indicates the actual matrix for the solution phase
788 . -mat_mkl_cpardiso_68 - Message level information, use 1 to get detailed information on the solver options
789 . -mat_mkl_cpardiso_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
790 . -mat_mkl_cpardiso_1  - Use default values
791 . -mat_mkl_cpardiso_2  - Fill-in reducing ordering for the input matrix
792 . -mat_mkl_cpardiso_4  - Preconditioned CGS/CG
793 . -mat_mkl_cpardiso_5  - User permutation
794 . -mat_mkl_cpardiso_6  - Write solution on x
795 . -mat_mkl_cpardiso_8  - Iterative refinement step
796 . -mat_mkl_cpardiso_10 - Pivoting perturbation
797 . -mat_mkl_cpardiso_11 - Scaling vectors
798 . -mat_mkl_cpardiso_12 - Solve with transposed or conjugate transposed matrix A
799 . -mat_mkl_cpardiso_13 - Improved accuracy using (non-) symmetric weighted matching
800 . -mat_mkl_cpardiso_18 - Numbers of non-zero elements
801 . -mat_mkl_cpardiso_19 - Report number of floating point operations
802 . -mat_mkl_cpardiso_21 - Pivoting for symmetric indefinite matrices
803 . -mat_mkl_cpardiso_24 - Parallel factorization control
804 . -mat_mkl_cpardiso_25 - Parallel forward/backward solve control
805 . -mat_mkl_cpardiso_27 - Matrix checker
806 . -mat_mkl_cpardiso_31 - Partial solve and computing selected components of the solution vectors
807 . -mat_mkl_cpardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
808 - -mat_mkl_cpardiso_60 - Intel MKL Cluster PARDISO mode
809 
810   Level: beginner
811 
812   Notes:
813   Use `-mat_mkl_cpardiso_68 1` to display the number of threads the solver is using. MKL does not provide a way to directly access this
814   information.
815 
816   For more information on the options check
817   <https://www.intel.com/content/www/us/en/docs/onemkl/developer-reference-c/2023-2/onemkl-pardiso-parallel-direct-sparse-solver-iface.html>
818 
819 .seealso: [](ch_matrices), `Mat`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatMkl_CPardisoSetCntl()`, `MatGetFactor()`, `MATSOLVERMKL_PARDISO`
820 M*/
821 
822 static PetscErrorCode MatFactorGetSolverType_mkl_cpardiso(Mat A, MatSolverType *type)
823 {
824   PetscFunctionBegin;
825   *type = MATSOLVERMKL_CPARDISO;
826   PetscFunctionReturn(PETSC_SUCCESS);
827 }
828 
829 /* MatGetFactor for MPI AIJ matrices */
830 static PetscErrorCode MatGetFactor_mpiaij_mkl_cpardiso(Mat A, MatFactorType ftype, Mat *F)
831 {
832   Mat               B;
833   Mat_MKL_CPARDISO *mat_mkl_cpardiso;
834   PetscBool         isSeqAIJ, isMPIBAIJ, isMPISBAIJ;
835 
836   PetscFunctionBegin;
837   /* Create the factorization matrix */
838 
839   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSeqAIJ));
840   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIBAIJ, &isMPIBAIJ));
841   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isMPISBAIJ));
842   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
843   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
844   PetscCall(PetscStrallocpy("mkl_cpardiso", &((PetscObject)B)->type_name));
845   PetscCall(MatSetUp(B));
846 
847   PetscCall(PetscNew(&mat_mkl_cpardiso));
848 
849   if (isSeqAIJ) mat_mkl_cpardiso->ConvertToTriples = MatCopy_seqaij_seqaij_MKL_CPARDISO;
850   else if (isMPIBAIJ) mat_mkl_cpardiso->ConvertToTriples = MatConvertToTriples_mpibaij_mpibaij_MKL_CPARDISO;
851   else if (isMPISBAIJ) mat_mkl_cpardiso->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij_MKL_CPARDISO;
852   else mat_mkl_cpardiso->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij_MKL_CPARDISO;
853 
854   if (ftype == MAT_FACTOR_LU) B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_CPARDISO;
855   else B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_CPARDISO;
856   B->ops->destroy = MatDestroy_MKL_CPARDISO;
857 
858   B->ops->view    = MatView_MKL_CPARDISO;
859   B->ops->getinfo = MatGetInfo_MKL_CPARDISO;
860 
861   B->factortype = ftype;
862   B->assembled  = PETSC_TRUE; /* required by -ksp_view */
863 
864   B->data = mat_mkl_cpardiso;
865 
866   /* set solvertype */
867   PetscCall(PetscFree(B->solvertype));
868   PetscCall(PetscStrallocpy(MATSOLVERMKL_CPARDISO, &B->solvertype));
869 
870   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mkl_cpardiso));
871   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMkl_CPardisoSetCntl_C", MatMkl_CPardisoSetCntl_MKL_CPARDISO));
872   PetscCall(PetscInitialize_MKL_CPARDISO(A, mat_mkl_cpardiso));
873 
874   *F = B;
875   PetscFunctionReturn(PETSC_SUCCESS);
876 }
877 
878 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MKL_CPardiso(void)
879 {
880   PetscFunctionBegin;
881   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_mkl_cpardiso));
882   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_mkl_cpardiso));
883   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATMPIBAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_mkl_cpardiso));
884   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATMPISBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_mpiaij_mkl_cpardiso));
885   PetscFunctionReturn(PETSC_SUCCESS);
886 }
887