Actual source code: mkl_pardiso.c

  1: #include <../src/mat/impls/aij/seq/aij.h>
  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:   PetscBLASInt iparm_copy[IPARM_SIZE], mtype_copy;

 45:   PetscCallVoid(PetscBLASIntCast(*mtype, &mtype_copy));
 46:   pardisoinit(pt, &mtype_copy, iparm_copy);
 47:   for (PetscInt 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

 56: #define PetscCallPardiso(f) PetscCallExternalVoid("MKL_PARDISO", f);

 58: /*
 59:    Internal data structure.
 60:  */
 61: typedef struct {
 62:   /* Configuration vector*/
 63:   INT_TYPE iparm[IPARM_SIZE];

 65:   /*
 66:      Internal MKL PARDISO memory location.
 67:      After the first call to MKL PARDISO do not modify pt, as that could cause a serious memory leak.
 68:    */
 69:   void *pt[IPARM_SIZE];

 71:   /* Basic MKL PARDISO info */
 72:   INT_TYPE phase, maxfct, mnum, mtype, n, nrhs, msglvl, err;

 74:   /* Matrix nonzero structure and values */
 75:   void     *a;
 76:   INT_TYPE *ia, *ja;

 78:   /* Number of non-zero elements */
 79:   INT_TYPE nz;

 81:   /* Row permutaton vector */
 82:   INT_TYPE *perm;

 84:   /* Define if matrix preserves sparse structure. */
 85:   MatStructure matstruc;

 87:   PetscBool needsym;
 88:   PetscBool freeaij;

 90:   /* Schur complement */
 91:   PetscScalar *schur;
 92:   PetscInt     schur_size;
 93:   PetscInt    *schur_idxs;
 94:   PetscScalar *schur_work;
 95:   PetscBLASInt schur_work_size;
 96:   PetscBool    solve_interior;

 98:   /* True if MKL PARDISO function have been used. */
 99:   PetscBool CleanUp;

101:   /* Conversion to a format suitable for MKL */
102:   PetscErrorCode (*Convert)(Mat, PetscBool, MatReuse, PetscBool *, INT_TYPE *, INT_TYPE **, INT_TYPE **, PetscScalar **);
103: } Mat_MKL_PARDISO;

105: static PetscErrorCode MatMKLPardiso_Convert_seqsbaij(Mat A, PetscBool sym, MatReuse reuse, PetscBool *free, INT_TYPE *nnz, INT_TYPE **r, INT_TYPE **c, PetscScalar **v)
106: {
107:   Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ *)A->data;
108:   PetscInt      bs = A->rmap->bs, i;

110:   PetscFunctionBegin;
111:   PetscCheck(sym, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "This should not happen");
112:   *v = aa->a;
113:   if (bs == 1) { /* already in the correct format */
114:     /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
115:     *r    = (INT_TYPE *)aa->i;
116:     *c    = (INT_TYPE *)aa->j;
117:     *nnz  = (INT_TYPE)aa->nz;
118:     *free = PETSC_FALSE;
119:   } else if (reuse == MAT_INITIAL_MATRIX) {
120:     PetscInt  m = A->rmap->n, nz = aa->nz;
121:     PetscInt *row, *col;
122:     PetscCall(PetscMalloc2(m + 1, &row, nz, &col));
123:     for (i = 0; i < m + 1; i++) row[i] = aa->i[i] + 1;
124:     for (i = 0; i < nz; i++) col[i] = aa->j[i] + 1;
125:     *r    = (INT_TYPE *)row;
126:     *c    = (INT_TYPE *)col;
127:     *nnz  = (INT_TYPE)nz;
128:     *free = PETSC_TRUE;
129:   }
130:   PetscFunctionReturn(PETSC_SUCCESS);
131: }

133: static PetscErrorCode MatMKLPardiso_Convert_seqbaij(Mat A, PetscBool sym, MatReuse reuse, PetscBool *free, INT_TYPE *nnz, INT_TYPE **r, INT_TYPE **c, PetscScalar **v)
134: {
135:   Mat_SeqBAIJ *aa = (Mat_SeqBAIJ *)A->data;
136:   PetscInt     bs = A->rmap->bs, i;

138:   PetscFunctionBegin;
139:   if (!sym) {
140:     *v = aa->a;
141:     if (bs == 1) { /* already in the correct format */
142:       /* though PetscInt and INT_TYPE are of the same size since they are defined differently the Intel compiler requires a cast */
143:       *r    = (INT_TYPE *)aa->i;
144:       *c    = (INT_TYPE *)aa->j;
145:       *nnz  = (INT_TYPE)aa->nz;
146:       *free = PETSC_FALSE;
147:       PetscFunctionReturn(PETSC_SUCCESS);
148:     } else if (reuse == MAT_INITIAL_MATRIX) {
149:       PetscInt  m = A->rmap->n, nz = aa->nz;
150:       PetscInt *row, *col;
151:       PetscCall(PetscMalloc2(m + 1, &row, nz, &col));
152:       for (i = 0; i < m + 1; i++) row[i] = aa->i[i] + 1;
153:       for (i = 0; i < nz; i++) col[i] = aa->j[i] + 1;
154:       *r   = (INT_TYPE *)row;
155:       *c   = (INT_TYPE *)col;
156:       *nnz = (INT_TYPE)nz;
157:     }
158:     *free = PETSC_TRUE;
159:   } else {
160:     SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "This should not happen");
161:   }
162:   PetscFunctionReturn(PETSC_SUCCESS);
163: }

