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