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