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