165: static PetscErrorCode MatMKLPardiso_Convert_seqaij(Mat A, PetscBool sym, MatReuse reuse, PetscBool *free, INT_TYPE *nnz, INT_TYPE **r, INT_TYPE **c, PetscScalar **v)
166: {
167:   Mat_SeqAIJ  *aa = (Mat_SeqAIJ *)A->data;
168:   PetscScalar *aav;

170:   PetscFunctionBegin;
171:   PetscCall(MatSeqAIJGetArrayRead(A, (const PetscScalar **)&aav));
172:   if (!sym) { /* already in the correct format */
173:     *v    = aav;
174:     *r    = (INT_TYPE *)aa->i;
175:     *c    = (INT_TYPE *)aa->j;
176:     *nnz  = (INT_TYPE)aa->nz;
177:     *free = PETSC_FALSE;
178:   } else if (reuse == MAT_INITIAL_MATRIX) { /* need to get the triangular part */
179:     PetscScalar    *vals, *vv;
180:     PetscInt       *row, *col, *jj;
181:     PetscInt        m = A->rmap->n, nz, i;
182:     const PetscInt *adiag;

184:     PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
185:     nz = 0;
186:     for (i = 0; i < m; i++) nz += aa->i[i + 1] - adiag[i];
187:     PetscCall(PetscMalloc2(m + 1, &row, nz, &col));
188:     PetscCall(PetscMalloc1(nz, &vals));
189:     jj = col;
190:     vv = vals;

192:     row[0] = 0;
193:     for (i = 0; i < m; i++) {
194:       PetscInt    *aj = aa->j + adiag[i];
195:       PetscScalar *av = aav + adiag[i];
196:       PetscInt     rl = aa->i[i + 1] - adiag[i], j;

198:       for (j = 0; j < rl; j++) {
199:         *jj = *aj;
200:         jj++;
201:         aj++;
202:         *vv = *av;
203:         vv++;
204:         av++;
205:       }
206:       row[i + 1] = row[i] + rl;
207:     }
208:     *v    = vals;
209:     *r    = (INT_TYPE *)row;
210:     *c    = (INT_TYPE *)col;
211:     *nnz  = (INT_TYPE)nz;
212:     *free = PETSC_TRUE;
213:   } else {
214:     PetscScalar    *vv;
215:     PetscInt        m = A->rmap->n, i;
216:     const PetscInt *adiag;

218:     PetscCall(MatGetDiagonalMarkers_SeqAIJ(A, &adiag, NULL));
219:     vv = *v;
220:     for (i = 0; i < m; i++) {
221:       PetscScalar *av = aav + adiag[i];
222:       PetscInt     rl = aa->i[i + 1] - adiag[i], j;
223:       for (j = 0; j < rl; j++) {
224:         *vv = *av;
225:         vv++;
226:         av++;
227:       }
228:     }
229:     *free = PETSC_TRUE;
230:   }
231:   PetscCall(MatSeqAIJRestoreArrayRead(A, (const PetscScalar **)&aav));
232:   PetscFunctionReturn(PETSC_SUCCESS);
233: }

235: static PetscErrorCode MatMKLPardisoSolveSchur_Private(Mat F, PetscScalar *B, PetscScalar *X)
236: {
237:   Mat_MKL_PARDISO     *mpardiso = (Mat_MKL_PARDISO *)F->data;
238:   Mat                  S, Xmat, Bmat;
239:   MatFactorSchurStatus schurstatus;

241:   PetscFunctionBegin;
242:   PetscCall(MatFactorGetSchurComplement(F, &S, &schurstatus));
243:   PetscCheck(X != B || schurstatus != MAT_FACTOR_SCHUR_INVERTED, PETSC_COMM_SELF, PETSC_ERR_SUP, "X and B cannot point to the same address");
244:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mpardiso->schur_size, mpardiso->nrhs, B, &Bmat));
245:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, mpardiso->schur_size, mpardiso->nrhs, X, &Xmat));
246:   PetscCall(MatSetType(Bmat, ((PetscObject)S)->type_name));
247:   PetscCall(MatSetType(Xmat, ((PetscObject)S)->type_name));
248: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
249:   PetscCall(MatBindToCPU(Xmat, S->boundtocpu));
250:   PetscCall(MatBindToCPU(Bmat, S->boundtocpu));
251: #endif

253: #if defined(PETSC_USE_COMPLEX)
254:   PetscCheck(mpardiso->iparm[12 - 1] != 1, PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Hermitian solve not implemented yet");
255: #endif

257:   switch (schurstatus) {
258:   case MAT_FACTOR_SCHUR_FACTORED:
259:     if (!mpardiso->iparm[12 - 1]) {
260:       PetscCall(MatMatSolve(S, Bmat, Xmat));
261:     } else { /* transpose solve */
262:       PetscCall(MatMatSolveTranspose(S, Bmat, Xmat));
263:     }
264:     break;
265:   case MAT_FACTOR_SCHUR_INVERTED:
266:     PetscCall(MatProductCreateWithMat(S, Bmat, NULL, Xmat));
267:     if (!mpardiso->iparm[12 - 1]) {
268:       PetscCall(MatProductSetType(Xmat, MATPRODUCT_AB));
269:     } else { /* transpose solve */
270:       PetscCall(MatProductSetType(Xmat, MATPRODUCT_AtB));
271:     }
272:     PetscCall(MatProductSetFromOptions(Xmat));
273:     PetscCall(MatProductSymbolic(Xmat));
274:     PetscCall(MatProductNumeric(Xmat));
275:     PetscCall(MatProductClear(Xmat));
276:     break;
277:   default:
278:     SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", (int)F->schur_status);
279:     break;
280:   }
281:   PetscCall(MatFactorRestoreSchurComplement(F, &S, schurstatus));
282:   PetscCall(MatDestroy(&Bmat));
283:   PetscCall(MatDestroy(&Xmat));
284:   PetscFunctionReturn(PETSC_SUCCESS);
285: }

287: static 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;

