Actual source code: mkl_pardiso.c

petsc-3.12.0 2019-09-29
Report Typos and Errors
  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) {
114:     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
115:   }
116:   *v      = aa->a;
117:   if (bs == 1) { /* already in the correct format */
118:     /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
119:     *r    = (INT_TYPE*)aa->i;
120:     *c    = (INT_TYPE*)aa->j;
121:     *nnz  = (INT_TYPE)aa->nz;
122:     *free = PETSC_FALSE;
123:   } else if (reuse == MAT_INITIAL_MATRIX) {
124:     PetscInt m = A->rmap->n,nz = aa->nz;
125:     PetscInt *row,*col;
126:     PetscMalloc2(m+1,&row,nz,&col);
127:     for (i=0; i<m+1; i++) {
128:       row[i] = aa->i[i]+1;
129:     }
130:     for (i=0; i<nz; i++) {
131:       col[i] = aa->j[i]+1;
132:     }
133:     *r    = (INT_TYPE*)row;
134:     *c    = (INT_TYPE*)col;
135:     *nnz  = (INT_TYPE)nz;
136:     *free = PETSC_TRUE;
137:   }
138:   return(0);
139: }

141: PetscErrorCode MatMKLPardiso_Convert_seqbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
142: {
143:   Mat_SeqBAIJ    *aa = (Mat_SeqBAIJ*)A->data;
144:   PetscInt       bs  = A->rmap->bs,i;

148:   if (!sym) {
149:     *v      = aa->a;
150:     if (bs == 1) { /* already in the correct format */
151:       /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
152:       *r    = (INT_TYPE*)aa->i;
153:       *c    = (INT_TYPE*)aa->j;
154:       *nnz  = (INT_TYPE)aa->nz;
155:       *free = PETSC_FALSE;
156:       return(0);
157:     } else if (reuse == MAT_INITIAL_MATRIX) {
158:       PetscInt m = A->rmap->n,nz = aa->nz;
159:       PetscInt *row,*col;
160:       PetscMalloc2(m+1,&row,nz,&col);
161:       for (i=0; i<m+1; i++) {
162:         row[i] = aa->i[i]+1;
163:       }
164:       for (i=0; i<nz; i++) {
165:         col[i] = aa->j[i]+1;
166:       }
167:       *r    = (INT_TYPE*)row;
168:       *c    = (INT_TYPE*)col;
169:       *nnz  = (INT_TYPE)nz;
170:     }
171:     *free = PETSC_TRUE;
172:   } else {
173:     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen");
174:   }
175:   return(0);
176: }

178: PetscErrorCode MatMKLPardiso_Convert_seqaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,PetscScalar **v)
179: {
180:   Mat_SeqAIJ     *aa = (Mat_SeqAIJ*)A->data;
181:   PetscScalar    *aav;

185:   MatSeqAIJGetArrayRead(A,(const PetscScalar**)&aav);
186:   if (!sym) { /* already in the correct format */
187:     *v    = aav;
188:     *r    = (INT_TYPE*)aa->i;
189:     *c    = (INT_TYPE*)aa->j;
190:     *nnz  = (INT_TYPE)aa->nz;
191:     *free = PETSC_FALSE;
192:   } else if (reuse == MAT_INITIAL_MATRIX) { /* need to get the triangular part */
193:     PetscScalar *vals,*vv;
194:     PetscInt    *row,*col,*jj;
195:     PetscInt    m = A->rmap->n,nz,i;

197:     nz = 0;
198:     for (i=0; i<m; i++) nz += aa->i[i+1] - aa->diag[i];
199:     PetscMalloc2(m+1,&row,nz,&col);
200:     PetscMalloc1(nz,&vals);
201:     jj = col;
202:     vv = vals;

204:     row[0] = 0;
205:     for (i=0; i<m; i++) {
206:       PetscInt    *aj = aa->j + aa->diag[i];
207:       PetscScalar *av = aav + aa->diag[i];
208:       PetscInt    rl  = aa->i[i+1] - aa->diag[i],j;

210:       for (j=0; j<rl; j++) {
211:         *jj = *aj; jj++; aj++;
212:         *vv = *av; vv++; av++;
213:       }
214:       row[i+1] = row[i] + rl;
215:     }
216:     *v    = vals;
217:     *r    = (INT_TYPE*)row;
218:     *c    = (INT_TYPE*)col;
219:     *nnz  = (INT_TYPE)nz;
220:     *free = PETSC_TRUE;
221:   } else {
222:     PetscScalar *vv;
223:     PetscInt    m = A->rmap->n,i;

225:     vv = *v;
226:     for (i=0; i<m; i++) {
227:       PetscScalar *av = aav + aa->diag[i];
228:       PetscInt    rl  = aa->i[i+1] - aa->diag[i],j;
229:       for (j=0; j<rl; j++) {
230:         *vv = *av; vv++; av++;
231:       }
232:     }
233:     *free = PETSC_TRUE;
234:   }
235:   MatSeqAIJRestoreArrayRead(A,(const PetscScalar**)&aav);
236:   return(0);
237: }


