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