xref: /petsc/src/mat/impls/aij/seq/mkl_pardiso/mkl_pardiso.c (revision 39bb8a47cdcfd31b9e015b972652b0f43fe71e8b)
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/dense/seq/dense.h>
7 
8 #include <stdio.h>
9 #include <stdlib.h>
10 #include <math.h>
11 #include <mkl.h>
12 
13 /*
14  *  Possible mkl_pardiso phases that controls the execution of the solver.
15  *  For more information check mkl_pardiso manual.
16  */
17 #define JOB_ANALYSIS 11
18 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION 12
19 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13
20 #define JOB_NUMERICAL_FACTORIZATION 22
21 #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 23
22 #define JOB_SOLVE_ITERATIVE_REFINEMENT 33
23 #define JOB_SOLVE_FORWARD_SUBSTITUTION 331
24 #define JOB_SOLVE_DIAGONAL_SUBSTITUTION 332
25 #define JOB_SOLVE_BACKWARD_SUBSTITUTION 333
26 #define JOB_RELEASE_OF_LU_MEMORY 0
27 #define JOB_RELEASE_OF_ALL_MEMORY -1
28 
29 #define IPARM_SIZE 64
30 
31 #if defined(PETSC_USE_64BIT_INDICES)
32  #if defined(PETSC_HAVE_LIBMKL_INTEL_ILP64)
33   /* sizeof(MKL_INT) == sizeof(long long int) if ilp64*/
34   #define INT_TYPE long long int
35   #define MKL_PARDISO pardiso
36   #define MKL_PARDISO_INIT pardisoinit
37  #else
38   #define INT_TYPE long long int
39   #define MKL_PARDISO pardiso_64
40   #define MKL_PARDISO_INIT pardiso_64init
41  #endif
42 #else
43  #define INT_TYPE int
44  #define MKL_PARDISO pardiso
45  #define MKL_PARDISO_INIT pardisoinit
46 #endif
47 
48 
49 /*
50  *  Internal data structure.
51  *  For more information check mkl_pardiso manual.
52  */
53 typedef struct {
54 
55   /* Configuration vector*/
56   INT_TYPE     iparm[IPARM_SIZE];
57 
58   /*
59    * Internal mkl_pardiso memory location.
60    * After the first call to mkl_pardiso do not modify pt, as that could cause a serious memory leak.
61    */
62   void         *pt[IPARM_SIZE];
63 
64   /* Basic mkl_pardiso info*/
65   INT_TYPE     phase, maxfct, mnum, mtype, n, nrhs, msglvl, err;
66 
67   /* Matrix structure*/
68   void         *a;
69   INT_TYPE     *ia, *ja;
70 
71   /* Number of non-zero elements*/
72   INT_TYPE     nz;
73 
74   /* Row permutaton vector*/
75   INT_TYPE     *perm;
76 
77   /* Define if matrix preserves sparse structure.*/
78   MatStructure matstruc;
79 
80   /* True if mkl_pardiso function have been used.*/
81   PetscBool CleanUp;
82 
83 } Mat_MKL_PARDISO;
84 
85 
86 void pardiso_64init(void *pt, INT_TYPE *mtype, INT_TYPE iparm []){
87   int     iparm_copy[IPARM_SIZE], mtype_copy, i;
88   mtype_copy = *mtype;
89   pardisoinit(pt, &mtype_copy, iparm_copy);
90   for(i = 0; i < IPARM_SIZE; i++)
91     iparm[i] = iparm_copy[i];
92 }
93 
94 
95 /*
96  * Copy the elements of matrix A.
97  * Input:
98  *   - Mat A: MATSEQAIJ matrix
99  *   - int shift: matrix index.
100  *     - 0 for c representation
101  *     - 1 for fortran representation
102  *   - MatReuse reuse:
103  *     - MAT_INITIAL_MATRIX: Create a new aij representation
104  *     - MAT_REUSE_MATRIX: Reuse all aij representation and just change values
105  * Output:
106  *   - int *nnz: Number of nonzero-elements.
107  *   - int **r pointer to i index
108  *   - int **c pointer to j elements
109  *   - MATRIXTYPE **v: Non-zero elements
110  */
111 #undef __FUNCT__
112 #define __FUNCT__ "MatCopy_MKL_PARDISO"
113 PetscErrorCode MatCopy_MKL_PARDISO(Mat A, MatReuse reuse, INT_TYPE *nnz, INT_TYPE **r, INT_TYPE **c, void **v){
114 
115   Mat_SeqAIJ     *aa=(Mat_SeqAIJ*)A->data;
116 
117   PetscFunctionBegin;
118   *v=aa->a;
119   if (reuse == MAT_INITIAL_MATRIX) {
120     *r   = (INT_TYPE*)aa->i;
121     *c   = (INT_TYPE*)aa->j;
122     *nnz = aa->nz;
123   }
124   PetscFunctionReturn(0);
125 }
126 
127 
128 /*
129  * Free memory for Mat_MKL_PARDISO structure and pointers to objects.
130  */
131 #undef __FUNCT__
132 #define __FUNCT__ "MatDestroy_MKL_PARDISO"
133 PetscErrorCode MatDestroy_MKL_PARDISO(Mat A){
134   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr;
135   PetscErrorCode ierr;
136 
137   PetscFunctionBegin;
138   /* Terminate instance, deallocate memories */
139   if (mat_mkl_pardiso->CleanUp) {
140     mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;
141 
142 
143     MKL_PARDISO (mat_mkl_pardiso->pt,
144       &mat_mkl_pardiso->maxfct,
145       &mat_mkl_pardiso->mnum,
146       &mat_mkl_pardiso->mtype,
147       &mat_mkl_pardiso->phase,
148       &mat_mkl_pardiso->n,
149       NULL,
150       NULL,
151       NULL,
152       mat_mkl_pardiso->perm,
153       &mat_mkl_pardiso->nrhs,
154       mat_mkl_pardiso->iparm,
155       &mat_mkl_pardiso->msglvl,
156       NULL,
157       NULL,
158       &mat_mkl_pardiso->err);
159   }
160   ierr = PetscFree(mat_mkl_pardiso->perm);CHKERRQ(ierr);
161   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
162 
163   /* clear composed functions */
164   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
165   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);CHKERRQ(ierr);
166   PetscFunctionReturn(0);
167 }
168 
169 
170 /*
171  * Computes Ax = b
172  */
173 #undef __FUNCT__
174 #define __FUNCT__ "MatSolve_MKL_PARDISO"
175 PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x)
176 {
177   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr;
178   PetscErrorCode    ierr;
179   PetscScalar       *xarray;
180   const PetscScalar *barray;
181 
182   PetscFunctionBegin;
183 
184 
185   mat_mkl_pardiso->nrhs = 1;
186   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
187   ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
188 
189   /* solve phase */
190   /*-------------*/
191   mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
192   MKL_PARDISO (mat_mkl_pardiso->pt,
193     &mat_mkl_pardiso->maxfct,
194     &mat_mkl_pardiso->mnum,
195     &mat_mkl_pardiso->mtype,
196     &mat_mkl_pardiso->phase,
197     &mat_mkl_pardiso->n,
198     mat_mkl_pardiso->a,
199     mat_mkl_pardiso->ia,
200     mat_mkl_pardiso->ja,
201     mat_mkl_pardiso->perm,
202     &mat_mkl_pardiso->nrhs,
203     mat_mkl_pardiso->iparm,
204     &mat_mkl_pardiso->msglvl,
205     (void*)barray,
206     (void*)xarray,
207     &mat_mkl_pardiso->err);
208 
209 
210   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);
211   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
212   ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
213   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
214   PetscFunctionReturn(0);
215 }
216 
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "MatSolveTranspose_MKL_PARDISO"
220 PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x)
221 {
222   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr;
223   PetscErrorCode ierr;
224 
225   PetscFunctionBegin;
226 #if defined(PETSC_USE_COMPLEX)
227   mat_mkl_pardiso->iparm[12 - 1] = 1;
228 #else
229   mat_mkl_pardiso->iparm[12 - 1] = 2;
230 #endif
231   ierr = MatSolve_MKL_PARDISO(A,b,x);CHKERRQ(ierr);
232   mat_mkl_pardiso->iparm[12 - 1] = 0;
233   PetscFunctionReturn(0);
234 }
235 
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "MatMatSolve_MKL_PARDISO"
239 PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X)
240 {
241   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr;
242   PetscErrorCode    ierr;
243   PetscScalar       *barray, *xarray;
244   PetscBool         flg;
245 
246   PetscFunctionBegin;
247 
248   ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr);
249   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix");
250   ierr = PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr);
251   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix");
252 
253   ierr = MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);CHKERRQ(ierr);
254 
255   if(mat_mkl_pardiso->nrhs > 0){
256     ierr = MatDenseGetArray(B,&barray);
257     ierr = MatDenseGetArray(X,&xarray);
258 
259     /* solve phase */
260     /*-------------*/
261     mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
262     MKL_PARDISO (mat_mkl_pardiso->pt,
263       &mat_mkl_pardiso->maxfct,
264       &mat_mkl_pardiso->mnum,
265       &mat_mkl_pardiso->mtype,
266       &mat_mkl_pardiso->phase,
267       &mat_mkl_pardiso->n,
268       mat_mkl_pardiso->a,
269       mat_mkl_pardiso->ia,
270       mat_mkl_pardiso->ja,
271       mat_mkl_pardiso->perm,
272       &mat_mkl_pardiso->nrhs,
273       mat_mkl_pardiso->iparm,
274       &mat_mkl_pardiso->msglvl,
275       (void*)barray,
276       (void*)xarray,
277       &mat_mkl_pardiso->err);
278     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);
279   }
280   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
281   PetscFunctionReturn(0);
282 
283 }
284 
285 /*
286  * LU Decomposition
287  */
288 #undef __FUNCT__
289 #define __FUNCT__ "MatFactorNumeric_MKL_PARDISO"
290 PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info){
291   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->spptr;
292   PetscErrorCode ierr;
293 
294   /* numerical factorization phase */
295   /*-------------------------------*/
296 
297   PetscFunctionBegin;
298 
299   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
300   ierr = MatCopy_MKL_PARDISO(A, MAT_REUSE_MATRIX, &mat_mkl_pardiso->nz, &mat_mkl_pardiso->ia, &mat_mkl_pardiso->ja, &mat_mkl_pardiso->a);CHKERRQ(ierr);
301 
302   /* numerical factorization phase */
303   /*-------------------------------*/
304   mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION;
305   MKL_PARDISO (mat_mkl_pardiso->pt,
306     &mat_mkl_pardiso->maxfct,
307     &mat_mkl_pardiso->mnum,
308     &mat_mkl_pardiso->mtype,
309     &mat_mkl_pardiso->phase,
310     &mat_mkl_pardiso->n,
311     mat_mkl_pardiso->a,
312     mat_mkl_pardiso->ia,
313     mat_mkl_pardiso->ja,
314     mat_mkl_pardiso->perm,
315     &mat_mkl_pardiso->nrhs,
316     mat_mkl_pardiso->iparm,
317     &mat_mkl_pardiso->msglvl,
318     NULL,
319     NULL,
320     &mat_mkl_pardiso->err);
321   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);
322 
323   mat_mkl_pardiso->matstruc     = SAME_NONZERO_PATTERN;
324   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
325   PetscFunctionReturn(0);
326 }
327 
328 /* Sets mkl_pardiso options from the options database */
329 #undef __FUNCT__
330 #define __FUNCT__ "PetscSetMKL_PARDISOFromOptions"
331 PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A){
332   Mat_MKL_PARDISO     *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr;
333   PetscErrorCode      ierr;
334   PetscInt            icntl;
335   PetscBool           flg;
336   int                 pt[IPARM_SIZE], threads = 1;
337 
338   PetscFunctionBegin;
339   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");CHKERRQ(ierr);
340   ierr = PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use","None",threads,&threads,&flg);CHKERRQ(ierr);
341   if (flg) mkl_set_num_threads(threads);
342 
343   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);
344   if (flg) mat_mkl_pardiso->maxfct = icntl;
345 
346   ierr = PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);CHKERRQ(ierr);
347   if (flg) mat_mkl_pardiso->mnum = icntl;
348 
349   ierr = PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);CHKERRQ(ierr);
350   if (flg) mat_mkl_pardiso->msglvl = icntl;
351 
352   ierr = PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);CHKERRQ(ierr);
353   if(flg){
354    mat_mkl_pardiso->mtype = icntl;
355    MKL_PARDISO_INIT(&pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
356 #if defined(PETSC_USE_REAL_SINGLE)
357     mat_mkl_pardiso->iparm[27] = 1;
358 #else
359     mat_mkl_pardiso->iparm[27] = 0;
360 #endif
361     mat_mkl_pardiso->iparm[34] = 1;
362   }
363   ierr = PetscOptionsInt("-mat_mkl_pardiso_1","Use default values","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);CHKERRQ(ierr);
364 
365   if(flg && icntl != 0){
366     ierr = PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);CHKERRQ(ierr);
367     if (flg) mat_mkl_pardiso->iparm[1] = icntl;
368 
369     ierr = PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);CHKERRQ(ierr);
370     if (flg) mat_mkl_pardiso->iparm[3] = icntl;
371 
372     ierr = PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);CHKERRQ(ierr);
373     if (flg) mat_mkl_pardiso->iparm[4] = icntl;
374 
375     ierr = PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);CHKERRQ(ierr);
376     if (flg) mat_mkl_pardiso->iparm[5] = icntl;
377 
378     ierr = PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);CHKERRQ(ierr);
379     if (flg) mat_mkl_pardiso->iparm[7] = icntl;
380 
381     ierr = PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);CHKERRQ(ierr);
382     if (flg) mat_mkl_pardiso->iparm[9] = icntl;
383 
384     ierr = PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);CHKERRQ(ierr);
385     if (flg) mat_mkl_pardiso->iparm[10] = icntl;
386 
387     ierr = PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);CHKERRQ(ierr);
388     if (flg) mat_mkl_pardiso->iparm[11] = icntl;
389 
390     ierr = PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);CHKERRQ(ierr);
391     if (flg) mat_mkl_pardiso->iparm[12] = icntl;
392 
393     ierr = PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);CHKERRQ(ierr);
394     if (flg) mat_mkl_pardiso->iparm[17] = icntl;
395 
396     ierr = PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);CHKERRQ(ierr);
397     if (flg) mat_mkl_pardiso->iparm[18] = icntl;
398 
399     ierr = PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);CHKERRQ(ierr);
400     if (flg) mat_mkl_pardiso->iparm[20] = icntl;
401 
402     ierr = PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);CHKERRQ(ierr);
403     if (flg) mat_mkl_pardiso->iparm[23] = icntl;
404 
405     ierr = PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);CHKERRQ(ierr);
406     if (flg) mat_mkl_pardiso->iparm[24] = icntl;
407 
408     ierr = PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);CHKERRQ(ierr);
409     if (flg) mat_mkl_pardiso->iparm[26] = icntl;
410 
411     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);
412     if (flg) mat_mkl_pardiso->iparm[30] = icntl;
413 
414     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);
415     if (flg) mat_mkl_pardiso->iparm[33] = icntl;
416 
417     ierr = PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);CHKERRQ(ierr);
418     if (flg) mat_mkl_pardiso->iparm[59] = icntl;
419   }
420   PetscOptionsEnd();
421   PetscFunctionReturn(0);
422 }
423 
424 
425 #undef __FUNCT__
426 #define __FUNCT__ "PetscInitializeMKL_PARDISO"
427 PetscErrorCode PetscInitializeMKL_PARDISO(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso){
428   PetscErrorCode ierr;
429   PetscInt       i;
430 
431   PetscFunctionBegin;
432 
433   for ( i = 0; i < IPARM_SIZE; i++ )
434   {
435     mat_mkl_pardiso->iparm[i] = 0;
436   }
437 
438   for ( i = 0; i < IPARM_SIZE; i++ )
439   {
440     mat_mkl_pardiso->pt[i] = 0;
441   }
442 
443   /*Default options for both sym and unsym */
444   mat_mkl_pardiso->iparm[ 0] =  1; /* Solver default parameters overriden with provided by iparm */
445   mat_mkl_pardiso->iparm[ 1] =  2; /* Metis reordering */
446   mat_mkl_pardiso->iparm[ 5] =  0; /* Write solution into x */
447   mat_mkl_pardiso->iparm[ 7] =  2; /* Max number of iterative refinement steps */
448   mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
449   mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
450 #if 0
451  mat_mkl_pardiso->iparm[23] =  1; /* Parallel factorization control*/
452 #endif
453   mat_mkl_pardiso->iparm[34] =  1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */
454   mat_mkl_pardiso->iparm[39] =  0; /* Input: matrix/rhs/solution stored on master
455                                       */
456 
457   mat_mkl_pardiso->CleanUp = PETSC_FALSE;
458   mat_mkl_pardiso->maxfct  = 1; /* Maximum number of numerical factorizations. */
459   mat_mkl_pardiso->mnum    = 1; /* Which factorization to use. */
460   mat_mkl_pardiso->msglvl  = 0; /* 0: do not print 1: Print statistical information in file */
461   mat_mkl_pardiso->phase = -1;
462   mat_mkl_pardiso->err = 0;
463 
464   mat_mkl_pardiso->n = A->rmap->N;
465   mat_mkl_pardiso->nrhs = 1;
466   mat_mkl_pardiso->err = 0;
467   mat_mkl_pardiso->phase = -1;
468   if(ftype == MAT_FACTOR_LU)
469   {
470     /*Default type for non-sym*/
471 #if defined(PETSC_USE_COMPLEX)
472     mat_mkl_pardiso->mtype = 13;
473 #else
474     mat_mkl_pardiso->mtype = 11;
475 #endif
476 
477     mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */
478     mat_mkl_pardiso->iparm[10] =  1; /* Use nonsymmetric permutation and scaling MPS */
479     mat_mkl_pardiso->iparm[12] =  1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
480 
481   }
482   else
483   {
484     /*Default type for sym*/
485 #if defined(PETSC_USE_COMPLEX)
486     mat_mkl_pardiso ->mtype  = 3;
487 #else
488     mat_mkl_pardiso ->mtype  = -2;
489 #endif
490 
491     mat_mkl_pardiso->iparm[ 9] =  13; /* Perturb the pivot elements with 1E-13 */
492     mat_mkl_pardiso->iparm[10] =  0; /* Use nonsymmetric permutation and scaling MPS */
493     mat_mkl_pardiso->iparm[12] =  1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
494 
495 /*    mat_mkl_pardiso->iparm[20] =  1; */ /* Apply 1x1 and 2x2 Bunch-Kaufman pivoting during the factorization process */
496 
497 #if defined(PETSC_USE_DEBUG)
498     mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */
499 #endif
500   }
501 
502   ierr = PetscMalloc(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);CHKERRQ(ierr);
503   for(i = 0; i < A->rmap->N; i++)
504   {
505     mat_mkl_pardiso->perm[i] = 0;
506   }
507 
508   PetscFunctionReturn(0);
509 }
510 
511 
512 /*
513  * Symbolic decomposition. Mkl_Pardiso analysis phase.
514  */
515 #undef __FUNCT__
516 #define __FUNCT__ "MatFactorSymbolic_AIJMKL_PARDISO_Private"
517 PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info){
518 
519   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr;
520   PetscErrorCode ierr;
521 
522   PetscFunctionBegin;
523 
524   mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
525 
526   /* Set MKL_PARDISO options from the options database */
527   ierr = PetscSetMKL_PARDISOFromOptions(F,A);CHKERRQ(ierr);
528 
529   ierr = MatCopy_MKL_PARDISO(A, MAT_INITIAL_MATRIX, &mat_mkl_pardiso->nz, &mat_mkl_pardiso->ia, &mat_mkl_pardiso->ja, &mat_mkl_pardiso->a);CHKERRQ(ierr);
530   mat_mkl_pardiso->n = A->rmap->N;
531 
532   /* analysis phase */
533   /*----------------*/
534 
535   mat_mkl_pardiso->phase = JOB_ANALYSIS;
536 
537   MKL_PARDISO (mat_mkl_pardiso->pt,
538     &mat_mkl_pardiso->maxfct,
539     &mat_mkl_pardiso->mnum,
540     &mat_mkl_pardiso->mtype,
541     &mat_mkl_pardiso->phase,
542     &mat_mkl_pardiso->n,
543     mat_mkl_pardiso->a,
544     mat_mkl_pardiso->ia,
545     mat_mkl_pardiso->ja,
546     mat_mkl_pardiso->perm,
547     &mat_mkl_pardiso->nrhs,
548     mat_mkl_pardiso->iparm,
549     &mat_mkl_pardiso->msglvl,
550     NULL,
551     NULL,
552     &mat_mkl_pardiso->err);
553 
554   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);
555 
556   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
557 
558   if(F->factortype == MAT_FACTOR_LU)
559   {
560     F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO;
561   }
562   else
563   {
564     F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO;
565   }
566   F->ops->solve           = MatSolve_MKL_PARDISO;
567   F->ops->solvetranspose  = MatSolveTranspose_MKL_PARDISO;
568   F->ops->matsolve        = MatMatSolve_MKL_PARDISO;
569   PetscFunctionReturn(0);
570 }
571 
572 #undef __FUNCT__
573 #define __FUNCT__ "MatLUFactorSymbolic_AIJMKL_PARDISO"
574 PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info){
575 
576   PetscErrorCode ierr;
577 
578   PetscFunctionBegin;
579 
580   ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);
581   CHKERRQ(ierr);
582 
583   PetscFunctionReturn(0);
584 }
585 
586 #undef __FUNCT__
587 #define __FUNCT__ "MatCholeskyFactorSymbolic_AIJMKL_PARDISO"
588 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info){
589 
590   PetscErrorCode ierr;
591 
592   PetscFunctionBegin;
593 
594   ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);
595   CHKERRQ(ierr);
596 
597   PetscFunctionReturn(0);
598 }
599 
600 #undef __FUNCT__
601 #define __FUNCT__ "MatView_MKL_PARDISO"
602 PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer){
603   PetscErrorCode    ierr;
604   PetscBool         iascii;
605   PetscViewerFormat format;
606   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr;
607   PetscInt          i;
608 
609   PetscFunctionBegin;
610   /* check if matrix is mkl_pardiso type */
611   if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(0);
612 
613   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
614   if (iascii) {
615     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
616     if (format == PETSC_VIEWER_ASCII_INFO) {
617       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");CHKERRQ(ierr);
618       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase:             %d \n",mat_mkl_pardiso->phase);CHKERRQ(ierr);
619       for(i = 1; i <= 64; i++){
620         ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]:     %d \n",i, mat_mkl_pardiso->iparm[i - 1]);CHKERRQ(ierr);
621       }
622       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct:     %d \n", mat_mkl_pardiso->maxfct);CHKERRQ(ierr);
623       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum:     %d \n", mat_mkl_pardiso->mnum);CHKERRQ(ierr);
624       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype:     %d \n", mat_mkl_pardiso->mtype);CHKERRQ(ierr);
625       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n:     %d \n", mat_mkl_pardiso->n);CHKERRQ(ierr);
626       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs:     %d \n", mat_mkl_pardiso->nrhs);CHKERRQ(ierr);
627       ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl:     %d \n", mat_mkl_pardiso->msglvl);CHKERRQ(ierr);
628     }
629   }
630   PetscFunctionReturn(0);
631 }
632 
633 
634 #undef __FUNCT__
635 #define __FUNCT__ "MatGetInfo_MKL_PARDISO"
636 PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info){
637   Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)A->spptr;
638 
639   PetscFunctionBegin;
640   info->block_size        = 1.0;
641   info->nz_allocated      = mat_mkl_pardiso->nz + 0.0;
642   info->nz_unneeded       = 0.0;
643   info->assemblies        = 0.0;
644   info->mallocs           = 0.0;
645   info->memory            = 0.0;
646   info->fill_ratio_given  = 0;
647   info->fill_ratio_needed = 0;
648   info->factor_mallocs    = 0;
649   PetscFunctionReturn(0);
650 }
651 
652 #undef __FUNCT__
653 #define __FUNCT__ "MatMkl_PardisoSetCntl_MKL_PARDISO"
654 PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival){
655   Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)F->spptr;
656   PetscFunctionBegin;
657   if(icntl <= 64){
658     mat_mkl_pardiso->iparm[icntl - 1] = ival;
659   } else {
660     if(icntl == 65)
661       mkl_set_num_threads((int)ival);
662     else if(icntl == 66)
663       mat_mkl_pardiso->maxfct = ival;
664     else if(icntl == 67)
665       mat_mkl_pardiso->mnum = ival;
666     else if(icntl == 68)
667       mat_mkl_pardiso->msglvl = ival;
668     else if(icntl == 69){
669       int pt[IPARM_SIZE];
670       mat_mkl_pardiso->mtype = ival;
671       MKL_PARDISO_INIT(&pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
672 #if defined(PETSC_USE_REAL_SINGLE)
673       mat_mkl_pardiso->iparm[27] = 1;
674 #else
675       mat_mkl_pardiso->iparm[27] = 0;
676 #endif
677       mat_mkl_pardiso->iparm[34] = 1;
678     }
679   }
680   PetscFunctionReturn(0);
681 }
682 
683 #undef __FUNCT__
684 #define __FUNCT__ "MatMkl_PardisoSetCntl"
685 /*@
686   MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters
687 
688    Logically Collective on Mat
689 
690    Input Parameters:
691 +  F - the factored matrix obtained by calling MatGetFactor()
692 .  icntl - index of Mkl_Pardiso parameter
693 -  ival - value of Mkl_Pardiso parameter
694 
695   Options Database:
696 .   -mat_mkl_pardiso_<icntl> <ival>
697 
698    Level: beginner
699 
700    References: Mkl_Pardiso Users' Guide
701 
702 .seealso: MatGetFactor()
703 @*/
704 PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival)
705 {
706   PetscErrorCode ierr;
707 
708   PetscFunctionBegin;
709   ierr = PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr);
710   PetscFunctionReturn(0);
711 }
712 
713 
714 /*MC
715   MATSOLVERMKL_PARDISO -  A matrix type providing direct solvers (LU) for
716   sequential matrices via the external package MKL_PARDISO.
717 
718   Works with MATSEQAIJ matrices
719 
720   Options Database Keys:
721 + -mat_mkl_pardiso_65 - Number of thrads to use
722 . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
723 . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase
724 . -mat_mkl_pardiso_68 - Message level information
725 . -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
726 . -mat_mkl_pardiso_1 - Use default values
727 . -mat_mkl_pardiso_2 - Fill-in reducing ordering for the input matrix
728 . -mat_mkl_pardiso_4 - Preconditioned CGS/CG
729 . -mat_mkl_pardiso_5 - User permutation
730 . -mat_mkl_pardiso_6 - Write solution on x
731 . -mat_mkl_pardiso_8 - Iterative refinement step
732 . -mat_mkl_pardiso_10 - Pivoting perturbation
733 . -mat_mkl_pardiso_11 - Scaling vectors
734 . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A
735 . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching
736 . -mat_mkl_pardiso_18 - Numbers of non-zero elements
737 . -mat_mkl_pardiso_19 - Report number of floating point operations
738 . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices
739 . -mat_mkl_pardiso_24 - Parallel factorization control
740 . -mat_mkl_pardiso_25 - Parallel forward/backward solve control
741 . -mat_mkl_pardiso_27 - Matrix checker
742 . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors
743 . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
744 - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode
745 
746   Level: beginner
747 
748   For more information please check  mkl_pardiso manual
749 
750 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
751 
752 M*/
753 
754 
755 #undef __FUNCT__
756 #define __FUNCT__ "MatFactorGetSolverPackage_mkl_pardiso"
757 static PetscErrorCode MatFactorGetSolverPackage_mkl_pardiso(Mat A, const MatSolverPackage *type){
758   PetscFunctionBegin;
759   *type = MATSOLVERMKL_PARDISO;
760   PetscFunctionReturn(0);
761 }
762 
763 /* MatGetFactor for Seq sbAIJ matrices */
764 #undef __FUNCT__
765 #define __FUNCT__ "MatGetFactor_sbaij_mkl_pardiso"
766 PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F){
767   Mat            B;
768   PetscErrorCode ierr;
769   Mat_MKL_PARDISO *mat_mkl_pardiso;
770   PetscBool      isSeqSBAIJ;
771   PetscInt       bs;
772 
773   PetscFunctionBegin;
774   /* Create the factorization matrix */
775 
776 
777   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
778   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
779   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
780   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
781   if (isSeqSBAIJ) {
782     ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr);
783   } else {
784     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Is not allowed other types of matrices apart from MATSEQSBAIJ.");
785   }
786 
787   ierr = MatGetBlockSize(A,&bs); CHKERRQ(ierr);
788 
789   if(bs != 1)
790   {
791     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrice MATSEQSBAIJ with block size other than 1 is not supported by Pardiso");
792   }
793 
794   if(ftype != MAT_FACTOR_CHOLESKY)
795   {
796     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrice MATSEQAIJ should be used only with MAT_FACTOR_CHOLESKY.");
797   }
798 
799   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO;
800   B->factortype                  = MAT_FACTOR_CHOLESKY;
801 
802   B->ops->destroy = MatDestroy_MKL_PARDISO;
803   B->ops->view    = MatView_MKL_PARDISO;
804   B->factortype   = ftype;
805   B->ops->getinfo = MatGetInfo_MKL_PARDISO;
806   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
807 
808   ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr);
809   B->spptr = mat_mkl_pardiso;
810   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr);
811   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr);
812   ierr = PetscInitializeMKL_PARDISO(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr);
813 
814   *F = B;
815   PetscFunctionReturn(0);
816 }
817 
818 /* MatGetFactor for Seq AIJ matrices */
819 #undef __FUNCT__
820 #define __FUNCT__ "MatGetFactor_aij_mkl_pardiso"
821 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F){
822   Mat            B;
823   PetscErrorCode ierr;
824   Mat_MKL_PARDISO *mat_mkl_pardiso;
825   PetscBool      isSeqAIJ;
826 
827   PetscFunctionBegin;
828   /* Create the factorization matrix */
829 
830 
831   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
832   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
833   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
834   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
835   if (isSeqAIJ) {
836     ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr);
837   } else {
838     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Is not allowed other types of matrices apart from MATSEQAIJ.");
839   }
840 
841   if(ftype != MAT_FACTOR_LU)
842   {
843     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrice MATSEQAIJ should be used only with MAT_FACTOR_LU.");
844   }
845 
846   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO;
847   B->factortype            = MAT_FACTOR_LU;
848 
849   B->ops->destroy = MatDestroy_MKL_PARDISO;
850   B->ops->view    = MatView_MKL_PARDISO;
851   B->factortype   = ftype;
852   B->ops->getinfo = MatGetInfo_MKL_PARDISO;
853   B->assembled    = PETSC_TRUE;           /* required by -ksp_view */
854 
855   ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr);
856   B->spptr = mat_mkl_pardiso;
857   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr);
858   ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr);
859   ierr = PetscInitializeMKL_PARDISO(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr);
860 
861   *F = B;
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "MatSolverPackageRegister_MKL_Pardiso"
867 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MKL_Pardiso(void)
868 {
869   PetscErrorCode ierr;
870 
871   PetscFunctionBegin;
872   ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,   MAT_FACTOR_LU,      MatGetFactor_aij_mkl_pardiso  );CHKERRQ(ierr);
873   ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mkl_pardiso);CHKERRQ(ierr);
874 
875   PetscFunctionReturn(0);
876 }
877 
878