240: static PetscErrorCode MatMKLPardisoSolveSchur_Private(Mat F, PetscScalar *B, PetscScalar *X)
241: {
242:   Mat_MKL_PARDISO      *mpardiso = (Mat_MKL_PARDISO*)F->data;
243:   Mat                  S,Xmat,Bmat;
244:   MatFactorSchurStatus schurstatus;
245:   PetscErrorCode       ierr;

248:   MatFactorFactorizeSchurComplement(F);
249:   MatFactorGetSchurComplement(F,&S,&schurstatus);
250:   if (X == B && schurstatus == MAT_FACTOR_SCHUR_INVERTED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"X and B cannot point to the same address");
251:   MatCreateSeqDense(PETSC_COMM_SELF,mpardiso->schur_size,mpardiso->nrhs,B,&Bmat);
252:   MatCreateSeqDense(PETSC_COMM_SELF,mpardiso->schur_size,mpardiso->nrhs,X,&Xmat);
253:   MatSetType(Bmat,((PetscObject)S)->type_name);
254:   MatSetType(Xmat,((PetscObject)S)->type_name);
255: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
256:   MatPinToCPU(Xmat,S->pinnedtocpu);
257:   MatPinToCPU(Bmat,S->pinnedtocpu);
258: #endif

260: #if defined(PETSC_USE_COMPLEX)
261:   if (mpardiso->iparm[12-1] == 1) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Hermitian solve not implemented yet");
262: #endif

264:   switch (schurstatus) {
265:   case MAT_FACTOR_SCHUR_FACTORED:
266:     if (!mpardiso->iparm[12-1]) {
267:       MatMatSolve(S,Bmat,Xmat);
268:     } else { /* transpose solve */
269:       MatMatSolveTranspose(S,Bmat,Xmat);
270:     }
271:     break;
272:   case MAT_FACTOR_SCHUR_INVERTED:
273:     if (!mpardiso->iparm[12-1]) {
274:       MatMatMult(S,Bmat,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Xmat);
275:     } else { /* transpose solve */
276:       MatTransposeMatMult(S,Bmat,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Xmat);
277:     }
278:     break;
279:   default:
280:     SETERRQ1(PetscObjectComm((PetscObject)F),PETSC_ERR_SUP,"Unhandled MatFactorSchurStatus %D",F->schur_status);
281:     break;
282:   }
283:   MatFactorRestoreSchurComplement(F,&S,schurstatus);
284:   MatDestroy(&Bmat);
285:   MatDestroy(&Xmat);
286:   return(0);
287: }

289: PetscErrorCode MatFactorSetSchurIS_MKL_PARDISO(Mat F, IS is)
290: {
291:   Mat_MKL_PARDISO   *mpardiso = (Mat_MKL_PARDISO*)F->data;
292:   const PetscScalar *arr;
293:   const PetscInt    *idxs;
294:   PetscInt          size,i;
295:   PetscMPIInt       csize;
296:   PetscBool         sorted;
297:   PetscErrorCode    ierr;

300:   MPI_Comm_size(PetscObjectComm((PetscObject)F),&csize);
301:   if (csize > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MKL_PARDISO parallel Schur complements not yet supported from PETSc");
302:   ISSorted(is,&sorted);
303:   if (!sorted) {
304:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"IS for MKL_PARDISO Schur complements needs to be sorted");
305:   }
306:   ISGetLocalSize(is,&size);
307:   PetscFree(mpardiso->schur_work);
308:   PetscBLASIntCast(PetscMax(mpardiso->n,2*size),&mpardiso->schur_work_size);
309:   PetscMalloc1(mpardiso->schur_work_size,&mpardiso->schur_work);
310:   MatDestroy(&F->schur);
311:   MatCreateSeqDense(PETSC_COMM_SELF,size,size,NULL,&F->schur);
312:   MatDenseGetArrayRead(F->schur,&arr);
313:   mpardiso->schur      = (PetscScalar*)arr;
314:   mpardiso->schur_size = size;
315:   MatDenseRestoreArrayRead(F->schur,&arr);
316:   if (mpardiso->mtype == 2) {
317:     MatSetOption(F->schur,MAT_SPD,PETSC_TRUE);
318:   }

320:   PetscFree(mpardiso->schur_idxs);
321:   PetscMalloc1(size,&mpardiso->schur_idxs);
322:   PetscArrayzero(mpardiso->perm,mpardiso->n);
323:   ISGetIndices(is,&idxs);
324:   PetscArraycpy(mpardiso->schur_idxs,idxs,size);
325:   for (i=0;i<size;i++) mpardiso->perm[idxs[i]] = 1;
326:   ISRestoreIndices(is,&idxs);
327:   if (size) { /* turn on Schur switch if the set of indices is not empty */
328:     mpardiso->iparm[36-1] = 2;
329:   }
330:   return(0);
331: }

