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