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