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 PetscBool CleanUpPastix; /* Boolean indicating if we call PaStiX clean step */ 39 VecScatter scat_rhs; 40 VecScatter scat_sol; 41 Vec b_seq; 42 PetscBool 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,PetscBool 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 PetscBool isSBAIJ; 78 PetscBool isSeqSBAIJ; 79 PetscBool isMpiSBAIJ; 80 PetscBool 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_COMM_SELF,PETSC_ERR_PLIB,"overflow on k %D",k); 160 } 161 } 162 } 163 } 164 { 165 PetscScalar *tmpvalues; 166 PetscInt *tmprows,*tmpcolptr; 167 tmpvalues = malloc(nnz*sizeof(PetscScalar)); if (!tmpvalues) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory"); 168 tmprows = malloc(nnz*sizeof(PetscInt));if (!tmprows) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory"); 169 tmpcolptr = malloc((*n+1)*sizeof(PetscInt));if (!tmpcolptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory"); 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; 327 328 PetscFunctionBegin; 329 330 ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 331 ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 332 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 333 (F)->ops->solve = MatSolve_PaStiX; 334 335 /* Initialize a PASTIX instance */ 336 ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr); 337 ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank); CHKERRQ(ierr); 338 ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize); CHKERRQ(ierr); 339 340 /* Set pastix options */ 341 lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO; 342 lu->iparm[IPARM_START_TASK] = API_TASK_INIT; 343 lu->iparm[IPARM_END_TASK] = API_TASK_INIT; 344 lu->rhsnbr = 1; 345 346 /* Call to set default pastix options */ 347 pastix((pastix_data_t **)&(lu->pastix_data), 348 (MPI_Comm) lu->pastix_comm, 349 (pastix_int_t) lu->n, 350 (pastix_int_t*) lu->colptr, 351 (pastix_int_t*) lu->row, 352 (pastix_float_t*) lu->val, 353 (pastix_int_t*) lu->perm, 354 (pastix_int_t*) lu->invp, 355 (pastix_float_t*) lu->rhs, 356 (pastix_int_t) lu->rhsnbr, 357 (pastix_int_t*) lu->iparm, 358 (double*) lu->dparm); 359 360 ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr); 361 362 icntl=-1; 363 lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT; 364 ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr); 365 if ((flg && icntl > 0) || PetscLogPrintInfo) { 366 lu->iparm[IPARM_VERBOSE] = icntl; 367 } 368 icntl=-1; 369 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); 370 if ((flg && icntl > 0)) { 371 lu->iparm[IPARM_THREAD_NBR] = icntl; 372 } 373 PetscOptionsEnd(); 374 valOnly = PETSC_FALSE; 375 } else { 376 valOnly = PETSC_TRUE; 377 } 378 379 lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES; 380 381 /* convert mpi A to seq mat A */ 382 ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 383 ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 384 ierr = ISDestroy(isrow);CHKERRQ(ierr); 385 386 ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val);CHKERRQ(ierr); 387 ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr); 388 ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr); 389 390 ierr = PetscMalloc((lu->n)*sizeof(PetscInt) ,&(lu->perm));CHKERRQ(ierr); 391 ierr = PetscMalloc((lu->n)*sizeof(PetscInt) ,&(lu->invp));CHKERRQ(ierr); 392 393 if (isSym) { 394 /* On symmetric matrix, LLT */ 395 lu->iparm[IPARM_SYM] = API_SYM_YES; 396 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT; 397 } else { 398 /* On unsymmetric matrix, LU */ 399 lu->iparm[IPARM_SYM] = API_SYM_NO; 400 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 401 } 402 403 /*----------------*/ 404 if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 405 if (!(isSeqAIJ || isSeqSBAIJ)) { 406 /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 407 ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 408 ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 409 ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 410 ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 411 ierr = VecSetFromOptions(b);CHKERRQ(ierr); 412 413 ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 414 ierr = VecScatterCreate(lu->b_seq,is_iden,b,is_iden,&lu->scat_sol);CHKERRQ(ierr); 415 ierr = ISDestroy(is_iden);CHKERRQ(ierr); 416 ierr = VecDestroy(b);CHKERRQ(ierr); 417 } 418 lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING; 419 lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 420 421 pastix((pastix_data_t **)&(lu->pastix_data), 422 (MPI_Comm) lu->pastix_comm, 423 (pastix_int_t) lu->n, 424 (pastix_int_t*) lu->colptr, 425 (pastix_int_t*) lu->row, 426 (pastix_float_t*) lu->val, 427 (pastix_int_t*) lu->perm, 428 (pastix_int_t*) lu->invp, 429 (pastix_float_t*) lu->rhs, 430 (pastix_int_t) lu->rhsnbr, 431 (pastix_int_t*) lu->iparm, 432 (double*) lu->dparm); 433 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]); 434 } else { 435 lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT; 436 lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 437 pastix((pastix_data_t **)&(lu->pastix_data), 438 (MPI_Comm) lu->pastix_comm, 439 (pastix_int_t) lu->n, 440 (pastix_int_t*) lu->colptr, 441 (pastix_int_t*) lu->row, 442 (pastix_float_t*) lu->val, 443 (pastix_int_t*) lu->perm, 444 (pastix_int_t*) lu->invp, 445 (pastix_float_t*) lu->rhs, 446 (pastix_int_t) lu->rhsnbr, 447 (pastix_int_t*) lu->iparm, 448 (double*) lu->dparm); 449 450 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]); 451 } 452 453 if (lu->commSize > 1){ 454 if ((F)->factortype == MAT_FACTOR_LU){ 455 F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 456 } else { 457 F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 458 } 459 F_diag->assembled = PETSC_TRUE; 460 } 461 (F)->assembled = PETSC_TRUE; 462 lu->matstruc = SAME_NONZERO_PATTERN; 463 lu->CleanUpPastix = PETSC_TRUE; 464 PetscFunctionReturn(0); 465 } 466 467 /* Note the Petsc r and c permutations are ignored */ 468 #undef __FUNCT__ 469 #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX" 470 PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 471 { 472 Mat_Pastix *lu = (Mat_Pastix*)F->spptr; 473 474 PetscFunctionBegin; 475 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 476 lu->iparm[IPARM_SYM] = API_SYM_YES; 477 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 478 F->ops->lufactornumeric = MatFactorNumeric_PaStiX; 479 PetscFunctionReturn(0); 480 } 481 482 483 /* Note the Petsc r permutation is ignored */ 484 #undef __FUNCT__ 485 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX" 486 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info) 487 { 488 Mat_Pastix *lu = (Mat_Pastix*)(F)->spptr; 489 490 PetscFunctionBegin; 491 lu->iparm[IPARM_FACTORIZATION] = API_FACT_LLT; 492 lu->iparm[IPARM_SYM] = API_SYM_NO; 493 lu->matstruc = DIFFERENT_NONZERO_PATTERN; 494 (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX; 495 PetscFunctionReturn(0); 496 } 497 498 #undef __FUNCT__ 499 #define __FUNCT__ "MatView_PaStiX" 500 PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer) 501 { 502 PetscErrorCode ierr; 503 PetscBool iascii; 504 PetscViewerFormat format; 505 506 PetscFunctionBegin; 507 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 508 if (iascii) { 509 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 510 if (format == PETSC_VIEWER_ASCII_INFO){ 511 Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 512 513 ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr); 514 ierr = PetscViewerASCIIPrintf(viewer," Matrix type : %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES)?"Symmetric":"Unsymmetric"));CHKERRQ(ierr); 515 ierr = PetscViewerASCIIPrintf(viewer," Level of printing (0,1,2): %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr); 516 ierr = PetscViewerASCIIPrintf(viewer," Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr); 517 ierr = PetscPrintf(PETSC_COMM_SELF," Error : %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr); 518 } 519 } 520 PetscFunctionReturn(0); 521 } 522 523 524 /*MC 525 MATSOLVERPASTIX - A solver package providing direct solvers (LU) for distributed 526 and sequential matrices via the external package PaStiX. 527 528 Use ./configure --download-pastix to have PETSc installed with PaStiX 529 530 Options Database Keys: 531 + -mat_pastix_verbose <0,1,2> - print level 532 - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task. 533 534 Level: beginner 535 536 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 537 538 M*/ 539 540 541 #undef __FUNCT__ 542 #define __FUNCT__ "MatGetInfo_PaStiX" 543 PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info) 544 { 545 Mat_Pastix *lu =(Mat_Pastix*)A->spptr; 546 547 PetscFunctionBegin; 548 info->block_size = 1.0; 549 info->nz_allocated = lu->iparm[IPARM_NNZEROS]; 550 info->nz_used = lu->iparm[IPARM_NNZEROS]; 551 info->nz_unneeded = 0.0; 552 info->assemblies = 0.0; 553 info->mallocs = 0.0; 554 info->memory = 0.0; 555 info->fill_ratio_given = 0; 556 info->fill_ratio_needed = 0; 557 info->factor_mallocs = 0; 558 PetscFunctionReturn(0); 559 } 560 561 EXTERN_C_BEGIN 562 #undef __FUNCT__ 563 #define __FUNCT__ "MatFactorGetSolverPackage_pastix" 564 PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type) 565 { 566 PetscFunctionBegin; 567 *type = MATSOLVERPASTIX; 568 PetscFunctionReturn(0); 569 } 570 EXTERN_C_END 571 572 EXTERN_C_BEGIN 573 /* 574 The seq and mpi versions of this function are the same 575 */ 576 #undef __FUNCT__ 577 #define __FUNCT__ "MatGetFactor_seqaij_pastix" 578 PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F) 579 { 580 Mat B; 581 PetscErrorCode ierr; 582 Mat_Pastix *pastix; 583 584 PetscFunctionBegin; 585 if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 586 /* Create the factorization matrix */ 587 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 588 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 589 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 590 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 591 592 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 593 B->ops->view = MatView_PaStiX; 594 B->ops->getinfo = MatGetInfo_PaStiX; 595 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix", MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 596 B->factortype = MAT_FACTOR_LU; 597 598 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 599 pastix->CleanUpPastix = PETSC_FALSE; 600 pastix->isAIJ = PETSC_TRUE; 601 pastix->scat_rhs = PETSC_NULL; 602 pastix->scat_sol = PETSC_NULL; 603 pastix->MatDestroy = B->ops->destroy; 604 B->ops->destroy = MatDestroy_Pastix; 605 B->spptr = (void*)pastix; 606 607 *F = B; 608 PetscFunctionReturn(0); 609 } 610 EXTERN_C_END 611 612 613 EXTERN_C_BEGIN 614 #undef __FUNCT__ 615 #define __FUNCT__ "MatGetFactor_mpiaij_pastix" 616 PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F) 617 { 618 Mat B; 619 PetscErrorCode ierr; 620 Mat_Pastix *pastix; 621 622 PetscFunctionBegin; 623 if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 624 /* Create the factorization matrix */ 625 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 626 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 627 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 628 ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 629 ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 630 631 B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 632 B->ops->view = MatView_PaStiX; 633 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 634 B->factortype = MAT_FACTOR_LU; 635 636 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 637 pastix->CleanUpPastix = PETSC_FALSE; 638 pastix->isAIJ = PETSC_TRUE; 639 pastix->scat_rhs = PETSC_NULL; 640 pastix->scat_sol = PETSC_NULL; 641 pastix->MatDestroy = B->ops->destroy; 642 B->ops->destroy = MatDestroy_Pastix; 643 B->spptr = (void*)pastix; 644 645 *F = B; 646 PetscFunctionReturn(0); 647 } 648 EXTERN_C_END 649 650 EXTERN_C_BEGIN 651 #undef __FUNCT__ 652 #define __FUNCT__ "MatGetFactor_seqsbaij_pastix" 653 PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 654 { 655 Mat B; 656 PetscErrorCode ierr; 657 Mat_Pastix *pastix; 658 659 PetscFunctionBegin; 660 if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 661 /* Create the factorization matrix */ 662 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 663 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 664 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 665 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 666 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 667 668 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 669 B->ops->view = MatView_PaStiX; 670 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 671 B->factortype = MAT_FACTOR_CHOLESKY; 672 673 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 674 pastix->CleanUpPastix = PETSC_FALSE; 675 pastix->isAIJ = PETSC_TRUE; 676 pastix->scat_rhs = PETSC_NULL; 677 pastix->scat_sol = PETSC_NULL; 678 pastix->MatDestroy = B->ops->destroy; 679 B->ops->destroy = MatDestroy_Pastix; 680 B->spptr = (void*)pastix; 681 682 *F = B; 683 PetscFunctionReturn(0); 684 } 685 EXTERN_C_END 686 687 EXTERN_C_BEGIN 688 #undef __FUNCT__ 689 #define __FUNCT__ "MatGetFactor_mpisbaij_pastix" 690 PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 691 { 692 Mat B; 693 PetscErrorCode ierr; 694 Mat_Pastix *pastix; 695 696 PetscFunctionBegin; 697 if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 698 699 /* Create the factorization matrix */ 700 ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 701 ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 702 ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 703 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 704 ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 705 706 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 707 B->ops->view = MatView_PaStiX; 708 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 709 B->factortype = MAT_FACTOR_CHOLESKY; 710 711 ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 712 pastix->CleanUpPastix = PETSC_FALSE; 713 pastix->isAIJ = PETSC_TRUE; 714 pastix->scat_rhs = PETSC_NULL; 715 pastix->scat_sol = PETSC_NULL; 716 pastix->MatDestroy = B->ops->destroy; 717 B->ops->destroy = MatDestroy_Pastix; 718 B->spptr = (void*)pastix; 719 720 *F = B; 721 PetscFunctionReturn(0); 722 } 723 EXTERN_C_END 724