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