296:   PetscFunctionBegin;
297:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)F), &csize));
298:   PetscCheck(csize <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "MKL PARDISO parallel Schur complements not yet supported from PETSc");
299:   PetscCall(ISSorted(is, &sorted));
300:   PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_SUP, "IS for MKL PARDISO Schur complements needs to be sorted");
301:   PetscCall(ISGetLocalSize(is, &size));
302:   PetscCall(PetscFree(mpardiso->schur_work));
303:   PetscCall(PetscBLASIntCast(PetscMax(mpardiso->n, 2 * size), &mpardiso->schur_work_size));
304:   PetscCall(PetscMalloc1(mpardiso->schur_work_size, &mpardiso->schur_work));
305:   PetscCall(MatDestroy(&F->schur));
306:   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, size, size, NULL, &F->schur));
307:   PetscCall(MatDenseGetArrayRead(F->schur, &arr));
308:   mpardiso->schur      = (PetscScalar *)arr;
309:   mpardiso->schur_size = size;
310:   PetscCall(MatDenseRestoreArrayRead(F->schur, &arr));
311:   if (mpardiso->mtype == 2) PetscCall(MatSetOption(F->schur, MAT_SPD, PETSC_TRUE));

313:   PetscCall(PetscFree(mpardiso->schur_idxs));
314:   PetscCall(PetscMalloc1(size, &mpardiso->schur_idxs));
315:   PetscCall(PetscArrayzero(mpardiso->perm, mpardiso->n));
316:   PetscCall(ISGetIndices(is, &idxs));
317:   PetscCall(PetscArraycpy(mpardiso->schur_idxs, idxs, size));
318:   for (i = 0; i < size; i++) mpardiso->perm[idxs[i]] = 1;
319:   PetscCall(ISRestoreIndices(is, &idxs));
320:   if (size) { /* turn on Schur switch if the set of indices is not empty */
321:     mpardiso->iparm[36 - 1] = 2;
322:   }
323:   PetscFunctionReturn(PETSC_SUCCESS);
324: }

326: static PetscErrorCode MatDestroy_MKL_PARDISO(Mat A)
327: {
328:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;

330:   PetscFunctionBegin;
331:   if (mat_mkl_pardiso->CleanUp) {
332:     mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY;

334:     PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, NULL, NULL, NULL, NULL, &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm,
335:                                  &mat_mkl_pardiso->msglvl, NULL, NULL, &mat_mkl_pardiso->err));
336:   }
337:   PetscCall(PetscFree(mat_mkl_pardiso->perm));
338:   PetscCall(PetscFree(mat_mkl_pardiso->schur_work));
339:   PetscCall(PetscFree(mat_mkl_pardiso->schur_idxs));
340:   if (mat_mkl_pardiso->freeaij) {
341:     PetscCall(PetscFree2(mat_mkl_pardiso->ia, mat_mkl_pardiso->ja));
342:     if (mat_mkl_pardiso->iparm[34] == 1) PetscCall(PetscFree(mat_mkl_pardiso->a));
343:   }
344:   PetscCall(PetscFree(A->data));

346:   /* clear composed functions */
347:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
348:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorSetSchurIS_C", NULL));
349:   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMkl_PardisoSetCntl_C", NULL));
350:   PetscFunctionReturn(PETSC_SUCCESS);
351: }

353: static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce)
354: {
355:   PetscFunctionBegin;
356:   if (reduce) { /* data given for the whole matrix */
357:     PetscInt i, m = 0, p = 0;
358:     for (i = 0; i < mpardiso->nrhs; i++) {
359:       for (PetscInt j = 0; j < mpardiso->schur_size; j++) schur[p + j] = whole[m + mpardiso->schur_idxs[j]];
360:       m += mpardiso->n;
361:       p += mpardiso->schur_size;
362:     }
363:   } else { /* from Schur to whole */
364:     PetscInt i, m = 0, p = 0;
365:     for (i = 0; i < mpardiso->nrhs; i++) {
366:       for (PetscInt j = 0; j < mpardiso->schur_size; j++) whole[m + mpardiso->schur_idxs[j]] = schur[p + j];
367:       m += mpardiso->n;
368:       p += mpardiso->schur_size;
369:     }
370:   }
371:   PetscFunctionReturn(PETSC_SUCCESS);
372: }

374: static PetscErrorCode MatSolve_MKL_PARDISO(Mat A, Vec b, Vec x)
375: {
376:   Mat_MKL_PARDISO   *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;
377:   PetscScalar       *xarray;
378:   const PetscScalar *barray;

380:   PetscFunctionBegin;
381:   mat_mkl_pardiso->nrhs = 1;
382:   PetscCall(VecGetArrayWrite(x, &xarray));
383:   PetscCall(VecGetArrayRead(b, &barray));

385:   if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
386:   else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;

388:   if (barray == xarray) { /* if the two vectors share the same memory */
389:     PetscScalar *work;
390:     if (!mat_mkl_pardiso->schur_work) {
391:       PetscCall(PetscMalloc1(mat_mkl_pardiso->n, &work));
392:     } else {
393:       work = mat_mkl_pardiso->schur_work;
394:     }
395:     mat_mkl_pardiso->iparm[6 - 1] = 1;
396:     PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja, NULL,
397:                                  &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)xarray, (void *)work, &mat_mkl_pardiso->err));
398:     if (!mat_mkl_pardiso->schur_work) PetscCall(PetscFree(work));
399:   } else {
400:     mat_mkl_pardiso->iparm[6 - 1] = 0;
401:     PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja,
402:                                  mat_mkl_pardiso->perm, &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_pardiso->err));
403:   }
404:   PetscCall(VecRestoreArrayRead(b, &barray));

406:   PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);

408:   if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
409:     if (!mat_mkl_pardiso->solve_interior) {
410:       PetscInt shift = mat_mkl_pardiso->schur_size;

412:       PetscCall(MatFactorFactorizeSchurComplement(A));
413:       /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
414:       if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;

416:       /* solve Schur complement */
417:       PetscCall(MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso, xarray, mat_mkl_pardiso->schur_work, PETSC_TRUE));
418:       PetscCall(MatMKLPardisoSolveSchur_Private(A, mat_mkl_pardiso->schur_work, mat_mkl_pardiso->schur_work + shift));
419:       PetscCall(MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso, xarray, mat_mkl_pardiso->schur_work + shift, PETSC_FALSE));
420:     } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substituted to xarray[schur] will be neglected */
421:       PetscInt i;
422:       for (i = 0; i < mat_mkl_pardiso->schur_size; i++) xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.;
423:     }

