13bf14a46SMatthew Knepley #define PETSCMAT_DLL 23bf14a46SMatthew Knepley 33bf14a46SMatthew Knepley /* 43bf14a46SMatthew Knepley Provides an interface to the PaStiX sparse solver 53bf14a46SMatthew Knepley */ 63bf14a46SMatthew Knepley #include "../src/mat/impls/aij/seq/aij.h" 73bf14a46SMatthew Knepley #include "../src/mat/impls/aij/mpi/mpiaij.h" 83bf14a46SMatthew Knepley #include "../src/mat/impls/sbaij/seq/sbaij.h" 93bf14a46SMatthew Knepley #include "../src/mat/impls/sbaij/mpi/mpisbaij.h" 103bf14a46SMatthew Knepley 113bf14a46SMatthew Knepley EXTERN_C_BEGIN 123bf14a46SMatthew Knepley #include "mpi.h" 133bf14a46SMatthew Knepley #include "pastix.h" 143bf14a46SMatthew Knepley EXTERN_C_END 153bf14a46SMatthew Knepley 163bf14a46SMatthew Knepley typedef struct Mat_Pastix_ { 173bf14a46SMatthew Knepley pastix_data_t *pastix_data; /* Pastix data storage structure */ 183bf14a46SMatthew Knepley MatStructure matstruc; 193bf14a46SMatthew Knepley PetscInt n; /* Number of columns in the matrix */ 203bf14a46SMatthew Knepley PetscInt *colptr; /* Index of first element of each column in row and val */ 213bf14a46SMatthew Knepley PetscInt *row; /* Row of each element of the matrix */ 223bf14a46SMatthew Knepley PetscScalar *val; /* Value of each element of the matrix */ 233bf14a46SMatthew Knepley PetscInt *perm; /* Permutation tabular */ 243bf14a46SMatthew Knepley PetscInt *invp; /* Reverse permutation tabular */ 253bf14a46SMatthew Knepley PetscScalar *rhs; /* Rhight-hand-side member */ 263bf14a46SMatthew Knepley PetscInt rhsnbr; /* Rhight-hand-side number (must be 1) */ 273bf14a46SMatthew Knepley PetscInt iparm[64]; /* Integer parameters */ 283bf14a46SMatthew Knepley double dparm[64]; /* Floating point parameters */ 293bf14a46SMatthew Knepley MPI_Comm pastix_comm; /* PaStiX MPI communicator */ 303bf14a46SMatthew Knepley PetscMPIInt commRank; /* MPI rank */ 313bf14a46SMatthew Knepley PetscMPIInt commSize; /* MPI communicator size */ 323bf14a46SMatthew Knepley PetscTruth CleanUpPastix; /* Boolean indicating if we call PaStiX clean step */ 333bf14a46SMatthew Knepley VecScatter scat_rhs; 343bf14a46SMatthew Knepley VecScatter scat_sol; 353bf14a46SMatthew Knepley Vec b_seq,x_seq; 363bf14a46SMatthew Knepley PetscTruth isAIJ; 373bf14a46SMatthew Knepley PetscInt nSolve; /* Number of consecutive solve */ 383bf14a46SMatthew Knepley PetscErrorCode (*MatDestroy)(Mat); 393bf14a46SMatthew Knepley } Mat_Pastix; 403bf14a46SMatthew Knepley 413bf14a46SMatthew Knepley EXTERN PetscErrorCode MatDuplicate_Pastix(Mat,MatDuplicateOption,Mat*); 423bf14a46SMatthew Knepley 433bf14a46SMatthew Knepley /* 443bf14a46SMatthew Knepley convert Petsc seqaij matrix to CSC: colptr[n], row[nz], val[nz] 453bf14a46SMatthew Knepley 463bf14a46SMatthew Knepley input: 473bf14a46SMatthew Knepley A - matrix in seqaij or mpisbaij (bs=1) format 483bf14a46SMatthew Knepley valOnly - FALSE: spaces are allocated and values are set for the CSC 493bf14a46SMatthew Knepley TRUE: Only fill values 503bf14a46SMatthew Knepley output: 513bf14a46SMatthew Knepley n - Size of the matrix 523bf14a46SMatthew Knepley colptr - Index of first element of each column in row and val 533bf14a46SMatthew Knepley row - Row of each element of the matrix 543bf14a46SMatthew Knepley values - Value of each element of the matrix 553bf14a46SMatthew Knepley */ 56*41c8de11SBarry Smith PetscErrorCode MatConvertToCSC(Mat A,PetscTruth valOnly,PetscInt *n,PetscInt **colptr,PetscInt **row,PetscScalar **values) 57*41c8de11SBarry Smith { 583bf14a46SMatthew Knepley Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 593bf14a46SMatthew Knepley PetscInt *rowptr = aa->i; 603bf14a46SMatthew Knepley PetscInt *col = aa->j; 613bf14a46SMatthew Knepley PetscScalar *rvalues = aa->a; 623bf14a46SMatthew Knepley PetscInt m = A->rmap->N; 63745c78f7SBarry Smith PetscInt nnz; 643bf14a46SMatthew Knepley PetscInt i,j, k; 653bf14a46SMatthew Knepley PetscInt base = 1; 663bf14a46SMatthew Knepley PetscInt idx; 673bf14a46SMatthew Knepley PetscErrorCode ierr; 683bf14a46SMatthew Knepley PetscInt colidx; 693bf14a46SMatthew Knepley PetscInt *colcount; 70*41c8de11SBarry Smith PetscTruth isSBAIJ; 71*41c8de11SBarry Smith PetscTruth isSeqSBAIJ; 72*41c8de11SBarry Smith PetscTruth isMpiSBAIJ; 73745c78f7SBarry Smith PetscTruth isSym; 743bf14a46SMatthew Knepley 753bf14a46SMatthew Knepley 763bf14a46SMatthew Knepley PetscFunctionBegin; 773bf14a46SMatthew Knepley /* Allocate the CSC */ 783bf14a46SMatthew Knepley 79*41c8de11SBarry Smith ierr = MatIsSymmetric(A,0.0,&isSym);CHKERRQ(ierr); 80*41c8de11SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSBAIJ,&isSBAIJ);CHKERRQ(ierr); 81*41c8de11SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 82*41c8de11SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPISBAIJ,&isMpiSBAIJ);CHKERRQ(ierr); 833bf14a46SMatthew Knepley 84745c78f7SBarry Smith *n = A->cmap->N; 85745c78f7SBarry Smith 86745c78f7SBarry Smith /* PaStiX only needs triangular matrix if matrix is symmetric 87745c78f7SBarry Smith */ 88*41c8de11SBarry Smith if (isSym && !(isSBAIJ || isSeqSBAIJ || isMpiSBAIJ)) { 89745c78f7SBarry Smith nnz = (aa->nz - *n)/2 + *n; 90745c78f7SBarry Smith } 91*41c8de11SBarry Smith else { 92745c78f7SBarry Smith nnz = aa->nz; 933bf14a46SMatthew Knepley } 943bf14a46SMatthew Knepley 953bf14a46SMatthew Knepley if (!valOnly){ 963bf14a46SMatthew Knepley ierr = PetscMalloc(((*n)+1) *sizeof(PetscInt) ,colptr );CHKERRQ(ierr); 973bf14a46SMatthew Knepley ierr = PetscMalloc( nnz *sizeof(PetscInt) ,row);CHKERRQ(ierr); 983bf14a46SMatthew Knepley ierr = PetscMalloc( nnz *sizeof(PetscScalar),values);CHKERRQ(ierr); 993bf14a46SMatthew Knepley 100*41c8de11SBarry Smith if (isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) { 101*41c8de11SBarry Smith fprintf(stdout, "prout\n"); 102*41c8de11SBarry Smith ierr = PetscMemcpy (*colptr, rowptr, ((*n)+1)*sizeof(PetscInt));CHKERRQ(ierr); 103*41c8de11SBarry Smith for (i = 0; i < *n+1; i++) 104*41c8de11SBarry Smith (*colptr)[i] += base; 105*41c8de11SBarry Smith ierr = PetscMemcpy (*row, col, (nnz)*sizeof(PetscInt));CHKERRQ(ierr); 106*41c8de11SBarry Smith for (i = 0; i < nnz; i++) 107*41c8de11SBarry Smith (*row)[i] += base; 108*41c8de11SBarry Smith ierr = PetscMemcpy (*values, rvalues, (nnz)*sizeof(PetscScalar));CHKERRQ(ierr); 109*41c8de11SBarry Smith } else { 110*41c8de11SBarry Smith ierr = PetscMalloc((*n)*sizeof(PetscInt) ,&colcount);CHKERRQ(ierr); 111*41c8de11SBarry Smith 1123bf14a46SMatthew Knepley for (i = 0; i < m; i++) 1133bf14a46SMatthew Knepley colcount[i] = 0; 1143bf14a46SMatthew Knepley /* Fill-in colptr */ 1153bf14a46SMatthew Knepley for (i = 0; i < m; i++) 1163bf14a46SMatthew Knepley for (j = rowptr[i]; j < rowptr[i+1]; j++) 117745c78f7SBarry Smith if (!isSym || col[j] <= i) 1183bf14a46SMatthew Knepley colcount[col[j]]++; 119745c78f7SBarry Smith 1203bf14a46SMatthew Knepley (*colptr)[0] = base; 1213bf14a46SMatthew Knepley for (j = 0; j < *n; j++) { 1223bf14a46SMatthew Knepley (*colptr)[j+1] = (*colptr)[j] + colcount[j]; 123745c78f7SBarry Smith /* in next loop we fill starting from (*colptr)[colidx] - base */ 1243bf14a46SMatthew Knepley colcount[j] = -base; 1253bf14a46SMatthew Knepley } 1263bf14a46SMatthew Knepley 1273bf14a46SMatthew Knepley /* Fill-in rows and values */ 1283bf14a46SMatthew Knepley for (i = 0; i < m; i++) { 1293bf14a46SMatthew Knepley for (j = rowptr[i]; j < rowptr[i+1]; j++) { 130*41c8de11SBarry Smith if (!isSym || col[j] <= i) { 1313bf14a46SMatthew Knepley colidx = col[j]; 1323bf14a46SMatthew Knepley idx = (*colptr)[colidx] + colcount[colidx]; 1333bf14a46SMatthew Knepley (*row)[idx] = i + base; 1343bf14a46SMatthew Knepley (*values)[idx] = rvalues[j]; 1353bf14a46SMatthew Knepley colcount[colidx]++; 1363bf14a46SMatthew Knepley } 1373bf14a46SMatthew Knepley } 1383bf14a46SMatthew Knepley } 139*41c8de11SBarry Smith ierr = PetscFree(colcount);CHKERRQ(ierr); 140745c78f7SBarry Smith } 141*41c8de11SBarry Smith } else { 142745c78f7SBarry Smith /* Fill-in only values */ 1433bf14a46SMatthew Knepley for (i = 0; i < m; i++) { 1443bf14a46SMatthew Knepley for (j = rowptr[i]; j < rowptr[i+1]; j++) { 1453bf14a46SMatthew Knepley colidx = col[j]; 146*41c8de11SBarry Smith if ((isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) ||!isSym || col[j] <= i) 147745c78f7SBarry Smith { 148745c78f7SBarry Smith /* look for the value to fill */ 149745c78f7SBarry Smith for (k = (*colptr)[colidx] - base; 150745c78f7SBarry Smith k < (*colptr)[colidx + 1] - base; 151745c78f7SBarry Smith k++) { 1523bf14a46SMatthew Knepley if ((*row)[k] == i) { 1533bf14a46SMatthew Knepley (*values)[k] = rvalues[j]; 1543bf14a46SMatthew Knepley break; 1553bf14a46SMatthew Knepley } 1563bf14a46SMatthew Knepley } 157745c78f7SBarry Smith /* shouldn't happen, overflow */ 1583bf14a46SMatthew Knepley if (k == (*colptr)[colidx + 1] - base) 1593bf14a46SMatthew Knepley PetscFunctionReturn(1); 1603bf14a46SMatthew Knepley } 1613bf14a46SMatthew Knepley } 1623bf14a46SMatthew Knepley } 163745c78f7SBarry Smith } 1643bf14a46SMatthew Knepley PetscFunctionReturn(0); 1653bf14a46SMatthew Knepley } 1663bf14a46SMatthew Knepley 1673bf14a46SMatthew Knepley 1683bf14a46SMatthew Knepley 1693bf14a46SMatthew Knepley #undef __FUNCT__ 1703bf14a46SMatthew Knepley #define __FUNCT__ "MatDestroy_Pastix" 1713bf14a46SMatthew Knepley /* 1723bf14a46SMatthew Knepley Call clean step of PaStiX if lu->CleanUpPastix == true. 1733bf14a46SMatthew Knepley Free the CSC matrix. 1743bf14a46SMatthew Knepley */ 1753bf14a46SMatthew Knepley PetscErrorCode MatDestroy_Pastix(Mat A) 1763bf14a46SMatthew Knepley { 1773bf14a46SMatthew Knepley Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 1783bf14a46SMatthew Knepley PetscErrorCode ierr; 1793bf14a46SMatthew Knepley PetscMPIInt size=lu->commSize; 180745c78f7SBarry Smith 1813bf14a46SMatthew Knepley PetscFunctionBegin; 1823bf14a46SMatthew Knepley if (lu->CleanUpPastix) { 1833bf14a46SMatthew Knepley /* Terminate instance, deallocate memories */ 1843bf14a46SMatthew Knepley if (size > 1){ 1853bf14a46SMatthew Knepley ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr); 1863bf14a46SMatthew Knepley ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr); 1873bf14a46SMatthew Knepley if (lu->nSolve && lu->scat_sol){ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);} 1883bf14a46SMatthew Knepley if (lu->nSolve && lu->x_seq){ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr);} 1893bf14a46SMatthew Knepley } 1903bf14a46SMatthew Knepley 1913bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK]=API_TASK_CLEAN; 1923bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK]=API_TASK_CLEAN; 1933bf14a46SMatthew Knepley 1943bf14a46SMatthew Knepley pastix((pastix_data_t **)&(lu->pastix_data), 1953bf14a46SMatthew Knepley lu->pastix_comm, 1963bf14a46SMatthew Knepley (pastix_int_t) lu->n, 1973bf14a46SMatthew Knepley (pastix_int_t*) lu->colptr, 1983bf14a46SMatthew Knepley (pastix_int_t*) lu->row, 1993bf14a46SMatthew Knepley (pastix_float_t*) lu->val, 2003bf14a46SMatthew Knepley (pastix_int_t*) lu->perm, 2013bf14a46SMatthew Knepley (pastix_int_t*) lu->invp, 2023bf14a46SMatthew Knepley (pastix_float_t*) lu->rhs, 2033bf14a46SMatthew Knepley (pastix_int_t) lu->rhsnbr, 2043bf14a46SMatthew Knepley (pastix_int_t*) lu->iparm, 2053bf14a46SMatthew Knepley lu->dparm); 2063bf14a46SMatthew Knepley 2073bf14a46SMatthew Knepley ierr = PetscFree(lu->colptr);CHKERRQ(ierr); 2083bf14a46SMatthew Knepley ierr = PetscFree(lu->row); CHKERRQ(ierr); 2093bf14a46SMatthew Knepley ierr = PetscFree(lu->val); CHKERRQ(ierr); 2103bf14a46SMatthew Knepley ierr = PetscFree(lu->perm); CHKERRQ(ierr); 2113bf14a46SMatthew Knepley ierr = PetscFree(lu->invp); CHKERRQ(ierr); 2123bf14a46SMatthew Knepley /* ierr = PetscFree(lu->rhs); CHKERRQ(ierr); */ 2133bf14a46SMatthew Knepley ierr = MPI_Comm_free(&(lu->pastix_comm));CHKERRQ(ierr); 2143bf14a46SMatthew Knepley } 2153bf14a46SMatthew Knepley ierr = (lu->MatDestroy)(A);CHKERRQ(ierr); 2163bf14a46SMatthew Knepley PetscFunctionReturn(0); 2173bf14a46SMatthew Knepley } 2183bf14a46SMatthew Knepley 2193bf14a46SMatthew Knepley #undef __FUNCT__ 2203bf14a46SMatthew Knepley #define __FUNCT__ "MatSolve_PaStiX" 2213bf14a46SMatthew Knepley /* 2223bf14a46SMatthew Knepley Gather right-hand-side. 2233bf14a46SMatthew Knepley Call for Solve step. 2243bf14a46SMatthew Knepley Scatter solution. 2253bf14a46SMatthew Knepley */ 2263bf14a46SMatthew Knepley PetscErrorCode MatSolve_PaStiX(Mat A,Vec b,Vec x) 2273bf14a46SMatthew Knepley { 2283bf14a46SMatthew Knepley Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 2293bf14a46SMatthew Knepley PetscScalar *array; 2303bf14a46SMatthew Knepley Vec x_seq; 2313bf14a46SMatthew Knepley PetscErrorCode ierr; 2323bf14a46SMatthew Knepley 2333bf14a46SMatthew Knepley PetscFunctionBegin; 2343bf14a46SMatthew Knepley lu->rhsnbr = 1; 2353bf14a46SMatthew Knepley x_seq = lu->b_seq; 2363bf14a46SMatthew Knepley if (lu->commSize > 1){ 2373bf14a46SMatthew Knepley /* PaStiX only supports centralized rhs. Scatter b into a seqential rhs vector */ 2383bf14a46SMatthew Knepley ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2393bf14a46SMatthew Knepley ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 240b5e56a35SBarry Smith ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr); 241*41c8de11SBarry Smith } else { /* size == 1 */ 2423bf14a46SMatthew Knepley ierr = VecCopy(b,x);CHKERRQ(ierr); 2433bf14a46SMatthew Knepley ierr = VecGetArray(x,&array);CHKERRQ(ierr); 2443bf14a46SMatthew Knepley } 2453bf14a46SMatthew Knepley lu->rhs = array; 2463bf14a46SMatthew Knepley if (lu->commSize == 1){ 2473bf14a46SMatthew Knepley ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 2483bf14a46SMatthew Knepley } else { 2493bf14a46SMatthew Knepley ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 2503bf14a46SMatthew Knepley } 2513bf14a46SMatthew Knepley 2523bf14a46SMatthew Knepley /* solve phase */ 2533bf14a46SMatthew Knepley /*-------------*/ 2543bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE; 2553bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_REFINE; 256745c78f7SBarry Smith lu->iparm[IPARM_RHS_MAKING] = API_RHS_B; 2573bf14a46SMatthew Knepley 2583bf14a46SMatthew Knepley pastix((pastix_data_t **)&(lu->pastix_data), 2593bf14a46SMatthew Knepley (MPI_Comm) lu->pastix_comm, 2603bf14a46SMatthew Knepley (pastix_int_t) lu->n, 2613bf14a46SMatthew Knepley (pastix_int_t*) lu->colptr, 2623bf14a46SMatthew Knepley (pastix_int_t*) lu->row, 2633bf14a46SMatthew Knepley (pastix_float_t*) lu->val, 2643bf14a46SMatthew Knepley (pastix_int_t*) lu->perm, 2653bf14a46SMatthew Knepley (pastix_int_t*) lu->invp, 2663bf14a46SMatthew Knepley (pastix_float_t*) lu->rhs, 2673bf14a46SMatthew Knepley (pastix_int_t) lu->rhsnbr, 2683bf14a46SMatthew Knepley (pastix_int_t*) lu->iparm, 2693bf14a46SMatthew Knepley (double*) lu->dparm); 2703bf14a46SMatthew Knepley 2713bf14a46SMatthew Knepley if (lu->iparm[IPARM_ERROR_NUMBER] < 0) { 272b5e56a35SBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in solve phase: lu->iparm[IPARM_ERROR_NUMBER] = %d\n",lu->iparm[IPARM_ERROR_NUMBER] ); 2733bf14a46SMatthew Knepley } 2743bf14a46SMatthew Knepley 2753bf14a46SMatthew Knepley if (lu->commSize == 1){ 2763bf14a46SMatthew Knepley ierr = VecRestoreArray(x,&(lu->rhs));CHKERRQ(ierr); 2773bf14a46SMatthew Knepley } else { 2783bf14a46SMatthew Knepley ierr = VecRestoreArray(x_seq,&(lu->rhs));CHKERRQ(ierr); 2793bf14a46SMatthew Knepley } 2803bf14a46SMatthew Knepley 2813bf14a46SMatthew Knepley if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */ 2823bf14a46SMatthew Knepley ierr = VecScatterBegin(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2833bf14a46SMatthew Knepley ierr = VecScatterEnd(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2843bf14a46SMatthew Knepley } 2853bf14a46SMatthew Knepley lu->nSolve++; 2863bf14a46SMatthew Knepley PetscFunctionReturn(0); 2873bf14a46SMatthew Knepley } 2883bf14a46SMatthew Knepley 2893bf14a46SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) 2903bf14a46SMatthew Knepley /* 2913bf14a46SMatthew Knepley TODO: Fill this function 2923bf14a46SMatthew Knepley I didn't fill this function 2933bf14a46SMatthew Knepley because I didn't understood its goal. 2943bf14a46SMatthew Knepley */ 2953bf14a46SMatthew Knepley 2963bf14a46SMatthew Knepley /* 2973bf14a46SMatthew Knepley input: 2983bf14a46SMatthew Knepley F: numeric factor 2993bf14a46SMatthew Knepley output: 3003bf14a46SMatthew Knepley nneg: total number of pivots 3013bf14a46SMatthew Knepley nzero: 0 3023bf14a46SMatthew Knepley npos: (global dimension of F) - nneg 3033bf14a46SMatthew Knepley */ 3043bf14a46SMatthew Knepley 3053bf14a46SMatthew Knepley #undef __FUNCT__ 3063bf14a46SMatthew Knepley #define __FUNCT__ "MatGetInertia_SBAIJPASTIX" 3073bf14a46SMatthew Knepley PetscErrorCode MatGetInertia_SBAIJPASTIX(Mat F,int *nneg,int *nzero,int *npos) 3083bf14a46SMatthew Knepley { 3093bf14a46SMatthew Knepley PetscFunctionBegin; 3103bf14a46SMatthew Knepley /* ierr = MPI_Comm_size(((PetscObject)F)->comm,&size);CHKERRQ(ierr); */ 3113bf14a46SMatthew Knepley /* /\* PASTIX 4.3.1 calls ScaLAPACK when ICNTL(13)=0 (default), which does not offer the possibility to compute the inertia of a dense matrix. Set ICNTL(13)=1 to skip ScaLAPACK *\/ */ 3123bf14a46SMatthew Knepley /* if (size > 1 && lu->id.ICNTL(13) != 1){ */ 3133bf14a46SMatthew Knepley /* SETERRQ1(PETSC_ERR_ARG_WRONG,"ICNTL(13)=%d. -mat_pastix_icntl_13 must be set as 1 for correct global matrix inertia\n",lu->id.INFOG(13)); */ 3143bf14a46SMatthew Knepley /* } */ 3153bf14a46SMatthew Knepley /* if (nneg){ */ 3163bf14a46SMatthew Knepley /* if (!lu->commSize){ */ 3173bf14a46SMatthew Knepley /* *nneg = lu->id.INFOG(12); */ 3183bf14a46SMatthew Knepley /* } */ 3193bf14a46SMatthew Knepley /* ierr = MPI_Bcast(nneg,1,MPI_INT,0,lu->comm_pastix);CHKERRQ(ierr); */ 3203bf14a46SMatthew Knepley /* } */ 3213bf14a46SMatthew Knepley /* if (nzero) *nzero = lu->iparm[IPARM_NNZEROS]; */ 3223bf14a46SMatthew Knepley /* if (npos) *npos = F->rmap->N - (*nneg); */ 3233bf14a46SMatthew Knepley PetscFunctionReturn(0); 3243bf14a46SMatthew Knepley } 3253bf14a46SMatthew Knepley #endif /* !defined(PETSC_USE_COMPLEX) */ 3263bf14a46SMatthew Knepley 3273bf14a46SMatthew Knepley /* 3283bf14a46SMatthew Knepley Numeric factorisation using PaStiX solver. 3293bf14a46SMatthew Knepley 3303bf14a46SMatthew Knepley */ 3313bf14a46SMatthew Knepley #undef __FUNCT__ 3323bf14a46SMatthew Knepley #define __FUNCT__ "MatFactorNumeric_PASTIX" 3333bf14a46SMatthew Knepley PetscErrorCode MatFactorNumeric_PaStiX(Mat F,Mat A,const MatFactorInfo *info) 3343bf14a46SMatthew Knepley { 3353bf14a46SMatthew Knepley Mat_Pastix *lu =(Mat_Pastix*)(F)->spptr; 336*41c8de11SBarry Smith Mat *tseq; 3373bf14a46SMatthew Knepley PetscErrorCode ierr = 0; 338b5e56a35SBarry Smith PetscInt icntl; 339b5e56a35SBarry Smith PetscInt M=A->rmap->N; 3403bf14a46SMatthew Knepley PetscTruth valOnly,flg, isSym; 3413bf14a46SMatthew Knepley Mat F_diag; 3423bf14a46SMatthew Knepley IS is_iden; 3433bf14a46SMatthew Knepley Vec b; 3443bf14a46SMatthew Knepley IS isrow; 3453bf14a46SMatthew Knepley PetscTruth isSeqAIJ,isSeqSBAIJ; 3463bf14a46SMatthew Knepley 3473bf14a46SMatthew Knepley PetscFunctionBegin; 348*41c8de11SBarry Smith 3493bf14a46SMatthew Knepley ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 3503bf14a46SMatthew Knepley ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 3513bf14a46SMatthew Knepley if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 3523bf14a46SMatthew Knepley (F)->ops->solve = MatSolve_PaStiX; 3533bf14a46SMatthew Knepley 3543bf14a46SMatthew Knepley /* Initialize a PASTIX instance */ 3553bf14a46SMatthew Knepley ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr); 3563bf14a46SMatthew Knepley ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank); CHKERRQ(ierr); 3573bf14a46SMatthew Knepley ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize); CHKERRQ(ierr); 3583bf14a46SMatthew Knepley 3593bf14a46SMatthew Knepley /* Set pastix options */ 3603bf14a46SMatthew Knepley lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO; 3613bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_INIT; 3623bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_INIT; 3633bf14a46SMatthew Knepley lu->rhsnbr = 1; 3643bf14a46SMatthew Knepley 3653bf14a46SMatthew Knepley /* Call to set default pastix options */ 3663bf14a46SMatthew Knepley pastix((pastix_data_t **)&(lu->pastix_data), 3673bf14a46SMatthew Knepley (MPI_Comm) lu->pastix_comm, 3683bf14a46SMatthew Knepley (pastix_int_t) lu->n, 3693bf14a46SMatthew Knepley (pastix_int_t*) lu->colptr, 3703bf14a46SMatthew Knepley (pastix_int_t*) lu->row, 3713bf14a46SMatthew Knepley (pastix_float_t*) lu->val, 3723bf14a46SMatthew Knepley (pastix_int_t*) lu->perm, 3733bf14a46SMatthew Knepley (pastix_int_t*) lu->invp, 3743bf14a46SMatthew Knepley (pastix_float_t*) lu->rhs, 3753bf14a46SMatthew Knepley (pastix_int_t) lu->rhsnbr, 3763bf14a46SMatthew Knepley (pastix_int_t*) lu->iparm, 3773bf14a46SMatthew Knepley (double*) lu->dparm); 3783bf14a46SMatthew Knepley 3793bf14a46SMatthew Knepley ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr); 3803bf14a46SMatthew Knepley 3813bf14a46SMatthew Knepley icntl=-1; 382*41c8de11SBarry Smith lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT; 383*41c8de11SBarry Smith ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr); 3843bf14a46SMatthew Knepley if ((flg && icntl > 0) || PetscLogPrintInfo) { 3853bf14a46SMatthew Knepley lu->iparm[IPARM_VERBOSE] = icntl; 3863bf14a46SMatthew Knepley } 3873bf14a46SMatthew Knepley icntl=-1; 388*41c8de11SBarry Smith 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); 3893bf14a46SMatthew Knepley if ((flg && icntl > 0)) { 3903bf14a46SMatthew Knepley lu->iparm[IPARM_THREAD_NBR] = icntl; 3913bf14a46SMatthew Knepley } 3923bf14a46SMatthew Knepley PetscOptionsEnd(); 3933bf14a46SMatthew Knepley valOnly = PETSC_FALSE; 394*41c8de11SBarry Smith } else { 3953bf14a46SMatthew Knepley valOnly = PETSC_TRUE; 3963bf14a46SMatthew Knepley } 3973bf14a46SMatthew Knepley 3983bf14a46SMatthew Knepley lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES; 3993bf14a46SMatthew Knepley 4003bf14a46SMatthew Knepley /* convert mpi A to seq mat A */ 4013bf14a46SMatthew Knepley ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 4023bf14a46SMatthew Knepley ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 4033bf14a46SMatthew Knepley ierr = ISDestroy(isrow);CHKERRQ(ierr); 4043bf14a46SMatthew Knepley 405*41c8de11SBarry Smith ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val); CHKERRQ(ierr); 406*41c8de11SBarry Smith ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr); 407*41c8de11SBarry Smith ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr); 408*41c8de11SBarry Smith 4093bf14a46SMatthew Knepley ierr = PetscMalloc((lu->n)*sizeof(PetscInt) ,&(lu->perm));CHKERRQ(ierr); 4103bf14a46SMatthew Knepley ierr = PetscMalloc((lu->n)*sizeof(PetscInt) ,&(lu->invp));CHKERRQ(ierr); 4113bf14a46SMatthew Knepley 4123bf14a46SMatthew Knepley 4133bf14a46SMatthew Knepley if (isSym) { 414745c78f7SBarry Smith /* On symmetric matrix, LLT */ 4153bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_YES; 416*41c8de11SBarry Smith lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT; 4173bf14a46SMatthew Knepley } 4183bf14a46SMatthew Knepley else { 419745c78f7SBarry Smith /* On unsymmetric matrix, LU */ 4203bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_NO; 4213bf14a46SMatthew Knepley lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 4223bf14a46SMatthew Knepley } 4233bf14a46SMatthew Knepley 4243bf14a46SMatthew Knepley /*----------------*/ 4253bf14a46SMatthew Knepley if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){ 4263bf14a46SMatthew Knepley if (!(isSeqAIJ || isSeqSBAIJ)) { 4273bf14a46SMatthew Knepley /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 4283bf14a46SMatthew Knepley ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 4293bf14a46SMatthew Knepley ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 4303bf14a46SMatthew Knepley ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 4313bf14a46SMatthew Knepley ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 4323bf14a46SMatthew Knepley ierr = VecSetFromOptions(b);CHKERRQ(ierr); 4333bf14a46SMatthew Knepley 4343bf14a46SMatthew Knepley ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 4353bf14a46SMatthew Knepley ierr = VecScatterCreate(lu->b_seq,is_iden,b,is_iden,&lu->scat_sol);CHKERRQ(ierr); 4363bf14a46SMatthew Knepley ierr = ISDestroy(is_iden);CHKERRQ(ierr); 4373bf14a46SMatthew Knepley ierr = VecDestroy(b);CHKERRQ(ierr); 4383bf14a46SMatthew Knepley } 4393bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING; 4403bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 4413bf14a46SMatthew Knepley 4423bf14a46SMatthew Knepley pastix((pastix_data_t **)&(lu->pastix_data), 4433bf14a46SMatthew Knepley (MPI_Comm) lu->pastix_comm, 4443bf14a46SMatthew Knepley (pastix_int_t) lu->n, 4453bf14a46SMatthew Knepley (pastix_int_t*) lu->colptr, 4463bf14a46SMatthew Knepley (pastix_int_t*) lu->row, 4473bf14a46SMatthew Knepley (pastix_float_t*) lu->val, 4483bf14a46SMatthew Knepley (pastix_int_t*) lu->perm, 4493bf14a46SMatthew Knepley (pastix_int_t*) lu->invp, 4503bf14a46SMatthew Knepley (pastix_float_t*) lu->rhs, 4513bf14a46SMatthew Knepley (pastix_int_t) lu->rhsnbr, 4523bf14a46SMatthew Knepley (pastix_int_t*) lu->iparm, 4533bf14a46SMatthew Knepley (double*) lu->dparm); 4543bf14a46SMatthew Knepley if (lu->iparm[IPARM_ERROR_NUMBER] < 0) { 455*41c8de11SBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%d\n",lu->iparm[IPARM_ERROR_NUMBER]); 4563bf14a46SMatthew Knepley } 457*41c8de11SBarry Smith } else { 4583bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT; 4593bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 4603bf14a46SMatthew Knepley pastix((pastix_data_t **)&(lu->pastix_data), 4613bf14a46SMatthew Knepley (MPI_Comm) lu->pastix_comm, 4623bf14a46SMatthew Knepley (pastix_int_t) lu->n, 4633bf14a46SMatthew Knepley (pastix_int_t*) lu->colptr, 4643bf14a46SMatthew Knepley (pastix_int_t*) lu->row, 4653bf14a46SMatthew Knepley (pastix_float_t*) lu->val, 4663bf14a46SMatthew Knepley (pastix_int_t*) lu->perm, 4673bf14a46SMatthew Knepley (pastix_int_t*) lu->invp, 4683bf14a46SMatthew Knepley (pastix_float_t*) lu->rhs, 4693bf14a46SMatthew Knepley (pastix_int_t) lu->rhsnbr, 4703bf14a46SMatthew Knepley (pastix_int_t*) lu->iparm, 4713bf14a46SMatthew Knepley (double*) lu->dparm); 4723bf14a46SMatthew Knepley 4733bf14a46SMatthew Knepley if (lu->iparm[IPARM_ERROR_NUMBER] < 0) { 474*41c8de11SBarry Smith SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%d\n",lu->iparm[IPARM_ERROR_NUMBER]); 4753bf14a46SMatthew Knepley } 4763bf14a46SMatthew Knepley } 4773bf14a46SMatthew Knepley 4783bf14a46SMatthew Knepley if (lu->commSize > 1){ 4793bf14a46SMatthew Knepley if ((F)->factor == MAT_FACTOR_LU){ 4803bf14a46SMatthew Knepley F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 4813bf14a46SMatthew Knepley } else { 4823bf14a46SMatthew Knepley F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 4833bf14a46SMatthew Knepley } 4843bf14a46SMatthew Knepley F_diag->assembled = PETSC_TRUE; 4853bf14a46SMatthew Knepley if (lu->nSolve){ 4863bf14a46SMatthew Knepley ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr); 4873bf14a46SMatthew Knepley ierr = VecDestroy(lu->x_seq);CHKERRQ(ierr); 4883bf14a46SMatthew Knepley } 4893bf14a46SMatthew Knepley } 4903bf14a46SMatthew Knepley (F)->assembled = PETSC_TRUE; 4913bf14a46SMatthew Knepley lu->matstruc = SAME_NONZERO_PATTERN; 4923bf14a46SMatthew Knepley lu->CleanUpPastix = PETSC_TRUE; 4933bf14a46SMatthew Knepley lu->nSolve = 0; 4943bf14a46SMatthew Knepley PetscFunctionReturn(0); 4953bf14a46SMatthew Knepley } 4963bf14a46SMatthew Knepley 4973bf14a46SMatthew Knepley 4983bf14a46SMatthew Knepley /* Note the Petsc r and c permutations are ignored */ 4993bf14a46SMatthew Knepley #undef __FUNCT__ 5003bf14a46SMatthew Knepley #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX" 5013bf14a46SMatthew Knepley PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 5023bf14a46SMatthew Knepley { 5033bf14a46SMatthew Knepley Mat_Pastix *lu = (Mat_Pastix*)F->spptr; 5043bf14a46SMatthew Knepley 5053bf14a46SMatthew Knepley PetscFunctionBegin; 5063bf14a46SMatthew Knepley lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 5073bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_YES; 5083bf14a46SMatthew Knepley lu->matstruc = DIFFERENT_NONZERO_PATTERN; 5093bf14a46SMatthew Knepley F->ops->lufactornumeric = MatFactorNumeric_PaStiX; 5103bf14a46SMatthew Knepley PetscFunctionReturn(0); 5113bf14a46SMatthew Knepley } 5123bf14a46SMatthew Knepley 5133bf14a46SMatthew Knepley 5143bf14a46SMatthew Knepley /* Note the Petsc r permutation is ignored */ 5153bf14a46SMatthew Knepley #undef __FUNCT__ 5163bf14a46SMatthew Knepley #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX" 5173bf14a46SMatthew Knepley PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info) 5183bf14a46SMatthew Knepley { 5193bf14a46SMatthew Knepley Mat_Pastix *lu = (Mat_Pastix*)(F)->spptr; 5203bf14a46SMatthew Knepley 5213bf14a46SMatthew Knepley PetscFunctionBegin; 5223bf14a46SMatthew Knepley lu->iparm[IPARM_FACTORIZATION] = API_FACT_LLT; 5233bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_NO; 5243bf14a46SMatthew Knepley lu->matstruc = DIFFERENT_NONZERO_PATTERN; 5253bf14a46SMatthew Knepley (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX; 5263bf14a46SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) 5273bf14a46SMatthew Knepley (F)->ops->getinertia = MatGetInertia_SBAIJPASTIX; 5283bf14a46SMatthew Knepley #endif 5293bf14a46SMatthew Knepley PetscFunctionReturn(0); 5303bf14a46SMatthew Knepley } 5313bf14a46SMatthew Knepley 5323bf14a46SMatthew Knepley #undef __FUNCT__ 5333bf14a46SMatthew Knepley #define __FUNCT__ "MatView_PaStiX" 5343bf14a46SMatthew Knepley PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer) 5353bf14a46SMatthew Knepley { 5363bf14a46SMatthew Knepley PetscErrorCode ierr; 5373bf14a46SMatthew Knepley PetscTruth iascii; 5383bf14a46SMatthew Knepley PetscViewerFormat format; 5393bf14a46SMatthew Knepley 5403bf14a46SMatthew Knepley PetscFunctionBegin; 5413bf14a46SMatthew Knepley ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 5423bf14a46SMatthew Knepley if (iascii) { 5433bf14a46SMatthew Knepley ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 5443bf14a46SMatthew Knepley if (format == PETSC_VIEWER_ASCII_INFO){ 545b5e56a35SBarry Smith Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 546b5e56a35SBarry Smith 547b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr); 548b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Matrix type : %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES)?"Symmetric":"Unsymmetric"));CHKERRQ(ierr); 549b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Level of printing (0,1,2): %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr); 550b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr); 551b5e56a35SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Error : %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr); 5523bf14a46SMatthew Knepley } 5533bf14a46SMatthew Knepley } 5543bf14a46SMatthew Knepley PetscFunctionReturn(0); 5553bf14a46SMatthew Knepley } 5563bf14a46SMatthew Knepley 5573bf14a46SMatthew Knepley 5583bf14a46SMatthew Knepley /*MC 559b5e56a35SBarry Smith MAT_SOLVER_PASTIX - A solver package providing direct solvers (LU) for distributed 5603bf14a46SMatthew Knepley and sequential matrices via the external package PaStiX. 5613bf14a46SMatthew Knepley 562b5e56a35SBarry Smith Use config/configure.py --download-pastix to have PETSc installed with PaStiX 5633bf14a46SMatthew Knepley 5643bf14a46SMatthew Knepley Options Database Keys: 565b5e56a35SBarry Smith + -mat_pastix_verbose <0,1,2> - print level 566b5e56a35SBarry Smith - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task. 5673bf14a46SMatthew Knepley 5683bf14a46SMatthew Knepley Level: beginner 5693bf14a46SMatthew Knepley 570*41c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 571*41c8de11SBarry Smith 5723bf14a46SMatthew Knepley M*/ 5733bf14a46SMatthew Knepley 5743bf14a46SMatthew Knepley 5753bf14a46SMatthew Knepley #undef __FUNCT__ 5763bf14a46SMatthew Knepley #define __FUNCT__ "MatGetInfo_PaStiX" 5773bf14a46SMatthew Knepley PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info) 5783bf14a46SMatthew Knepley { 5793bf14a46SMatthew Knepley Mat_Pastix *lu =(Mat_Pastix*)A->spptr; 5803bf14a46SMatthew Knepley 5813bf14a46SMatthew Knepley PetscFunctionBegin; 5823bf14a46SMatthew Knepley info->block_size = 1.0; 5833bf14a46SMatthew Knepley info->nz_allocated = lu->iparm[IPARM_NNZEROS]; 5843bf14a46SMatthew Knepley info->nz_used = lu->iparm[IPARM_NNZEROS]; 5853bf14a46SMatthew Knepley info->nz_unneeded = 0.0; 5863bf14a46SMatthew Knepley info->assemblies = 0.0; 5873bf14a46SMatthew Knepley info->mallocs = 0.0; 5883bf14a46SMatthew Knepley info->memory = 0.0; 5893bf14a46SMatthew Knepley info->fill_ratio_given = 0; 5903bf14a46SMatthew Knepley info->fill_ratio_needed = 0; 5913bf14a46SMatthew Knepley info->factor_mallocs = 0; 5923bf14a46SMatthew Knepley PetscFunctionReturn(0); 5933bf14a46SMatthew Knepley } 5943bf14a46SMatthew Knepley 5953bf14a46SMatthew Knepley EXTERN_C_BEGIN 5963bf14a46SMatthew Knepley #undef __FUNCT__ 5973bf14a46SMatthew Knepley #define __FUNCT__ "MatFactorGetSolverPackage_pastix" 5983bf14a46SMatthew Knepley PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type) 5993bf14a46SMatthew Knepley { 6003bf14a46SMatthew Knepley PetscFunctionBegin; 6013bf14a46SMatthew Knepley *type = MAT_SOLVER_PASTIX; 6023bf14a46SMatthew Knepley PetscFunctionReturn(0); 6033bf14a46SMatthew Knepley } 6043bf14a46SMatthew Knepley EXTERN_C_END 6053bf14a46SMatthew Knepley 6063bf14a46SMatthew Knepley EXTERN_C_BEGIN 6073bf14a46SMatthew Knepley /* 6083bf14a46SMatthew Knepley The seq and mpi versions of this function are the same 6093bf14a46SMatthew Knepley */ 6103bf14a46SMatthew Knepley #undef __FUNCT__ 6113bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_seqaij_pastix" 6123bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F) 6133bf14a46SMatthew Knepley { 6143bf14a46SMatthew Knepley Mat B; 6153bf14a46SMatthew Knepley PetscErrorCode ierr; 6163bf14a46SMatthew Knepley Mat_Pastix *pastix; 6173bf14a46SMatthew Knepley 6183bf14a46SMatthew Knepley PetscFunctionBegin; 6193bf14a46SMatthew Knepley if (ftype != MAT_FACTOR_LU) { 6203bf14a46SMatthew Knepley SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 6213bf14a46SMatthew Knepley } 6223bf14a46SMatthew Knepley /* Create the factorization matrix */ 6233bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 6243bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 6253bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 6263bf14a46SMatthew Knepley ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 6273bf14a46SMatthew Knepley 6283bf14a46SMatthew Knepley B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 6293bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 6303bf14a46SMatthew Knepley B->ops->getinfo = MatGetInfo_PaStiX; 6313bf14a46SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C", 6323bf14a46SMatthew Knepley "MatFactorGetSolverPackage_pastix", 6333bf14a46SMatthew Knepley MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 6343bf14a46SMatthew Knepley B->factor = MAT_FACTOR_LU; 6353bf14a46SMatthew Knepley 6363bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 6373bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 6383bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 6393bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 6403bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 6413bf14a46SMatthew Knepley pastix->nSolve = 0; 6423bf14a46SMatthew Knepley pastix->MatDestroy = B->ops->destroy; 6433bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 6443bf14a46SMatthew Knepley B->spptr = (void*)pastix; 6453bf14a46SMatthew Knepley 6463bf14a46SMatthew Knepley *F = B; 6473bf14a46SMatthew Knepley PetscFunctionReturn(0); 6483bf14a46SMatthew Knepley } 6493bf14a46SMatthew Knepley EXTERN_C_END 6503bf14a46SMatthew Knepley 651b5e56a35SBarry Smith 6523bf14a46SMatthew Knepley EXTERN_C_BEGIN 6533bf14a46SMatthew Knepley #undef __FUNCT__ 6543bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_mpiaij_pastix" 6553bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F) 6563bf14a46SMatthew Knepley { 6573bf14a46SMatthew Knepley Mat B; 6583bf14a46SMatthew Knepley PetscErrorCode ierr; 6593bf14a46SMatthew Knepley Mat_Pastix *pastix; 6603bf14a46SMatthew Knepley 6613bf14a46SMatthew Knepley PetscFunctionBegin; 662*41c8de11SBarry Smith if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 6633bf14a46SMatthew Knepley /* Create the factorization matrix */ 6643bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 6653bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 6663bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 6673bf14a46SMatthew Knepley ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 6683bf14a46SMatthew Knepley ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 6693bf14a46SMatthew Knepley 6703bf14a46SMatthew Knepley B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 6713bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 6723bf14a46SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)B, 6733bf14a46SMatthew Knepley "MatFactorGetSolverPackage_C", 6743bf14a46SMatthew Knepley "MatFactorGetSolverPackage_pastix", 6753bf14a46SMatthew Knepley MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 6763bf14a46SMatthew Knepley B->factor = MAT_FACTOR_LU; 6773bf14a46SMatthew Knepley 6783bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 6793bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 6803bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 6813bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 6823bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 6833bf14a46SMatthew Knepley pastix->nSolve = 0; 6843bf14a46SMatthew Knepley pastix->MatDestroy = B->ops->destroy; 6853bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 6863bf14a46SMatthew Knepley B->spptr = (void*)pastix; 6873bf14a46SMatthew Knepley 6883bf14a46SMatthew Knepley *F = B; 6893bf14a46SMatthew Knepley PetscFunctionReturn(0); 6903bf14a46SMatthew Knepley } 6913bf14a46SMatthew Knepley EXTERN_C_END 6923bf14a46SMatthew Knepley 6933bf14a46SMatthew Knepley EXTERN_C_BEGIN 6943bf14a46SMatthew Knepley #undef __FUNCT__ 6953bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_seqsbaij_pastix" 6963bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 6973bf14a46SMatthew Knepley { 6983bf14a46SMatthew Knepley Mat B; 6993bf14a46SMatthew Knepley PetscErrorCode ierr; 7003bf14a46SMatthew Knepley Mat_Pastix *pastix; 7013bf14a46SMatthew Knepley 7023bf14a46SMatthew Knepley PetscFunctionBegin; 7033bf14a46SMatthew Knepley if (ftype != MAT_FACTOR_CHOLESKY) { 7043bf14a46SMatthew Knepley SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 7053bf14a46SMatthew Knepley } 7063bf14a46SMatthew Knepley /* Create the factorization matrix */ 7073bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 7083bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 7093bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 7103bf14a46SMatthew Knepley ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 7113bf14a46SMatthew Knepley ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 7123bf14a46SMatthew Knepley 7133bf14a46SMatthew Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 7143bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 7153bf14a46SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)B, 7163bf14a46SMatthew Knepley "MatFactorGetSolverPackage_C", 7173bf14a46SMatthew Knepley "MatFactorGetSolverPackage_pastix", 7183bf14a46SMatthew Knepley MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 7193bf14a46SMatthew Knepley 7203bf14a46SMatthew Knepley B->factor = MAT_FACTOR_CHOLESKY; 7213bf14a46SMatthew Knepley 7223bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 7233bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 7243bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 7253bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 7263bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 7273bf14a46SMatthew Knepley pastix->nSolve = 0; 7283bf14a46SMatthew Knepley pastix->MatDestroy = B->ops->destroy; 7293bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 7303bf14a46SMatthew Knepley B->spptr = (void*)pastix; 7313bf14a46SMatthew Knepley 7323bf14a46SMatthew Knepley *F = B; 7333bf14a46SMatthew Knepley PetscFunctionReturn(0); 7343bf14a46SMatthew Knepley } 7353bf14a46SMatthew Knepley EXTERN_C_END 7363bf14a46SMatthew Knepley 7373bf14a46SMatthew Knepley EXTERN_C_BEGIN 7383bf14a46SMatthew Knepley #undef __FUNCT__ 7393bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_mpisbaij_pastix" 7403bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 7413bf14a46SMatthew Knepley { 7423bf14a46SMatthew Knepley Mat B; 7433bf14a46SMatthew Knepley PetscErrorCode ierr; 7443bf14a46SMatthew Knepley Mat_Pastix *pastix; 7453bf14a46SMatthew Knepley 7463bf14a46SMatthew Knepley PetscFunctionBegin; 747*41c8de11SBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 748*41c8de11SBarry Smith 7493bf14a46SMatthew Knepley /* Create the factorization matrix */ 7503bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 7513bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 7523bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 7533bf14a46SMatthew Knepley ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 7543bf14a46SMatthew Knepley ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 7553bf14a46SMatthew Knepley 7563bf14a46SMatthew Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 7573bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 7583bf14a46SMatthew Knepley ierr = PetscObjectComposeFunctionDynamic((PetscObject)B, 7593bf14a46SMatthew Knepley "MatFactorGetSolverPackage_C", 7603bf14a46SMatthew Knepley "MatFactorGetSolverPackage_pastix", 7613bf14a46SMatthew Knepley MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 7623bf14a46SMatthew Knepley B->factor = MAT_FACTOR_CHOLESKY; 7633bf14a46SMatthew Knepley 7643bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 7653bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 7663bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 7673bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 7683bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 7693bf14a46SMatthew Knepley pastix->nSolve = 0; 7703bf14a46SMatthew Knepley pastix->MatDestroy = B->ops->destroy; 7713bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 7723bf14a46SMatthew Knepley B->spptr = (void*)pastix; 7733bf14a46SMatthew Knepley 7743bf14a46SMatthew Knepley *F = B; 7753bf14a46SMatthew Knepley PetscFunctionReturn(0); 7763bf14a46SMatthew Knepley } 7773bf14a46SMatthew Knepley EXTERN_C_END 778