1 #define PETSCMAT_DLL 2 3 /* 4 Provides an interface to the PaStiX sparse solver 5 */ 6 #include "../src/mat/impls/aij/seq/aij.h" 7 #include "../src/mat/impls/aij/mpi/mpiaij.h" 8 #include "../src/mat/impls/sbaij/seq/sbaij.h" 9 #include "../src/mat/impls/sbaij/mpi/mpisbaij.h" 10 11 #if defined(PETSC_HAVE_STDLIB_H) 12 #include <stdlib.h> 13 #endif 14 #if defined(PETSC_HAVE_STRING_H) 15 #include <string.h> 16 #endif 17 18 EXTERN_C_BEGIN 19 #include "pastix.h" 20 EXTERN_C_END 21 22 typedef struct Mat_Pastix_ { 23 pastix_data_t *pastix_data; /* Pastix data storage structure */ 24 MatStructure matstruc; 25 PetscInt n; /* Number of columns in the matrix */ 26 PetscInt *colptr; /* Index of first element of each column in row and val */ 27 PetscInt *row; /* Row of each element of the matrix */ 28 PetscScalar *val; /* Value of each element of the matrix */ 29 PetscInt *perm; /* Permutation tabular */ 30 PetscInt *invp; /* Reverse permutation tabular */ 31 PetscScalar *rhs; /* Rhight-hand-side member */ 32 PetscInt rhsnbr; /* Rhight-hand-side number (must be 1) */ 33 PetscInt iparm[64]; /* Integer parameters */ 34 double dparm[64]; /* Floating point parameters */ 35 MPI_Comm pastix_comm; /* PaStiX MPI communicator */ 36 PetscMPIInt commRank; /* MPI rank */ 37 PetscMPIInt commSize; /* MPI communicator size */ 38 PetscTruth CleanUpPastix; /* Boolean indicating if we call PaStiX clean step */ 39 VecScatter scat_rhs; 40 VecScatter scat_sol; 41 Vec b_seq; 42 PetscTruth isAIJ; 43 PetscErrorCode (*MatDestroy)(Mat); 44 } Mat_Pastix; 45 46 EXTERN PetscErrorCode MatDuplicate_Pastix(Mat,MatDuplicateOption,Mat*); 47 48 #undef __FUNCT__ 49 #define __FUNCT__ "MatConvertToCSC" 50 /* 51 convert Petsc seqaij matrix to CSC: colptr[n], row[nz], val[nz] 52 53 input: 54 A - matrix in seqaij or mpisbaij (bs=1) format 55 valOnly - FALSE: spaces are allocated and values are set for the CSC 56 TRUE: Only fill values 57 output: 58 n - Size of the matrix 59 colptr - Index of first element of each column in row and val 60 row - Row of each element of the matrix 61 values - Value of each element of the matrix 62 */ 63 PetscErrorCode MatConvertToCSC(Mat A,PetscTruth valOnly,PetscInt *n,PetscInt **colptr,PetscInt **row,PetscScalar **values) 64 { 65 Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 66 PetscInt *rowptr = aa->i; 67 PetscInt *col = aa->j; 68 PetscScalar *rvalues = aa->a; 69 PetscInt m = A->rmap->N; 70 PetscInt nnz; 71 PetscInt i,j, k; 72 PetscInt base = 1; 73 PetscInt idx; 74 PetscErrorCode ierr; 75 PetscInt colidx; 76 PetscInt *colcount; 77 PetscTruth isSBAIJ; 78 PetscTruth isSeqSBAIJ; 79 PetscTruth isMpiSBAIJ; 80 PetscTruth isSym; 81 82 PetscFunctionBegin; 83 ierr = MatIsSymmetric(A,0.0,&isSym);CHKERRQ(ierr); 84 ierr = PetscTypeCompare((PetscObject)A,MATSBAIJ,&isSBAIJ);CHKERRQ(ierr); 85 ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 86 ierr = PetscTypeCompare((PetscObject)A,MATMPISBAIJ,&isMpiSBAIJ);CHKERRQ(ierr); 87 88 *n = A->cmap->N; 89 90 /* PaStiX only needs triangular matrix if matrix is symmetric 91 */ 92 if (isSym && !(isSBAIJ || isSeqSBAIJ || isMpiSBAIJ)) { 93 nnz = (aa->nz - *n)/2 + *n; 94 } 95 else { 96 nnz = aa->nz; 97 } 98 99 if (!valOnly){ 100 ierr = PetscMalloc(((*n)+1) *sizeof(PetscInt) ,colptr );CHKERRQ(ierr); 101 ierr = PetscMalloc( nnz *sizeof(PetscInt) ,row);CHKERRQ(ierr); 102 ierr = PetscMalloc( nnz *sizeof(PetscScalar),values);CHKERRQ(ierr); 103 104 if (isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) { 105 ierr = PetscMemcpy (*colptr, rowptr, ((*n)+1)*sizeof(PetscInt));CHKERRQ(ierr); 106 for (i = 0; i < *n+1; i++) 107 (*colptr)[i] += base; 108 ierr = PetscMemcpy (*row, col, (nnz)*sizeof(PetscInt));CHKERRQ(ierr); 109 for (i = 0; i < nnz; i++) 110 (*row)[i] += base; 111 ierr = PetscMemcpy (*values, rvalues, (nnz)*sizeof(PetscScalar));CHKERRQ(ierr); 112 } else { 113 ierr = PetscMalloc((*n)*sizeof(PetscInt) ,&colcount);CHKERRQ(ierr); 114 115 for (i = 0; i < m; i++) colcount[i] = 0; 116 /* Fill-in colptr */ 117 for (i = 0; i < m; i++) { 118 for (j = rowptr[i]; j < rowptr[i+1]; j++) { 119 if (!isSym || col[j] <= i) colcount[col[j]]++; 120 } 121 } 122 123 (*colptr)[0] = base; 124 for (j = 0; j < *n; j++) { 125 (*colptr)[j+1] = (*colptr)[j] + colcount[j]; 126 /* in next loop we fill starting from (*colptr)[colidx] - base */ 127 colcount[j] = -base; 128 } 129 130 /* Fill-in rows and values */ 131 for (i = 0; i < m; i++) { 132 for (j = rowptr[i]; j < rowptr[i+1]; j++) { 133 if (!isSym || col[j] <= i) { 134 colidx = col[j]; 135 idx = (*colptr)[colidx] + colcount[colidx]; 136 (*row)[idx] = i + base; 137 (*values)[idx] = rvalues[j]; 138 colcount[colidx]++; 139 } 140 } 141 } 142 ierr = PetscFree(colcount);CHKERRQ(ierr); 143 } 144 } else { 145 /* Fill-in only values */ 146 for (i = 0; i < m; i++) { 147 for (j = rowptr[i]; j < rowptr[i+1]; j++) { 148 colidx = col[j]; 149 if ((isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) ||!isSym || col[j] <= i) 150 { 151 /* look for the value to fill */ 152 for (k = (*colptr)[colidx] - base; k < (*colptr)[colidx + 1] - base; k++) { 153 if (((*row)[k]-base) == i) { 154 (*values)[k] = rvalues[j]; 155 break; 156 } 157 } 158 /* data structure of sparse matrix has changed */ 159 if (k == (*colptr)[colidx + 1] - base) SETERRQ1(PETSC_ERR_PLIB,"overflow on k %D",k); 160 } 161 } 162 } 163 } 164 { 165 PetscScalar *tmpvalues; 166 PetscInt *tmprows,*tmpcolptr; 167 ierr = PetscMalloc3(nnz,PetscScalar,&tmpvalues,nnz,PetscInt,&tmprows,(*n+1),PetscInt,&tmpcolptr);CHKERRQ(ierr); 168 if (sizeof(PetscScalar) != sizeof(pastix_float_t)) { 169 SETERRQ2(PETSC_ERR_SUP,"sizeof(PetscScalar) %d != sizeof(pastix_float_t) %d",sizeof(PetscScalar),sizeof(pastix_float_t)); 170 } 171 if (sizeof(PetscInt) != sizeof(pastix_int_t)) { 172 SETERRQ2(PETSC_ERR_SUP,"sizeof(PetscInt) %d != sizeof(pastix_int_t) %d",sizeof(PetscInt),sizeof(pastix_int_t)); 173 } 174 175 ierr = PetscMemcpy(tmpcolptr,*colptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr); 176 ierr = PetscMemcpy(tmprows,*row,nnz*sizeof(PetscInt));CHKERRQ(ierr); 177 ierr = PetscMemcpy(tmpvalues,*values,nnz*sizeof(PetscScalar));CHKERRQ(ierr); 178 ierr = PetscFree(*row);CHKERRQ(ierr); 179 ierr = PetscFree(*values);CHKERRQ(ierr); 180 181 pastix_checkMatrix(MPI_COMM_WORLD,API_VERBOSE_NO,((isSym != 0) ? API_SYM_YES : API_SYM_NO),API_YES,*n,&tmpcolptr,&tmprows,&tmpvalues,NULL,1); 182 183 ierr = PetscMemcpy(*colptr,tmpcolptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr); 184 ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscInt),row);CHKERRQ(ierr); 185 ierr = PetscMemcpy(*row,tmprows,((*colptr)[*n]-1)*sizeof(PetscInt));CHKERRQ(ierr); 186 ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscScalar),values);CHKERRQ(ierr); 187 ierr = PetscMemcpy(*values,tmpvalues,((*colptr)[*n]-1)*sizeof(PetscScalar));CHKERRQ(ierr); 188 ierr = PetscFree3(tmpvalues,tmprows,tmpcolptr);CHKERRQ(ierr); 189 } 190 PetscFunctionReturn(0); 191 } 192 193 194 195 #undef __FUNCT__ 196 #define __FUNCT__ "MatDestroy_Pastix" 197 /* 198 Call clean step of PaStiX if lu->CleanUpPastix == true. 199 Free the CSC matrix. 200 */ 201 PetscErrorCode MatDestroy_Pastix(Mat A) 202 { 203 Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 204 PetscErrorCode ierr; 205 PetscMPIInt size=lu->commSize; 206 207 PetscFunctionBegin; 208 if (lu->CleanUpPastix) { 209 /* Terminate instance, deallocate memories */ 210 if (size > 1){ 211 ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 212 ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 213 ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 214 } 215 216 lu->iparm[IPARM_START_TASK]=API_TASK_CLEAN; 217 lu->iparm[IPARM_END_TASK]=API_TASK_CLEAN; 218 219 pastix((pastix_data_t **)&(lu->pastix_data), 220 lu->pastix_comm, 221 (pastix_int_t) lu->n, 222 (pastix_int_t*) lu->colptr, 223 (pastix_int_t*) lu->row, 224 (pastix_float_t*) lu->val, 225 (pastix_int_t*) lu->perm, 226 (pastix_int_t*) lu->invp, 227 (pastix_float_t*) lu->rhs, 228 (pastix_int_t) lu->rhsnbr, 229 (pastix_int_t*) lu->iparm, 230 lu->dparm); 231 232 ierr = PetscFree(lu->colptr);CHKERRQ(ierr); 233 ierr = PetscFree(lu->row); CHKERRQ(ierr); 234 ierr = PetscFree(lu->val); CHKERRQ(ierr); 235 ierr = PetscFree(lu->perm); CHKERRQ(ierr); 236 ierr = PetscFree(lu->invp); CHKERRQ(ierr); 237 ierr = MPI_Comm_free(&(lu->pastix_comm));CHKERRQ(ierr); 238 } 239 ierr = (lu->MatDestroy)(A);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 #undef __FUNCT__ 244 #define __FUNCT__ "MatSolve_PaStiX" 245 /* 246 Gather right-hand-side. 247 Call for Solve step. 248 Scatter solution. 249 */ 250 PetscErrorCode MatSolve_PaStiX(Mat A,Vec b,Vec x) 251 { 252 Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 253 PetscScalar *array; 254 Vec x_seq; 255 PetscErrorCode ierr; 256 257 PetscFunctionBegin; 258 lu->rhsnbr = 1; 259 x_seq = lu->b_seq; 260 if (lu->commSize > 1){ 261 /* PaStiX only supports centralized rhs. Scatter b into a seqential rhs vector */ 262 ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 263 ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 264 ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr); 265 } else { /* size == 1 */ 266 ierr = VecCopy(b,x);CHKERRQ(ierr); 267 ierr = VecGetArray(x,&array);CHKERRQ(ierr); 268 } 269 lu->rhs = array; 270 if (lu->commSize == 1){ 271 ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 272 } else { 273 ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 274 } 275 276 /* solve phase */ 277 /*-------------*/ 278 lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE; 279 lu->iparm[IPARM_END_TASK] = API_TASK_REFINE; 280 lu->iparm[IPARM_RHS_MAKING] = API_RHS_B; 281 282 pastix((pastix_data_t **)&(lu->pastix_data), 283 (MPI_Comm) lu->pastix_comm, 284 (pastix_int_t) lu->n, 285 (pastix_int_t*) lu->colptr, 286 (pastix_int_t*) lu->row, 287 (pastix_float_t*) lu->val, 288 (pastix_int_t*) lu->perm, 289 (pastix_int_t*) lu->invp, 290 (pastix_float_t*) lu->rhs, 291 (pastix_int_t) lu->rhsnbr, 292 (pastix_int_t*) lu->iparm, 293 (double*) lu->dparm); 294 295 if (lu->iparm[IPARM_ERROR_NUMBER] < 0) { 296 SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in solve phase: lu->iparm[IPARM_ERROR_NUMBER] = %d\n",lu->iparm[IPARM_ERROR_NUMBER] ); 297 } 298 299 if (lu->commSize == 1){ 300 ierr = VecRestoreArray(x,&(lu->rhs));CHKERRQ(ierr); 301 } else { 302 ierr = VecRestoreArray(x_seq,&(lu->rhs));CHKERRQ(ierr); 303 } 304 305 if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */ 306 ierr = VecScatterBegin(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 307 ierr = VecScatterEnd(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 308 } 309 PetscFunctionReturn(0); 310 } 311 312 /* 313 Numeric factorisation using PaStiX solver. 314 315 */ 316 #undef __FUNCT__ 317 #define __FUNCT__ "MatFactorNumeric_PASTIX" 318 PetscErrorCode MatFactorNumeric_PaStiX(Mat F,Mat A,const MatFactorInfo *info) 319 { 320 Mat_Pastix *lu =(Mat_Pastix*)(F)->spptr; 321 Mat *tseq; 322 PetscErrorCode ierr = 0; 323 PetscInt icntl; 324 PetscInt M=A->rmap->N; 325 PetscTruth valOnly,flg, isSym; 326 Mat F_diag; 327 IS is_iden; 328 Vec b; 329 IS isrow; 330 PetscTruth isSeqAIJ,isSeqSBAIJ; 331 332 PetscFunctionBegin; 333 334 ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 335 ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 336 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 337 (F)->ops->solve = MatSolve_PaStiX; 338 339 /* Initialize a PASTIX instance */ 340 ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr); 341 ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank); CHKERRQ(ierr); 342 ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize); CHKERRQ(ierr); 343 344 /* Set pastix options */ 345 lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO; 346 lu->iparm[IPARM_START_TASK] = API_TASK_INIT; 347 lu->iparm[IPARM_END_TASK] = API_TASK_INIT; 348 lu->rhsnbr = 1; 349 350 /* Call to set default pastix options */ 351 pastix((pastix_data_t **)&(lu->pastix_data), 352 (MPI_Comm) lu->pastix_comm, 353 (pastix_int_t) lu->n, 354 (pastix_int_t*) lu->colptr, 355 (pastix_int_t*) lu->row, 356 (pastix_float_t*) lu->val, 357 (pastix_int_t*) lu->perm, 358 (pastix_int_t*) lu->invp, 359 (pastix_float_t*) lu->rhs, 360 (pastix_int_t) lu->rhsnbr, 361 (pastix_int_t*) lu->iparm, 362 (double*) lu->dparm); 363 364 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr); 365 366 icntl=-1; 367 lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT; 368 ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr); 369 if ((flg && icntl > 0) || PetscLogPrintInfo) { 370 lu->iparm[IPARM_VERBOSE] = icntl; 371 } 372 icntl=-1; 373 ierr = PetscOptionsInt("-mat_pastix_threadnbr","iparm[IPARM_THREAD_NBR] : Number of thread by MPI node","None",lu->iparm[IPARM_THREAD_NBR],&icntl,PETSC_NULL);CHKERRQ(ierr); 374 if ((flg && icntl > 0)) { 375 lu->iparm[IPARM_THREAD_NBR] = icntl; 376 } 377 PetscOptionsEnd(); 378 valOnly = PETSC_FALSE; 379 } else { 380 valOnly = PETSC_TRUE; 381 } 382 383 lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES; 384 385 /* convert mpi A to seq mat A */ 386 ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 387 ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 388 ierr = ISDestroy(isrow);CHKERRQ(ierr); 389 390 ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val);CHKERRQ(ierr); 391 ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr); 392 ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr); 393 394 ierr = PetscMalloc((lu->n)*sizeof(PetscInt) ,&(lu->perm));CHKERRQ(ierr); 395 ierr = PetscMalloc((lu->n)*sizeof(PetscInt) ,&(lu->invp));CHKERRQ(ierr); 396 397 if (isSym) { 398 /* On symmetric matrix, LLT */ 399 lu->iparm[IPARM_SYM] = API_SYM_YES; 400 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT; 401 } else { 402 /* On unsymmetric matrix, LU */ 403 lu->iparm[IPARM_SYM] = API_SYM_NO; 404 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 405 } 406 407 /*----------------*/ 408 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 409 if (!(isSeqAIJ || isSeqSBAIJ)) { 410 /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 411 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 412 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 413 ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 414 ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 415 ierr = VecSetFromOptions(b);CHKERRQ(ierr); 416 417 ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 418 ierr = VecScatterCreate(lu->b_seq,is_iden,b,is_iden,&lu->scat_sol);CHKERRQ(ierr); 419 ierr = ISDestroy(is_iden);CHKERRQ(ierr); 420 ierr = VecDestroy(b);CHKERRQ(ierr); 421 } 422 lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING; 423 lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 424 425 pastix((pastix_data_t **)&(lu->pastix_data), 426 (MPI_Comm) lu->pastix_comm, 427 (pastix_int_t) lu->n, 428 (pastix_int_t*) lu->colptr, 429 (pastix_int_t*) lu->row, 430 (pastix_float_t*) lu->val, 431 (pastix_int_t*) lu->perm, 432 (pastix_int_t*) lu->invp, 433 (pastix_float_t*) lu->rhs, 434 (pastix_int_t) lu->rhsnbr, 435 (pastix_int_t*) lu->iparm, 436 (double*) lu->dparm); 437 if (lu->iparm[IPARM_ERROR_NUMBER] < 0) { 438 SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%d\n",lu->iparm[IPARM_ERROR_NUMBER]); 439 } 440 } else { 441 lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT; 442 lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 443 pastix((pastix_data_t **)&(lu->pastix_data), 444 (MPI_Comm) lu->pastix_comm, 445 (pastix_int_t) lu->n, 446 (pastix_int_t*) lu->colptr, 447 (pastix_int_t*) lu->row, 448 (pastix_float_t*) lu->val, 449 (pastix_int_t*) lu->perm, 450 (pastix_int_t*) lu->invp, 451 (pastix_float_t*) lu->rhs, 452 (pastix_int_t) lu->rhsnbr, 453 (pastix_int_t*) lu->iparm, 454 (double*) lu->dparm); 455 456 if (lu->iparm[IPARM_ERROR_NUMBER] < 0) { 457 SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%d\n",lu->iparm[IPARM_ERROR_NUMBER]); 458 } 459 } 460 461 if (lu->commSize > 1){ 462 if ((F)->factor == MAT_FACTOR_LU){ 463 F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 464 } else { 465 F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 466 } 467 F_diag->assembled = PETSC_TRUE; 468 } 469 (F)->assembled = PETSC_TRUE; 470 lu->matstruc = SAME_NONZERO_PATTERN; 471 lu->CleanUpPastix = PETSC_TRUE; 472 PetscFunctionReturn(0); 473 } 474 475 /* Note the Petsc r and c permutations are ignored */ 476 #undef __FUNCT__ 477 #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX" 478 PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 479 { 480 Mat_Pastix *lu = (Mat_Pastix*)F->spptr; 481 482 PetscFunctionBegin; 483 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 484 lu->iparm[IPARM_SYM] = API_SYM_YES; 485 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 486 F->ops->lufactornumeric = MatFactorNumeric_PaStiX; 487 PetscFunctionReturn(0); 488 } 489 490 491 /* Note the Petsc r permutation is ignored */ 492 #undef __FUNCT__ 493 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX" 494 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info) 495 { 496 Mat_Pastix *lu = (Mat_Pastix*)(F)->spptr; 497 498 PetscFunctionBegin; 499 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LLT; 500 lu->iparm[IPARM_SYM] = API_SYM_NO; 501 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 502 (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX; 503 PetscFunctionReturn(0); 504 } 505 506 #undef __FUNCT__ 507 #define __FUNCT__ "MatView_PaStiX" 508 PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer) 509 { 510 PetscErrorCode ierr; 511 PetscTruth iascii; 512 PetscViewerFormat format; 513 514 PetscFunctionBegin; 515 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 516 if (iascii) { 517 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 518 if (format == PETSC_VIEWER_ASCII_INFO){ 519 Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 520 521 ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr); 522 ierr = PetscViewerASCIIPrintf(viewer," Matrix type : %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES)?"Symmetric":"Unsymmetric"));CHKERRQ(ierr); 523 ierr = PetscViewerASCIIPrintf(viewer," Level of printing (0,1,2): %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr); 524 ierr = PetscViewerASCIIPrintf(viewer," Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr); 525 ierr = PetscPrintf(PETSC_COMM_SELF," Error : %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr); 526 } 527 } 528 PetscFunctionReturn(0); 529 } 530 531 532 /*MC 533 MAT_SOLVER_PASTIX - A solver package providing direct solvers (LU) for distributed 534 and sequential matrices via the external package PaStiX. 535 536 Use config/configure.py --download-pastix to have PETSc installed with PaStiX 537 538 Options Database Keys: 539 + -mat_pastix_verbose <0,1,2> - print level 540 - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task. 541 542 Level: beginner 543 544 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 545 546 M*/ 547 548 549 #undef __FUNCT__ 550 #define __FUNCT__ "MatGetInfo_PaStiX" 551 PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info) 552 { 553 Mat_Pastix *lu =(Mat_Pastix*)A->spptr; 554 555 PetscFunctionBegin; 556 info->block_size = 1.0; 557 info->nz_allocated = lu->iparm[IPARM_NNZEROS]; 558 info->nz_used = lu->iparm[IPARM_NNZEROS]; 559 info->nz_unneeded = 0.0; 560 info->assemblies = 0.0; 561 info->mallocs = 0.0; 562 info->memory = 0.0; 563 info->fill_ratio_given = 0; 564 info->fill_ratio_needed = 0; 565 info->factor_mallocs = 0; 566 PetscFunctionReturn(0); 567 } 568 569 EXTERN_C_BEGIN 570 #undef __FUNCT__ 571 #define __FUNCT__ "MatFactorGetSolverPackage_pastix" 572 PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type) 573 { 574 PetscFunctionBegin; 575 *type = MAT_SOLVER_PASTIX; 576 PetscFunctionReturn(0); 577 } 578 EXTERN_C_END 579 580 EXTERN_C_BEGIN 581 /* 582 The seq and mpi versions of this function are the same 583 */ 584 #undef __FUNCT__ 585 #define __FUNCT__ "MatGetFactor_seqaij_pastix" 586 PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F) 587 { 588 Mat B; 589 PetscErrorCode ierr; 590 Mat_Pastix *pastix; 591 592 PetscFunctionBegin; 593 if (ftype != MAT_FACTOR_LU) { 594 SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 595 } 596 /* Create the factorization matrix */ 597 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 598 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 599 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 600 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 601 602 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 603 B->ops->view = MatView_PaStiX; 604 B->ops->getinfo = MatGetInfo_PaStiX; 605 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix", MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 606 B->factor = MAT_FACTOR_LU; 607 608 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 609 pastix->CleanUpPastix = PETSC_FALSE; 610 pastix->isAIJ = PETSC_TRUE; 611 pastix->scat_rhs = PETSC_NULL; 612 pastix->scat_sol = PETSC_NULL; 613 pastix->MatDestroy = B->ops->destroy; 614 B->ops->destroy = MatDestroy_Pastix; 615 B->spptr = (void*)pastix; 616 617 *F = B; 618 PetscFunctionReturn(0); 619 } 620 EXTERN_C_END 621 622 623 EXTERN_C_BEGIN 624 #undef __FUNCT__ 625 #define __FUNCT__ "MatGetFactor_mpiaij_pastix" 626 PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F) 627 { 628 Mat B; 629 PetscErrorCode ierr; 630 Mat_Pastix *pastix; 631 632 PetscFunctionBegin; 633 if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 634 /* Create the factorization matrix */ 635 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 636 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 637 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 638 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 639 ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 640 641 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 642 B->ops->view = MatView_PaStiX; 643 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 644 B->factor = MAT_FACTOR_LU; 645 646 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 647 pastix->CleanUpPastix = PETSC_FALSE; 648 pastix->isAIJ = PETSC_TRUE; 649 pastix->scat_rhs = PETSC_NULL; 650 pastix->scat_sol = PETSC_NULL; 651 pastix->MatDestroy = B->ops->destroy; 652 B->ops->destroy = MatDestroy_Pastix; 653 B->spptr = (void*)pastix; 654 655 *F = B; 656 PetscFunctionReturn(0); 657 } 658 EXTERN_C_END 659 660 EXTERN_C_BEGIN 661 #undef __FUNCT__ 662 #define __FUNCT__ "MatGetFactor_seqsbaij_pastix" 663 PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 664 { 665 Mat B; 666 PetscErrorCode ierr; 667 Mat_Pastix *pastix; 668 669 PetscFunctionBegin; 670 if (ftype != MAT_FACTOR_CHOLESKY) { 671 SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 672 } 673 /* Create the factorization matrix */ 674 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 675 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 676 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 677 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 678 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 679 680 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 681 B->ops->view = MatView_PaStiX; 682 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 683 684 B->factor = MAT_FACTOR_CHOLESKY; 685 686 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 687 pastix->CleanUpPastix = PETSC_FALSE; 688 pastix->isAIJ = PETSC_TRUE; 689 pastix->scat_rhs = PETSC_NULL; 690 pastix->scat_sol = PETSC_NULL; 691 pastix->MatDestroy = B->ops->destroy; 692 B->ops->destroy = MatDestroy_Pastix; 693 B->spptr = (void*)pastix; 694 695 *F = B; 696 PetscFunctionReturn(0); 697 } 698 EXTERN_C_END 699 700 EXTERN_C_BEGIN 701 #undef __FUNCT__ 702 #define __FUNCT__ "MatGetFactor_mpisbaij_pastix" 703 PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 704 { 705 Mat B; 706 PetscErrorCode ierr; 707 Mat_Pastix *pastix; 708 709 PetscFunctionBegin; 710 if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 711 712 /* Create the factorization matrix */ 713 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 714 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 715 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 716 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 717 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 718 719 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 720 B->ops->view = MatView_PaStiX; 721 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 722 B->factor = MAT_FACTOR_CHOLESKY; 723 724 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 725 pastix->CleanUpPastix = PETSC_FALSE; 726 pastix->isAIJ = PETSC_TRUE; 727 pastix->scat_rhs = PETSC_NULL; 728 pastix->scat_sol = PETSC_NULL; 729 pastix->MatDestroy = B->ops->destroy; 730 B->ops->destroy = MatDestroy_Pastix; 731 B->spptr = (void*)pastix; 732 733 *F = B; 734 PetscFunctionReturn(0); 735 } 736 EXTERN_C_END 737