425:     /* expansion phase */
426:     mat_mkl_pardiso->iparm[6 - 1] = 1;
427:     mat_mkl_pardiso->phase        = JOB_SOLVE_BACKWARD_SUBSTITUTION;
428:     PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja,
429:                                  mat_mkl_pardiso->perm, &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)xarray, (void *)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
430:                                  &mat_mkl_pardiso->err));
431:     PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);
432:     mat_mkl_pardiso->iparm[6 - 1] = 0;
433:   }
434:   PetscCall(VecRestoreArrayWrite(x, &xarray));
435:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
436:   PetscFunctionReturn(PETSC_SUCCESS);
437: }

439: static PetscErrorCode MatForwardSolve_MKL_PARDISO(Mat A, Vec b, Vec x)
440: {
441:   Mat_MKL_PARDISO   *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;
442:   PetscScalar       *xarray;
443:   const PetscScalar *barray;

445:   PetscFunctionBegin;
446:   PetscCheck(!mat_mkl_pardiso->schur, PETSC_COMM_SELF, PETSC_ERR_SUP, "Forward substitution not supported with Schur complement");

448:   mat_mkl_pardiso->nrhs = 1;
449:   PetscCall(VecGetArrayWrite(x, &xarray));
450:   PetscCall(VecGetArrayRead(b, &barray));

452:   mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;

454:   PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja, mat_mkl_pardiso->perm,
455:                                &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_pardiso->err));
456:   PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);

458:   PetscCall(VecRestoreArrayRead(b, &barray));
459:   PetscCall(VecRestoreArrayWrite(x, &xarray));
460:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
461:   PetscFunctionReturn(PETSC_SUCCESS);
462: }

464: static PetscErrorCode MatBackwardSolve_MKL_PARDISO(Mat A, Vec b, Vec x)
465: {
466:   Mat_MKL_PARDISO   *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;
467:   PetscScalar       *xarray;
468:   const PetscScalar *barray;

470:   PetscFunctionBegin;
471:   PetscCheck(!mat_mkl_pardiso->schur, PETSC_COMM_SELF, PETSC_ERR_SUP, "Backward substitution not supported with Schur complement");

473:   mat_mkl_pardiso->nrhs = 1;
474:   PetscCall(VecGetArrayWrite(x, &xarray));
475:   PetscCall(VecGetArrayRead(b, &barray));

477:   mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION;

479:   PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja, mat_mkl_pardiso->perm,
480:                                &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_pardiso->err));
481:   PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);

483:   PetscCall(VecRestoreArrayRead(b, &barray));
484:   PetscCall(VecRestoreArrayWrite(x, &xarray));
485:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }

489: static PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A, Vec b, Vec x)
490: {
491:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;
492:   PetscInt         oiparm12;

494:   PetscFunctionBegin;
495:   oiparm12                       = mat_mkl_pardiso->iparm[12 - 1];
496:   mat_mkl_pardiso->iparm[12 - 1] = 2;
497:   PetscCall(MatSolve_MKL_PARDISO(A, b, x));
498:   mat_mkl_pardiso->iparm[12 - 1] = oiparm12;
499:   PetscFunctionReturn(PETSC_SUCCESS);
500: }

502: static PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A, Mat B, Mat X)
503: {
504:   Mat_MKL_PARDISO   *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;
505:   const PetscScalar *barray;
506:   PetscScalar       *xarray;
507:   PetscBool          flg;

509:   PetscFunctionBegin;
510:   PetscCall(PetscObjectBaseTypeCompare((PetscObject)B, MATSEQDENSE, &flg));
511:   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix B must be MATSEQDENSE matrix");
512:   if (X != B) {
513:     PetscCall(PetscObjectBaseTypeCompare((PetscObject)X, MATSEQDENSE, &flg));
514:     PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Matrix X must be MATSEQDENSE matrix");
515:   }

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

519:   if (mat_mkl_pardiso->nrhs > 0) {
520:     PetscCall(MatDenseGetArrayRead(B, &barray));
521:     PetscCall(MatDenseGetArrayWrite(X, &xarray));

523:     PetscCheck(barray != xarray, PETSC_COMM_SELF, PETSC_ERR_SUP, "B and X cannot share the same memory location");
524:     if (!mat_mkl_pardiso->schur) mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT;
525:     else mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION;

527:     PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja,
528:                                  mat_mkl_pardiso->perm, &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_pardiso->err));
529:     PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);

531:     PetscCall(MatDenseRestoreArrayRead(B, &barray));
532:     if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */
533:       PetscScalar *o_schur_work = NULL;

535:       /* solve Schur complement */
536:       if (!mat_mkl_pardiso->solve_interior) {
537:         PetscInt shift = mat_mkl_pardiso->schur_size * mat_mkl_pardiso->nrhs, scale;
538:         PetscInt mem   = mat_mkl_pardiso->n * mat_mkl_pardiso->nrhs;

540:         PetscCall(MatFactorFactorizeSchurComplement(A));
541:         /* allocate extra memory if it is needed */
542:         scale = 1;
543:         if (A->schur_status == MAT_FACTOR_SCHUR_INVERTED) scale = 2;
544:         mem *= scale;
545:         if (mem > mat_mkl_pardiso->schur_work_size) {
546:           o_schur_work = mat_mkl_pardiso->schur_work;
547:           PetscCall(PetscMalloc1(mem, &mat_mkl_pardiso->schur_work));
548:         }
549:         /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */
550:         if (A->schur_status != MAT_FACTOR_SCHUR_INVERTED) shift = 0;
551:         PetscCall(MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso, xarray, mat_mkl_pardiso->schur_work, PETSC_TRUE));
552:         PetscCall(MatMKLPardisoSolveSchur_Private(A, mat_mkl_pardiso->schur_work, mat_mkl_pardiso->schur_work + shift));
553:         PetscCall(MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso, xarray, mat_mkl_pardiso->schur_work + shift, PETSC_FALSE));
554:       } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substituted to xarray[schur,n] will be neglected */
555:         PetscInt i, n, m = 0;
556:         for (n = 0; n < mat_mkl_pardiso->nrhs; n++) {
557:           for (i = 0; i < mat_mkl_pardiso->schur_size; i++) xarray[mat_mkl_pardiso->schur_idxs[i] + m] = 0.;
558:           m += mat_mkl_pardiso->n;
559:         }
560:       }