333: PetscErrorCode MatDestroy_MKL_PARDISO(Mat A)
334: {
335:   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
336:   PetscErrorCode  ierr;

339:   if (mat_mkl_pardiso->CleanUp) {
340:     mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;

342:     MKL_PARDISO (mat_mkl_pardiso->pt,
343:       &mat_mkl_pardiso->maxfct,
344:       &mat_mkl_pardiso->mnum,
345:       &mat_mkl_pardiso->mtype,
346:       &mat_mkl_pardiso->phase,
347:       &mat_mkl_pardiso->n,
348:       NULL,
349:       NULL,
350:       NULL,
351:       NULL,
352:       &mat_mkl_pardiso->nrhs,
353:       mat_mkl_pardiso->iparm,
354:       &mat_mkl_pardiso->msglvl,
355:       NULL,
356:       NULL,
357:       &mat_mkl_pardiso->err);
358:   }
359:   PetscFree(mat_mkl_pardiso->perm);
360:   PetscFree(mat_mkl_pardiso->schur_work);
361:   PetscFree(mat_mkl_pardiso->schur_idxs);
362:   if (mat_mkl_pardiso->freeaij) {
363:     PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);
364:     if (mat_mkl_pardiso->iparm[34] == 1) {
365:       PetscFree(mat_mkl_pardiso->a);
366:     }
367:   }
368:   PetscFree(A->data);

370:   /* clear composed functions */
371:   PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverType_C",NULL);
372:   PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);
373:   PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);
374:   return(0);
375: }

377: static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce)
378: {
380:   if (reduce) { /* data given for the whole matrix */
381:     PetscInt i,m=0,p=0;
382:     for (i=0;i<mpardiso->nrhs;i++) {
383:       PetscInt j;
384:       for (j=0;j<mpardiso->schur_size;j++) {
385:         schur[p+j] = whole[m+mpardiso->schur_idxs[j]];
386:       }
387:       m += mpardiso->n;
388:       p += mpardiso->schur_size;
389:     }
390:   } else { /* from Schur to whole */
391:     PetscInt i,m=0,p=0;
392:     for (i=0;i<mpardiso->nrhs;i++) {
393:       PetscInt j;
394:       for (j=0;j<mpardiso->schur_size;j++) {
395:         whole[m+mpardiso->schur_idxs[j]] = schur[p+j];
396:       }
397:       m += mpardiso->n;
398:       p += mpardiso->schur_size;
399:     }
400:   }
401:   return(0);
402: }

404: PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x)
405: {
406:   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
407:   PetscErrorCode    ierr;
408:   PetscScalar       *xarray;
409:   const PetscScalar *barray;

412:   mat_mkl_pardiso->nrhs = 1;
413:   VecGetArray(x,&xarray);
414:   VecGetArrayRead(b,&barray);

416:   if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
417:   else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;

419:   if (barray == xarray) { /* if the two vectors share the same memory */
420:     PetscScalar *work;
421:     if (!mat_mkl_pardiso->schur_work) {
422:       PetscMalloc1(mat_mkl_pardiso->n,&work);
423:     } else {
424:       work = mat_mkl_pardiso->schur_work;
425:     }
426:     mat_mkl_pardiso->iparm[6-1] = 1;
427:     MKL_PARDISO (mat_mkl_pardiso->pt,
428:       &mat_mkl_pardiso->maxfct,
429:       &mat_mkl_pardiso->mnum,
430:       &mat_mkl_pardiso->mtype,
431:       &mat_mkl_pardiso->phase,
432:       &mat_mkl_pardiso->n,
433:       mat_mkl_pardiso->a,
434:       mat_mkl_pardiso->ia,
435:       mat_mkl_pardiso->ja,
436:       NULL,
437:       &mat_mkl_pardiso->nrhs,
438:       mat_mkl_pardiso->iparm,
439:       &mat_mkl_pardiso->msglvl,
440:       (void*)xarray,
441:       (void*)work,
442:       &mat_mkl_pardiso->err);
443:     if (!mat_mkl_pardiso->schur_work) {
444:       PetscFree(work);
445:     }
446:   } else {
447:     mat_mkl_pardiso->iparm[6-1] = 0;
448:     MKL_PARDISO (mat_mkl_pardiso->pt,
449:       &mat_mkl_pardiso->maxfct,
450:       &mat_mkl_pardiso->mnum,
451:       &mat_mkl_pardiso->mtype,
452:       &mat_mkl_pardiso->phase,
453:       &mat_mkl_pardiso->n,
454:       mat_mkl_pardiso->a,
455:       mat_mkl_pardiso->ia,
456:       mat_mkl_pardiso->ja,
457:       mat_mkl_pardiso->perm,
458:       &mat_mkl_pardiso->nrhs,
459:       mat_mkl_pardiso->iparm,
460:       &mat_mkl_pardiso->msglvl,
461:       (void*)barray,
462:       (void*)xarray,
463:       &mat_mkl_pardiso->err);
464:   }
465:   VecRestoreArrayRead(b,&barray);

467:   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);

469:   if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
470:     PetscInt shift = mat_mkl_pardiso->schur_size;

472:     /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
473:     if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;

