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