xref: /petsc/src/mat/impls/aij/mpi/pastix/pastix.c (revision 07cda0dce865a0d4ab78ff237336a1a9351bce15)
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   PetscTruth     CleanUpPastix;            /* Boolean indicating if we call PaStiX clean step      */
39   VecScatter     scat_rhs;
40   VecScatter     scat_sol;
41   Vec            b_seq;
42   PetscTruth     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,PetscTruth 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   PetscTruth      isSBAIJ;
78   PetscTruth      isSeqSBAIJ;
79   PetscTruth      isMpiSBAIJ;
80   PetscTruth      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_ERR_PLIB,"overflow on k %D",k);
160 	  }
161       }
162     }
163   }
164   {
165     PetscScalar *tmpvalues;
166     PetscInt    *tmprows,*tmpcolptr;
167     tmpvalues = malloc(nnz*sizeof(PetscScalar)); if (!tmpvalues) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory");
168     tmprows   = malloc(nnz*sizeof(PetscInt));if (!tmprows) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory");
169     tmpcolptr = malloc((*n+1)*sizeof(PetscInt));if (!tmpcolptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MEM,"Unable to allocate memory");
170 
171     if (sizeof(PetscScalar) != sizeof(pastix_float_t)) {
172       SETERRQ2(PETSC_ERR_SUP,"sizeof(PetscScalar) %d != sizeof(pastix_float_t) %d",sizeof(PetscScalar),sizeof(pastix_float_t));
173     }
174     if (sizeof(PetscInt) != sizeof(pastix_int_t)) {
175       SETERRQ2(PETSC_ERR_SUP,"sizeof(PetscInt) %d != sizeof(pastix_int_t) %d",sizeof(PetscInt),sizeof(pastix_int_t));
176     }
177 
178     ierr = PetscMemcpy(tmpcolptr,*colptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr);
179     ierr = PetscMemcpy(tmprows,*row,nnz*sizeof(PetscInt));CHKERRQ(ierr);
180     ierr = PetscMemcpy(tmpvalues,*values,nnz*sizeof(PetscScalar));CHKERRQ(ierr);
181     ierr = PetscFree(*row);CHKERRQ(ierr);
182     ierr = PetscFree(*values);CHKERRQ(ierr);
183 
184     pastix_checkMatrix(MPI_COMM_WORLD,API_VERBOSE_NO,((isSym != 0) ? API_SYM_YES : API_SYM_NO),API_YES,*n,&tmpcolptr,&tmprows,&tmpvalues,NULL,1);
185 
186     ierr = PetscMemcpy(*colptr,tmpcolptr,(*n+1)*sizeof(PetscInt));CHKERRQ(ierr);
187     ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscInt),row);CHKERRQ(ierr);
188     ierr = PetscMemcpy(*row,tmprows,((*colptr)[*n]-1)*sizeof(PetscInt));CHKERRQ(ierr);
189     ierr = PetscMalloc(((*colptr)[*n]-1)*sizeof(PetscScalar),values);CHKERRQ(ierr);
190     ierr = PetscMemcpy(*values,tmpvalues,((*colptr)[*n]-1)*sizeof(PetscScalar));CHKERRQ(ierr);
191     free(tmpvalues);
192     free(tmprows);
193     free(tmpcolptr);
194   }
195   PetscFunctionReturn(0);
196 }
197 
198 
199 
200 #undef __FUNCT__
201 #define __FUNCT__ "MatDestroy_Pastix"
202 /*
203   Call clean step of PaStiX if lu->CleanUpPastix == true.
204   Free the CSC matrix.
205  */
206 PetscErrorCode MatDestroy_Pastix(Mat A)
207 {
208   Mat_Pastix      *lu=(Mat_Pastix*)A->spptr;
209   PetscErrorCode   ierr;
210   PetscMPIInt      size=lu->commSize;
211 
212   PetscFunctionBegin;
213   if (lu->CleanUpPastix) {
214     /* Terminate instance, deallocate memories */
215     if (size > 1){
216       ierr = VecScatterDestroy(lu->scat_rhs);CHKERRQ(ierr);
217       ierr = VecDestroy(lu->b_seq);CHKERRQ(ierr);
218       ierr = VecScatterDestroy(lu->scat_sol);CHKERRQ(ierr);
219     }
220 
221     lu->iparm[IPARM_START_TASK]=API_TASK_CLEAN;
222     lu->iparm[IPARM_END_TASK]=API_TASK_CLEAN;
223 
224     pastix((pastix_data_t **)&(lu->pastix_data),
225 	                      lu->pastix_comm,
226 	   (pastix_int_t)     lu->n,
227 	   (pastix_int_t*)    lu->colptr,
228 	   (pastix_int_t*)    lu->row,
229 	   (pastix_float_t*)  lu->val,
230 	   (pastix_int_t*)    lu->perm,
231 	   (pastix_int_t*)    lu->invp,
232 	   (pastix_float_t*)  lu->rhs,
233 	   (pastix_int_t)     lu->rhsnbr,
234 	   (pastix_int_t*)    lu->iparm,
235 	                      lu->dparm);
236 
237     ierr = PetscFree(lu->colptr);CHKERRQ(ierr);
238     ierr = PetscFree(lu->row);  CHKERRQ(ierr);
239     ierr = PetscFree(lu->val);  CHKERRQ(ierr);
240     ierr = PetscFree(lu->perm); CHKERRQ(ierr);
241     ierr = PetscFree(lu->invp); CHKERRQ(ierr);
242     ierr = MPI_Comm_free(&(lu->pastix_comm));CHKERRQ(ierr);
243   }
244   ierr = (lu->MatDestroy)(A);CHKERRQ(ierr);
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "MatSolve_PaStiX"
250 /*
251   Gather right-hand-side.
252   Call for Solve step.
253   Scatter solution.
254  */
255 PetscErrorCode MatSolve_PaStiX(Mat A,Vec b,Vec x)
256 {
257   Mat_Pastix     *lu=(Mat_Pastix*)A->spptr;
258   PetscScalar    *array;
259   Vec             x_seq;
260   PetscErrorCode  ierr;
261 
262   PetscFunctionBegin;
263   lu->rhsnbr = 1;
264   x_seq = lu->b_seq;
265   if (lu->commSize > 1){
266     /* PaStiX only supports centralized rhs. Scatter b into a seqential rhs vector */
267     ierr = VecScatterBegin(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
268     ierr = VecScatterEnd(lu->scat_rhs,b,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
269     ierr = VecGetArray(x_seq,&array);CHKERRQ(ierr);
270   } else {  /* size == 1 */
271     ierr = VecCopy(b,x);CHKERRQ(ierr);
272     ierr = VecGetArray(x,&array);CHKERRQ(ierr);
273   }
274   lu->rhs = array;
275   if (lu->commSize == 1){
276     ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
277   } else {
278     ierr = VecRestoreArray(x_seq,&array);CHKERRQ(ierr);
279   }
280 
281   /* solve phase */
282   /*-------------*/
283   lu->iparm[IPARM_START_TASK] = API_TASK_SOLVE;
284   lu->iparm[IPARM_END_TASK]   = API_TASK_REFINE;
285   lu->iparm[IPARM_RHS_MAKING] = API_RHS_B;
286 
287   pastix((pastix_data_t **)&(lu->pastix_data),
288 	 (MPI_Comm)         lu->pastix_comm,
289 	 (pastix_int_t)     lu->n,
290 	 (pastix_int_t*)    lu->colptr,
291 	 (pastix_int_t*)    lu->row,
292 	 (pastix_float_t*)  lu->val,
293 	 (pastix_int_t*)    lu->perm,
294 	 (pastix_int_t*)    lu->invp,
295 	 (pastix_float_t*)  lu->rhs,
296 	 (pastix_int_t)     lu->rhsnbr,
297 	 (pastix_int_t*)    lu->iparm,
298 	 (double*)          lu->dparm);
299 
300   if (lu->iparm[IPARM_ERROR_NUMBER] < 0) {
301     SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in solve phase: lu->iparm[IPARM_ERROR_NUMBER] = %d\n",lu->iparm[IPARM_ERROR_NUMBER] );
302   }
303 
304   if (lu->commSize == 1){
305     ierr = VecRestoreArray(x,&(lu->rhs));CHKERRQ(ierr);
306   } else {
307     ierr = VecRestoreArray(x_seq,&(lu->rhs));CHKERRQ(ierr);
308   }
309 
310   if (lu->commSize > 1) { /* convert PaStiX centralized solution to petsc mpi x */
311     ierr = VecScatterBegin(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
312     ierr = VecScatterEnd(lu->scat_sol,x_seq,x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
313   }
314   PetscFunctionReturn(0);
315 }
316 
317 /*
318   Numeric factorisation using PaStiX solver.
319 
320  */
321 #undef __FUNCT__
322 #define __FUNCT__ "MatFactorNumeric_PASTIX"
323 PetscErrorCode MatFactorNumeric_PaStiX(Mat F,Mat A,const MatFactorInfo *info)
324 {
325   Mat_Pastix    *lu =(Mat_Pastix*)(F)->spptr;
326   Mat           *tseq;
327   PetscErrorCode ierr = 0;
328   PetscInt       icntl;
329   PetscInt       M=A->rmap->N;
330   PetscTruth     valOnly,flg, isSym;
331   Mat            F_diag;
332   IS             is_iden;
333   Vec            b;
334   IS             isrow;
335   PetscTruth     isSeqAIJ,isSeqSBAIJ;
336 
337   PetscFunctionBegin;
338 
339   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&isSeqAIJ);CHKERRQ(ierr);
340   ierr = PetscTypeCompare((PetscObject)A,MATSEQSBAIJ,&isSeqSBAIJ);CHKERRQ(ierr);
341   if (lu->matstruc == DIFFERENT_NONZERO_PATTERN){
342     (F)->ops->solve   = MatSolve_PaStiX;
343 
344     /* Initialize a PASTIX instance */
345     ierr = MPI_Comm_dup(((PetscObject)A)->comm,&(lu->pastix_comm));CHKERRQ(ierr);
346     ierr = MPI_Comm_rank(lu->pastix_comm, &lu->commRank);         CHKERRQ(ierr);
347     ierr = MPI_Comm_size(lu->pastix_comm, &lu->commSize);         CHKERRQ(ierr);
348 
349     /* Set pastix options */
350     lu->iparm[IPARM_MODIFY_PARAMETER] = API_NO;
351     lu->iparm[IPARM_START_TASK]       = API_TASK_INIT;
352     lu->iparm[IPARM_END_TASK]         = API_TASK_INIT;
353     lu->rhsnbr = 1;
354 
355     /* Call to set default pastix options */
356     pastix((pastix_data_t **)&(lu->pastix_data),
357 	   (MPI_Comm)         lu->pastix_comm,
358 	   (pastix_int_t)     lu->n,
359 	   (pastix_int_t*)    lu->colptr,
360 	   (pastix_int_t*)    lu->row,
361 	   (pastix_float_t*)  lu->val,
362 	   (pastix_int_t*)    lu->perm,
363 	   (pastix_int_t*)    lu->invp,
364 	   (pastix_float_t*)  lu->rhs,
365 	   (pastix_int_t)     lu->rhsnbr,
366 	   (pastix_int_t*)    lu->iparm,
367 	   (double*)          lu->dparm);
368 
369     ierr = PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"PaStiX Options","Mat");CHKERRQ(ierr);
370 
371     icntl=-1;
372     lu->iparm[IPARM_VERBOSE] = API_VERBOSE_NOT;
373     ierr = PetscOptionsInt("-mat_pastix_verbose","iparm[IPARM_VERBOSE] : level of printing (0 to 2)","None",lu->iparm[IPARM_VERBOSE],&icntl,&flg);CHKERRQ(ierr);
374     if ((flg && icntl > 0) || PetscLogPrintInfo) {
375       lu->iparm[IPARM_VERBOSE] =  icntl;
376     }
377     icntl=-1;
378     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);
379     if ((flg && icntl > 0)) {
380       lu->iparm[IPARM_THREAD_NBR] = icntl;
381     }
382     PetscOptionsEnd();
383     valOnly = PETSC_FALSE;
384   }  else {
385     valOnly = PETSC_TRUE;
386   }
387 
388   lu->iparm[IPARM_MATRIX_VERIFICATION] = API_YES;
389 
390   /* convert mpi A to seq mat A */
391   ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr);
392   ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr);
393   ierr = ISDestroy(isrow);CHKERRQ(ierr);
394 
395   ierr = MatConvertToCSC(*tseq,valOnly, &lu->n, &lu->colptr, &lu->row, &lu->val);CHKERRQ(ierr);
396   ierr = MatIsSymmetric(*tseq,0.0,&isSym);CHKERRQ(ierr);
397   ierr = MatDestroyMatrices(1,&tseq);CHKERRQ(ierr);
398 
399   ierr = PetscMalloc((lu->n)*sizeof(PetscInt)   ,&(lu->perm));CHKERRQ(ierr);
400   ierr = PetscMalloc((lu->n)*sizeof(PetscInt)   ,&(lu->invp));CHKERRQ(ierr);
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) {
443       SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%d\n",lu->iparm[IPARM_ERROR_NUMBER]);
444     }
445   } else {
446     lu->iparm[IPARM_START_TASK] = API_TASK_NUMFACT;
447     lu->iparm[IPARM_END_TASK]   = API_TASK_NUMFACT;
448     pastix((pastix_data_t **)&(lu->pastix_data),
449 	   (MPI_Comm)         lu->pastix_comm,
450 	   (pastix_int_t)     lu->n,
451 	   (pastix_int_t*)    lu->colptr,
452 	   (pastix_int_t*)    lu->row,
453 	   (pastix_float_t*)  lu->val,
454 	   (pastix_int_t*)    lu->perm,
455 	   (pastix_int_t*)    lu->invp,
456 	   (pastix_float_t*)  lu->rhs,
457 	   (pastix_int_t)     lu->rhsnbr,
458 	   (pastix_int_t*)    lu->iparm,
459 	   (double*)          lu->dparm);
460 
461     if (lu->iparm[IPARM_ERROR_NUMBER] < 0) {
462       SETERRQ1(PETSC_ERR_LIB,"Error reported by PaStiX in analysis phase: iparm(IPARM_ERROR_NUMBER)=%d\n",lu->iparm[IPARM_ERROR_NUMBER]);
463     }
464   }
465 
466   if (lu->commSize > 1){
467     if ((F)->factor == MAT_FACTOR_LU){
468       F_diag = ((Mat_MPIAIJ *)(F)->data)->A;
469     } else {
470       F_diag = ((Mat_MPISBAIJ *)(F)->data)->A;
471     }
472     F_diag->assembled = PETSC_TRUE;
473   }
474   (F)->assembled     = PETSC_TRUE;
475   lu->matstruc       = SAME_NONZERO_PATTERN;
476   lu->CleanUpPastix  = PETSC_TRUE;
477   PetscFunctionReturn(0);
478 }
479 
480 /* Note the Petsc r and c permutations are ignored */
481 #undef __FUNCT__
482 #define __FUNCT__ "MatLUFactorSymbolic_AIJPASTIX"
483 PetscErrorCode MatLUFactorSymbolic_AIJPASTIX(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
484 {
485   Mat_Pastix      *lu = (Mat_Pastix*)F->spptr;
486 
487   PetscFunctionBegin;
488   lu->iparm[IPARM_FACTORIZATION] = API_FACT_LU;
489   lu->iparm[IPARM_SYM]           = API_SYM_YES;
490   lu->matstruc                   = DIFFERENT_NONZERO_PATTERN;
491   F->ops->lufactornumeric        = MatFactorNumeric_PaStiX;
492   PetscFunctionReturn(0);
493 }
494 
495 
496 /* Note the Petsc r permutation is ignored */
497 #undef __FUNCT__
498 #define __FUNCT__ "MatCholeskyFactorSymbolic_SBAIJPASTIX"
499 PetscErrorCode MatCholeskyFactorSymbolic_SBAIJPASTIX(Mat F,Mat A,IS r,const MatFactorInfo *info)
500 {
501   Mat_Pastix      *lu = (Mat_Pastix*)(F)->spptr;
502 
503   PetscFunctionBegin;
504   lu->iparm[IPARM_FACTORIZATION]  = API_FACT_LLT;
505   lu->iparm[IPARM_SYM]            = API_SYM_NO;
506   lu->matstruc                    = DIFFERENT_NONZERO_PATTERN;
507   (F)->ops->choleskyfactornumeric = MatFactorNumeric_PaStiX;
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "MatView_PaStiX"
513 PetscErrorCode MatView_PaStiX(Mat A,PetscViewer viewer)
514 {
515   PetscErrorCode    ierr;
516   PetscTruth        iascii;
517   PetscViewerFormat format;
518 
519   PetscFunctionBegin;
520   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
521   if (iascii) {
522     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
523     if (format == PETSC_VIEWER_ASCII_INFO){
524       Mat_Pastix      *lu=(Mat_Pastix*)A->spptr;
525 
526       ierr = PetscViewerASCIIPrintf(viewer,"PaStiX run parameters:\n");CHKERRQ(ierr);
527       ierr = PetscViewerASCIIPrintf(viewer,"  Matrix type :                      %s \n",((lu->iparm[IPARM_SYM] == API_SYM_YES)?"Symmetric":"Unsymmetric"));CHKERRQ(ierr);
528       ierr = PetscViewerASCIIPrintf(viewer,"  Level of printing (0,1,2):         %d \n",lu->iparm[IPARM_VERBOSE]);CHKERRQ(ierr);
529       ierr = PetscViewerASCIIPrintf(viewer,"  Number of refinements iterations : %d \n",lu->iparm[IPARM_NBITER]);CHKERRQ(ierr);
530       ierr = PetscPrintf(PETSC_COMM_SELF,"  Error :                        %g \n",lu->dparm[DPARM_RELATIVE_ERROR]);CHKERRQ(ierr);
531     }
532   }
533   PetscFunctionReturn(0);
534 }
535 
536 
537 /*MC
538      MAT_SOLVER_PASTIX  - A solver package providing direct solvers (LU) for distributed
539   and sequential matrices via the external package PaStiX.
540 
541   Use config/configure.py --download-pastix to have PETSc installed with PaStiX
542 
543   Options Database Keys:
544 + -mat_pastix_verbose   <0,1,2>   - print level
545 - -mat_pastix_threadnbr <integer> - Set the thread number by MPI task.
546 
547   Level: beginner
548 
549 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
550 
551 M*/
552 
553 
554 #undef __FUNCT__
555 #define __FUNCT__ "MatGetInfo_PaStiX"
556 PetscErrorCode MatGetInfo_PaStiX(Mat A,MatInfoType flag,MatInfo *info)
557 {
558     Mat_Pastix  *lu =(Mat_Pastix*)A->spptr;
559 
560     PetscFunctionBegin;
561     info->block_size        = 1.0;
562     info->nz_allocated      = lu->iparm[IPARM_NNZEROS];
563     info->nz_used           = lu->iparm[IPARM_NNZEROS];
564     info->nz_unneeded       = 0.0;
565     info->assemblies        = 0.0;
566     info->mallocs           = 0.0;
567     info->memory            = 0.0;
568     info->fill_ratio_given  = 0;
569     info->fill_ratio_needed = 0;
570     info->factor_mallocs    = 0;
571     PetscFunctionReturn(0);
572 }
573 
574 EXTERN_C_BEGIN
575 #undef __FUNCT__
576 #define __FUNCT__ "MatFactorGetSolverPackage_pastix"
577 PetscErrorCode MatFactorGetSolverPackage_pastix(Mat A,const MatSolverPackage *type)
578 {
579   PetscFunctionBegin;
580   *type = MAT_SOLVER_PASTIX;
581   PetscFunctionReturn(0);
582 }
583 EXTERN_C_END
584 
585 EXTERN_C_BEGIN
586 /*
587     The seq and mpi versions of this function are the same
588 */
589 #undef __FUNCT__
590 #define __FUNCT__ "MatGetFactor_seqaij_pastix"
591 PetscErrorCode MatGetFactor_seqaij_pastix(Mat A,MatFactorType ftype,Mat *F)
592 {
593   Mat            B;
594   PetscErrorCode ierr;
595   Mat_Pastix    *pastix;
596 
597   PetscFunctionBegin;
598   if (ftype != MAT_FACTOR_LU) {
599     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
600   }
601   /* Create the factorization matrix */
602   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
603   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
604   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
605   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
606 
607   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
608   B->ops->view             = MatView_PaStiX;
609   B->ops->getinfo          = MatGetInfo_PaStiX;
610   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix", MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
611   B->factor                = MAT_FACTOR_LU;
612 
613   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
614   pastix->CleanUpPastix             = PETSC_FALSE;
615   pastix->isAIJ                     = PETSC_TRUE;
616   pastix->scat_rhs                  = PETSC_NULL;
617   pastix->scat_sol                  = PETSC_NULL;
618   pastix->MatDestroy                = B->ops->destroy;
619   B->ops->destroy                   = MatDestroy_Pastix;
620   B->spptr                          = (void*)pastix;
621 
622   *F = B;
623   PetscFunctionReturn(0);
624 }
625 EXTERN_C_END
626 
627 
628 EXTERN_C_BEGIN
629 #undef __FUNCT__
630 #define __FUNCT__ "MatGetFactor_mpiaij_pastix"
631 PetscErrorCode MatGetFactor_mpiaij_pastix(Mat A,MatFactorType ftype,Mat *F)
632 {
633   Mat            B;
634   PetscErrorCode ierr;
635   Mat_Pastix    *pastix;
636 
637   PetscFunctionBegin;
638   if (ftype != MAT_FACTOR_LU) SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc AIJ matrices with PaStiX Cholesky, use SBAIJ matrix");
639   /* Create the factorization matrix */
640   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
641   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
642   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
643   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
644   ierr = MatMPIAIJSetPreallocation(B,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
645 
646   B->ops->lufactorsymbolic = MatLUFactorSymbolic_AIJPASTIX;
647   B->ops->view             = MatView_PaStiX;
648   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
649   B->factor                = MAT_FACTOR_LU;
650 
651   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
652   pastix->CleanUpPastix             = PETSC_FALSE;
653   pastix->isAIJ                     = PETSC_TRUE;
654   pastix->scat_rhs                  = PETSC_NULL;
655   pastix->scat_sol                  = PETSC_NULL;
656   pastix->MatDestroy                = B->ops->destroy;
657   B->ops->destroy                  = MatDestroy_Pastix;
658   B->spptr                         = (void*)pastix;
659 
660   *F = B;
661   PetscFunctionReturn(0);
662 }
663 EXTERN_C_END
664 
665 EXTERN_C_BEGIN
666 #undef __FUNCT__
667 #define __FUNCT__ "MatGetFactor_seqsbaij_pastix"
668 PetscErrorCode MatGetFactor_seqsbaij_pastix(Mat A,MatFactorType ftype,Mat *F)
669 {
670   Mat            B;
671   PetscErrorCode ierr;
672   Mat_Pastix    *pastix;
673 
674   PetscFunctionBegin;
675   if (ftype != MAT_FACTOR_CHOLESKY) {
676     SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
677   }
678   /* Create the factorization matrix */
679   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
680   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
681   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
682   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
683   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
684 
685   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
686   B->ops->view                   = MatView_PaStiX;
687   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
688 
689   B->factor                      = MAT_FACTOR_CHOLESKY;
690 
691   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
692   pastix->CleanUpPastix             = PETSC_FALSE;
693   pastix->isAIJ                     = PETSC_TRUE;
694   pastix->scat_rhs                  = PETSC_NULL;
695   pastix->scat_sol                  = PETSC_NULL;
696   pastix->MatDestroy                = B->ops->destroy;
697   B->ops->destroy                  = MatDestroy_Pastix;
698   B->spptr                         = (void*)pastix;
699 
700   *F = B;
701   PetscFunctionReturn(0);
702 }
703 EXTERN_C_END
704 
705 EXTERN_C_BEGIN
706 #undef __FUNCT__
707 #define __FUNCT__ "MatGetFactor_mpisbaij_pastix"
708 PetscErrorCode MatGetFactor_mpisbaij_pastix(Mat A,MatFactorType ftype,Mat *F)
709 {
710   Mat            B;
711   PetscErrorCode ierr;
712   Mat_Pastix    *pastix;
713 
714   PetscFunctionBegin;
715   if (ftype != MAT_FACTOR_CHOLESKY) SETERRQ(PETSC_ERR_SUP,"Cannot use PETSc SBAIJ matrices with PaStiX LU, use AIJ matrix");
716 
717   /* Create the factorization matrix */
718   ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr);
719   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
720   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
721   ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr);
722   ierr = MatMPISBAIJSetPreallocation(B,1,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr);
723 
724   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SBAIJPASTIX;
725   B->ops->view                   = MatView_PaStiX;
726   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_pastix",MatFactorGetSolverPackage_pastix);CHKERRQ(ierr);
727   B->factor                      = MAT_FACTOR_CHOLESKY;
728 
729   ierr = PetscNewLog(B,Mat_Pastix,&pastix);CHKERRQ(ierr);
730   pastix->CleanUpPastix             = PETSC_FALSE;
731   pastix->isAIJ                     = PETSC_TRUE;
732   pastix->scat_rhs                  = PETSC_NULL;
733   pastix->scat_sol                  = PETSC_NULL;
734   pastix->MatDestroy                = B->ops->destroy;
735   B->ops->destroy                   = MatDestroy_Pastix;
736   B->spptr                          = (void*)pastix;
737 
738   *F = B;
739   PetscFunctionReturn(0);
740 }
741 EXTERN_C_END
742