475:     if (!mat_mkl_pardiso->solve_interior) {
476:       /* solve Schur complement */
477:       MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);
478:       MatMKLPardisoSolveSchur_Private(A,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);
479:       MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);
480:     } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substituted to xarray[schur] will be neglected */
481:       PetscInt i;
482:       for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
483:         xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.;
484:       }
485:     }

487:     /* expansion phase */
488:     mat_mkl_pardiso->iparm[6-1] = 1;
489:     mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
490:     MKL_PARDISO (mat_mkl_pardiso->pt,
491:       &mat_mkl_pardiso->maxfct,
492:       &mat_mkl_pardiso->mnum,
493:       &mat_mkl_pardiso->mtype,
494:       &mat_mkl_pardiso->phase,
495:       &mat_mkl_pardiso->n,
496:       mat_mkl_pardiso->a,
497:       mat_mkl_pardiso->ia,
498:       mat_mkl_pardiso->ja,
499:       mat_mkl_pardiso->perm,
500:       &mat_mkl_pardiso->nrhs,
501:       mat_mkl_pardiso->iparm,
502:       &mat_mkl_pardiso->msglvl,
503:       (void*)xarray,
504:       (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
505:       &mat_mkl_pardiso->err);

507:     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);
508:     mat_mkl_pardiso->iparm[6-1] = 0;
509:   }
510:   VecRestoreArray(x,&xarray);
511:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
512:   return(0);
513: }

515: PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x)
516: {
517:   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
518:   PetscInt        oiparm12;
519:   PetscErrorCode  ierr;

522:   oiparm12 = mat_mkl_pardiso->iparm[12 - 1];
523:   mat_mkl_pardiso->iparm[12 - 1] = 2;
524:   MatSolve_MKL_PARDISO(A,b,x);
525:   mat_mkl_pardiso->iparm[12 - 1] = oiparm12;
526:   return(0);
527: }

529: PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X)
530: {
531:   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->data;
532:   PetscErrorCode    ierr;
533:   const PetscScalar *barray;
534:   PetscScalar       *xarray;
535:   PetscBool         flg;

538:   PetscObjectBaseTypeCompare((PetscObject)B,MATSEQDENSE,&flg);
539:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix");
540:   PetscObjectBaseTypeCompare((PetscObject)X,MATSEQDENSE,&flg);
541:   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix");

543:   MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);

545:   if (mat_mkl_pardiso->nrhs > 0) {
546:     MatDenseGetArrayRead(B,&barray);
547:     MatDenseGetArray(X,&xarray);

549:     if (barray == xarray) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"B and X cannot share the same memory location");
550:     if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
551:     else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;
552:     mat_mkl_pardiso->iparm[6-1] = 0;

554:     MKL_PARDISO (mat_mkl_pardiso->pt,
555:       &mat_mkl_pardiso->maxfct,
556:       &mat_mkl_pardiso->mnum,
557:       &mat_mkl_pardiso->mtype,
558:       &mat_mkl_pardiso->phase,
559:       &mat_mkl_pardiso->n,
560:       mat_mkl_pardiso->a,
561:       mat_mkl_pardiso->ia,
562:       mat_mkl_pardiso->ja,
563:       mat_mkl_pardiso->perm,
564:       &mat_mkl_pardiso->nrhs,
565:       mat_mkl_pardiso->iparm,
566:       &mat_mkl_pardiso->msglvl,
567:       (void*)barray,
568:       (void*)xarray,
569:       &mat_mkl_pardiso->err);
570:     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);

572:     MatDenseRestoreArrayRead(B,&barray);
573:     if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
574:       PetscScalar *o_schur_work = NULL;
575:       PetscInt    shift = mat_mkl_pardiso->schur_size*mat_mkl_pardiso->nrhs,scale;
576:       PetscInt    mem = mat_mkl_pardiso->n*mat_mkl_pardiso->nrhs;

578:       /* allocate extra memory if it is needed */
579:       scale = 1;
580:       if (A->schur_status == MAT_FACTOR_SCHUR_INVERTED) scale = 2;

582:       mem *= scale;
583:       if (mem > mat_mkl_pardiso->schur_work_size) {
584:         o_schur_work = mat_mkl_pardiso->schur_work;
585:         PetscMalloc1(mem,&mat_mkl_pardiso->schur_work);
586:       }

588:       /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
589:       if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;

591:       /* solve Schur complement */
592:       if (!mat_mkl_pardiso->solve_interior) {
593:         MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);
594:         MatMKLPardisoSolveSchur_Private(A,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);
595:         MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);
596:       } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substituted to xarray[schur,n] will be neglected */
597:         PetscInt i,n,m=0;
598:         for (n=0;n<mat_mkl_pardiso->nrhs;n++) {
599:           for (i=0;i<mat_mkl_pardiso->schur_size;i++) {
600:             xarray[mat_mkl_pardiso->schur_idxs[i]+m] = 0.;
601:           }
602:           m += mat_mkl_pardiso->n;
603:         }
604:       }