562:       /* expansion phase */
563:       mat_mkl_pardiso->iparm[6 - 1] = 1;
564:       mat_mkl_pardiso->phase        = JOB_SOLVE_BACKWARD_SUBSTITUTION;
565:       PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja,
566:                                    mat_mkl_pardiso->perm, &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, (void *)xarray, (void *)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */
567:                                    &mat_mkl_pardiso->err));
568:       if (o_schur_work) { /* restore original Schur_work (minimal size) */
569:         PetscCall(PetscFree(mat_mkl_pardiso->schur_work));
570:         mat_mkl_pardiso->schur_work = o_schur_work;
571:       }
572:       PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);
573:       mat_mkl_pardiso->iparm[6 - 1] = 0;
574:     }
575:     PetscCall(MatDenseRestoreArrayWrite(X, &xarray));
576:   }
577:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;
578:   PetscFunctionReturn(PETSC_SUCCESS);
579: }

581: static PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F, Mat A, const MatFactorInfo *info)
582: {
583:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)F->data;

585:   PetscFunctionBegin;
586:   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
587:   PetscCall((*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));

589:   mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION;
590:   PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja, mat_mkl_pardiso->perm,
591:                                &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, NULL, (void *)mat_mkl_pardiso->schur, &mat_mkl_pardiso->err));
592:   PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);

594:   /* report flops */
595:   if (mat_mkl_pardiso->iparm[18] > 0) PetscCall(PetscLogFlops(PetscPowRealInt(10., 6) * mat_mkl_pardiso->iparm[18]));

597:   if (F->schur) { /* schur output from pardiso is in row major format */
598: #if defined(PETSC_HAVE_CUDA)
599:     F->schur->offloadmask = PETSC_OFFLOAD_CPU;
600: #endif
601:     PetscCall(MatFactorRestoreSchurComplement(F, NULL, MAT_FACTOR_SCHUR_UNFACTORED));
602:     PetscCall(MatTranspose(F->schur, MAT_INPLACE_MATRIX, &F->schur));
603:   }
604:   mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN;
605:   mat_mkl_pardiso->CleanUp  = PETSC_TRUE;
606:   PetscFunctionReturn(PETSC_SUCCESS);
607: }

609: static PetscErrorCode MatSetFromOptions_MKL_PARDISO(Mat F, Mat A)
610: {
611:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)F->data;
612:   PetscInt         icntl, bs, threads = 1;
613:   PetscBool        flg;

615:   PetscFunctionBegin;
616:   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "MKL_PARDISO Options", "Mat");

618:   PetscCall(PetscOptionsInt("-mat_mkl_pardiso_65", "Suggested number of threads to use within MKL PARDISO", "None", threads, &threads, &flg));
619:   if (flg) PetscSetMKL_PARDISOThreads((int)threads);

621:   PetscCall(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));
622:   if (flg) mat_mkl_pardiso->maxfct = icntl;

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

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

630:   PetscCall(PetscOptionsInt("-mat_mkl_pardiso_69", "Defines the matrix type", "None", mat_mkl_pardiso->mtype, &icntl, &flg));
631:   if (flg) {
632:     void *pt[IPARM_SIZE];
633:     mat_mkl_pardiso->mtype = icntl;
634:     icntl                  = mat_mkl_pardiso->iparm[34];
635:     bs                     = mat_mkl_pardiso->iparm[36];
636:     MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
637: #if defined(PETSC_USE_REAL_SINGLE)
638:     mat_mkl_pardiso->iparm[27] = 1;
639: #else
640:     mat_mkl_pardiso->iparm[27] = 0;
641: #endif
642:     mat_mkl_pardiso->iparm[34] = icntl;
643:     mat_mkl_pardiso->iparm[36] = bs;
644:   }

646:   PetscCall(PetscOptionsInt("-mat_mkl_pardiso_1", "Use default values (if 0)", "None", mat_mkl_pardiso->iparm[0], &icntl, &flg));
647:   if (flg) mat_mkl_pardiso->iparm[0] = icntl;

649:   PetscCall(PetscOptionsInt("-mat_mkl_pardiso_2", "Fill-in reducing ordering for the input matrix", "None", mat_mkl_pardiso->iparm[1], &icntl, &flg));
650:   if (flg) mat_mkl_pardiso->iparm[1] = icntl;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

