xref: /petsc/src/mat/impls/aij/mpi/pastix/pastix.c (revision 1cc06b555e92f8ec64db10330b8bbd830e5bc876)
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 
9dbbbd53dSSatish Balay #if defined(PETSC_USE_COMPLEX)
105ec454cfSSatish Balay   #define _H_COMPLEX
115ec454cfSSatish Balay #endif
125ec454cfSSatish Balay 
133bf14a46SMatthew Knepley EXTERN_C_BEGIN
14c6db04a5SJed Brown #include <pastix.h>
153bf14a46SMatthew Knepley EXTERN_C_END
163bf14a46SMatthew Knepley 
17519f805aSKarl Rupp #if defined(PETSC_USE_COMPLEX)
18519f805aSKarl Rupp   #if defined(PETSC_USE_REAL_SINGLE)
195ec454cfSSatish Balay     #define PASTIX_CALL c_pastix
205ec454cfSSatish Balay   #else
215ec454cfSSatish Balay     #define PASTIX_CALL z_pastix
225ec454cfSSatish Balay   #endif
23d41469e0Sxavier lacoste 
24d41469e0Sxavier lacoste #else /* PETSC_USE_COMPLEX */
25d41469e0Sxavier lacoste 
26519f805aSKarl Rupp   #if defined(PETSC_USE_REAL_SINGLE)
275ec454cfSSatish Balay     #define PASTIX_CALL s_pastix
285ec454cfSSatish Balay   #else
295ec454cfSSatish Balay     #define PASTIX_CALL d_pastix
305ec454cfSSatish Balay   #endif
31d41469e0Sxavier lacoste 
32d41469e0Sxavier lacoste #endif /* PETSC_USE_COMPLEX */
33d41469e0Sxavier lacoste 
34dc4f8bb9SStefano Zampini #define PetscPastixCall(...) \
35dc4f8bb9SStefano Zampini   do { \
36dc4f8bb9SStefano Zampini     PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); \
37dc4f8bb9SStefano Zampini     PetscStackCallExternalVoid(PetscStringize(PASTIX_CALL), PASTIX_CALL(__VA_ARGS__)); \
38dc4f8bb9SStefano Zampini     PetscCall(PetscFPTrapPop()); \
39dc4f8bb9SStefano Zampini   } while (0)
40dc4f8bb9SStefano Zampini 
41dbbbd53dSSatish Balay typedef PetscScalar PastixScalar;
42dbbbd53dSSatish Balay 
433bf14a46SMatthew Knepley typedef struct Mat_Pastix_ {
443bf14a46SMatthew Knepley   pastix_data_t *pastix_data; /* Pastix data storage structure                        */
453bf14a46SMatthew Knepley   MatStructure   matstruc;
463bf14a46SMatthew Knepley   PetscInt       n;                 /* Number of columns in the matrix                      */
473bf14a46SMatthew Knepley   PetscInt      *colptr;            /* Index of first element of each column in row and val */
483bf14a46SMatthew Knepley   PetscInt      *row;               /* Row of each element of the matrix                    */
493bf14a46SMatthew Knepley   PetscScalar   *val;               /* Value of each element of the matrix                  */
503bf14a46SMatthew Knepley   PetscInt      *perm;              /* Permutation tabular                                  */
513bf14a46SMatthew Knepley   PetscInt      *invp;              /* Reverse permutation tabular                          */
523bf14a46SMatthew Knepley   PetscScalar   *rhs;               /* Rhight-hand-side member                              */
533bf14a46SMatthew Knepley   PetscInt       rhsnbr;            /* Rhight-hand-side number (must be 1)                  */
54baed9decSBarry Smith   PetscInt       iparm[IPARM_SIZE]; /* Integer parameters                                   */
55baed9decSBarry Smith   double         dparm[DPARM_SIZE]; /* Floating point parameters                            */
563bf14a46SMatthew Knepley   MPI_Comm       pastix_comm;       /* PaStiX MPI communicator                              */
573bf14a46SMatthew Knepley   PetscMPIInt    commRank;          /* MPI rank                                             */
583bf14a46SMatthew Knepley   PetscMPIInt    commSize;          /* MPI communicator size                                */
59ace3abfcSBarry Smith   PetscBool      CleanUpPastix;     /* Boolean indicating if we call PaStiX clean step      */
603bf14a46SMatthew Knepley   VecScatter     scat_rhs;
613bf14a46SMatthew Knepley   VecScatter     scat_sol;
62f31ce8a6SBarry Smith   Vec            b_seq;
633bf14a46SMatthew Knepley } Mat_Pastix;
643bf14a46SMatthew Knepley 
6509573ac7SBarry Smith extern PetscErrorCode MatDuplicate_Pastix(Mat, MatDuplicateOption, Mat *);
663bf14a46SMatthew Knepley 
673bf14a46SMatthew Knepley /*
683bf14a46SMatthew Knepley    convert Petsc seqaij matrix to CSC: colptr[n], row[nz], val[nz]
693bf14a46SMatthew Knepley 
703bf14a46SMatthew Knepley   input:
713bf14a46SMatthew Knepley     A       - matrix in seqaij or mpisbaij (bs=1) format
723bf14a46SMatthew Knepley     valOnly - FALSE: spaces are allocated and values are set for the CSC
733bf14a46SMatthew Knepley               TRUE:  Only fill values
743bf14a46SMatthew Knepley   output:
753bf14a46SMatthew Knepley     n       - Size of the matrix
763bf14a46SMatthew Knepley     colptr  - Index of first element of each column in row and val
773bf14a46SMatthew Knepley     row     - Row of each element of the matrix
783bf14a46SMatthew Knepley     values  - Value of each element of the matrix
793bf14a46SMatthew Knepley  */
80d71ae5a4SJacob Faibussowitsch PetscErrorCode MatConvertToCSC(Mat A, PetscBool valOnly, PetscInt *n, PetscInt **colptr, PetscInt **row, PetscScalar **values)
81d71ae5a4SJacob Faibussowitsch {
823bf14a46SMatthew Knepley   Mat_SeqAIJ  *aa      = (Mat_SeqAIJ *)A->data;
833bf14a46SMatthew Knepley   PetscInt    *rowptr  = aa->i;
843bf14a46SMatthew Knepley   PetscInt    *col     = aa->j;
853bf14a46SMatthew Knepley   PetscScalar *rvalues = aa->a;
863bf14a46SMatthew Knepley   PetscInt     m       = A->rmap->N;
87745c78f7SBarry Smith   PetscInt     nnz;
883bf14a46SMatthew Knepley   PetscInt     i, j, k;
893bf14a46SMatthew Knepley   PetscInt     base = 1;
903bf14a46SMatthew Knepley   PetscInt     idx;
913bf14a46SMatthew Knepley   PetscInt     colidx;
923bf14a46SMatthew Knepley   PetscInt    *colcount;
93ace3abfcSBarry Smith   PetscBool    isSBAIJ;
94ace3abfcSBarry Smith   PetscBool    isSeqSBAIJ;
95ace3abfcSBarry Smith   PetscBool    isMpiSBAIJ;
96ace3abfcSBarry Smith   PetscBool    isSym;
973bf14a46SMatthew Knepley 
983bf14a46SMatthew Knepley   PetscFunctionBegin;
999566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetric(A, 0.0, &isSym));
1009566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSBAIJ, &isSBAIJ));
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSeqSBAIJ));
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPISBAIJ, &isMpiSBAIJ));
1033bf14a46SMatthew Knepley 
104745c78f7SBarry Smith   *n = A->cmap->N;
105745c78f7SBarry Smith 
106745c78f7SBarry Smith   /* PaStiX only needs triangular matrix if matrix is symmetric
107745c78f7SBarry Smith    */
1082205254eSKarl Rupp   if (isSym && !(isSBAIJ || isSeqSBAIJ || isMpiSBAIJ)) nnz = (aa->nz - *n) / 2 + *n;
1092205254eSKarl Rupp   else nnz = aa->nz;
1103bf14a46SMatthew Knepley 
1113bf14a46SMatthew Knepley   if (!valOnly) {
1129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1((*n) + 1, colptr));
1139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nnz, row));
1149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nnz, values));
1153bf14a46SMatthew Knepley 
11641c8de11SBarry Smith     if (isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) {
1179566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(*colptr, rowptr, (*n) + 1));
1182205254eSKarl Rupp       for (i = 0; i < *n + 1; i++) (*colptr)[i] += base;
1199566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(*row, col, nnz));
1202205254eSKarl Rupp       for (i = 0; i < nnz; i++) (*row)[i] += base;
1219566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(*values, rvalues, nnz));
12241c8de11SBarry Smith     } else {
1239566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(*n, &colcount));
12441c8de11SBarry Smith 
125f31ce8a6SBarry Smith       for (i = 0; i < m; i++) colcount[i] = 0;
1263bf14a46SMatthew Knepley       /* Fill-in colptr */
127f31ce8a6SBarry Smith       for (i = 0; i < m; i++) {
128f31ce8a6SBarry Smith         for (j = rowptr[i]; j < rowptr[i + 1]; j++) {
129f31ce8a6SBarry Smith           if (!isSym || col[j] <= i) colcount[col[j]]++;
130f31ce8a6SBarry Smith         }
131f31ce8a6SBarry Smith       }
132745c78f7SBarry Smith 
1333bf14a46SMatthew Knepley       (*colptr)[0] = base;
1343bf14a46SMatthew Knepley       for (j = 0; j < *n; j++) {
1353bf14a46SMatthew Knepley         (*colptr)[j + 1] = (*colptr)[j] + colcount[j];
136745c78f7SBarry Smith         /* in next loop we fill starting from (*colptr)[colidx] - base */
1373bf14a46SMatthew Knepley         colcount[j] = -base;
1383bf14a46SMatthew Knepley       }
1393bf14a46SMatthew Knepley 
1403bf14a46SMatthew Knepley       /* Fill-in rows and values */
1413bf14a46SMatthew Knepley       for (i = 0; i < m; i++) {
1423bf14a46SMatthew Knepley         for (j = rowptr[i]; j < rowptr[i + 1]; j++) {
14341c8de11SBarry Smith           if (!isSym || col[j] <= i) {
1443bf14a46SMatthew Knepley             colidx         = col[j];
1453bf14a46SMatthew Knepley             idx            = (*colptr)[colidx] + colcount[colidx];
1463bf14a46SMatthew Knepley             (*row)[idx]    = i + base;
1473bf14a46SMatthew Knepley             (*values)[idx] = rvalues[j];
1483bf14a46SMatthew Knepley             colcount[colidx]++;
1493bf14a46SMatthew Knepley           }
1503bf14a46SMatthew Knepley         }
1513bf14a46SMatthew Knepley       }
1529566063dSJacob Faibussowitsch       PetscCall(PetscFree(colcount));
153745c78f7SBarry Smith     }
15441c8de11SBarry Smith   } else {
155745c78f7SBarry Smith     /* Fill-in only values */
1563bf14a46SMatthew Knepley     for (i = 0; i < m; i++) {
1573bf14a46SMatthew Knepley       for (j = rowptr[i]; j < rowptr[i + 1]; j++) {
1583bf14a46SMatthew Knepley         colidx = col[j];
1592205254eSKarl Rupp         if ((isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) || !isSym || col[j] <= i) {
160745c78f7SBarry Smith           /* look for the value to fill */
161f31ce8a6SBarry Smith           for (k = (*colptr)[colidx] - base; k < (*colptr)[colidx + 1] - base; k++) {
162eb1f6c34SBarry Smith             if (((*row)[k] - base) == i) {
1633bf14a46SMatthew Knepley               (*values)[k] = rvalues[j];
1643bf14a46SMatthew Knepley               break;
1653bf14a46SMatthew Knepley             }
1663bf14a46SMatthew Knepley           }
167f31ce8a6SBarry Smith           /* data structure of sparse matrix has changed */
168aed4548fSBarry Smith           PetscCheck(k != (*colptr)[colidx + 1] - base, PETSC_COMM_SELF, PETSC_ERR_PLIB, "overflow on k %" PetscInt_FMT, k);
1693bf14a46SMatthew Knepley         }
1703bf14a46SMatthew Knepley       }
1713bf14a46SMatthew Knepley     }
172745c78f7SBarry Smith   }
1733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1743bf14a46SMatthew Knepley }
1753bf14a46SMatthew Knepley 
1763bf14a46SMatthew Knepley /*
1773bf14a46SMatthew Knepley   Call clean step of PaStiX if lu->CleanUpPastix == true.
1783bf14a46SMatthew Knepley   Free the CSC matrix.
1793bf14a46SMatthew Knepley  */
180d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDestroy_Pastix(Mat A)
181d71ae5a4SJacob Faibussowitsch {
18258c95f1bSBarry Smith   Mat_Pastix *lu = (Mat_Pastix *)A->data;
183745c78f7SBarry Smith 
1843bf14a46SMatthew Knepley   PetscFunctionBegin;
18558c95f1bSBarry Smith   if (lu->CleanUpPastix) {
1863bf14a46SMatthew Knepley     /* Terminate instance, deallocate memories */
1879566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&lu->scat_rhs));
1889566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lu->b_seq));
1899566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&lu->scat_sol));
1903bf14a46SMatthew Knepley 
1913bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK] = API_TASK_CLEAN;
1923bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]   = API_TASK_CLEAN;
1933bf14a46SMatthew Knepley 
194dc4f8bb9SStefano Zampini     PetscPastixCall(&(lu->pastix_data), lu->pastix_comm, lu->n, lu->colptr, lu->row, (PastixScalar *)lu->val, lu->perm, lu->invp, (PastixScalar *)lu->rhs, lu->rhsnbr, lu->iparm, lu->dparm);
19508401ef6SPierre Jolivet     PetscCheck(lu->iparm[IPARM_ERROR_NUMBER] == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by PaStiX in destroy: iparm(IPARM_ERROR_NUMBER)=%" PetscInt_FMT, lu->iparm[IPARM_ERROR_NUMBER]);
1969566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->colptr));
1979566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->row));
1989566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->val));
1999566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->perm));
2009566063dSJacob Faibussowitsch     PetscCall(PetscFree(lu->invp));
2019566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_free(&(lu->pastix_comm)));
2023bf14a46SMatthew Knepley   }
2032e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatFactorGetSolverType_C", NULL));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2063bf14a46SMatthew Knepley }
2073bf14a46SMatthew Knepley 
2083bf14a46SMatthew Knepley /*
2093bf14a46SMatthew Knepley   Gather right-hand-side.
2103bf14a46SMatthew Knepley   Call for Solve step.
2113bf14a46SMatthew Knepley   Scatter solution.
2123bf14a46SMatthew Knepley  */
213d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_PaStiX(Mat A, Vec b, Vec x)
214d71ae5a4SJacob Faibussowitsch {
21558c95f1bSBarry Smith   Mat_Pastix  *lu = (Mat_Pastix *)A->data;
2163bf14a46SMatthew Knepley   PetscScalar *array;
2173bf14a46SMatthew Knepley   Vec          x_seq;
2183bf14a46SMatthew Knepley 
2193bf14a46SMatthew Knepley   PetscFunctionBegin;
2203bf14a46SMatthew Knepley   lu->rhsnbr = 1;
2213bf14a46SMatthew Knepley   x_seq      = lu->b_seq;
2223bf14a46SMatthew Knepley   if (lu->commSize > 1) {
22341ffd417SStefano Zampini     /* PaStiX only supports centralized rhs. Scatter b into a sequential rhs vector */
2249566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(lu->scat_rhs, b, x_seq, INSERT_VALUES, SCATTER_FORWARD));
2259566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(lu->scat_rhs, b, x_seq, INSERT_VALUES, SCATTER_FORWARD));
2269566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x_seq, &array));
22741c8de11SBarry Smith   } else { /* size == 1 */
2289566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, x));
2299566063dSJacob Faibussowitsch     PetscCall(VecGetArray(x, &array));
2303bf14a46SMatthew Knepley   }
2313bf14a46SMatthew Knepley   lu->rhs = array;
2323bf14a46SMatthew Knepley   if (lu->commSize == 1) {
2339566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &array));
2343bf14a46SMatthew Knepley   } else {
2359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x_seq, &array));
2363bf14a46SMatthew Knepley   }
2373bf14a46SMatthew Knepley 
2383bf14a46SMatthew Knepley   /* solve phase */
2393bf14a46SMatthew Knepley   lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE;
2403bf14a46SMatthew Knepley   lu->iparm[IPARM_END_TASK]   = API_TASK_REFINE;
241745c78f7SBarry Smith   lu->iparm[IPARM_RHS_MAKING] = API_RHS_B;
2423bf14a46SMatthew Knepley 
243dc4f8bb9SStefano Zampini   PetscPastixCall(&(lu->pastix_data), lu->pastix_comm, lu->n, lu->colptr, lu->row, (PastixScalar *)lu->val, lu->perm, lu->invp, (PastixScalar *)lu->rhs, lu->rhsnbr, lu->iparm, lu->dparm);
24408401ef6SPierre Jolivet   PetscCheck(lu->iparm[IPARM_ERROR_NUMBER] == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by PaStiX in solve phase: lu->iparm[IPARM_ERROR_NUMBER] = %" PetscInt_FMT, lu->iparm[IPARM_ERROR_NUMBER]);
2453bf14a46SMatthew Knepley 
2463bf14a46SMatthew Knepley   if (lu->commSize == 1) {
2479566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x, &(lu->rhs)));
2483bf14a46SMatthew Knepley   } else {
2499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(x_seq, &(lu->rhs)));
2503bf14a46SMatthew Knepley   }
2513bf14a46SMatthew Knepley 
2523bf14a46SMatthew Knepley   if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */
2539566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(lu->scat_sol, x_seq, x, INSERT_VALUES, SCATTER_FORWARD));
2549566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(lu->scat_sol, x_seq, x, INSERT_VALUES, SCATTER_FORWARD));
2553bf14a46SMatthew Knepley   }
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2573bf14a46SMatthew Knepley }
2583bf14a46SMatthew Knepley 
2593bf14a46SMatthew Knepley /*
2603bf14a46SMatthew Knepley   Numeric factorisation using PaStiX solver.
2613bf14a46SMatthew Knepley 
2623bf14a46SMatthew Knepley  */
263d71ae5a4SJacob Faibussowitsch PetscErrorCode MatFactorNumeric_PaStiX(Mat F, Mat A, const MatFactorInfo *info)
264d71ae5a4SJacob Faibussowitsch {
26558c95f1bSBarry Smith   Mat_Pastix *lu = (Mat_Pastix *)(F)->data;
26641c8de11SBarry Smith   Mat        *tseq;
267b5e56a35SBarry Smith   PetscInt    icntl;
268b5e56a35SBarry Smith   PetscInt    M = A->rmap->N;
269ace3abfcSBarry Smith   PetscBool   valOnly, flg, isSym;
2703bf14a46SMatthew Knepley   IS          is_iden;
2713bf14a46SMatthew Knepley   Vec         b;
2723bf14a46SMatthew Knepley   IS          isrow;
27351a30905SBarry Smith   PetscBool   isSeqAIJ, isSeqSBAIJ, isMPIAIJ;
2743bf14a46SMatthew Knepley 
2753bf14a46SMatthew Knepley   PetscFunctionBegin;
2769566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQAIJ, &isSeqAIJ));
2779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATMPIAIJ, &isMPIAIJ));
2789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATSEQSBAIJ, &isSeqSBAIJ));
2793bf14a46SMatthew Knepley   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN) {
2803bf14a46SMatthew Knepley     (F)->ops->solve = MatSolve_PaStiX;
2813bf14a46SMatthew Knepley 
2823bf14a46SMatthew Knepley     /* Initialize a PASTIX instance */
2839566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)A), &(lu->pastix_comm)));
2849566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(lu->pastix_comm, &lu->commRank));
2859566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(lu->pastix_comm, &lu->commSize));
2863bf14a46SMatthew Knepley 
2873bf14a46SMatthew Knepley     /* Set pastix options */
2883bf14a46SMatthew Knepley     lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO;
2893bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK]       = API_TASK_INIT;
2903bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]         = API_TASK_INIT;
2912205254eSKarl Rupp 
2923bf14a46SMatthew Knepley     lu->rhsnbr = 1;
2933bf14a46SMatthew Knepley 
2943bf14a46SMatthew Knepley     /* Call to set default pastix options */
295dc4f8bb9SStefano Zampini     PetscPastixCall(&(lu->pastix_data), lu->pastix_comm, lu->n, lu->colptr, lu->row, (PastixScalar *)lu->val, lu->perm, lu->invp, (PastixScalar *)lu->rhs, lu->rhsnbr, lu->iparm, lu->dparm);
29608401ef6SPierre Jolivet     PetscCheck(lu->iparm[IPARM_ERROR_NUMBER] == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by PaStiX in MatFactorNumeric: iparm(IPARM_ERROR_NUMBER)=%" PetscInt_FMT, lu->iparm[IPARM_ERROR_NUMBER]);
2973bf14a46SMatthew Knepley 
29826cc229bSBarry Smith     PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "PaStiX Options", "Mat");
2993bf14a46SMatthew Knepley     icntl                    = -1;
30041c8de11SBarry Smith     lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT;
3019566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_pastix_verbose", "iparm[IPARM_VERBOSE] : level of printing (0 to 2)", "None", lu->iparm[IPARM_VERBOSE], &icntl, &flg));
302ad540459SPierre Jolivet     if ((flg && icntl >= 0) || PetscLogPrintInfo) lu->iparm[IPARM_VERBOSE] = icntl;
3033bf14a46SMatthew Knepley     icntl = -1;
3049566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-mat_pastix_threadnbr", "iparm[IPARM_THREAD_NBR] : Number of thread by MPI node", "None", lu->iparm[IPARM_THREAD_NBR], &icntl, &flg));
305ad540459SPierre Jolivet     if ((flg && icntl > 0)) lu->iparm[IPARM_THREAD_NBR] = icntl;
306d0609cedSBarry Smith     PetscOptionsEnd();
3073bf14a46SMatthew Knepley     valOnly = PETSC_FALSE;
30841c8de11SBarry Smith   } else {
3095d6241c8SBarry Smith     if (isSeqAIJ || isMPIAIJ) {
3109566063dSJacob Faibussowitsch       PetscCall(PetscFree(lu->colptr));
3119566063dSJacob Faibussowitsch       PetscCall(PetscFree(lu->row));
3129566063dSJacob Faibussowitsch       PetscCall(PetscFree(lu->val));
3135d6241c8SBarry Smith       valOnly = PETSC_FALSE;
3145d6241c8SBarry Smith     } else valOnly = PETSC_TRUE;
3153bf14a46SMatthew Knepley   }
3163bf14a46SMatthew Knepley 
3173bf14a46SMatthew Knepley   lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES;
3183bf14a46SMatthew Knepley 
3193bf14a46SMatthew Knepley   /* convert mpi A to seq mat A */
3209566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &isrow));
3219566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A, 1, &isrow, &isrow, MAT_INITIAL_MATRIX, &tseq));
3229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
3233bf14a46SMatthew Knepley 
3249566063dSJacob Faibussowitsch   PetscCall(MatConvertToCSC(*tseq, valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val));
3259566063dSJacob Faibussowitsch   PetscCall(MatIsSymmetric(*tseq, 0.0, &isSym));
3269566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &tseq));
32741c8de11SBarry Smith 
3285d6241c8SBarry Smith   if (!lu->perm) {
3299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lu->n, &(lu->perm)));
3309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(lu->n, &(lu->invp)));
3315d6241c8SBarry Smith   }
3323bf14a46SMatthew Knepley 
3333bf14a46SMatthew Knepley   if (isSym) {
334745c78f7SBarry Smith     /* On symmetric matrix, LLT */
3353bf14a46SMatthew Knepley     lu->iparm[IPARM_SYM]           = API_SYM_YES;
33641c8de11SBarry Smith     lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT;
337f31ce8a6SBarry Smith   } else {
338745c78f7SBarry Smith     /* On unsymmetric matrix, LU */
3393bf14a46SMatthew Knepley     lu->iparm[IPARM_SYM]           = API_SYM_NO;
3403bf14a46SMatthew Knepley     lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
3413bf14a46SMatthew Knepley   }
3423bf14a46SMatthew Knepley 
3433bf14a46SMatthew Knepley   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN) {
34458c95f1bSBarry Smith     if (!(isSeqAIJ || isSeqSBAIJ) && !lu->b_seq) {
3453bf14a46SMatthew Knepley       /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
3469566063dSJacob Faibussowitsch       PetscCall(VecCreateSeq(PETSC_COMM_SELF, A->cmap->N, &lu->b_seq));
3479566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, A->cmap->N, 0, 1, &is_iden));
3489566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(A, NULL, &b));
3499566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(b, is_iden, lu->b_seq, is_iden, &lu->scat_rhs));
3509566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(lu->b_seq, is_iden, b, is_iden, &lu->scat_sol));
3519566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is_iden));
3529566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&b));
3533bf14a46SMatthew Knepley     }
3543bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING;
3553bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
3563bf14a46SMatthew Knepley 
357dc4f8bb9SStefano Zampini     PetscPastixCall(&(lu->pastix_data), lu->pastix_comm, lu->n, lu->colptr, lu->row, (PastixScalar *)lu->val, lu->perm, lu->invp, (PastixScalar *)lu->rhs, lu->rhsnbr, lu->iparm, lu->dparm);
35808401ef6SPierre Jolivet     PetscCheck(lu->iparm[IPARM_ERROR_NUMBER] == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%" PetscInt_FMT, lu->iparm[IPARM_ERROR_NUMBER]);
35941c8de11SBarry Smith   } else {
3603bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT;
3613bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
362dc4f8bb9SStefano Zampini     PetscPastixCall(&(lu->pastix_data), lu->pastix_comm, lu->n, lu->colptr, lu->row, (PastixScalar *)lu->val, lu->perm, lu->invp, (PastixScalar *)lu->rhs, lu->rhsnbr, lu->iparm, lu->dparm);
36308401ef6SPierre Jolivet     PetscCheck(lu->iparm[IPARM_ERROR_NUMBER] == 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%" PetscInt_FMT, lu->iparm[IPARM_ERROR_NUMBER]);
3643bf14a46SMatthew Knepley   }
3653bf14a46SMatthew Knepley 
3663bf14a46SMatthew Knepley   (F)->assembled    = PETSC_TRUE;
3673bf14a46SMatthew Knepley   lu->matstruc      = SAME_NONZERO_PATTERN;
3683bf14a46SMatthew Knepley   lu->CleanUpPastix = PETSC_TRUE;
3693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3703bf14a46SMatthew Knepley }
3713bf14a46SMatthew Knepley 
3723bf14a46SMatthew Knepley /* Note the Petsc r and c permutations are ignored */
373d71ae5a4SJacob Faibussowitsch PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F, Mat A, IS r, IS c, const MatFactorInfo *info)
374d71ae5a4SJacob Faibussowitsch {
37558c95f1bSBarry Smith   Mat_Pastix *lu = (Mat_Pastix *)F->data;
3763bf14a46SMatthew Knepley 
3773bf14a46SMatthew Knepley   PetscFunctionBegin;
3783bf14a46SMatthew Knepley   lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
3793bf14a46SMatthew Knepley   lu->iparm[IPARM_SYM]           = API_SYM_YES;
3803bf14a46SMatthew Knepley   lu->matstruc                   = DIFFERENT_NONZERO_PATTERN;
3813bf14a46SMatthew Knepley   F->ops->lufactornumeric        = MatFactorNumeric_PaStiX;
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3833bf14a46SMatthew Knepley }
3843bf14a46SMatthew Knepley 
385d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F, Mat A, IS r, const MatFactorInfo *info)
386d71ae5a4SJacob Faibussowitsch {
38758c95f1bSBarry Smith   Mat_Pastix *lu = (Mat_Pastix *)(F)->data;
3883bf14a46SMatthew Knepley 
3893bf14a46SMatthew Knepley   PetscFunctionBegin;
3903bf14a46SMatthew Knepley   lu->iparm[IPARM_FACTORIZATION]  = API_FACT_LLT;
3913bf14a46SMatthew Knepley   lu->iparm[IPARM_SYM]            = API_SYM_NO;
3923bf14a46SMatthew Knepley   lu->matstruc                    = DIFFERENT_NONZERO_PATTERN;
3933bf14a46SMatthew Knepley   (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX;
3943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3953bf14a46SMatthew Knepley }
3963bf14a46SMatthew Knepley 
397d71ae5a4SJacob Faibussowitsch PetscErrorCode MatView_PaStiX(Mat A, PetscViewer viewer)
398d71ae5a4SJacob Faibussowitsch {
399ace3abfcSBarry Smith   PetscBool         iascii;
4003bf14a46SMatthew Knepley   PetscViewerFormat format;
4013bf14a46SMatthew Knepley 
4023bf14a46SMatthew Knepley   PetscFunctionBegin;
4039566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
4043bf14a46SMatthew Knepley   if (iascii) {
4059566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
4063bf14a46SMatthew Knepley     if (format == PETSC_VIEWER_ASCII_INFO) {
40758c95f1bSBarry Smith       Mat_Pastix *lu = (Mat_Pastix *)A->data;
408b5e56a35SBarry Smith 
4099566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "PaStiX run parameters:\n"));
4109566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Matrix type :                      %s \n", ((lu->iparm[IPARM_SYM] == API_SYM_YES) ? "Symmetric" : "Unsymmetric")));
4119566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Level of printing (0,1,2):         %" PetscInt_FMT " \n", lu->iparm[IPARM_VERBOSE]));
4129566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of refinements iterations : %" PetscInt_FMT " \n", lu->iparm[IPARM_NBITER]));
4139566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Error :                        %g \n", lu->dparm[DPARM_RELATIVE_ERROR]));
4143bf14a46SMatthew Knepley     }
4153bf14a46SMatthew Knepley   }
4163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4173bf14a46SMatthew Knepley }
4183bf14a46SMatthew Knepley 
4193bf14a46SMatthew Knepley /*MC
4202692d6eeSBarry Smith      MATSOLVERPASTIX  - A solver package providing direct solvers (LU) for distributed
4213bf14a46SMatthew Knepley   and sequential matrices via the external package PaStiX.
4223bf14a46SMatthew Knepley 
4232ef1f0ffSBarry Smith   Use `./configure` `--download-pastix` `--download-ptscotch`  to have PETSc installed with PasTiX
424c2b89b5dSBarry Smith 
4252ef1f0ffSBarry Smith   Use `-pc_type lu` `-pc_factor_mat_solver_type pastix` to use this direct solver
4263bf14a46SMatthew Knepley 
4273bf14a46SMatthew Knepley   Options Database Keys:
4282ef1f0ffSBarry Smith + -mat_pastix_verbose   <0,1,2>   - print level of information messages from PaStiX
4292ef1f0ffSBarry Smith - -mat_pastix_threadnbr <integer> - Set the number of threads for each MPI process
4303bf14a46SMatthew Knepley 
43195452b02SPatrick Sanan   Notes:
43295452b02SPatrick Sanan     This only works for matrices with symmetric nonzero structure, if you pass it a matrix with
4332ef1f0ffSBarry Smith    nonsymmetric structure PasTiX, and hence, PETSc return with an error.
434baed9decSBarry Smith 
4353bf14a46SMatthew Knepley   Level: beginner
4363bf14a46SMatthew Knepley 
437*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCFactorSetMatSolverType()`, `MatSolverType`, `MatGetFactor()`
4383bf14a46SMatthew Knepley M*/
4393bf14a46SMatthew Knepley 
440d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_PaStiX(Mat A, MatInfoType flag, MatInfo *info)
441d71ae5a4SJacob Faibussowitsch {
44258c95f1bSBarry Smith   Mat_Pastix *lu = (Mat_Pastix *)A->data;
4433bf14a46SMatthew Knepley 
4443bf14a46SMatthew Knepley   PetscFunctionBegin;
4453bf14a46SMatthew Knepley   info->block_size        = 1.0;
4463bf14a46SMatthew Knepley   info->nz_allocated      = lu->iparm[IPARM_NNZEROS];
4473bf14a46SMatthew Knepley   info->nz_used           = lu->iparm[IPARM_NNZEROS];
4483bf14a46SMatthew Knepley   info->nz_unneeded       = 0.0;
4493bf14a46SMatthew Knepley   info->assemblies        = 0.0;
4503bf14a46SMatthew Knepley   info->mallocs           = 0.0;
4513bf14a46SMatthew Knepley   info->memory            = 0.0;
4523bf14a46SMatthew Knepley   info->fill_ratio_given  = 0;
4533bf14a46SMatthew Knepley   info->fill_ratio_needed = 0;
4543bf14a46SMatthew Knepley   info->factor_mallocs    = 0;
4553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4563bf14a46SMatthew Knepley }
4573bf14a46SMatthew Knepley 
458d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_pastix(Mat A, MatSolverType *type)
459d71ae5a4SJacob Faibussowitsch {
4603bf14a46SMatthew Knepley   PetscFunctionBegin;
4612692d6eeSBarry Smith   *type = MATSOLVERPASTIX;
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4633bf14a46SMatthew Knepley }
4643bf14a46SMatthew Knepley 
4653bf14a46SMatthew Knepley /*
4663bf14a46SMatthew Knepley     The seq and mpi versions of this function are the same
4673bf14a46SMatthew Knepley */
468d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_seqaij_pastix(Mat A, MatFactorType ftype, Mat *F)
469d71ae5a4SJacob Faibussowitsch {
4703bf14a46SMatthew Knepley   Mat         B;
4713bf14a46SMatthew Knepley   Mat_Pastix *pastix;
4723bf14a46SMatthew Knepley 
4733bf14a46SMatthew Knepley   PetscFunctionBegin;
4745f80ce2aSJacob Faibussowitsch   PetscCheck(ftype == MAT_FACTOR_LU, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
4753bf14a46SMatthew Knepley   /* Create the factorization matrix */
4769566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
4779566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
4789566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("pastix", &((PetscObject)B)->type_name));
4799566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
4803bf14a46SMatthew Knepley 
48166e17bc3SBarry Smith   B->trivialsymbolic       = PETSC_TRUE;
4823bf14a46SMatthew Knepley   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
4833bf14a46SMatthew Knepley   B->ops->view             = MatView_PaStiX;
4843bf14a46SMatthew Knepley   B->ops->getinfo          = MatGetInfo_PaStiX;
4852205254eSKarl Rupp 
4869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_pastix));
4872205254eSKarl Rupp 
488d5f3da31SBarry Smith   B->factortype = MAT_FACTOR_LU;
4893bf14a46SMatthew Knepley 
49000c67f3bSHong Zhang   /* set solvertype */
4919566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
4929566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPASTIX, &B->solvertype));
49300c67f3bSHong Zhang 
4944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pastix));
4952205254eSKarl Rupp 
4963bf14a46SMatthew Knepley   pastix->CleanUpPastix = PETSC_FALSE;
4970298fd71SBarry Smith   pastix->scat_rhs      = NULL;
4980298fd71SBarry Smith   pastix->scat_sol      = NULL;
49958c95f1bSBarry Smith   B->ops->getinfo       = MatGetInfo_External;
5003bf14a46SMatthew Knepley   B->ops->destroy       = MatDestroy_Pastix;
50158c95f1bSBarry Smith   B->data               = (void *)pastix;
5023bf14a46SMatthew Knepley 
5033bf14a46SMatthew Knepley   *F = B;
5043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5053bf14a46SMatthew Knepley }
5063bf14a46SMatthew Knepley 
507d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A, MatFactorType ftype, Mat *F)
508d71ae5a4SJacob Faibussowitsch {
5093bf14a46SMatthew Knepley   Mat         B;
5103bf14a46SMatthew Knepley   Mat_Pastix *pastix;
5113bf14a46SMatthew Knepley 
5123bf14a46SMatthew Knepley   PetscFunctionBegin;
5135f80ce2aSJacob Faibussowitsch   PetscCheck(ftype == MAT_FACTOR_LU, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
5143bf14a46SMatthew Knepley   /* Create the factorization matrix */
5159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
5179566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("pastix", &((PetscObject)B)->type_name));
5189566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
5193bf14a46SMatthew Knepley 
52066e17bc3SBarry Smith   B->trivialsymbolic       = PETSC_TRUE;
5213bf14a46SMatthew Knepley   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
5223bf14a46SMatthew Knepley   B->ops->view             = MatView_PaStiX;
52358c95f1bSBarry Smith   B->ops->getinfo          = MatGetInfo_PaStiX;
5249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_pastix));
5252205254eSKarl Rupp 
526d5f3da31SBarry Smith   B->factortype = MAT_FACTOR_LU;
5273bf14a46SMatthew Knepley 
52800c67f3bSHong Zhang   /* set solvertype */
5299566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
5309566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPASTIX, &B->solvertype));
53100c67f3bSHong Zhang 
5324dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pastix));
5332205254eSKarl Rupp 
5343bf14a46SMatthew Knepley   pastix->CleanUpPastix = PETSC_FALSE;
5350298fd71SBarry Smith   pastix->scat_rhs      = NULL;
5360298fd71SBarry Smith   pastix->scat_sol      = NULL;
53758c95f1bSBarry Smith   B->ops->getinfo       = MatGetInfo_External;
5383bf14a46SMatthew Knepley   B->ops->destroy       = MatDestroy_Pastix;
53958c95f1bSBarry Smith   B->data               = (void *)pastix;
5403bf14a46SMatthew Knepley 
5413bf14a46SMatthew Knepley   *F = B;
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5433bf14a46SMatthew Knepley }
5443bf14a46SMatthew Knepley 
545d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A, MatFactorType ftype, Mat *F)
546d71ae5a4SJacob Faibussowitsch {
5473bf14a46SMatthew Knepley   Mat         B;
5483bf14a46SMatthew Knepley   Mat_Pastix *pastix;
5493bf14a46SMatthew Knepley 
5503bf14a46SMatthew Knepley   PetscFunctionBegin;
5515f80ce2aSJacob Faibussowitsch   PetscCheck(ftype == MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
5523bf14a46SMatthew Knepley   /* Create the factorization matrix */
5539566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5549566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
5559566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("pastix", &((PetscObject)B)->type_name));
5569566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
5573bf14a46SMatthew Knepley 
55866e17bc3SBarry Smith   B->trivialsymbolic             = PETSC_TRUE;
5593bf14a46SMatthew Knepley   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
5603bf14a46SMatthew Knepley   B->ops->view                   = MatView_PaStiX;
56158c95f1bSBarry Smith   B->ops->getinfo                = MatGetInfo_PaStiX;
5629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_pastix));
5632205254eSKarl Rupp 
564d5f3da31SBarry Smith   B->factortype = MAT_FACTOR_CHOLESKY;
5653bf14a46SMatthew Knepley 
56600c67f3bSHong Zhang   /* set solvertype */
5679566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
5689566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPASTIX, &B->solvertype));
56900c67f3bSHong Zhang 
5704dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pastix));
5712205254eSKarl Rupp 
5723bf14a46SMatthew Knepley   pastix->CleanUpPastix = PETSC_FALSE;
5730298fd71SBarry Smith   pastix->scat_rhs      = NULL;
5740298fd71SBarry Smith   pastix->scat_sol      = NULL;
57558c95f1bSBarry Smith   B->ops->getinfo       = MatGetInfo_External;
5763bf14a46SMatthew Knepley   B->ops->destroy       = MatDestroy_Pastix;
57758c95f1bSBarry Smith   B->data               = (void *)pastix;
5783bf14a46SMatthew Knepley   *F                    = B;
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5803bf14a46SMatthew Knepley }
5813bf14a46SMatthew Knepley 
582d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A, MatFactorType ftype, Mat *F)
583d71ae5a4SJacob Faibussowitsch {
5843bf14a46SMatthew Knepley   Mat         B;
5853bf14a46SMatthew Knepley   Mat_Pastix *pastix;
5863bf14a46SMatthew Knepley 
5873bf14a46SMatthew Knepley   PetscFunctionBegin;
5885f80ce2aSJacob Faibussowitsch   PetscCheck(ftype == MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
58941c8de11SBarry Smith 
5903bf14a46SMatthew Knepley   /* Create the factorization matrix */
5919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
5939566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("pastix", &((PetscObject)B)->type_name));
5949566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
5953bf14a46SMatthew Knepley 
5963bf14a46SMatthew Knepley   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
5973bf14a46SMatthew Knepley   B->ops->view                   = MatView_PaStiX;
59858c95f1bSBarry Smith   B->ops->getinfo                = MatGetInfo_PaStiX;
59958c95f1bSBarry Smith   B->ops->destroy                = MatDestroy_Pastix;
6009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_pastix));
6012205254eSKarl Rupp 
602d5f3da31SBarry Smith   B->factortype = MAT_FACTOR_CHOLESKY;
6033bf14a46SMatthew Knepley 
60400c67f3bSHong Zhang   /* set solvertype */
6059566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
6069566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPASTIX, &B->solvertype));
60700c67f3bSHong Zhang 
6084dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&pastix));
6092205254eSKarl Rupp 
6103bf14a46SMatthew Knepley   pastix->CleanUpPastix = PETSC_FALSE;
6110298fd71SBarry Smith   pastix->scat_rhs      = NULL;
6120298fd71SBarry Smith   pastix->scat_sol      = NULL;
61358c95f1bSBarry Smith   B->data               = (void *)pastix;
6143bf14a46SMatthew Knepley 
6153bf14a46SMatthew Knepley   *F = B;
6163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6173bf14a46SMatthew Knepley }
618f7a08781SBarry Smith 
619d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatSolverTypeRegister_Pastix(void)
620d71ae5a4SJacob Faibussowitsch {
62142c9c57cSBarry Smith   PetscFunctionBegin;
6229566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATMPIAIJ, MAT_FACTOR_LU, MatGetFactor_mpiaij_pastix));
6239566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATSEQAIJ, MAT_FACTOR_LU, MatGetFactor_seqaij_pastix));
6249566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATMPISBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_mpisbaij_pastix));
6259566063dSJacob Faibussowitsch   PetscCall(MatSolverTypeRegister(MATSOLVERPASTIX, MATSEQSBAIJ, MAT_FACTOR_CHOLESKY, MatGetFactor_seqsbaij_pastix));
6263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62742c9c57cSBarry Smith }
628