xref: /petsc/src/mat/impls/aij/seq/mkl_pardiso/mkl_pardiso.c (revision 4f02bc6a7df4e4b03c783003a8e0899ee35fcb05)
1 #if defined(PETSC_HAVE_LIBMKL_INTEL_ILP64)
2 #define MKL_ILP64
3 #endif
4 
5 #include <../src/mat/impls/aij/seq/aij.h>        /*I "petscmat.h" I*/
6 #include <../src/mat/impls/sbaij/seq/sbaij.h>
7 #include <../src/mat/impls/dense/seq/dense.h>
8 #include <petscblaslapack.h>
9 
10 #include <stdio.h>
11 #include <stdlib.h>
12 #include <math.h>
13 #include <mkl_pardiso.h>
14 
15 PETSC_EXTERN void PetscSetMKL_PARDISOThreads(int);
16 
17 /*
18  *  Possible mkl_pardiso phases that controls the execution of the solver.
19  *  For more information check mkl_pardiso manual.
20  */
21 #define JOB_ANALYSIS 11
22 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION 12
23 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13
24 #define JOB_NUMERICAL_FACTORIZATION 22
25 #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 23
26 #define JOB_SOLVE_ITERATIVE_REFINEMENT 33
27 #define JOB_SOLVE_FORWARD_SUBSTITUTION 331
28 #define JOB_SOLVE_DIAGONAL_SUBSTITUTION 332
29 #define JOB_SOLVE_BACKWARD_SUBSTITUTION 333
30 #define JOB_RELEASE_OF_LU_MEMORY 0
31 #define JOB_RELEASE_OF_ALL_MEMORY -1
32 
33 #define IPARM_SIZE 64
34 
35 #if defined(PETSC_USE_64BIT_INDICES)
36  #if defined(PETSC_HAVE_LIBMKL_INTEL_ILP64)
37   /* sizeof(MKL_INT) == sizeof(long long int) if ilp64*/
38   #define INT_TYPE long long int
39   #define MKL_PARDISO pardiso
40   #define MKL_PARDISO_INIT pardisoinit
41  #else
42   #define INT_TYPE long long int
43   #define MKL_PARDISO pardiso_64
44   #define MKL_PARDISO_INIT pardiso_64init
45  #endif
46 #else
47  #define INT_TYPE int
48  #define MKL_PARDISO pardiso
49  #define MKL_PARDISO_INIT pardisoinit
50 #endif
51 
52 
53 /*
54  *  Internal data structure.
55  *  For more information check mkl_pardiso manual.
56  */
57 typedef struct {
58 
59   /* Configuration vector*/
60   INT_TYPE     iparm[IPARM_SIZE];
61 
62   /*
63    * Internal mkl_pardiso memory location.
64    * After the first call to mkl_pardiso do not modify pt, as that could cause a serious memory leak.
65    */
66   void         *pt[IPARM_SIZE];
67 
68   /* Basic mkl_pardiso info*/
69   INT_TYPE     phase, maxfct, mnum, mtype, n, nrhs, msglvl, err;
70 
71   /* Matrix structure*/
72   void         *a;
73   INT_TYPE     *ia, *ja;
74 
75   /* Number of non-zero elements*/
76   INT_TYPE     nz;
77 
78   /* Row permutaton vector*/
79   INT_TYPE     *perm;
80 
81   /* Define if matrix preserves sparse structure.*/
82   MatStructure matstruc;
83 
84   PetscBool    needsym;
85   PetscBool    freeaij;
86 
87   /* Schur complement */
88   PetscScalar  *schur;
89   PetscInt     schur_size;
90   PetscInt     *schur_idxs;
91   PetscScalar  *schur_work;
92   PetscBLASInt schur_work_size;
93   PetscInt     schur_solver_type;
94   PetscInt     *schur_pivots;
95   PetscBool    schur_factored;
96   PetscBool    schur_inverted;
97   PetscBool    solve_interior;
98 
99   /* True if mkl_pardiso function have been used.*/
100   PetscBool CleanUp;
101 
102   /* Conversion to a format suitable for MKL */
103   PetscErrorCode (*Convert)(Mat, PetscBool, MatReuse, PetscBool*, INT_TYPE*, INT_TYPE**, INT_TYPE**, void**);
104   PetscErrorCode (*Destroy)(Mat);
105 } Mat_MKL_PARDISO;
106 
107 #undef __FUNCT__
108 #define __FUNCT__ "MatMKLPardiso_Convert_seqsbaij"
109 PetscErrorCode MatMKLPardiso_Convert_seqsbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,void **v)
110 {
111   Mat_SeqSBAIJ   *aa = (Mat_SeqSBAIJ*)A->data;
112   PetscInt       bs  = A->rmap->bs;
113 
114   PetscFunctionBegin;
115   if (!sym) {
116     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
117   }
118   if (bs == 1) { /* already in the correct format */
119     *v    = aa->a;
120     *r    = aa->i;
121     *c    = aa->j;
122     *nnz  = aa->nz;
123     *free = PETSC_FALSE;
124   } else {
125     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Conversion from SeqSBAIJ to MKL Pardiso format still need to be implemented");
126   }
127   PetscFunctionReturn(0);
128 }
129 
130 #undef __FUNCT__
131 #define __FUNCT__ "MatMKLPardiso_Convert_seqbaij"
132 PetscErrorCode MatMKLPardiso_Convert_seqbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,void **v)
133 {
134   PetscFunctionBegin;
135   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Conversion from SeqBAIJ to MKL Pardiso format still need to be implemented");
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "MatMKLPardiso_Convert_seqaij"
141 PetscErrorCode MatMKLPardiso_Convert_seqaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,void **v)
142 {
143   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)A->data;
144   PetscErrorCode ierr;
145 
146   PetscFunctionBegin;
147   if (!sym) { /* already in the correct format */
148     *v    = aa->a;
149     *r    = aa->i;
150     *c    = aa->j;
151     *nnz  = aa->nz;
152     *free = PETSC_FALSE;
153     PetscFunctionReturn(0);
154   }
155   /* need to get the triangular part */
156   if (reuse == MAT_INITIAL_MATRIX) {
157     PetscScalar *vals,*vv;
158     PetscInt    *row,*col,*jj;
159     PetscInt    m = A->rmap->n,nz,i;
160 
161     nz = 0;
162     for (i=0; i<m; i++) {
163       nz += aa->i[i+1] - aa->diag[i];
164     }
165     ierr = PetscMalloc2(m+1,&row,nz,&col);CHKERRQ(ierr);
166     ierr = PetscMalloc1(nz,&vals);CHKERRQ(ierr);
167     jj = col;
168     vv = vals;
169 
170     row[0] = 0;
171     for (i=0; i<m; i++) {
172       PetscInt    *aj = aa->j + aa->diag[i];
173       PetscScalar *av = aa->a + aa->diag[i];
174       PetscInt    rl = aa->i[i+1] - aa->diag[i],j;
175       for (j=0; j<rl; j++) {
176         *jj = *aj; jj++; aj++;
177         *vv = *av; vv++; av++;
178       }
179       row[i+1]    = row[i] + rl;
180     }
181     *v    = vals;
182     *r    = row;
183     *c    = col;
184     *nnz  = nz;
185     *free = PETSC_TRUE;
186   } else {
187     PetscScalar *vv;
188     PetscInt    m = A->rmap->n,i;
189 
190     vv = *v;
191     for (i=0; i<m; i++) {
192       PetscScalar *av = aa->a + aa->diag[i];
193       PetscInt    rl = aa->i[i+1] - aa->diag[i],j;
194       for (j=0; j<rl; j++) {
195         *vv = *av; vv++; av++;
196       }
197     }
198     *free = PETSC_TRUE;
199   }
200   PetscFunctionReturn(0);
201 }
202 
203 void pardiso_64init(void *pt, INT_TYPE *mtype, INT_TYPE iparm [])
204 {
205   int iparm_copy[IPARM_SIZE], mtype_copy, i;
206 
207   mtype_copy = *mtype;
208   pardisoinit(pt, &mtype_copy, iparm_copy);
209   for(i = 0; i < IPARM_SIZE; i++){
210     iparm[i] = iparm_copy[i];
211   }
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "MatMKLPardisoFactorSchur_Private"
216 static PetscErrorCode MatMKLPardisoFactorSchur_Private(Mat_MKL_PARDISO* mpardiso)
217 {
218   PetscBLASInt   B_N,B_ierr;
219   PetscScalar    *work,val;
220   PetscBLASInt   lwork = -1;
221   PetscErrorCode ierr;
222 
223   PetscFunctionBegin;
224   if (mpardiso->schur_factored) {
225     PetscFunctionReturn(0);
226   }
227   ierr = PetscBLASIntCast(mpardiso->schur_size,&B_N);CHKERRQ(ierr);
228   switch (mpardiso->schur_solver_type) {
229     case 1: /* hermitian solver */
230       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
231       PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,mpardiso->schur,&B_N,&B_ierr));
232       ierr = PetscFPTrapPop();CHKERRQ(ierr);
233       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr);
234       break;
235     case 2: /* symmetric */
236       if (!mpardiso->schur_pivots) {
237         ierr = PetscMalloc1(B_N,&mpardiso->schur_pivots);CHKERRQ(ierr);
238       }
239       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
240       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,&val,&lwork,&B_ierr));
241       ierr = PetscBLASIntCast((PetscInt)val,&lwork);CHKERRQ(ierr);
242       ierr = PetscMalloc1(lwork,&work);CHKERRQ(ierr);
243       PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,work,&lwork,&B_ierr));
244       ierr = PetscFree(work);CHKERRQ(ierr);
245       ierr = PetscFPTrapPop();CHKERRQ(ierr);
246       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr);
247       break;
248     default: /* general */
249       if (!mpardiso->schur_pivots) {
250         ierr = PetscMalloc1(B_N,&mpardiso->schur_pivots);CHKERRQ(ierr);
251       }
252       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
253       PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,&B_ierr));
254       ierr = PetscFPTrapPop();CHKERRQ(ierr);
255       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr);
256       break;
257   }
258   mpardiso->schur_factored = PETSC_TRUE;
259   PetscFunctionReturn(0);
260 }
261 
262 #undef __FUNCT__
263 #define __FUNCT__ "MatMKLPardisoInvertSchur_Private"
264 static PetscErrorCode MatMKLPardisoInvertSchur_Private(Mat_MKL_PARDISO* mpardiso)
265 {
266   PetscBLASInt   B_N,B_ierr;
267   PetscErrorCode ierr;
268 
269   PetscFunctionBegin;
270   ierr = MatMKLPardisoFactorSchur_Private(mpardiso);CHKERRQ(ierr);
271   ierr = PetscBLASIntCast(mpardiso->schur_size,&B_N);CHKERRQ(ierr);
272   switch (mpardiso->schur_solver_type) {
273     case 1: /* hermitian solver */
274       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
275       PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,mpardiso->schur,&B_N,&B_ierr));
276       ierr = PetscFPTrapPop();CHKERRQ(ierr);
277       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr);
278       break;
279     case 2: /* symmetric */
280       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
281       PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,mpardiso->schur_work,&B_ierr));
282       ierr = PetscFPTrapPop();CHKERRQ(ierr);
283       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr);
284       break;
285     default: /* general */
286       ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
287       PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,mpardiso->schur_work,&mpardiso->schur_work_size,&B_ierr));
288       ierr = PetscFPTrapPop();CHKERRQ(ierr);
289       if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr);
290       break;
291   }
292   mpardiso->schur_inverted = PETSC_TRUE;
293   PetscFunctionReturn(0);
294 }
295 
296 #undef __FUNCT__
297 #define __FUNCT__ "MatMKLPardisoSolveSchur_Private"
298 static PetscErrorCode MatMKLPardisoSolveSchur_Private(Mat_MKL_PARDISO* mpardiso, PetscScalar *B, PetscScalar *X)
299 {
300   PetscScalar    one=1.,zero=0.,*schur_rhs,*schur_sol;
301   PetscBLASInt   B_N,B_Nrhs,B_ierr;
302   char           type[2];
303   PetscErrorCode ierr;
304 
305   PetscFunctionBegin;
306   ierr = MatMKLPardisoFactorSchur_Private(mpardiso);CHKERRQ(ierr);
307   ierr = PetscBLASIntCast(mpardiso->schur_size,&B_N);CHKERRQ(ierr);
308   ierr = PetscBLASIntCast(mpardiso->nrhs,&B_Nrhs);CHKERRQ(ierr);
309   if (X == B && mpardiso->schur_inverted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"X and B cannot point to the same address");
310   if (X != B) { /* using LAPACK *TRS subroutines */
311     ierr = PetscMemcpy(X,B,B_N*B_Nrhs*sizeof(PetscScalar));CHKERRQ(ierr);
312   }
313   schur_rhs = B;
314   schur_sol = X;
315   switch (mpardiso->schur_solver_type) {
316     case 1: /* hermitian solver */
317       if (mpardiso->schur_inverted) { /* BLAShemm should go here */
318         PetscStackCallBLAS("BLASsymm",BLASsymm_("L","L",&B_N,&B_Nrhs,&one,mpardiso->schur,&B_N,schur_rhs,&B_N,&zero,schur_sol,&B_N));
319       } else {
320         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
321         PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&B_N,&B_Nrhs,mpardiso->schur,&B_N,schur_sol,&B_N,&B_ierr));
322         ierr = PetscFPTrapPop();CHKERRQ(ierr);
323         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr);
324       }
325       break;
326     case 2: /* symmetric solver */
327       if (mpardiso->schur_inverted) {
328         PetscStackCallBLAS("BLASsymm",BLASsymm_("L","L",&B_N,&B_Nrhs,&one,mpardiso->schur,&B_N,schur_rhs,&B_N,&zero,schur_sol,&B_N));
329       } else {
330         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
331         PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&B_N,&B_Nrhs,mpardiso->schur,&B_N,mpardiso->schur_pivots,schur_sol,&B_N,&B_ierr));
332         ierr = PetscFPTrapPop();CHKERRQ(ierr);
333         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr);
334       }
335       break;
336     default: /* general */
337       switch (mpardiso->iparm[12-1]) {
338         case 1:
339           sprintf(type,"C");
340           break;
341         case 2:
342           sprintf(type,"T");
343           break;
344         default:
345           sprintf(type,"N");
346           break;
347       }
348       if (mpardiso->schur_inverted) {
349         PetscStackCallBLAS("BLASgemm",BLASgemm_(type,"N",&B_N,&B_Nrhs,&B_N,&one,mpardiso->schur,&B_N,schur_rhs,&B_N,&zero,schur_sol,&B_N));
350       } else {
351         ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
352         PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_(type,&B_N,&B_Nrhs,mpardiso->schur,&B_N,mpardiso->schur_pivots,schur_sol,&B_N,&B_ierr));
353         ierr = PetscFPTrapPop();CHKERRQ(ierr);
354         if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr);
355       }
356       break;
357   }
358   PetscFunctionReturn(0);
359 }
360 
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "MatFactorSetSchurIS_MKL_PARDISO"
364 PetscErrorCode MatFactorSetSchurIS_MKL_PARDISO(Mat F, IS is)
365 {
366   Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr;
367   const PetscInt  *idxs;
368   PetscInt        size,i;
369   PetscMPIInt     csize;
370   PetscBool       sorted;
371   PetscErrorCode  ierr;
372 
373   PetscFunctionBegin;
374   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&csize);CHKERRQ(ierr);
375   if (csize > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MKL_PARDISO parallel Schur complements not yet supported from PETSc\n");
376   ierr = ISSorted(is,&sorted);CHKERRQ(ierr);
377   if (!sorted) {
378     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IS for MKL_PARDISO Schur complements needs to be sorted\n");
379   }
380   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
381   if (mpardiso->schur_size != size) {
382     mpardiso->schur_size = size;
383     ierr = PetscFree2(mpardiso->schur,mpardiso->schur_work);CHKERRQ(ierr);
384     ierr = PetscFree(mpardiso->schur_idxs);CHKERRQ(ierr);
385     ierr = PetscFree(mpardiso->schur_pivots);CHKERRQ(ierr);
386     ierr = PetscBLASIntCast(PetscMax(mpardiso->n,2*size),&mpardiso->schur_work_size);CHKERRQ(ierr);
387     ierr = PetscMalloc2(size*size,&mpardiso->schur,mpardiso->schur_work_size,&mpardiso->schur_work);CHKERRQ(ierr);
388     ierr = PetscMalloc1(size,&mpardiso->schur_idxs);CHKERRQ(ierr);
389   }
390   ierr = PetscMemzero(mpardiso->perm,mpardiso->n*sizeof(INT_TYPE));CHKERRQ(ierr);
391   ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr);
392   ierr = PetscMemcpy(mpardiso->schur_idxs,idxs,size*sizeof(PetscInt));CHKERRQ(ierr);
393   for (i=0;i<size;i++) mpardiso->perm[idxs[i]] = 1;
394   ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr);
395   if (size) { /* turn on Schur switch if we the set of indices is not empty */
396     mpardiso->iparm[36-1] = 2;
397   }
398   mpardiso->schur_factored = PETSC_FALSE;
399   mpardiso->schur_inverted = PETSC_FALSE;
400   PetscFunctionReturn(0);
401 }
402 
403 #undef __FUNCT__
404 #define __FUNCT__ "MatFactorCreateSchurComplement_MKL_PARDISO"
405 PetscErrorCode MatFactorCreateSchurComplement_MKL_PARDISO(Mat F,Mat* S)
406 {
407   Mat             St;
408   Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr;
409   PetscScalar     *array;
410   PetscErrorCode  ierr;
411 
412   PetscFunctionBegin;
413   if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
414   else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
415 
416   ierr = MatCreate(PetscObjectComm((PetscObject)F),&St);CHKERRQ(ierr);
417   ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mpardiso->schur_size,mpardiso->schur_size);CHKERRQ(ierr);
418   ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr);
419   ierr = MatSetUp(St);CHKERRQ(ierr);
420   ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr);
421   ierr = PetscMemcpy(array,mpardiso->schur,mpardiso->schur_size*mpardiso->schur_size*sizeof(PetscScalar));CHKERRQ(ierr);
422   ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr);
423   *S = St;
424   PetscFunctionReturn(0);
425 }
426 
427 #undef __FUNCT__
428 #define __FUNCT__ "MatFactorGetSchurComplement_MKL_PARDISO"
429 PetscErrorCode MatFactorGetSchurComplement_MKL_PARDISO(Mat F,Mat* S)
430 {
431   Mat             St;
432   Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr;
433   PetscErrorCode  ierr;
434 
435   PetscFunctionBegin;
436   if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
437   else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
438 
439   ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)F),mpardiso->schur_size,mpardiso->schur_size,mpardiso->schur,&St);CHKERRQ(ierr);
440   *S = St;
441   PetscFunctionReturn(0);
442 }
443 
444 #undef __FUNCT__
445 #define __FUNCT__ "MatFactorInvertSchurComplement_MKL_PARDISO"
446 PetscErrorCode MatFactorInvertSchurComplement_MKL_PARDISO(Mat F)
447 {
448   Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr;
449   PetscErrorCode  ierr;
450 
451   PetscFunctionBegin;
452   if (!mpardiso->iparm[36-1]) { /* do nothing */
453     PetscFunctionReturn(0);
454   }
455   if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
456   ierr = MatMKLPardisoInvertSchur_Private(mpardiso);CHKERRQ(ierr);
457   PetscFunctionReturn(0);
458 }
459 
460 #undef __FUNCT__
461 #define __FUNCT__ "MatFactorSolveSchurComplement_MKL_PARDISO"
462 PetscErrorCode MatFactorSolveSchurComplement_MKL_PARDISO(Mat F, Vec rhs, Vec sol)
463 {
464   Mat_MKL_PARDISO   *mpardiso =(Mat_MKL_PARDISO*)F->spptr;
465   PetscScalar       *asol,*arhs;
466   PetscErrorCode ierr;
467 
468   PetscFunctionBegin;
469   if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
470   else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
471 
472   mpardiso->nrhs = 1;
473   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr);
474   ierr = VecGetArray(sol,&asol);CHKERRQ(ierr);
475   ierr = MatMKLPardisoSolveSchur_Private(mpardiso,arhs,asol);CHKERRQ(ierr);
476   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr);
477   ierr = VecRestoreArray(sol,&asol);CHKERRQ(ierr);
478   PetscFunctionReturn(0);
479 }
480 
481 #undef __FUNCT__
482 #define __FUNCT__ "MatFactorSolveSchurComplementTranspose_MKL_PARDISO"
483 PetscErrorCode MatFactorSolveSchurComplementTranspose_MKL_PARDISO(Mat F, Vec rhs, Vec sol)
484 {
485   Mat_MKL_PARDISO   *mpardiso =(Mat_MKL_PARDISO*)F->spptr;
486   PetscScalar       *asol,*arhs;
487   PetscErrorCode ierr;
488 
489   PetscFunctionBegin;
490   if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it");
491   else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before");
492 
493   mpardiso->nrhs = 1;
494   ierr = VecGetArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr);
495   ierr = VecGetArray(sol,&asol);CHKERRQ(ierr);
496   mpardiso->iparm[12 - 1] = 2;
497   ierr = MatMKLPardisoSolveSchur_Private(mpardiso,arhs,asol);CHKERRQ(ierr);
498   mpardiso->iparm[12 - 1] = 0;
499   ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr);
500   ierr = VecRestoreArray(sol,&asol);CHKERRQ(ierr);
501   PetscFunctionReturn(0);
502 }
503 
504 #undef __FUNCT__
505 #define __FUNCT__ "MatDestroy_MKL_PARDISO"
506 PetscErrorCode MatDestroy_MKL_PARDISO(Mat A)
507 {
508   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr;
509   PetscErrorCode  ierr;
510 
511   PetscFunctionBegin;
512   if (mat_mkl_pardiso->CleanUp) {
513     mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;
514 
515     MKL_PARDISO (mat_mkl_pardiso->pt,
516       &mat_mkl_pardiso->maxfct,
517       &mat_mkl_pardiso->mnum,
518       &mat_mkl_pardiso->mtype,
519       &mat_mkl_pardiso->phase,
520       &mat_mkl_pardiso->n,
521       NULL,
522       NULL,
523       NULL,
524       mat_mkl_pardiso->perm,
525       &mat_mkl_pardiso->nrhs,
526       mat_mkl_pardiso->iparm,
527       &mat_mkl_pardiso->msglvl,
528       NULL,
529       NULL,
530       &mat_mkl_pardiso->err);
531   }
532   ierr = PetscFree(mat_mkl_pardiso->perm);CHKERRQ(ierr);
533   ierr = PetscFree2(mat_mkl_pardiso->schur,mat_mkl_pardiso->schur_work);CHKERRQ(ierr);
534   ierr = PetscFree(mat_mkl_pardiso->schur_idxs);CHKERRQ(ierr);
535   ierr = PetscFree(mat_mkl_pardiso->schur_pivots);CHKERRQ(ierr);
536   if (mat_mkl_pardiso->freeaij) {
537     ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr);
538     ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr);
539   }
540   if (mat_mkl_pardiso->Destroy) {
541     ierr = (mat_mkl_pardiso->Destroy)(A);CHKERRQ(ierr);
542   }
543   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
544 
545   /* clear composed functions */
546   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
547   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr);
548   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr);
549   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);CHKERRQ(ierr);
550   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);CHKERRQ(ierr);
551   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);CHKERRQ(ierr);
552   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);CHKERRQ(ierr);
553   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);CHKERRQ(ierr);
554   PetscFunctionReturn(0);
555 }
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "MatMKLPardisoScatterSchur_Private"
559 static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce)
560 {
561   PetscFunctionBegin;
562   if (reduce) { /* data given for the whole matrix */
563     PetscInt i,m=0,p=0;
564     for (i=0;i<mpardiso->nrhs;i++) {
565       PetscInt j;
566       for (j=0;j<mpardiso->schur_size;j++) {
567         schur[p+j] = whole[m+mpardiso->schur_idxs[j]];
568       }
569       m += mpardiso->n;
570       p += mpardiso->schur_size;
571     }
572   } else { /* from Schur to whole */
573     PetscInt i,m=0,p=0;
574     for (i=0;i<mpardiso->nrhs;i++) {
575       PetscInt j;
576       for (j=0;j<mpardiso->schur_size;j++) {
577         whole[m+mpardiso->schur_idxs[j]] = schur[p+j];
578       }
579       m += mpardiso->n;
580       p += mpardiso->schur_size;
581     }
582   }
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "MatSolve_MKL_PARDISO"
588 PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x)
589 {
590   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr;
591   PetscErrorCode    ierr;
592   PetscScalar       *xarray;
593   const PetscScalar *barray;
594 
595   PetscFunctionBegin;
596   mat_mkl_pardiso->nrhs = 1;
597   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
598   ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
599 
600   if (!mat_mkl_pardiso->schur) {
601     mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
602   } else {
603     mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
604   }
605   mat_mkl_pardiso->iparm[6-1] = 0;
606 
607   MKL_PARDISO (mat_mkl_pardiso->pt,
608     &mat_mkl_pardiso->maxfct,
609     &mat_mkl_pardiso->mnum,
610     &mat_mkl_pardiso->mtype,
611     &mat_mkl_pardiso->phase,
612     &mat_mkl_pardiso->n,
613     mat_mkl_pardiso->a,
614     mat_mkl_pardiso->ia,
615     mat_mkl_pardiso->ja,
616     mat_mkl_pardiso->perm,
617     &mat_mkl_pardiso->nrhs,
618     mat_mkl_pardiso->iparm,
619     &mat_mkl_pardiso->msglvl,
620     (void*)barray,
621     (void*)xarray,
622     &mat_mkl_pardiso->err);
623 
624   if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err);
625 
626   if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
627     PetscInt shift = mat_mkl_pardiso->schur_size;
628 
629     /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
630     if (!mat_mkl_pardiso->schur_inverted) {
631       shift = 0;
632     }
633 
634     if (!mat_mkl_pardiso->solve_interior) {
635       /* solve Schur complement */
636       ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr);
637       ierr = MatMKLPardisoSolveSchur_Private(mat_mkl_pardiso,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr);
638       ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr);
639     } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substitued to xarray[schur] will be neglected */
640       PetscInt i;
641       for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
642         xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.;
643       }
644     }
645     /* expansion phase */
646     mat_mkl_pardiso->iparm[6-1] = 1;
647     mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
648     MKL_PARDISO (mat_mkl_pardiso->pt,
649       &mat_mkl_pardiso->maxfct,
650       &mat_mkl_pardiso->mnum,
651       &mat_mkl_pardiso->mtype,
652       &mat_mkl_pardiso->phase,
653       &mat_mkl_pardiso->n,
654       mat_mkl_pardiso->a,
655       mat_mkl_pardiso->ia,
656       mat_mkl_pardiso->ja,
657       mat_mkl_pardiso->perm,
658       &mat_mkl_pardiso->nrhs,
659       mat_mkl_pardiso->iparm,
660       &mat_mkl_pardiso->msglvl,
661       (void*)xarray,
662       (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
663       &mat_mkl_pardiso->err);
664 
665     if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err);
666     mat_mkl_pardiso->iparm[6-1] = 0;
667   }
668   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
669   ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
670   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "MatSolveTranspose_MKL_PARDISO"
676 PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x)
677 {
678   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr;
679   PetscErrorCode  ierr;
680 
681   PetscFunctionBegin;
682 #if defined(PETSC_USE_COMPLEX)
683   mat_mkl_pardiso->iparm[12 - 1] = 1;
684 #else
685   mat_mkl_pardiso->iparm[12 - 1] = 2;
686 #endif
687   ierr = MatSolve_MKL_PARDISO(A,b,x);CHKERRQ(ierr);
688   mat_mkl_pardiso->iparm[12 - 1] = 0;
689   PetscFunctionReturn(0);
690 }
691 
692 #undef __FUNCT__
693 #define __FUNCT__ "MatMatSolve_MKL_PARDISO"
694 PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X)
695 {
696   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr;
697   PetscErrorCode    ierr;
698   PetscScalar       *barray, *xarray;
699   PetscBool         flg;
700 
701   PetscFunctionBegin;
702   ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr);
703   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix");
704   ierr = PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr);
705   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix");
706 
707   ierr = MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);CHKERRQ(ierr);
708 
709   if(mat_mkl_pardiso->nrhs > 0){
710     ierr = MatDenseGetArray(B,&barray);
711     ierr = MatDenseGetArray(X,&xarray);
712 
713     if (!mat_mkl_pardiso->schur) {
714       mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
715     } else {
716       mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
717     }
718     mat_mkl_pardiso->iparm[6-1] = 0;
719 
720     MKL_PARDISO (mat_mkl_pardiso->pt,
721       &mat_mkl_pardiso->maxfct,
722       &mat_mkl_pardiso->mnum,
723       &mat_mkl_pardiso->mtype,
724       &mat_mkl_pardiso->phase,
725       &mat_mkl_pardiso->n,
726       mat_mkl_pardiso->a,
727       mat_mkl_pardiso->ia,
728       mat_mkl_pardiso->ja,
729       mat_mkl_pardiso->perm,
730       &mat_mkl_pardiso->nrhs,
731       mat_mkl_pardiso->iparm,
732       &mat_mkl_pardiso->msglvl,
733       (void*)barray,
734       (void*)xarray,
735       &mat_mkl_pardiso->err);
736     if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err);
737 
738     if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
739       PetscScalar *o_schur_work = NULL;
740       PetscInt    shift = mat_mkl_pardiso->schur_size*mat_mkl_pardiso->nrhs,scale;
741       PetscInt    mem = mat_mkl_pardiso->n*mat_mkl_pardiso->nrhs;
742 
743       /* allocate extra memory if it is needed */
744       scale = 1;
745       if (mat_mkl_pardiso->schur_inverted) {
746         scale = 2;
747       }
748       mem *= scale;
749       if (mem > mat_mkl_pardiso->schur_work_size) {
750         o_schur_work = mat_mkl_pardiso->schur_work;
751         ierr = PetscMalloc1(mem,&mat_mkl_pardiso->schur_work);CHKERRQ(ierr);
752       }
753 
754       /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
755       if (!mat_mkl_pardiso->schur_inverted) {
756         shift = 0;
757       }
758 
759       /* solve Schur complement */
760       if (!mat_mkl_pardiso->solve_interior) {
761         ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr);
762         ierr = MatMKLPardisoSolveSchur_Private(mat_mkl_pardiso,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr);
763         ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr);
764       } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substitued to xarray[schur,n] will be neglected */
765         PetscInt i,n,m=0;
766         for (n=0;n<mat_mkl_pardiso->nrhs;n++) {
767           for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
768             xarray[mat_mkl_pardiso->schur_idxs[i]+m] = 0.;
769           }
770           m += mat_mkl_pardiso->n;
771         }
772       }
773       /* expansion phase */
774       mat_mkl_pardiso->iparm[6-1] = 1;
775       mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
776       MKL_PARDISO (mat_mkl_pardiso->pt,
777         &mat_mkl_pardiso->maxfct,
778         &mat_mkl_pardiso->mnum,
779         &mat_mkl_pardiso->mtype,
780         &mat_mkl_pardiso->phase,
781         &mat_mkl_pardiso->n,
782         mat_mkl_pardiso->a,
783         mat_mkl_pardiso->ia,
784         mat_mkl_pardiso->ja,
785         mat_mkl_pardiso->perm,
786         &mat_mkl_pardiso->nrhs,
787         mat_mkl_pardiso->iparm,
788         &mat_mkl_pardiso->msglvl,
789         (void*)xarray,
790         (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
791         &mat_mkl_pardiso->err);
792       if (o_schur_work) { /* restore original schur_work (minimal size) */
793         ierr = PetscFree(mat_mkl_pardiso->schur_work);CHKERRQ(ierr);
794         mat_mkl_pardiso->schur_work = o_schur_work;
795       }
796       if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err);
797       mat_mkl_pardiso->iparm[6-1] = 0;
798     }
799   }
800   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
801   PetscFunctionReturn(0);
802 }
803 
804 #undef __FUNCT__
805 #define __FUNCT__ "MatFactorNumeric_MKL_PARDISO"
806 PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info)
807 {
808   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->spptr;
809   PetscErrorCode  ierr;
810 
811   PetscFunctionBegin;
812   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
813   ierr = (*mat_mkl_pardiso->Convert)(A,mat_mkl_pardiso->needsym,MAT_REUSE_MATRIX,&mat_mkl_pardiso->freeaij,&mat_mkl_pardiso->nz,&mat_mkl_pardiso->ia,&mat_mkl_pardiso->ja,&mat_mkl_pardiso->a);CHKERRQ(ierr);
814 
815   mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION;
816   MKL_PARDISO (mat_mkl_pardiso->pt,
817     &mat_mkl_pardiso->maxfct,
818     &mat_mkl_pardiso->mnum,
819     &mat_mkl_pardiso->mtype,
820     &mat_mkl_pardiso->phase,
821     &mat_mkl_pardiso->n,
822     mat_mkl_pardiso->a,
823     mat_mkl_pardiso->ia,
824     mat_mkl_pardiso->ja,
825     mat_mkl_pardiso->perm,
826     &mat_mkl_pardiso->nrhs,
827     mat_mkl_pardiso->iparm,
828     &mat_mkl_pardiso->msglvl,
829     NULL,
830     (void*)mat_mkl_pardiso->schur,
831     &mat_mkl_pardiso->err);
832   if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err);
833 
834   if (mat_mkl_pardiso->schur) { /* schur output from pardiso is in row major format */
835     PetscInt j,k,n=mat_mkl_pardiso->schur_size;
836     if (!mat_mkl_pardiso->schur_solver_type) {
837       for (j=0; j<n; j++) {
838         for (k=0; k<j; k++) {
839           PetscScalar tmp = mat_mkl_pardiso->schur[j + k*n];
840           mat_mkl_pardiso->schur[j + k*n] = mat_mkl_pardiso->schur[k + j*n];
841           mat_mkl_pardiso->schur[k + j*n] = tmp;
842         }
843       }
844     } else { /* we could use row-major in LAPACK routines (e.g. use 'U' instead of 'L'; instead, I prefer consistency between data structures and swap to column major */
845       for (j=0; j<n; j++) {
846         for (k=0; k<j; k++) {
847           mat_mkl_pardiso->schur[j + k*n] = mat_mkl_pardiso->schur[k + j*n];
848         }
849       }
850     }
851   }
852   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
853   mat_mkl_pardiso->CleanUp  = PETSC_TRUE;
854   mat_mkl_pardiso->schur_factored = PETSC_FALSE;
855   mat_mkl_pardiso->schur_inverted = PETSC_FALSE;
856   PetscFunctionReturn(0);
857 }
858 
859 #undef __FUNCT__
860 #define __FUNCT__ "PetscSetMKL_PARDISOFromOptions"
861 PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A)
862 {
863   Mat_MKL_PARDISO     *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr;
864   PetscErrorCode      ierr;
865   PetscInt            icntl,threads=1;
866   PetscBool           flg;
867 
868   PetscFunctionBegin;
869   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");CHKERRQ(ierr);
870 
871   ierr = PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use within PARDISO","None",threads,&threads,&flg);CHKERRQ(ierr);
872   if (flg) PetscSetMKL_PARDISOThreads((int)threads);
873 
874   ierr = PetscOptionsInt("-mat_mkl_pardiso_66","Maximum number of factors with identical sparsity structure that must be kept in memory at the same time","None",mat_mkl_pardiso->maxfct,&icntl,&flg);CHKERRQ(ierr);
875   if (flg) mat_mkl_pardiso->maxfct = icntl;
876 
877   ierr = PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);CHKERRQ(ierr);
878   if (flg) mat_mkl_pardiso->mnum = icntl;
879 
880   ierr = PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);CHKERRQ(ierr);
881   if (flg) mat_mkl_pardiso->msglvl = icntl;
882 
883   ierr = PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);CHKERRQ(ierr);
884   if(flg){
885     void *pt[IPARM_SIZE];
886     mat_mkl_pardiso->mtype = icntl;
887     MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
888 #if defined(PETSC_USE_REAL_SINGLE)
889     mat_mkl_pardiso->iparm[27] = 1;
890 #else
891     mat_mkl_pardiso->iparm[27] = 0;
892 #endif
893     mat_mkl_pardiso->iparm[34] = 1; /* use 0-based indexing */
894   }
895   ierr = PetscOptionsInt("-mat_mkl_pardiso_1","Use default values","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);CHKERRQ(ierr);
896 
897   if(flg && icntl != 0){
898     ierr = PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);CHKERRQ(ierr);
899     if (flg) mat_mkl_pardiso->iparm[1] = icntl;
900 
901     ierr = PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);CHKERRQ(ierr);
902     if (flg) mat_mkl_pardiso->iparm[3] = icntl;
903 
904     ierr = PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);CHKERRQ(ierr);
905     if (flg) mat_mkl_pardiso->iparm[4] = icntl;
906 
907     ierr = PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);CHKERRQ(ierr);
908     if (flg) mat_mkl_pardiso->iparm[5] = icntl;
909 
910     ierr = PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);CHKERRQ(ierr);
911     if (flg) mat_mkl_pardiso->iparm[7] = icntl;
912 
913     ierr = PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);CHKERRQ(ierr);
914     if (flg) mat_mkl_pardiso->iparm[9] = icntl;
915 
916     ierr = PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);CHKERRQ(ierr);
917     if (flg) mat_mkl_pardiso->iparm[10] = icntl;
918 
919     ierr = PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);CHKERRQ(ierr);
920     if (flg) mat_mkl_pardiso->iparm[11] = icntl;
921 
922     ierr = PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);CHKERRQ(ierr);
923     if (flg) mat_mkl_pardiso->iparm[12] = icntl;
924 
925     ierr = PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);CHKERRQ(ierr);
926     if (flg) mat_mkl_pardiso->iparm[17] = icntl;
927 
928     ierr = PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);CHKERRQ(ierr);
929     if (flg) mat_mkl_pardiso->iparm[18] = icntl;
930 
931     ierr = PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);CHKERRQ(ierr);
932     if (flg) mat_mkl_pardiso->iparm[20] = icntl;
933 
934     ierr = PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);CHKERRQ(ierr);
935     if (flg) mat_mkl_pardiso->iparm[23] = icntl;
936 
937     ierr = PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);CHKERRQ(ierr);
938     if (flg) mat_mkl_pardiso->iparm[24] = icntl;
939 
940     ierr = PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);CHKERRQ(ierr);
941     if (flg) mat_mkl_pardiso->iparm[26] = icntl;
942 
943     ierr = PetscOptionsInt("-mat_mkl_pardiso_31","Partial solve and computing selected components of the solution vectors","None",mat_mkl_pardiso->iparm[30],&icntl,&flg);CHKERRQ(ierr);
944     if (flg) mat_mkl_pardiso->iparm[30] = icntl;
945 
946     ierr = PetscOptionsInt("-mat_mkl_pardiso_34","Optimal number of threads for conditional numerical reproducibility (CNR) mode","None",mat_mkl_pardiso->iparm[33],&icntl,&flg);CHKERRQ(ierr);
947     if (flg) mat_mkl_pardiso->iparm[33] = icntl;
948 
949     ierr = PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);CHKERRQ(ierr);
950     if (flg) mat_mkl_pardiso->iparm[59] = icntl;
951   }
952   PetscOptionsEnd();
953   PetscFunctionReturn(0);
954 }
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "MatFactorMKL_PARDISOInitialize_Private"
958 PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso)
959 {
960   PetscErrorCode ierr;
961   PetscInt       i;
962 
963   PetscFunctionBegin;
964   for ( i = 0; i < IPARM_SIZE; i++ ){
965     mat_mkl_pardiso->iparm[i] = 0;
966   }
967   for ( i = 0; i < IPARM_SIZE; i++ ){
968     mat_mkl_pardiso->pt[i] = 0;
969   }
970   /* Default options for both sym and unsym */
971   mat_mkl_pardiso->iparm[ 0] =  1; /* Solver default parameters overriden with provided by iparm */
972   mat_mkl_pardiso->iparm[ 1] =  2; /* Metis reordering */
973   mat_mkl_pardiso->iparm[ 5] =  0; /* Write solution into x */
974   mat_mkl_pardiso->iparm[ 7] =  0; /* Max number of iterative refinement steps */
975   mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
976   mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
977 #if 0
978   mat_mkl_pardiso->iparm[23] =  1; /* Parallel factorization control*/
979 #endif
980   mat_mkl_pardiso->iparm[34] =  1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */
981   mat_mkl_pardiso->iparm[39] =  0; /* Input: matrix/rhs/solution stored on master */
982 
983   mat_mkl_pardiso->CleanUp   = PETSC_FALSE;
984   mat_mkl_pardiso->maxfct    = 1; /* Maximum number of numerical factorizations. */
985   mat_mkl_pardiso->mnum      = 1; /* Which factorization to use. */
986   mat_mkl_pardiso->msglvl    = 0; /* 0: do not print 1: Print statistical information in file */
987   mat_mkl_pardiso->phase     = -1;
988   mat_mkl_pardiso->err       = 0;
989 
990   mat_mkl_pardiso->n         = A->rmap->N;
991   mat_mkl_pardiso->nrhs      = 1;
992   mat_mkl_pardiso->err       = 0;
993   mat_mkl_pardiso->phase     = -1;
994 
995   if(ftype == MAT_FACTOR_LU){
996     mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */
997     mat_mkl_pardiso->iparm[10] =  1; /* Use nonsymmetric permutation and scaling MPS */
998     mat_mkl_pardiso->iparm[12] =  1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
999 
1000   } else {
1001     mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */
1002     mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */
1003     mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
1004 /*    mat_mkl_pardiso->iparm[20] =  1; */ /* Apply 1x1 and 2x2 Bunch-Kaufman pivoting during the factorization process */
1005 #if defined(PETSC_USE_DEBUG)
1006     mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */
1007 #endif
1008   }
1009   ierr = PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);CHKERRQ(ierr);
1010   for(i = 0; i < A->rmap->N; i++){
1011     mat_mkl_pardiso->perm[i] = 0;
1012   }
1013   mat_mkl_pardiso->schur_size = 0;
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 #undef __FUNCT__
1018 #define __FUNCT__ "MatFactorSymbolic_AIJMKL_PARDISO_Private"
1019 PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info)
1020 {
1021   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr;
1022   PetscErrorCode  ierr;
1023 
1024   PetscFunctionBegin;
1025   mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
1026   ierr = PetscSetMKL_PARDISOFromOptions(F,A);CHKERRQ(ierr);
1027 
1028   /* throw away any previously computed structure */
1029   if (mat_mkl_pardiso->freeaij) {
1030     ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr);
1031     ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr);
1032   }
1033   ierr = (*mat_mkl_pardiso->Convert)(A,mat_mkl_pardiso->needsym,MAT_INITIAL_MATRIX,&mat_mkl_pardiso->freeaij,&mat_mkl_pardiso->nz,&mat_mkl_pardiso->ia,&mat_mkl_pardiso->ja,&mat_mkl_pardiso->a);CHKERRQ(ierr);
1034   mat_mkl_pardiso->n = A->rmap->N;
1035 
1036   mat_mkl_pardiso->phase = JOB_ANALYSIS;
1037 
1038   MKL_PARDISO (mat_mkl_pardiso->pt,
1039     &mat_mkl_pardiso->maxfct,
1040     &mat_mkl_pardiso->mnum,
1041     &mat_mkl_pardiso->mtype,
1042     &mat_mkl_pardiso->phase,
1043     &mat_mkl_pardiso->n,
1044     mat_mkl_pardiso->a,
1045     mat_mkl_pardiso->ia,
1046     mat_mkl_pardiso->ja,
1047     mat_mkl_pardiso->perm,
1048     &mat_mkl_pardiso->nrhs,
1049     mat_mkl_pardiso->iparm,
1050     &mat_mkl_pardiso->msglvl,
1051     NULL,
1052     NULL,
1053     &mat_mkl_pardiso->err);
1054   if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d\n. Please check manual",mat_mkl_pardiso->err);
1055 
1056   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
1057 
1058   if (F->factortype == MAT_FACTOR_LU) {
1059     F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO;
1060   } else {
1061     F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO;
1062   }
1063   F->ops->solve           = MatSolve_MKL_PARDISO;
1064   F->ops->solvetranspose  = MatSolveTranspose_MKL_PARDISO;
1065   F->ops->matsolve        = MatMatSolve_MKL_PARDISO;
1066   PetscFunctionReturn(0);
1067 }
1068 
1069 #undef __FUNCT__
1070 #define __FUNCT__ "MatLUFactorSymbolic_AIJMKL_PARDISO"
1071 PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
1072 {
1073   PetscErrorCode ierr;
1074 
1075   PetscFunctionBegin;
1076   ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr);
1077   PetscFunctionReturn(0);
1078 }
1079 
1080 #undef __FUNCT__
1081 #define __FUNCT__ "MatCholeskyFactorSymbolic_AIJMKL_PARDISO"
1082 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info)
1083 {
1084   PetscErrorCode ierr;
1085 
1086   PetscFunctionBegin;
1087   ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr);
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 #undef __FUNCT__
1092 #define __FUNCT__ "MatView_MKL_PARDISO"
1093 PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer)
1094 {
1095   PetscErrorCode    ierr;
1096   PetscBool         iascii;
1097   PetscViewerFormat format;
1098   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr;
1099   PetscInt          i;
1100 
1101   PetscFunctionBegin;
1102   if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(0);
1103 
1104   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1105   if (iascii) {
1106     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1107     if (format == PETSC_VIEWER_ASCII_INFO) {
1108       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");CHKERRQ(ierr);
1109       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase:             %d \n",mat_mkl_pardiso->phase);CHKERRQ(ierr);
1110       for(i = 1; i <= 64; i++){
1111         ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]:     %d \n",i, mat_mkl_pardiso->iparm[i - 1]);CHKERRQ(ierr);
1112       }
1113       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct:     %d \n", mat_mkl_pardiso->maxfct);CHKERRQ(ierr);
1114       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum:     %d \n", mat_mkl_pardiso->mnum);CHKERRQ(ierr);
1115       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype:     %d \n", mat_mkl_pardiso->mtype);CHKERRQ(ierr);
1116       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n:     %d \n", mat_mkl_pardiso->n);CHKERRQ(ierr);
1117       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs:     %d \n", mat_mkl_pardiso->nrhs);CHKERRQ(ierr);
1118       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl:     %d \n", mat_mkl_pardiso->msglvl);CHKERRQ(ierr);
1119     }
1120   }
1121   PetscFunctionReturn(0);
1122 }
1123 
1124 
1125 #undef __FUNCT__
1126 #define __FUNCT__ "MatGetInfo_MKL_PARDISO"
1127 PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info)
1128 {
1129   Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)A->spptr;
1130 
1131   PetscFunctionBegin;
1132   info->block_size        = 1.0;
1133   info->nz_used           = mat_mkl_pardiso->nz;
1134   info->nz_allocated      = mat_mkl_pardiso->nz;
1135   info->nz_unneeded       = 0.0;
1136   info->assemblies        = 0.0;
1137   info->mallocs           = 0.0;
1138   info->memory            = 0.0;
1139   info->fill_ratio_given  = 0;
1140   info->fill_ratio_needed = 0;
1141   info->factor_mallocs    = 0;
1142   PetscFunctionReturn(0);
1143 }
1144 
1145 #undef __FUNCT__
1146 #define __FUNCT__ "MatMkl_PardisoSetCntl_MKL_PARDISO"
1147 PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival)
1148 {
1149   Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)F->spptr;
1150 
1151   PetscFunctionBegin;
1152   if(icntl <= 64){
1153     mat_mkl_pardiso->iparm[icntl - 1] = ival;
1154   } else {
1155     if(icntl == 65)
1156       PetscSetMKL_PARDISOThreads(ival);
1157     else if(icntl == 66)
1158       mat_mkl_pardiso->maxfct = ival;
1159     else if(icntl == 67)
1160       mat_mkl_pardiso->mnum = ival;
1161     else if(icntl == 68)
1162       mat_mkl_pardiso->msglvl = ival;
1163     else if(icntl == 69){
1164       void *pt[IPARM_SIZE];
1165       mat_mkl_pardiso->mtype = ival;
1166       MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
1167 #if defined(PETSC_USE_REAL_SINGLE)
1168       mat_mkl_pardiso->iparm[27] = 1;
1169 #else
1170       mat_mkl_pardiso->iparm[27] = 0;
1171 #endif
1172       mat_mkl_pardiso->iparm[34] = 1;
1173     } else if(icntl==70) {
1174       mat_mkl_pardiso->solve_interior = !!ival;
1175     }
1176   }
1177   PetscFunctionReturn(0);
1178 }
1179 
1180 #undef __FUNCT__
1181 #define __FUNCT__ "MatMkl_PardisoSetCntl"
1182 /*@
1183   MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters
1184 
1185    Logically Collective on Mat
1186 
1187    Input Parameters:
1188 +  F - the factored matrix obtained by calling MatGetFactor()
1189 .  icntl - index of Mkl_Pardiso parameter
1190 -  ival - value of Mkl_Pardiso parameter
1191 
1192   Options Database:
1193 .   -mat_mkl_pardiso_<icntl> <ival>
1194 
1195    Level: beginner
1196 
1197    References: Mkl_Pardiso Users' Guide
1198 
1199 .seealso: MatGetFactor()
1200 @*/
1201 PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival)
1202 {
1203   PetscErrorCode ierr;
1204 
1205   PetscFunctionBegin;
1206   ierr = PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
1207   PetscFunctionReturn(0);
1208 }
1209 
1210 /*MC
1211   MATSOLVERMKL_PARDISO -  A matrix type providing direct solvers (LU) for
1212   sequential matrices via the external package MKL_PARDISO.
1213 
1214   Works with MATSEQAIJ matrices
1215 
1216   Use -pc_type lu -pc_factor_mat_solver_package mkl_pardiso to us this direct solver
1217 
1218   Options Database Keys:
1219 + -mat_mkl_pardiso_65 - Number of threads to use within MKL_PARDISO
1220 . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
1221 . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase
1222 . -mat_mkl_pardiso_68 - Message level information
1223 . -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
1224 . -mat_mkl_pardiso_1  - Use default values
1225 . -mat_mkl_pardiso_2  - Fill-in reducing ordering for the input matrix
1226 . -mat_mkl_pardiso_4  - Preconditioned CGS/CG
1227 . -mat_mkl_pardiso_5  - User permutation
1228 . -mat_mkl_pardiso_6  - Write solution on x
1229 . -mat_mkl_pardiso_8  - Iterative refinement step
1230 . -mat_mkl_pardiso_10 - Pivoting perturbation
1231 . -mat_mkl_pardiso_11 - Scaling vectors
1232 . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A
1233 . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching
1234 . -mat_mkl_pardiso_18 - Numbers of non-zero elements
1235 . -mat_mkl_pardiso_19 - Report number of floating point operations
1236 . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices
1237 . -mat_mkl_pardiso_24 - Parallel factorization control
1238 . -mat_mkl_pardiso_25 - Parallel forward/backward solve control
1239 . -mat_mkl_pardiso_27 - Matrix checker
1240 . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors
1241 . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
1242 - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode
1243 
1244   Level: beginner
1245 
1246   For more information please check  mkl_pardiso manual
1247 
1248 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
1249 
1250 M*/
1251 #undef __FUNCT__
1252 #define __FUNCT__ "MatFactorGetSolverPackage_mkl_pardiso"
1253 static PetscErrorCode MatFactorGetSolverPackage_mkl_pardiso(Mat A, const MatSolverPackage *type)
1254 {
1255   PetscFunctionBegin;
1256   *type = MATSOLVERMKL_PARDISO;
1257   PetscFunctionReturn(0);
1258 }
1259 
1260 #undef __FUNCT__
1261 #define __FUNCT__ "MatGetFactor_aij_mkl_pardiso"
1262 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F)
1263 {
1264   Mat             B;
1265   PetscErrorCode  ierr;
1266   Mat_MKL_PARDISO *mat_mkl_pardiso;
1267   PetscBool       isSeqAIJ,isSeqBAIJ,isSeqSBAIJ;
1268 
1269   PetscFunctionBegin;
1270   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
1271   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr);
1272   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
1273   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
1274   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
1275   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
1276   ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
1277   ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1278   ierr = MatSeqSBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr);
1279 
1280   ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr);
1281   B->spptr = mat_mkl_pardiso;
1282 
1283   ierr = MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr);
1284   if (ftype == MAT_FACTOR_LU) {
1285     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO;
1286     B->factortype            = MAT_FACTOR_LU;
1287     mat_mkl_pardiso->needsym = PETSC_FALSE;
1288     if (isSeqAIJ) {
1289       mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1290     } else if (isSeqBAIJ) {
1291       mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1292     } else if (isSeqSBAIJ) {
1293       SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead");
1294     } else {
1295       SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU with %s format",((PetscObject)A)->type_name);
1296     }
1297 #if defined(PETSC_USE_COMPLEX)
1298     mat_mkl_pardiso->schur_solver_type = 0; /* use a general solver for the moment */
1299     mat_mkl_pardiso->mtype = 13;
1300 #else
1301     if (A->structurally_symmetric) {
1302       mat_mkl_pardiso->mtype = 1;
1303     } else {
1304       mat_mkl_pardiso->mtype = 11;
1305     }
1306 #endif
1307     mat_mkl_pardiso->schur_solver_type = 0;
1308   } else {
1309     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO;
1310     B->factortype                  = MAT_FACTOR_CHOLESKY;
1311     if (isSeqAIJ) {
1312       mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1313     } else if (isSeqBAIJ) {
1314       mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1315     } else if (isSeqSBAIJ) {
1316       mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij;
1317     } else {
1318       SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with %s format",((PetscObject)A)->type_name);
1319     }
1320     mat_mkl_pardiso->needsym = PETSC_TRUE;
1321 #if defined(PETSC_USE_COMPLEX)
1322     mat_mkl_pardiso->schur_solver_type = 0; /* use a general solver for the moment */
1323     mat_mkl_pardiso->mtype = 13;
1324 #else
1325     if (A->spd_set && A->spd) {
1326       mat_mkl_pardiso->schur_solver_type = 1;
1327       mat_mkl_pardiso->mtype = 2;
1328     } else {
1329       mat_mkl_pardiso->schur_solver_type = 2;
1330       mat_mkl_pardiso->mtype = -2;
1331     }
1332 #endif
1333   }
1334   mat_mkl_pardiso->Destroy = B->ops->destroy;
1335   B->ops->destroy          = MatDestroy_MKL_PARDISO;
1336   B->ops->view             = MatView_MKL_PARDISO;
1337   B->factortype            = ftype;
1338   B->ops->getinfo          = MatGetInfo_MKL_PARDISO;
1339   B->assembled             = PETSC_TRUE;
1340 
1341   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr);
1342   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MKL_PARDISO);CHKERRQ(ierr);
1343   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MKL_PARDISO);CHKERRQ(ierr);
1344   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MKL_PARDISO);CHKERRQ(ierr);
1345   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MKL_PARDISO);CHKERRQ(ierr);
1346   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MKL_PARDISO);CHKERRQ(ierr);
1347   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MKL_PARDISO);CHKERRQ(ierr);
1348   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr);
1349 
1350   *F = B;
1351   PetscFunctionReturn(0);
1352 }
1353 
1354 #undef __FUNCT__
1355 #define __FUNCT__ "MatSolverPackageRegister_MKL_Pardiso"
1356 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MKL_Pardiso(void)
1357 {
1358   PetscErrorCode ierr;
1359 
1360   PetscFunctionBegin;
1361   ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1362   ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1363   ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr);
1364   PetscFunctionReturn(0);
1365 }
1366 
1367