700:   PetscCall(PetscOptionsInt("-mat_mkl_pardiso_60", "Intel MKL PARDISO mode", "None", mat_mkl_pardiso->iparm[59], &icntl, &flg));
701:   if (flg) mat_mkl_pardiso->iparm[59] = icntl;
702:   PetscOptionsEnd();
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: static PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso)
707: {
708:   PetscInt  bs;
709:   PetscBool match;

711:   PetscFunctionBegin;
712:   for (PetscInt i = 0; i < IPARM_SIZE; i++) mat_mkl_pardiso->iparm[i] = 0;
713:   for (PetscInt i = 0; i < IPARM_SIZE; i++) mat_mkl_pardiso->pt[i] = 0;
714: #if defined(PETSC_USE_REAL_SINGLE)
715:   mat_mkl_pardiso->iparm[27] = 1;
716: #else
717:   mat_mkl_pardiso->iparm[27] = 0;
718: #endif
719:   /* Default options for both sym and unsym */
720:   mat_mkl_pardiso->iparm[0]  = 1;  /* Solver default parameters overridden with provided by iparm */
721:   mat_mkl_pardiso->iparm[1]  = 2;  /* Metis reordering */
722:   mat_mkl_pardiso->iparm[5]  = 0;  /* Write solution into x */
723:   mat_mkl_pardiso->iparm[7]  = 0;  /* Max number of iterative refinement steps */
724:   mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */
725:   mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */
726: #if 0
727:   mat_mkl_pardiso->iparm[23] =  1; /* Parallel factorization control*/
728: #endif
729:   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &match, MATSEQBAIJ, MATSEQSBAIJ, ""));
730:   PetscCall(MatGetBlockSize(A, &bs));
731:   if (!match || bs == 1) {
732:     mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */
733:     mat_mkl_pardiso->n         = A->rmap->N;
734:   } else {
735:     mat_mkl_pardiso->iparm[34] = 0; /* Cluster Sparse Solver use Fortran-style indexing for ia and ja arrays */
736:     mat_mkl_pardiso->iparm[36] = bs;
737:     mat_mkl_pardiso->n         = A->rmap->N / bs;
738:   }
739:   mat_mkl_pardiso->iparm[39] = 0; /* Input: matrix/rhs/solution stored on rank-0 */

741:   mat_mkl_pardiso->CleanUp = PETSC_FALSE;
742:   mat_mkl_pardiso->maxfct  = 1; /* Maximum number of numerical factorizations. */
743:   mat_mkl_pardiso->mnum    = 1; /* Which factorization to use. */
744:   mat_mkl_pardiso->msglvl  = 0; /* 0: do not print 1: Print statistical information in file */
745:   mat_mkl_pardiso->phase   = -1;
746:   mat_mkl_pardiso->err     = 0;

748:   mat_mkl_pardiso->nrhs  = 1;
749:   mat_mkl_pardiso->err   = 0;
750:   mat_mkl_pardiso->phase = -1;

752:   if (ftype == MAT_FACTOR_LU) {
753:     mat_mkl_pardiso->iparm[9]  = 13; /* Perturb the pivot elements with 1E-13 */
754:     mat_mkl_pardiso->iparm[10] = 1;  /* Use nonsymmetric permutation and scaling MPS */
755:     mat_mkl_pardiso->iparm[12] = 1;  /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
756:   } else {
757:     mat_mkl_pardiso->iparm[9]  = 8; /* Perturb the pivot elements with 1E-8 */
758:     mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */
759:     mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */
760: #if defined(PETSC_USE_DEBUG)
761:     mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */
762: #endif
763:   }
764:   PetscCall(PetscCalloc1(A->rmap->N * sizeof(INT_TYPE), &mat_mkl_pardiso->perm));
765:   mat_mkl_pardiso->schur_size = 0;
766:   PetscFunctionReturn(PETSC_SUCCESS);
767: }

769: static PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F, Mat A, const MatFactorInfo *info)
770: {
771:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)F->data;

773:   PetscFunctionBegin;
774:   mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN;
775:   PetscCall(MatSetFromOptions_MKL_PARDISO(F, A));
776:   /* throw away any previously computed structure */
777:   if (mat_mkl_pardiso->freeaij) {
778:     PetscCall(PetscFree2(mat_mkl_pardiso->ia, mat_mkl_pardiso->ja));
779:     if (mat_mkl_pardiso->iparm[34] == 1) PetscCall(PetscFree(mat_mkl_pardiso->a));
780:   }
781:   PetscCall((*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));
782:   if (mat_mkl_pardiso->iparm[34] == 1) mat_mkl_pardiso->n = A->rmap->N;
783:   else mat_mkl_pardiso->n = A->rmap->N / A->rmap->bs;

785:   mat_mkl_pardiso->phase = JOB_ANALYSIS;

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

790:   PetscCallPardiso(MKL_PARDISO(mat_mkl_pardiso->pt, &mat_mkl_pardiso->maxfct, &mat_mkl_pardiso->mnum, &mat_mkl_pardiso->mtype, &mat_mkl_pardiso->phase, &mat_mkl_pardiso->n, mat_mkl_pardiso->a, mat_mkl_pardiso->ia, mat_mkl_pardiso->ja, mat_mkl_pardiso->perm,
791:                                &mat_mkl_pardiso->nrhs, mat_mkl_pardiso->iparm, &mat_mkl_pardiso->msglvl, NULL, NULL, &mat_mkl_pardiso->err));
792:   PetscCheck(mat_mkl_pardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL PARDISO: err=%" PetscInt_FMT ". Please check manual", (PetscInt)mat_mkl_pardiso->err);

794:   mat_mkl_pardiso->CleanUp = PETSC_TRUE;

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

799:   F->ops->solve          = MatSolve_MKL_PARDISO;
800:   F->ops->solvetranspose = MatSolveTranspose_MKL_PARDISO;
801:   F->ops->matsolve       = MatMatSolve_MKL_PARDISO;
802:   if (F->factortype == MAT_FACTOR_LU || (!PetscDefined(USE_COMPLEX) && F->factortype == MAT_FACTOR_CHOLESKY && A->spd == PETSC_BOOL3_TRUE)) {
803:     F->ops->backwardsolve = MatBackwardSolve_MKL_PARDISO;
804:     F->ops->forwardsolve  = MatForwardSolve_MKL_PARDISO;
805:   }
806:   PetscFunctionReturn(PETSC_SUCCESS);
807: }

809: static PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
810: {
811:   PetscFunctionBegin;
812:   PetscCall(MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info));
813:   PetscFunctionReturn(PETSC_SUCCESS);
814: }

816: #if !defined(PETSC_USE_COMPLEX)
817: static PetscErrorCode MatGetInertia_MKL_PARDISO(Mat F, PetscInt *nneg, PetscInt *nzero, PetscInt *npos)
818: {
819:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)F->data;