606:       /* expansion phase */
607:       mat_mkl_pardiso->iparm[6-1] = 1;
608:       mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;
609:       MKL_PARDISO (mat_mkl_pardiso->pt,
610:         &mat_mkl_pardiso->maxfct,
611:         &mat_mkl_pardiso->mnum,
612:         &mat_mkl_pardiso->mtype,
613:         &mat_mkl_pardiso->phase,
614:         &mat_mkl_pardiso->n,
615:         mat_mkl_pardiso->a,
616:         mat_mkl_pardiso->ia,
617:         mat_mkl_pardiso->ja,
618:         mat_mkl_pardiso->perm,
619:         &mat_mkl_pardiso->nrhs,
620:         mat_mkl_pardiso->iparm,
621:         &mat_mkl_pardiso->msglvl,
622:         (void*)xarray,
623:         (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
624:         &mat_mkl_pardiso->err);
625:       if (o_schur_work) { /* restore original schur_work (minimal size) */
626:         PetscFree(mat_mkl_pardiso->schur_work);
627:         mat_mkl_pardiso->schur_work = o_schur_work;
628:       }
629:       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);
630:       mat_mkl_pardiso->iparm[6-1] = 0;
631:     }
632:   }
633:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
634:   return(0);
635: }

637: PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info)
638: {
639:   Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->data;
640:   PetscErrorCode  ierr;

643:   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
644:   (*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);

646:   mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION;
647:   MKL_PARDISO (mat_mkl_pardiso->pt,
648:     &mat_mkl_pardiso->maxfct,
649:     &mat_mkl_pardiso->mnum,
650:     &mat_mkl_pardiso->mtype,
651:     &mat_mkl_pardiso->phase,
652:     &mat_mkl_pardiso->n,
653:     mat_mkl_pardiso->a,
654:     mat_mkl_pardiso->ia,
655:     mat_mkl_pardiso->ja,
656:     mat_mkl_pardiso->perm,
657:     &mat_mkl_pardiso->nrhs,
658:     mat_mkl_pardiso->iparm,
659:     &mat_mkl_pardiso->msglvl,
660:     NULL,
661:     (void*)mat_mkl_pardiso->schur,
662:     &mat_mkl_pardiso->err);
663:   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);

665:   /* report flops */
666:   if (mat_mkl_pardiso->iparm[18] > 0) {
667:     PetscLogFlops(PetscPowRealInt(10.,6)*mat_mkl_pardiso->iparm[18]);
668:   }

670:   if (F->schur) { /* schur output from pardiso is in row major format */
671: #if defined(PETSC_HAVE_CUDA)
672:     F->schur->valid_GPU_matrix = PETSC_OFFLOAD_CPU;
673: #endif
674:     MatFactorRestoreSchurComplement(F,NULL,MAT_FACTOR_SCHUR_UNFACTORED);
675:     MatTranspose(F->schur,MAT_INPLACE_MATRIX,&F->schur);
676:   }
677:   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
678:   mat_mkl_pardiso->CleanUp  = PETSC_TRUE;
679:   return(0);
680: }

682: PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A)
683: {
684:   Mat_MKL_PARDISO     *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
685:   PetscErrorCode      ierr;
686:   PetscInt            icntl,threads=1;
687:   PetscBool           flg;

690:   PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");

692:   PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use within PARDISO","None",threads,&threads,&flg);
693:   if (flg) PetscSetMKL_PARDISOThreads((int)threads);

695:   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);
696:   if (flg) mat_mkl_pardiso->maxfct = icntl;

698:   PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);
699:   if (flg) mat_mkl_pardiso->mnum = icntl;

701:   PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);
702:   if (flg) mat_mkl_pardiso->msglvl = icntl;

704:   PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);
705:   if (flg) {
706:     void *pt[IPARM_SIZE];
707:     mat_mkl_pardiso->mtype = icntl;
708:     MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
709: #if defined(PETSC_USE_REAL_SINGLE)
710:     mat_mkl_pardiso->iparm[27] = 1;
711: #else
712:     mat_mkl_pardiso->iparm[27] = 0;
713: #endif
714:     mat_mkl_pardiso->iparm[34] = 1; /* use 0-based indexing */
715:   }
716:   PetscOptionsInt("-mat_mkl_pardiso_1","Use default values (if 0)","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);

718:   if (flg && icntl != 0) {
719:     PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);
720:     if (flg) mat_mkl_pardiso->iparm[1] = icntl;

722:     PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);
723:     if (flg) mat_mkl_pardiso->iparm[3] = icntl;

725:     PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);
726:     if (flg) mat_mkl_pardiso->iparm[4] = icntl;

728:     PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);
729:     if (flg) mat_mkl_pardiso->iparm[5] = icntl;

731:     PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);
732:     if (flg) mat_mkl_pardiso->iparm[7] = icntl;

734:     PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);
735:     if (flg) mat_mkl_pardiso->iparm[9] = icntl;

737:     PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);
738:     if (flg) mat_mkl_pardiso->iparm[10] = icntl;

740:     PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);
741:     if (flg) mat_mkl_pardiso->iparm[11] = icntl;

743:     PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);
744:     if (flg) mat_mkl_pardiso->iparm[12] = icntl;

