xref: /petsc/src/mat/impls/aij/mpi/pastix/pastix.c (revision 01a79839fc82a7dabb7a87cd2a8bb532c6bfa88d)
1 #define PETSCMAT_DLL
2 
3 /*
4     Provides an interface to the PaStiX sparse solver
5 */
6 #include "../src/mat/impls/aij/seq/aij.h"
7 #include "../src/mat/impls/aij/mpi/mpiaij.h"
8 #include "../src/mat/impls/sbaij/seq/sbaij.h"
9 #include "../src/mat/impls/sbaij/mpi/mpisbaij.h"
10 
11 #if defined(PETSC_HAVE_STDLIB_H)
12 #include <stdlib.h>
13 #endif
14 #if defined(PETSC_HAVE_STRING_H)
15 #include <string.h>
16 #endif
17 
18 EXTERN_C_BEGIN
19 #include "pastix.h"
20 EXTERN_C_END
21 
22 typedef struct Mat_Pastix_ {
23   pastix_data_t *pastix_data;              /* Pastix data storage structure                        */
24   MatStructure   matstruc;
25   PetscInt       n;                        /* Number of columns in the matrix                      */
26   PetscInt       *colptr;                  /* Index of first element of each column in row and val */
27   PetscInt       *row;                     /* Row of each element of the matrix                    */
28   PetscScalar    *val;                     /* Value of each element of the matrix                  */
29   PetscInt       *perm;                    /* Permutation tabular                                  */
30   PetscInt       *invp;                    /* Reverse permutation tabular                          */
31   PetscScalar    *rhs;                     /* Rhight-hand-side member                              */
32   PetscInt       rhsnbr;                   /* Rhight-hand-side number (must be 1)                  */
33   PetscInt       iparm[64];                /* Integer parameters                                   */
34   double         dparm[64];                /* Floating point parameters                            */
35   MPI_Comm       pastix_comm;              /* PaStiX MPI communicator                              */
36   PetscMPIInt    commRank;                 /* MPI rank                                             */
37   PetscMPIInt    commSize;                 /* MPI communicator size                                */
38   PetscBool      CleanUpPastix;            /* Boolean indicating if we call PaStiX clean step      */
39   VecScatter     scat_rhs;
40   VecScatter     scat_sol;
41   Vec            b_seq;
42   PetscBool      isAIJ;
43   PetscErrorCode (*MatDestroy)(Mat);
44 } Mat_Pastix;
45 
46 extern PetscErrorCode MatDuplicate_Pastix(Mat,MatDuplicateOption,Mat*);
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "MatConvertToCSC"
50 /*
51    convert Petsc seqaij matrix to CSC: colptr[n], row[nz], val[nz]
52 
53   input:
54     A       - matrix in seqaij or mpisbaij (bs=1) format
55     valOnly - FALSE: spaces are allocated and values are set for the CSC
56               TRUE:  Only fill values
57   output:
58     n       - Size of the matrix
59     colptr  - Index of first element of each column in row and val
60     row     - Row of each element of the matrix
61     values  - Value of each element of the matrix
62  */
63 PetscErrorCode MatConvertToCSC(Mat A,PetscBool  valOnly,PetscInt *n,PetscInt **colptr,PetscInt **row,PetscScalar **values)
64 {
65   Mat_SeqAIJ     *aa      = (Mat_SeqAIJ*)A->data;
66   PetscInt       *rowptr  = aa->i;
67   PetscInt       *col     = aa->j;
68   PetscScalar    *rvalues = aa->a;
69   PetscInt        m       = A->rmap->N;
70   PetscInt        nnz;
71   PetscInt        i,j, k;
72   PetscInt        base = 1;
73   PetscInt        idx;
74   PetscErrorCode  ierr;
75   PetscInt        colidx;
76   PetscInt       *colcount;
77   PetscBool       isSBAIJ;
78   PetscBool       isSeqSBAIJ;
79   PetscBool       isMpiSBAIJ;
80   PetscBool       isSym;
81 
82   PetscFunctionBegin;
83   ierr = MatIsSymmetric(A,0.0,&isSym);CHKERRQ(ierr);
84   ierr = PetscTypeCompare((PetscObject)A,MATSBAIJ,&isSBAIJ);CHKERRQ(ierr);
85   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
86   ierr = PetscTypeCompare((PetscObject)A,MATMPISBAIJ,&isMpiSBAIJ);CHKERRQ(ierr);
87 
88   *n = A->cmap->N;
89 
90   /* PaStiX only needs triangular matrix if matrix is symmetric
91    */
92   if (isSym && !(isSBAIJ || isSeqSBAIJ || isMpiSBAIJ)) {
93     nnz = (aa->nz - *n)/2 + *n;
94   }
95   else {
96     nnz     = aa->nz;
97   }
98 
99   if (!valOnly){
100     ierr = PetscMalloc(((*n)+1) *sizeof(PetscInt)   ,colptr );CHKERRQ(ierr);
101     ierr = PetscMalloc( nnz     *sizeof(PetscInt)   ,row);CHKERRQ(ierr);
102     ierr = PetscMalloc( nnz     *sizeof(PetscScalar),values);CHKERRQ(ierr);
103 
104     if (isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) {
105 	ierr = PetscMemcpy (*colptr, rowptr, ((*n)+1)*sizeof(PetscInt));CHKERRQ(ierr);
106 	for (i = 0; i < *n+1; i++)
107 	  (*colptr)[i] += base;
108 	ierr = PetscMemcpy (*row, col, (nnz)*sizeof(PetscInt));CHKERRQ(ierr);
109 	for (i = 0; i < nnz; i++)
110 	  (*row)[i] += base;
111 	ierr = PetscMemcpy (*values, rvalues, (nnz)*sizeof(PetscScalar));CHKERRQ(ierr);
112     } else {
113       ierr = PetscMalloc((*n)*sizeof(PetscInt)   ,&colcount);CHKERRQ(ierr);
114 
115       for (i = 0; i < m; i++) colcount[i] = 0;
116       /* Fill-in colptr */
117       for (i = 0; i < m; i++) {
118 	for (j = rowptr[i]; j < rowptr[i+1]; j++) {
119 	  if (!isSym || col[j] <= i)  colcount[col[j]]++;
120         }
121       }
122 
123       (*colptr)[0] = base;
124       for (j = 0; j < *n; j++) {
125 	(*colptr)[j+1] = (*colptr)[j] + colcount[j];
126 	/* in next loop we fill starting from (*colptr)[colidx] - base */
127 	colcount[j] = -base;
128       }
129 
130       /* Fill-in rows and values */
131       for (i = 0; i < m; i++) {
132 	for (j = rowptr[i]; j < rowptr[i+1]; j++) {
133 	  if (!isSym || col[j] <= i) {
134 	    colidx = col[j];
135 	    idx    = (*colptr)[colidx] + colcount[colidx];
136 	    (*row)[idx]    = i + base;
137 	    (*values)[idx] = rvalues[j];
138 	    colcount[colidx]++;
139 	  }
140 	}
141       }
142       ierr = PetscFree(colcount);CHKERRQ(ierr);
143     }
144   } else {
145     /* Fill-in only values */
146     for (i = 0; i < m; i++) {
147       for (j = rowptr[i]; j < rowptr[i+1]; j++) {
148 	colidx = col[j];
149 	if ((isSBAIJ || isSeqSBAIJ || isMpiSBAIJ) ||!isSym || col[j] <= i)
150 	  {
151 	    /* look for the value to fill */
152 	    for (k = (*colptr)[colidx] - base; k < (*colptr)[colidx + 1] - base; k++) {
153 	      if (((*row)[k]-base) == i) {
154 		(*values)[k] = rvalues[j];
155 		break;
156 	      }
157 	    }
158 	    /* data structure of sparse matrix has changed */
159 	    if (k == (*colptr)[colidx + 1] - base) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"overflow on k %D",k);
160 	  }
161       }
162     }
163   }
164   {
165     PetscScalar *tmpvalues;
166     PetscInt    *tmprows,*tmpcolptr;
167     tmpvalues = (PetscScalar*)malloc(nnz*sizeof(PetscScalar)); if (!tmpvalues) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory");
168     tmprows   = (PetscInt*)malloc(nnz*sizeof(PetscInt));if (!tmprows) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory");
169     tmpcolptr = (PetscInt*)malloc((*n+1)*sizeof(PetscInt));if (!tmpcolptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory");
170 
171     ierr = PetscMemcpy(tmpcolptr,*colptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr);
172     ierr = PetscMemcpy(tmprows,*row,nnz*sizeof(PetscInt));CHKERRQ(ierr);
173     ierr = PetscMemcpy(tmpvalues,*values,nnz*sizeof(PetscScalar));CHKERRQ(ierr);
174     ierr = PetscFree(*row);CHKERRQ(ierr);
175     ierr = PetscFree(*values);CHKERRQ(ierr);
176 
177     pastix_checkMatrix(MPI_COMM_WORLD,API_VERBOSE_NO,((isSym != 0) ? API_SYM_YES : API_SYM_NO),API_YES,*n,&tmpcolptr,&tmprows,&tmpvalues,NULL,1);
178 
179     ierr = PetscMemcpy(*colptr,tmpcolptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr);
180     ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscInt),row);CHKERRQ(ierr);
181     ierr = PetscMemcpy(*row,tmprows,((*colptr)[*n]-1)*sizeof(PetscInt));CHKERRQ(ierr);
182     ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscScalar),values);CHKERRQ(ierr);
183     ierr = PetscMemcpy(*values,tmpvalues,((*colptr)[*n]-1)*sizeof(PetscScalar));CHKERRQ(ierr);
184     free(tmpvalues);
185     free(tmprows);
186     free(tmpcolptr);
187 
188   }
189   PetscFunctionReturn(0);
190 }
191 
192 
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "MatDestroy_Pastix"
196 /*
197   Call clean step of PaStiX if lu->CleanUpPastix == true.
198   Free the CSC matrix.
199  */
200 PetscErrorCode MatDestroy_Pastix(Mat A)
201 {
202   Mat_Pastix      *lu=(Mat_Pastix*)A->spptr;
203   PetscErrorCode   ierr;
204   PetscMPIInt      size=lu->commSize;
205 
206   PetscFunctionBegin;
207   if (lu->CleanUpPastix) {
208     /* Terminate instance, deallocate memories */
209     if (size > 1){
210       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
211       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
212       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
213     }
214 
215     lu->iparm[IPARM_START_TASK]=API_TASK_CLEAN;
216     lu->iparm[IPARM_END_TASK]=API_TASK_CLEAN;
217 
218     pastix((pastix_data_t **)&(lu->pastix_data),
219 	                      lu->pastix_comm,
220 	   (pastix_int_t)     lu->n,
221 	   (pastix_int_t*)    lu->colptr,
222 	   (pastix_int_t*)    lu->row,
223 	   (pastix_float_t*)  lu->val,
224 	   (pastix_int_t*)    lu->perm,
225 	   (pastix_int_t*)    lu->invp,
226 	   (pastix_float_t*)  lu->rhs,
227 	   (pastix_int_t)     lu->rhsnbr,
228 	   (pastix_int_t*)    lu->iparm,
229 	                      lu->dparm);
230 
231     ierr = PetscFree(lu->colptr);CHKERRQ(ierr);
232     ierr = PetscFree(lu->row);  CHKERRQ(ierr);
233     ierr = PetscFree(lu->val);  CHKERRQ(ierr);
234     ierr = PetscFree(lu->perm); CHKERRQ(ierr);
235     ierr = PetscFree(lu->invp); CHKERRQ(ierr);
236     ierr = MPI_Comm_free(&(lu->pastix_comm));CHKERRQ(ierr);
237   }
238   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
239   PetscFunctionReturn(0);
240 }
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "MatSolve_PaStiX"
244 /*
245   Gather right-hand-side.
246   Call for Solve step.
247   Scatter solution.
248  */
249 PetscErrorCode MatSolve_PaStiX(Mat A,Vec b,Vec x)
250 {
251   Mat_Pastix     *lu=(Mat_Pastix*)A->spptr;
252   PetscScalar    *array;
253   Vec             x_seq;
254   PetscErrorCode  ierr;
255 
256   PetscFunctionBegin;
257   lu->rhsnbr = 1;
258   x_seq = lu->b_seq;
259   if (lu->commSize > 1){
260     /* PaStiX only supports centralized rhs. Scatter b into a seqential rhs vector */
261     ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
262     ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
263     ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);
264   } else {  /* size == 1 */
265     ierr = VecCopy(b,x);CHKERRQ(ierr);
266     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
267   }
268   lu->rhs = array;
269   if (lu->commSize == 1){
270     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
271   } else {
272     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
273   }
274 
275   /* solve phase */
276   /*-------------*/
277   lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE;
278   lu->iparm[IPARM_END_TASK]   = API_TASK_REFINE;
279   lu->iparm[IPARM_RHS_MAKING] = API_RHS_B;
280 
281   pastix((pastix_data_t **)&(lu->pastix_data),
282 	 (MPI_Comm)         lu->pastix_comm,
283 	 (pastix_int_t)     lu->n,
284 	 (pastix_int_t*)    lu->colptr,
285 	 (pastix_int_t*)    lu->row,
286 	 (pastix_float_t*)  lu->val,
287 	 (pastix_int_t*)    lu->perm,
288 	 (pastix_int_t*)    lu->invp,
289 	 (pastix_float_t*)  lu->rhs,
290 	 (pastix_int_t)     lu->rhsnbr,
291 	 (pastix_int_t*)    lu->iparm,
292 	 (double*)          lu->dparm);
293 
294   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] );
295 
296   if (lu->commSize == 1){
297     ierr = VecRestoreArray(x,&(lu->rhs));CHKERRQ(ierr);
298   } else {
299     ierr = VecRestoreArray(x_seq,&(lu->rhs));CHKERRQ(ierr);
300   }
301 
302   if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */
303     ierr = VecScatterBegin(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
304     ierr = VecScatterEnd(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
305   }
306   PetscFunctionReturn(0);
307 }
308 
309 /*
310   Numeric factorisation using PaStiX solver.
311 
312  */
313 #undef __FUNCT__
314 #define __FUNCT__ "MatFactorNumeric_PaStiX"
315 PetscErrorCode MatFactorNumeric_PaStiX(Mat F,Mat A,const MatFactorInfo *info)
316 {
317   Mat_Pastix    *lu =(Mat_Pastix*)(F)->spptr;
318   Mat           *tseq;
319   PetscErrorCode ierr = 0;
320   PetscInt       icntl;
321   PetscInt       M=A->rmap->N;
322   PetscBool      valOnly,flg, isSym;
323   Mat            F_diag;
324   IS             is_iden;
325   Vec            b;
326   IS             isrow;
327   PetscBool      isSeqAIJ,isSeqSBAIJ,isMPIAIJ;
328 
329   PetscFunctionBegin;
330 
331   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
332   ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&isMPIAIJ);CHKERRQ(ierr);
333   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
334   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
335     (F)->ops->solve   = MatSolve_PaStiX;
336 
337     /* Initialize a PASTIX instance */
338     ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr);
339     ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank);         CHKERRQ(ierr);
340     ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize);         CHKERRQ(ierr);
341 
342     /* Set pastix options */
343     lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO;
344     lu->iparm[IPARM_START_TASK]       = API_TASK_INIT;
345     lu->iparm[IPARM_END_TASK]         = API_TASK_INIT;
346     lu->rhsnbr = 1;
347 
348     /* Call to set default pastix options */
349     pastix((pastix_data_t **)&(lu->pastix_data),
350 	   (MPI_Comm)         lu->pastix_comm,
351 	   (pastix_int_t)     lu->n,
352 	   (pastix_int_t*)    lu->colptr,
353 	   (pastix_int_t*)    lu->row,
354 	   (pastix_float_t*)  lu->val,
355 	   (pastix_int_t*)    lu->perm,
356 	   (pastix_int_t*)    lu->invp,
357 	   (pastix_float_t*)  lu->rhs,
358 	   (pastix_int_t)     lu->rhsnbr,
359 	   (pastix_int_t*)    lu->iparm,
360 	   (double*)          lu->dparm);
361 
362     ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr);
363 
364     icntl=-1;
365     lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT;
366     ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr);
367     if ((flg && icntl > 0) || PetscLogPrintInfo) {
368       lu->iparm[IPARM_VERBOSE] =  icntl;
369     }
370     icntl=-1;
371     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);
372     if ((flg && icntl > 0)) {
373       lu->iparm[IPARM_THREAD_NBR] = icntl;
374     }
375     PetscOptionsEnd();
376     valOnly = PETSC_FALSE;
377   }  else {
378     if (isSeqAIJ || isMPIAIJ)  {
379       ierr = PetscFree(lu->colptr);CHKERRQ(ierr);
380       ierr = PetscFree(lu->row);CHKERRQ(ierr);
381       ierr = PetscFree(lu->val);CHKERRQ(ierr);
382       valOnly = PETSC_FALSE;
383     } else valOnly = PETSC_TRUE;
384   }
385 
386   lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES;
387 
388   /* convert mpi A to seq mat A */
389   ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr);
390   ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr);
391   ierr = ISDestroy(isrow);CHKERRQ(ierr);
392 
393   ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val);CHKERRQ(ierr);
394   ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr);
395   ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr);
396 
397   if (!lu->perm) {
398     ierr = PetscMalloc((lu->n)*sizeof(PetscInt)   ,&(lu->perm));CHKERRQ(ierr);
399     ierr = PetscMalloc((lu->n)*sizeof(PetscInt)   ,&(lu->invp));CHKERRQ(ierr);
400   }
401 
402   if (isSym) {
403     /* On symmetric matrix, LLT */
404     lu->iparm[IPARM_SYM] = API_SYM_YES;
405     lu->iparm[IPARM_FACTORIZATION] = API_FACT_LDLT;
406   } else {
407     /* On unsymmetric matrix, LU */
408     lu->iparm[IPARM_SYM] = API_SYM_NO;
409     lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
410   }
411 
412   /*----------------*/
413   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
414     if (!(isSeqAIJ || isSeqSBAIJ)) {
415       /* PaStiX only supports centralized rhs. Create scatter scat_rhs for repeated use in MatSolve() */
416 	ierr = VecCreateSeq(PETSC_COMM_SELF,A->cmap->N,&lu->b_seq);CHKERRQ(ierr);
417 	ierr = ISCreateStride(PETSC_COMM_SELF,A->cmap->N,0,1,&is_iden);CHKERRQ(ierr);
418 	ierr = VecCreate(((PetscObject)A)->comm,&b);CHKERRQ(ierr);
419 	ierr = VecSetSizes(b,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr);
420 	ierr = VecSetFromOptions(b);CHKERRQ(ierr);
421 
422 	ierr = VecScatterCreate(b,is_iden,lu->b_seq,is_iden,&lu->scat_rhs);CHKERRQ(ierr);
423 	ierr = VecScatterCreate(lu->b_seq,is_iden,b,is_iden,&lu->scat_sol);CHKERRQ(ierr);
424 	ierr = ISDestroy(is_iden);CHKERRQ(ierr);
425 	ierr = VecDestroy(b);CHKERRQ(ierr);
426     }
427     lu->iparm[IPARM_START_TASK] = API_TASK_ORDERING;
428     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
429 
430     pastix((pastix_data_t **)&(lu->pastix_data),
431 	   (MPI_Comm)         lu->pastix_comm,
432 	   (pastix_int_t)     lu->n,
433 	   (pastix_int_t*)    lu->colptr,
434 	   (pastix_int_t*)    lu->row,
435 	   (pastix_float_t*)  lu->val,
436 	   (pastix_int_t*)    lu->perm,
437 	   (pastix_int_t*)    lu->invp,
438 	   (pastix_float_t*)  lu->rhs,
439 	   (pastix_int_t)     lu->rhsnbr,
440 	   (pastix_int_t*)    lu->iparm,
441 	   (double*)          lu->dparm);
442     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]);
443   } else {
444     lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT;
445     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
446     pastix((pastix_data_t **)&(lu->pastix_data),
447 	   (MPI_Comm)         lu->pastix_comm,
448 	   (pastix_int_t)     lu->n,
449 	   (pastix_int_t*)    lu->colptr,
450 	   (pastix_int_t*)    lu->row,
451 	   (pastix_float_t*)  lu->val,
452 	   (pastix_int_t*)    lu->perm,
453 	   (pastix_int_t*)    lu->invp,
454 	   (pastix_float_t*)  lu->rhs,
455 	   (pastix_int_t)     lu->rhsnbr,
456 	   (pastix_int_t*)    lu->iparm,
457 	   (double*)          lu->dparm);
458 
459     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]);
460   }
461 
462   if (lu->commSize > 1){
463     if ((F)->factortype == MAT_FACTOR_LU){
464       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
465     } else {
466       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
467     }
468     F_diag->assembled = PETSC_TRUE;
469   }
470   (F)->assembled     = PETSC_TRUE;
471   lu->matstruc       = SAME_NONZERO_PATTERN;
472   lu->CleanUpPastix  = PETSC_TRUE;
473   PetscFunctionReturn(0);
474 }
475 
476 /* Note the Petsc r and c permutations are ignored */
477 #undef __FUNCT__
478 #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX"
479 PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
480 {
481   Mat_Pastix      *lu = (Mat_Pastix*)F->spptr;
482 
483   PetscFunctionBegin;
484   lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
485   lu->iparm[IPARM_SYM]           = API_SYM_YES;
486   lu->matstruc                   = DIFFERENT_NONZERO_PATTERN;
487   F->ops->lufactornumeric        = MatFactorNumeric_PaStiX;
488   PetscFunctionReturn(0);
489 }
490 
491 
492 /* Note the Petsc r permutation is ignored */
493 #undef __FUNCT__
494 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX"
495 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info)
496 {
497   Mat_Pastix      *lu = (Mat_Pastix*)(F)->spptr;
498 
499   PetscFunctionBegin;
500   lu->iparm[IPARM_FACTORIZATION]  = API_FACT_LLT;
501   lu->iparm[IPARM_SYM]            = API_SYM_NO;
502   lu->matstruc                    = DIFFERENT_NONZERO_PATTERN;
503   (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX;
504   PetscFunctionReturn(0);
505 }
506 
507 #undef __FUNCT__
508 #define __FUNCT__ "MatView_PaStiX"
509 PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer)
510 {
511   PetscErrorCode    ierr;
512   PetscBool         iascii;
513   PetscViewerFormat format;
514 
515   PetscFunctionBegin;
516   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
517   if (iascii) {
518     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
519     if (format == PETSC_VIEWER_ASCII_INFO){
520       Mat_Pastix      *lu=(Mat_Pastix*)A->spptr;
521 
522       ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr);
523       ierr = PetscViewerASCIIPrintf(viewer,"  Matrix type :                      %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES)?"Symmetric":"Unsymmetric"));CHKERRQ(ierr);
524       ierr = PetscViewerASCIIPrintf(viewer,"  Level of printing (0,1,2):         %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr);
525       ierr = PetscViewerASCIIPrintf(viewer,"  Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr);
526       ierr = PetscPrintf(PETSC_COMM_SELF,"  Error :                        %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr);
527     }
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 
533 /*MC
534      MATSOLVERPASTIX  - A solver package providing direct solvers (LU) for distributed
535   and sequential matrices via the external package PaStiX.
536 
537   Use ./configure --download-pastix to have PETSc installed with PaStiX
538 
539   Options Database Keys:
540 + -mat_pastix_verbose   <0,1,2>   - print level
541 - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task.
542 
543   Level: beginner
544 
545 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
546 
547 M*/
548 
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "MatGetInfo_PaStiX"
552 PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info)
553 {
554     Mat_Pastix  *lu =(Mat_Pastix*)A->spptr;
555 
556     PetscFunctionBegin;
557     info->block_size        = 1.0;
558     info->nz_allocated      = lu->iparm[IPARM_NNZEROS];
559     info->nz_used           = lu->iparm[IPARM_NNZEROS];
560     info->nz_unneeded       = 0.0;
561     info->assemblies        = 0.0;
562     info->mallocs           = 0.0;
563     info->memory            = 0.0;
564     info->fill_ratio_given  = 0;
565     info->fill_ratio_needed = 0;
566     info->factor_mallocs    = 0;
567     PetscFunctionReturn(0);
568 }
569 
570 EXTERN_C_BEGIN
571 #undef __FUNCT__
572 #define __FUNCT__ "MatFactorGetSolverPackage_pastix"
573 PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type)
574 {
575   PetscFunctionBegin;
576   *type = MATSOLVERPASTIX;
577   PetscFunctionReturn(0);
578 }
579 EXTERN_C_END
580 
581 EXTERN_C_BEGIN
582 /*
583     The seq and mpi versions of this function are the same
584 */
585 #undef __FUNCT__
586 #define __FUNCT__ "MatGetFactor_seqaij_pastix"
587 PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F)
588 {
589   Mat            B;
590   PetscErrorCode ierr;
591   Mat_Pastix    *pastix;
592 
593   PetscFunctionBegin;
594   if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
595   /* Create the factorization matrix */
596   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
597   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
598   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
599   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
600 
601   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
602   B->ops->view             = MatView_PaStiX;
603   B->ops->getinfo          = MatGetInfo_PaStiX;
604   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix", MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
605   B->factortype            = MAT_FACTOR_LU;
606 
607   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
608   pastix->CleanUpPastix             = PETSC_FALSE;
609   pastix->isAIJ                     = PETSC_TRUE;
610   pastix->scat_rhs                  = PETSC_NULL;
611   pastix->scat_sol                  = PETSC_NULL;
612   pastix->MatDestroy                = B->ops->destroy;
613   B->ops->destroy                   = MatDestroy_Pastix;
614   B->spptr                          = (void*)pastix;
615 
616   *F = B;
617   PetscFunctionReturn(0);
618 }
619 EXTERN_C_END
620 
621 
622 EXTERN_C_BEGIN
623 #undef __FUNCT__
624 #define __FUNCT__ "MatGetFactor_mpiaij_pastix"
625 PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F)
626 {
627   Mat            B;
628   PetscErrorCode ierr;
629   Mat_Pastix    *pastix;
630 
631   PetscFunctionBegin;
632   if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
633   /* Create the factorization matrix */
634   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
635   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
636   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
637   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
638   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
639 
640   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
641   B->ops->view             = MatView_PaStiX;
642   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
643   B->factortype            = MAT_FACTOR_LU;
644 
645   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
646   pastix->CleanUpPastix             = PETSC_FALSE;
647   pastix->isAIJ                     = PETSC_TRUE;
648   pastix->scat_rhs                  = PETSC_NULL;
649   pastix->scat_sol                  = PETSC_NULL;
650   pastix->MatDestroy                = B->ops->destroy;
651   B->ops->destroy                  = MatDestroy_Pastix;
652   B->spptr                         = (void*)pastix;
653 
654   *F = B;
655   PetscFunctionReturn(0);
656 }
657 EXTERN_C_END
658 
659 EXTERN_C_BEGIN
660 #undef __FUNCT__
661 #define __FUNCT__ "MatGetFactor_seqsbaij_pastix"
662 PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F)
663 {
664   Mat            B;
665   PetscErrorCode ierr;
666   Mat_Pastix    *pastix;
667 
668   PetscFunctionBegin;
669   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
670   /* Create the factorization matrix */
671   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
672   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
673   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
674   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
675   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
676 
677   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
678   B->ops->view                   = MatView_PaStiX;
679   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
680   B->factortype                  = MAT_FACTOR_CHOLESKY;
681 
682   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
683   pastix->CleanUpPastix             = PETSC_FALSE;
684   pastix->isAIJ                     = PETSC_TRUE;
685   pastix->scat_rhs                  = PETSC_NULL;
686   pastix->scat_sol                  = PETSC_NULL;
687   pastix->MatDestroy                = B->ops->destroy;
688   B->ops->destroy                  = MatDestroy_Pastix;
689   B->spptr                         = (void*)pastix;
690 
691   *F = B;
692   PetscFunctionReturn(0);
693 }
694 EXTERN_C_END
695 
696 EXTERN_C_BEGIN
697 #undef __FUNCT__
698 #define __FUNCT__ "MatGetFactor_mpisbaij_pastix"
699 PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F)
700 {
701   Mat            B;
702   PetscErrorCode ierr;
703   Mat_Pastix    *pastix;
704 
705   PetscFunctionBegin;
706   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
707 
708   /* Create the factorization matrix */
709   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
710   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
711   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
712   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
713   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
714 
715   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
716   B->ops->view                   = MatView_PaStiX;
717   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
718   B->factortype                  = MAT_FACTOR_CHOLESKY;
719 
720   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
721   pastix->CleanUpPastix             = PETSC_FALSE;
722   pastix->isAIJ                     = PETSC_TRUE;
723   pastix->scat_rhs                  = PETSC_NULL;
724   pastix->scat_sol                  = PETSC_NULL;
725   pastix->MatDestroy                = B->ops->destroy;
726   B->ops->destroy                   = MatDestroy_Pastix;
727   B->spptr                          = (void*)pastix;
728 
729   *F = B;
730   PetscFunctionReturn(0);
731 }
732 EXTERN_C_END
733