1 2 #include <petscsys.h> 3 #include <../src/mat/impls/aij/mpi/mpiaij.h> /*I "petscmat.h" I*/ 4 #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 5 6 #if defined(PETSC_HAVE_MKL_INTEL_ILP64) 7 #define MKL_ILP64 8 #endif 9 #include <mkl.h> 10 #include <mkl_cluster_sparse_solver.h> 11 12 /* 13 * Possible mkl_cpardiso phases that controls the execution of the solver. 14 * For more information check mkl_cpardiso 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 27 28 #define IPARM_SIZE 64 29 #define INT_TYPE MKL_INT 30 31 static const char *Err_MSG_CPardiso(int errNo) 32 { 33 switch (errNo) { 34 case -1: 35 return "input inconsistent"; 36 break; 37 case -2: 38 return "not enough memory"; 39 break; 40 case -3: 41 return "reordering problem"; 42 break; 43 case -4: 44 return "zero pivot, numerical factorization or iterative refinement problem"; 45 break; 46 case -5: 47 return "unclassified (internal) error"; 48 break; 49 case -6: 50 return "preordering failed (matrix types 11, 13 only)"; 51 break; 52 case -7: 53 return "diagonal matrix problem"; 54 break; 55 case -8: 56 return "32-bit integer overflow problem"; 57 break; 58 case -9: 59 return "not enough memory for OOC"; 60 break; 61 case -10: 62 return "problems with opening OOC temporary files"; 63 break; 64 case -11: 65 return "read/write problems with the OOC data file"; 66 break; 67 default: 68 return "unknown error"; 69 } 70 } 71 72 /* 73 * Internal data structure. 74 * For more information check mkl_cpardiso manual. 75 */ 76 77 typedef struct { 78 /* Configuration vector */ 79 INT_TYPE iparm[IPARM_SIZE]; 80 81 /* 82 * Internal mkl_cpardiso memory location. 83 * After the first call to mkl_cpardiso do not modify pt, as that could cause a serious memory leak. 84 */ 85 void *pt[IPARM_SIZE]; 86 87 MPI_Fint comm_mkl_cpardiso; 88 89 /* Basic mkl_cpardiso info*/ 90 INT_TYPE phase, maxfct, mnum, mtype, n, nrhs, msglvl, err; 91 92 /* Matrix structure */ 93 PetscScalar *a; 94 95 INT_TYPE *ia, *ja; 96 97 /* Number of non-zero elements */ 98 INT_TYPE nz; 99 100 /* Row permutaton vector*/ 101 INT_TYPE *perm; 102 103 /* Define is matrix preserve sparce structure. */ 104 MatStructure matstruc; 105 106 PetscErrorCode (*ConvertToTriples)(Mat, MatReuse, PetscInt *, PetscInt **, PetscInt **, PetscScalar **); 107 108 /* True if mkl_cpardiso function have been used. */ 109 PetscBool CleanUp; 110 } Mat_MKL_CPARDISO; 111 112 /* 113 * Copy the elements of matrix A. 114 * Input: 115 * - Mat A: MATSEQAIJ matrix 116 * - int shift: matrix index. 117 * - 0 for c representation 118 * - 1 for fortran representation 119 * - MatReuse reuse: 120 * - MAT_INITIAL_MATRIX: Create a new aij representation 121 * - MAT_REUSE_MATRIX: Reuse all aij representation and just change values 122 * Output: 123 * - int *nnz: Number of nonzero-elements. 124 * - int **r pointer to i index 125 * - int **c pointer to j elements 126 * - MATRIXTYPE **v: Non-zero elements 127 */ 128 PetscErrorCode MatCopy_seqaij_seqaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v) 129 { 130 Mat_SeqAIJ *aa = (Mat_SeqAIJ *)A->data; 131 132 PetscFunctionBegin; 133 *v = aa->a; 134 if (reuse == MAT_INITIAL_MATRIX) { 135 *r = (INT_TYPE *)aa->i; 136 *c = (INT_TYPE *)aa->j; 137 *nnz = aa->nz; 138 } 139 PetscFunctionReturn(0); 140 } 141 142 PetscErrorCode MatConvertToTriples_mpiaij_mpiaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v) 143 { 144 const PetscInt *ai, *aj, *bi, *bj, *garray, m = A->rmap->n, *ajj, *bjj; 145 PetscInt rstart, nz, i, j, countA, countB; 146 PetscInt *row, *col; 147 const PetscScalar *av, *bv; 148 PetscScalar *val; 149 Mat_MPIAIJ *mat = (Mat_MPIAIJ *)A->data; 150 Mat_SeqAIJ *aa = (Mat_SeqAIJ *)(mat->A)->data; 151 Mat_SeqAIJ *bb = (Mat_SeqAIJ *)(mat->B)->data; 152 PetscInt colA_start, jB, jcol; 153 154 PetscFunctionBegin; 155 ai = aa->i; 156 aj = aa->j; 157 bi = bb->i; 158 bj = bb->j; 159 rstart = A->rmap->rstart; 160 av = aa->a; 161 bv = bb->a; 162 163 garray = mat->garray; 164 165 if (reuse == MAT_INITIAL_MATRIX) { 166 nz = aa->nz + bb->nz; 167 *nnz = nz; 168 PetscCall(PetscMalloc3(m + 1, &row, nz, &col, nz, &val)); 169 *r = row; 170 *c = col; 171 *v = val; 172 } else { 173 row = *r; 174 col = *c; 175 val = *v; 176 } 177 178 nz = 0; 179 for (i = 0; i < m; i++) { 180 row[i] = nz; 181 countA = ai[i + 1] - ai[i]; 182 countB = bi[i + 1] - bi[i]; 183 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 184 bjj = bj + bi[i]; 185 186 /* B part, smaller col index */ 187 colA_start = rstart + ajj[0]; /* the smallest global col index of A */ 188 jB = 0; 189 for (j = 0; j < countB; j++) { 190 jcol = garray[bjj[j]]; 191 if (jcol > colA_start) break; 192 col[nz] = jcol; 193 val[nz++] = *bv++; 194 } 195 jB = j; 196 197 /* A part */ 198 for (j = 0; j < countA; j++) { 199 col[nz] = rstart + ajj[j]; 200 val[nz++] = *av++; 201 } 202 203 /* B part, larger col index */ 204 for (j = jB; j < countB; j++) { 205 col[nz] = garray[bjj[j]]; 206 val[nz++] = *bv++; 207 } 208 } 209 row[m] = nz; 210 211 PetscFunctionReturn(0); 212 } 213 214 PetscErrorCode MatConvertToTriples_mpibaij_mpibaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v) 215 { 216 const PetscInt *ai, *aj, *bi, *bj, *garray, bs = A->rmap->bs, bs2 = bs * bs, m = A->rmap->n / bs, *ajj, *bjj; 217 PetscInt rstart, nz, i, j, countA, countB; 218 PetscInt *row, *col; 219 const PetscScalar *av, *bv; 220 PetscScalar *val; 221 Mat_MPIBAIJ *mat = (Mat_MPIBAIJ *)A->data; 222 Mat_SeqBAIJ *aa = (Mat_SeqBAIJ *)(mat->A)->data; 223 Mat_SeqBAIJ *bb = (Mat_SeqBAIJ *)(mat->B)->data; 224 PetscInt colA_start, jB, jcol; 225 226 PetscFunctionBegin; 227 ai = aa->i; 228 aj = aa->j; 229 bi = bb->i; 230 bj = bb->j; 231 rstart = A->rmap->rstart / bs; 232 av = aa->a; 233 bv = bb->a; 234 235 garray = mat->garray; 236 237 if (reuse == MAT_INITIAL_MATRIX) { 238 nz = aa->nz + bb->nz; 239 *nnz = nz; 240 PetscCall(PetscMalloc3(m + 1, &row, nz, &col, nz * bs2, &val)); 241 *r = row; 242 *c = col; 243 *v = val; 244 } else { 245 row = *r; 246 col = *c; 247 val = *v; 248 } 249 250 nz = 0; 251 for (i = 0; i < m; i++) { 252 row[i] = nz + 1; 253 countA = ai[i + 1] - ai[i]; 254 countB = bi[i + 1] - bi[i]; 255 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 256 bjj = bj + bi[i]; 257 258 /* B part, smaller col index */ 259 colA_start = rstart + (countA > 0 ? ajj[0] : 0); /* the smallest global col index of A */ 260 jB = 0; 261 for (j = 0; j < countB; j++) { 262 jcol = garray[bjj[j]]; 263 if (jcol > colA_start) break; 264 col[nz++] = jcol + 1; 265 } 266 jB = j; 267 PetscCall(PetscArraycpy(val, bv, jB * bs2)); 268 val += jB * bs2; 269 bv += jB * bs2; 270 271 /* A part */ 272 for (j = 0; j < countA; j++) col[nz++] = rstart + ajj[j] + 1; 273 PetscCall(PetscArraycpy(val, av, countA * bs2)); 274 val += countA * bs2; 275 av += countA * bs2; 276 277 /* B part, larger col index */ 278 for (j = jB; j < countB; j++) col[nz++] = garray[bjj[j]] + 1; 279 PetscCall(PetscArraycpy(val, bv, (countB - jB) * bs2)); 280 val += (countB - jB) * bs2; 281 bv += (countB - jB) * bs2; 282 } 283 row[m] = nz + 1; 284 285 PetscFunctionReturn(0); 286 } 287 288 PetscErrorCode MatConvertToTriples_mpisbaij_mpisbaij_MKL_CPARDISO(Mat A, MatReuse reuse, PetscInt *nnz, PetscInt **r, PetscInt **c, PetscScalar **v) 289 { 290 const PetscInt *ai, *aj, *bi, *bj, *garray, bs = A->rmap->bs, bs2 = bs * bs, m = A->rmap->n / bs, *ajj, *bjj; 291 PetscInt rstart, nz, i, j, countA, countB; 292 PetscInt *row, *col; 293 const PetscScalar *av, *bv; 294 PetscScalar *val; 295 Mat_MPISBAIJ *mat = (Mat_MPISBAIJ *)A->data; 296 Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ *)(mat->A)->data; 297 Mat_SeqBAIJ *bb = (Mat_SeqBAIJ *)(mat->B)->data; 298 299 PetscFunctionBegin; 300 ai = aa->i; 301 aj = aa->j; 302 bi = bb->i; 303 bj = bb->j; 304 rstart = A->rmap->rstart / bs; 305 av = aa->a; 306 bv = bb->a; 307 308 garray = mat->garray; 309 310 if (reuse == MAT_INITIAL_MATRIX) { 311 nz = aa->nz + bb->nz; 312 *nnz = nz; 313 PetscCall(PetscMalloc3(m + 1, &row, nz, &col, nz * bs2, &val)); 314 *r = row; 315 *c = col; 316 *v = val; 317 } else { 318 row = *r; 319 col = *c; 320 val = *v; 321 } 322 323 nz = 0; 324 for (i = 0; i < m; i++) { 325 row[i] = nz + 1; 326 countA = ai[i + 1] - ai[i]; 327 countB = bi[i + 1] - bi[i]; 328 ajj = aj + ai[i]; /* ptr to the beginning of this row */ 329 bjj = bj + bi[i]; 330 331 /* A part */ 332 for (j = 0; j < countA; j++) col[nz++] = rstart + ajj[j] + 1; 333 PetscCall(PetscArraycpy(val, av, countA * bs2)); 334 val += countA * bs2; 335 av += countA * bs2; 336 337 /* B part, larger col index */ 338 for (j = 0; j < countB; j++) col[nz++] = garray[bjj[j]] + 1; 339 PetscCall(PetscArraycpy(val, bv, countB * bs2)); 340 val += countB * bs2; 341 bv += countB * bs2; 342 } 343 row[m] = nz + 1; 344 345 PetscFunctionReturn(0); 346 } 347 348 /* 349 * Free memory for Mat_MKL_CPARDISO structure and pointers to objects. 350 */ 351 PetscErrorCode MatDestroy_MKL_CPARDISO(Mat A) 352 { 353 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data; 354 MPI_Comm comm; 355 356 PetscFunctionBegin; 357 /* Terminate instance, deallocate memories */ 358 if (mat_mkl_cpardiso->CleanUp) { 359 mat_mkl_cpardiso->phase = JOB_RELEASE_OF_ALL_MEMORY; 360 361 cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, NULL, NULL, NULL, mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, 362 mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err); 363 } 364 365 if (mat_mkl_cpardiso->ConvertToTriples != MatCopy_seqaij_seqaij_MKL_CPARDISO) PetscCall(PetscFree3(mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, mat_mkl_cpardiso->a)); 366 comm = MPI_Comm_f2c(mat_mkl_cpardiso->comm_mkl_cpardiso); 367 PetscCallMPI(MPI_Comm_free(&comm)); 368 PetscCall(PetscFree(A->data)); 369 370 /* clear composed functions */ 371 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL)); 372 PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatMkl_CPardisoSetCntl_C", NULL)); 373 PetscFunctionReturn(0); 374 } 375 376 /* 377 * Computes Ax = b 378 */ 379 PetscErrorCode MatSolve_MKL_CPARDISO(Mat A, Vec b, Vec x) 380 { 381 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)(A)->data; 382 PetscScalar *xarray; 383 const PetscScalar *barray; 384 385 PetscFunctionBegin; 386 mat_mkl_cpardiso->nrhs = 1; 387 PetscCall(VecGetArray(x, &xarray)); 388 PetscCall(VecGetArrayRead(b, &barray)); 389 390 /* solve phase */ 391 /*-------------*/ 392 mat_mkl_cpardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 393 cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, 394 mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err); 395 396 PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL_CPARDISO: err=%d, msg = \"%s\". Please check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err)); 397 398 PetscCall(VecRestoreArray(x, &xarray)); 399 PetscCall(VecRestoreArrayRead(b, &barray)); 400 mat_mkl_cpardiso->CleanUp = PETSC_TRUE; 401 PetscFunctionReturn(0); 402 } 403 404 PetscErrorCode MatSolveTranspose_MKL_CPARDISO(Mat A, Vec b, Vec x) 405 { 406 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data; 407 408 PetscFunctionBegin; 409 #if defined(PETSC_USE_COMPLEX) 410 mat_mkl_cpardiso->iparm[12 - 1] = 1; 411 #else 412 mat_mkl_cpardiso->iparm[12 - 1] = 2; 413 #endif 414 PetscCall(MatSolve_MKL_CPARDISO(A, b, x)); 415 mat_mkl_cpardiso->iparm[12 - 1] = 0; 416 PetscFunctionReturn(0); 417 } 418 419 PetscErrorCode MatMatSolve_MKL_CPARDISO(Mat A, Mat B, Mat X) 420 { 421 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)(A)->data; 422 PetscScalar *xarray; 423 const PetscScalar *barray; 424 425 PetscFunctionBegin; 426 PetscCall(MatGetSize(B, NULL, (PetscInt *)&mat_mkl_cpardiso->nrhs)); 427 428 if (mat_mkl_cpardiso->nrhs > 0) { 429 PetscCall(MatDenseGetArrayRead(B, &barray)); 430 PetscCall(MatDenseGetArray(X, &xarray)); 431 432 PetscCheck(barray != xarray, PETSC_COMM_SELF, PETSC_ERR_SUP, "B and X cannot share the same memory location"); 433 434 /* solve phase */ 435 /*-------------*/ 436 mat_mkl_cpardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 437 cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, 438 mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, (void *)barray, (void *)xarray, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err); 439 PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL_CPARDISO: err=%d, msg = \"%s\". Please check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err)); 440 PetscCall(MatDenseRestoreArrayRead(B, &barray)); 441 PetscCall(MatDenseRestoreArray(X, &xarray)); 442 } 443 mat_mkl_cpardiso->CleanUp = PETSC_TRUE; 444 PetscFunctionReturn(0); 445 } 446 447 /* 448 * LU Decomposition 449 */ 450 PetscErrorCode MatFactorNumeric_MKL_CPARDISO(Mat F, Mat A, const MatFactorInfo *info) 451 { 452 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)(F)->data; 453 454 PetscFunctionBegin; 455 mat_mkl_cpardiso->matstruc = SAME_NONZERO_PATTERN; 456 PetscCall((*mat_mkl_cpardiso->ConvertToTriples)(A, MAT_REUSE_MATRIX, &mat_mkl_cpardiso->nz, &mat_mkl_cpardiso->ia, &mat_mkl_cpardiso->ja, &mat_mkl_cpardiso->a)); 457 458 mat_mkl_cpardiso->phase = JOB_NUMERICAL_FACTORIZATION; 459 cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, 460 mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, &mat_mkl_cpardiso->err); 461 PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL_CPARDISO: err=%d, msg = \"%s\". Please check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err)); 462 463 mat_mkl_cpardiso->matstruc = SAME_NONZERO_PATTERN; 464 mat_mkl_cpardiso->CleanUp = PETSC_TRUE; 465 PetscFunctionReturn(0); 466 } 467 468 /* Sets mkl_cpardiso options from the options database */ 469 PetscErrorCode MatSetFromOptions_MKL_CPARDISO(Mat F, Mat A) 470 { 471 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data; 472 PetscInt icntl, threads; 473 PetscBool flg; 474 475 PetscFunctionBegin; 476 PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "MKL_CPARDISO Options", "Mat"); 477 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_65", "Suggested number of threads to use within MKL_CPARDISO", "None", threads, &threads, &flg)); 478 if (flg) mkl_set_num_threads((int)threads); 479 480 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_66", "Maximum number of factors with identical sparsity structure that must be kept in memory at the same time", "None", mat_mkl_cpardiso->maxfct, &icntl, &flg)); 481 if (flg) mat_mkl_cpardiso->maxfct = icntl; 482 483 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_67", "Indicates the actual matrix for the solution phase", "None", mat_mkl_cpardiso->mnum, &icntl, &flg)); 484 if (flg) mat_mkl_cpardiso->mnum = icntl; 485 486 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_68", "Message level information", "None", mat_mkl_cpardiso->msglvl, &icntl, &flg)); 487 if (flg) mat_mkl_cpardiso->msglvl = icntl; 488 489 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_69", "Defines the matrix type", "None", mat_mkl_cpardiso->mtype, &icntl, &flg)); 490 if (flg) mat_mkl_cpardiso->mtype = icntl; 491 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_1", "Use default values", "None", mat_mkl_cpardiso->iparm[0], &icntl, &flg)); 492 493 if (flg && icntl != 0) { 494 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_2", "Fill-in reducing ordering for the input matrix", "None", mat_mkl_cpardiso->iparm[1], &icntl, &flg)); 495 if (flg) mat_mkl_cpardiso->iparm[1] = icntl; 496 497 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_4", "Preconditioned CGS/CG", "None", mat_mkl_cpardiso->iparm[3], &icntl, &flg)); 498 if (flg) mat_mkl_cpardiso->iparm[3] = icntl; 499 500 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_5", "User permutation", "None", mat_mkl_cpardiso->iparm[4], &icntl, &flg)); 501 if (flg) mat_mkl_cpardiso->iparm[4] = icntl; 502 503 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_6", "Write solution on x", "None", mat_mkl_cpardiso->iparm[5], &icntl, &flg)); 504 if (flg) mat_mkl_cpardiso->iparm[5] = icntl; 505 506 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_8", "Iterative refinement step", "None", mat_mkl_cpardiso->iparm[7], &icntl, &flg)); 507 if (flg) mat_mkl_cpardiso->iparm[7] = icntl; 508 509 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_10", "Pivoting perturbation", "None", mat_mkl_cpardiso->iparm[9], &icntl, &flg)); 510 if (flg) mat_mkl_cpardiso->iparm[9] = icntl; 511 512 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_11", "Scaling vectors", "None", mat_mkl_cpardiso->iparm[10], &icntl, &flg)); 513 if (flg) mat_mkl_cpardiso->iparm[10] = icntl; 514 515 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_12", "Solve with transposed or conjugate transposed matrix A", "None", mat_mkl_cpardiso->iparm[11], &icntl, &flg)); 516 if (flg) mat_mkl_cpardiso->iparm[11] = icntl; 517 518 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_13", "Improved accuracy using (non-) symmetric weighted matching", "None", mat_mkl_cpardiso->iparm[12], &icntl, &flg)); 519 if (flg) mat_mkl_cpardiso->iparm[12] = icntl; 520 521 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_18", "Numbers of non-zero elements", "None", mat_mkl_cpardiso->iparm[17], &icntl, &flg)); 522 if (flg) mat_mkl_cpardiso->iparm[17] = icntl; 523 524 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_19", "Report number of floating point operations", "None", mat_mkl_cpardiso->iparm[18], &icntl, &flg)); 525 if (flg) mat_mkl_cpardiso->iparm[18] = icntl; 526 527 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_21", "Pivoting for symmetric indefinite matrices", "None", mat_mkl_cpardiso->iparm[20], &icntl, &flg)); 528 if (flg) mat_mkl_cpardiso->iparm[20] = icntl; 529 530 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_24", "Parallel factorization control", "None", mat_mkl_cpardiso->iparm[23], &icntl, &flg)); 531 if (flg) mat_mkl_cpardiso->iparm[23] = icntl; 532 533 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_25", "Parallel forward/backward solve control", "None", mat_mkl_cpardiso->iparm[24], &icntl, &flg)); 534 if (flg) mat_mkl_cpardiso->iparm[24] = icntl; 535 536 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_27", "Matrix checker", "None", mat_mkl_cpardiso->iparm[26], &icntl, &flg)); 537 if (flg) mat_mkl_cpardiso->iparm[26] = icntl; 538 539 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_31", "Partial solve and computing selected components of the solution vectors", "None", mat_mkl_cpardiso->iparm[30], &icntl, &flg)); 540 if (flg) mat_mkl_cpardiso->iparm[30] = icntl; 541 542 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_34", "Optimal number of threads for conditional numerical reproducibility (CNR) mode", "None", mat_mkl_cpardiso->iparm[33], &icntl, &flg)); 543 if (flg) mat_mkl_cpardiso->iparm[33] = icntl; 544 545 PetscCall(PetscOptionsInt("-mat_mkl_cpardiso_60", "Intel MKL_CPARDISO mode", "None", mat_mkl_cpardiso->iparm[59], &icntl, &flg)); 546 if (flg) mat_mkl_cpardiso->iparm[59] = icntl; 547 } 548 549 PetscOptionsEnd(); 550 PetscFunctionReturn(0); 551 } 552 553 PetscErrorCode PetscInitialize_MKL_CPARDISO(Mat A, Mat_MKL_CPARDISO *mat_mkl_cpardiso) 554 { 555 PetscInt bs; 556 PetscBool match; 557 PetscMPIInt size; 558 MPI_Comm comm; 559 560 PetscFunctionBegin; 561 562 PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)A), &comm)); 563 PetscCallMPI(MPI_Comm_size(comm, &size)); 564 mat_mkl_cpardiso->comm_mkl_cpardiso = MPI_Comm_c2f(comm); 565 566 mat_mkl_cpardiso->CleanUp = PETSC_FALSE; 567 mat_mkl_cpardiso->maxfct = 1; 568 mat_mkl_cpardiso->mnum = 1; 569 mat_mkl_cpardiso->n = A->rmap->N; 570 if (mat_mkl_cpardiso->iparm[36]) mat_mkl_cpardiso->n /= mat_mkl_cpardiso->iparm[36]; 571 mat_mkl_cpardiso->msglvl = 0; 572 mat_mkl_cpardiso->nrhs = 1; 573 mat_mkl_cpardiso->err = 0; 574 mat_mkl_cpardiso->phase = -1; 575 #if defined(PETSC_USE_COMPLEX) 576 mat_mkl_cpardiso->mtype = 13; 577 #else 578 mat_mkl_cpardiso->mtype = 11; 579 #endif 580 581 #if defined(PETSC_USE_REAL_SINGLE) 582 mat_mkl_cpardiso->iparm[27] = 1; 583 #else 584 mat_mkl_cpardiso->iparm[27] = 0; 585 #endif 586 587 mat_mkl_cpardiso->iparm[0] = 1; /* Solver default parameters overriden with provided by iparm */ 588 mat_mkl_cpardiso->iparm[1] = 2; /* Use METIS for fill-in reordering */ 589 mat_mkl_cpardiso->iparm[5] = 0; /* Write solution into x */ 590 mat_mkl_cpardiso->iparm[7] = 2; /* Max number of iterative refinement steps */ 591 mat_mkl_cpardiso->iparm[9] = 13; /* Perturb the pivot elements with 1E-13 */ 592 mat_mkl_cpardiso->iparm[10] = 1; /* Use nonsymmetric permutation and scaling MPS */ 593 mat_mkl_cpardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 594 mat_mkl_cpardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */ 595 mat_mkl_cpardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */ 596 mat_mkl_cpardiso->iparm[26] = 1; /* Check input data for correctness */ 597 598 mat_mkl_cpardiso->iparm[39] = 0; 599 if (size > 1) { 600 mat_mkl_cpardiso->iparm[39] = 2; 601 mat_mkl_cpardiso->iparm[40] = A->rmap->rstart; 602 mat_mkl_cpardiso->iparm[41] = A->rmap->rend - 1; 603 } 604 PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &match, MATMPIBAIJ, MATMPISBAIJ, "")); 605 if (match) { 606 PetscCall(MatGetBlockSize(A, &bs)); 607 mat_mkl_cpardiso->iparm[36] = bs; 608 mat_mkl_cpardiso->iparm[40] /= bs; 609 mat_mkl_cpardiso->iparm[41] /= bs; 610 mat_mkl_cpardiso->iparm[40]++; 611 mat_mkl_cpardiso->iparm[41]++; 612 mat_mkl_cpardiso->iparm[34] = 0; /* Fortran style */ 613 } else { 614 mat_mkl_cpardiso->iparm[34] = 1; /* C style */ 615 } 616 617 mat_mkl_cpardiso->perm = 0; 618 PetscFunctionReturn(0); 619 } 620 621 /* 622 * Symbolic decomposition. Mkl_Pardiso analysis phase. 623 */ 624 PetscErrorCode MatLUFactorSymbolic_AIJMKL_CPARDISO(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info) 625 { 626 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data; 627 628 PetscFunctionBegin; 629 mat_mkl_cpardiso->matstruc = DIFFERENT_NONZERO_PATTERN; 630 631 /* Set MKL_CPARDISO options from the options database */ 632 PetscCall(MatSetFromOptions_MKL_CPARDISO(F, A)); 633 PetscCall((*mat_mkl_cpardiso->ConvertToTriples)(A, MAT_INITIAL_MATRIX, &mat_mkl_cpardiso->nz, &mat_mkl_cpardiso->ia, &mat_mkl_cpardiso->ja, &mat_mkl_cpardiso->a)); 634 635 mat_mkl_cpardiso->n = A->rmap->N; 636 if (mat_mkl_cpardiso->iparm[36]) mat_mkl_cpardiso->n /= mat_mkl_cpardiso->iparm[36]; 637 638 /* analysis phase */ 639 /*----------------*/ 640 mat_mkl_cpardiso->phase = JOB_ANALYSIS; 641 642 cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, 643 mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err); 644 645 PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL_CPARDISO: err=%d, msg = \"%s\".Check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err)); 646 647 mat_mkl_cpardiso->CleanUp = PETSC_TRUE; 648 F->ops->lufactornumeric = MatFactorNumeric_MKL_CPARDISO; 649 F->ops->solve = MatSolve_MKL_CPARDISO; 650 F->ops->solvetranspose = MatSolveTranspose_MKL_CPARDISO; 651 F->ops->matsolve = MatMatSolve_MKL_CPARDISO; 652 PetscFunctionReturn(0); 653 } 654 655 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_CPARDISO(Mat F, Mat A, IS perm, const MatFactorInfo *info) 656 { 657 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data; 658 659 PetscFunctionBegin; 660 mat_mkl_cpardiso->matstruc = DIFFERENT_NONZERO_PATTERN; 661 662 /* Set MKL_CPARDISO options from the options database */ 663 PetscCall(MatSetFromOptions_MKL_CPARDISO(F, A)); 664 PetscCall((*mat_mkl_cpardiso->ConvertToTriples)(A, MAT_INITIAL_MATRIX, &mat_mkl_cpardiso->nz, &mat_mkl_cpardiso->ia, &mat_mkl_cpardiso->ja, &mat_mkl_cpardiso->a)); 665 666 mat_mkl_cpardiso->n = A->rmap->N; 667 if (mat_mkl_cpardiso->iparm[36]) mat_mkl_cpardiso->n /= mat_mkl_cpardiso->iparm[36]; 668 PetscCheck(!PetscDefined(USE_COMPLEX), PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for PARDISO CHOLESKY with complex scalars! Use MAT_FACTOR_LU instead"); 669 if (A->spd == PETSC_BOOL3_TRUE) mat_mkl_cpardiso->mtype = 2; 670 else mat_mkl_cpardiso->mtype = -2; 671 672 /* analysis phase */ 673 /*----------------*/ 674 mat_mkl_cpardiso->phase = JOB_ANALYSIS; 675 676 cluster_sparse_solver(mat_mkl_cpardiso->pt, &mat_mkl_cpardiso->maxfct, &mat_mkl_cpardiso->mnum, &mat_mkl_cpardiso->mtype, &mat_mkl_cpardiso->phase, &mat_mkl_cpardiso->n, mat_mkl_cpardiso->a, mat_mkl_cpardiso->ia, mat_mkl_cpardiso->ja, 677 mat_mkl_cpardiso->perm, &mat_mkl_cpardiso->nrhs, mat_mkl_cpardiso->iparm, &mat_mkl_cpardiso->msglvl, NULL, NULL, &mat_mkl_cpardiso->comm_mkl_cpardiso, (PetscInt *)&mat_mkl_cpardiso->err); 678 679 PetscCheck(mat_mkl_cpardiso->err >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by MKL_CPARDISO: err=%d, msg = \"%s\".Check manual", mat_mkl_cpardiso->err, Err_MSG_CPardiso(mat_mkl_cpardiso->err)); 680 681 mat_mkl_cpardiso->CleanUp = PETSC_TRUE; 682 F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_CPARDISO; 683 F->ops->solve = MatSolve_MKL_CPARDISO; 684 F->ops->solvetranspose = MatSolveTranspose_MKL_CPARDISO; 685 F->ops->matsolve = MatMatSolve_MKL_CPARDISO; 686 PetscFunctionReturn(0); 687 } 688 689 PetscErrorCode MatView_MKL_CPARDISO(Mat A, PetscViewer viewer) 690 { 691 PetscBool iascii; 692 PetscViewerFormat format; 693 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data; 694 PetscInt i; 695 696 PetscFunctionBegin; 697 /* check if matrix is mkl_cpardiso type */ 698 if (A->ops->solve != MatSolve_MKL_CPARDISO) PetscFunctionReturn(0); 699 700 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 701 if (iascii) { 702 PetscCall(PetscViewerGetFormat(viewer, &format)); 703 if (format == PETSC_VIEWER_ASCII_INFO) { 704 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO run parameters:\n")); 705 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO phase: %d \n", mat_mkl_cpardiso->phase)); 706 for (i = 1; i <= 64; i++) PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO iparm[%d]: %d \n", i, mat_mkl_cpardiso->iparm[i - 1])); 707 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO maxfct: %d \n", mat_mkl_cpardiso->maxfct)); 708 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO mnum: %d \n", mat_mkl_cpardiso->mnum)); 709 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO mtype: %d \n", mat_mkl_cpardiso->mtype)); 710 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO n: %d \n", mat_mkl_cpardiso->n)); 711 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO nrhs: %d \n", mat_mkl_cpardiso->nrhs)); 712 PetscCall(PetscViewerASCIIPrintf(viewer, "MKL_CPARDISO msglvl: %d \n", mat_mkl_cpardiso->msglvl)); 713 } 714 } 715 PetscFunctionReturn(0); 716 } 717 718 PetscErrorCode MatGetInfo_MKL_CPARDISO(Mat A, MatInfoType flag, MatInfo *info) 719 { 720 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)A->data; 721 722 PetscFunctionBegin; 723 info->block_size = 1.0; 724 info->nz_allocated = mat_mkl_cpardiso->nz + 0.0; 725 info->nz_unneeded = 0.0; 726 info->assemblies = 0.0; 727 info->mallocs = 0.0; 728 info->memory = 0.0; 729 info->fill_ratio_given = 0; 730 info->fill_ratio_needed = 0; 731 info->factor_mallocs = 0; 732 PetscFunctionReturn(0); 733 } 734 735 PetscErrorCode MatMkl_CPardisoSetCntl_MKL_CPARDISO(Mat F, PetscInt icntl, PetscInt ival) 736 { 737 Mat_MKL_CPARDISO *mat_mkl_cpardiso = (Mat_MKL_CPARDISO *)F->data; 738 739 PetscFunctionBegin; 740 if (icntl <= 64) { 741 mat_mkl_cpardiso->iparm[icntl - 1] = ival; 742 } else { 743 if (icntl == 65) mkl_set_num_threads((int)ival); 744 else if (icntl == 66) mat_mkl_cpardiso->maxfct = ival; 745 else if (icntl == 67) mat_mkl_cpardiso->mnum = ival; 746 else if (icntl == 68) mat_mkl_cpardiso->msglvl = ival; 747 else if (icntl == 69) mat_mkl_cpardiso->mtype = ival; 748 } 749 PetscFunctionReturn(0); 750 } 751 752 /*@ 753 MatMkl_CPardisoSetCntl - Set Mkl_Pardiso parameters 754 755 Logically Collective 756 757 Input Parameters: 758 + F - the factored matrix obtained by calling `MatGetFactor()` 759 . icntl - index of Mkl_Pardiso parameter 760 - ival - value of Mkl_Pardiso parameter 761 762 Options Database Key: 763 . -mat_mkl_cpardiso_<icntl> <ival> - set the option numbered icntl to ival 764 765 Level: Intermediate 766 767 Note: 768 This routine cannot be used if you are solving the linear system with `TS`, `SNES`, or `KSP`, only if you directly call `MatGetFactor()` so use the options 769 database approach when working with `TS`, `SNES`, or `KSP`. 770 771 References: 772 . * - Mkl_Pardiso Users' Guide 773 774 .seealso: `MatGetFactor()`, `MATMPIAIJ`, `MATSOLVERMKL_CPARDISO` 775 @*/ 776 PetscErrorCode MatMkl_CPardisoSetCntl(Mat F, PetscInt icntl, PetscInt ival) 777 { 778 PetscFunctionBegin; 779 PetscTryMethod(F, "MatMkl_CPardisoSetCntl_C", (Mat, PetscInt, PetscInt), (F, icntl, ival)); 780 PetscFunctionReturn(0); 781 } 782 783 /*MC 784 MATSOLVERMKL_CPARDISO - A matrix type providing direct solvers (LU) for parallel matrices via the external package MKL_CPARDISO. 785 786 Works with `MATMPIAIJ` matrices 787 788 Use -pc_type lu -pc_factor_mat_solver_type mkl_cpardiso to use this direct solver 789 790 Options Database Keys: 791 + -mat_mkl_cpardiso_65 - Suggested number of threads to use within MKL_CPARDISO 792 . -mat_mkl_cpardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time 793 . -mat_mkl_cpardiso_67 - Indicates the actual matrix for the solution phase 794 . -mat_mkl_cpardiso_68 - Message level information, use 1 to get detailed information on the solver options 795 . -mat_mkl_cpardiso_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 796 . -mat_mkl_cpardiso_1 - Use default values 797 . -mat_mkl_cpardiso_2 - Fill-in reducing ordering for the input matrix 798 . -mat_mkl_cpardiso_4 - Preconditioned CGS/CG 799 . -mat_mkl_cpardiso_5 - User permutation 800 . -mat_mkl_cpardiso_6 - Write solution on x 801 . -mat_mkl_cpardiso_8 - Iterative refinement step 802 . -mat_mkl_cpardiso_10 - Pivoting perturbation 803 . -mat_mkl_cpardiso_11 - Scaling vectors 804 . -mat_mkl_cpardiso_12 - Solve with transposed or conjugate transposed matrix A 805 . -mat_mkl_cpardiso_13 - Improved accuracy using (non-) symmetric weighted matching 806 . -mat_mkl_cpardiso_18 - Numbers of non-zero elements 807 . -mat_mkl_cpardiso_19 - Report number of floating point operations 808 . -mat_mkl_cpardiso_21 - Pivoting for symmetric indefinite matrices 809 . -mat_mkl_cpardiso_24 - Parallel factorization control 810 . -mat_mkl_cpardiso_25 - Parallel forward/backward solve control 811 . -mat_mkl_cpardiso_27 - Matrix checker 812 . -mat_mkl_cpardiso_31 - Partial solve and computing selected components of the solution vectors 813 . -mat_mkl_cpardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode 814 - -mat_mkl_cpardiso_60 - Intel MKL_CPARDISO mode 815 816 Level: beginner 817 818 Notes: 819 Use -mat_mkl_cpardiso_68 1 to display the number of threads the solver is using. MKL does not provide a way to directly access this 820 information. 821 822 For more information on the options check the MKL_CPARDISO manual 823 824 .seealso: `PCFactorSetMatSolverType()`, `MatSolverType`, `MatMkl_CPardisoSetCntl()`, `MatGetFactor()` 825 M*/ 826 827 static PetscErrorCode MatFactorGetSolverType_mkl_cpardiso(Mat A, MatSolverType *type) 828 { 829 PetscFunctionBegin; 830 *type = MATSOLVERMKL_CPARDISO; 831 PetscFunctionReturn(0); 832 } 833 834 /* MatGetFactor for MPI AIJ matrices */ 835 static PetscErrorCode MatGetFactor_mpiaij_mkl_cpardiso(Mat A, MatFactorType ftype, Mat *F) 836 { 837 Mat B; 838 Mat_MKL_CPARDISO *mat_mkl_cpardiso; 839 PetscBool isSeqAIJ, isMPIBAIJ, isMPISBAIJ; 840 841 PetscFunctionBegin; 842 /* Create the factorization matrix */ 843 844 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSeqAIJ)); 845 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIBAIJ, &isMPIBAIJ)); 846 PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isMPISBAIJ)); 847 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 848 PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 849 PetscCall(PetscStrallocpy("mkl_cpardiso", &((PetscObject)B)->type_name)); 850 PetscCall(MatSetUp(B)); 851 852 PetscCall(PetscNew(&mat_mkl_cpardiso)); 853 854 if (isSeqAIJ) mat_mkl_cpardiso->ConvertToTriples = MatCopy_seqaij_seqaij_MKL_CPARDISO; 855 else if (isMPIBAIJ) mat_mkl_cpardiso->ConvertToTriples = MatConvertToTriples_mpibaij_mpibaij_MKL_CPARDISO; 856 else if (isMPISBAIJ) mat_mkl_cpardiso->ConvertToTriples = MatConvertToTriples_mpisbaij_mpisbaij_MKL_CPARDISO; 857 else mat_mkl_cpardiso->ConvertToTriples = MatConvertToTriples_mpiaij_mpiaij_MKL_CPARDISO; 858 859 if (ftype == MAT_FACTOR_LU) B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_CPARDISO; 860 else B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_CPARDISO; 861 B->ops->destroy = MatDestroy_MKL_CPARDISO; 862 863 B->ops->view = MatView_MKL_CPARDISO; 864 B->ops->getinfo = MatGetInfo_MKL_CPARDISO; 865 866 B->factortype = ftype; 867 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 868 869 B->data = mat_mkl_cpardiso; 870 871 /* set solvertype */ 872 PetscCall(PetscFree(B->solvertype)); 873 PetscCall(PetscStrallocpy(MATSOLVERMKL_CPARDISO, &B->solvertype)); 874 875 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_mkl_cpardiso)); 876 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatMkl_CPardisoSetCntl_C", MatMkl_CPardisoSetCntl_MKL_CPARDISO)); 877 PetscCall(PetscInitialize_MKL_CPARDISO(A, mat_mkl_cpardiso)); 878 879 *F = B; 880 PetscFunctionReturn(0); 881 } 882 883 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_MKL_CPardiso(void) 884 { 885 PetscFunctionBegin; 886 PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_mkl_cpardiso)); 887 PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_mkl_cpardiso)); 888 PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATMPIBAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_mkl_cpardiso)); 889 PetscCall(MatSolverTypeRegister(MATSOLVERMKL_CPARDISO, MATMPISBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_mpiaij_mkl_cpardiso)); 890 PetscFunctionReturn(0); 891 } 892