821:   PetscFunctionBegin;
822:   if (nneg) *nneg = mat_mkl_pardiso->iparm[22];
823:   if (npos) *npos = mat_mkl_pardiso->iparm[21];
824:   if (nzero) *nzero = F->rmap->N - (mat_mkl_pardiso->iparm[22] + mat_mkl_pardiso->iparm[21]);
825:   PetscFunctionReturn(PETSC_SUCCESS);
826: }
827: #endif

829: static PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F, Mat A, IS r, const MatFactorInfo *info)
830: {
831:   PetscFunctionBegin;
832:   PetscCall(MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info));
833:   F->ops->getinertia = NULL;
834: #if !defined(PETSC_USE_COMPLEX)
835:   F->ops->getinertia = MatGetInertia_MKL_PARDISO;
836: #endif
837:   PetscFunctionReturn(PETSC_SUCCESS);
838: }

840: static PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer)
841: {
842:   PetscBool         isascii;
843:   PetscViewerFormat format;
844:   Mat_MKL_PARDISO  *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;

846:   PetscFunctionBegin;
847:   if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(PETSC_SUCCESS);

849:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
850:   if (isascii) {
851:     PetscCall(PetscViewerGetFormat(viewer, &format));
852:     if (format == PETSC_VIEWER_ASCII_INFO) {
853:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO run parameters:\n"));
854:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO phase:             %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->phase));
855:       for (PetscInt i = 1; i <= 64; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO iparm[%" PetscInt_FMT "]:     %" PetscInt_FMT "\n", i, (PetscInt)mat_mkl_pardiso->iparm[i - 1]));
856:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO maxfct:     %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->maxfct));
857:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO mnum:     %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->mnum));
858:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO mtype:     %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->mtype));
859:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO n:     %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->n));
860:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO nrhs:     %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->nrhs));
861:       PetscCall(PetscViewerASCIIPrintf(viewer, "MKL PARDISO msglvl:     %" PetscInt_FMT "\n", (PetscInt)mat_mkl_pardiso->msglvl));
862:     }
863:   }
864:   PetscFunctionReturn(PETSC_SUCCESS);
865: }

867: static PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info)
868: {
869:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)A->data;

871:   PetscFunctionBegin;
872:   info->block_size        = 1.0;
873:   info->nz_used           = mat_mkl_pardiso->iparm[17];
874:   info->nz_allocated      = mat_mkl_pardiso->iparm[17];
875:   info->nz_unneeded       = 0.0;
876:   info->assemblies        = 0.0;
877:   info->mallocs           = 0.0;
878:   info->memory            = 0.0;
879:   info->fill_ratio_given  = 0;
880:   info->fill_ratio_needed = 0;
881:   info->factor_mallocs    = 0;
882:   PetscFunctionReturn(PETSC_SUCCESS);
883: }

885: static PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F, PetscInt icntl, PetscInt ival)
886: {
887:   PetscInt         backup, bs;
888:   Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO *)F->data;

890:   PetscFunctionBegin;
891:   if (icntl <= 64) {
892:     mat_mkl_pardiso->iparm[icntl - 1] = ival;
893:   } else {
894:     if (icntl == 65) PetscSetMKL_PARDISOThreads((int)ival);
895:     else if (icntl == 66) mat_mkl_pardiso->maxfct = ival;
896:     else if (icntl == 67) mat_mkl_pardiso->mnum = ival;
897:     else if (icntl == 68) mat_mkl_pardiso->msglvl = ival;
898:     else if (icntl == 69) {
899:       void *pt[IPARM_SIZE];
900:       backup                 = mat_mkl_pardiso->iparm[34];
901:       bs                     = mat_mkl_pardiso->iparm[36];
902:       mat_mkl_pardiso->mtype = ival;
903:       MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm);
904: #if defined(PETSC_USE_REAL_SINGLE)
905:       mat_mkl_pardiso->iparm[27] = 1;
906: #else
907:       mat_mkl_pardiso->iparm[27] = 0;
908: #endif
909:       mat_mkl_pardiso->iparm[34] = backup;
910:       mat_mkl_pardiso->iparm[36] = bs;
911:     } else if (icntl == 70) mat_mkl_pardiso->solve_interior = (PetscBool)!!ival;
912:   }
913:   PetscFunctionReturn(PETSC_SUCCESS);
914: }

916: /*@
917:   MatMkl_PardisoSetCntl - Set MKL PARDISO <https://www.intel.com/content/www/us/en/docs/onemkl/developer-reference-c/2023-2/onemkl-pardiso-parallel-direct-sparse-solver-iface.html> parameters

919:   Logically Collective

921:   Input Parameters:
922: + F     - the factored matrix obtained by calling `MatGetFactor()`
923: . icntl - index of MKL PARDISO parameter
924: - ival  - value of MKL PARDISO parameter

926:   Options Database Key:
927: . -mat_mkl_pardiso_ICNTL ival - change the option numbered ICNTL to the value `ival`

929:   Level: beginner

931: .seealso: [](ch_matrices), `Mat`, `MATSOLVERMKL_PARDISO`, `MatGetFactor()`
932: @*/
933: PetscErrorCode MatMkl_PardisoSetCntl(Mat F, PetscInt icntl, PetscInt ival)
934: {
935:   PetscFunctionBegin;
936:   PetscTryMethod(F, "MatMkl_PardisoSetCntl_C", (Mat, PetscInt, PetscInt), (F, icntl, ival));
937:   PetscFunctionReturn(PETSC_SUCCESS);
938: }

