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 PetscBool isSeqSBAIJ; 499 PetscErrorCode ierr; 500 501 PetscFunctionBegin; 502 if (mat_mkl_pardiso->CleanUp) { 503 mat_mkl_pardiso->phase = JOB_RELEASE_OF_ALL_MEMORY; 504 505 MKL_PARDISO (mat_mkl_pardiso->pt, 506 &mat_mkl_pardiso->maxfct, 507 &mat_mkl_pardiso->mnum, 508 &mat_mkl_pardiso->mtype, 509 &mat_mkl_pardiso->phase, 510 &mat_mkl_pardiso->n, 511 NULL, 512 NULL, 513 NULL, 514 mat_mkl_pardiso->perm, 515 &mat_mkl_pardiso->nrhs, 516 mat_mkl_pardiso->iparm, 517 &mat_mkl_pardiso->msglvl, 518 NULL, 519 NULL, 520 &mat_mkl_pardiso->err); 521 } 522 ierr = PetscFree(mat_mkl_pardiso->perm);CHKERRQ(ierr); 523 ierr = PetscFree2(mat_mkl_pardiso->schur,mat_mkl_pardiso->schur_work);CHKERRQ(ierr); 524 ierr = PetscFree(mat_mkl_pardiso->schur_idxs);CHKERRQ(ierr); 525 ierr = PetscFree(mat_mkl_pardiso->schur_pivots);CHKERRQ(ierr); 526 if (mat_mkl_pardiso->freeaij) { 527 ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr); 528 ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr); 529 } 530 if (mat_mkl_pardiso->Destroy) { 531 ierr = (mat_mkl_pardiso->Destroy)(A);CHKERRQ(ierr); 532 } 533 ierr = PetscFree(A->spptr);CHKERRQ(ierr); 534 535 /* clear composed functions */ 536 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr); 537 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSetSchurIS_C",NULL);CHKERRQ(ierr); 538 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorCreateSchurComplement_C",NULL);CHKERRQ(ierr); 539 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSchurComplement_C",NULL);CHKERRQ(ierr); 540 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorInvertSchurComplement_C",NULL);CHKERRQ(ierr); 541 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplement_C",NULL);CHKERRQ(ierr); 542 ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorSolveSchurComplementTranspose_C",NULL);CHKERRQ(ierr); 543 ierr = PetscObjectComposeFunction((PetscObject)A,"MatMkl_PardisoSetCntl_C",NULL);CHKERRQ(ierr); 544 PetscFunctionReturn(0); 545 } 546 547 #undef __FUNCT__ 548 #define __FUNCT__ "MatMKLPardisoScatterSchur_Private" 549 static PetscErrorCode MatMKLPardisoScatterSchur_Private(Mat_MKL_PARDISO *mpardiso, PetscScalar *whole, PetscScalar *schur, PetscBool reduce) 550 { 551 PetscFunctionBegin; 552 if (reduce) { /* data given for the whole matrix */ 553 PetscInt i,m=0,p=0; 554 for (i=0;i<mpardiso->nrhs;i++) { 555 PetscInt j; 556 for (j=0;j<mpardiso->schur_size;j++) { 557 schur[p+j] = whole[m+mpardiso->schur_idxs[j]]; 558 } 559 m += mpardiso->n; 560 p += mpardiso->schur_size; 561 } 562 } else { /* from Schur to whole */ 563 PetscInt i,m=0,p=0; 564 for (i=0;i<mpardiso->nrhs;i++) { 565 PetscInt j; 566 for (j=0;j<mpardiso->schur_size;j++) { 567 whole[m+mpardiso->schur_idxs[j]] = schur[p+j]; 568 } 569 m += mpardiso->n; 570 p += mpardiso->schur_size; 571 } 572 } 573 PetscFunctionReturn(0); 574 } 575 576 #undef __FUNCT__ 577 #define __FUNCT__ "MatSolve_MKL_PARDISO" 578 PetscErrorCode MatSolve_MKL_PARDISO(Mat A,Vec b,Vec x) 579 { 580 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr; 581 PetscErrorCode ierr; 582 PetscScalar *xarray; 583 const PetscScalar *barray; 584 585 PetscFunctionBegin; 586 mat_mkl_pardiso->nrhs = 1; 587 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 588 ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr); 589 590 if (!mat_mkl_pardiso->schur) { 591 mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 592 } else { 593 mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION; 594 } 595 mat_mkl_pardiso->iparm[6-1] = 0; 596 597 MKL_PARDISO (mat_mkl_pardiso->pt, 598 &mat_mkl_pardiso->maxfct, 599 &mat_mkl_pardiso->mnum, 600 &mat_mkl_pardiso->mtype, 601 &mat_mkl_pardiso->phase, 602 &mat_mkl_pardiso->n, 603 mat_mkl_pardiso->a, 604 mat_mkl_pardiso->ia, 605 mat_mkl_pardiso->ja, 606 mat_mkl_pardiso->perm, 607 &mat_mkl_pardiso->nrhs, 608 mat_mkl_pardiso->iparm, 609 &mat_mkl_pardiso->msglvl, 610 (void*)barray, 611 (void*)xarray, 612 &mat_mkl_pardiso->err); 613 614 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); 615 616 if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */ 617 PetscInt shift = mat_mkl_pardiso->schur_size; 618 619 /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */ 620 if (!mat_mkl_pardiso->schur_inverted) { 621 shift = 0; 622 } 623 624 if (!mat_mkl_pardiso->solve_interior) { 625 /* solve Schur complement */ 626 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr); 627 ierr = MatMKLPardisoSolveSchur_Private(mat_mkl_pardiso,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr); 628 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr); 629 } else { /* if we are solving for the interior problem, any value in barray[schur] forward-substitued to xarray[schur] will be neglected */ 630 PetscInt i; 631 for (i=0;i<mat_mkl_pardiso->schur_size;i++) { 632 xarray[mat_mkl_pardiso->schur_idxs[i]] = 0.; 633 } 634 } 635 /* expansion phase */ 636 mat_mkl_pardiso->iparm[6-1] = 1; 637 mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION; 638 MKL_PARDISO (mat_mkl_pardiso->pt, 639 &mat_mkl_pardiso->maxfct, 640 &mat_mkl_pardiso->mnum, 641 &mat_mkl_pardiso->mtype, 642 &mat_mkl_pardiso->phase, 643 &mat_mkl_pardiso->n, 644 mat_mkl_pardiso->a, 645 mat_mkl_pardiso->ia, 646 mat_mkl_pardiso->ja, 647 mat_mkl_pardiso->perm, 648 &mat_mkl_pardiso->nrhs, 649 mat_mkl_pardiso->iparm, 650 &mat_mkl_pardiso->msglvl, 651 (void*)xarray, 652 (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */ 653 &mat_mkl_pardiso->err); 654 655 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); 656 mat_mkl_pardiso->iparm[6-1] = 0; 657 } 658 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 659 ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr); 660 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 661 PetscFunctionReturn(0); 662 } 663 664 #undef __FUNCT__ 665 #define __FUNCT__ "MatSolveTranspose_MKL_PARDISO" 666 PetscErrorCode MatSolveTranspose_MKL_PARDISO(Mat A,Vec b,Vec x) 667 { 668 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 669 PetscErrorCode ierr; 670 671 PetscFunctionBegin; 672 #if defined(PETSC_USE_COMPLEX) 673 mat_mkl_pardiso->iparm[12 - 1] = 1; 674 #else 675 mat_mkl_pardiso->iparm[12 - 1] = 2; 676 #endif 677 ierr = MatSolve_MKL_PARDISO(A,b,x);CHKERRQ(ierr); 678 mat_mkl_pardiso->iparm[12 - 1] = 0; 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "MatMatSolve_MKL_PARDISO" 684 PetscErrorCode MatMatSolve_MKL_PARDISO(Mat A,Mat B,Mat X) 685 { 686 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(A)->spptr; 687 PetscErrorCode ierr; 688 PetscScalar *barray, *xarray; 689 PetscBool flg; 690 691 PetscFunctionBegin; 692 ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&flg);CHKERRQ(ierr); 693 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATSEQDENSE matrix"); 694 ierr = PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&flg);CHKERRQ(ierr); 695 if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATSEQDENSE matrix"); 696 697 ierr = MatGetSize(B,NULL,(PetscInt*)&mat_mkl_pardiso->nrhs);CHKERRQ(ierr); 698 699 if(mat_mkl_pardiso->nrhs > 0){ 700 ierr = MatDenseGetArray(B,&barray); 701 ierr = MatDenseGetArray(X,&xarray); 702 703 if (!mat_mkl_pardiso->schur) { 704 mat_mkl_pardiso->phase = JOB_SOLVE_ITERATIVE_REFINEMENT; 705 } else { 706 mat_mkl_pardiso->phase = JOB_SOLVE_FORWARD_SUBSTITUTION; 707 } 708 mat_mkl_pardiso->iparm[6-1] = 0; 709 710 MKL_PARDISO (mat_mkl_pardiso->pt, 711 &mat_mkl_pardiso->maxfct, 712 &mat_mkl_pardiso->mnum, 713 &mat_mkl_pardiso->mtype, 714 &mat_mkl_pardiso->phase, 715 &mat_mkl_pardiso->n, 716 mat_mkl_pardiso->a, 717 mat_mkl_pardiso->ia, 718 mat_mkl_pardiso->ja, 719 mat_mkl_pardiso->perm, 720 &mat_mkl_pardiso->nrhs, 721 mat_mkl_pardiso->iparm, 722 &mat_mkl_pardiso->msglvl, 723 (void*)barray, 724 (void*)xarray, 725 &mat_mkl_pardiso->err); 726 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); 727 728 if (mat_mkl_pardiso->schur) { /* solve Schur complement and expand solution */ 729 PetscScalar *o_schur_work = NULL; 730 PetscInt shift = mat_mkl_pardiso->schur_size*mat_mkl_pardiso->nrhs,scale; 731 PetscInt mem = mat_mkl_pardiso->n*mat_mkl_pardiso->nrhs; 732 733 /* allocate extra memory if it is needed */ 734 scale = 1; 735 if (mat_mkl_pardiso->schur_inverted) { 736 scale = 2; 737 } 738 mem *= scale; 739 if (mem > mat_mkl_pardiso->schur_work_size) { 740 o_schur_work = mat_mkl_pardiso->schur_work; 741 ierr = PetscMalloc1(mem,&mat_mkl_pardiso->schur_work);CHKERRQ(ierr); 742 } 743 744 /* if inverted, uses BLAS *MM subroutines, otherwise LAPACK *TRS */ 745 if (!mat_mkl_pardiso->schur_inverted) { 746 shift = 0; 747 } 748 749 /* solve Schur complement */ 750 if (!mat_mkl_pardiso->solve_interior) { 751 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work,PETSC_TRUE);CHKERRQ(ierr); 752 ierr = MatMKLPardisoSolveSchur_Private(mat_mkl_pardiso,mat_mkl_pardiso->schur_work,mat_mkl_pardiso->schur_work+shift);CHKERRQ(ierr); 753 ierr = MatMKLPardisoScatterSchur_Private(mat_mkl_pardiso,xarray,mat_mkl_pardiso->schur_work+shift,PETSC_FALSE);CHKERRQ(ierr); 754 } else { /* if we are solving for the interior problem, any value in barray[schur,n] forward-substitued to xarray[schur,n] will be neglected */ 755 PetscInt i,n,m=0; 756 for (n=0;n<mat_mkl_pardiso->nrhs;n++) { 757 for (i=0;i<mat_mkl_pardiso->schur_size;i++) { 758 xarray[mat_mkl_pardiso->schur_idxs[i]+m] = 0.; 759 } 760 m += mat_mkl_pardiso->n; 761 } 762 } 763 /* expansion phase */ 764 mat_mkl_pardiso->iparm[6-1] = 1; 765 mat_mkl_pardiso->phase = JOB_SOLVE_BACKWARD_SUBSTITUTION; 766 MKL_PARDISO (mat_mkl_pardiso->pt, 767 &mat_mkl_pardiso->maxfct, 768 &mat_mkl_pardiso->mnum, 769 &mat_mkl_pardiso->mtype, 770 &mat_mkl_pardiso->phase, 771 &mat_mkl_pardiso->n, 772 mat_mkl_pardiso->a, 773 mat_mkl_pardiso->ia, 774 mat_mkl_pardiso->ja, 775 mat_mkl_pardiso->perm, 776 &mat_mkl_pardiso->nrhs, 777 mat_mkl_pardiso->iparm, 778 &mat_mkl_pardiso->msglvl, 779 (void*)xarray, 780 (void*)mat_mkl_pardiso->schur_work, /* according to the specs, the solution vector is always used */ 781 &mat_mkl_pardiso->err); 782 if (o_schur_work) { /* restore original schur_work (minimal size) */ 783 ierr = PetscFree(mat_mkl_pardiso->schur_work);CHKERRQ(ierr); 784 mat_mkl_pardiso->schur_work = o_schur_work; 785 } 786 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); 787 mat_mkl_pardiso->iparm[6-1] = 0; 788 } 789 } 790 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 791 PetscFunctionReturn(0); 792 } 793 794 #undef __FUNCT__ 795 #define __FUNCT__ "MatFactorNumeric_MKL_PARDISO" 796 PetscErrorCode MatFactorNumeric_MKL_PARDISO(Mat F,Mat A,const MatFactorInfo *info) 797 { 798 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)(F)->spptr; 799 PetscErrorCode ierr; 800 801 PetscFunctionBegin; 802 mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN; 803 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); 804 805 mat_mkl_pardiso->phase = JOB_NUMERICAL_FACTORIZATION; 806 MKL_PARDISO (mat_mkl_pardiso->pt, 807 &mat_mkl_pardiso->maxfct, 808 &mat_mkl_pardiso->mnum, 809 &mat_mkl_pardiso->mtype, 810 &mat_mkl_pardiso->phase, 811 &mat_mkl_pardiso->n, 812 mat_mkl_pardiso->a, 813 mat_mkl_pardiso->ia, 814 mat_mkl_pardiso->ja, 815 mat_mkl_pardiso->perm, 816 &mat_mkl_pardiso->nrhs, 817 mat_mkl_pardiso->iparm, 818 &mat_mkl_pardiso->msglvl, 819 NULL, 820 (void*)mat_mkl_pardiso->schur, 821 &mat_mkl_pardiso->err); 822 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); 823 824 if (mat_mkl_pardiso->schur) { /* schur output from pardiso is in row major format */ 825 PetscInt j,k,n=mat_mkl_pardiso->schur_size; 826 if (!mat_mkl_pardiso->schur_solver_type) { 827 for (j=0; j<n; j++) { 828 for (k=0; k<j; k++) { 829 PetscScalar tmp = mat_mkl_pardiso->schur[j + k*n]; 830 mat_mkl_pardiso->schur[j + k*n] = mat_mkl_pardiso->schur[k + j*n]; 831 mat_mkl_pardiso->schur[k + j*n] = tmp; 832 } 833 } 834 } 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 */ 835 for (j=0; j<n; j++) { 836 for (k=0; k<j; k++) { 837 mat_mkl_pardiso->schur[j + k*n] = mat_mkl_pardiso->schur[k + j*n]; 838 } 839 } 840 } 841 } 842 mat_mkl_pardiso->matstruc = SAME_NONZERO_PATTERN; 843 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 844 mat_mkl_pardiso->schur_factored = PETSC_FALSE; 845 mat_mkl_pardiso->schur_inverted = PETSC_FALSE; 846 mat_mkl_pardiso->schur_solver_type = 0; 847 PetscFunctionReturn(0); 848 } 849 850 #undef __FUNCT__ 851 #define __FUNCT__ "PetscSetMKL_PARDISOFromOptions" 852 PetscErrorCode PetscSetMKL_PARDISOFromOptions(Mat F, Mat A) 853 { 854 Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr; 855 PetscErrorCode ierr; 856 PetscInt icntl, threads = 1; 857 PetscBool flg; 858 859 PetscFunctionBegin; 860 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MKL_PARDISO Options","Mat");CHKERRQ(ierr); 861 ierr = PetscOptionsInt("-mat_mkl_pardiso_65","Number of threads to use","None",threads,&threads,&flg);CHKERRQ(ierr); 862 if (flg) mkl_set_num_threads((int)threads); 863 864 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); 865 if (flg) mat_mkl_pardiso->maxfct = icntl; 866 867 ierr = PetscOptionsInt("-mat_mkl_pardiso_67","Indicates the actual matrix for the solution phase","None",mat_mkl_pardiso->mnum,&icntl,&flg);CHKERRQ(ierr); 868 if (flg) mat_mkl_pardiso->mnum = icntl; 869 870 ierr = PetscOptionsInt("-mat_mkl_pardiso_68","Message level information","None",mat_mkl_pardiso->msglvl,&icntl,&flg);CHKERRQ(ierr); 871 if (flg) mat_mkl_pardiso->msglvl = icntl; 872 873 ierr = PetscOptionsInt("-mat_mkl_pardiso_69","Defines the matrix type","None",mat_mkl_pardiso->mtype,&icntl,&flg);CHKERRQ(ierr); 874 if(flg){ 875 void *pt[IPARM_SIZE]; 876 mat_mkl_pardiso->mtype = icntl; 877 MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm); 878 #if defined(PETSC_USE_REAL_SINGLE) 879 mat_mkl_pardiso->iparm[27] = 1; 880 #else 881 mat_mkl_pardiso->iparm[27] = 0; 882 #endif 883 mat_mkl_pardiso->iparm[34] = 1; 884 } 885 ierr = PetscOptionsInt("-mat_mkl_pardiso_1","Use default values","None",mat_mkl_pardiso->iparm[0],&icntl,&flg);CHKERRQ(ierr); 886 887 if(flg && icntl != 0){ 888 ierr = PetscOptionsInt("-mat_mkl_pardiso_2","Fill-in reducing ordering for the input matrix","None",mat_mkl_pardiso->iparm[1],&icntl,&flg);CHKERRQ(ierr); 889 if (flg) mat_mkl_pardiso->iparm[1] = icntl; 890 891 ierr = PetscOptionsInt("-mat_mkl_pardiso_4","Preconditioned CGS/CG","None",mat_mkl_pardiso->iparm[3],&icntl,&flg);CHKERRQ(ierr); 892 if (flg) mat_mkl_pardiso->iparm[3] = icntl; 893 894 ierr = PetscOptionsInt("-mat_mkl_pardiso_5","User permutation","None",mat_mkl_pardiso->iparm[4],&icntl,&flg);CHKERRQ(ierr); 895 if (flg) mat_mkl_pardiso->iparm[4] = icntl; 896 897 ierr = PetscOptionsInt("-mat_mkl_pardiso_6","Write solution on x","None",mat_mkl_pardiso->iparm[5],&icntl,&flg);CHKERRQ(ierr); 898 if (flg) mat_mkl_pardiso->iparm[5] = icntl; 899 900 ierr = PetscOptionsInt("-mat_mkl_pardiso_8","Iterative refinement step","None",mat_mkl_pardiso->iparm[7],&icntl,&flg);CHKERRQ(ierr); 901 if (flg) mat_mkl_pardiso->iparm[7] = icntl; 902 903 ierr = PetscOptionsInt("-mat_mkl_pardiso_10","Pivoting perturbation","None",mat_mkl_pardiso->iparm[9],&icntl,&flg);CHKERRQ(ierr); 904 if (flg) mat_mkl_pardiso->iparm[9] = icntl; 905 906 ierr = PetscOptionsInt("-mat_mkl_pardiso_11","Scaling vectors","None",mat_mkl_pardiso->iparm[10],&icntl,&flg);CHKERRQ(ierr); 907 if (flg) mat_mkl_pardiso->iparm[10] = icntl; 908 909 ierr = PetscOptionsInt("-mat_mkl_pardiso_12","Solve with transposed or conjugate transposed matrix A","None",mat_mkl_pardiso->iparm[11],&icntl,&flg);CHKERRQ(ierr); 910 if (flg) mat_mkl_pardiso->iparm[11] = icntl; 911 912 ierr = PetscOptionsInt("-mat_mkl_pardiso_13","Improved accuracy using (non-) symmetric weighted matching","None",mat_mkl_pardiso->iparm[12],&icntl,&flg);CHKERRQ(ierr); 913 if (flg) mat_mkl_pardiso->iparm[12] = icntl; 914 915 ierr = PetscOptionsInt("-mat_mkl_pardiso_18","Numbers of non-zero elements","None",mat_mkl_pardiso->iparm[17],&icntl,&flg);CHKERRQ(ierr); 916 if (flg) mat_mkl_pardiso->iparm[17] = icntl; 917 918 ierr = PetscOptionsInt("-mat_mkl_pardiso_19","Report number of floating point operations","None",mat_mkl_pardiso->iparm[18],&icntl,&flg);CHKERRQ(ierr); 919 if (flg) mat_mkl_pardiso->iparm[18] = icntl; 920 921 ierr = PetscOptionsInt("-mat_mkl_pardiso_21","Pivoting for symmetric indefinite matrices","None",mat_mkl_pardiso->iparm[20],&icntl,&flg);CHKERRQ(ierr); 922 if (flg) mat_mkl_pardiso->iparm[20] = icntl; 923 924 ierr = PetscOptionsInt("-mat_mkl_pardiso_24","Parallel factorization control","None",mat_mkl_pardiso->iparm[23],&icntl,&flg);CHKERRQ(ierr); 925 if (flg) mat_mkl_pardiso->iparm[23] = icntl; 926 927 ierr = PetscOptionsInt("-mat_mkl_pardiso_25","Parallel forward/backward solve control","None",mat_mkl_pardiso->iparm[24],&icntl,&flg);CHKERRQ(ierr); 928 if (flg) mat_mkl_pardiso->iparm[24] = icntl; 929 930 ierr = PetscOptionsInt("-mat_mkl_pardiso_27","Matrix checker","None",mat_mkl_pardiso->iparm[26],&icntl,&flg);CHKERRQ(ierr); 931 if (flg) mat_mkl_pardiso->iparm[26] = icntl; 932 933 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); 934 if (flg) mat_mkl_pardiso->iparm[30] = icntl; 935 936 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); 937 if (flg) mat_mkl_pardiso->iparm[33] = icntl; 938 939 ierr = PetscOptionsInt("-mat_mkl_pardiso_60","Intel MKL_PARDISO mode","None",mat_mkl_pardiso->iparm[59],&icntl,&flg);CHKERRQ(ierr); 940 if (flg) mat_mkl_pardiso->iparm[59] = icntl; 941 } 942 PetscOptionsEnd(); 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "MatFactorMKL_PARDISOInitialize_Private" 948 PetscErrorCode MatFactorMKL_PARDISOInitialize_Private(Mat A, MatFactorType ftype, Mat_MKL_PARDISO *mat_mkl_pardiso) 949 { 950 PetscErrorCode ierr; 951 PetscInt i; 952 953 PetscFunctionBegin; 954 for ( i = 0; i < IPARM_SIZE; i++ ){ 955 mat_mkl_pardiso->iparm[i] = 0; 956 } 957 for ( i = 0; i < IPARM_SIZE; i++ ){ 958 mat_mkl_pardiso->pt[i] = 0; 959 } 960 /* Default options for both sym and unsym */ 961 mat_mkl_pardiso->iparm[ 0] = 1; /* Solver default parameters overriden with provided by iparm */ 962 mat_mkl_pardiso->iparm[ 1] = 2; /* Metis reordering */ 963 mat_mkl_pardiso->iparm[ 5] = 0; /* Write solution into x */ 964 mat_mkl_pardiso->iparm[ 7] = 0; /* Max number of iterative refinement steps */ 965 mat_mkl_pardiso->iparm[17] = -1; /* Output: Number of nonzeros in the factor LU */ 966 mat_mkl_pardiso->iparm[18] = -1; /* Output: Mflops for LU factorization */ 967 #if 0 968 mat_mkl_pardiso->iparm[23] = 1; /* Parallel factorization control*/ 969 #endif 970 mat_mkl_pardiso->iparm[34] = 1; /* Cluster Sparse Solver use C-style indexing for ia and ja arrays */ 971 mat_mkl_pardiso->iparm[39] = 0; /* Input: matrix/rhs/solution stored on master */ 972 973 mat_mkl_pardiso->CleanUp = PETSC_FALSE; 974 mat_mkl_pardiso->maxfct = 1; /* Maximum number of numerical factorizations. */ 975 mat_mkl_pardiso->mnum = 1; /* Which factorization to use. */ 976 mat_mkl_pardiso->msglvl = 0; /* 0: do not print 1: Print statistical information in file */ 977 mat_mkl_pardiso->phase = -1; 978 mat_mkl_pardiso->err = 0; 979 980 mat_mkl_pardiso->n = A->rmap->N; 981 mat_mkl_pardiso->nrhs = 1; 982 mat_mkl_pardiso->err = 0; 983 mat_mkl_pardiso->phase = -1; 984 985 if(ftype == MAT_FACTOR_LU){ 986 mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */ 987 mat_mkl_pardiso->iparm[10] = 1; /* Use nonsymmetric permutation and scaling MPS */ 988 mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 989 990 } else { 991 mat_mkl_pardiso->iparm[ 9] = 13; /* Perturb the pivot elements with 1E-13 */ 992 mat_mkl_pardiso->iparm[10] = 0; /* Use nonsymmetric permutation and scaling MPS */ 993 mat_mkl_pardiso->iparm[12] = 1; /* Switch on Maximum Weighted Matching algorithm (default for non-symmetric) */ 994 /* mat_mkl_pardiso->iparm[20] = 1; */ /* Apply 1x1 and 2x2 Bunch-Kaufman pivoting during the factorization process */ 995 #if defined(PETSC_USE_DEBUG) 996 mat_mkl_pardiso->iparm[26] = 1; /* Matrix checker */ 997 #endif 998 } 999 ierr = PetscMalloc1(A->rmap->N*sizeof(INT_TYPE), &mat_mkl_pardiso->perm);CHKERRQ(ierr); 1000 for(i = 0; i < A->rmap->N; i++){ 1001 mat_mkl_pardiso->perm[i] = 0; 1002 } 1003 mat_mkl_pardiso->schur_size = 0; 1004 PetscFunctionReturn(0); 1005 } 1006 1007 #undef __FUNCT__ 1008 #define __FUNCT__ "MatFactorSymbolic_AIJMKL_PARDISO_Private" 1009 PetscErrorCode MatFactorSymbolic_AIJMKL_PARDISO_Private(Mat F,Mat A,const MatFactorInfo *info) 1010 { 1011 Mat_MKL_PARDISO *mat_mkl_pardiso = (Mat_MKL_PARDISO*)F->spptr; 1012 PetscErrorCode ierr; 1013 1014 PetscFunctionBegin; 1015 mat_mkl_pardiso->matstruc = DIFFERENT_NONZERO_PATTERN; 1016 ierr = PetscSetMKL_PARDISOFromOptions(F,A);CHKERRQ(ierr); 1017 1018 /* throw away any previously computed structure */ 1019 if (mat_mkl_pardiso->freeaij) { 1020 ierr = PetscFree2(mat_mkl_pardiso->ia,mat_mkl_pardiso->ja);CHKERRQ(ierr); 1021 ierr = PetscFree(mat_mkl_pardiso->a);CHKERRQ(ierr); 1022 } 1023 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); 1024 mat_mkl_pardiso->n = A->rmap->N; 1025 1026 mat_mkl_pardiso->phase = JOB_ANALYSIS; 1027 1028 MKL_PARDISO (mat_mkl_pardiso->pt, 1029 &mat_mkl_pardiso->maxfct, 1030 &mat_mkl_pardiso->mnum, 1031 &mat_mkl_pardiso->mtype, 1032 &mat_mkl_pardiso->phase, 1033 &mat_mkl_pardiso->n, 1034 mat_mkl_pardiso->a, 1035 mat_mkl_pardiso->ia, 1036 mat_mkl_pardiso->ja, 1037 mat_mkl_pardiso->perm, 1038 &mat_mkl_pardiso->nrhs, 1039 mat_mkl_pardiso->iparm, 1040 &mat_mkl_pardiso->msglvl, 1041 NULL, 1042 NULL, 1043 &mat_mkl_pardiso->err); 1044 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); 1045 1046 mat_mkl_pardiso->CleanUp = PETSC_TRUE; 1047 1048 if (F->factortype == MAT_FACTOR_LU) { 1049 F->ops->lufactornumeric = MatFactorNumeric_MKL_PARDISO; 1050 } else { 1051 F->ops->choleskyfactornumeric = MatFactorNumeric_MKL_PARDISO; 1052 } 1053 F->ops->solve = MatSolve_MKL_PARDISO; 1054 F->ops->solvetranspose = MatSolveTranspose_MKL_PARDISO; 1055 F->ops->matsolve = MatMatSolve_MKL_PARDISO; 1056 PetscFunctionReturn(0); 1057 } 1058 1059 #undef __FUNCT__ 1060 #define __FUNCT__ "MatLUFactorSymbolic_AIJMKL_PARDISO" 1061 PetscErrorCode MatLUFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 1062 { 1063 PetscErrorCode ierr; 1064 1065 PetscFunctionBegin; 1066 ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr); 1067 PetscFunctionReturn(0); 1068 } 1069 1070 #undef __FUNCT__ 1071 #define __FUNCT__ "MatCholeskyFactorSymbolic_AIJMKL_PARDISO" 1072 PetscErrorCode MatCholeskyFactorSymbolic_AIJMKL_PARDISO(Mat F,Mat A,IS r,const MatFactorInfo *info) 1073 { 1074 PetscErrorCode ierr; 1075 1076 PetscFunctionBegin; 1077 ierr = MatFactorSymbolic_AIJMKL_PARDISO_Private(F, A, info);CHKERRQ(ierr); 1078 PetscFunctionReturn(0); 1079 } 1080 1081 #undef __FUNCT__ 1082 #define __FUNCT__ "MatView_MKL_PARDISO" 1083 PetscErrorCode MatView_MKL_PARDISO(Mat A, PetscViewer viewer) 1084 { 1085 PetscErrorCode ierr; 1086 PetscBool iascii; 1087 PetscViewerFormat format; 1088 Mat_MKL_PARDISO *mat_mkl_pardiso=(Mat_MKL_PARDISO*)A->spptr; 1089 PetscInt i; 1090 1091 PetscFunctionBegin; 1092 if (A->ops->solve != MatSolve_MKL_PARDISO) PetscFunctionReturn(0); 1093 1094 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1095 if (iascii) { 1096 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1097 if (format == PETSC_VIEWER_ASCII_INFO) { 1098 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO run parameters:\n");CHKERRQ(ierr); 1099 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO phase: %d \n",mat_mkl_pardiso->phase);CHKERRQ(ierr); 1100 for(i = 1; i <= 64; i++){ 1101 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO iparm[%d]: %d \n",i, mat_mkl_pardiso->iparm[i - 1]);CHKERRQ(ierr); 1102 } 1103 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO maxfct: %d \n", mat_mkl_pardiso->maxfct);CHKERRQ(ierr); 1104 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mnum: %d \n", mat_mkl_pardiso->mnum);CHKERRQ(ierr); 1105 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO mtype: %d \n", mat_mkl_pardiso->mtype);CHKERRQ(ierr); 1106 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO n: %d \n", mat_mkl_pardiso->n);CHKERRQ(ierr); 1107 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO nrhs: %d \n", mat_mkl_pardiso->nrhs);CHKERRQ(ierr); 1108 ierr = PetscViewerASCIIPrintf(viewer,"MKL_PARDISO msglvl: %d \n", mat_mkl_pardiso->msglvl);CHKERRQ(ierr); 1109 } 1110 } 1111 PetscFunctionReturn(0); 1112 } 1113 1114 1115 #undef __FUNCT__ 1116 #define __FUNCT__ "MatGetInfo_MKL_PARDISO" 1117 PetscErrorCode MatGetInfo_MKL_PARDISO(Mat A, MatInfoType flag, MatInfo *info) 1118 { 1119 Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)A->spptr; 1120 1121 PetscFunctionBegin; 1122 info->block_size = 1.0; 1123 info->nz_allocated = mat_mkl_pardiso->nz + 0.0; 1124 info->nz_unneeded = 0.0; 1125 info->assemblies = 0.0; 1126 info->mallocs = 0.0; 1127 info->memory = 0.0; 1128 info->fill_ratio_given = 0; 1129 info->fill_ratio_needed = 0; 1130 info->factor_mallocs = 0; 1131 PetscFunctionReturn(0); 1132 } 1133 1134 #undef __FUNCT__ 1135 #define __FUNCT__ "MatMkl_PardisoSetCntl_MKL_PARDISO" 1136 PetscErrorCode MatMkl_PardisoSetCntl_MKL_PARDISO(Mat F,PetscInt icntl,PetscInt ival) 1137 { 1138 Mat_MKL_PARDISO *mat_mkl_pardiso =(Mat_MKL_PARDISO*)F->spptr; 1139 1140 PetscFunctionBegin; 1141 if(icntl <= 64){ 1142 mat_mkl_pardiso->iparm[icntl - 1] = ival; 1143 } else { 1144 if(icntl == 65) 1145 mkl_set_num_threads((int)ival); 1146 else if(icntl == 66) 1147 mat_mkl_pardiso->maxfct = ival; 1148 else if(icntl == 67) 1149 mat_mkl_pardiso->mnum = ival; 1150 else if(icntl == 68) 1151 mat_mkl_pardiso->msglvl = ival; 1152 else if(icntl == 69){ 1153 void *pt[IPARM_SIZE]; 1154 mat_mkl_pardiso->mtype = ival; 1155 MKL_PARDISO_INIT(pt, &mat_mkl_pardiso->mtype, mat_mkl_pardiso->iparm); 1156 #if defined(PETSC_USE_REAL_SINGLE) 1157 mat_mkl_pardiso->iparm[27] = 1; 1158 #else 1159 mat_mkl_pardiso->iparm[27] = 0; 1160 #endif 1161 mat_mkl_pardiso->iparm[34] = 1; 1162 } else if(icntl==70) { 1163 mat_mkl_pardiso->solve_interior = !!ival; 1164 } 1165 } 1166 PetscFunctionReturn(0); 1167 } 1168 1169 #undef __FUNCT__ 1170 #define __FUNCT__ "MatMkl_PardisoSetCntl" 1171 /*@ 1172 MatMkl_PardisoSetCntl - Set Mkl_Pardiso parameters 1173 1174 Logically Collective on Mat 1175 1176 Input Parameters: 1177 + F - the factored matrix obtained by calling MatGetFactor() 1178 . icntl - index of Mkl_Pardiso parameter 1179 - ival - value of Mkl_Pardiso parameter 1180 1181 Options Database: 1182 . -mat_mkl_pardiso_<icntl> <ival> 1183 1184 Level: beginner 1185 1186 References: Mkl_Pardiso Users' Guide 1187 1188 .seealso: MatGetFactor() 1189 @*/ 1190 PetscErrorCode MatMkl_PardisoSetCntl(Mat F,PetscInt icntl,PetscInt ival) 1191 { 1192 PetscErrorCode ierr; 1193 1194 PetscFunctionBegin; 1195 ierr = PetscTryMethod(F,"MatMkl_PardisoSetCntl_C",(Mat,PetscInt,PetscInt),(F,icntl,ival));CHKERRQ(ierr); 1196 PetscFunctionReturn(0); 1197 } 1198 1199 /*MC 1200 MATSOLVERMKL_PARDISO - A matrix type providing direct solvers (LU) for 1201 sequential matrices via the external package MKL_PARDISO. 1202 1203 Works with MATSEQAIJ matrices 1204 1205 Use -pc_type lu -pc_factor_mat_solver_package mkl_pardiso to us this direct solver 1206 1207 Options Database Keys: 1208 + -mat_mkl_pardiso_65 - Number of threads to use 1209 . -mat_mkl_pardiso_66 - Maximum number of factors with identical sparsity structure that must be kept in memory at the same time 1210 . -mat_mkl_pardiso_67 - Indicates the actual matrix for the solution phase 1211 . -mat_mkl_pardiso_68 - Message level information 1212 . -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 1213 . -mat_mkl_pardiso_1 - Use default values 1214 . -mat_mkl_pardiso_2 - Fill-in reducing ordering for the input matrix 1215 . -mat_mkl_pardiso_4 - Preconditioned CGS/CG 1216 . -mat_mkl_pardiso_5 - User permutation 1217 . -mat_mkl_pardiso_6 - Write solution on x 1218 . -mat_mkl_pardiso_8 - Iterative refinement step 1219 . -mat_mkl_pardiso_10 - Pivoting perturbation 1220 . -mat_mkl_pardiso_11 - Scaling vectors 1221 . -mat_mkl_pardiso_12 - Solve with transposed or conjugate transposed matrix A 1222 . -mat_mkl_pardiso_13 - Improved accuracy using (non-) symmetric weighted matching 1223 . -mat_mkl_pardiso_18 - Numbers of non-zero elements 1224 . -mat_mkl_pardiso_19 - Report number of floating point operations 1225 . -mat_mkl_pardiso_21 - Pivoting for symmetric indefinite matrices 1226 . -mat_mkl_pardiso_24 - Parallel factorization control 1227 . -mat_mkl_pardiso_25 - Parallel forward/backward solve control 1228 . -mat_mkl_pardiso_27 - Matrix checker 1229 . -mat_mkl_pardiso_31 - Partial solve and computing selected components of the solution vectors 1230 . -mat_mkl_pardiso_34 - Optimal number of threads for conditional numerical reproducibility (CNR) mode 1231 - -mat_mkl_pardiso_60 - Intel MKL_PARDISO mode 1232 1233 Level: beginner 1234 1235 For more information please check mkl_pardiso manual 1236 1237 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 1238 1239 M*/ 1240 #undef __FUNCT__ 1241 #define __FUNCT__ "MatFactorGetSolverPackage_mkl_pardiso" 1242 static PetscErrorCode MatFactorGetSolverPackage_mkl_pardiso(Mat A, const MatSolverPackage *type) 1243 { 1244 PetscFunctionBegin; 1245 *type = MATSOLVERMKL_PARDISO; 1246 PetscFunctionReturn(0); 1247 } 1248 1249 #undef __FUNCT__ 1250 #define __FUNCT__ "MatGetFactor_aij_mkl_pardiso" 1251 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_mkl_pardiso(Mat A,MatFactorType ftype,Mat *F) 1252 { 1253 Mat B; 1254 PetscErrorCode ierr; 1255 Mat_MKL_PARDISO *mat_mkl_pardiso; 1256 PetscBool isSeqAIJ,isSeqBAIJ,isSeqSBAIJ; 1257 1258 PetscFunctionBegin; 1259 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 1260 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&isSeqBAIJ);CHKERRQ(ierr); 1261 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 1262 ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 1263 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 1264 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 1265 ierr = MatSeqAIJSetPreallocation(B,0,NULL);CHKERRQ(ierr); 1266 ierr = MatSeqBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 1267 ierr = MatSeqSBAIJSetPreallocation(B,A->rmap->bs,0,NULL);CHKERRQ(ierr); 1268 1269 ierr = PetscNewLog(B,&mat_mkl_pardiso);CHKERRQ(ierr); 1270 B->spptr = mat_mkl_pardiso; 1271 1272 ierr = MatFactorMKL_PARDISOInitialize_Private(A, ftype, mat_mkl_pardiso);CHKERRQ(ierr); 1273 if (ftype == MAT_FACTOR_LU) { 1274 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJMKL_PARDISO; 1275 B->factortype = MAT_FACTOR_LU; 1276 mat_mkl_pardiso->needsym = PETSC_FALSE; 1277 if (isSeqAIJ) { 1278 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij; 1279 } else if (isSeqBAIJ) { 1280 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij; 1281 } else if (isSeqSBAIJ) { 1282 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU factor with SEQSBAIJ format! Use MAT_FACTOR_CHOLESKY instead"); 1283 } else { 1284 SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO LU with %s format",((PetscObject)A)->type_name); 1285 } 1286 #if defined(PETSC_USE_COMPLEX) 1287 mat_mkl_pardiso->schur_solver_type = 0; /* use a general solver for the moment */ 1288 mat_mkl_pardiso->mtype = 13; 1289 #else 1290 if (A->structurally_symmetric) { 1291 mat_mkl_pardiso->mtype = 1; 1292 } else { 1293 mat_mkl_pardiso->mtype = 11; 1294 } 1295 #endif 1296 mat_mkl_pardiso->schur_solver_type = 0; 1297 } else { 1298 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_AIJMKL_PARDISO; 1299 B->factortype = MAT_FACTOR_CHOLESKY; 1300 if (isSeqAIJ) { 1301 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqaij; 1302 } else if (isSeqBAIJ) { 1303 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqbaij; 1304 } else if (isSeqSBAIJ) { 1305 mat_mkl_pardiso->Convert = MatMKLPardiso_Convert_seqsbaij; 1306 } else { 1307 SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for PARDISO CHOLESKY with %s format",((PetscObject)A)->type_name); 1308 } 1309 mat_mkl_pardiso->needsym = PETSC_TRUE; 1310 #if defined(PETSC_USE_COMPLEX) 1311 mat_mkl_pardiso->schur_solver_type = 0; /* use a general solver for the moment */ 1312 mat_mkl_pardiso->mtype = 13; 1313 #else 1314 if (A->spd_set && A->spd) { 1315 mat_mkl_pardiso->schur_solver_type = 2; 1316 mat_mkl_pardiso->mtype = 2; 1317 } else { 1318 mat_mkl_pardiso->schur_solver_type = 1; 1319 mat_mkl_pardiso->mtype = -2; 1320 } 1321 #endif 1322 } 1323 mat_mkl_pardiso->Destroy = B->ops->destroy; 1324 B->ops->destroy = MatDestroy_MKL_PARDISO; 1325 B->ops->view = MatView_MKL_PARDISO; 1326 B->factortype = ftype; 1327 B->ops->getinfo = MatGetInfo_MKL_PARDISO; 1328 B->assembled = PETSC_TRUE; 1329 1330 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_mkl_pardiso);CHKERRQ(ierr); 1331 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSetSchurIS_C",MatFactorSetSchurIS_MKL_PARDISO);CHKERRQ(ierr); 1332 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorCreateSchurComplement_C",MatFactorCreateSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1333 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSchurComplement_C",MatFactorGetSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1334 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorInvertSchurComplement_C",MatFactorInvertSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1335 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplement_C",MatFactorSolveSchurComplement_MKL_PARDISO);CHKERRQ(ierr); 1336 ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorSolveSchurComplementTranspose_C",MatFactorSolveSchurComplementTranspose_MKL_PARDISO);CHKERRQ(ierr); 1337 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMkl_PardisoSetCntl_C",MatMkl_PardisoSetCntl_MKL_PARDISO);CHKERRQ(ierr); 1338 1339 *F = B; 1340 PetscFunctionReturn(0); 1341 } 1342 1343 #undef __FUNCT__ 1344 #define __FUNCT__ "MatSolverPackageRegister_MKL_Pardiso" 1345 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_MKL_Pardiso(void) 1346 { 1347 PetscErrorCode ierr; 1348 1349 PetscFunctionBegin; 1350 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_LU,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr); 1351 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr); 1352 ierr = MatSolverPackageRegister(MATSOLVERMKL_PARDISO,MATSEQSBAIJ,MAT_FACTOR_CHOLESKY,MatGetFactor_aij_mkl_pardiso);CHKERRQ(ierr); 1353 PetscFunctionReturn(0); 1354 } 1355 1356