1 #if defined(PETSC_HAVE_LIBMKL_INTEL_ILP64) 2 #define MKL_ILP64 3 #endif 4 5 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 6 #include <../src/mat/impls/sbaij/seq/sbaij.h> /*I "petscmat.h" I*/ 7 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 8 9 #include <stdio.h> 10 #include <stdlib.h> 11 #include <math.h> 12 #include <mkl.h> 13 14 /* 15 * Possible mkl_pardiso phases that controls the execution of the solver. 16 * For more information check mkl_pardiso manual. 17 */ 18 #define JOB_ANALYSIS 11 19 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION 12 20 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13 21 #define JOB_NUMERICAL_FACTORIZATION 22 22 #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 23 23 #define JOB_SOLVE_ITERATIVE_REFINEMENT 33 24 #define JOB_SOLVE_FORWARD_SUBSTITUTION 331 25 #define JOB_SOLVE_DIAGONAL_SUBSTITUTION 332 26 #define JOB_SOLVE_BACKWARD_SUBSTITUTION 333 27 #define JOB_RELEASE_OF_LU_MEMORY 0 28 #define JOB_RELEASE_OF_ALL_MEMORY -1 29 30 #define IPARM_SIZE 64 31 32 #if defined(PETSC_USE_64BIT_INDICES) 33 #if defined(PETSC_HAVE_LIBMKL_INTEL_ILP64) 34 /* sizeof(MKL_INT) == sizeof(long long int) if ilp64*/ 35 #define INT_TYPE long long int 36 #define MKL_PARDISO pardiso 37 #define MKL_PARDISO_INIT pardisoinit 38 #else 39 #define INT_TYPE long long int 40 #define MKL_PARDISO pardiso_64 41 #define MKL_PARDISO_INIT pardiso_64init 42 #endif 43 #else 44 #define INT_TYPE int 45 #define MKL_PARDISO pardiso 46 #define MKL_PARDISO_INIT pardisoinit 47 #endif 48 49 50 /* 51 * Internal data structure. 52 * For more information check mkl_pardiso manual. 53 */ 54 typedef struct { 55 56 /* Configuration vector*/ 57 INT_TYPE iparm[IPARM_SIZE]; 58 59 /* 60 * Internal mkl_pardiso memory location. 61 * After the first call to mkl_pardiso do not modify pt, as that could cause a serious memory leak. 62 */ 63 void *pt[IPARM_SIZE]; 64 65 /* Basic mkl_pardiso info*/ 66 INT_TYPE phase, maxfct, mnum, mtype, n, nrhs, msglvl, err; 67 68 /* Matrix structure*/ 69 void *a; 70 INT_TYPE *ia, *ja; 71 72 /* Number of non-zero elements*/ 73 INT_TYPE nz; 74 75 /* Row permutaton vector*/ 76 INT_TYPE *perm; 77 78 /* Define if matrix preserves sparse structure.*/ 79 MatStructure matstruc; 80 81 /* True if mkl_pardiso function have been used.*/ 82 PetscBool CleanUp; 83 } Mat_MKL_PARDISO; 84 85 86 void pardiso_64init(void *pt, INT_TYPE *mtype, INT_TYPE iparm []) 87 { 88 int iparm_copy[IPARM_SIZE], mtype_copy, i; 89 90 mtype_copy = *mtype; 91 pardisoinit(pt, &mtype_copy, iparm_copy); 92 for(i = 0; i < IPARM_SIZE; i++){ 93 iparm[i] = iparm_copy[i]; 94 } 95 } 96 97 98 /* 99 * Copy the elements of matrix A. 100 * Input: 101 * - Mat A: MATSEQAIJ matrix 102 * - int shift: matrix index. 103 * - 0 for c representation 104 * - 1 for fortran representation 105 * - MatReuse reuse: 106 * - MAT_INITIAL_MATRIX: Create a new aij representation 107 * - MAT_REUSE_MATRIX: Reuse all aij representation and just change values 108 * Output: 109 * - int *nnz: Number of nonzero-elements. 110 * - int **r pointer to i index 111 * - int **c pointer to j elements 112 * - MATRIXTYPE **v: Non-zero elements 113 */ 114 #undef __FUNCT__ 115 #define __FUNCT__ "MatCopy_MKL_PARDISO" 116 PetscErrorCode MatCopy_MKL_PARDISO(Mat A, MatReuse reuse, INT_TYPE *nnz, INT_TYPE **r, INT_TYPE **c, void **v) 117 { 118 Mat_SeqAIJ *aa=(Mat_SeqAIJ*)A->data; 119 120 PetscFunctionBegin; 121 *v=aa->a; 122 if (reuse == MAT_INITIAL_MATRIX) { 123 *r = (INT_TYPE*)aa->i; 124 *c = (INT_TYPE*)aa->j; 125 *nnz = aa->nz; 126 } 127 PetscFunctionReturn(0); 128 } 129 130 /* 131 * Free memory for Mat_MKL_PARDISO structure and pointers to objects. 132 */ 133 #undef __FUNCT__ 134 #define __FUNCT__ "MatDestroy_MKL_PARDISO" 135 PetscErrorCode MatDestroy_MKL_PARDISO(Mat A) 136 { 137 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 138 PetscBool isSeqSBAIJ; 139 PetscErrorCode ierr; 140 141 PetscFunctionBegin; 142 /* Terminate instance, deallocate memories */ 143 if (mat_mkl_pardiso->CleanUp) { 144 mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY; 145 146 MKL_PARDISO (mat_mkl_pardiso->pt, 147 &mat_mkl_pardiso->maxfct, 148 &mat_mkl_pardiso->mnum, 149 &mat_mkl_pardiso->mtype, 150 &mat_mkl_pardiso->phase, 151 &mat_mkl_pardiso->n, 152 NULL, 153 NULL, 154 NULL, 155 mat_mkl_pardiso->perm, 156 &mat_mkl_pardiso->nrhs, 157 mat_mkl_pardiso->iparm, 158 &mat_mkl_pardiso->msglvl, 159 NULL, 160 NULL, 161 &mat_mkl_pardiso->err); 162 } 163 ierr = PetscFree(mat_mkl_pardiso->perm);CHKERRQ(ierr); 164 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 165 166 /* clear composed functions */ 167 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 168 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);CHKERRQ(ierr); 169 170 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 171 if (isSeqSBAIJ) {ierr = MatDestroy_SeqSBAIJ(A);CHKERRQ(ierr);} 172 else {ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);} 173 PetscFunctionReturn(0); 174 } 175 176 /* 177 * Computes Ax = b 178 */ 179 #undef __FUNCT__ 180 #define __FUNCT__ "MatSolve_MKL_PARDISO" 181 PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x) 182 { 183 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr; 184 PetscErrorCode ierr; 185 PetscScalar *xarray; 186 const PetscScalar *barray; 187 188 PetscFunctionBegin; 189 mat_mkl_pardiso->nrhs = 1; 190 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 191 ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 192 193 /* solve phase */ 194 /*-------------*/ 195 mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 196 MKL_PARDISO (mat_mkl_pardiso->pt, 197 &mat_mkl_pardiso->maxfct, 198 &mat_mkl_pardiso->mnum, 199 &mat_mkl_pardiso->mtype, 200 &mat_mkl_pardiso->phase, 201 &mat_mkl_pardiso->n, 202 mat_mkl_pardiso->a, 203 mat_mkl_pardiso->ia, 204 mat_mkl_pardiso->ja, 205 mat_mkl_pardiso->perm, 206 &mat_mkl_pardiso->nrhs, 207 mat_mkl_pardiso->iparm, 208 &mat_mkl_pardiso->msglvl, 209 (void*)barray, 210 (void*)xarray, 211 &mat_mkl_pardiso->err); 212 213 if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err); 214 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 215 ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 216 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 217 PetscFunctionReturn(0); 218 } 219 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "MatSolveTranspose_MKL_PARDISO" 223 PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x) 224 { 225 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 226 PetscErrorCode ierr; 227 228 PetscFunctionBegin; 229 #if defined(PETSC_USE_COMPLEX) 230 mat_mkl_pardiso->iparm[12 - 1] = 1; 231 #else 232 mat_mkl_pardiso->iparm[12 - 1] = 2; 233 #endif 234 ierr = MatSolve_MKL_PARDISO(A,b,x);CHKERRQ(ierr); 235 mat_mkl_pardiso->iparm[12 - 1] = 0; 236 PetscFunctionReturn(0); 237 } 238 239 240 #undef __FUNCT__ 241 #define __FUNCT__ "MatMatSolve_MKL_PARDISO" 242 PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X) 243 { 244 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr; 245 PetscErrorCode ierr; 246 PetscScalar *barray, *xarray; 247 PetscBool flg; 248 249 PetscFunctionBegin; 250 ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr); 251 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix"); 252 ierr = PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr); 253 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix"); 254 255 ierr = MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);CHKERRQ(ierr); 256 257 if(mat_mkl_pardiso->nrhs > 0){ 258 ierr = MatDenseGetArray(B,&barray); 259 ierr = MatDenseGetArray(X,&xarray); 260 261 /* solve phase */ 262 /*-------------*/ 263 mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 264 MKL_PARDISO (mat_mkl_pardiso->pt, 265 &mat_mkl_pardiso->maxfct, 266 &mat_mkl_pardiso->mnum, 267 &mat_mkl_pardiso->mtype, 268 &mat_mkl_pardiso->phase, 269 &mat_mkl_pardiso->n, 270 mat_mkl_pardiso->a, 271 mat_mkl_pardiso->ia, 272 mat_mkl_pardiso->ja, 273 mat_mkl_pardiso->perm, 274 &mat_mkl_pardiso->nrhs, 275 mat_mkl_pardiso->iparm, 276 &mat_mkl_pardiso->msglvl, 277 (void*)barray, 278 (void*)xarray, 279 &mat_mkl_pardiso->err); 280 if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err); 281 } 282 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 283 PetscFunctionReturn(0); 284 } 285 286 /* 287 * LU Decomposition 288 */ 289 #undef __FUNCT__ 290 #define __FUNCT__ "MatFactorNumeric_MKL_PARDISO" 291 PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info) 292 { 293 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->spptr; 294 PetscErrorCode ierr; 295 296 /* numerical factorization phase */ 297 /*-------------------------------*/ 298 PetscFunctionBegin; 299 mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN; 300 ierr = MatCopy_MKL_PARDISO(A, MAT_REUSE_MATRIX, &mat_mkl_pardiso->nz, &mat_mkl_pardiso->ia, &mat_mkl_pardiso->ja, &mat_mkl_pardiso->a);CHKERRQ(ierr); 301 302 /* numerical factorization phase */ 303 /*-------------------------------*/ 304 mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION; 305 MKL_PARDISO (mat_mkl_pardiso->pt, 306 &mat_mkl_pardiso->maxfct, 307 &mat_mkl_pardiso->mnum, 308 &mat_mkl_pardiso->mtype, 309 &mat_mkl_pardiso->phase, 310 &mat_mkl_pardiso->n, 311 mat_mkl_pardiso->a, 312 mat_mkl_pardiso->ia, 313 mat_mkl_pardiso->ja, 314 mat_mkl_pardiso->perm, 315 &mat_mkl_pardiso->nrhs, 316 mat_mkl_pardiso->iparm, 317 &mat_mkl_pardiso->msglvl, 318 NULL, 319 NULL, 320 &mat_mkl_pardiso->err); 321 if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d. Please check manual\n",mat_mkl_pardiso->err); 322 323 mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN; 324 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 325 PetscFunctionReturn(0); 326 } 327 328 /* Sets mkl_pardiso options from the options database */ 329 #undef __FUNCT__ 330 #define __FUNCT__ "PetscSetMKL_PARDISOFromOptions" 331 PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A) 332 { 333 Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr; 334 PetscErrorCode ierr; 335 PetscInt icntl, threads = 1; 336 PetscBool flg; 337 338 PetscFunctionBegin; 339 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");CHKERRQ(ierr); 340 ierr = PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use","None",threads,&threads,&flg);CHKERRQ(ierr); 341 if (flg) mkl_set_num_threads((int)threads); 342 343 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); 344 if (flg) mat_mkl_pardiso->maxfct = icntl; 345 346 ierr = PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);CHKERRQ(ierr); 347 if (flg) mat_mkl_pardiso->mnum = icntl; 348 349 ierr = PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);CHKERRQ(ierr); 350 if (flg) mat_mkl_pardiso->msglvl = icntl; 351 352 ierr = PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);CHKERRQ(ierr); 353 if(flg){ 354 void *pt[IPARM_SIZE]; 355 mat_mkl_pardiso->mtype = icntl; 356 MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm); 357 #if defined(PETSC_USE_REAL_SINGLE) 358 mat_mkl_pardiso->iparm[27] = 1; 359 #else 360 mat_mkl_pardiso->iparm[27] = 0; 361 #endif 362 mat_mkl_pardiso->iparm[34] = 1; 363 } 364 ierr = PetscOptionsInt("-mat_mkl_pardiso_1","Use default values","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);CHKERRQ(ierr); 365 366 if(flg && icntl != 0){ 367 ierr = PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);CHKERRQ(ierr); 368 if (flg) mat_mkl_pardiso->iparm[1] = icntl; 369 370 ierr = PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);CHKERRQ(ierr); 371 if (flg) mat_mkl_pardiso->iparm[3] = icntl; 372 373 ierr = PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);CHKERRQ(ierr); 374 if (flg) mat_mkl_pardiso->iparm[4] = icntl; 375 376 ierr = PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);CHKERRQ(ierr); 377 if (flg) mat_mkl_pardiso->iparm[5] = icntl; 378 379 ierr = PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);CHKERRQ(ierr); 380 if (flg) mat_mkl_pardiso->iparm[7] = icntl; 381 382 ierr = PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);CHKERRQ(ierr); 383 if (flg) mat_mkl_pardiso->iparm[9] = icntl; 384 385 ierr = PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);CHKERRQ(ierr); 386 if (flg) mat_mkl_pardiso->iparm[10] = icntl; 387 388 ierr = PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);CHKERRQ(ierr); 389 if (flg) mat_mkl_pardiso->iparm[11] = icntl; 390 391 ierr = PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);CHKERRQ(ierr); 392 if (flg) mat_mkl_pardiso->iparm[12] = icntl; 393 394 ierr = PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);CHKERRQ(ierr); 395 if (flg) mat_mkl_pardiso->iparm[17] = icntl; 396 397 ierr = PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);CHKERRQ(ierr); 398 if (flg) mat_mkl_pardiso->iparm[18] = icntl; 399 400 ierr = PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);CHKERRQ(ierr); 401 if (flg) mat_mkl_pardiso->iparm[20] = icntl; 402 403 ierr = PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);CHKERRQ(ierr); 404 if (flg) mat_mkl_pardiso->iparm[23] = icntl; 405 406 ierr = PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);CHKERRQ(ierr); 407 if (flg) mat_mkl_pardiso->iparm[24] = icntl; 408 409 ierr = PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);CHKERRQ(ierr); 410 if (flg) mat_mkl_pardiso->iparm[26] = icntl; 411 412 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); 413 if (flg) mat_mkl_pardiso->iparm[30] = icntl; 414 415 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); 416 if (flg) mat_mkl_pardiso->iparm[33] = icntl; 417 418 ierr = PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);CHKERRQ(ierr); 419 if (flg) mat_mkl_pardiso->iparm[59] = icntl; 420 } 421 PetscOptionsEnd(); 422 PetscFunctionReturn(0); 423 } 424 425 #undef __FUNCT__ 426 #define __FUNCT__ "MatFactorMKL_PARDISOInitialize_Private" 427 PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso) 428 { 429 PetscErrorCode ierr; 430 PetscInt i; 431 432 PetscFunctionBegin; 433 for ( i = 0; i < IPARM_SIZE; i++ ){ 434 mat_mkl_pardiso->iparm[i] = 0; 435 } 436 437 for ( i = 0; i < IPARM_SIZE; i++ ){ 438 mat_mkl_pardiso->pt[i] = 0; 439 } 440 441 /* Default options for both sym and unsym */ 442 mat_mkl_pardiso->iparm[ 0] = 1; /* Solver default parameters overriden with provided by iparm */ 443 mat_mkl_pardiso->iparm[ 1] = 2; /* Metis reordering */ 444 mat_mkl_pardiso->iparm[ 5] = 0; /* Write solution into x */ 445 mat_mkl_pardiso->iparm[ 7] = 2; /* Max number of iterative refinement steps */ 446 mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */ 447 mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */ 448 #if 0 449 mat_mkl_pardiso->iparm[23] = 1; /* Parallel factorization control*/ 450 #endif 451 mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */ 452 mat_mkl_pardiso->iparm[39] = 0; /* Input: matrix/rhs/solution stored on master */ 453 454 mat_mkl_pardiso->CleanUp = PETSC_FALSE; 455 mat_mkl_pardiso->maxfct = 1; /* Maximum number of numerical factorizations. */ 456 mat_mkl_pardiso->mnum = 1; /* Which factorization to use. */ 457 mat_mkl_pardiso->msglvl = 0; /* 0: do not print 1: Print statistical information in file */ 458 mat_mkl_pardiso->phase = -1; 459 mat_mkl_pardiso->err = 0; 460 461 mat_mkl_pardiso->n = A->rmap->N; 462 mat_mkl_pardiso->nrhs = 1; 463 mat_mkl_pardiso->err = 0; 464 mat_mkl_pardiso->phase = -1; 465 466 if(ftype == MAT_FACTOR_LU){ 467 /* Default type for non-sym */ 468 #if defined(PETSC_USE_COMPLEX) 469 mat_mkl_pardiso->mtype = 13; 470 #else 471 mat_mkl_pardiso->mtype = 11; 472 #endif 473 474 mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */ 475 mat_mkl_pardiso->iparm[10] = 1; /* Use nonsymmetric permutation and scaling MPS */ 476 mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 477 478 } else { 479 /* Default type for sym */ 480 #if defined(PETSC_USE_COMPLEX) 481 mat_mkl_pardiso ->mtype = 3; 482 #else 483 mat_mkl_pardiso ->mtype = -2; 484 #endif 485 mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */ 486 mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */ 487 mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 488 /* mat_mkl_pardiso->iparm[20] = 1; */ /* Apply 1x1 and 2x2 Bunch-Kaufman pivoting during the factorization process */ 489 #if defined(PETSC_USE_DEBUG) 490 mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */ 491 #endif 492 } 493 ierr = PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);CHKERRQ(ierr); 494 for(i = 0; i < A->rmap->N; i++){ 495 mat_mkl_pardiso->perm[i] = 0; 496 } 497 PetscFunctionReturn(0); 498 } 499 500 /* 501 * Symbolic decomposition. Mkl_Pardiso analysis phase. 502 */ 503 #undef __FUNCT__ 504 #define __FUNCT__ "MatFactorSymbolic_AIJMKL_PARDISO_Private" 505 PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info) 506 { 507 Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr; 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN; 512 513 /* Set MKL_PARDISO options from the options database */ 514 ierr = PetscSetMKL_PARDISOFromOptions(F,A);CHKERRQ(ierr); 515 516 ierr = MatCopy_MKL_PARDISO(A, MAT_INITIAL_MATRIX, &mat_mkl_pardiso->nz, &mat_mkl_pardiso->ia, &mat_mkl_pardiso->ja, &mat_mkl_pardiso->a);CHKERRQ(ierr); 517 mat_mkl_pardiso->n = A->rmap->N; 518 519 /* analysis phase */ 520 /*----------------*/ 521 mat_mkl_pardiso->phase = JOB_ANALYSIS; 522 523 MKL_PARDISO (mat_mkl_pardiso->pt, 524 &mat_mkl_pardiso->maxfct, 525 &mat_mkl_pardiso->mnum, 526 &mat_mkl_pardiso->mtype, 527 &mat_mkl_pardiso->phase, 528 &mat_mkl_pardiso->n, 529 mat_mkl_pardiso->a, 530 mat_mkl_pardiso->ia, 531 mat_mkl_pardiso->ja, 532 mat_mkl_pardiso->perm, 533 &mat_mkl_pardiso->nrhs, 534 mat_mkl_pardiso->iparm, 535 &mat_mkl_pardiso->msglvl, 536 NULL, 537 NULL, 538 &mat_mkl_pardiso->err); 539 if (mat_mkl_pardiso->err < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error reported by MKL_PARDISO: err=%d\n. Please check manual",mat_mkl_pardiso->err); 540 541 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 542 543 if(F->factortype == MAT_FACTOR_LU){ 544 F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO; 545 } else { 546 F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO; 547 } 548 F->ops->solve = MatSolve_MKL_PARDISO; 549 F->ops->solvetranspose = MatSolveTranspose_MKL_PARDISO; 550 F->ops->matsolve = MatMatSolve_MKL_PARDISO; 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "MatLUFactorSymbolic_AIJMKL_PARDISO" 556 PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 557 { 558 PetscErrorCode ierr; 559 560 PetscFunctionBegin; 561 ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr); 562 PetscFunctionReturn(0); 563 } 564 565 #undef __FUNCT__ 566 #define __FUNCT__ "MatCholeskyFactorSymbolic_AIJMKL_PARDISO" 567 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info) 568 { 569 PetscErrorCode ierr; 570 571 PetscFunctionBegin; 572 ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr); 573 PetscFunctionReturn(0); 574 } 575 576 #undef __FUNCT__ 577 #define __FUNCT__ "MatView_MKL_PARDISO" 578 PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer) 579 { 580 PetscErrorCode ierr; 581 PetscBool iascii; 582 PetscViewerFormat format; 583 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 584 PetscInt i; 585 586 PetscFunctionBegin; 587 /* check if matrix is mkl_pardiso type */ 588 if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(0); 589 590 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 591 if (iascii) { 592 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 593 if (format == PETSC_VIEWER_ASCII_INFO) { 594 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");CHKERRQ(ierr); 595 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase: %d \n",mat_mkl_pardiso->phase);CHKERRQ(ierr); 596 for(i = 1; i <= 64; i++){ 597 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]: %d \n",i, mat_mkl_pardiso->iparm[i - 1]);CHKERRQ(ierr); 598 } 599 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct: %d \n", mat_mkl_pardiso->maxfct);CHKERRQ(ierr); 600 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum: %d \n", mat_mkl_pardiso->mnum);CHKERRQ(ierr); 601 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype: %d \n", mat_mkl_pardiso->mtype);CHKERRQ(ierr); 602 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n: %d \n", mat_mkl_pardiso->n);CHKERRQ(ierr); 603 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs: %d \n", mat_mkl_pardiso->nrhs);CHKERRQ(ierr); 604 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl: %d \n", mat_mkl_pardiso->msglvl);CHKERRQ(ierr); 605 } 606 } 607 PetscFunctionReturn(0); 608 } 609 610 611 #undef __FUNCT__ 612 #define __FUNCT__ "MatGetInfo_MKL_PARDISO" 613 PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info) 614 { 615 Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)A->spptr; 616 617 PetscFunctionBegin; 618 info->block_size = 1.0; 619 info->nz_allocated = mat_mkl_pardiso->nz + 0.0; 620 info->nz_unneeded = 0.0; 621 info->assemblies = 0.0; 622 info->mallocs = 0.0; 623 info->memory = 0.0; 624 info->fill_ratio_given = 0; 625 info->fill_ratio_needed = 0; 626 info->factor_mallocs = 0; 627 PetscFunctionReturn(0); 628 } 629 630 #undef __FUNCT__ 631 #define __FUNCT__ "MatMkl_PardisoSetCntl_MKL_PARDISO" 632 PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival) 633 { 634 Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)F->spptr; 635 636 PetscFunctionBegin; 637 if(icntl <= 64){ 638 mat_mkl_pardiso->iparm[icntl - 1] = ival; 639 } else { 640 if(icntl == 65) 641 mkl_set_num_threads((int)ival); 642 else if(icntl == 66) 643 mat_mkl_pardiso->maxfct = ival; 644 else if(icntl == 67) 645 mat_mkl_pardiso->mnum = ival; 646 else if(icntl == 68) 647 mat_mkl_pardiso->msglvl = ival; 648 else if(icntl == 69){ 649 void *pt[IPARM_SIZE]; 650 mat_mkl_pardiso->mtype = ival; 651 MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm); 652 #if defined(PETSC_USE_REAL_SINGLE) 653 mat_mkl_pardiso->iparm[27] = 1; 654 #else 655 mat_mkl_pardiso->iparm[27] = 0; 656 #endif 657 mat_mkl_pardiso->iparm[34] = 1; 658 } 659 } 660 PetscFunctionReturn(0); 661 } 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "MatMkl_PardisoSetCntl" 665 /*@ 666 MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters 667 668 Logically Collective on Mat 669 670 Input Parameters: 671 + F - the factored matrix obtained by calling MatGetFactor() 672 . icntl - index of Mkl_Pardiso parameter 673 - ival - value of Mkl_Pardiso parameter 674 675 Options Database: 676 . -mat_mkl_pardiso_<icntl> <ival> 677 678 Level: beginner 679 680 References: Mkl_Pardiso Users' Guide 681 682 .seealso: MatGetFactor() 683 @*/ 684 PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival) 685 { 686 PetscErrorCode ierr; 687 688 PetscFunctionBegin; 689 ierr = PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 690 PetscFunctionReturn(0); 691 } 692 693 /*MC 694 MATSOLVERMKL_PARDISO - A matrix type providing direct solvers (LU) for 695 sequential matrices via the external package MKL_PARDISO. 696 697 Works with MATSEQAIJ matrices 698 699 Use -pc_type lu -pc_factor_mat_solver_package mkl_pardiso to us this direct solver 700 701 Options Database Keys: 702 + -mat_mkl_pardiso_65 - Number of threads to use 703 . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time 704 . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase 705 . -mat_mkl_pardiso_68 - Message level information 706 . -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 707 . -mat_mkl_pardiso_1 - Use default values 708 . -mat_mkl_pardiso_2 - Fill-in reducing ordering for the input matrix 709 . -mat_mkl_pardiso_4 - Preconditioned CGS/CG 710 . -mat_mkl_pardiso_5 - User permutation 711 . -mat_mkl_pardiso_6 - Write solution on x 712 . -mat_mkl_pardiso_8 - Iterative refinement step 713 . -mat_mkl_pardiso_10 - Pivoting perturbation 714 . -mat_mkl_pardiso_11 - Scaling vectors 715 . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A 716 . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching 717 . -mat_mkl_pardiso_18 - Numbers of non-zero elements 718 . -mat_mkl_pardiso_19 - Report number of floating point operations 719 . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices 720 . -mat_mkl_pardiso_24 - Parallel factorization control 721 . -mat_mkl_pardiso_25 - Parallel forward/backward solve control 722 . -mat_mkl_pardiso_27 - Matrix checker 723 . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors 724 . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode 725 - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode 726 727 Level: beginner 728 729 For more information please check mkl_pardiso manual 730 731 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 732 733 M*/ 734 #undef __FUNCT__ 735 #define __FUNCT__ "MatFactorGetSolverPackage_mkl_pardiso" 736 static PetscErrorCode MatFactorGetSolverPackage_mkl_pardiso(Mat A, const MatSolverPackage *type) 737 { 738 PetscFunctionBegin; 739 *type = MATSOLVERMKL_PARDISO; 740 PetscFunctionReturn(0); 741 } 742 743 /* MatGetFactor for Seq sbAIJ matrices */ 744 #undef __FUNCT__ 745 #define __FUNCT__ "MatGetFactor_sbaij_mkl_pardiso" 746 PETSC_EXTERN PetscErrorCode MatGetFactor_sbaij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F) 747 { 748 Mat B; 749 PetscErrorCode ierr; 750 Mat_MKL_PARDISO *mat_mkl_pardiso; 751 PetscBool isSeqSBAIJ; 752 PetscInt bs; 753 754 PetscFunctionBegin; 755 /* Create the factorization matrix */ 756 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 757 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 758 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 759 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 760 ierr = MatSeqSBAIJSetPreallocation(B,1,0,NULL);CHKERRQ(ierr); 761 ierr = MatGetBlockSize(A,&bs); CHKERRQ(ierr); 762 763 if(bs != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrice MATSEQSBAIJ with block size other than 1 is not supported by Pardiso"); 764 if(ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrice MATSEQAIJ should be used only with MAT_FACTOR_CHOLESKY."); 765 766 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO; 767 B->factortype = MAT_FACTOR_CHOLESKY; 768 B->ops->destroy = MatDestroy_MKL_PARDISO; 769 B->ops->view = MatView_MKL_PARDISO; 770 B->factortype = ftype; 771 B->ops->getinfo = MatGetInfo_MKL_PARDISO; 772 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 773 774 ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr); 775 B->spptr = mat_mkl_pardiso; 776 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr); 777 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr); 778 ierr = MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr); 779 *F = B; 780 PetscFunctionReturn(0); 781 } 782 783 /* MatGetFactor for Seq AIJ matrices */ 784 #undef __FUNCT__ 785 #define __FUNCT__ "MatGetFactor_aij_mkl_pardiso" 786 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F) 787 { 788 Mat B; 789 PetscErrorCode ierr; 790 Mat_MKL_PARDISO *mat_mkl_pardiso; 791 PetscBool isSeqAIJ; 792 793 PetscFunctionBegin; 794 /* Create the factorization matrix */ 795 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 796 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 797 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 798 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 799 ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 800 801 if(ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrice MATSEQAIJ should be used only with MAT_FACTOR_LU."); 802 803 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO; 804 B->factortype = MAT_FACTOR_LU; 805 B->ops->destroy = MatDestroy_MKL_PARDISO; 806 B->ops->view = MatView_MKL_PARDISO; 807 B->factortype = ftype; 808 B->ops->getinfo = MatGetInfo_MKL_PARDISO; 809 B->assembled = PETSC_TRUE; /* required by -ksp_view */ 810 811 ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr); 812 B->spptr = mat_mkl_pardiso; 813 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr); 814 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr); 815 ierr = MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr); 816 817 *F = B; 818 PetscFunctionReturn(0); 819 } 820 821 #undef __FUNCT__ 822 #define __FUNCT__ "MatSolverPackageRegister_MKL_Pardiso" 823 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MKL_Pardiso(void) 824 { 825 PetscErrorCode ierr; 826 827 PetscFunctionBegin; 828 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_aij_mkl_pardiso );CHKERRQ(ierr); 829 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ, MAT_FACTOR_CHOLESKY,MatGetFactor_sbaij_mkl_pardiso);CHKERRQ(ierr); 830 PetscFunctionReturn(0); 831 } 832 833