940: /*MC
941:   MATSOLVERMKL_PARDISO -  A matrix type providing direct solvers, LU, for
942:   `MATSEQAIJ` matrices via the external package MKL PARDISO
943:   <https://www.intel.com/content/www/us/en/docs/onemkl/developer-reference-c/2024-0/onemkl-pardiso-parallel-direct-sparse-solver-iface.html>.

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

947:   Options Database Keys:
948: + -mat_mkl_pardiso_65 - Suggested number of threads to use within MKL PARDISO
949: . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time
950: . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase
951: . -mat_mkl_pardiso_68 - Message level information, use 1 to get detailed information on the solver options
952: . -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
953: . -mat_mkl_pardiso_1  - Use default values
954: . -mat_mkl_pardiso_2  - Fill-in reducing ordering for the input matrix
955: . -mat_mkl_pardiso_4  - Preconditioned CGS/CG
956: . -mat_mkl_pardiso_5  - User permutation
957: . -mat_mkl_pardiso_6  - Write solution on x
958: . -mat_mkl_pardiso_8  - Iterative refinement step
959: . -mat_mkl_pardiso_10 - Pivoting perturbation
960: . -mat_mkl_pardiso_11 - Scaling vectors
961: . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A
962: . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching
963: . -mat_mkl_pardiso_18 - Numbers of non-zero elements
964: . -mat_mkl_pardiso_19 - Report number of floating point operations
965: . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices
966: . -mat_mkl_pardiso_24 - Parallel factorization control
967: . -mat_mkl_pardiso_25 - Parallel forward/backward solve control
968: . -mat_mkl_pardiso_27 - Matrix checker
969: . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors
970: . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode
971: - -mat_mkl_pardiso_60 - Intel MKL PARDISO mode

973:   Level: beginner

975:   Notes:
976:   Use `-mat_mkl_pardiso_68 1` to display the number of threads the solver is using. MKL does not provide a way to directly access this
977:   information.

979:   For more information on the options check the MKL PARDISO manual

981: .seealso: [](ch_matrices), `Mat`, `MATSEQAIJ`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatMkl_PardisoSetCntl()`, `MATSOLVERMKL_CPARDISO`
982: M*/
983: static PetscErrorCode MatFactorGetSolverType_mkl_pardiso(Mat A, MatSolverType *type)
984: {
985:   PetscFunctionBegin;
986:   *type = MATSOLVERMKL_PARDISO;
987:   PetscFunctionReturn(PETSC_SUCCESS);
988: }

990: PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A, MatFactorType ftype, Mat *F)
991: {
992:   Mat              B;
993:   Mat_MKL_PARDISO *mat_mkl_pardiso;
994:   PetscBool        isSeqAIJ, isSeqBAIJ, isSeqSBAIJ;

996:   PetscFunctionBegin;
997:   PetscCall(PetscObjectBaseTypeCompare((PetscObject)A, MATSEQAIJ, &isSeqAIJ));
998:   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQBAIJ, &isSeqBAIJ));
999:   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSeqSBAIJ));
1000:   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
1001:   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1002:   PetscCall(PetscStrallocpy("mkl_pardiso", &((PetscObject)B)->type_name));
1003:   PetscCall(MatSetUp(B));

1005:   PetscCall(PetscNew(&mat_mkl_pardiso));
1006:   B->data = mat_mkl_pardiso;

1008:   PetscCall(MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso));
1009:   if (ftype == MAT_FACTOR_LU) {
1010:     B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO;
1011:     B->factortype            = MAT_FACTOR_LU;
1012:     mat_mkl_pardiso->needsym = PETSC_FALSE;
1013:     if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1014:     else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1015:     else {
1016:       PetscCheck(!isSeqSBAIJ, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for MKL PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead");
1017:       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for MKL PARDISO LU with %s format", ((PetscObject)A)->type_name);
1018:     }
1019: #if defined(PETSC_USE_COMPLEX)
1020:     mat_mkl_pardiso->mtype = 13;
1021: #else
1022:     mat_mkl_pardiso->mtype = 11;
1023: #endif
1024:   } else {
1025:     B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO;
1026:     B->factortype                  = MAT_FACTOR_CHOLESKY;
1027:     if (isSeqAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij;
1028:     else if (isSeqBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij;
1029:     else if (isSeqSBAIJ) mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij;
1030:     else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for PARDISO CHOLESKY with %s format", ((PetscObject)A)->type_name);

1032:     mat_mkl_pardiso->needsym = PETSC_TRUE;
1033: #if !defined(PETSC_USE_COMPLEX)
1034:     if (A->spd == PETSC_BOOL3_TRUE) mat_mkl_pardiso->mtype = 2;
1035:     else mat_mkl_pardiso->mtype = -2;
1036: #else
1037:     mat_mkl_pardiso->mtype = 6;
1038:     PetscCheck(A->hermitian != PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for MKL PARDISO CHOLESKY with Hermitian matrices! Use MAT_FACTOR_LU instead");
1039: #endif
1040:   }
1041:   B->ops->destroy = MatDestroy_MKL_PARDISO;
1042:   B->ops->view    = MatView_MKL_PARDISO;
1043:   B->ops->getinfo = MatGetInfo_MKL_PARDISO;
1044:   B->factortype   = ftype;
1045:   B->assembled    = PETSC_TRUE;

1047:   PetscCall(PetscFree(B->solvertype));
1048:   PetscCall(PetscStrallocpy(MATSOLVERMKL_PARDISO, &B->solvertype));

1050:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mkl_pardiso));
1051:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorSetSchurIS_C", MatFactorSetSchurIS_MKL_PARDISO));
1052:   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMkl_PardisoSetCntl_C", MatMkl_PardisoSetCntl_MKL_PARDISO));

1054:   *F = B;
1055:   PetscFunctionReturn(PETSC_SUCCESS);
1056: }

1058: PETSC_INTERN PetscErrorCode MatSolverTypeRegister_MKL_Pardiso(void)
1059: {
1060:   PetscFunctionBegin;
1061:   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_PARDISO, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mkl_pardiso));
1062:   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_PARDISO, MATSEQAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mkl_pardiso));
1063:   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_PARDISO, MATSEQBAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mkl_pardiso));
1064:   PetscCall(MatSolverTypeRegister(MATSOLVERMKL_PARDISO, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_aij_mkl_pardiso));
1065:   PetscFunctionReturn(PETSC_SUCCESS);
1066: }