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