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> 7 #include <../src/mat/impls/dense/seq/dense.h> 8 #include <petscblaslapack.h> 9 10 #include <stdio.h> 11 #include <stdlib.h> 12 #include <math.h> 13 #include <mkl.h> 14 15 /* 16 * Possible mkl_pardiso phases that controls the execution of the solver. 17 * For more information check mkl_pardiso manual. 18 */ 19 #define JOB_ANALYSIS 11 20 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION 12 21 #define JOB_ANALYSIS_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 13 22 #define JOB_NUMERICAL_FACTORIZATION 22 23 #define JOB_NUMERICAL_FACTORIZATION_SOLVE_ITERATIVE_REFINEMENT 23 24 #define JOB_SOLVE_ITERATIVE_REFINEMENT 33 25 #define JOB_SOLVE_FORWARD_SUBSTITUTION 331 26 #define JOB_SOLVE_DIAGONAL_SUBSTITUTION 332 27 #define JOB_SOLVE_BACKWARD_SUBSTITUTION 333 28 #define JOB_RELEASE_OF_LU_MEMORY 0 29 #define JOB_RELEASE_OF_ALL_MEMORY -1 30 31 #define IPARM_SIZE 64 32 33 #if defined(PETSC_USE_64BIT_INDICES) 34 #if defined(PETSC_HAVE_LIBMKL_INTEL_ILP64) 35 /* sizeof(MKL_INT) == sizeof(long long int) if ilp64*/ 36 #define INT_TYPE long long int 37 #define MKL_PARDISO pardiso 38 #define MKL_PARDISO_INIT pardisoinit 39 #else 40 #define INT_TYPE long long int 41 #define MKL_PARDISO pardiso_64 42 #define MKL_PARDISO_INIT pardiso_64init 43 #endif 44 #else 45 #define INT_TYPE int 46 #define MKL_PARDISO pardiso 47 #define MKL_PARDISO_INIT pardisoinit 48 #endif 49 50 51 /* 52 * Internal data structure. 53 * For more information check mkl_pardiso manual. 54 */ 55 typedef struct { 56 57 /* Configuration vector*/ 58 INT_TYPE iparm[IPARM_SIZE]; 59 60 /* 61 * Internal mkl_pardiso memory location. 62 * After the first call to mkl_pardiso do not modify pt, as that could cause a serious memory leak. 63 */ 64 void *pt[IPARM_SIZE]; 65 66 /* Basic mkl_pardiso info*/ 67 INT_TYPE phase, maxfct, mnum, mtype, n, nrhs, msglvl, err; 68 69 /* Matrix structure*/ 70 void *a; 71 INT_TYPE *ia, *ja; 72 73 /* Number of non-zero elements*/ 74 INT_TYPE nz; 75 76 /* Row permutaton vector*/ 77 INT_TYPE *perm; 78 79 /* Define if matrix preserves sparse structure.*/ 80 MatStructure matstruc; 81 82 PetscBool needsym; 83 PetscBool freeaij; 84 85 /* Schur complement */ 86 PetscScalar *schur; 87 PetscInt schur_size; 88 PetscInt *schur_idxs; 89 PetscScalar *schur_work; 90 PetscBLASInt schur_work_size; 91 PetscInt schur_solver_type; 92 PetscInt *schur_pivots; 93 PetscBool schur_factored; 94 PetscBool schur_inverted; 95 PetscBool solve_interior; 96 97 /* True if mkl_pardiso function have been used.*/ 98 PetscBool CleanUp; 99 100 /* Conversion to a format suitable for MKL */ 101 PetscErrorCode (*Convert)(Mat, PetscBool, MatReuse, PetscBool*, INT_TYPE*, INT_TYPE**, INT_TYPE**, void**); 102 PetscErrorCode (*Destroy)(Mat); 103 } Mat_MKL_PARDISO; 104 105 #undef __FUNCT__ 106 #define __FUNCT__ "MatMKLPardiso_Convert_seqsbaij" 107 PetscErrorCode MatMKLPardiso_Convert_seqsbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,void **v) 108 { 109 Mat_SeqSBAIJ *aa = (Mat_SeqSBAIJ*)A->data; 110 PetscInt bs = A->rmap->bs; 111 112 PetscFunctionBegin; 113 if (!sym) { 114 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"This should not happen"); 115 } 116 if (bs == 1) { /* already in the correct format */ 117 *v = aa->a; 118 *r = aa->i; 119 *c = aa->j; 120 *nnz = aa->nz; 121 *free = PETSC_FALSE; 122 } else { 123 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Conversion from SeqSBAIJ to MKL Pardiso format still need to be implemented"); 124 } 125 PetscFunctionReturn(0); 126 } 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "MatMKLPardiso_Convert_seqbaij" 130 PetscErrorCode MatMKLPardiso_Convert_seqbaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,void **v) 131 { 132 Mat_SeqBAIJ *aa=(Mat_SeqBAIJ*)A->data; 133 PetscFunctionBegin; 134 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Conversion from SeqBAIJ to MKL Pardiso format still need to be implemented"); 135 PetscFunctionReturn(0); 136 } 137 138 #undef __FUNCT__ 139 #define __FUNCT__ "MatMKLPardiso_Convert_seqaij" 140 PetscErrorCode MatMKLPardiso_Convert_seqaij(Mat A,PetscBool sym,MatReuse reuse,PetscBool *free,INT_TYPE *nnz,INT_TYPE **r,INT_TYPE **c,void **v) 141 { 142 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 143 PetscErrorCode ierr; 144 145 PetscFunctionBegin; 146 if (!sym) { /* already in the correct format */ 147 *v = aa->a; 148 *r = aa->i; 149 *c = aa->j; 150 *nnz = aa->nz; 151 *free = PETSC_FALSE; 152 PetscFunctionReturn(0); 153 } 154 /* need to get the triangular part */ 155 if (reuse == MAT_INITIAL_MATRIX) { 156 PetscScalar *vals,*vv; 157 PetscInt *row,*col,*jj; 158 PetscInt m = A->rmap->n,nz,i; 159 160 nz = 0; 161 for (i=0; i<m; i++) { 162 nz += aa->i[i+1] - aa->diag[i]; 163 } 164 ierr = PetscMalloc2(m+1,&row,nz,&col);CHKERRQ(ierr); 165 ierr = PetscMalloc1(nz,&vals);CHKERRQ(ierr); 166 jj = col; 167 vv = vals; 168 169 row[0] = 0; 170 for (i=0; i<m; i++) { 171 PetscInt *aj = aa->j + aa->diag[i]; 172 PetscScalar *av = aa->a + aa->diag[i]; 173 PetscInt rl = aa->i[i+1] - aa->diag[i],j; 174 for (j=0; j<rl; j++) { 175 *jj = *aj; jj++; aj++; 176 *vv = *av; vv++; av++; 177 } 178 row[i+1] = row[i] + rl; 179 } 180 *v = vals; 181 *r = row; 182 *c = col; 183 *nnz = nz; 184 *free = PETSC_TRUE; 185 } else { 186 PetscScalar *vv; 187 PetscInt m = A->rmap->n,i; 188 189 vv = *v; 190 for (i=0; i<m; i++) { 191 PetscInt *aj = aa->j + aa->diag[i]; 192 PetscScalar *av = aa->a + aa->diag[i]; 193 PetscInt rl = aa->i[i+1] - aa->diag[i],j; 194 for (j=0; j<rl; j++) { 195 *vv = *av; vv++; av++; 196 } 197 } 198 *free = PETSC_TRUE; 199 } 200 PetscFunctionReturn(0); 201 } 202 203 void pardiso_64init(void *pt, INT_TYPE *mtype, INT_TYPE iparm []) 204 { 205 int iparm_copy[IPARM_SIZE], mtype_copy, i; 206 207 mtype_copy = *mtype; 208 pardisoinit(pt, &mtype_copy, iparm_copy); 209 for(i = 0; i < IPARM_SIZE; i++){ 210 iparm[i] = iparm_copy[i]; 211 } 212 } 213 214 #undef __FUNCT__ 215 #define __FUNCT__ "MatMKLPardisoFactorSchur_Private" 216 static PetscErrorCode MatMKLPardisoFactorSchur_Private(Mat_MKL_PARDISO* mpardiso) 217 { 218 PetscBLASInt B_N,B_ierr; 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 if (mpardiso->schur_factored) { 223 PetscFunctionReturn(0); 224 } 225 ierr = PetscBLASIntCast(mpardiso->schur_size,&B_N);CHKERRQ(ierr); 226 switch (mpardiso->schur_solver_type) { 227 case 1: /* symmetric */ 228 if (!mpardiso->schur_pivots) { 229 ierr = PetscMalloc1(B_N,&mpardiso->schur_pivots);CHKERRQ(ierr); 230 } 231 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 232 PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,mpardiso->schur_work,&mpardiso->schur_work_size,&B_ierr)); 233 ierr = PetscFPTrapPop();CHKERRQ(ierr); 234 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRF Lapack routine %d",(int)B_ierr); 235 break; 236 case 2: /* hermitian solver */ 237 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 238 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&B_N,mpardiso->schur,&B_N,&B_ierr)); 239 ierr = PetscFPTrapPop();CHKERRQ(ierr); 240 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRF Lapack routine %d",(int)B_ierr); 241 break; 242 default: /* general */ 243 if (!mpardiso->schur_pivots) { 244 ierr = PetscMalloc1(B_N,&mpardiso->schur_pivots);CHKERRQ(ierr); 245 } 246 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 247 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&B_N,&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,&B_ierr)); 248 ierr = PetscFPTrapPop();CHKERRQ(ierr); 249 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRF Lapack routine %d",(int)B_ierr); 250 break; 251 } 252 mpardiso->schur_factored = PETSC_TRUE; 253 PetscFunctionReturn(0); 254 } 255 256 #undef __FUNCT__ 257 #define __FUNCT__ "MatMKLPardisoInvertSchur_Private" 258 static PetscErrorCode MatMKLPardisoInvertSchur_Private(Mat_MKL_PARDISO* mpardiso) 259 { 260 PetscBLASInt B_N,B_ierr; 261 PetscErrorCode ierr; 262 263 PetscFunctionBegin; 264 ierr = MatMKLPardisoFactorSchur_Private(mpardiso);CHKERRQ(ierr); 265 ierr = PetscBLASIntCast(mpardiso->schur_size,&B_N);CHKERRQ(ierr); 266 switch (mpardiso->schur_solver_type) { 267 case 1: /* symmetric */ 268 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 269 PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,mpardiso->schur_work,&B_ierr)); 270 ierr = PetscFPTrapPop();CHKERRQ(ierr); 271 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRI Lapack routine %d",(int)B_ierr); 272 break; 273 case 2: /* hermitian solver */ 274 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 275 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&B_N,mpardiso->schur,&B_N,&B_ierr)); 276 ierr = PetscFPTrapPop();CHKERRQ(ierr); 277 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRI Lapack routine %d",(int)B_ierr); 278 break; 279 default: /* general */ 280 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 281 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&B_N,mpardiso->schur,&B_N,mpardiso->schur_pivots,mpardiso->schur_work,&mpardiso->schur_work_size,&B_ierr)); 282 ierr = PetscFPTrapPop();CHKERRQ(ierr); 283 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRI Lapack routine %d",(int)B_ierr); 284 break; 285 } 286 mpardiso->schur_inverted = PETSC_TRUE; 287 PetscFunctionReturn(0); 288 } 289 290 #undef __FUNCT__ 291 #define __FUNCT__ "MatMKLPardisoSolveSchur_Private" 292 static PetscErrorCode MatMKLPardisoSolveSchur_Private(Mat_MKL_PARDISO* mpardiso, PetscScalar *B, PetscScalar *X) 293 { 294 PetscScalar one=1.,zero=0.,*schur_rhs,*schur_sol; 295 PetscBLASInt B_N,B_Nrhs,B_ierr; 296 char type[2]; 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 ierr = MatMKLPardisoFactorSchur_Private(mpardiso);CHKERRQ(ierr); 301 ierr = PetscBLASIntCast(mpardiso->schur_size,&B_N);CHKERRQ(ierr); 302 ierr = PetscBLASIntCast(mpardiso->nrhs,&B_Nrhs);CHKERRQ(ierr); 303 if (X == B && mpardiso->schur_inverted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"X and B cannot point to the same address"); 304 if (X != B) { /* using LAPACK *TRS subroutines */ 305 ierr = PetscMemcpy(X,B,B_N*B_Nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 306 } 307 schur_rhs = B; 308 schur_sol = X; 309 switch (mpardiso->schur_solver_type) { 310 case 1: /* symmetric solver */ 311 if (mpardiso->schur_inverted) { 312 PetscStackCallBLAS("BLASsymm",BLASsymm_("L","L",&B_N,&B_Nrhs,&one,mpardiso->schur,&B_N,schur_rhs,&B_N,&zero,schur_sol,&B_N)); 313 } else { 314 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 315 PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&B_N,&B_Nrhs,mpardiso->schur,&B_N,mpardiso->schur_pivots,schur_sol,&B_N,&B_ierr)); 316 ierr = PetscFPTrapPop();CHKERRQ(ierr); 317 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in SYTRS Lapack routine %d",(int)B_ierr); 318 } 319 break; 320 case 2: /* hermitian solver */ 321 if (mpardiso->schur_inverted) { /* BLAShemm should go here */ 322 PetscStackCallBLAS("BLASsymm",BLASsymm_("L","L",&B_N,&B_Nrhs,&one,mpardiso->schur,&B_N,schur_rhs,&B_N,&zero,schur_sol,&B_N)); 323 } else { 324 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 325 PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&B_N,&B_Nrhs,mpardiso->schur,&B_N,schur_sol,&B_N,&B_ierr)); 326 ierr = PetscFPTrapPop();CHKERRQ(ierr); 327 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in POTRS Lapack routine %d",(int)B_ierr); 328 } 329 break; 330 default: /* general */ 331 switch (mpardiso->iparm[12-1]) { 332 case 1: 333 sprintf(type,"C"); 334 break; 335 case 2: 336 sprintf(type,"T"); 337 break; 338 default: 339 sprintf(type,"N"); 340 break; 341 } 342 if (mpardiso->schur_inverted) { 343 PetscStackCallBLAS("BLASgemm",BLASgemm_(type,"N",&B_N,&B_Nrhs,&B_N,&one,mpardiso->schur,&B_N,schur_rhs,&B_N,&zero,schur_sol,&B_N)); 344 } else { 345 ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 346 PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_(type,&B_N,&B_Nrhs,mpardiso->schur,&B_N,mpardiso->schur_pivots,schur_sol,&B_N,&B_ierr)); 347 ierr = PetscFPTrapPop();CHKERRQ(ierr); 348 if (B_ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in GETRS Lapack routine %d",(int)B_ierr); 349 } 350 break; 351 } 352 PetscFunctionReturn(0); 353 } 354 355 356 #undef __FUNCT__ 357 #define __FUNCT__ "MatFactorSetSchurIS_MKL_PARDISO" 358 PetscErrorCode MatFactorSetSchurIS_MKL_PARDISO(Mat F, IS is) 359 { 360 Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr; 361 const PetscInt *idxs; 362 PetscInt size,i; 363 PetscMPIInt csize; 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)F),&csize);CHKERRQ(ierr); 368 if (csize > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MKL_PARDISO parallel Schur complements not yet supported from PETSc\n"); 369 ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 370 if (mpardiso->schur_size != size) { 371 mpardiso->schur_size = size; 372 ierr = PetscFree2(mpardiso->schur,mpardiso->schur_work);CHKERRQ(ierr); 373 ierr = PetscFree(mpardiso->schur_idxs);CHKERRQ(ierr); 374 ierr = PetscFree(mpardiso->schur_pivots);CHKERRQ(ierr); 375 ierr = PetscBLASIntCast(PetscMax(mpardiso->n,2*size),&mpardiso->schur_work_size);CHKERRQ(ierr); 376 ierr = PetscMalloc2(size*size,&mpardiso->schur,mpardiso->schur_work_size,&mpardiso->schur_work);CHKERRQ(ierr); 377 ierr = PetscMalloc1(size,&mpardiso->schur_idxs);CHKERRQ(ierr); 378 } 379 ierr = PetscMemzero(mpardiso->perm,mpardiso->n*sizeof(INT_TYPE));CHKERRQ(ierr); 380 ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 381 ierr = PetscMemcpy(mpardiso->schur_idxs,idxs,size*sizeof(PetscInt));CHKERRQ(ierr); 382 for (i=0;i<size;i++) mpardiso->perm[idxs[i]] = 1; 383 ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 384 if (size) { /* turn on Schur switch if we the set of indices is not empty */ 385 mpardiso->iparm[36-1] = 2; 386 } 387 mpardiso->schur_factored = PETSC_FALSE; 388 mpardiso->schur_inverted = PETSC_FALSE; 389 PetscFunctionReturn(0); 390 } 391 392 #undef __FUNCT__ 393 #define __FUNCT__ "MatFactorCreateSchurComplement_MKL_PARDISO" 394 PetscErrorCode MatFactorCreateSchurComplement_MKL_PARDISO(Mat F,Mat* S) 395 { 396 Mat St; 397 Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr; 398 PetscScalar *array; 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 403 else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 404 405 ierr = MatCreate(PetscObjectComm((PetscObject)F),&St);CHKERRQ(ierr); 406 ierr = MatSetSizes(St,PETSC_DECIDE,PETSC_DECIDE,mpardiso->schur_size,mpardiso->schur_size);CHKERRQ(ierr); 407 ierr = MatSetType(St,MATDENSE);CHKERRQ(ierr); 408 ierr = MatSetUp(St);CHKERRQ(ierr); 409 ierr = MatDenseGetArray(St,&array);CHKERRQ(ierr); 410 ierr = PetscMemcpy(array,mpardiso->schur,mpardiso->schur_size*mpardiso->schur_size*sizeof(PetscScalar));CHKERRQ(ierr); 411 ierr = MatDenseRestoreArray(St,&array);CHKERRQ(ierr); 412 *S = St; 413 PetscFunctionReturn(0); 414 } 415 416 #undef __FUNCT__ 417 #define __FUNCT__ "MatFactorGetSchurComplement_MKL_PARDISO" 418 PetscErrorCode MatFactorGetSchurComplement_MKL_PARDISO(Mat F,Mat* S) 419 { 420 Mat St; 421 Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr; 422 PetscErrorCode ierr; 423 424 PetscFunctionBegin; 425 if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 426 else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 427 428 ierr = MatCreateSeqDense(PetscObjectComm((PetscObject)F),mpardiso->schur_size,mpardiso->schur_size,mpardiso->schur,&St);CHKERRQ(ierr); 429 *S = St; 430 PetscFunctionReturn(0); 431 } 432 433 #undef __FUNCT__ 434 #define __FUNCT__ "MatFactorInvertSchurComplement_MKL_PARDISO" 435 PetscErrorCode MatFactorInvertSchurComplement_MKL_PARDISO(Mat F) 436 { 437 Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr; 438 PetscErrorCode ierr; 439 440 PetscFunctionBegin; 441 if (!mpardiso->iparm[36-1]) { /* do nothing */ 442 PetscFunctionReturn(0); 443 } 444 if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 445 ierr = MatMKLPardisoInvertSchur_Private(mpardiso);CHKERRQ(ierr); 446 PetscFunctionReturn(0); 447 } 448 449 #undef __FUNCT__ 450 #define __FUNCT__ "MatFactorSolveSchurComplement_MKL_PARDISO" 451 PetscErrorCode MatFactorSolveSchurComplement_MKL_PARDISO(Mat F, Vec rhs, Vec sol) 452 { 453 Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr; 454 PetscScalar *asol,*arhs; 455 PetscErrorCode ierr; 456 457 PetscFunctionBegin; 458 if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 459 else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 460 461 mpardiso->nrhs = 1; 462 ierr = VecGetArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr); 463 ierr = VecGetArray(sol,&asol);CHKERRQ(ierr); 464 ierr = MatMKLPardisoSolveSchur_Private(mpardiso,arhs,asol);CHKERRQ(ierr); 465 ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr); 466 ierr = VecRestoreArray(sol,&asol);CHKERRQ(ierr); 467 PetscFunctionReturn(0); 468 } 469 470 #undef __FUNCT__ 471 #define __FUNCT__ "MatFactorSolveSchurComplementTranspose_MKL_PARDISO" 472 PetscErrorCode MatFactorSolveSchurComplementTranspose_MKL_PARDISO(Mat F, Vec rhs, Vec sol) 473 { 474 Mat_MKL_PARDISO *mpardiso =(Mat_MKL_PARDISO*)F->spptr; 475 PetscScalar *asol,*arhs; 476 PetscErrorCode ierr; 477 478 PetscFunctionBegin; 479 if (!mpardiso->iparm[36-1]) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur complement mode not selected! You should call MatFactorSetSchurIS to enable it"); 480 else if (!mpardiso->schur_size) SETERRQ(PetscObjectComm((PetscObject)F),PETSC_ERR_ORDER,"Schur indices not set! You should call MatFactorSetSchurIS before"); 481 482 mpardiso->nrhs = 1; 483 ierr = VecGetArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr); 484 ierr = VecGetArray(sol,&asol);CHKERRQ(ierr); 485 mpardiso->iparm[12 - 1] = 2; 486 ierr = MatMKLPardisoSolveSchur_Private(mpardiso,arhs,asol);CHKERRQ(ierr); 487 mpardiso->iparm[12 - 1] = 0; 488 ierr = VecRestoreArrayRead(rhs,(const PetscScalar**)&arhs);CHKERRQ(ierr); 489 ierr = VecRestoreArray(sol,&asol);CHKERRQ(ierr); 490 PetscFunctionReturn(0); 491 } 492 493 #undef __FUNCT__ 494 #define __FUNCT__ "MatDestroy_MKL_PARDISO" 495 PetscErrorCode MatDestroy_MKL_PARDISO(Mat A) 496 { 497 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 498 PetscErrorCode ierr; 499 500 PetscFunctionBegin; 501 if (mat_mkl_pardiso->CleanUp) { 502 mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY; 503 504 MKL_PARDISO (mat_mkl_pardiso->pt, 505 &mat_mkl_pardiso->maxfct, 506 &mat_mkl_pardiso->mnum, 507 &mat_mkl_pardiso->mtype, 508 &mat_mkl_pardiso->phase, 509 &mat_mkl_pardiso->n, 510 NULL, 511 NULL, 512 NULL, 513 mat_mkl_pardiso->perm, 514 &mat_mkl_pardiso->nrhs, 515 mat_mkl_pardiso->iparm, 516 &mat_mkl_pardiso->msglvl, 517 NULL, 518 NULL, 519 &mat_mkl_pardiso->err); 520 } 521 ierr = PetscFree(mat_mkl_pardiso->perm);CHKERRQ(ierr); 522 ierr = PetscFree2(mat_mkl_pardiso->schur,mat_mkl_pardiso->schur_work);CHKERRQ(ierr); 523 ierr = PetscFree(mat_mkl_pardiso->schur_idxs);CHKERRQ(ierr); 524 ierr = PetscFree(mat_mkl_pardiso->schur_pivots);CHKERRQ(ierr); 525 if (mat_mkl_pardiso->freeaij) { 526 ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr); 527 ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr); 528 } 529 if (mat_mkl_pardiso->Destroy) { 530 ierr = (mat_mkl_pardiso->Destroy)(A);CHKERRQ(ierr); 531 } 532 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 533 534 /* clear composed functions */ 535 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 536 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr); 537 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr); 538 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);CHKERRQ(ierr); 539 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);CHKERRQ(ierr); 540 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);CHKERRQ(ierr); 541 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);CHKERRQ(ierr); 542 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);CHKERRQ(ierr); 543 PetscFunctionReturn(0); 544 } 545 546 #undef __FUNCT__ 547 #define __FUNCT__ "MatMKLPardisoScatterSchur_Private" 548 static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce) 549 { 550 PetscFunctionBegin; 551 if (reduce) { /* data given for the whole matrix */ 552 PetscInt i,m=0,p=0; 553 for (i=0;i<mpardiso->nrhs;i++) { 554 PetscInt j; 555 for (j=0;j<mpardiso->schur_size;j++) { 556 schur[p+j] = whole[m+mpardiso->schur_idxs[j]]; 557 } 558 m += mpardiso->n; 559 p += mpardiso->schur_size; 560 } 561 } else { /* from Schur to whole */ 562 PetscInt i,m=0,p=0; 563 for (i=0;i<mpardiso->nrhs;i++) { 564 PetscInt j; 565 for (j=0;j<mpardiso->schur_size;j++) { 566 whole[m+mpardiso->schur_idxs[j]] = schur[p+j]; 567 } 568 m += mpardiso->n; 569 p += mpardiso->schur_size; 570 } 571 } 572 PetscFunctionReturn(0); 573 } 574 575 #undef __FUNCT__ 576 #define __FUNCT__ "MatSolve_MKL_PARDISO" 577 PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x) 578 { 579 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr; 580 PetscErrorCode ierr; 581 PetscScalar *xarray; 582 const PetscScalar *barray; 583 584 PetscFunctionBegin; 585 mat_mkl_pardiso->nrhs = 1; 586 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 587 ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 588 589 if (!mat_mkl_pardiso->schur) { 590 mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 591 } else { 592 mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION; 593 } 594 mat_mkl_pardiso->iparm[6-1] = 0; 595 596 MKL_PARDISO (mat_mkl_pardiso->pt, 597 &mat_mkl_pardiso->maxfct, 598 &mat_mkl_pardiso->mnum, 599 &mat_mkl_pardiso->mtype, 600 &mat_mkl_pardiso->phase, 601 &mat_mkl_pardiso->n, 602 mat_mkl_pardiso->a, 603 mat_mkl_pardiso->ia, 604 mat_mkl_pardiso->ja, 605 mat_mkl_pardiso->perm, 606 &mat_mkl_pardiso->nrhs, 607 mat_mkl_pardiso->iparm, 608 &mat_mkl_pardiso->msglvl, 609 (void*)barray, 610 (void*)xarray, 611 &mat_mkl_pardiso->err); 612 613 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); 614 615 if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */ 616 PetscInt shift = mat_mkl_pardiso->schur_size; 617 618 /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */ 619 if (!mat_mkl_pardiso->schur_inverted) { 620 shift = 0; 621 } 622 623 if (!mat_mkl_pardiso->solve_interior) { 624 /* solve Schur complement */ 625 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr); 626 ierr = MatMKLPardisoSolveSchur_Private(mat_mkl_pardiso,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr); 627 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr); 628 } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substitued to xarray[schur] will be neglected */ 629 PetscInt i; 630 for (i=0;i<mat_mkl_pardiso->schur_size;i++) { 631 xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.; 632 } 633 } 634 /* expansion phase */ 635 mat_mkl_pardiso->iparm[6-1] = 1; 636 mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION; 637 MKL_PARDISO (mat_mkl_pardiso->pt, 638 &mat_mkl_pardiso->maxfct, 639 &mat_mkl_pardiso->mnum, 640 &mat_mkl_pardiso->mtype, 641 &mat_mkl_pardiso->phase, 642 &mat_mkl_pardiso->n, 643 mat_mkl_pardiso->a, 644 mat_mkl_pardiso->ia, 645 mat_mkl_pardiso->ja, 646 mat_mkl_pardiso->perm, 647 &mat_mkl_pardiso->nrhs, 648 mat_mkl_pardiso->iparm, 649 &mat_mkl_pardiso->msglvl, 650 (void*)xarray, 651 (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */ 652 &mat_mkl_pardiso->err); 653 654 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); 655 mat_mkl_pardiso->iparm[6-1] = 0; 656 } 657 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 658 ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 659 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 660 PetscFunctionReturn(0); 661 } 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "MatSolveTranspose_MKL_PARDISO" 665 PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x) 666 { 667 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 668 PetscErrorCode ierr; 669 670 PetscFunctionBegin; 671 #if defined(PETSC_USE_COMPLEX) 672 mat_mkl_pardiso->iparm[12 - 1] = 1; 673 #else 674 mat_mkl_pardiso->iparm[12 - 1] = 2; 675 #endif 676 ierr = MatSolve_MKL_PARDISO(A,b,x);CHKERRQ(ierr); 677 mat_mkl_pardiso->iparm[12 - 1] = 0; 678 PetscFunctionReturn(0); 679 } 680 681 #undef __FUNCT__ 682 #define __FUNCT__ "MatMatSolve_MKL_PARDISO" 683 PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X) 684 { 685 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr; 686 PetscErrorCode ierr; 687 PetscScalar *barray, *xarray; 688 PetscBool flg; 689 690 PetscFunctionBegin; 691 ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr); 692 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix"); 693 ierr = PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr); 694 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix"); 695 696 ierr = MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);CHKERRQ(ierr); 697 698 if(mat_mkl_pardiso->nrhs > 0){ 699 ierr = MatDenseGetArray(B,&barray); 700 ierr = MatDenseGetArray(X,&xarray); 701 702 if (!mat_mkl_pardiso->schur) { 703 mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 704 } else { 705 mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION; 706 } 707 mat_mkl_pardiso->iparm[6-1] = 0; 708 709 MKL_PARDISO (mat_mkl_pardiso->pt, 710 &mat_mkl_pardiso->maxfct, 711 &mat_mkl_pardiso->mnum, 712 &mat_mkl_pardiso->mtype, 713 &mat_mkl_pardiso->phase, 714 &mat_mkl_pardiso->n, 715 mat_mkl_pardiso->a, 716 mat_mkl_pardiso->ia, 717 mat_mkl_pardiso->ja, 718 mat_mkl_pardiso->perm, 719 &mat_mkl_pardiso->nrhs, 720 mat_mkl_pardiso->iparm, 721 &mat_mkl_pardiso->msglvl, 722 (void*)barray, 723 (void*)xarray, 724 &mat_mkl_pardiso->err); 725 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); 726 727 if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */ 728 PetscScalar *o_schur_work = NULL; 729 PetscInt shift = mat_mkl_pardiso->schur_size*mat_mkl_pardiso->nrhs,scale; 730 PetscInt mem = mat_mkl_pardiso->n*mat_mkl_pardiso->nrhs; 731 732 /* allocate extra memory if it is needed */ 733 scale = 1; 734 if (mat_mkl_pardiso->schur_inverted) { 735 scale = 2; 736 } 737 mem *= scale; 738 if (mem > mat_mkl_pardiso->schur_work_size) { 739 o_schur_work = mat_mkl_pardiso->schur_work; 740 ierr = PetscMalloc1(mem,&mat_mkl_pardiso->schur_work);CHKERRQ(ierr); 741 } 742 743 /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */ 744 if (!mat_mkl_pardiso->schur_inverted) { 745 shift = 0; 746 } 747 748 /* solve Schur complement */ 749 if (!mat_mkl_pardiso->solve_interior) { 750 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr); 751 ierr = MatMKLPardisoSolveSchur_Private(mat_mkl_pardiso,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr); 752 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr); 753 } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substitued to xarray[schur,n] will be neglected */ 754 PetscInt i,n,m=0; 755 for (n=0;n<mat_mkl_pardiso->nrhs;n++) { 756 for (i=0;i<mat_mkl_pardiso->schur_size;i++) { 757 xarray[mat_mkl_pardiso->schur_idxs[i]+m] = 0.; 758 } 759 m += mat_mkl_pardiso->n; 760 } 761 } 762 /* expansion phase */ 763 mat_mkl_pardiso->iparm[6-1] = 1; 764 mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION; 765 MKL_PARDISO (mat_mkl_pardiso->pt, 766 &mat_mkl_pardiso->maxfct, 767 &mat_mkl_pardiso->mnum, 768 &mat_mkl_pardiso->mtype, 769 &mat_mkl_pardiso->phase, 770 &mat_mkl_pardiso->n, 771 mat_mkl_pardiso->a, 772 mat_mkl_pardiso->ia, 773 mat_mkl_pardiso->ja, 774 mat_mkl_pardiso->perm, 775 &mat_mkl_pardiso->nrhs, 776 mat_mkl_pardiso->iparm, 777 &mat_mkl_pardiso->msglvl, 778 (void*)xarray, 779 (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */ 780 &mat_mkl_pardiso->err); 781 if (o_schur_work) { /* restore original schur_work (minimal size) */ 782 ierr = PetscFree(mat_mkl_pardiso->schur_work);CHKERRQ(ierr); 783 mat_mkl_pardiso->schur_work = o_schur_work; 784 } 785 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); 786 mat_mkl_pardiso->iparm[6-1] = 0; 787 } 788 } 789 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNCT__ 794 #define __FUNCT__ "MatFactorNumeric_MKL_PARDISO" 795 PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info) 796 { 797 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->spptr; 798 PetscErrorCode ierr; 799 800 PetscFunctionBegin; 801 mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN; 802 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,&mat_mkl_pardiso->a);CHKERRQ(ierr); 803 804 mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION; 805 MKL_PARDISO (mat_mkl_pardiso->pt, 806 &mat_mkl_pardiso->maxfct, 807 &mat_mkl_pardiso->mnum, 808 &mat_mkl_pardiso->mtype, 809 &mat_mkl_pardiso->phase, 810 &mat_mkl_pardiso->n, 811 mat_mkl_pardiso->a, 812 mat_mkl_pardiso->ia, 813 mat_mkl_pardiso->ja, 814 mat_mkl_pardiso->perm, 815 &mat_mkl_pardiso->nrhs, 816 mat_mkl_pardiso->iparm, 817 &mat_mkl_pardiso->msglvl, 818 NULL, 819 (void*)mat_mkl_pardiso->schur, 820 &mat_mkl_pardiso->err); 821 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); 822 823 if (mat_mkl_pardiso->schur) { /* schur output from pardiso is in row major format */ 824 PetscInt j,k,n=mat_mkl_pardiso->schur_size; 825 if (!mat_mkl_pardiso->schur_solver_type) { 826 for (j=0; j<n; j++) { 827 for (k=0; k<j; k++) { 828 PetscScalar tmp = mat_mkl_pardiso->schur[j + k*n]; 829 mat_mkl_pardiso->schur[j + k*n] = mat_mkl_pardiso->schur[k + j*n]; 830 mat_mkl_pardiso->schur[k + j*n] = tmp; 831 } 832 } 833 } else { /* we could use row-major in LAPACK routines (e.g. use 'U' instead of 'L'; instead, I prefer consistency between data structures and swap to column major */ 834 for (j=0; j<n; j++) { 835 for (k=0; k<j; k++) { 836 mat_mkl_pardiso->schur[j + k*n] = mat_mkl_pardiso->schur[k + j*n]; 837 } 838 } 839 } 840 } 841 mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN; 842 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 843 mat_mkl_pardiso->schur_factored = PETSC_FALSE; 844 mat_mkl_pardiso->schur_inverted = PETSC_FALSE; 845 mat_mkl_pardiso->schur_solver_type = 0; 846 PetscFunctionReturn(0); 847 } 848 849 #undef __FUNCT__ 850 #define __FUNCT__ "PetscSetMKL_PARDISOFromOptions" 851 PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A) 852 { 853 Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr; 854 PetscErrorCode ierr; 855 PetscInt icntl, threads = 1; 856 PetscBool flg; 857 858 PetscFunctionBegin; 859 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");CHKERRQ(ierr); 860 ierr = PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use","None",threads,&threads,&flg);CHKERRQ(ierr); 861 if (flg) mkl_set_num_threads((int)threads); 862 863 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); 864 if (flg) mat_mkl_pardiso->maxfct = icntl; 865 866 ierr = PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);CHKERRQ(ierr); 867 if (flg) mat_mkl_pardiso->mnum = icntl; 868 869 ierr = PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);CHKERRQ(ierr); 870 if (flg) mat_mkl_pardiso->msglvl = icntl; 871 872 ierr = PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);CHKERRQ(ierr); 873 if(flg){ 874 void *pt[IPARM_SIZE]; 875 mat_mkl_pardiso->mtype = icntl; 876 MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm); 877 #if defined(PETSC_USE_REAL_SINGLE) 878 mat_mkl_pardiso->iparm[27] = 1; 879 #else 880 mat_mkl_pardiso->iparm[27] = 0; 881 #endif 882 mat_mkl_pardiso->iparm[34] = 1; 883 } 884 ierr = PetscOptionsInt("-mat_mkl_pardiso_1","Use default values","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);CHKERRQ(ierr); 885 886 if(flg && icntl != 0){ 887 ierr = PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);CHKERRQ(ierr); 888 if (flg) mat_mkl_pardiso->iparm[1] = icntl; 889 890 ierr = PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);CHKERRQ(ierr); 891 if (flg) mat_mkl_pardiso->iparm[3] = icntl; 892 893 ierr = PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);CHKERRQ(ierr); 894 if (flg) mat_mkl_pardiso->iparm[4] = icntl; 895 896 ierr = PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);CHKERRQ(ierr); 897 if (flg) mat_mkl_pardiso->iparm[5] = icntl; 898 899 ierr = PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);CHKERRQ(ierr); 900 if (flg) mat_mkl_pardiso->iparm[7] = icntl; 901 902 ierr = PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);CHKERRQ(ierr); 903 if (flg) mat_mkl_pardiso->iparm[9] = icntl; 904 905 ierr = PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);CHKERRQ(ierr); 906 if (flg) mat_mkl_pardiso->iparm[10] = icntl; 907 908 ierr = PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);CHKERRQ(ierr); 909 if (flg) mat_mkl_pardiso->iparm[11] = icntl; 910 911 ierr = PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);CHKERRQ(ierr); 912 if (flg) mat_mkl_pardiso->iparm[12] = icntl; 913 914 ierr = PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);CHKERRQ(ierr); 915 if (flg) mat_mkl_pardiso->iparm[17] = icntl; 916 917 ierr = PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);CHKERRQ(ierr); 918 if (flg) mat_mkl_pardiso->iparm[18] = icntl; 919 920 ierr = PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);CHKERRQ(ierr); 921 if (flg) mat_mkl_pardiso->iparm[20] = icntl; 922 923 ierr = PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);CHKERRQ(ierr); 924 if (flg) mat_mkl_pardiso->iparm[23] = icntl; 925 926 ierr = PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);CHKERRQ(ierr); 927 if (flg) mat_mkl_pardiso->iparm[24] = icntl; 928 929 ierr = PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);CHKERRQ(ierr); 930 if (flg) mat_mkl_pardiso->iparm[26] = icntl; 931 932 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); 933 if (flg) mat_mkl_pardiso->iparm[30] = icntl; 934 935 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); 936 if (flg) mat_mkl_pardiso->iparm[33] = icntl; 937 938 ierr = PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);CHKERRQ(ierr); 939 if (flg) mat_mkl_pardiso->iparm[59] = icntl; 940 } 941 PetscOptionsEnd(); 942 PetscFunctionReturn(0); 943 } 944 945 #undef __FUNCT__ 946 #define __FUNCT__ "MatFactorMKL_PARDISOInitialize_Private" 947 PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso) 948 { 949 PetscErrorCode ierr; 950 PetscInt i; 951 952 PetscFunctionBegin; 953 for ( i = 0; i < IPARM_SIZE; i++ ){ 954 mat_mkl_pardiso->iparm[i] = 0; 955 } 956 for ( i = 0; i < IPARM_SIZE; i++ ){ 957 mat_mkl_pardiso->pt[i] = 0; 958 } 959 /* Default options for both sym and unsym */ 960 mat_mkl_pardiso->iparm[ 0] = 1; /* Solver default parameters overriden with provided by iparm */ 961 mat_mkl_pardiso->iparm[ 1] = 2; /* Metis reordering */ 962 mat_mkl_pardiso->iparm[ 5] = 0; /* Write solution into x */ 963 mat_mkl_pardiso->iparm[ 7] = 0; /* Max number of iterative refinement steps */ 964 mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */ 965 mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */ 966 #if 0 967 mat_mkl_pardiso->iparm[23] = 1; /* Parallel factorization control*/ 968 #endif 969 mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */ 970 mat_mkl_pardiso->iparm[39] = 0; /* Input: matrix/rhs/solution stored on master */ 971 972 mat_mkl_pardiso->CleanUp = PETSC_FALSE; 973 mat_mkl_pardiso->maxfct = 1; /* Maximum number of numerical factorizations. */ 974 mat_mkl_pardiso->mnum = 1; /* Which factorization to use. */ 975 mat_mkl_pardiso->msglvl = 0; /* 0: do not print 1: Print statistical information in file */ 976 mat_mkl_pardiso->phase = -1; 977 mat_mkl_pardiso->err = 0; 978 979 mat_mkl_pardiso->n = A->rmap->N; 980 mat_mkl_pardiso->nrhs = 1; 981 mat_mkl_pardiso->err = 0; 982 mat_mkl_pardiso->phase = -1; 983 984 if(ftype == MAT_FACTOR_LU){ 985 mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */ 986 mat_mkl_pardiso->iparm[10] = 1; /* Use nonsymmetric permutation and scaling MPS */ 987 mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 988 989 } else { 990 mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */ 991 mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */ 992 mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 993 /* mat_mkl_pardiso->iparm[20] = 1; */ /* Apply 1x1 and 2x2 Bunch-Kaufman pivoting during the factorization process */ 994 #if defined(PETSC_USE_DEBUG) 995 mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */ 996 #endif 997 } 998 ierr = PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);CHKERRQ(ierr); 999 for(i = 0; i < A->rmap->N; i++){ 1000 mat_mkl_pardiso->perm[i] = 0; 1001 } 1002 mat_mkl_pardiso->schur_size = 0; 1003 PetscFunctionReturn(0); 1004 } 1005 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "MatFactorSymbolic_AIJMKL_PARDISO_Private" 1008 PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info) 1009 { 1010 Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr; 1011 PetscErrorCode ierr; 1012 1013 PetscFunctionBegin; 1014 mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN; 1015 ierr = PetscSetMKL_PARDISOFromOptions(F,A);CHKERRQ(ierr); 1016 1017 /* throw away any previously computed structure */ 1018 if (mat_mkl_pardiso->freeaij) { 1019 ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr); 1020 ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr); 1021 } 1022 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,&mat_mkl_pardiso->a);CHKERRQ(ierr); 1023 mat_mkl_pardiso->n = A->rmap->N; 1024 1025 mat_mkl_pardiso->phase = JOB_ANALYSIS; 1026 1027 MKL_PARDISO (mat_mkl_pardiso->pt, 1028 &mat_mkl_pardiso->maxfct, 1029 &mat_mkl_pardiso->mnum, 1030 &mat_mkl_pardiso->mtype, 1031 &mat_mkl_pardiso->phase, 1032 &mat_mkl_pardiso->n, 1033 mat_mkl_pardiso->a, 1034 mat_mkl_pardiso->ia, 1035 mat_mkl_pardiso->ja, 1036 mat_mkl_pardiso->perm, 1037 &mat_mkl_pardiso->nrhs, 1038 mat_mkl_pardiso->iparm, 1039 &mat_mkl_pardiso->msglvl, 1040 NULL, 1041 NULL, 1042 &mat_mkl_pardiso->err); 1043 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); 1044 1045 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 1046 1047 if (F->factortype == MAT_FACTOR_LU) { 1048 F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO; 1049 } else { 1050 F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO; 1051 } 1052 F->ops->solve = MatSolve_MKL_PARDISO; 1053 F->ops->solvetranspose = MatSolveTranspose_MKL_PARDISO; 1054 F->ops->matsolve = MatMatSolve_MKL_PARDISO; 1055 PetscFunctionReturn(0); 1056 } 1057 1058 #undef __FUNCT__ 1059 #define __FUNCT__ "MatLUFactorSymbolic_AIJMKL_PARDISO" 1060 PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1061 { 1062 PetscErrorCode ierr; 1063 1064 PetscFunctionBegin; 1065 ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr); 1066 PetscFunctionReturn(0); 1067 } 1068 1069 #undef __FUNCT__ 1070 #define __FUNCT__ "MatCholeskyFactorSymbolic_AIJMKL_PARDISO" 1071 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info) 1072 { 1073 PetscErrorCode ierr; 1074 1075 PetscFunctionBegin; 1076 ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr); 1077 PetscFunctionReturn(0); 1078 } 1079 1080 #undef __FUNCT__ 1081 #define __FUNCT__ "MatView_MKL_PARDISO" 1082 PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer) 1083 { 1084 PetscErrorCode ierr; 1085 PetscBool iascii; 1086 PetscViewerFormat format; 1087 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 1088 PetscInt i; 1089 1090 PetscFunctionBegin; 1091 if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(0); 1092 1093 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1094 if (iascii) { 1095 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1096 if (format == PETSC_VIEWER_ASCII_INFO) { 1097 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");CHKERRQ(ierr); 1098 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase: %d \n",mat_mkl_pardiso->phase);CHKERRQ(ierr); 1099 for(i = 1; i <= 64; i++){ 1100 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]: %d \n",i, mat_mkl_pardiso->iparm[i - 1]);CHKERRQ(ierr); 1101 } 1102 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct: %d \n", mat_mkl_pardiso->maxfct);CHKERRQ(ierr); 1103 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum: %d \n", mat_mkl_pardiso->mnum);CHKERRQ(ierr); 1104 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype: %d \n", mat_mkl_pardiso->mtype);CHKERRQ(ierr); 1105 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n: %d \n", mat_mkl_pardiso->n);CHKERRQ(ierr); 1106 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs: %d \n", mat_mkl_pardiso->nrhs);CHKERRQ(ierr); 1107 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl: %d \n", mat_mkl_pardiso->msglvl);CHKERRQ(ierr); 1108 } 1109 } 1110 PetscFunctionReturn(0); 1111 } 1112 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "MatGetInfo_MKL_PARDISO" 1116 PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info) 1117 { 1118 Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)A->spptr; 1119 1120 PetscFunctionBegin; 1121 info->block_size = 1.0; 1122 info->nz_allocated = mat_mkl_pardiso->nz + 0.0; 1123 info->nz_unneeded = 0.0; 1124 info->assemblies = 0.0; 1125 info->mallocs = 0.0; 1126 info->memory = 0.0; 1127 info->fill_ratio_given = 0; 1128 info->fill_ratio_needed = 0; 1129 info->factor_mallocs = 0; 1130 PetscFunctionReturn(0); 1131 } 1132 1133 #undef __FUNCT__ 1134 #define __FUNCT__ "MatMkl_PardisoSetCntl_MKL_PARDISO" 1135 PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival) 1136 { 1137 Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)F->spptr; 1138 1139 PetscFunctionBegin; 1140 if(icntl <= 64){ 1141 mat_mkl_pardiso->iparm[icntl - 1] = ival; 1142 } else { 1143 if(icntl == 65) 1144 mkl_set_num_threads((int)ival); 1145 else if(icntl == 66) 1146 mat_mkl_pardiso->maxfct = ival; 1147 else if(icntl == 67) 1148 mat_mkl_pardiso->mnum = ival; 1149 else if(icntl == 68) 1150 mat_mkl_pardiso->msglvl = ival; 1151 else if(icntl == 69){ 1152 void *pt[IPARM_SIZE]; 1153 mat_mkl_pardiso->mtype = ival; 1154 MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm); 1155 #if defined(PETSC_USE_REAL_SINGLE) 1156 mat_mkl_pardiso->iparm[27] = 1; 1157 #else 1158 mat_mkl_pardiso->iparm[27] = 0; 1159 #endif 1160 mat_mkl_pardiso->iparm[34] = 1; 1161 } else if(icntl==70) { 1162 mat_mkl_pardiso->solve_interior = !!ival; 1163 } 1164 } 1165 PetscFunctionReturn(0); 1166 } 1167 1168 #undef __FUNCT__ 1169 #define __FUNCT__ "MatMkl_PardisoSetCntl" 1170 /*@ 1171 MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters 1172 1173 Logically Collective on Mat 1174 1175 Input Parameters: 1176 + F - the factored matrix obtained by calling MatGetFactor() 1177 . icntl - index of Mkl_Pardiso parameter 1178 - ival - value of Mkl_Pardiso parameter 1179 1180 Options Database: 1181 . -mat_mkl_pardiso_<icntl> <ival> 1182 1183 Level: beginner 1184 1185 References: Mkl_Pardiso Users' Guide 1186 1187 .seealso: MatGetFactor() 1188 @*/ 1189 PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival) 1190 { 1191 PetscErrorCode ierr; 1192 1193 PetscFunctionBegin; 1194 ierr = PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 1195 PetscFunctionReturn(0); 1196 } 1197 1198 /*MC 1199 MATSOLVERMKL_PARDISO - A matrix type providing direct solvers (LU) for 1200 sequential matrices via the external package MKL_PARDISO. 1201 1202 Works with MATSEQAIJ matrices 1203 1204 Use -pc_type lu -pc_factor_mat_solver_package mkl_pardiso to us this direct solver 1205 1206 Options Database Keys: 1207 + -mat_mkl_pardiso_65 - Number of threads to use 1208 . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time 1209 . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase 1210 . -mat_mkl_pardiso_68 - Message level information 1211 . -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 1212 . -mat_mkl_pardiso_1 - Use default values 1213 . -mat_mkl_pardiso_2 - Fill-in reducing ordering for the input matrix 1214 . -mat_mkl_pardiso_4 - Preconditioned CGS/CG 1215 . -mat_mkl_pardiso_5 - User permutation 1216 . -mat_mkl_pardiso_6 - Write solution on x 1217 . -mat_mkl_pardiso_8 - Iterative refinement step 1218 . -mat_mkl_pardiso_10 - Pivoting perturbation 1219 . -mat_mkl_pardiso_11 - Scaling vectors 1220 . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A 1221 . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching 1222 . -mat_mkl_pardiso_18 - Numbers of non-zero elements 1223 . -mat_mkl_pardiso_19 - Report number of floating point operations 1224 . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices 1225 . -mat_mkl_pardiso_24 - Parallel factorization control 1226 . -mat_mkl_pardiso_25 - Parallel forward/backward solve control 1227 . -mat_mkl_pardiso_27 - Matrix checker 1228 . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors 1229 . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode 1230 - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode 1231 1232 Level: beginner 1233 1234 For more information please check mkl_pardiso manual 1235 1236 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 1237 1238 M*/ 1239 #undef __FUNCT__ 1240 #define __FUNCT__ "MatFactorGetSolverPackage_mkl_pardiso" 1241 static PetscErrorCode MatFactorGetSolverPackage_mkl_pardiso(Mat A, const MatSolverPackage *type) 1242 { 1243 PetscFunctionBegin; 1244 *type = MATSOLVERMKL_PARDISO; 1245 PetscFunctionReturn(0); 1246 } 1247 1248 #undef __FUNCT__ 1249 #define __FUNCT__ "MatGetFactor_aij_mkl_pardiso" 1250 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F) 1251 { 1252 Mat B; 1253 PetscErrorCode ierr; 1254 Mat_MKL_PARDISO *mat_mkl_pardiso; 1255 PetscBool isSeqAIJ,isSeqBAIJ,isSeqSBAIJ; 1256 1257 PetscFunctionBegin; 1258 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 1259 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 1260 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 1261 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 1262 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1263 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1264 ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 1265 ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 1266 ierr = MatSeqSBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 1267 1268 ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr); 1269 B->spptr = mat_mkl_pardiso; 1270 1271 ierr = MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr); 1272 if (ftype == MAT_FACTOR_LU) { 1273 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO; 1274 B->factortype = MAT_FACTOR_LU; 1275 mat_mkl_pardiso->needsym = PETSC_FALSE; 1276 if (isSeqAIJ) { 1277 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij; 1278 } else if (isSeqBAIJ) { 1279 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij; 1280 } else if (isSeqSBAIJ) { 1281 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead"); 1282 } else { 1283 SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU with %s format",((PetscObject)A)->type_name); 1284 } 1285 #if defined(PETSC_USE_COMPLEX) 1286 mat_mkl_pardiso->schur_solver_type = 0; /* use a general solver for the moment */ 1287 mat_mkl_pardiso->mtype = 13; 1288 #else 1289 if (A->structurally_symmetric) { 1290 mat_mkl_pardiso->mtype = 1; 1291 } else { 1292 mat_mkl_pardiso->mtype = 11; 1293 } 1294 #endif 1295 mat_mkl_pardiso->schur_solver_type = 0; 1296 } else { 1297 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO; 1298 B->factortype = MAT_FACTOR_CHOLESKY; 1299 if (isSeqAIJ) { 1300 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij; 1301 } else if (isSeqBAIJ) { 1302 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij; 1303 } else if (isSeqSBAIJ) { 1304 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij; 1305 } else { 1306 SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with %s format",((PetscObject)A)->type_name); 1307 } 1308 mat_mkl_pardiso->needsym = PETSC_TRUE; 1309 #if defined(PETSC_USE_COMPLEX) 1310 mat_mkl_pardiso->schur_solver_type = 0; /* use a general solver for the moment */ 1311 mat_mkl_pardiso->mtype = 13; 1312 #else 1313 if (A->spd_set && A->spd) { 1314 mat_mkl_pardiso->schur_solver_type = 2; 1315 mat_mkl_pardiso->mtype = 2; 1316 } else { 1317 mat_mkl_pardiso->schur_solver_type = 1; 1318 mat_mkl_pardiso->mtype = -2; 1319 } 1320 #endif 1321 } 1322 mat_mkl_pardiso->Destroy = B->ops->destroy; 1323 B->ops->destroy = MatDestroy_MKL_PARDISO; 1324 B->ops->view = MatView_MKL_PARDISO; 1325 B->factortype = ftype; 1326 B->ops->getinfo = MatGetInfo_MKL_PARDISO; 1327 B->assembled = PETSC_TRUE; 1328 1329 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr); 1330 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MKL_PARDISO);CHKERRQ(ierr); 1331 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1332 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1333 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1334 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1335 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MKL_PARDISO);CHKERRQ(ierr); 1336 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr); 1337 1338 *F = B; 1339 PetscFunctionReturn(0); 1340 } 1341 1342 #undef __FUNCT__ 1343 #define __FUNCT__ "MatSolverPackageRegister_MKL_Pardiso" 1344 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MKL_Pardiso(void) 1345 { 1346 PetscErrorCode ierr; 1347 1348 PetscFunctionBegin; 1349 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr); 1350 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr); 1351 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr); 1352 PetscFunctionReturn(0); 1353 } 1354 1355