746:     PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);
747:     if (flg) mat_mkl_pardiso->iparm[17] = icntl;

749:     PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations (0 to disable)","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);
750:     if (flg) mat_mkl_pardiso->iparm[18] = icntl;

752:     PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);
753:     if (flg) mat_mkl_pardiso->iparm[20] = icntl;

755:     PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);
756:     if (flg) mat_mkl_pardiso->iparm[23] = icntl;

758:     PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);
759:     if (flg) mat_mkl_pardiso->iparm[24] = icntl;

761:     PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);
762:     if (flg) mat_mkl_pardiso->iparm[26] = icntl;

764:     PetscOptionsInt("-mat_mkl_pardiso_31","Partial solve and computing selected components of the solution vectors","None",mat_mkl_pardiso->iparm[30],&icntl,&flg);
765:     if (flg) mat_mkl_pardiso->iparm[30] = icntl;

767:     PetscOptionsInt("-mat_mkl_pardiso_34","Optimal number of threads for conditional numerical reproducibility (CNR) mode","None",mat_mkl_pardiso->iparm[33],&icntl,&flg);
768:     if (flg) mat_mkl_pardiso->iparm[33] = icntl;

770:     PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);
771:     if (flg) mat_mkl_pardiso->iparm[59] = icntl;
772:   }
773:   PetscOptionsEnd();
774:   return(0);
775: }

777: PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso)
778: {
780:   PetscInt       i,bs;
781:   PetscBool      match;

784:   for (i=0; i<IPARM_SIZE; i++) mat_mkl_pardiso->iparm[i] = 0;
785:   for (i=0; i<IPARM_SIZE; i++) mat_mkl_pardiso->pt[i] = 0;
786:   /* Default options for both sym and unsym */
787:   mat_mkl_pardiso->iparm[ 0] =  1; /* Solver default parameters overriden with provided by iparm */
788:   mat_mkl_pardiso->iparm[ 1] =  2; /* Metis reordering */
789:   mat_mkl_pardiso->iparm[ 5] =  0; /* Write solution into x */
790:   mat_mkl_pardiso->iparm[ 7] =  0; /* Max number of iterative refinement steps */
791:   mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
792:   mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
793: #if 0
794:   mat_mkl_pardiso->iparm[23] =  1; /* Parallel factorization control*/
795: #endif
796:   PetscObjectTypeCompareAny((PetscObject)A,&match,MATSEQBAIJ,MATSEQSBAIJ,"");
797:   MatGetBlockSize(A,&bs);
798:   if (!match || bs == 1) {
799:     mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */
800:     mat_mkl_pardiso->n         = A->rmap->N;
801:   } else {
802:     mat_mkl_pardiso->iparm[34] = 0; /* Cluster Sparse Solver use Fortran-style indexing for ia and ja arrays */
803:     mat_mkl_pardiso->iparm[36] = bs;
804:     mat_mkl_pardiso->n         = A->rmap->N/bs;
805:   }
806:   mat_mkl_pardiso->iparm[39] =  0; /* Input: matrix/rhs/solution stored on master */

808:   mat_mkl_pardiso->CleanUp   = PETSC_FALSE;
809:   mat_mkl_pardiso->maxfct    = 1; /* Maximum number of numerical factorizations. */
810:   mat_mkl_pardiso->mnum      = 1; /* Which factorization to use. */
811:   mat_mkl_pardiso->msglvl    = 0; /* 0: do not print 1: Print statistical information in file */
812:   mat_mkl_pardiso->phase     = -1;
813:   mat_mkl_pardiso->err       = 0;

815:   mat_mkl_pardiso->nrhs      = 1;
816:   mat_mkl_pardiso->err       = 0;
817:   mat_mkl_pardiso->phase     = -1;

819:   if (ftype == MAT_FACTOR_LU) {
820:     mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */
821:     mat_mkl_pardiso->iparm[10] =  1; /* Use nonsymmetric permutation and scaling MPS */
822:     mat_mkl_pardiso->iparm[12] =  1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
823:   } else {
824:     mat_mkl_pardiso->iparm[ 9] = 8; /* Perturb the pivot elements with 1E-8 */
825:     mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */
826:     mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
827: #if defined(PETSC_USE_DEBUG)
828:     mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */
829: #endif
830:   }
831:   PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);
832:   for (i=0; i<A->rmap->N; i++) {
833:     mat_mkl_pardiso->perm[i] = 0;
834:   }
835:   mat_mkl_pardiso->schur_size = 0;
836:   return(0);
837: }

839: PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info)
840: {
841:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;
842:   PetscErrorCode  ierr;

845:   mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
846:   PetscSetMKL_PARDISOFromOptions(F,A);

848:   /* throw away any previously computed structure */
849:   if (mat_mkl_pardiso->freeaij) {
850:     PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);
851:     if (mat_mkl_pardiso->iparm[34] == 1) {
852:       PetscFree(mat_mkl_pardiso->a);
853:     }
854:   }
855:   (*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);
856:   if (mat_mkl_pardiso->iparm[34] == 1) mat_mkl_pardiso->n = A->rmap->N;
857:   else mat_mkl_pardiso->n = A->rmap->N/A->rmap->bs;

