Actual source code: mkl_pardiso.c
petsc-3.12.5 2020-03-29
1: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
3: #include <../src/mat/impls/dense/seq/dense.h>
5: #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
6: #define MKL_ILP64
7: #endif
8: #include <mkl_pardiso.h>
10: PETSC_EXTERN void PetscSetMKL_PARDISOThreads(int);
12: /*
13: * Possible mkl_pardiso phases that controls the execution of the solver.
14: * For more information check mkl_pardiso manual.
15: */
16: #define JOB_ANALYSIS 11
17: #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION 12
18: #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13
19: #define JOB_NUMERICAL_FACTORIZATION 22
20: #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 23
21: #define JOB_SOLVE_ITERATIVE_REFINEMENT 33
22: #define JOB_SOLVE_FORWARD_SUBSTITUTION 331
23: #define JOB_SOLVE_DIAGONAL_SUBSTITUTION 332
24: #define JOB_SOLVE_BACKWARD_SUBSTITUTION 333
25: #define JOB_RELEASE_OF_LU_MEMORY 0
26: #define JOB_RELEASE_OF_ALL_MEMORY -1
28: #define IPARM_SIZE 64
30: #if defined(PETSC_USE_64BIT_INDICES)
31: #if defined(PETSC_HAVE_MKL_INTEL_ILP64)
32: #define INT_TYPE long long int
33: #define MKL_PARDISO pardiso
34: #define MKL_PARDISO_INIT pardisoinit
35: #else
36: /* this is the case where the MKL BLAS/LAPACK are 32 bit integers but the 64 bit integer version of
37: of Pardiso code is used; hence the need for the 64 below*/
38: #define INT_TYPE long long int
39: #define MKL_PARDISO pardiso_64
40: #define MKL_PARDISO_INIT pardiso_64init
41: void pardiso_64init(void *pt, INT_TYPE *mtype, INT_TYPE iparm [])
42: {
43: int iparm_copy[IPARM_SIZE], mtype_copy, i;
45: mtype_copy = *mtype;
46: pardisoinit(pt, &mtype_copy, iparm_copy);
47: for (i=0; i<IPARM_SIZE; i++) iparm[i] = iparm_copy[i];
48: }
49: #endif
50: #else
51: #define INT_TYPE int
52: #define MKL_PARDISO pardiso
53: #define MKL_PARDISO_INIT pardisoinit
54: #endif
57: /*
58: * Internal data structure.
59: * For more information check mkl_pardiso manual.
60: */
61: typedef struct {
63: /* Configuration vector*/
64: INT_TYPE iparm[IPARM_SIZE];
66: /*
67: * Internal mkl_pardiso memory location.
68: * After the first call to mkl_pardiso do not modify pt, as that could cause a serious memory leak.
69: */
70: void *pt[IPARM_SIZE];
72: /* Basic mkl_pardiso info*/
73: INT_TYPE phase, maxfct, mnum, mtype, n, nrhs, msglvl, err;
75: /* Matrix structure*/
76: void *a;
77: INT_TYPE *ia, *ja;
79: /* Number of non-zero elements*/
80: INT_TYPE nz;
82: /* Row permutaton vector*/
83: INT_TYPE *perm;
85: /* Define if matrix preserves sparse structure.*/
86: MatStructure matstruc;
88: PetscBool needsym;
89: PetscBool freeaij;
91: /* Schur complement */
92: PetscScalar *schur;
93: PetscInt schur_size;
94: PetscInt *schur_idxs;
95: PetscScalar *schur_work;
96: PetscBLASInt schur_work_size;
97: PetscBool solve_interior;
99: /* True if mkl_pardiso function have been used.*/
100: PetscBool CleanUp;
102: /* Conversion to a format suitable for MKL */
103: PetscErrorCode (*Convert)(Mat, PetscBool, MatReuse, PetscBool*, INT_TYPE*, INT_TYPE**, INT_TYPE**, PetscScalar**);
104: } Mat_MKL_PARDISO;
106: PetscErrorCode MatMKLPardiso_Convert_seqsbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
107: {
108: Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data;
109: PetscInt bs = A->rmap->bs,i;
113: if (!sym) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
114: *v = aa->a;
115: if (bs == 1) { /* already in the correct format */
116: /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
117: *r = (INT_TYPE*)aa->i;
118: *c = (INT_TYPE*)aa->j;
119: *nnz = (INT_TYPE)aa->nz;
120: *free = PETSC_FALSE;
121: } else if (reuse == MAT_INITIAL_MATRIX) {
122: PetscInt m = A->rmap->n,nz = aa->nz;
123: PetscInt *row,*col;
124: PetscMalloc2(m+1,&row,nz,&col);
125: for (i=0; i<m+1; i++) {
126: row[i] = aa->i[i]+1;
127: }
128: for (i=0; i<nz; i++) {
129: col[i] = aa->j[i]+1;
130: }
131: *r = (INT_TYPE*)row;
132: *c = (INT_TYPE*)col;
133: *nnz = (INT_TYPE)nz;
134: *free = PETSC_TRUE;
135: }
136: return(0);
137: }
139: PetscErrorCode MatMKLPardiso_Convert_seqbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
140: {
141: Mat_SeqBAIJ *aa = (Mat_SeqBAIJ*)A->data;
142: PetscInt bs = A->rmap->bs,i;
146: if (!sym) {
147: *v = aa->a;
148: if (bs == 1) { /* already in the correct format */
149: /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
150: *r = (INT_TYPE*)aa->i;
151: *c = (INT_TYPE*)aa->j;
152: *nnz = (INT_TYPE)aa->nz;
153: *free = PETSC_FALSE;
154: return(0);
155: } else if (reuse == MAT_INITIAL_MATRIX) {
156: PetscInt m = A->rmap->n,nz = aa->nz;
157: PetscInt *row,*col;
158: PetscMalloc2(m+1,&row,nz,&col);
159: for (i=0; i<m+1; i++) {
160: row[i] = aa->i[i]+1;
161: }
162: for (i=0; i<nz; i++) {
163: col[i] = aa->j[i]+1;
164: }
165: *r = (INT_TYPE*)row;
166: *c = (INT_TYPE*)col;
167: *nnz = (INT_TYPE)nz;
168: }
169: *free = PETSC_TRUE;
170: } else {
171: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
172: }
173: return(0);
174: }
176: PetscErrorCode MatMKLPardiso_Convert_seqaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
177: {
178: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data;
179: PetscScalar *aav;
183: MatSeqAIJGetArrayRead(A,(const PetscScalar**)&aav);
184: if (!sym) { /* already in the correct format */
185: *v = aav;
186: *r = (INT_TYPE*)aa->i;
187: *c = (INT_TYPE*)aa->j;
188: *nnz = (INT_TYPE)aa->nz;
189: *free = PETSC_FALSE;
190: } else if (reuse == MAT_INITIAL_MATRIX) { /* need to get the triangular part */
191: PetscScalar *vals,*vv;
192: PetscInt *row,*col,*jj;
193: PetscInt m = A->rmap->n,nz,i;
195: nz = 0;
196: for (i=0; i<m; i++) nz += aa->i[i+1] - aa->diag[i];
197: PetscMalloc2(m+1,&row,nz,&col);
198: PetscMalloc1(nz,&vals);
199: jj = col;
200: vv = vals;
202: row[0] = 0;
203: for (i=0; i<m; i++) {
204: PetscInt *aj = aa->j + aa->diag[i];
205: PetscScalar *av = aav + aa->diag[i];
206: PetscInt rl = aa->i[i+1] - aa->diag[i],j;
208: for (j=0; j<rl; j++) {
209: *jj = *aj; jj++; aj++;
210: *vv = *av; vv++; av++;
211: }
212: row[i+1] = row[i] + rl;
213: }
214: *v = vals;
215: *r = (INT_TYPE*)row;
216: *c = (INT_TYPE*)col;
217: *nnz = (INT_TYPE)nz;
218: *free = PETSC_TRUE;
219: } else {
220: PetscScalar *vv;
221: PetscInt m = A->rmap->n,i;
223: vv = *v;
224: for (i=0; i<m; i++) {
225: PetscScalar *av = aav + aa->diag[i];
226: PetscInt rl = aa->i[i+1] - aa->diag[i],j;
227: for (j=0; j<rl; j++) {
228: *vv = *av; vv++; av++;
229: }
230: }
231: *free = PETSC_TRUE;
232: }
233: MatSeqAIJRestoreArrayRead(A,(const PetscScalar**)&aav);
234: return(0);
235: }
238: static PetscErrorCode MatMKLPardisoSolveSchur_Private(Mat F, PetscScalar *B, PetscScalar *X)
239: {
240: Mat_MKL_PARDISO *mpardiso = (Mat_MKL_PARDISO*)F->data;
241: Mat S,Xmat,Bmat;
242: MatFactorSchurStatus schurstatus;
243: PetscErrorCode ierr;
246: MatFactorFactorizeSchurComplement(F);
247: MatFactorGetSchurComplement(F,&S,&schurstatus);
248: if (X == B && schurstatus == MAT_FACTOR_SCHUR_INVERTED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"X and B cannot point to the same address");
249: MatCreateSeqDense(PETSC_COMM_SELF,mpardiso->schur_size,mpardiso->nrhs,B,&Bmat);
250: MatCreateSeqDense(PETSC_COMM_SELF,mpardiso->schur_size,mpardiso->nrhs,X,&Xmat);
251: MatSetType(Bmat,((PetscObject)S)->type_name);
252: MatSetType(Xmat,((PetscObject)S)->type_name);
253: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
254: MatPinToCPU(Xmat,S->pinnedtocpu);
255: MatPinToCPU(Bmat,S->pinnedtocpu);
256: #endif
258: #if defined(PETSC_USE_COMPLEX)
259: if (mpardiso->iparm[12-1] == 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Hermitian solve not implemented yet");
260: #endif
262: switch (schurstatus) {
263: case MAT_FACTOR_SCHUR_FACTORED:
264: if (!mpardiso->iparm[12-1]) {
265: MatMatSolve(S,Bmat,Xmat);
266: } else { /* transpose solve */
267: MatMatSolveTranspose(S,Bmat,Xmat);
268: }
269: break;
270: case MAT_FACTOR_SCHUR_INVERTED:
271: if (!mpardiso->iparm[12-1]) {
272: MatMatMult(S,Bmat,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Xmat);
273: } else { /* transpose solve */
274: MatTransposeMatMult(S,Bmat,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Xmat);
275: }
276: break;
277: default:
278: SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
279: break;
280: }
281: MatFactorRestoreSchurComplement(F,&S,schurstatus);
282: MatDestroy(&Bmat);
283: MatDestroy(&Xmat);
284: return(0);
285: }
287: PetscErrorCode MatFactorSetSchurIS_MKL_PARDISO(Mat F, IS is)
288: {
289: Mat_MKL_PARDISO *mpardiso = (Mat_MKL_PARDISO*)F->data;
290: const PetscScalar *arr;
291: const PetscInt *idxs;
292: PetscInt size,i;
293: PetscMPIInt csize;
294: PetscBool sorted;
295: PetscErrorCode ierr;
298: MPI_Comm_size(PetscObjectComm((PetscObject)F),&csize);
299: if (csize > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MKL_PARDISO parallel Schur complements not yet supported from PETSc");
300: ISSorted(is,&sorted);
301: if (!sorted) {
302: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IS for MKL_PARDISO Schur complements needs to be sorted");
303: }
304: ISGetLocalSize(is,&size);
305: PetscFree(mpardiso->schur_work);
306: PetscBLASIntCast(PetscMax(mpardiso->n,2*size),&mpardiso->schur_work_size);
307: PetscMalloc1(mpardiso->schur_work_size,&mpardiso->schur_work);
308: MatDestroy(&F->schur);
309: MatCreateSeqDense(PETSC_COMM_SELF,size,size,NULL,&F->schur);
310: MatDenseGetArrayRead(F->schur,&arr);
311: mpardiso->schur = (PetscScalar*)arr;
312: mpardiso->schur_size = size;
313: MatDenseRestoreArrayRead(F->schur,&arr);
314: if (mpardiso->mtype == 2) {
315: MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);
316: }
318: PetscFree(mpardiso->schur_idxs);
319: PetscMalloc1(size,&mpardiso->schur_idxs);
320: PetscArrayzero(mpardiso->perm,mpardiso->n);
321: ISGetIndices(is,&idxs);
322: PetscArraycpy(mpardiso->schur_idxs,idxs,size);
323: for (i=0;i<size;i++) mpardiso->perm[idxs[i]] = 1;
324: ISRestoreIndices(is,&idxs);
325: if (size) { /* turn on Schur switch if the set of indices is not empty */
326: mpardiso->iparm[36-1] = 2;
327: }
328: return(0);
329: }
331: PetscErrorCode MatDestroy_MKL_PARDISO(Mat A)
332: {
333: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
334: PetscErrorCode ierr;
337: if (mat_mkl_pardiso->CleanUp) {
338: mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;
340: MKL_PARDISO (mat_mkl_pardiso->pt,
341: &mat_mkl_pardiso->maxfct,
342: &mat_mkl_pardiso->mnum,
343: &mat_mkl_pardiso->mtype,
344: &mat_mkl_pardiso->phase,
345: &mat_mkl_pardiso->n,
346: NULL,
347: NULL,
348: NULL,
349: NULL,
350: &mat_mkl_pardiso->nrhs,
351: mat_mkl_pardiso->iparm,
352: &mat_mkl_pardiso->msglvl,
353: NULL,
354: NULL,
355: &mat_mkl_pardiso->err);
356: }
357: PetscFree(mat_mkl_pardiso->perm);
358: PetscFree(mat_mkl_pardiso->schur_work);
359: PetscFree(mat_mkl_pardiso->schur_idxs);
360: if (mat_mkl_pardiso->freeaij) {
361: PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);
362: if (mat_mkl_pardiso->iparm[34] == 1) {
363: PetscFree(mat_mkl_pardiso->a);
364: }
365: }
366: PetscFree(A->data);
368: /* clear composed functions */
369: PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
370: PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
371: PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);
372: return(0);
373: }
375: static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce)
376: {
378: if (reduce) { /* data given for the whole matrix */
379: PetscInt i,m=0,p=0;
380: for (i=0;i<mpardiso->nrhs;i++) {
381: PetscInt j;
382: for (j=0;j<mpardiso->schur_size;j++) {
383: schur[p+j] = whole[m+mpardiso->schur_idxs[j]];
384: }
385: m += mpardiso->n;
386: p += mpardiso->schur_size;
387: }
388: } else { /* from Schur to whole */
389: PetscInt i,m=0,p=0;
390: for (i=0;i<mpardiso->nrhs;i++) {
391: PetscInt j;
392: for (j=0;j<mpardiso->schur_size;j++) {
393: whole[m+mpardiso->schur_idxs[j]] = schur[p+j];
394: }
395: m += mpardiso->n;
396: p += mpardiso->schur_size;
397: }
398: }
399: return(0);
400: }
402: PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x)
403: {
404: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
405: PetscErrorCode ierr;
406: PetscScalar *xarray;
407: const PetscScalar *barray;
410: mat_mkl_pardiso->nrhs = 1;
411: VecGetArray(x,&xarray);
412: VecGetArrayRead(b,&barray);
414: if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
415: else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
417: if (barray == xarray) { /* if the two vectors share the same memory */
418: PetscScalar *work;
419: if (!mat_mkl_pardiso->schur_work) {
420: PetscMalloc1(mat_mkl_pardiso->n,&work);
421: } else {
422: work = mat_mkl_pardiso->schur_work;
423: }
424: mat_mkl_pardiso->iparm[6-1] = 1;
425: MKL_PARDISO (mat_mkl_pardiso->pt,
426: &mat_mkl_pardiso->maxfct,
427: &mat_mkl_pardiso->mnum,
428: &mat_mkl_pardiso->mtype,
429: &mat_mkl_pardiso->phase,
430: &mat_mkl_pardiso->n,
431: mat_mkl_pardiso->a,
432: mat_mkl_pardiso->ia,
433: mat_mkl_pardiso->ja,
434: NULL,
435: &mat_mkl_pardiso->nrhs,
436: mat_mkl_pardiso->iparm,
437: &mat_mkl_pardiso->msglvl,
438: (void*)xarray,
439: (void*)work,
440: &mat_mkl_pardiso->err);
441: if (!mat_mkl_pardiso->schur_work) {
442: PetscFree(work);
443: }
444: } else {
445: mat_mkl_pardiso->iparm[6-1] = 0;
446: MKL_PARDISO (mat_mkl_pardiso->pt,
447: &mat_mkl_pardiso->maxfct,
448: &mat_mkl_pardiso->mnum,
449: &mat_mkl_pardiso->mtype,
450: &mat_mkl_pardiso->phase,
451: &mat_mkl_pardiso->n,
452: mat_mkl_pardiso->a,
453: mat_mkl_pardiso->ia,
454: mat_mkl_pardiso->ja,
455: mat_mkl_pardiso->perm,
456: &mat_mkl_pardiso->nrhs,
457: mat_mkl_pardiso->iparm,
458: &mat_mkl_pardiso->msglvl,
459: (void*)barray,
460: (void*)xarray,
461: &mat_mkl_pardiso->err);
462: }
463: VecRestoreArrayRead(b,&barray);
465: if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
467: if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
468: PetscInt shift = mat_mkl_pardiso->schur_size;
470: /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
471: if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;
473: if (!mat_mkl_pardiso->solve_interior) {
474: /* solve Schur complement */
475: MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);
476: MatMKLPardisoSolveSchur_Private(A,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);
477: MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);
478: } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substituted to xarray[schur] will be neglected */
479: PetscInt i;
480: for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
481: xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.;
482: }
483: }
485: /* expansion phase */
486: mat_mkl_pardiso->iparm[6-1] = 1;
487: mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
488: MKL_PARDISO (mat_mkl_pardiso->pt,
489: &mat_mkl_pardiso->maxfct,
490: &mat_mkl_pardiso->mnum,
491: &mat_mkl_pardiso->mtype,
492: &mat_mkl_pardiso->phase,
493: &mat_mkl_pardiso->n,
494: mat_mkl_pardiso->a,
495: mat_mkl_pardiso->ia,
496: mat_mkl_pardiso->ja,
497: mat_mkl_pardiso->perm,
498: &mat_mkl_pardiso->nrhs,
499: mat_mkl_pardiso->iparm,
500: &mat_mkl_pardiso->msglvl,
501: (void*)xarray,
502: (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
503: &mat_mkl_pardiso->err);
505: if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
506: mat_mkl_pardiso->iparm[6-1] = 0;
507: }
508: VecRestoreArray(x,&xarray);
509: mat_mkl_pardiso->CleanUp = PETSC_TRUE;
510: return(0);
511: }
513: PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x)
514: {
515: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
516: PetscInt oiparm12;
517: PetscErrorCode ierr;
520: oiparm12 = mat_mkl_pardiso->iparm[12 - 1];
521: mat_mkl_pardiso->iparm[12 - 1] = 2;
522: MatSolve_MKL_PARDISO(A,b,x);
523: mat_mkl_pardiso->iparm[12 - 1] = oiparm12;
524: return(0);
525: }
527: PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X)
528: {
529: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->data;
530: PetscErrorCode ierr;
531: const PetscScalar *barray;
532: PetscScalar *xarray;
533: PetscBool flg;
536: PetscObjectBaseTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
537: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix");
538: PetscObjectBaseTypeCompare((PetscObject)X,MATSEQDENSE,&flg);
539: if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix");
541: MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);
543: if (mat_mkl_pardiso->nrhs > 0) {
544: MatDenseGetArrayRead(B,&barray);
545: MatDenseGetArray(X,&xarray);
547: if (barray == xarray) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"B and X cannot share the same memory location");
548: if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
549: else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
550: mat_mkl_pardiso->iparm[6-1] = 0;
552: MKL_PARDISO (mat_mkl_pardiso->pt,
553: &mat_mkl_pardiso->maxfct,
554: &mat_mkl_pardiso->mnum,
555: &mat_mkl_pardiso->mtype,
556: &mat_mkl_pardiso->phase,
557: &mat_mkl_pardiso->n,
558: mat_mkl_pardiso->a,
559: mat_mkl_pardiso->ia,
560: mat_mkl_pardiso->ja,
561: mat_mkl_pardiso->perm,
562: &mat_mkl_pardiso->nrhs,
563: mat_mkl_pardiso->iparm,
564: &mat_mkl_pardiso->msglvl,
565: (void*)barray,
566: (void*)xarray,
567: &mat_mkl_pardiso->err);
568: if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
570: MatDenseRestoreArrayRead(B,&barray);
571: if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
572: PetscScalar *o_schur_work = NULL;
573: PetscInt shift = mat_mkl_pardiso->schur_size*mat_mkl_pardiso->nrhs,scale;
574: PetscInt mem = mat_mkl_pardiso->n*mat_mkl_pardiso->nrhs;
576: /* allocate extra memory if it is needed */
577: scale = 1;
578: if (A->schur_status == MAT_FACTOR_SCHUR_INVERTED) scale = 2;
580: mem *= scale;
581: if (mem > mat_mkl_pardiso->schur_work_size) {
582: o_schur_work = mat_mkl_pardiso->schur_work;
583: PetscMalloc1(mem,&mat_mkl_pardiso->schur_work);
584: }
586: /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
587: if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;
589: /* solve Schur complement */
590: if (!mat_mkl_pardiso->solve_interior) {
591: MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);
592: MatMKLPardisoSolveSchur_Private(A,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);
593: MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);
594: } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substituted to xarray[schur,n] will be neglected */
595: PetscInt i,n,m=0;
596: for (n=0;n<mat_mkl_pardiso->nrhs;n++) {
597: for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
598: xarray[mat_mkl_pardiso->schur_idxs[i]+m] = 0.;
599: }
600: m += mat_mkl_pardiso->n;
601: }
602: }
604: /* expansion phase */
605: mat_mkl_pardiso->iparm[6-1] = 1;
606: mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
607: MKL_PARDISO (mat_mkl_pardiso->pt,
608: &mat_mkl_pardiso->maxfct,
609: &mat_mkl_pardiso->mnum,
610: &mat_mkl_pardiso->mtype,
611: &mat_mkl_pardiso->phase,
612: &mat_mkl_pardiso->n,
613: mat_mkl_pardiso->a,
614: mat_mkl_pardiso->ia,
615: mat_mkl_pardiso->ja,
616: mat_mkl_pardiso->perm,
617: &mat_mkl_pardiso->nrhs,
618: mat_mkl_pardiso->iparm,
619: &mat_mkl_pardiso->msglvl,
620: (void*)xarray,
621: (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
622: &mat_mkl_pardiso->err);
623: if (o_schur_work) { /* restore original schur_work (minimal size) */
624: PetscFree(mat_mkl_pardiso->schur_work);
625: mat_mkl_pardiso->schur_work = o_schur_work;
626: }
627: if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
628: mat_mkl_pardiso->iparm[6-1] = 0;
629: }
630: }
631: mat_mkl_pardiso->CleanUp = PETSC_TRUE;
632: return(0);
633: }
635: PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info)
636: {
637: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->data;
638: PetscErrorCode ierr;
641: mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
642: (*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,(PetscScalar**)&mat_mkl_pardiso->a);
644: mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION;
645: MKL_PARDISO (mat_mkl_pardiso->pt,
646: &mat_mkl_pardiso->maxfct,
647: &mat_mkl_pardiso->mnum,
648: &mat_mkl_pardiso->mtype,
649: &mat_mkl_pardiso->phase,
650: &mat_mkl_pardiso->n,
651: mat_mkl_pardiso->a,
652: mat_mkl_pardiso->ia,
653: mat_mkl_pardiso->ja,
654: mat_mkl_pardiso->perm,
655: &mat_mkl_pardiso->nrhs,
656: mat_mkl_pardiso->iparm,
657: &mat_mkl_pardiso->msglvl,
658: NULL,
659: (void*)mat_mkl_pardiso->schur,
660: &mat_mkl_pardiso->err);
661: if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
663: /* report flops */
664: if (mat_mkl_pardiso->iparm[18] > 0) {
665: PetscLogFlops(PetscPowRealInt(10.,6)*mat_mkl_pardiso->iparm[18]);
666: }
668: if (F->schur) { /* schur output from pardiso is in row major format */
669: #if defined(PETSC_HAVE_CUDA)
670: F->schur->offloadmask = PETSC_OFFLOAD_CPU;
671: #endif
672: MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);
673: MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);
674: }
675: mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
676: mat_mkl_pardiso->CleanUp = PETSC_TRUE;
677: return(0);
678: }
680: PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A)
681: {
682: Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
683: PetscErrorCode ierr;
684: PetscInt icntl,threads=1;
685: PetscBool flg;
688: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");
690: PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use within PARDISO","None",threads,&threads,&flg);
691: if (flg) PetscSetMKL_PARDISOThreads((int)threads);
693: 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);
694: if (flg) mat_mkl_pardiso->maxfct = icntl;
696: PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);
697: if (flg) mat_mkl_pardiso->mnum = icntl;
699: PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);
700: if (flg) mat_mkl_pardiso->msglvl = icntl;
702: PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);
703: if (flg) {
704: void *pt[IPARM_SIZE];
705: mat_mkl_pardiso->mtype = icntl;
706: MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
707: #if defined(PETSC_USE_REAL_SINGLE)
708: mat_mkl_pardiso->iparm[27] = 1;
709: #else
710: mat_mkl_pardiso->iparm[27] = 0;
711: #endif
712: mat_mkl_pardiso->iparm[34] = 1; /* use 0-based indexing */
713: }
714: PetscOptionsInt("-mat_mkl_pardiso_1","Use default values (if 0)","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);
716: if (flg && icntl != 0) {
717: PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);
718: if (flg) mat_mkl_pardiso->iparm[1] = icntl;
720: PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);
721: if (flg) mat_mkl_pardiso->iparm[3] = icntl;
723: PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);
724: if (flg) mat_mkl_pardiso->iparm[4] = icntl;
726: PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);
727: if (flg) mat_mkl_pardiso->iparm[5] = icntl;
729: PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);
730: if (flg) mat_mkl_pardiso->iparm[7] = icntl;
732: PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);
733: if (flg) mat_mkl_pardiso->iparm[9] = icntl;
735: PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);
736: if (flg) mat_mkl_pardiso->iparm[10] = icntl;
738: PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);
739: if (flg) mat_mkl_pardiso->iparm[11] = icntl;
741: PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);
742: if (flg) mat_mkl_pardiso->iparm[12] = icntl;
744: PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);
745: if (flg) mat_mkl_pardiso->iparm[17] = icntl;
747: PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations (0 to disable)","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);
748: if (flg) mat_mkl_pardiso->iparm[18] = icntl;
750: PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);
751: if (flg) mat_mkl_pardiso->iparm[20] = icntl;
753: PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);
754: if (flg) mat_mkl_pardiso->iparm[23] = icntl;
756: PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);
757: if (flg) mat_mkl_pardiso->iparm[24] = icntl;
759: PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);
760: if (flg) mat_mkl_pardiso->iparm[26] = icntl;
762: PetscOptionsInt("-mat_mkl_pardiso_31","Partial solve and computing selected components of the solution vectors","None",mat_mkl_pardiso->iparm[30],&icntl,&flg);
763: if (flg) mat_mkl_pardiso->iparm[30] = icntl;
765: PetscOptionsInt("-mat_mkl_pardiso_34","Optimal number of threads for conditional numerical reproducibility (CNR) mode","None",mat_mkl_pardiso->iparm[33],&icntl,&flg);
766: if (flg) mat_mkl_pardiso->iparm[33] = icntl;
768: PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);
769: if (flg) mat_mkl_pardiso->iparm[59] = icntl;
770: }
771: PetscOptionsEnd();
772: return(0);
773: }
775: PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso)
776: {
778: PetscInt i,bs;
779: PetscBool match;
782: for (i=0; i<IPARM_SIZE; i++) mat_mkl_pardiso->iparm[i] = 0;
783: for (i=0; i<IPARM_SIZE; i++) mat_mkl_pardiso->pt[i] = 0;
784: /* Default options for both sym and unsym */
785: mat_mkl_pardiso->iparm[ 0] = 1; /* Solver default parameters overriden with provided by iparm */
786: mat_mkl_pardiso->iparm[ 1] = 2; /* Metis reordering */
787: mat_mkl_pardiso->iparm[ 5] = 0; /* Write solution into x */
788: mat_mkl_pardiso->iparm[ 7] = 0; /* Max number of iterative refinement steps */
789: mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
790: mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
791: #if 0
792: mat_mkl_pardiso->iparm[23] = 1; /* Parallel factorization control*/
793: #endif
794: PetscObjectTypeCompareAny((PetscObject)A,&match,MATSEQBAIJ,MATSEQSBAIJ,"");
795: MatGetBlockSize(A,&bs);
796: if (!match || bs == 1) {
797: mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */
798: mat_mkl_pardiso->n = A->rmap->N;
799: } else {
800: mat_mkl_pardiso->iparm[34] = 0; /* Cluster Sparse Solver use Fortran-style indexing for ia and ja arrays */
801: mat_mkl_pardiso->iparm[36] = bs;
802: mat_mkl_pardiso->n = A->rmap->N/bs;
803: }
804: mat_mkl_pardiso->iparm[39] = 0; /* Input: matrix/rhs/solution stored on master */
806: mat_mkl_pardiso->CleanUp = PETSC_FALSE;
807: mat_mkl_pardiso->maxfct = 1; /* Maximum number of numerical factorizations. */
808: mat_mkl_pardiso->mnum = 1; /* Which factorization to use. */
809: mat_mkl_pardiso->msglvl = 0; /* 0: do not print 1: Print statistical information in file */
810: mat_mkl_pardiso->phase = -1;
811: mat_mkl_pardiso->err = 0;
813: mat_mkl_pardiso->nrhs = 1;
814: mat_mkl_pardiso->err = 0;
815: mat_mkl_pardiso->phase = -1;
817: if (ftype == MAT_FACTOR_LU) {
818: mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */
819: mat_mkl_pardiso->iparm[10] = 1; /* Use nonsymmetric permutation and scaling MPS */
820: mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
821: } else {
822: mat_mkl_pardiso->iparm[ 9] = 8; /* Perturb the pivot elements with 1E-8 */
823: mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */
824: mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
825: #if defined(PETSC_USE_DEBUG)
826: mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */
827: #endif
828: }
829: PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);
830: for (i=0; i<A->rmap->N; i++) {
831: mat_mkl_pardiso->perm[i] = 0;
832: }
833: mat_mkl_pardiso->schur_size = 0;
834: return(0);
835: }
837: PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info)
838: {
839: Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
840: PetscErrorCode ierr;
843: mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
844: PetscSetMKL_PARDISOFromOptions(F,A);
846: /* throw away any previously computed structure */
847: if (mat_mkl_pardiso->freeaij) {
848: PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);
849: if (mat_mkl_pardiso->iparm[34] == 1) {
850: PetscFree(mat_mkl_pardiso->a);
851: }
852: }
853: (*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,(PetscScalar**)&mat_mkl_pardiso->a);
854: if (mat_mkl_pardiso->iparm[34] == 1) mat_mkl_pardiso->n = A->rmap->N;
855: else mat_mkl_pardiso->n = A->rmap->N/A->rmap->bs;
857: mat_mkl_pardiso->phase = JOB_ANALYSIS;
859: /* reset flops counting if requested */
860: if (mat_mkl_pardiso->iparm[18]) mat_mkl_pardiso->iparm[18] = -1;
862: MKL_PARDISO (mat_mkl_pardiso->pt,
863: &mat_mkl_pardiso->maxfct,
864: &mat_mkl_pardiso->mnum,
865: &mat_mkl_pardiso->mtype,
866: &mat_mkl_pardiso->phase,
867: &mat_mkl_pardiso->n,
868: mat_mkl_pardiso->a,
869: mat_mkl_pardiso->ia,
870: mat_mkl_pardiso->ja,
871: mat_mkl_pardiso->perm,
872: &mat_mkl_pardiso->nrhs,
873: mat_mkl_pardiso->iparm,
874: &mat_mkl_pardiso->msglvl,
875: NULL,
876: NULL,
877: &mat_mkl_pardiso->err);
878: if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual",mat_mkl_pardiso->err);
880: mat_mkl_pardiso->CleanUp = PETSC_TRUE;
882: if (F->factortype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO;
883: else F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO;
885: F->ops->solve = MatSolve_MKL_PARDISO;
886: F->ops->solvetranspose = MatSolveTranspose_MKL_PARDISO;
887: F->ops->matsolve = MatMatSolve_MKL_PARDISO;
888: return(0);
889: }
891: PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
892: {
896: MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);
897: return(0);
898: }
900: #if !defined(PETSC_USE_COMPLEX)
901: PetscErrorCode MatGetInertia_MKL_PARDISO(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
902: {
903: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)F->data;
906: if (nneg) *nneg = mat_mkl_pardiso->iparm[22];
907: if (npos) *npos = mat_mkl_pardiso->iparm[21];
908: if (nzero) *nzero = F->rmap->N - (mat_mkl_pardiso->iparm[22] + mat_mkl_pardiso->iparm[21]);
909: return(0);
910: }
911: #endif
913: PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info)
914: {
918: MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);
919: #if defined(PETSC_USE_COMPLEX)
920: F->ops->getinertia = NULL;
921: #else
922: F->ops->getinertia = MatGetInertia_MKL_PARDISO;
923: #endif
924: return(0);
925: }
927: PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer)
928: {
929: PetscErrorCode ierr;
930: PetscBool iascii;
931: PetscViewerFormat format;
932: Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
933: PetscInt i;
936: if (A->ops->solve != MatSolve_MKL_PARDISO) return(0);
938: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
939: if (iascii) {
940: PetscViewerGetFormat(viewer,&format);
941: if (format == PETSC_VIEWER_ASCII_INFO) {
942: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");
943: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase: %d \n",mat_mkl_pardiso->phase);
944: for (i=1; i<=64; i++) {
945: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]: %d \n",i, mat_mkl_pardiso->iparm[i - 1]);
946: }
947: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct: %d \n", mat_mkl_pardiso->maxfct);
948: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum: %d \n", mat_mkl_pardiso->mnum);
949: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype: %d \n", mat_mkl_pardiso->mtype);
950: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n: %d \n", mat_mkl_pardiso->n);
951: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs: %d \n", mat_mkl_pardiso->nrhs);
952: PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl: %d \n", mat_mkl_pardiso->msglvl);
953: }
954: }
955: return(0);
956: }
959: PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info)
960: {
961: Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)A->data;
964: info->block_size = 1.0;
965: info->nz_used = mat_mkl_pardiso->iparm[17];
966: info->nz_allocated = mat_mkl_pardiso->iparm[17];
967: info->nz_unneeded = 0.0;
968: info->assemblies = 0.0;
969: info->mallocs = 0.0;
970: info->memory = 0.0;
971: info->fill_ratio_given = 0;
972: info->fill_ratio_needed = 0;
973: info->factor_mallocs = 0;
974: return(0);
975: }
977: PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival)
978: {
979: Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
982: if (icntl <= 64) {
983: mat_mkl_pardiso->iparm[icntl - 1] = ival;
984: } else {
985: if (icntl == 65) PetscSetMKL_PARDISOThreads(ival);
986: else if (icntl == 66) mat_mkl_pardiso->maxfct = ival;
987: else if (icntl == 67) mat_mkl_pardiso->mnum = ival;
988: else if (icntl == 68) mat_mkl_pardiso->msglvl = ival;
989: else if (icntl == 69) {
990: void *pt[IPARM_SIZE];
991: mat_mkl_pardiso->mtype = ival;
992: MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
993: #if defined(PETSC_USE_REAL_SINGLE)
994: mat_mkl_pardiso->iparm[27] = 1;
995: #else
996: mat_mkl_pardiso->iparm[27] = 0;
997: #endif
998: mat_mkl_pardiso->iparm[34] = 1;
999: } else if (icntl==70) mat_mkl_pardiso->solve_interior = (PetscBool)!!ival;
1000: }
1001: return(0);
1002: }
1004: /*@
1005: MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters
1007: Logically Collective on Mat
1009: Input Parameters:
1010: + F - the factored matrix obtained by calling MatGetFactor()
1011: . icntl - index of Mkl_Pardiso parameter
1012: - ival - value of Mkl_Pardiso parameter
1014: Options Database:
1015: . -mat_mkl_pardiso_<icntl> <ival>
1017: Level: beginner
1019: References:
1020: . Mkl_Pardiso Users' Guide
1022: .seealso: MatGetFactor()
1023: @*/
1024: PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival)
1025: {
1029: PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
1030: return(0);
1031: }
1033: /*MC
1034: MATSOLVERMKL_PARDISO - A matrix type providing direct solvers (LU) for
1035: sequential matrices via the external package MKL_PARDISO.
1037: Works with MATSEQAIJ matrices
1039: Use -pc_type lu -pc_factor_mat_solver_type mkl_pardiso to use this direct solver
1041: Options Database Keys:
1042: + -mat_mkl_pardiso_65 - Number of threads to use within MKL_PARDISO
1043: . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
1044: . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase
1045: . -mat_mkl_pardiso_68 - Message level information
1046: . -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
1047: . -mat_mkl_pardiso_1 - Use default values
1048: . -mat_mkl_pardiso_2 - Fill-in reducing ordering for the input matrix
1049: . -mat_mkl_pardiso_4 - Preconditioned CGS/CG
1050: . -mat_mkl_pardiso_5 - User permutation
1051: . -mat_mkl_pardiso_6 - Write solution on x
1052: . -mat_mkl_pardiso_8 - Iterative refinement step
1053: . -mat_mkl_pardiso_10 - Pivoting perturbation
1054: . -mat_mkl_pardiso_11 - Scaling vectors
1055: . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A
1056: . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching
1057: . -mat_mkl_pardiso_18 - Numbers of non-zero elements
1058: . -mat_mkl_pardiso_19 - Report number of floating point operations
1059: . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices
1060: . -mat_mkl_pardiso_24 - Parallel factorization control
1061: . -mat_mkl_pardiso_25 - Parallel forward/backward solve control
1062: . -mat_mkl_pardiso_27 - Matrix checker
1063: . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors
1064: . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
1065: - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode
1067: Level: beginner
1069: For more information please check mkl_pardiso manual
1071: .seealso: PCFactorSetMatSolverType(), MatSolverType
1073: M*/
1074: static PetscErrorCode MatFactorGetSolverType_mkl_pardiso(Mat A, MatSolverType *type)
1075: {
1077: *type = MATSOLVERMKL_PARDISO;
1078: return(0);
1079: }
1081: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F)
1082: {
1083: Mat B;
1084: PetscErrorCode ierr;
1085: Mat_MKL_PARDISO *mat_mkl_pardiso;
1086: PetscBool isSeqAIJ,isSeqBAIJ,isSeqSBAIJ;
1089: PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
1090: PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
1091: PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
1092: MatCreate(PetscObjectComm((PetscObject)A),&B);
1093: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1094: PetscStrallocpy("mkl_pardiso",&((PetscObject)B)->type_name);
1095: MatSetUp(B);
1097: PetscNewLog(B,&mat_mkl_pardiso);
1098: B->data = mat_mkl_pardiso;
1100: MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);
1101: if (ftype == MAT_FACTOR_LU) {
1102: B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO;
1103: B->factortype = MAT_FACTOR_LU;
1104: mat_mkl_pardiso->needsym = PETSC_FALSE;
1105: if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1106: else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1107: else if (isSeqSBAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead");
1108: else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU with %s format",((PetscObject)A)->type_name);
1109: #if defined(PETSC_USE_COMPLEX)
1110: mat_mkl_pardiso->mtype = 13;
1111: #else
1112: mat_mkl_pardiso->mtype = 11;
1113: #endif
1114: } else {
1115: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO;
1116: B->factortype = MAT_FACTOR_CHOLESKY;
1117: if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1118: else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1119: else if (isSeqSBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij;
1120: else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with %s format",((PetscObject)A)->type_name);
1122: mat_mkl_pardiso->needsym = PETSC_TRUE;
1123: #if !defined(PETSC_USE_COMPLEX)
1124: if (A->spd_set && A->spd) mat_mkl_pardiso->mtype = 2;
1125: else mat_mkl_pardiso->mtype = -2;
1126: #else
1127: mat_mkl_pardiso->mtype = 6;
1128: if (A->hermitian) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with Hermitian matrices! Use MAT_FACTOR_LU instead");
1129: #endif
1130: }
1131: B->ops->destroy = MatDestroy_MKL_PARDISO;
1132: B->ops->view = MatView_MKL_PARDISO;
1133: B->ops->getinfo = MatGetInfo_MKL_PARDISO;
1134: B->factortype = ftype;
1135: B->assembled = PETSC_TRUE;
1137: PetscFree(B->solvertype);
1138: PetscStrallocpy(MATSOLVERMKL_PARDISO,&B->solvertype);
1140: PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverType_C",MatFactorGetSolverType_mkl_pardiso);
1141: PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MKL_PARDISO);
1142: PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);
1144: *F = B;
1145: return(0);
1146: }
1148: PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MKL_Pardiso(void)
1149: {
1153: MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);
1154: MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);
1155: MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQBAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);
1156: MatSolverTypeRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);
1157: return(0);
1158: }