13bf14a46SMatthew Knepley /* 23bf14a46SMatthew Knepley Provides an interface to the PaStiX sparse solver 33bf14a46SMatthew Knepley */ 4c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 5c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 6c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 7c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h> 83bf14a46SMatthew Knepley 943801a69SSatish Balay #if defined(PETSC_HAVE_STDLIB_H) 1043801a69SSatish Balay #include <stdlib.h> 1143801a69SSatish Balay #endif 1243801a69SSatish Balay #if defined(PETSC_HAVE_STRING_H) 1343801a69SSatish Balay #include <string.h> 1443801a69SSatish Balay #endif 1543801a69SSatish Balay 165ec454cfSSatish Balay #if defined(PETSC_USE_COMPLEX) && defined(__cplusplus) 175ec454cfSSatish Balay #define _H_COMPLEX 185ec454cfSSatish Balay #endif 195ec454cfSSatish Balay 203bf14a46SMatthew Knepley EXTERN_C_BEGIN 21c6db04a5SJed Brown #include <pastix.h> 223bf14a46SMatthew Knepley EXTERN_C_END 233bf14a46SMatthew Knepley 24519f805aSKarl Rupp #if defined(PETSC_USE_COMPLEX) 25519f805aSKarl Rupp #if defined(PETSC_USE_REAL_SINGLE) 265ec454cfSSatish Balay #define PASTIX_CALL c_pastix 275ec454cfSSatish Balay #define PASTIX_CHECKMATRIX c_pastix_checkMatrix 285ec454cfSSatish Balay #define PastixScalar COMPLEX 295ec454cfSSatish Balay #else 305ec454cfSSatish Balay #define PASTIX_CALL z_pastix 315ec454cfSSatish Balay #define PASTIX_CHECKMATRIX z_pastix_checkMatrix 325ec454cfSSatish Balay #define PastixScalar DCOMPLEX 335ec454cfSSatish Balay #endif 34d41469e0Sxavier lacoste 35d41469e0Sxavier lacoste #else /* PETSC_USE_COMPLEX */ 36d41469e0Sxavier lacoste 37519f805aSKarl Rupp #if defined(PETSC_USE_REAL_SINGLE) 385ec454cfSSatish Balay #define PASTIX_CALL s_pastix 395ec454cfSSatish Balay #define PASTIX_CHECKMATRIX s_pastix_checkMatrix 405ec454cfSSatish Balay #define PastixScalar float 415ec454cfSSatish Balay #else 425ec454cfSSatish Balay #define PASTIX_CALL d_pastix 435ec454cfSSatish Balay #define PASTIX_CHECKMATRIX d_pastix_checkMatrix 445ec454cfSSatish Balay #define PastixScalar double 455ec454cfSSatish Balay #endif 46d41469e0Sxavier lacoste 47d41469e0Sxavier lacoste #endif /* PETSC_USE_COMPLEX */ 48d41469e0Sxavier lacoste 493bf14a46SMatthew Knepley typedef struct Mat_Pastix_ { 503bf14a46SMatthew Knepley pastix_data_t *pastix_data; /* Pastix data storage structure */ 513bf14a46SMatthew Knepley MatStructure matstruc; 523bf14a46SMatthew Knepley PetscInt n; /* Number of columns in the matrix */ 533bf14a46SMatthew Knepley PetscInt *colptr; /* Index of first element of each column in row and val */ 543bf14a46SMatthew Knepley PetscInt *row; /* Row of each element of the matrix */ 553bf14a46SMatthew Knepley PetscScalar *val; /* Value of each element of the matrix */ 563bf14a46SMatthew Knepley PetscInt *perm; /* Permutation tabular */ 573bf14a46SMatthew Knepley PetscInt *invp; /* Reverse permutation tabular */ 583bf14a46SMatthew Knepley PetscScalar *rhs; /* Rhight-hand-side member */ 593bf14a46SMatthew Knepley PetscInt rhsnbr; /* Rhight-hand-side number (must be 1) */ 603bf14a46SMatthew Knepley PetscInt iparm[64]; /* Integer parameters */ 613bf14a46SMatthew Knepley double dparm[64]; /* Floating point parameters */ 623bf14a46SMatthew Knepley MPI_Comm pastix_comm; /* PaStiX MPI communicator */ 633bf14a46SMatthew Knepley PetscMPIInt commRank; /* MPI rank */ 643bf14a46SMatthew Knepley PetscMPIInt commSize; /* MPI communicator size */ 65ace3abfcSBarry Smith PetscBool CleanUpPastix; /* Boolean indicating if we call PaStiX clean step */ 663bf14a46SMatthew Knepley VecScatter scat_rhs; 673bf14a46SMatthew Knepley VecScatter scat_sol; 68f31ce8a6SBarry Smith Vec b_seq; 69ace3abfcSBarry Smith PetscBool isAIJ; 70bf0cc555SLisandro Dalcin PetscErrorCode (*Destroy)(Mat); 713bf14a46SMatthew Knepley } Mat_Pastix; 723bf14a46SMatthew Knepley 7309573ac7SBarry Smith extern PetscErrorCode MatDuplicate_Pastix(Mat,MatDuplicateOption,Mat*); 743bf14a46SMatthew Knepley 75eb1f6c34SBarry Smith #undef __FUNCT__ 76eb1f6c34SBarry Smith #define __FUNCT__ "MatConvertToCSC" 773bf14a46SMatthew Knepley /* 783bf14a46SMatthew Knepley convert Petsc seqaij matrix to CSC: colptr[n], row[nz], val[nz] 793bf14a46SMatthew Knepley 803bf14a46SMatthew Knepley input: 813bf14a46SMatthew Knepley A - matrix in seqaij or mpisbaij (bs=1) format 823bf14a46SMatthew Knepley valOnly - FALSE: spaces are allocated and values are set for the CSC 833bf14a46SMatthew Knepley TRUE: Only fill values 843bf14a46SMatthew Knepley output: 853bf14a46SMatthew Knepley n - Size of the matrix 863bf14a46SMatthew Knepley colptr - Index of first element of each column in row and val 873bf14a46SMatthew Knepley row - Row of each element of the matrix 883bf14a46SMatthew Knepley values - Value of each element of the matrix 893bf14a46SMatthew Knepley */ 90ace3abfcSBarry Smith PetscErrorCode MatConvertToCSC(Mat A,PetscBool valOnly,PetscInt *n,PetscInt **colptr,PetscInt **row,PetscScalar **values) 9141c8de11SBarry Smith { 923bf14a46SMatthew Knepley Mat_SeqAIJ *aa = (Mat_SeqAIJ*)A->data; 933bf14a46SMatthew Knepley PetscInt *rowptr = aa->i; 943bf14a46SMatthew Knepley PetscInt *col = aa->j; 953bf14a46SMatthew Knepley PetscScalar *rvalues = aa->a; 963bf14a46SMatthew Knepley PetscInt m = A->rmap->N; 97745c78f7SBarry Smith PetscInt nnz; 983bf14a46SMatthew Knepley PetscInt i,j, k; 993bf14a46SMatthew Knepley PetscInt base = 1; 1003bf14a46SMatthew Knepley PetscInt idx; 1013bf14a46SMatthew Knepley PetscErrorCode ierr; 1023bf14a46SMatthew Knepley PetscInt colidx; 1033bf14a46SMatthew Knepley PetscInt *colcount; 104ace3abfcSBarry Smith PetscBool isSBAIJ; 105ace3abfcSBarry Smith PetscBool isSeqSBAIJ; 106ace3abfcSBarry Smith PetscBool isMpiSBAIJ; 107ace3abfcSBarry Smith PetscBool isSym; 108d41469e0Sxavier lacoste PetscBool flg; 109d41469e0Sxavier lacoste PetscInt icntl; 110d41469e0Sxavier lacoste PetscInt verb; 111d41469e0Sxavier lacoste PetscInt check; 1123bf14a46SMatthew Knepley 1133bf14a46SMatthew Knepley PetscFunctionBegin; 11441c8de11SBarry Smith ierr = MatIsSymmetric(A,0.0,&isSym);CHKERRQ(ierr); 115251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSBAIJ,&isSBAIJ);CHKERRQ(ierr); 116251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 117251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPISBAIJ,&isMpiSBAIJ);CHKERRQ(ierr); 1183bf14a46SMatthew Knepley 119745c78f7SBarry Smith *n = A->cmap->N; 120745c78f7SBarry Smith 121745c78f7SBarry Smith /* PaStiX only needs triangular matrix if matrix is symmetric 122745c78f7SBarry Smith */ 123*2205254eSKarl Rupp if (isSym && !(isSBAIJ || isSeqSBAIJ || isMpiSBAIJ)) nnz = (aa->nz - *n)/2 + *n; 124*2205254eSKarl Rupp else nnz = aa->nz; 1253bf14a46SMatthew Knepley 1263bf14a46SMatthew Knepley if (!valOnly) { 1273bf14a46SMatthew Knepley ierr = PetscMalloc(((*n)+1) *sizeof(PetscInt) ,colptr);CHKERRQ(ierr); 1283bf14a46SMatthew Knepley ierr = PetscMalloc(nnz *sizeof(PetscInt) ,row);CHKERRQ(ierr); 1293bf14a46SMatthew Knepley ierr = PetscMalloc(nnz *sizeof(PetscScalar),values);CHKERRQ(ierr); 1303bf14a46SMatthew Knepley 13141c8de11SBarry Smith if (isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) { 13241c8de11SBarry Smith ierr = PetscMemcpy (*colptr, rowptr, ((*n)+1)*sizeof(PetscInt));CHKERRQ(ierr); 133*2205254eSKarl Rupp for (i = 0; i < *n+1; i++) (*colptr)[i] += base; 13441c8de11SBarry Smith ierr = PetscMemcpy (*row, col, (nnz)*sizeof(PetscInt));CHKERRQ(ierr); 135*2205254eSKarl Rupp for (i = 0; i < nnz; i++) (*row)[i] += base; 13641c8de11SBarry Smith ierr = PetscMemcpy (*values, rvalues, (nnz)*sizeof(PetscScalar));CHKERRQ(ierr); 13741c8de11SBarry Smith } else { 13841c8de11SBarry Smith ierr = PetscMalloc((*n)*sizeof(PetscInt),&colcount);CHKERRQ(ierr); 13941c8de11SBarry Smith 140f31ce8a6SBarry Smith for (i = 0; i < m; i++) colcount[i] = 0; 1413bf14a46SMatthew Knepley /* Fill-in colptr */ 142f31ce8a6SBarry Smith for (i = 0; i < m; i++) { 143f31ce8a6SBarry Smith for (j = rowptr[i]; j < rowptr[i+1]; j++) { 144f31ce8a6SBarry Smith if (!isSym || col[j] <= i) colcount[col[j]]++; 145f31ce8a6SBarry Smith } 146f31ce8a6SBarry Smith } 147745c78f7SBarry Smith 1483bf14a46SMatthew Knepley (*colptr)[0] = base; 1493bf14a46SMatthew Knepley for (j = 0; j < *n; j++) { 1503bf14a46SMatthew Knepley (*colptr)[j+1] = (*colptr)[j] + colcount[j]; 151745c78f7SBarry Smith /* in next loop we fill starting from (*colptr)[colidx] - base */ 1523bf14a46SMatthew Knepley colcount[j] = -base; 1533bf14a46SMatthew Knepley } 1543bf14a46SMatthew Knepley 1553bf14a46SMatthew Knepley /* Fill-in rows and values */ 1563bf14a46SMatthew Knepley for (i = 0; i < m; i++) { 1573bf14a46SMatthew Knepley for (j = rowptr[i]; j < rowptr[i+1]; j++) { 15841c8de11SBarry Smith if (!isSym || col[j] <= i) { 1593bf14a46SMatthew Knepley colidx = col[j]; 1603bf14a46SMatthew Knepley idx = (*colptr)[colidx] + colcount[colidx]; 1613bf14a46SMatthew Knepley (*row)[idx] = i + base; 1623bf14a46SMatthew Knepley (*values)[idx] = rvalues[j]; 1633bf14a46SMatthew Knepley colcount[colidx]++; 1643bf14a46SMatthew Knepley } 1653bf14a46SMatthew Knepley } 1663bf14a46SMatthew Knepley } 16741c8de11SBarry Smith ierr = PetscFree(colcount);CHKERRQ(ierr); 168745c78f7SBarry Smith } 16941c8de11SBarry Smith } else { 170745c78f7SBarry Smith /* Fill-in only values */ 1713bf14a46SMatthew Knepley for (i = 0; i < m; i++) { 1723bf14a46SMatthew Knepley for (j = rowptr[i]; j < rowptr[i+1]; j++) { 1733bf14a46SMatthew Knepley colidx = col[j]; 174*2205254eSKarl Rupp if ((isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) ||!isSym || col[j] <= i) { 175745c78f7SBarry Smith /* look for the value to fill */ 176f31ce8a6SBarry Smith for (k = (*colptr)[colidx] - base; k < (*colptr)[colidx + 1] - base; k++) { 177eb1f6c34SBarry Smith if (((*row)[k]-base) == i) { 1783bf14a46SMatthew Knepley (*values)[k] = rvalues[j]; 1793bf14a46SMatthew Knepley break; 1803bf14a46SMatthew Knepley } 1813bf14a46SMatthew Knepley } 182f31ce8a6SBarry Smith /* data structure of sparse matrix has changed */ 183e32f2f54SBarry Smith if (k == (*colptr)[colidx + 1] - base) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"overflow on k %D",k); 1843bf14a46SMatthew Knepley } 1853bf14a46SMatthew Knepley } 1863bf14a46SMatthew Knepley } 187745c78f7SBarry Smith } 188d41469e0Sxavier lacoste 189d41469e0Sxavier lacoste icntl =-1; 190d41469e0Sxavier lacoste check = 0; 191d41469e0Sxavier lacoste ierr = PetscOptionsInt("-mat_pastix_check","Check the matrix 0 : no, 1 : yes)","None",check,&icntl,&flg);CHKERRQ(ierr); 192*2205254eSKarl Rupp if ((flg && icntl >= 0) || PetscLogPrintInfo) check = icntl; 193*2205254eSKarl Rupp 194d41469e0Sxavier lacoste if (check == 1) { 19570fe17b1SSatish Balay PetscScalar *tmpvalues; 19670fe17b1SSatish Balay PetscInt *tmprows,*tmpcolptr; 1970f11a792SBarry Smith tmpvalues = (PetscScalar*)malloc(nnz*sizeof(PetscScalar)); if (!tmpvalues) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory"); 1980f11a792SBarry Smith tmprows = (PetscInt*) malloc(nnz*sizeof(PetscInt)); if (!tmprows) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory"); 1990f11a792SBarry Smith tmpcolptr = (PetscInt*) malloc((*n+1)*sizeof(PetscInt)); if (!tmpcolptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory"); 20043801a69SSatish Balay 20170fe17b1SSatish Balay ierr = PetscMemcpy(tmpcolptr,*colptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr); 20270fe17b1SSatish Balay ierr = PetscMemcpy(tmprows,*row,nnz*sizeof(PetscInt));CHKERRQ(ierr); 20370fe17b1SSatish Balay ierr = PetscMemcpy(tmpvalues,*values,nnz*sizeof(PetscScalar));CHKERRQ(ierr); 20443801a69SSatish Balay ierr = PetscFree(*row);CHKERRQ(ierr); 20543801a69SSatish Balay ierr = PetscFree(*values);CHKERRQ(ierr); 20643801a69SSatish Balay 207d41469e0Sxavier lacoste icntl=-1; 208d41469e0Sxavier lacoste verb = API_VERBOSE_NOT; 209d41469e0Sxavier lacoste ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",verb,&icntl,&flg);CHKERRQ(ierr); 210d41469e0Sxavier lacoste if ((flg && icntl >= 0) || PetscLogPrintInfo) { 211d41469e0Sxavier lacoste verb = icntl; 212d41469e0Sxavier lacoste } 2135ec454cfSSatish Balay PASTIX_CHECKMATRIX(MPI_COMM_WORLD,verb,((isSym != 0) ? API_SYM_YES : API_SYM_NO),API_YES,*n,&tmpcolptr,&tmprows,(PastixScalar**)&tmpvalues,NULL,1); 21443801a69SSatish Balay 21570fe17b1SSatish Balay ierr = PetscMemcpy(*colptr,tmpcolptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr); 21643801a69SSatish Balay ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscInt),row);CHKERRQ(ierr); 21770fe17b1SSatish Balay ierr = PetscMemcpy(*row,tmprows,((*colptr)[*n]-1)*sizeof(PetscInt));CHKERRQ(ierr); 21843801a69SSatish Balay ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscScalar),values);CHKERRQ(ierr); 21970fe17b1SSatish Balay ierr = PetscMemcpy(*values,tmpvalues,((*colptr)[*n]-1)*sizeof(PetscScalar));CHKERRQ(ierr); 220be76a908SBarry Smith free(tmpvalues); 221be76a908SBarry Smith free(tmprows); 222be76a908SBarry Smith free(tmpcolptr); 223be76a908SBarry Smith 22443801a69SSatish Balay } 2253bf14a46SMatthew Knepley PetscFunctionReturn(0); 2263bf14a46SMatthew Knepley } 2273bf14a46SMatthew Knepley 2283bf14a46SMatthew Knepley 2293bf14a46SMatthew Knepley 2303bf14a46SMatthew Knepley #undef __FUNCT__ 2313bf14a46SMatthew Knepley #define __FUNCT__ "MatDestroy_Pastix" 2323bf14a46SMatthew Knepley /* 2333bf14a46SMatthew Knepley Call clean step of PaStiX if lu->CleanUpPastix == true. 2343bf14a46SMatthew Knepley Free the CSC matrix. 2353bf14a46SMatthew Knepley */ 2363bf14a46SMatthew Knepley PetscErrorCode MatDestroy_Pastix(Mat A) 2373bf14a46SMatthew Knepley { 2383bf14a46SMatthew Knepley Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 2393bf14a46SMatthew Knepley PetscErrorCode ierr; 2403bf14a46SMatthew Knepley PetscMPIInt size=lu->commSize; 241745c78f7SBarry Smith 2423bf14a46SMatthew Knepley PetscFunctionBegin; 243bf0cc555SLisandro Dalcin if (lu && lu->CleanUpPastix) { 2443bf14a46SMatthew Knepley /* Terminate instance, deallocate memories */ 2453bf14a46SMatthew Knepley if (size > 1) { 2466bf464f9SBarry Smith ierr = VecScatterDestroy(&lu->scat_rhs);CHKERRQ(ierr); 2476bf464f9SBarry Smith ierr = VecDestroy(&lu->b_seq);CHKERRQ(ierr); 2486bf464f9SBarry Smith ierr = VecScatterDestroy(&lu->scat_sol);CHKERRQ(ierr); 2493bf14a46SMatthew Knepley } 2503bf14a46SMatthew Knepley 2513bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK]=API_TASK_CLEAN; 2523bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] =API_TASK_CLEAN; 2533bf14a46SMatthew Knepley 254d41469e0Sxavier lacoste PASTIX_CALL(&(lu->pastix_data), 2553bf14a46SMatthew Knepley lu->pastix_comm, 256d41469e0Sxavier lacoste lu->n, 257d41469e0Sxavier lacoste lu->colptr, 258d41469e0Sxavier lacoste lu->row, 2595ec454cfSSatish Balay (PastixScalar*)lu->val, 260d41469e0Sxavier lacoste lu->perm, 261d41469e0Sxavier lacoste lu->invp, 2625ec454cfSSatish Balay (PastixScalar*)lu->rhs, 263d41469e0Sxavier lacoste lu->rhsnbr, 264d41469e0Sxavier lacoste lu->iparm, 2653bf14a46SMatthew Knepley lu->dparm); 2663bf14a46SMatthew Knepley 2673bf14a46SMatthew Knepley ierr = PetscFree(lu->colptr);CHKERRQ(ierr); 2683bf14a46SMatthew Knepley ierr = PetscFree(lu->row);CHKERRQ(ierr); 2693bf14a46SMatthew Knepley ierr = PetscFree(lu->val);CHKERRQ(ierr); 2703bf14a46SMatthew Knepley ierr = PetscFree(lu->perm);CHKERRQ(ierr); 2713bf14a46SMatthew Knepley ierr = PetscFree(lu->invp);CHKERRQ(ierr); 2723bf14a46SMatthew Knepley ierr = MPI_Comm_free(&(lu->pastix_comm));CHKERRQ(ierr); 2733bf14a46SMatthew Knepley } 274bf0cc555SLisandro Dalcin if (lu && lu->Destroy) { 275bf0cc555SLisandro Dalcin ierr = (lu->Destroy)(A);CHKERRQ(ierr); 276bf0cc555SLisandro Dalcin } 277bf0cc555SLisandro Dalcin ierr = PetscFree(A->spptr);CHKERRQ(ierr); 2783bf14a46SMatthew Knepley PetscFunctionReturn(0); 2793bf14a46SMatthew Knepley } 2803bf14a46SMatthew Knepley 2813bf14a46SMatthew Knepley #undef __FUNCT__ 2823bf14a46SMatthew Knepley #define __FUNCT__ "MatSolve_PaStiX" 2833bf14a46SMatthew Knepley /* 2843bf14a46SMatthew Knepley Gather right-hand-side. 2853bf14a46SMatthew Knepley Call for Solve step. 2863bf14a46SMatthew Knepley Scatter solution. 2873bf14a46SMatthew Knepley */ 2883bf14a46SMatthew Knepley PetscErrorCode MatSolve_PaStiX(Mat A,Vec b,Vec x) 2893bf14a46SMatthew Knepley { 2903bf14a46SMatthew Knepley Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 2913bf14a46SMatthew Knepley PetscScalar *array; 2923bf14a46SMatthew Knepley Vec x_seq; 2933bf14a46SMatthew Knepley PetscErrorCode ierr; 2943bf14a46SMatthew Knepley 2953bf14a46SMatthew Knepley PetscFunctionBegin; 2963bf14a46SMatthew Knepley lu->rhsnbr = 1; 2973bf14a46SMatthew Knepley x_seq = lu->b_seq; 2983bf14a46SMatthew Knepley if (lu->commSize > 1) { 2993bf14a46SMatthew Knepley /* PaStiX only supports centralized rhs. Scatter b into a seqential rhs vector */ 3003bf14a46SMatthew Knepley ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3013bf14a46SMatthew Knepley ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 302b5e56a35SBarry Smith ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr); 30341c8de11SBarry Smith } else { /* size == 1 */ 3043bf14a46SMatthew Knepley ierr = VecCopy(b,x);CHKERRQ(ierr); 3053bf14a46SMatthew Knepley ierr = VecGetArray(x,&array);CHKERRQ(ierr); 3063bf14a46SMatthew Knepley } 3073bf14a46SMatthew Knepley lu->rhs = array; 3083bf14a46SMatthew Knepley if (lu->commSize == 1) { 3093bf14a46SMatthew Knepley ierr = VecRestoreArray(x,&array);CHKERRQ(ierr); 3103bf14a46SMatthew Knepley } else { 3113bf14a46SMatthew Knepley ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr); 3123bf14a46SMatthew Knepley } 3133bf14a46SMatthew Knepley 3143bf14a46SMatthew Knepley /* solve phase */ 3153bf14a46SMatthew Knepley /*-------------*/ 3163bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE; 3173bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_REFINE; 318745c78f7SBarry Smith lu->iparm[IPARM_RHS_MAKING] = API_RHS_B; 3193bf14a46SMatthew Knepley 320d41469e0Sxavier lacoste PASTIX_CALL(&(lu->pastix_data), 321d41469e0Sxavier lacoste lu->pastix_comm, 322d41469e0Sxavier lacoste lu->n, 323d41469e0Sxavier lacoste lu->colptr, 324d41469e0Sxavier lacoste lu->row, 3255ec454cfSSatish Balay (PastixScalar*)lu->val, 326d41469e0Sxavier lacoste lu->perm, 327d41469e0Sxavier lacoste lu->invp, 3285ec454cfSSatish Balay (PastixScalar*)lu->rhs, 329d41469e0Sxavier lacoste lu->rhsnbr, 330d41469e0Sxavier lacoste lu->iparm, 331d41469e0Sxavier lacoste lu->dparm); 3323bf14a46SMatthew Knepley 33365e19b50SBarry Smith 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]); 3343bf14a46SMatthew Knepley 3353bf14a46SMatthew Knepley if (lu->commSize == 1) { 3363bf14a46SMatthew Knepley ierr = VecRestoreArray(x,&(lu->rhs));CHKERRQ(ierr); 3373bf14a46SMatthew Knepley } else { 3383bf14a46SMatthew Knepley ierr = VecRestoreArray(x_seq,&(lu->rhs));CHKERRQ(ierr); 3393bf14a46SMatthew Knepley } 3403bf14a46SMatthew Knepley 3413bf14a46SMatthew Knepley if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */ 3423bf14a46SMatthew Knepley ierr = VecScatterBegin(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3433bf14a46SMatthew Knepley ierr = VecScatterEnd(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 3443bf14a46SMatthew Knepley } 3453bf14a46SMatthew Knepley PetscFunctionReturn(0); 3463bf14a46SMatthew Knepley } 3473bf14a46SMatthew Knepley 3483bf14a46SMatthew Knepley /* 3493bf14a46SMatthew Knepley Numeric factorisation using PaStiX solver. 3503bf14a46SMatthew Knepley 3513bf14a46SMatthew Knepley */ 3523bf14a46SMatthew Knepley #undef __FUNCT__ 35353c77d0aSJed Brown #define __FUNCT__ "MatFactorNumeric_PaStiX" 3543bf14a46SMatthew Knepley PetscErrorCode MatFactorNumeric_PaStiX(Mat F,Mat A,const MatFactorInfo *info) 3553bf14a46SMatthew Knepley { 3563bf14a46SMatthew Knepley Mat_Pastix *lu =(Mat_Pastix*)(F)->spptr; 35741c8de11SBarry Smith Mat *tseq; 3583bf14a46SMatthew Knepley PetscErrorCode ierr = 0; 359b5e56a35SBarry Smith PetscInt icntl; 360b5e56a35SBarry Smith PetscInt M=A->rmap->N; 361ace3abfcSBarry Smith PetscBool valOnly,flg, isSym; 3623bf14a46SMatthew Knepley Mat F_diag; 3633bf14a46SMatthew Knepley IS is_iden; 3643bf14a46SMatthew Knepley Vec b; 3653bf14a46SMatthew Knepley IS isrow; 36651a30905SBarry Smith PetscBool isSeqAIJ,isSeqSBAIJ,isMPIAIJ; 3673bf14a46SMatthew Knepley 3683bf14a46SMatthew Knepley PetscFunctionBegin; 369251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr); 370251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr); 371251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr); 3723bf14a46SMatthew Knepley if (lu->matstruc == DIFFERENT_NONZERO_PATTERN) { 3733bf14a46SMatthew Knepley (F)->ops->solve = MatSolve_PaStiX; 3743bf14a46SMatthew Knepley 3753bf14a46SMatthew Knepley /* Initialize a PASTIX instance */ 3763bf14a46SMatthew Knepley ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr); 3773bf14a46SMatthew Knepley ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank);CHKERRQ(ierr); 3783bf14a46SMatthew Knepley ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize);CHKERRQ(ierr); 3793bf14a46SMatthew Knepley 3803bf14a46SMatthew Knepley /* Set pastix options */ 3813bf14a46SMatthew Knepley lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO; 3823bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_INIT; 3833bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_INIT; 384*2205254eSKarl Rupp 3853bf14a46SMatthew Knepley lu->rhsnbr = 1; 3863bf14a46SMatthew Knepley 3873bf14a46SMatthew Knepley /* Call to set default pastix options */ 388d41469e0Sxavier lacoste PASTIX_CALL(&(lu->pastix_data), 389d41469e0Sxavier lacoste lu->pastix_comm, 390d41469e0Sxavier lacoste lu->n, 391d41469e0Sxavier lacoste lu->colptr, 392d41469e0Sxavier lacoste lu->row, 3935ec454cfSSatish Balay (PastixScalar*)lu->val, 394d41469e0Sxavier lacoste lu->perm, 395d41469e0Sxavier lacoste lu->invp, 3965ec454cfSSatish Balay (PastixScalar*)lu->rhs, 397d41469e0Sxavier lacoste lu->rhsnbr, 398d41469e0Sxavier lacoste lu->iparm, 399d41469e0Sxavier lacoste lu->dparm); 4003bf14a46SMatthew Knepley 4013bf14a46SMatthew Knepley ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr); 4023bf14a46SMatthew Knepley 4033bf14a46SMatthew Knepley icntl = -1; 404*2205254eSKarl Rupp 40541c8de11SBarry Smith lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT; 406*2205254eSKarl Rupp 40741c8de11SBarry Smith ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr); 408d41469e0Sxavier lacoste if ((flg && icntl >= 0) || PetscLogPrintInfo) { 4093bf14a46SMatthew Knepley lu->iparm[IPARM_VERBOSE] = icntl; 4103bf14a46SMatthew Knepley } 4113bf14a46SMatthew Knepley icntl=-1; 412e4e47003SBarry Smith ierr = PetscOptionsInt("-mat_pastix_threadnbr","iparm[IPARM_THREAD_NBR] : Number of thread by MPI node","None",lu->iparm[IPARM_THREAD_NBR],&icntl,&flg);CHKERRQ(ierr); 4133bf14a46SMatthew Knepley if ((flg && icntl > 0)) { 4143bf14a46SMatthew Knepley lu->iparm[IPARM_THREAD_NBR] = icntl; 4153bf14a46SMatthew Knepley } 4163bf14a46SMatthew Knepley PetscOptionsEnd(); 4173bf14a46SMatthew Knepley valOnly = PETSC_FALSE; 41841c8de11SBarry Smith } else { 4195d6241c8SBarry Smith if (isSeqAIJ || isMPIAIJ) { 4205d6241c8SBarry Smith ierr = PetscFree(lu->colptr);CHKERRQ(ierr); 4215d6241c8SBarry Smith ierr = PetscFree(lu->row);CHKERRQ(ierr); 4225d6241c8SBarry Smith ierr = PetscFree(lu->val);CHKERRQ(ierr); 4235d6241c8SBarry Smith valOnly = PETSC_FALSE; 4245d6241c8SBarry Smith } else valOnly = PETSC_TRUE; 4253bf14a46SMatthew Knepley } 4263bf14a46SMatthew Knepley 4273bf14a46SMatthew Knepley lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES; 4283bf14a46SMatthew Knepley 4293bf14a46SMatthew Knepley /* convert mpi A to seq mat A */ 4303bf14a46SMatthew Knepley ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr); 4313bf14a46SMatthew Knepley ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr); 4326bf464f9SBarry Smith ierr = ISDestroy(&isrow);CHKERRQ(ierr); 4333bf14a46SMatthew Knepley 43441c8de11SBarry Smith ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val);CHKERRQ(ierr); 43541c8de11SBarry Smith ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr); 43641c8de11SBarry Smith ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr); 43741c8de11SBarry Smith 4385d6241c8SBarry Smith if (!lu->perm) { 4393bf14a46SMatthew Knepley ierr = PetscMalloc((lu->n)*sizeof(PetscInt),&(lu->perm));CHKERRQ(ierr); 4403bf14a46SMatthew Knepley ierr = PetscMalloc((lu->n)*sizeof(PetscInt),&(lu->invp));CHKERRQ(ierr); 4415d6241c8SBarry Smith } 4423bf14a46SMatthew Knepley 4433bf14a46SMatthew Knepley if (isSym) { 444745c78f7SBarry Smith /* On symmetric matrix, LLT */ 4453bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_YES; 44641c8de11SBarry Smith lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT; 447f31ce8a6SBarry Smith } else { 448745c78f7SBarry Smith /* On unsymmetric matrix, LU */ 4493bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_NO; 4503bf14a46SMatthew Knepley lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 4513bf14a46SMatthew Knepley } 4523bf14a46SMatthew Knepley 4533bf14a46SMatthew Knepley /*----------------*/ 4543bf14a46SMatthew Knepley if (lu->matstruc == DIFFERENT_NONZERO_PATTERN) { 4553bf14a46SMatthew Knepley if (!(isSeqAIJ || isSeqSBAIJ)) { 4563bf14a46SMatthew Knepley /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */ 4573bf14a46SMatthew Knepley ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr); 4583bf14a46SMatthew Knepley ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr); 4593bf14a46SMatthew Knepley ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr); 4603bf14a46SMatthew Knepley ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 4613bf14a46SMatthew Knepley ierr = VecSetFromOptions(b);CHKERRQ(ierr); 4623bf14a46SMatthew Knepley 4633bf14a46SMatthew Knepley ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr); 4643bf14a46SMatthew Knepley ierr = VecScatterCreate(lu->b_seq,is_iden,b,is_iden,&lu->scat_sol);CHKERRQ(ierr); 4656bf464f9SBarry Smith ierr = ISDestroy(&is_iden);CHKERRQ(ierr); 4666bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 4673bf14a46SMatthew Knepley } 4683bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING; 4693bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 4703bf14a46SMatthew Knepley 471d41469e0Sxavier lacoste PASTIX_CALL(&(lu->pastix_data), 472d41469e0Sxavier lacoste lu->pastix_comm, 473d41469e0Sxavier lacoste lu->n, 474d41469e0Sxavier lacoste lu->colptr, 475d41469e0Sxavier lacoste lu->row, 4765ec454cfSSatish Balay (PastixScalar*)lu->val, 477d41469e0Sxavier lacoste lu->perm, 478d41469e0Sxavier lacoste lu->invp, 4795ec454cfSSatish Balay (PastixScalar*)lu->rhs, 480d41469e0Sxavier lacoste lu->rhsnbr, 481d41469e0Sxavier lacoste lu->iparm, 482d41469e0Sxavier lacoste lu->dparm); 48365e19b50SBarry Smith 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]); 48441c8de11SBarry Smith } else { 4853bf14a46SMatthew Knepley lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT; 4863bf14a46SMatthew Knepley lu->iparm[IPARM_END_TASK] = API_TASK_NUMFACT; 487d41469e0Sxavier lacoste PASTIX_CALL(&(lu->pastix_data), 488d41469e0Sxavier lacoste lu->pastix_comm, 489d41469e0Sxavier lacoste lu->n, 490d41469e0Sxavier lacoste lu->colptr, 491d41469e0Sxavier lacoste lu->row, 4925ec454cfSSatish Balay (PastixScalar*)lu->val, 493d41469e0Sxavier lacoste lu->perm, 494d41469e0Sxavier lacoste lu->invp, 4955ec454cfSSatish Balay (PastixScalar*)lu->rhs, 496d41469e0Sxavier lacoste lu->rhsnbr, 497d41469e0Sxavier lacoste lu->iparm, 498d41469e0Sxavier lacoste lu->dparm); 4993bf14a46SMatthew Knepley 50065e19b50SBarry Smith 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]); 5013bf14a46SMatthew Knepley } 5023bf14a46SMatthew Knepley 5033bf14a46SMatthew Knepley if (lu->commSize > 1) { 504d5f3da31SBarry Smith if ((F)->factortype == MAT_FACTOR_LU) { 5053bf14a46SMatthew Knepley F_diag = ((Mat_MPIAIJ *)(F)->data)->A; 5063bf14a46SMatthew Knepley } else { 5073bf14a46SMatthew Knepley F_diag = ((Mat_MPISBAIJ *)(F)->data)->A; 5083bf14a46SMatthew Knepley } 5093bf14a46SMatthew Knepley F_diag->assembled = PETSC_TRUE; 5103bf14a46SMatthew Knepley } 5113bf14a46SMatthew Knepley (F)->assembled = PETSC_TRUE; 5123bf14a46SMatthew Knepley lu->matstruc = SAME_NONZERO_PATTERN; 5133bf14a46SMatthew Knepley lu->CleanUpPastix = PETSC_TRUE; 5143bf14a46SMatthew Knepley PetscFunctionReturn(0); 5153bf14a46SMatthew Knepley } 5163bf14a46SMatthew Knepley 5173bf14a46SMatthew Knepley /* Note the Petsc r and c permutations are ignored */ 5183bf14a46SMatthew Knepley #undef __FUNCT__ 5193bf14a46SMatthew Knepley #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX" 5203bf14a46SMatthew Knepley PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info) 5213bf14a46SMatthew Knepley { 5223bf14a46SMatthew Knepley Mat_Pastix *lu = (Mat_Pastix*)F->spptr; 5233bf14a46SMatthew Knepley 5243bf14a46SMatthew Knepley PetscFunctionBegin; 5253bf14a46SMatthew Knepley lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU; 5263bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_YES; 5273bf14a46SMatthew Knepley lu->matstruc = DIFFERENT_NONZERO_PATTERN; 5283bf14a46SMatthew Knepley F->ops->lufactornumeric = MatFactorNumeric_PaStiX; 5293bf14a46SMatthew Knepley PetscFunctionReturn(0); 5303bf14a46SMatthew Knepley } 5313bf14a46SMatthew Knepley 5323bf14a46SMatthew Knepley 5333bf14a46SMatthew Knepley /* Note the Petsc r permutation is ignored */ 5343bf14a46SMatthew Knepley #undef __FUNCT__ 5353bf14a46SMatthew Knepley #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX" 5363bf14a46SMatthew Knepley PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info) 5373bf14a46SMatthew Knepley { 5383bf14a46SMatthew Knepley Mat_Pastix *lu = (Mat_Pastix*)(F)->spptr; 5393bf14a46SMatthew Knepley 5403bf14a46SMatthew Knepley PetscFunctionBegin; 5413bf14a46SMatthew Knepley lu->iparm[IPARM_FACTORIZATION] = API_FACT_LLT; 5423bf14a46SMatthew Knepley lu->iparm[IPARM_SYM] = API_SYM_NO; 5433bf14a46SMatthew Knepley lu->matstruc = DIFFERENT_NONZERO_PATTERN; 5443bf14a46SMatthew Knepley (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX; 5453bf14a46SMatthew Knepley PetscFunctionReturn(0); 5463bf14a46SMatthew Knepley } 5473bf14a46SMatthew Knepley 5483bf14a46SMatthew Knepley #undef __FUNCT__ 5493bf14a46SMatthew Knepley #define __FUNCT__ "MatView_PaStiX" 5503bf14a46SMatthew Knepley PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer) 5513bf14a46SMatthew Knepley { 5523bf14a46SMatthew Knepley PetscErrorCode ierr; 553ace3abfcSBarry Smith PetscBool iascii; 5543bf14a46SMatthew Knepley PetscViewerFormat format; 5553bf14a46SMatthew Knepley 5563bf14a46SMatthew Knepley PetscFunctionBegin; 557251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5583bf14a46SMatthew Knepley if (iascii) { 5593bf14a46SMatthew Knepley ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 5603bf14a46SMatthew Knepley if (format == PETSC_VIEWER_ASCII_INFO) { 561b5e56a35SBarry Smith Mat_Pastix *lu=(Mat_Pastix*)A->spptr; 562b5e56a35SBarry Smith 563b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr); 564b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Matrix type : %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES) ? "Symmetric" : "Unsymmetric"));CHKERRQ(ierr); 565b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Level of printing (0,1,2): %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr); 566b5e56a35SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr); 567b5e56a35SBarry Smith ierr = PetscPrintf(PETSC_COMM_SELF," Error : %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr); 5683bf14a46SMatthew Knepley } 5693bf14a46SMatthew Knepley } 5703bf14a46SMatthew Knepley PetscFunctionReturn(0); 5713bf14a46SMatthew Knepley } 5723bf14a46SMatthew Knepley 5733bf14a46SMatthew Knepley 5743bf14a46SMatthew Knepley /*MC 5752692d6eeSBarry Smith MATSOLVERPASTIX - A solver package providing direct solvers (LU) for distributed 5763bf14a46SMatthew Knepley and sequential matrices via the external package PaStiX. 5773bf14a46SMatthew Knepley 578e2e64c6bSBarry Smith Use ./configure --download-pastix to have PETSc installed with PaStiX 5793bf14a46SMatthew Knepley 5803bf14a46SMatthew Knepley Options Database Keys: 581b5e56a35SBarry Smith + -mat_pastix_verbose <0,1,2> - print level 582b5e56a35SBarry Smith - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task. 5833bf14a46SMatthew Knepley 5843bf14a46SMatthew Knepley Level: beginner 5853bf14a46SMatthew Knepley 58641c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage 58741c8de11SBarry Smith 5883bf14a46SMatthew Knepley M*/ 5893bf14a46SMatthew Knepley 5903bf14a46SMatthew Knepley 5913bf14a46SMatthew Knepley #undef __FUNCT__ 5923bf14a46SMatthew Knepley #define __FUNCT__ "MatGetInfo_PaStiX" 5933bf14a46SMatthew Knepley PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info) 5943bf14a46SMatthew Knepley { 5953bf14a46SMatthew Knepley Mat_Pastix *lu =(Mat_Pastix*)A->spptr; 5963bf14a46SMatthew Knepley 5973bf14a46SMatthew Knepley PetscFunctionBegin; 5983bf14a46SMatthew Knepley info->block_size = 1.0; 5993bf14a46SMatthew Knepley info->nz_allocated = lu->iparm[IPARM_NNZEROS]; 6003bf14a46SMatthew Knepley info->nz_used = lu->iparm[IPARM_NNZEROS]; 6013bf14a46SMatthew Knepley info->nz_unneeded = 0.0; 6023bf14a46SMatthew Knepley info->assemblies = 0.0; 6033bf14a46SMatthew Knepley info->mallocs = 0.0; 6043bf14a46SMatthew Knepley info->memory = 0.0; 6053bf14a46SMatthew Knepley info->fill_ratio_given = 0; 6063bf14a46SMatthew Knepley info->fill_ratio_needed = 0; 6073bf14a46SMatthew Knepley info->factor_mallocs = 0; 6083bf14a46SMatthew Knepley PetscFunctionReturn(0); 6093bf14a46SMatthew Knepley } 6103bf14a46SMatthew Knepley 6113bf14a46SMatthew Knepley EXTERN_C_BEGIN 6123bf14a46SMatthew Knepley #undef __FUNCT__ 6133bf14a46SMatthew Knepley #define __FUNCT__ "MatFactorGetSolverPackage_pastix" 6143bf14a46SMatthew Knepley PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type) 6153bf14a46SMatthew Knepley { 6163bf14a46SMatthew Knepley PetscFunctionBegin; 6172692d6eeSBarry Smith *type = MATSOLVERPASTIX; 6183bf14a46SMatthew Knepley PetscFunctionReturn(0); 6193bf14a46SMatthew Knepley } 6203bf14a46SMatthew Knepley EXTERN_C_END 6213bf14a46SMatthew Knepley 6223bf14a46SMatthew Knepley EXTERN_C_BEGIN 6233bf14a46SMatthew Knepley /* 6243bf14a46SMatthew Knepley The seq and mpi versions of this function are the same 6253bf14a46SMatthew Knepley */ 6263bf14a46SMatthew Knepley #undef __FUNCT__ 6273bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_seqaij_pastix" 6283bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F) 6293bf14a46SMatthew Knepley { 6303bf14a46SMatthew Knepley Mat B; 6313bf14a46SMatthew Knepley PetscErrorCode ierr; 6323bf14a46SMatthew Knepley Mat_Pastix *pastix; 6333bf14a46SMatthew Knepley 6343bf14a46SMatthew Knepley PetscFunctionBegin; 635e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 6363bf14a46SMatthew Knepley /* Create the factorization matrix */ 6373bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 6383bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 6393bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 6403bf14a46SMatthew Knepley ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 6413bf14a46SMatthew Knepley 6423bf14a46SMatthew Knepley B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 6433bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 6443bf14a46SMatthew Knepley B->ops->getinfo = MatGetInfo_PaStiX; 645*2205254eSKarl Rupp 64631e762f5SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix", MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 647*2205254eSKarl Rupp 648d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 6493bf14a46SMatthew Knepley 6503bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 651*2205254eSKarl Rupp 6523bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 6533bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 6543bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 6553bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 656bf0cc555SLisandro Dalcin pastix->Destroy = B->ops->destroy; 6573bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 6583bf14a46SMatthew Knepley B->spptr = (void*)pastix; 6593bf14a46SMatthew Knepley 6603bf14a46SMatthew Knepley *F = B; 6613bf14a46SMatthew Knepley PetscFunctionReturn(0); 6623bf14a46SMatthew Knepley } 6633bf14a46SMatthew Knepley EXTERN_C_END 6643bf14a46SMatthew Knepley 665b5e56a35SBarry Smith 6663bf14a46SMatthew Knepley EXTERN_C_BEGIN 6673bf14a46SMatthew Knepley #undef __FUNCT__ 6683bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_mpiaij_pastix" 6693bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F) 6703bf14a46SMatthew Knepley { 6713bf14a46SMatthew Knepley Mat B; 6723bf14a46SMatthew Knepley PetscErrorCode ierr; 6733bf14a46SMatthew Knepley Mat_Pastix *pastix; 6743bf14a46SMatthew Knepley 6753bf14a46SMatthew Knepley PetscFunctionBegin; 676e32f2f54SBarry Smith if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix"); 6773bf14a46SMatthew Knepley /* Create the factorization matrix */ 6783bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 6793bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 6803bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 6813bf14a46SMatthew Knepley ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr); 6823bf14a46SMatthew Knepley ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 6833bf14a46SMatthew Knepley 6843bf14a46SMatthew Knepley B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX; 6853bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 686*2205254eSKarl Rupp 68731e762f5SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 688*2205254eSKarl Rupp 689d5f3da31SBarry Smith B->factortype = MAT_FACTOR_LU; 6903bf14a46SMatthew Knepley 6913bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 692*2205254eSKarl Rupp 6933bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 6943bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 6953bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 6963bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 697bf0cc555SLisandro Dalcin pastix->Destroy = B->ops->destroy; 6983bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 6993bf14a46SMatthew Knepley B->spptr = (void*)pastix; 7003bf14a46SMatthew Knepley 7013bf14a46SMatthew Knepley *F = B; 7023bf14a46SMatthew Knepley PetscFunctionReturn(0); 7033bf14a46SMatthew Knepley } 7043bf14a46SMatthew Knepley EXTERN_C_END 7053bf14a46SMatthew Knepley 7063bf14a46SMatthew Knepley EXTERN_C_BEGIN 7073bf14a46SMatthew Knepley #undef __FUNCT__ 7083bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_seqsbaij_pastix" 7093bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 7103bf14a46SMatthew Knepley { 7113bf14a46SMatthew Knepley Mat B; 7123bf14a46SMatthew Knepley PetscErrorCode ierr; 7133bf14a46SMatthew Knepley Mat_Pastix *pastix; 7143bf14a46SMatthew Knepley 7153bf14a46SMatthew Knepley PetscFunctionBegin; 716e7e72b3dSBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 7173bf14a46SMatthew Knepley /* Create the factorization matrix */ 7183bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 7193bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 7203bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 7213bf14a46SMatthew Knepley ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 7223bf14a46SMatthew Knepley ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 7233bf14a46SMatthew Knepley 7243bf14a46SMatthew Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 7253bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 726*2205254eSKarl Rupp 72731e762f5SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 728*2205254eSKarl Rupp 729d5f3da31SBarry Smith B->factortype = MAT_FACTOR_CHOLESKY; 7303bf14a46SMatthew Knepley 7313bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 732*2205254eSKarl Rupp 7333bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 7343bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 7353bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 7363bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 737bf0cc555SLisandro Dalcin pastix->Destroy = B->ops->destroy; 7383bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 7393bf14a46SMatthew Knepley B->spptr = (void*)pastix; 7403bf14a46SMatthew Knepley 7413bf14a46SMatthew Knepley *F = B; 7423bf14a46SMatthew Knepley PetscFunctionReturn(0); 7433bf14a46SMatthew Knepley } 7443bf14a46SMatthew Knepley EXTERN_C_END 7453bf14a46SMatthew Knepley 7463bf14a46SMatthew Knepley EXTERN_C_BEGIN 7473bf14a46SMatthew Knepley #undef __FUNCT__ 7483bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_mpisbaij_pastix" 7493bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F) 7503bf14a46SMatthew Knepley { 7513bf14a46SMatthew Knepley Mat B; 7523bf14a46SMatthew Knepley PetscErrorCode ierr; 7533bf14a46SMatthew Knepley Mat_Pastix *pastix; 7543bf14a46SMatthew Knepley 7553bf14a46SMatthew Knepley PetscFunctionBegin; 756e32f2f54SBarry Smith if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix"); 75741c8de11SBarry Smith 7583bf14a46SMatthew Knepley /* Create the factorization matrix */ 7593bf14a46SMatthew Knepley ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 7603bf14a46SMatthew Knepley ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 7613bf14a46SMatthew Knepley ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr); 7623bf14a46SMatthew Knepley ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 7633bf14a46SMatthew Knepley ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 7643bf14a46SMatthew Knepley 7653bf14a46SMatthew Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX; 7663bf14a46SMatthew Knepley B->ops->view = MatView_PaStiX; 767*2205254eSKarl Rupp 76831e762f5SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr); 769*2205254eSKarl Rupp 770d5f3da31SBarry Smith B->factortype = MAT_FACTOR_CHOLESKY; 7713bf14a46SMatthew Knepley 7723bf14a46SMatthew Knepley ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr); 773*2205254eSKarl Rupp 7743bf14a46SMatthew Knepley pastix->CleanUpPastix = PETSC_FALSE; 7753bf14a46SMatthew Knepley pastix->isAIJ = PETSC_TRUE; 7763bf14a46SMatthew Knepley pastix->scat_rhs = PETSC_NULL; 7773bf14a46SMatthew Knepley pastix->scat_sol = PETSC_NULL; 778bf0cc555SLisandro Dalcin pastix->Destroy = B->ops->destroy; 7793bf14a46SMatthew Knepley B->ops->destroy = MatDestroy_Pastix; 7803bf14a46SMatthew Knepley B->spptr = (void*)pastix; 7813bf14a46SMatthew Knepley 7823bf14a46SMatthew Knepley *F = B; 7833bf14a46SMatthew Knepley PetscFunctionReturn(0); 7843bf14a46SMatthew Knepley } 7853bf14a46SMatthew Knepley EXTERN_C_END 786