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