859:   mat_mkl_pardiso->phase = JOB_ANALYSIS;

861:   /* reset flops counting if requested */
862:   if (mat_mkl_pardiso->iparm[18]) mat_mkl_pardiso->iparm[18] = -1;

864:   MKL_PARDISO (mat_mkl_pardiso->pt,
865:     &mat_mkl_pardiso->maxfct,
866:     &mat_mkl_pardiso->mnum,
867:     &mat_mkl_pardiso->mtype,
868:     &mat_mkl_pardiso->phase,
869:     &mat_mkl_pardiso->n,
870:     mat_mkl_pardiso->a,
871:     mat_mkl_pardiso->ia,
872:     mat_mkl_pardiso->ja,
873:     mat_mkl_pardiso->perm,
874:     &mat_mkl_pardiso->nrhs,
875:     mat_mkl_pardiso->iparm,
876:     &mat_mkl_pardiso->msglvl,
877:     NULL,
878:     NULL,
879:     &mat_mkl_pardiso->err);
880:   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);

882:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;

884:   if (F->factortype == MAT_FACTOR_LU) F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO;
885:   else F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO;

887:   F->ops->solve           = MatSolve_MKL_PARDISO;
888:   F->ops->solvetranspose  = MatSolveTranspose_MKL_PARDISO;
889:   F->ops->matsolve        = MatMatSolve_MKL_PARDISO;
890:   return(0);
891: }

893: PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
894: {

898:   MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);
899:   return(0);
900: }

902: #if !defined(PETSC_USE_COMPLEX)
903: PetscErrorCode MatGetInertia_MKL_PARDISO(Mat F,PetscInt *nneg,PetscInt *nzero,PetscInt *npos)
904: {
905:   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)F->data;

908:   if (nneg) *nneg = mat_mkl_pardiso->iparm[22];
909:   if (npos) *npos = mat_mkl_pardiso->iparm[21];
910:   if (nzero) *nzero = F->rmap->N - (mat_mkl_pardiso->iparm[22] + mat_mkl_pardiso->iparm[21]);
911:   return(0);
912: }
913: #endif

915: PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info)
916: {

920:   MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);
921: #if defined(PETSC_USE_COMPLEX)
922:   F->ops->getinertia = NULL;
923: #else
924:   F->ops->getinertia = MatGetInertia_MKL_PARDISO;
925: #endif
926:   return(0);
927: }

929: PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer)
930: {
931:   PetscErrorCode    ierr;
932:   PetscBool         iascii;
933:   PetscViewerFormat format;
934:   Mat_MKL_PARDISO   *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->data;
935:   PetscInt          i;

938:   if (A->ops->solve != MatSolve_MKL_PARDISO) return(0);

940:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
941:   if (iascii) {
942:     PetscViewerGetFormat(viewer,&format);
943:     if (format == PETSC_VIEWER_ASCII_INFO) {
944:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");
945:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase:             %d \n",mat_mkl_pardiso->phase);
946:       for (i=1; i<=64; i++) {
947:         PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]:     %d \n",i, mat_mkl_pardiso->iparm[i - 1]);
948:       }
949:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct:     %d \n", mat_mkl_pardiso->maxfct);
950:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum:     %d \n", mat_mkl_pardiso->mnum);
951:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype:     %d \n", mat_mkl_pardiso->mtype);
952:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n:     %d \n", mat_mkl_pardiso->n);
953:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs:     %d \n", mat_mkl_pardiso->nrhs);
954:       PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl:     %d \n", mat_mkl_pardiso->msglvl);
955:     }
956:   }
957:   return(0);
958: }


961: PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info)
962: {
963:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)A->data;

966:   info->block_size        = 1.0;
967:   info->nz_used           = mat_mkl_pardiso->iparm[17];
968:   info->nz_allocated      = mat_mkl_pardiso->iparm[17];
969:   info->nz_unneeded       = 0.0;
970:   info->assemblies        = 0.0;
971:   info->mallocs           = 0.0;
972:   info->memory            = 0.0;
973:   info->fill_ratio_given  = 0;
974:   info->fill_ratio_needed = 0;
975:   info->factor_mallocs    = 0;
976:   return(0);
977: }

979: PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival)
980: {
981:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->data;

984:   if (icntl <= 64) {
985:     mat_mkl_pardiso->iparm[icntl - 1] = ival;
986:   } else {
987:     if (icntl == 65) PetscSetMKL_PARDISOThreads(ival);
988:     else if (icntl == 66) mat_mkl_pardiso->maxfct = ival;
989:     else if (icntl == 67) mat_mkl_pardiso->mnum = ival;
990:     else if (icntl == 68) mat_mkl_pardiso->msglvl = ival;
991:     else if (icntl == 69) {
992:       void *pt[IPARM_SIZE];
993:       mat_mkl_pardiso->mtype = ival;
994:       MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
995: #if defined(PETSC_USE_REAL_SINGLE)
996:       mat_mkl_pardiso->iparm[27] = 1;
997: #else
998:       mat_mkl_pardiso->iparm[27] = 0;
999: #endif
1000:       mat_mkl_pardiso->iparm[34] = 1;
1001:     } else if (icntl==70) mat_mkl_pardiso->solve_interior = (PetscBool)!!ival;
1002:   }
1003:   return(0);
1004: }

