xref: /petsc/src/mat/impls/aij/mpi/pastix/pastix.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
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 
1143801a69SSatish Balay #if defined(PETSC_HAVE_STDLIB_H)
1243801a69SSatish Balay #include <stdlib.h>
1343801a69SSatish Balay #endif
1443801a69SSatish Balay #if defined(PETSC_HAVE_STRING_H)
1543801a69SSatish Balay #include <string.h>
1643801a69SSatish Balay #endif
1743801a69SSatish Balay 
183bf14a46SMatthew Knepley EXTERN_C_BEGIN
193bf14a46SMatthew Knepley #include "pastix.h"
203bf14a46SMatthew Knepley EXTERN_C_END
213bf14a46SMatthew Knepley 
223bf14a46SMatthew Knepley typedef struct Mat_Pastix_ {
233bf14a46SMatthew Knepley   pastix_data_t *pastix_data;              /* Pastix data storage structure                        */
243bf14a46SMatthew Knepley   MatStructure   matstruc;
253bf14a46SMatthew Knepley   PetscInt       n;                        /* Number of columns in the matrix                      */
263bf14a46SMatthew Knepley   PetscInt       *colptr;                  /* Index of first element of each column in row and val */
273bf14a46SMatthew Knepley   PetscInt       *row;                     /* Row of each element of the matrix                    */
283bf14a46SMatthew Knepley   PetscScalar    *val;                     /* Value of each element of the matrix                  */
293bf14a46SMatthew Knepley   PetscInt       *perm;                    /* Permutation tabular                                  */
303bf14a46SMatthew Knepley   PetscInt       *invp;                    /* Reverse permutation tabular                          */
313bf14a46SMatthew Knepley   PetscScalar    *rhs;                     /* Rhight-hand-side member                              */
323bf14a46SMatthew Knepley   PetscInt       rhsnbr;                   /* Rhight-hand-side number (must be 1)                  */
333bf14a46SMatthew Knepley   PetscInt       iparm[64];                /* Integer parameters                                   */
343bf14a46SMatthew Knepley   double         dparm[64];                /* Floating point parameters                            */
353bf14a46SMatthew Knepley   MPI_Comm       pastix_comm;              /* PaStiX MPI communicator                              */
363bf14a46SMatthew Knepley   PetscMPIInt    commRank;                 /* MPI rank                                             */
373bf14a46SMatthew Knepley   PetscMPIInt    commSize;                 /* MPI communicator size                                */
383bf14a46SMatthew Knepley   PetscTruth     CleanUpPastix;            /* Boolean indicating if we call PaStiX clean step      */
393bf14a46SMatthew Knepley   VecScatter     scat_rhs;
403bf14a46SMatthew Knepley   VecScatter     scat_sol;
41f31ce8a6SBarry Smith   Vec            b_seq;
423bf14a46SMatthew Knepley   PetscTruth     isAIJ;
433bf14a46SMatthew Knepley   PetscErrorCode (*MatDestroy)(Mat);
443bf14a46SMatthew Knepley } Mat_Pastix;
453bf14a46SMatthew Knepley 
463bf14a46SMatthew Knepley EXTERN PetscErrorCode MatDuplicate_Pastix(Mat,MatDuplicateOption,Mat*);
473bf14a46SMatthew Knepley 
48eb1f6c34SBarry Smith #undef __FUNCT__
49eb1f6c34SBarry Smith #define __FUNCT__ "MatConvertToCSC"
503bf14a46SMatthew Knepley /*
513bf14a46SMatthew Knepley    convert Petsc seqaij matrix to CSC: colptr[n], row[nz], val[nz]
523bf14a46SMatthew Knepley 
533bf14a46SMatthew Knepley   input:
543bf14a46SMatthew Knepley     A       - matrix in seqaij or mpisbaij (bs=1) format
553bf14a46SMatthew Knepley     valOnly - FALSE: spaces are allocated and values are set for the CSC
563bf14a46SMatthew Knepley               TRUE:  Only fill values
573bf14a46SMatthew Knepley   output:
583bf14a46SMatthew Knepley     n       - Size of the matrix
593bf14a46SMatthew Knepley     colptr  - Index of first element of each column in row and val
603bf14a46SMatthew Knepley     row     - Row of each element of the matrix
613bf14a46SMatthew Knepley     values  - Value of each element of the matrix
623bf14a46SMatthew Knepley  */
6341c8de11SBarry Smith PetscErrorCode MatConvertToCSC(Mat A,PetscTruth valOnly,PetscInt *n,PetscInt **colptr,PetscInt **row,PetscScalar **values)
6441c8de11SBarry Smith {
653bf14a46SMatthew Knepley   Mat_SeqAIJ     *aa      = (Mat_SeqAIJ*)A->data;
663bf14a46SMatthew Knepley   PetscInt       *rowptr  = aa->i;
673bf14a46SMatthew Knepley   PetscInt       *col     = aa->j;
683bf14a46SMatthew Knepley   PetscScalar    *rvalues = aa->a;
693bf14a46SMatthew Knepley   PetscInt        m       = A->rmap->N;
70745c78f7SBarry Smith   PetscInt        nnz;
713bf14a46SMatthew Knepley   PetscInt        i,j, k;
723bf14a46SMatthew Knepley   PetscInt        base = 1;
733bf14a46SMatthew Knepley   PetscInt        idx;
743bf14a46SMatthew Knepley   PetscErrorCode  ierr;
753bf14a46SMatthew Knepley   PetscInt        colidx;
763bf14a46SMatthew Knepley   PetscInt       *colcount;
7741c8de11SBarry Smith   PetscTruth      isSBAIJ;
7841c8de11SBarry Smith   PetscTruth      isSeqSBAIJ;
7941c8de11SBarry Smith   PetscTruth      isMpiSBAIJ;
80745c78f7SBarry Smith   PetscTruth      isSym;
813bf14a46SMatthew Knepley 
823bf14a46SMatthew Knepley   PetscFunctionBegin;
8341c8de11SBarry Smith   ierr = MatIsSymmetric(A,0.0,&isSym);CHKERRQ(ierr);
8441c8de11SBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATSBAIJ,&isSBAIJ);CHKERRQ(ierr);
8541c8de11SBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
8641c8de11SBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATMPISBAIJ,&isMpiSBAIJ);CHKERRQ(ierr);
873bf14a46SMatthew Knepley 
88745c78f7SBarry Smith   *n = A->cmap->N;
89745c78f7SBarry Smith 
90745c78f7SBarry Smith   /* PaStiX only needs triangular matrix if matrix is symmetric
91745c78f7SBarry Smith    */
9241c8de11SBarry Smith   if (isSym && !(isSBAIJ || isSeqSBAIJ || isMpiSBAIJ)) {
93745c78f7SBarry Smith     nnz = (aa->nz - *n)/2 + *n;
94745c78f7SBarry Smith   }
9541c8de11SBarry Smith   else {
96745c78f7SBarry Smith     nnz     = aa->nz;
973bf14a46SMatthew Knepley   }
983bf14a46SMatthew Knepley 
993bf14a46SMatthew Knepley   if (!valOnly){
1003bf14a46SMatthew Knepley     ierr = PetscMalloc(((*n)+1) *sizeof(PetscInt)   ,colptr );CHKERRQ(ierr);
1013bf14a46SMatthew Knepley     ierr = PetscMalloc( nnz     *sizeof(PetscInt)   ,row);CHKERRQ(ierr);
1023bf14a46SMatthew Knepley     ierr = PetscMalloc( nnz     *sizeof(PetscScalar),values);CHKERRQ(ierr);
1033bf14a46SMatthew Knepley 
10441c8de11SBarry Smith     if (isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) {
10541c8de11SBarry Smith 	ierr = PetscMemcpy (*colptr, rowptr, ((*n)+1)*sizeof(PetscInt));CHKERRQ(ierr);
10641c8de11SBarry Smith 	for (i = 0; i < *n+1; i++)
10741c8de11SBarry Smith 	  (*colptr)[i] += base;
10841c8de11SBarry Smith 	ierr = PetscMemcpy (*row, col, (nnz)*sizeof(PetscInt));CHKERRQ(ierr);
10941c8de11SBarry Smith 	for (i = 0; i < nnz; i++)
11041c8de11SBarry Smith 	  (*row)[i] += base;
11141c8de11SBarry Smith 	ierr = PetscMemcpy (*values, rvalues, (nnz)*sizeof(PetscScalar));CHKERRQ(ierr);
11241c8de11SBarry Smith     } else {
11341c8de11SBarry Smith       ierr = PetscMalloc((*n)*sizeof(PetscInt)   ,&colcount);CHKERRQ(ierr);
11441c8de11SBarry Smith 
115f31ce8a6SBarry Smith       for (i = 0; i < m; i++) colcount[i] = 0;
1163bf14a46SMatthew Knepley       /* Fill-in colptr */
117f31ce8a6SBarry Smith       for (i = 0; i < m; i++) {
118f31ce8a6SBarry Smith 	for (j = rowptr[i]; j < rowptr[i+1]; j++) {
119f31ce8a6SBarry Smith 	  if (!isSym || col[j] <= i)  colcount[col[j]]++;
120f31ce8a6SBarry Smith         }
121f31ce8a6SBarry Smith       }
122745c78f7SBarry Smith 
1233bf14a46SMatthew Knepley       (*colptr)[0] = base;
1243bf14a46SMatthew Knepley       for (j = 0; j < *n; j++) {
1253bf14a46SMatthew Knepley 	(*colptr)[j+1] = (*colptr)[j] + colcount[j];
126745c78f7SBarry Smith 	/* in next loop we fill starting from (*colptr)[colidx] - base */
1273bf14a46SMatthew Knepley 	colcount[j] = -base;
1283bf14a46SMatthew Knepley       }
1293bf14a46SMatthew Knepley 
1303bf14a46SMatthew Knepley       /* Fill-in rows and values */
1313bf14a46SMatthew Knepley       for (i = 0; i < m; i++) {
1323bf14a46SMatthew Knepley 	for (j = rowptr[i]; j < rowptr[i+1]; j++) {
13341c8de11SBarry Smith 	  if (!isSym || col[j] <= i) {
1343bf14a46SMatthew Knepley 	    colidx = col[j];
1353bf14a46SMatthew Knepley 	    idx    = (*colptr)[colidx] + colcount[colidx];
1363bf14a46SMatthew Knepley 	    (*row)[idx]    = i + base;
1373bf14a46SMatthew Knepley 	    (*values)[idx] = rvalues[j];
1383bf14a46SMatthew Knepley 	    colcount[colidx]++;
1393bf14a46SMatthew Knepley 	  }
1403bf14a46SMatthew Knepley 	}
1413bf14a46SMatthew Knepley       }
14241c8de11SBarry Smith       ierr = PetscFree(colcount);CHKERRQ(ierr);
143745c78f7SBarry Smith     }
14441c8de11SBarry Smith   } else {
145745c78f7SBarry Smith     /* Fill-in only values */
1463bf14a46SMatthew Knepley     for (i = 0; i < m; i++) {
1473bf14a46SMatthew Knepley       for (j = rowptr[i]; j < rowptr[i+1]; j++) {
1483bf14a46SMatthew Knepley 	colidx = col[j];
14941c8de11SBarry Smith 	if ((isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) ||!isSym || col[j] <= i)
150745c78f7SBarry Smith 	  {
151745c78f7SBarry Smith 	    /* look for the value to fill */
152f31ce8a6SBarry Smith 	    for (k = (*colptr)[colidx] - base; k < (*colptr)[colidx + 1] - base; k++) {
153eb1f6c34SBarry Smith 	      if (((*row)[k]-base) == i) {
1543bf14a46SMatthew Knepley 		(*values)[k] = rvalues[j];
1553bf14a46SMatthew Knepley 		break;
1563bf14a46SMatthew Knepley 	      }
1573bf14a46SMatthew Knepley 	    }
158f31ce8a6SBarry Smith 	    /* data structure of sparse matrix has changed */
159*e32f2f54SBarry Smith 	    if (k == (*colptr)[colidx + 1] - base) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"overflow on k %D",k);
1603bf14a46SMatthew Knepley 	  }
1613bf14a46SMatthew Knepley       }
1623bf14a46SMatthew Knepley     }
163745c78f7SBarry Smith   }
16443801a69SSatish Balay   {
16570fe17b1SSatish Balay     PetscScalar *tmpvalues;
16670fe17b1SSatish Balay     PetscInt    *tmprows,*tmpcolptr;
16770fe17b1SSatish Balay     ierr = PetscMalloc3(nnz,PetscScalar,&tmpvalues,nnz,PetscInt,&tmprows,(*n+1),PetscInt,&tmpcolptr);CHKERRQ(ierr);
16870fe17b1SSatish Balay     if (sizeof(PetscScalar) != sizeof(pastix_float_t)) {
169*e32f2f54SBarry Smith       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"sizeof(PetscScalar) %d != sizeof(pastix_float_t) %d",sizeof(PetscScalar),sizeof(pastix_float_t));
17043801a69SSatish Balay     }
17170fe17b1SSatish Balay     if (sizeof(PetscInt) != sizeof(pastix_int_t)) {
172*e32f2f54SBarry Smith       SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"sizeof(PetscInt) %d != sizeof(pastix_int_t) %d",sizeof(PetscInt),sizeof(pastix_int_t));
17370fe17b1SSatish Balay     }
17443801a69SSatish Balay 
17570fe17b1SSatish Balay     ierr = PetscMemcpy(tmpcolptr,*colptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr);
17670fe17b1SSatish Balay     ierr = PetscMemcpy(tmprows,*row,nnz*sizeof(PetscInt));CHKERRQ(ierr);
17770fe17b1SSatish Balay     ierr = PetscMemcpy(tmpvalues,*values,nnz*sizeof(PetscScalar));CHKERRQ(ierr);
17843801a69SSatish Balay     ierr = PetscFree(*row);CHKERRQ(ierr);
17943801a69SSatish Balay     ierr = PetscFree(*values);CHKERRQ(ierr);
18043801a69SSatish Balay 
181f31ce8a6SBarry Smith     pastix_checkMatrix(MPI_COMM_WORLD,API_VERBOSE_NO,((isSym != 0) ? API_SYM_YES : API_SYM_NO),API_YES,*n,&tmpcolptr,&tmprows,&tmpvalues,NULL,1);
18243801a69SSatish Balay 
18370fe17b1SSatish Balay     ierr = PetscMemcpy(*colptr,tmpcolptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr);
18443801a69SSatish Balay     ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscInt),row);CHKERRQ(ierr);
18570fe17b1SSatish Balay     ierr = PetscMemcpy(*row,tmprows,((*colptr)[*n]-1)*sizeof(PetscInt));CHKERRQ(ierr);
18643801a69SSatish Balay     ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscScalar),values);CHKERRQ(ierr);
18770fe17b1SSatish Balay     ierr = PetscMemcpy(*values,tmpvalues,((*colptr)[*n]-1)*sizeof(PetscScalar));CHKERRQ(ierr);
18870fe17b1SSatish Balay     ierr = PetscFree3(tmpvalues,tmprows,tmpcolptr);CHKERRQ(ierr);
18943801a69SSatish Balay   }
1903bf14a46SMatthew Knepley   PetscFunctionReturn(0);
1913bf14a46SMatthew Knepley }
1923bf14a46SMatthew Knepley 
1933bf14a46SMatthew Knepley 
1943bf14a46SMatthew Knepley 
1953bf14a46SMatthew Knepley #undef __FUNCT__
1963bf14a46SMatthew Knepley #define __FUNCT__ "MatDestroy_Pastix"
1973bf14a46SMatthew Knepley /*
1983bf14a46SMatthew Knepley   Call clean step of PaStiX if lu->CleanUpPastix == true.
1993bf14a46SMatthew Knepley   Free the CSC matrix.
2003bf14a46SMatthew Knepley  */
2013bf14a46SMatthew Knepley PetscErrorCode MatDestroy_Pastix(Mat A)
2023bf14a46SMatthew Knepley {
2033bf14a46SMatthew Knepley   Mat_Pastix      *lu=(Mat_Pastix*)A->spptr;
2043bf14a46SMatthew Knepley   PetscErrorCode   ierr;
2053bf14a46SMatthew Knepley   PetscMPIInt      size=lu->commSize;
206745c78f7SBarry Smith 
2073bf14a46SMatthew Knepley   PetscFunctionBegin;
2083bf14a46SMatthew Knepley   if (lu->CleanUpPastix) {
2093bf14a46SMatthew Knepley     /* Terminate instance, deallocate memories */
2103bf14a46SMatthew Knepley     if (size > 1){
2113bf14a46SMatthew Knepley       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
2123bf14a46SMatthew Knepley       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
213f31ce8a6SBarry Smith       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
2143bf14a46SMatthew Knepley     }
2153bf14a46SMatthew Knepley 
2163bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK]=API_TASK_CLEAN;
2173bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]=API_TASK_CLEAN;
2183bf14a46SMatthew Knepley 
2193bf14a46SMatthew Knepley     pastix((pastix_data_t **)&(lu->pastix_data),
2203bf14a46SMatthew Knepley 	                      lu->pastix_comm,
2213bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->n,
2223bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->colptr,
2233bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->row,
2243bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->val,
2253bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->perm,
2263bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->invp,
2273bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->rhs,
2283bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->rhsnbr,
2293bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->iparm,
2303bf14a46SMatthew Knepley 	                      lu->dparm);
2313bf14a46SMatthew Knepley 
2323bf14a46SMatthew Knepley     ierr = PetscFree(lu->colptr);CHKERRQ(ierr);
2333bf14a46SMatthew Knepley     ierr = PetscFree(lu->row);  CHKERRQ(ierr);
2343bf14a46SMatthew Knepley     ierr = PetscFree(lu->val);  CHKERRQ(ierr);
2353bf14a46SMatthew Knepley     ierr = PetscFree(lu->perm); CHKERRQ(ierr);
2363bf14a46SMatthew Knepley     ierr = PetscFree(lu->invp); CHKERRQ(ierr);
2373bf14a46SMatthew Knepley     ierr = MPI_Comm_free(&(lu->pastix_comm));CHKERRQ(ierr);
2383bf14a46SMatthew Knepley   }
2393bf14a46SMatthew Knepley   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
2403bf14a46SMatthew Knepley   PetscFunctionReturn(0);
2413bf14a46SMatthew Knepley }
2423bf14a46SMatthew Knepley 
2433bf14a46SMatthew Knepley #undef __FUNCT__
2443bf14a46SMatthew Knepley #define __FUNCT__ "MatSolve_PaStiX"
2453bf14a46SMatthew Knepley /*
2463bf14a46SMatthew Knepley   Gather right-hand-side.
2473bf14a46SMatthew Knepley   Call for Solve step.
2483bf14a46SMatthew Knepley   Scatter solution.
2493bf14a46SMatthew Knepley  */
2503bf14a46SMatthew Knepley PetscErrorCode MatSolve_PaStiX(Mat A,Vec b,Vec x)
2513bf14a46SMatthew Knepley {
2523bf14a46SMatthew Knepley   Mat_Pastix     *lu=(Mat_Pastix*)A->spptr;
2533bf14a46SMatthew Knepley   PetscScalar    *array;
2543bf14a46SMatthew Knepley   Vec             x_seq;
2553bf14a46SMatthew Knepley   PetscErrorCode  ierr;
2563bf14a46SMatthew Knepley 
2573bf14a46SMatthew Knepley   PetscFunctionBegin;
2583bf14a46SMatthew Knepley   lu->rhsnbr = 1;
2593bf14a46SMatthew Knepley   x_seq = lu->b_seq;
2603bf14a46SMatthew Knepley   if (lu->commSize > 1){
2613bf14a46SMatthew Knepley     /* PaStiX only supports centralized rhs. Scatter b into a seqential rhs vector */
2623bf14a46SMatthew Knepley     ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2633bf14a46SMatthew Knepley     ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
264b5e56a35SBarry Smith     ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);
26541c8de11SBarry Smith   } else {  /* size == 1 */
2663bf14a46SMatthew Knepley     ierr = VecCopy(b,x);CHKERRQ(ierr);
2673bf14a46SMatthew Knepley     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
2683bf14a46SMatthew Knepley   }
2693bf14a46SMatthew Knepley   lu->rhs = array;
2703bf14a46SMatthew Knepley   if (lu->commSize == 1){
2713bf14a46SMatthew Knepley     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
2723bf14a46SMatthew Knepley   } else {
2733bf14a46SMatthew Knepley     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
2743bf14a46SMatthew Knepley   }
2753bf14a46SMatthew Knepley 
2763bf14a46SMatthew Knepley   /* solve phase */
2773bf14a46SMatthew Knepley   /*-------------*/
2783bf14a46SMatthew Knepley   lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE;
2793bf14a46SMatthew Knepley   lu->iparm[IPARM_END_TASK]   = API_TASK_REFINE;
280745c78f7SBarry Smith   lu->iparm[IPARM_RHS_MAKING] = API_RHS_B;
2813bf14a46SMatthew Knepley 
2823bf14a46SMatthew Knepley   pastix((pastix_data_t **)&(lu->pastix_data),
2833bf14a46SMatthew Knepley 	 (MPI_Comm)         lu->pastix_comm,
2843bf14a46SMatthew Knepley 	 (pastix_int_t)     lu->n,
2853bf14a46SMatthew Knepley 	 (pastix_int_t*)    lu->colptr,
2863bf14a46SMatthew Knepley 	 (pastix_int_t*)    lu->row,
2873bf14a46SMatthew Knepley 	 (pastix_float_t*)  lu->val,
2883bf14a46SMatthew Knepley 	 (pastix_int_t*)    lu->perm,
2893bf14a46SMatthew Knepley 	 (pastix_int_t*)    lu->invp,
2903bf14a46SMatthew Knepley 	 (pastix_float_t*)  lu->rhs,
2913bf14a46SMatthew Knepley 	 (pastix_int_t)     lu->rhsnbr,
2923bf14a46SMatthew Knepley 	 (pastix_int_t*)    lu->iparm,
2933bf14a46SMatthew Knepley 	 (double*)          lu->dparm);
2943bf14a46SMatthew Knepley 
2953bf14a46SMatthew Knepley   if (lu->iparm[IPARM_ERROR_NUMBER] < 0) {
296*e32f2f54SBarry Smith     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] );
2973bf14a46SMatthew Knepley   }
2983bf14a46SMatthew Knepley 
2993bf14a46SMatthew Knepley   if (lu->commSize == 1){
3003bf14a46SMatthew Knepley     ierr = VecRestoreArray(x,&(lu->rhs));CHKERRQ(ierr);
3013bf14a46SMatthew Knepley   } else {
3023bf14a46SMatthew Knepley     ierr = VecRestoreArray(x_seq,&(lu->rhs));CHKERRQ(ierr);
3033bf14a46SMatthew Knepley   }
3043bf14a46SMatthew Knepley 
3053bf14a46SMatthew Knepley   if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */
3063bf14a46SMatthew Knepley     ierr = VecScatterBegin(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3073bf14a46SMatthew Knepley     ierr = VecScatterEnd(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3083bf14a46SMatthew Knepley   }
3093bf14a46SMatthew Knepley   PetscFunctionReturn(0);
3103bf14a46SMatthew Knepley }
3113bf14a46SMatthew Knepley 
3123bf14a46SMatthew Knepley /*
3133bf14a46SMatthew Knepley   Numeric factorisation using PaStiX solver.
3143bf14a46SMatthew Knepley 
3153bf14a46SMatthew Knepley  */
3163bf14a46SMatthew Knepley #undef __FUNCT__
3173bf14a46SMatthew Knepley #define __FUNCT__ "MatFactorNumeric_PASTIX"
3183bf14a46SMatthew Knepley PetscErrorCode MatFactorNumeric_PaStiX(Mat F,Mat A,const MatFactorInfo *info)
3193bf14a46SMatthew Knepley {
3203bf14a46SMatthew Knepley   Mat_Pastix    *lu =(Mat_Pastix*)(F)->spptr;
32141c8de11SBarry Smith   Mat           *tseq;
3223bf14a46SMatthew Knepley   PetscErrorCode ierr = 0;
323b5e56a35SBarry Smith   PetscInt       icntl;
324b5e56a35SBarry Smith   PetscInt       M=A->rmap->N;
3253bf14a46SMatthew Knepley   PetscTruth     valOnly,flg, isSym;
3263bf14a46SMatthew Knepley   Mat            F_diag;
3273bf14a46SMatthew Knepley   IS             is_iden;
3283bf14a46SMatthew Knepley   Vec            b;
3293bf14a46SMatthew Knepley   IS             isrow;
3303bf14a46SMatthew Knepley   PetscTruth     isSeqAIJ,isSeqSBAIJ;
3313bf14a46SMatthew Knepley 
3323bf14a46SMatthew Knepley   PetscFunctionBegin;
33341c8de11SBarry Smith 
3343bf14a46SMatthew Knepley   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
3353bf14a46SMatthew Knepley   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
3363bf14a46SMatthew Knepley   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
3373bf14a46SMatthew Knepley     (F)->ops->solve   = MatSolve_PaStiX;
3383bf14a46SMatthew Knepley 
3393bf14a46SMatthew Knepley     /* Initialize a PASTIX instance */
3403bf14a46SMatthew Knepley     ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr);
3413bf14a46SMatthew Knepley     ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank);         CHKERRQ(ierr);
3423bf14a46SMatthew Knepley     ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize);         CHKERRQ(ierr);
3433bf14a46SMatthew Knepley 
3443bf14a46SMatthew Knepley     /* Set pastix options */
3453bf14a46SMatthew Knepley     lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO;
3463bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK]       = API_TASK_INIT;
3473bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]         = API_TASK_INIT;
3483bf14a46SMatthew Knepley     lu->rhsnbr = 1;
3493bf14a46SMatthew Knepley 
3503bf14a46SMatthew Knepley     /* Call to set default pastix options */
3513bf14a46SMatthew Knepley     pastix((pastix_data_t **)&(lu->pastix_data),
3523bf14a46SMatthew Knepley 	   (MPI_Comm)         lu->pastix_comm,
3533bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->n,
3543bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->colptr,
3553bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->row,
3563bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->val,
3573bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->perm,
3583bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->invp,
3593bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->rhs,
3603bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->rhsnbr,
3613bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->iparm,
3623bf14a46SMatthew Knepley 	   (double*)          lu->dparm);
3633bf14a46SMatthew Knepley 
3643bf14a46SMatthew Knepley     ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr);
3653bf14a46SMatthew Knepley 
3663bf14a46SMatthew Knepley     icntl=-1;
36741c8de11SBarry Smith     lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT;
36841c8de11SBarry Smith     ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr);
3693bf14a46SMatthew Knepley     if ((flg && icntl > 0) || PetscLogPrintInfo) {
3703bf14a46SMatthew Knepley       lu->iparm[IPARM_VERBOSE] =  icntl;
3713bf14a46SMatthew Knepley     }
3723bf14a46SMatthew Knepley     icntl=-1;
37341c8de11SBarry 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);
3743bf14a46SMatthew Knepley     if ((flg && icntl > 0)) {
3753bf14a46SMatthew Knepley       lu->iparm[IPARM_THREAD_NBR] = icntl;
3763bf14a46SMatthew Knepley     }
3773bf14a46SMatthew Knepley     PetscOptionsEnd();
3783bf14a46SMatthew Knepley     valOnly = PETSC_FALSE;
37941c8de11SBarry Smith   }  else {
3803bf14a46SMatthew Knepley     valOnly = PETSC_TRUE;
3813bf14a46SMatthew Knepley   }
3823bf14a46SMatthew Knepley 
3833bf14a46SMatthew Knepley   lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES;
3843bf14a46SMatthew Knepley 
3853bf14a46SMatthew Knepley   /* convert mpi A to seq mat A */
3863bf14a46SMatthew Knepley   ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr);
3873bf14a46SMatthew Knepley   ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr);
3883bf14a46SMatthew Knepley   ierr = ISDestroy(isrow);CHKERRQ(ierr);
3893bf14a46SMatthew Knepley 
39041c8de11SBarry Smith   ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val);CHKERRQ(ierr);
39141c8de11SBarry Smith   ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr);
39241c8de11SBarry Smith   ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr);
39341c8de11SBarry Smith 
3943bf14a46SMatthew Knepley   ierr = PetscMalloc((lu->n)*sizeof(PetscInt)   ,&(lu->perm));CHKERRQ(ierr);
3953bf14a46SMatthew Knepley   ierr = PetscMalloc((lu->n)*sizeof(PetscInt)   ,&(lu->invp));CHKERRQ(ierr);
3963bf14a46SMatthew Knepley 
3973bf14a46SMatthew Knepley   if (isSym) {
398745c78f7SBarry Smith     /* On symmetric matrix, LLT */
3993bf14a46SMatthew Knepley     lu->iparm[IPARM_SYM] = API_SYM_YES;
40041c8de11SBarry Smith     lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT;
401f31ce8a6SBarry Smith   } else {
402745c78f7SBarry Smith     /* On unsymmetric matrix, LU */
4033bf14a46SMatthew Knepley     lu->iparm[IPARM_SYM] = API_SYM_NO;
4043bf14a46SMatthew Knepley     lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
4053bf14a46SMatthew Knepley   }
4063bf14a46SMatthew Knepley 
4073bf14a46SMatthew Knepley   /*----------------*/
4083bf14a46SMatthew Knepley   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
4093bf14a46SMatthew Knepley     if (!(isSeqAIJ || isSeqSBAIJ)) {
4103bf14a46SMatthew Knepley       /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
4113bf14a46SMatthew Knepley 	ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
4123bf14a46SMatthew Knepley 	ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
4133bf14a46SMatthew Knepley 	ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
4143bf14a46SMatthew Knepley 	ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
4153bf14a46SMatthew Knepley 	ierr = VecSetFromOptions(b);CHKERRQ(ierr);
4163bf14a46SMatthew Knepley 
4173bf14a46SMatthew Knepley 	ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
4183bf14a46SMatthew Knepley 	ierr = VecScatterCreate(lu->b_seq,is_iden,b,is_iden,&lu->scat_sol);CHKERRQ(ierr);
4193bf14a46SMatthew Knepley 	ierr = ISDestroy(is_iden);CHKERRQ(ierr);
4203bf14a46SMatthew Knepley 	ierr = VecDestroy(b);CHKERRQ(ierr);
4213bf14a46SMatthew Knepley     }
4223bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING;
4233bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
4243bf14a46SMatthew Knepley 
4253bf14a46SMatthew Knepley     pastix((pastix_data_t **)&(lu->pastix_data),
4263bf14a46SMatthew Knepley 	   (MPI_Comm)         lu->pastix_comm,
4273bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->n,
4283bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->colptr,
4293bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->row,
4303bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->val,
4313bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->perm,
4323bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->invp,
4333bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->rhs,
4343bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->rhsnbr,
4353bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->iparm,
4363bf14a46SMatthew Knepley 	   (double*)          lu->dparm);
4373bf14a46SMatthew Knepley     if (lu->iparm[IPARM_ERROR_NUMBER] < 0) {
438*e32f2f54SBarry Smith       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]);
4393bf14a46SMatthew Knepley     }
44041c8de11SBarry Smith   } else {
4413bf14a46SMatthew Knepley     lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT;
4423bf14a46SMatthew Knepley     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
4433bf14a46SMatthew Knepley     pastix((pastix_data_t **)&(lu->pastix_data),
4443bf14a46SMatthew Knepley 	   (MPI_Comm)         lu->pastix_comm,
4453bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->n,
4463bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->colptr,
4473bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->row,
4483bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->val,
4493bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->perm,
4503bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->invp,
4513bf14a46SMatthew Knepley 	   (pastix_float_t*)  lu->rhs,
4523bf14a46SMatthew Knepley 	   (pastix_int_t)     lu->rhsnbr,
4533bf14a46SMatthew Knepley 	   (pastix_int_t*)    lu->iparm,
4543bf14a46SMatthew Knepley 	   (double*)          lu->dparm);
4553bf14a46SMatthew Knepley 
4563bf14a46SMatthew Knepley     if (lu->iparm[IPARM_ERROR_NUMBER] < 0) {
457*e32f2f54SBarry Smith       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]);
4583bf14a46SMatthew Knepley     }
4593bf14a46SMatthew Knepley   }
4603bf14a46SMatthew Knepley 
4613bf14a46SMatthew Knepley   if (lu->commSize > 1){
462d5f3da31SBarry Smith     if ((F)->factortype == MAT_FACTOR_LU){
4633bf14a46SMatthew Knepley       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
4643bf14a46SMatthew Knepley     } else {
4653bf14a46SMatthew Knepley       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
4663bf14a46SMatthew Knepley     }
4673bf14a46SMatthew Knepley     F_diag->assembled = PETSC_TRUE;
4683bf14a46SMatthew Knepley   }
4693bf14a46SMatthew Knepley   (F)->assembled     = PETSC_TRUE;
4703bf14a46SMatthew Knepley   lu->matstruc       = SAME_NONZERO_PATTERN;
4713bf14a46SMatthew Knepley   lu->CleanUpPastix  = PETSC_TRUE;
4723bf14a46SMatthew Knepley   PetscFunctionReturn(0);
4733bf14a46SMatthew Knepley }
4743bf14a46SMatthew Knepley 
4753bf14a46SMatthew Knepley /* Note the Petsc r and c permutations are ignored */
4763bf14a46SMatthew Knepley #undef __FUNCT__
4773bf14a46SMatthew Knepley #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX"
4783bf14a46SMatthew Knepley PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
4793bf14a46SMatthew Knepley {
4803bf14a46SMatthew Knepley   Mat_Pastix      *lu = (Mat_Pastix*)F->spptr;
4813bf14a46SMatthew Knepley 
4823bf14a46SMatthew Knepley   PetscFunctionBegin;
4833bf14a46SMatthew Knepley   lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
4843bf14a46SMatthew Knepley   lu->iparm[IPARM_SYM]           = API_SYM_YES;
4853bf14a46SMatthew Knepley   lu->matstruc                   = DIFFERENT_NONZERO_PATTERN;
4863bf14a46SMatthew Knepley   F->ops->lufactornumeric        = MatFactorNumeric_PaStiX;
4873bf14a46SMatthew Knepley   PetscFunctionReturn(0);
4883bf14a46SMatthew Knepley }
4893bf14a46SMatthew Knepley 
4903bf14a46SMatthew Knepley 
4913bf14a46SMatthew Knepley /* Note the Petsc r permutation is ignored */
4923bf14a46SMatthew Knepley #undef __FUNCT__
4933bf14a46SMatthew Knepley #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX"
4943bf14a46SMatthew Knepley PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info)
4953bf14a46SMatthew Knepley {
4963bf14a46SMatthew Knepley   Mat_Pastix      *lu = (Mat_Pastix*)(F)->spptr;
4973bf14a46SMatthew Knepley 
4983bf14a46SMatthew Knepley   PetscFunctionBegin;
4993bf14a46SMatthew Knepley   lu->iparm[IPARM_FACTORIZATION]  = API_FACT_LLT;
5003bf14a46SMatthew Knepley   lu->iparm[IPARM_SYM]            = API_SYM_NO;
5013bf14a46SMatthew Knepley   lu->matstruc                    = DIFFERENT_NONZERO_PATTERN;
5023bf14a46SMatthew Knepley   (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX;
5033bf14a46SMatthew Knepley   PetscFunctionReturn(0);
5043bf14a46SMatthew Knepley }
5053bf14a46SMatthew Knepley 
5063bf14a46SMatthew Knepley #undef __FUNCT__
5073bf14a46SMatthew Knepley #define __FUNCT__ "MatView_PaStiX"
5083bf14a46SMatthew Knepley PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer)
5093bf14a46SMatthew Knepley {
5103bf14a46SMatthew Knepley   PetscErrorCode    ierr;
5113bf14a46SMatthew Knepley   PetscTruth        iascii;
5123bf14a46SMatthew Knepley   PetscViewerFormat format;
5133bf14a46SMatthew Knepley 
5143bf14a46SMatthew Knepley   PetscFunctionBegin;
5153bf14a46SMatthew Knepley   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
5163bf14a46SMatthew Knepley   if (iascii) {
5173bf14a46SMatthew Knepley     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
5183bf14a46SMatthew Knepley     if (format == PETSC_VIEWER_ASCII_INFO){
519b5e56a35SBarry Smith       Mat_Pastix      *lu=(Mat_Pastix*)A->spptr;
520b5e56a35SBarry Smith 
521b5e56a35SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr);
522b5e56a35SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Matrix type :                      %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES)?"Symmetric":"Unsymmetric"));CHKERRQ(ierr);
523b5e56a35SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Level of printing (0,1,2):         %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr);
524b5e56a35SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr);
525b5e56a35SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Error :                        %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr);
5263bf14a46SMatthew Knepley     }
5273bf14a46SMatthew Knepley   }
5283bf14a46SMatthew Knepley   PetscFunctionReturn(0);
5293bf14a46SMatthew Knepley }
5303bf14a46SMatthew Knepley 
5313bf14a46SMatthew Knepley 
5323bf14a46SMatthew Knepley /*MC
533b5e56a35SBarry Smith      MAT_SOLVER_PASTIX  - A solver package providing direct solvers (LU) for distributed
5343bf14a46SMatthew Knepley   and sequential matrices via the external package PaStiX.
5353bf14a46SMatthew Knepley 
536e2e64c6bSBarry Smith   Use ./configure --download-pastix to have PETSc installed with PaStiX
5373bf14a46SMatthew Knepley 
5383bf14a46SMatthew Knepley   Options Database Keys:
539b5e56a35SBarry Smith + -mat_pastix_verbose   <0,1,2>   - print level
540b5e56a35SBarry Smith - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task.
5413bf14a46SMatthew Knepley 
5423bf14a46SMatthew Knepley   Level: beginner
5433bf14a46SMatthew Knepley 
54441c8de11SBarry Smith .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
54541c8de11SBarry Smith 
5463bf14a46SMatthew Knepley M*/
5473bf14a46SMatthew Knepley 
5483bf14a46SMatthew Knepley 
5493bf14a46SMatthew Knepley #undef __FUNCT__
5503bf14a46SMatthew Knepley #define __FUNCT__ "MatGetInfo_PaStiX"
5513bf14a46SMatthew Knepley PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info)
5523bf14a46SMatthew Knepley {
5533bf14a46SMatthew Knepley     Mat_Pastix  *lu =(Mat_Pastix*)A->spptr;
5543bf14a46SMatthew Knepley 
5553bf14a46SMatthew Knepley     PetscFunctionBegin;
5563bf14a46SMatthew Knepley     info->block_size        = 1.0;
5573bf14a46SMatthew Knepley     info->nz_allocated      = lu->iparm[IPARM_NNZEROS];
5583bf14a46SMatthew Knepley     info->nz_used           = lu->iparm[IPARM_NNZEROS];
5593bf14a46SMatthew Knepley     info->nz_unneeded       = 0.0;
5603bf14a46SMatthew Knepley     info->assemblies        = 0.0;
5613bf14a46SMatthew Knepley     info->mallocs           = 0.0;
5623bf14a46SMatthew Knepley     info->memory            = 0.0;
5633bf14a46SMatthew Knepley     info->fill_ratio_given  = 0;
5643bf14a46SMatthew Knepley     info->fill_ratio_needed = 0;
5653bf14a46SMatthew Knepley     info->factor_mallocs    = 0;
5663bf14a46SMatthew Knepley     PetscFunctionReturn(0);
5673bf14a46SMatthew Knepley }
5683bf14a46SMatthew Knepley 
5693bf14a46SMatthew Knepley EXTERN_C_BEGIN
5703bf14a46SMatthew Knepley #undef __FUNCT__
5713bf14a46SMatthew Knepley #define __FUNCT__ "MatFactorGetSolverPackage_pastix"
5723bf14a46SMatthew Knepley PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type)
5733bf14a46SMatthew Knepley {
5743bf14a46SMatthew Knepley   PetscFunctionBegin;
5753bf14a46SMatthew Knepley   *type = MAT_SOLVER_PASTIX;
5763bf14a46SMatthew Knepley   PetscFunctionReturn(0);
5773bf14a46SMatthew Knepley }
5783bf14a46SMatthew Knepley EXTERN_C_END
5793bf14a46SMatthew Knepley 
5803bf14a46SMatthew Knepley EXTERN_C_BEGIN
5813bf14a46SMatthew Knepley /*
5823bf14a46SMatthew Knepley     The seq and mpi versions of this function are the same
5833bf14a46SMatthew Knepley */
5843bf14a46SMatthew Knepley #undef __FUNCT__
5853bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_seqaij_pastix"
5863bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F)
5873bf14a46SMatthew Knepley {
5883bf14a46SMatthew Knepley   Mat            B;
5893bf14a46SMatthew Knepley   PetscErrorCode ierr;
5903bf14a46SMatthew Knepley   Mat_Pastix    *pastix;
5913bf14a46SMatthew Knepley 
5923bf14a46SMatthew Knepley   PetscFunctionBegin;
5933bf14a46SMatthew Knepley   if (ftype != MAT_FACTOR_LU) {
594*e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
5953bf14a46SMatthew Knepley   }
5963bf14a46SMatthew Knepley   /* Create the factorization matrix */
5973bf14a46SMatthew Knepley   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
5983bf14a46SMatthew Knepley   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
5993bf14a46SMatthew Knepley   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
6003bf14a46SMatthew Knepley   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
6013bf14a46SMatthew Knepley 
6023bf14a46SMatthew Knepley   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
6033bf14a46SMatthew Knepley   B->ops->view             = MatView_PaStiX;
6043bf14a46SMatthew Knepley   B->ops->getinfo          = MatGetInfo_PaStiX;
60531e762f5SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix", MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
606d5f3da31SBarry Smith   B->factortype            = MAT_FACTOR_LU;
6073bf14a46SMatthew Knepley 
6083bf14a46SMatthew Knepley   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
6093bf14a46SMatthew Knepley   pastix->CleanUpPastix             = PETSC_FALSE;
6103bf14a46SMatthew Knepley   pastix->isAIJ                     = PETSC_TRUE;
6113bf14a46SMatthew Knepley   pastix->scat_rhs                  = PETSC_NULL;
6123bf14a46SMatthew Knepley   pastix->scat_sol                  = PETSC_NULL;
6133bf14a46SMatthew Knepley   pastix->MatDestroy                = B->ops->destroy;
6143bf14a46SMatthew Knepley   B->ops->destroy                   = MatDestroy_Pastix;
6153bf14a46SMatthew Knepley   B->spptr                          = (void*)pastix;
6163bf14a46SMatthew Knepley 
6173bf14a46SMatthew Knepley   *F = B;
6183bf14a46SMatthew Knepley   PetscFunctionReturn(0);
6193bf14a46SMatthew Knepley }
6203bf14a46SMatthew Knepley EXTERN_C_END
6213bf14a46SMatthew Knepley 
622b5e56a35SBarry Smith 
6233bf14a46SMatthew Knepley EXTERN_C_BEGIN
6243bf14a46SMatthew Knepley #undef __FUNCT__
6253bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_mpiaij_pastix"
6263bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F)
6273bf14a46SMatthew Knepley {
6283bf14a46SMatthew Knepley   Mat            B;
6293bf14a46SMatthew Knepley   PetscErrorCode ierr;
6303bf14a46SMatthew Knepley   Mat_Pastix    *pastix;
6313bf14a46SMatthew Knepley 
6323bf14a46SMatthew Knepley   PetscFunctionBegin;
633*e32f2f54SBarry Smith   if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
6343bf14a46SMatthew Knepley   /* Create the factorization matrix */
6353bf14a46SMatthew Knepley   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
6363bf14a46SMatthew Knepley   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
6373bf14a46SMatthew Knepley   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
6383bf14a46SMatthew Knepley   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
6393bf14a46SMatthew Knepley   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
6403bf14a46SMatthew Knepley 
6413bf14a46SMatthew Knepley   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
6423bf14a46SMatthew Knepley   B->ops->view             = MatView_PaStiX;
64331e762f5SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
644d5f3da31SBarry Smith   B->factortype            = MAT_FACTOR_LU;
6453bf14a46SMatthew Knepley 
6463bf14a46SMatthew Knepley   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
6473bf14a46SMatthew Knepley   pastix->CleanUpPastix             = PETSC_FALSE;
6483bf14a46SMatthew Knepley   pastix->isAIJ                     = PETSC_TRUE;
6493bf14a46SMatthew Knepley   pastix->scat_rhs                  = PETSC_NULL;
6503bf14a46SMatthew Knepley   pastix->scat_sol                  = PETSC_NULL;
6513bf14a46SMatthew Knepley   pastix->MatDestroy                = B->ops->destroy;
6523bf14a46SMatthew Knepley   B->ops->destroy                  = MatDestroy_Pastix;
6533bf14a46SMatthew Knepley   B->spptr                         = (void*)pastix;
6543bf14a46SMatthew Knepley 
6553bf14a46SMatthew Knepley   *F = B;
6563bf14a46SMatthew Knepley   PetscFunctionReturn(0);
6573bf14a46SMatthew Knepley }
6583bf14a46SMatthew Knepley EXTERN_C_END
6593bf14a46SMatthew Knepley 
6603bf14a46SMatthew Knepley EXTERN_C_BEGIN
6613bf14a46SMatthew Knepley #undef __FUNCT__
6623bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_seqsbaij_pastix"
6633bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F)
6643bf14a46SMatthew Knepley {
6653bf14a46SMatthew Knepley   Mat            B;
6663bf14a46SMatthew Knepley   PetscErrorCode ierr;
6673bf14a46SMatthew Knepley   Mat_Pastix    *pastix;
6683bf14a46SMatthew Knepley 
6693bf14a46SMatthew Knepley   PetscFunctionBegin;
6703bf14a46SMatthew Knepley   if (ftype != MAT_FACTOR_CHOLESKY) {
671*e32f2f54SBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
6723bf14a46SMatthew Knepley   }
6733bf14a46SMatthew Knepley   /* Create the factorization matrix */
6743bf14a46SMatthew Knepley   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
6753bf14a46SMatthew Knepley   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
6763bf14a46SMatthew Knepley   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
6773bf14a46SMatthew Knepley   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
6783bf14a46SMatthew Knepley   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
6793bf14a46SMatthew Knepley 
6803bf14a46SMatthew Knepley   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
6813bf14a46SMatthew Knepley   B->ops->view                   = MatView_PaStiX;
68231e762f5SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
683d5f3da31SBarry Smith   B->factortype                  = MAT_FACTOR_CHOLESKY;
6843bf14a46SMatthew Knepley 
6853bf14a46SMatthew Knepley   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
6863bf14a46SMatthew Knepley   pastix->CleanUpPastix             = PETSC_FALSE;
6873bf14a46SMatthew Knepley   pastix->isAIJ                     = PETSC_TRUE;
6883bf14a46SMatthew Knepley   pastix->scat_rhs                  = PETSC_NULL;
6893bf14a46SMatthew Knepley   pastix->scat_sol                  = PETSC_NULL;
6903bf14a46SMatthew Knepley   pastix->MatDestroy                = B->ops->destroy;
6913bf14a46SMatthew Knepley   B->ops->destroy                  = MatDestroy_Pastix;
6923bf14a46SMatthew Knepley   B->spptr                         = (void*)pastix;
6933bf14a46SMatthew Knepley 
6943bf14a46SMatthew Knepley   *F = B;
6953bf14a46SMatthew Knepley   PetscFunctionReturn(0);
6963bf14a46SMatthew Knepley }
6973bf14a46SMatthew Knepley EXTERN_C_END
6983bf14a46SMatthew Knepley 
6993bf14a46SMatthew Knepley EXTERN_C_BEGIN
7003bf14a46SMatthew Knepley #undef __FUNCT__
7013bf14a46SMatthew Knepley #define __FUNCT__ "MatGetFactor_mpisbaij_pastix"
7023bf14a46SMatthew Knepley PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F)
7033bf14a46SMatthew Knepley {
7043bf14a46SMatthew Knepley   Mat            B;
7053bf14a46SMatthew Knepley   PetscErrorCode ierr;
7063bf14a46SMatthew Knepley   Mat_Pastix    *pastix;
7073bf14a46SMatthew Knepley 
7083bf14a46SMatthew Knepley   PetscFunctionBegin;
709*e32f2f54SBarry Smith   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
71041c8de11SBarry Smith 
7113bf14a46SMatthew Knepley   /* Create the factorization matrix */
7123bf14a46SMatthew Knepley   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
7133bf14a46SMatthew Knepley   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
7143bf14a46SMatthew Knepley   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
7153bf14a46SMatthew Knepley   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
7163bf14a46SMatthew Knepley   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
7173bf14a46SMatthew Knepley 
7183bf14a46SMatthew Knepley   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
7193bf14a46SMatthew Knepley   B->ops->view                   = MatView_PaStiX;
72031e762f5SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
721d5f3da31SBarry Smith   B->factortype                  = MAT_FACTOR_CHOLESKY;
7223bf14a46SMatthew Knepley 
7233bf14a46SMatthew Knepley   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
7243bf14a46SMatthew Knepley   pastix->CleanUpPastix             = PETSC_FALSE;
7253bf14a46SMatthew Knepley   pastix->isAIJ                     = PETSC_TRUE;
7263bf14a46SMatthew Knepley   pastix->scat_rhs                  = PETSC_NULL;
7273bf14a46SMatthew Knepley   pastix->scat_sol                  = PETSC_NULL;
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
736