1006: /*@
1007:   MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters

1009:    Logically Collective on Mat

1011:    Input Parameters:
1012: +  F - the factored matrix obtained by calling MatGetFactor()
1013: .  icntl - index of Mkl_Pardiso parameter
1014: -  ival - value of Mkl_Pardiso parameter

1016:   Options Database:
1017: .   -mat_mkl_pardiso_<icntl> <ival>

1019:    Level: beginner

1021:    References:
1022: .      Mkl_Pardiso Users' Guide

1024: .seealso: MatGetFactor()
1025: @*/
1026: PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival)
1027: {

1031:   PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));
1032:   return(0);
1033: }

1035: /*MC
1036:   MATSOLVERMKL_PARDISO -  A matrix type providing direct solvers (LU) for
1037:   sequential matrices via the external package MKL_PARDISO.

1039:   Works with MATSEQAIJ matrices

1041:   Use -pc_type lu -pc_factor_mat_solver_type mkl_pardiso to use this direct solver

1043:   Options Database Keys:
1044: + -mat_mkl_pardiso_65 - Number of threads to use within MKL_PARDISO
1045: . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
1046: . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase
1047: . -mat_mkl_pardiso_68 - Message level information
1048: . -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
1049: . -mat_mkl_pardiso_1  - Use default values
1050: . -mat_mkl_pardiso_2  - Fill-in reducing ordering for the input matrix
1051: . -mat_mkl_pardiso_4  - Preconditioned CGS/CG
1052: . -mat_mkl_pardiso_5  - User permutation
1053: . -mat_mkl_pardiso_6  - Write solution on x
1054: . -mat_mkl_pardiso_8  - Iterative refinement step
1055: . -mat_mkl_pardiso_10 - Pivoting perturbation
1056: . -mat_mkl_pardiso_11 - Scaling vectors
1057: . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A
1058: . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching
1059: . -mat_mkl_pardiso_18 - Numbers of non-zero elements
1060: . -mat_mkl_pardiso_19 - Report number of floating point operations
1061: . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices
1062: . -mat_mkl_pardiso_24 - Parallel factorization control
1063: . -mat_mkl_pardiso_25 - Parallel forward/backward solve control
1064: . -mat_mkl_pardiso_27 - Matrix checker
1065: . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors
1066: . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
1067: - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode

1069:   Level: beginner

1071:   For more information please check  mkl_pardiso manual

1073: .seealso: PCFactorSetMatSolverType(), MatSolverType

1075: M*/
1076: static PetscErrorCode MatFactorGetSolverType_mkl_pardiso(Mat A, MatSolverType *type)
1077: {
1079:   *type = MATSOLVERMKL_PARDISO;
1080:   return(0);
1081: }

1083: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F)
1084: {
1085:   Mat             B;
1086:   PetscErrorCode  ierr;
1087:   Mat_MKL_PARDISO *mat_mkl_pardiso;
1088:   PetscBool       isSeqAIJ,isSeqBAIJ,isSeqSBAIJ;

1091:   PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);
1092:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);
1093:   PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);
1094:   MatCreate(PetscObjectComm((PetscObject)A),&B);
1095:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
1096:   PetscStrallocpy("mkl_pardiso",&((PetscObject)B)->type_name);
1097:   MatSetUp(B);

1099:   PetscNewLog(B,&mat_mkl_pardiso);
1100:   B->data = mat_mkl_pardiso;

1102:   MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);
1103:   if (ftype == MAT_FACTOR_LU) {
1104:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO;
1105:     B->factortype            = MAT_FACTOR_LU;
1106:     mat_mkl_pardiso->needsym = PETSC_FALSE;
1107:     if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1108:     else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1109:     else if (isSeqSBAIJ) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead");
1110:     else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU with %s format",((PetscObject)A)->type_name);
1111: #if defined(PETSC_USE_COMPLEX)
1112:     mat_mkl_pardiso->mtype = 13;
1113: #else
1114:     mat_mkl_pardiso->mtype = 11;
1115: #endif
1116:   } else {
1117: #if defined(PETSC_USE_COMPLEX)
1118:     SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with complex scalars! Use MAT_FACTOR_LU instead",((PetscObject)A)->type_name);
1119: #endif
1120:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO;
1121:     B->factortype                  = MAT_FACTOR_CHOLESKY;
1122:     if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1123:     else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1124:     else if (isSeqSBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij;
1125:     else SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with %s format",((PetscObject)A)->type_name);

1127:     mat_mkl_pardiso->needsym = PETSC_TRUE;
1128:     if (A->spd_set && A->spd) mat_mkl_pardiso->mtype = 2;
1129:     else                      mat_mkl_pardiso->mtype = -2;
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: }