xref: /petsc/src/mat/impls/aij/mpi/superlu_dist/superlu_dist.c (revision a7fac7c2b47dbcd8ece0cc2dfdfe0e63be1bb7b5)
1 
2 /*
3         Provides an interface to the SuperLU_DIST_2.2 sparse solver
4 */
5 
6 #include <../src/mat/impls/aij/seq/aij.h>
7 #include <../src/mat/impls/aij/mpi/mpiaij.h>
8 #if defined(PETSC_HAVE_STDLIB_H) /* This is to get around weird problem with SuperLU on cray */
9 #include <stdlib.h>
10 #endif
11 
12 #if defined(PETSC_USE_64BIT_INDICES)
13 /* ugly SuperLU_Dist variable telling it to use long long int */
14 #define _LONGINT
15 #endif
16 
17 EXTERN_C_BEGIN
18 #if defined(PETSC_USE_COMPLEX)
19 #include <superlu_zdefs.h>
20 #else
21 #include <superlu_ddefs.h>
22 #endif
23 EXTERN_C_END
24 
25 /*
26     GLOBAL - The sparse matrix and right hand side are all stored initially on process 0. Should be called centralized
27     DISTRIBUTED - The sparse matrix and right hand size are initially stored across the entire MPI communicator.
28 */
29 typedef enum {GLOBAL,DISTRIBUTED} SuperLU_MatInputMode;
30 const char *SuperLU_MatInputModes[] = {"GLOBAL","DISTRIBUTED","SuperLU_MatInputMode","PETSC_",0};
31 
32 typedef struct {
33   int_t                nprow,npcol,*row,*col;
34   gridinfo_t           grid;
35   superlu_options_t    options;
36   SuperMatrix          A_sup;
37   ScalePermstruct_t    ScalePermstruct;
38   LUstruct_t           LUstruct;
39   int                  StatPrint;
40   SuperLU_MatInputMode MatInputMode;
41   SOLVEstruct_t        SOLVEstruct;
42   fact_t               FactPattern;
43   MPI_Comm             comm_superlu;
44 #if defined(PETSC_USE_COMPLEX)
45   doublecomplex        *val;
46 #else
47   double               *val;
48 #endif
49   PetscBool            matsolve_iscalled,matmatsolve_iscalled;
50   PetscBool            CleanUpSuperLU_Dist;  /* Flag to clean up (non-global) SuperLU objects during Destroy */
51 } Mat_SuperLU_DIST;
52 
53 extern PetscErrorCode MatFactorInfo_SuperLU_DIST(Mat,PetscViewer);
54 extern PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat,Mat,const MatFactorInfo*);
55 extern PetscErrorCode MatDestroy_SuperLU_DIST(Mat);
56 extern PetscErrorCode MatView_SuperLU_DIST(Mat,PetscViewer);
57 extern PetscErrorCode MatSolve_SuperLU_DIST(Mat,Vec,Vec);
58 extern PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat,Mat,IS,IS,const MatFactorInfo*);
59 extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
60 
61 #undef __FUNCT__
62 #define __FUNCT__ "MatGetDiagonal_SuperLU_DIST"
63 PetscErrorCode MatGetDiagonal_SuperLU_DIST(Mat A,Vec v)
64 {
65   PetscFunctionBegin;
66   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Mat type: SuperLU_DIST factor");
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "MatDestroy_SuperLU_DIST"
72 PetscErrorCode MatDestroy_SuperLU_DIST(Mat A)
73 {
74   PetscErrorCode   ierr;
75   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->spptr;
76   PetscBool        flg;
77 
78   PetscFunctionBegin;
79   if (lu && lu->CleanUpSuperLU_Dist) {
80     /* Deallocate SuperLU_DIST storage */
81     if (lu->MatInputMode == GLOBAL) {
82       PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup));
83     } else {
84       PetscStackCall("SuperLU_DIST:Destroy_CompRowLoc_Matrix_dist",Destroy_CompRowLoc_Matrix_dist(&lu->A_sup));
85       if (lu->options.SolveInitialized) {
86 #if defined(PETSC_USE_COMPLEX)
87         PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
88 #else
89         PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
90 #endif
91       }
92     }
93     PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(A->cmap->N, &lu->grid, &lu->LUstruct));
94     PetscStackCall("SuperLU_DIST:ScalePermstructFree",ScalePermstructFree(&lu->ScalePermstruct));
95     PetscStackCall("SuperLU_DIST:LUstructFree",LUstructFree(&lu->LUstruct));
96 
97     /* Release the SuperLU_DIST process grid. */
98     PetscStackCall("SuperLU_DIST:superlu_gridexit",superlu_gridexit(&lu->grid));
99     ierr = MPI_Comm_free(&(lu->comm_superlu));CHKERRQ(ierr);
100   }
101   ierr = PetscFree(A->spptr);CHKERRQ(ierr);
102 
103   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
104   if (flg) {
105     ierr = MatDestroy_SeqAIJ(A);CHKERRQ(ierr);
106   } else {
107     ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
108   }
109   PetscFunctionReturn(0);
110 }
111 
112 #undef __FUNCT__
113 #define __FUNCT__ "MatSolve_SuperLU_DIST"
114 PetscErrorCode MatSolve_SuperLU_DIST(Mat A,Vec b_mpi,Vec x)
115 {
116   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->spptr;
117   PetscErrorCode   ierr;
118   PetscMPIInt      size;
119   PetscInt         m=A->rmap->n,M=A->rmap->N,N=A->cmap->N;
120   SuperLUStat_t    stat;
121   double           berr[1];
122   PetscScalar      *bptr;
123   PetscInt         nrhs=1;
124   Vec              x_seq;
125   IS               iden;
126   VecScatter       scat;
127   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
128   static PetscBool cite = PETSC_FALSE;
129 
130   PetscFunctionBegin;
131   ierr = PetscCitationsRegister("@article{lidemmel03,\n  author = {Xiaoye S. Li and James W. Demmel},\n  title = {{SuperLU_DIST}: A Scalable Distributed-Memory Sparse Direct\n           Solver for Unsymmetric Linear Systems},\n  journal = {ACM Trans. Mathematical Software},\n  volume = {29},\n  number = {2},\n  pages = {110-140},\n  year = 2003\n}\n",&cite);CHKERRQ(ierr);
132   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
133   if (size > 1 && lu->MatInputMode == GLOBAL) {
134     /* global mat input, convert b to x_seq */
135     ierr = VecCreateSeq(PETSC_COMM_SELF,N,&x_seq);CHKERRQ(ierr);
136     ierr = ISCreateStride(PETSC_COMM_SELF,N,0,1,&iden);CHKERRQ(ierr);
137     ierr = VecScatterCreate(b_mpi,iden,x_seq,iden,&scat);CHKERRQ(ierr);
138     ierr = ISDestroy(&iden);CHKERRQ(ierr);
139 
140     ierr = VecScatterBegin(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
141     ierr = VecScatterEnd(scat,b_mpi,x_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
142     ierr = VecGetArray(x_seq,&bptr);CHKERRQ(ierr);
143   } else { /* size==1 || distributed mat input */
144     if (lu->options.SolveInitialized && !lu->matsolve_iscalled) {
145       /* see comments in MatMatSolve() */
146 #if defined(PETSC_USE_COMPLEX)
147       PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
148 #else
149       PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
150 #endif
151       lu->options.SolveInitialized = NO;
152     }
153     ierr = VecCopy(b_mpi,x);CHKERRQ(ierr);
154     ierr = VecGetArray(x,&bptr);CHKERRQ(ierr);
155   }
156 
157   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED");
158 
159   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
160   if (lu->MatInputMode == GLOBAL) {
161 #if defined(PETSC_USE_COMPLEX)
162     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info));
163 #else
164     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr,M,nrhs,&lu->grid,&lu->LUstruct,berr,&stat,&info));
165 #endif
166   } else { /* distributed mat input */
167 #if defined(PETSC_USE_COMPLEX)
168     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
169 #else
170     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
171 #endif
172   }
173   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
174 
175   if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid);      /* Print the statistics. */
176   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
177 
178   if (size > 1 && lu->MatInputMode == GLOBAL) {
179     /* convert seq x to mpi x */
180     ierr = VecRestoreArray(x_seq,&bptr);CHKERRQ(ierr);
181     ierr = VecScatterBegin(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
182     ierr = VecScatterEnd(scat,x_seq,x,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
183     ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
184     ierr = VecDestroy(&x_seq);CHKERRQ(ierr);
185   } else {
186     ierr = VecRestoreArray(x,&bptr);CHKERRQ(ierr);
187 
188     lu->matsolve_iscalled    = PETSC_TRUE;
189     lu->matmatsolve_iscalled = PETSC_FALSE;
190   }
191   PetscFunctionReturn(0);
192 }
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "MatMatSolve_SuperLU_DIST"
196 PetscErrorCode MatMatSolve_SuperLU_DIST(Mat A,Mat B_mpi,Mat X)
197 {
198   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)A->spptr;
199   PetscErrorCode   ierr;
200   PetscMPIInt      size;
201   PetscInt         M=A->rmap->N,m=A->rmap->n,nrhs;
202   SuperLUStat_t    stat;
203   double           berr[1];
204   PetscScalar      *bptr;
205   int              info; /* SuperLU_Dist info code is ALWAYS an int, even with long long indices */
206   PetscBool        flg;
207 
208   PetscFunctionBegin;
209   if (lu->options.Fact != FACTORED) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"SuperLU_DIST options.Fact mush equal FACTORED");
210   ierr = PetscObjectTypeCompareAny((PetscObject)B_mpi,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
211   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
212   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
213   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
214 
215   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
216   if (size > 1 && lu->MatInputMode == GLOBAL) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatInputMode=GLOBAL for nproc %d>1 is not supported",size);
217   /* size==1 or distributed mat input */
218   if (lu->options.SolveInitialized && !lu->matmatsolve_iscalled) {
219     /* communication pattern of SOLVEstruct is unlikely created for matmatsolve,
220        thus destroy it and create a new SOLVEstruct.
221        Otherwise it may result in memory corruption or incorrect solution
222        See src/mat/examples/tests/ex125.c */
223 #if defined(PETSC_USE_COMPLEX)
224     PetscStackCall("SuperLU_DIST:zSolveFinalize",zSolveFinalize(&lu->options, &lu->SOLVEstruct));
225 #else
226     PetscStackCall("SuperLU_DIST:dSolveFinalize",dSolveFinalize(&lu->options, &lu->SOLVEstruct));
227 #endif
228     lu->options.SolveInitialized = NO;
229   }
230   ierr = MatCopy(B_mpi,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
231 
232   ierr = MatGetSize(B_mpi,NULL,&nrhs);CHKERRQ(ierr);
233 
234   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));        /* Initialize the statistics variables. */
235   ierr = MatDenseGetArray(X,&bptr);CHKERRQ(ierr);
236   if (lu->MatInputMode == GLOBAL) { /* size == 1 */
237 #if defined(PETSC_USE_COMPLEX)
238     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,(doublecomplex*)bptr, M, nrhs,&lu->grid, &lu->LUstruct, berr, &stat, &info));
239 #else
240     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct,bptr, M, nrhs, &lu->grid, &lu->LUstruct, berr, &stat, &info));
241 #endif
242   } else { /* distributed mat input */
243 #if defined(PETSC_USE_COMPLEX)
244     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,(doublecomplex*)bptr,m,nrhs,&lu->grid, &lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
245 #else
246     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options,&lu->A_sup,&lu->ScalePermstruct,bptr,m,nrhs,&lu->grid,&lu->LUstruct,&lu->SOLVEstruct,berr,&stat,&info));
247 #endif
248   }
249   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",info);
250   ierr = MatDenseRestoreArray(X,&bptr);CHKERRQ(ierr);
251 
252   if (lu->options.PrintStat) PStatPrint(&lu->options, &stat, &lu->grid); /* Print the statistics. */
253   PetscStackCall("SuperLU_DIST:PStatFree",PStatFree(&stat));
254   lu->matsolve_iscalled    = PETSC_FALSE;
255   lu->matmatsolve_iscalled = PETSC_TRUE;
256   PetscFunctionReturn(0);
257 }
258 
259 
260 #undef __FUNCT__
261 #define __FUNCT__ "MatLUFactorNumeric_SuperLU_DIST"
262 PetscErrorCode MatLUFactorNumeric_SuperLU_DIST(Mat F,Mat A,const MatFactorInfo *info)
263 {
264   Mat              *tseq,A_seq = NULL;
265   Mat_SeqAIJ       *aa,*bb;
266   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)(F)->spptr;
267   PetscErrorCode   ierr;
268   PetscInt         M=A->rmap->N,N=A->cmap->N,i,*ai,*aj,*bi,*bj,nz,rstart,*garray,
269                    m=A->rmap->n, colA_start,j,jcol,jB,countA,countB,*bjj,*ajj;
270   int              sinfo;   /* SuperLU_Dist info flag is always an int even with long long indices */
271   PetscMPIInt      size;
272   SuperLUStat_t    stat;
273   double           *berr=0;
274   IS               isrow;
275   Mat              F_diag=NULL;
276 #if defined(PETSC_USE_COMPLEX)
277   doublecomplex    *av, *bv;
278 #else
279   double           *av, *bv;
280 #endif
281 
282   PetscFunctionBegin;
283   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
284 
285   if (lu->MatInputMode == GLOBAL) { /* global mat input */
286     if (size > 1) { /* convert mpi A to seq mat A */
287       ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);CHKERRQ(ierr);
288       ierr = MatGetSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&tseq);CHKERRQ(ierr);
289       ierr = ISDestroy(&isrow);CHKERRQ(ierr);
290 
291       A_seq = *tseq;
292       ierr  = PetscFree(tseq);CHKERRQ(ierr);
293       aa    = (Mat_SeqAIJ*)A_seq->data;
294     } else {
295       PetscBool flg;
296       ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
297       if (flg) {
298         Mat_MPIAIJ *At = (Mat_MPIAIJ*)A->data;
299         A = At->A;
300       }
301       aa =  (Mat_SeqAIJ*)A->data;
302     }
303 
304     /* Convert Petsc NR matrix to SuperLU_DIST NC.
305        Note: memories of lu->val, col and row are allocated by CompRow_to_CompCol_dist()! */
306     if (lu->options.Fact != DOFACT) {/* successive numeric factorization, sparsity pattern is reused. */
307       PetscStackCall("SuperLU_DIST:Destroy_CompCol_Matrix_dist",Destroy_CompCol_Matrix_dist(&lu->A_sup));
308       if (lu->FactPattern == SamePattern_SameRowPerm) {
309         lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
310       } else { /* lu->FactPattern == SamePattern */
311         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct));
312         lu->options.Fact = SamePattern;
313       }
314     }
315 #if defined(PETSC_USE_COMPLEX)
316     PetscStackCall("SuperLU_DIST:zCompRow_to_CompCol_dist",zCompRow_to_CompCol_dist(M,N,aa->nz,(doublecomplex*)aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val,&lu->col, &lu->row));
317 #else
318     PetscStackCall("SuperLU_DIST:dCompRow_to_CompCol_dist",dCompRow_to_CompCol_dist(M,N,aa->nz,aa->a,(int_t*)aa->j,(int_t*)aa->i,&lu->val, &lu->col, &lu->row));
319 #endif
320 
321     /* Create compressed column matrix A_sup. */
322 #if defined(PETSC_USE_COMPLEX)
323     PetscStackCall("SuperLU_DIST:zCreate_CompCol_Matrix_dist",zCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_Z, SLU_GE));
324 #else
325     PetscStackCall("SuperLU_DIST:dCreate_CompCol_Matrix_dist",dCreate_CompCol_Matrix_dist(&lu->A_sup, M, N, aa->nz, lu->val, lu->col, lu->row, SLU_NC, SLU_D, SLU_GE));
326 #endif
327   } else { /* distributed mat input */
328     Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
329     aa=(Mat_SeqAIJ*)(mat->A)->data;
330     bb=(Mat_SeqAIJ*)(mat->B)->data;
331     ai=aa->i; aj=aa->j;
332     bi=bb->i; bj=bb->j;
333 #if defined(PETSC_USE_COMPLEX)
334     av=(doublecomplex*)aa->a;
335     bv=(doublecomplex*)bb->a;
336 #else
337     av=aa->a;
338     bv=bb->a;
339 #endif
340     rstart = A->rmap->rstart;
341     nz     = aa->nz + bb->nz;
342     garray = mat->garray;
343 
344     if (lu->options.Fact == DOFACT) { /* first numeric factorization */
345 #if defined(PETSC_USE_COMPLEX)
346       PetscStackCall("SuperLU_DIST:zallocateA_dist",zallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row));
347 #else
348       PetscStackCall("SuperLU_DIST:dallocateA_dist",dallocateA_dist(m, nz, &lu->val, &lu->col, &lu->row));
349 #endif
350     } else { /* successive numeric factorization, sparsity pattern and perm_c are reused. */
351       /* Destroy_CompRowLoc_Matrix_dist(&lu->A_sup); */ /* this leads to crash! However, see SuperLU_DIST_2.5/EXAMPLE/pzdrive2.c */
352       if (lu->FactPattern == SamePattern_SameRowPerm) {
353         lu->options.Fact = SamePattern_SameRowPerm; /* matrix has similar numerical values */
354       } else {
355         PetscStackCall("SuperLU_DIST:Destroy_LU",Destroy_LU(N, &lu->grid, &lu->LUstruct)); /* Deallocate storage associated with the L and U matrices. */
356         lu->options.Fact = SamePattern;
357       }
358     }
359     nz = 0;
360     for (i=0; i<m; i++) {
361       lu->row[i] = nz;
362       countA     = ai[i+1] - ai[i];
363       countB     = bi[i+1] - bi[i];
364       if (aj) {
365         ajj = aj + ai[i]; /* ptr to the beginning of this row */
366       } else {
367         ajj = NULL;
368       }
369       bjj = bj + bi[i];
370 
371       /* B part, smaller col index */
372       if (aj) {
373         colA_start = rstart + ajj[0]; /* the smallest global col index of A */
374       } else { /* superlu_dist does not require matrix has diagonal entries, thus aj=NULL would work */
375         colA_start = rstart;
376       }
377       jB         = 0;
378       for (j=0; j<countB; j++) {
379         jcol = garray[bjj[j]];
380         if (jcol > colA_start) {
381           jB = j;
382           break;
383         }
384         lu->col[nz]   = jcol;
385         lu->val[nz++] = *bv++;
386         if (j==countB-1) jB = countB;
387       }
388 
389       /* A part */
390       for (j=0; j<countA; j++) {
391         lu->col[nz]   = rstart + ajj[j];
392         lu->val[nz++] = *av++;
393       }
394 
395       /* B part, larger col index */
396       for (j=jB; j<countB; j++) {
397         lu->col[nz]   = garray[bjj[j]];
398         lu->val[nz++] = *bv++;
399       }
400     }
401     lu->row[m] = nz;
402 #if defined(PETSC_USE_COMPLEX)
403     PetscStackCall("SuperLU_DIST:zCreate_CompRowLoc_Matrix_dist",zCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_Z, SLU_GE));
404 #else
405     PetscStackCall("SuperLU_DIST:dCreate_CompRowLoc_Matrix_dist",dCreate_CompRowLoc_Matrix_dist(&lu->A_sup, M, N, nz, m, rstart,lu->val, lu->col, lu->row, SLU_NR_loc, SLU_D, SLU_GE));
406 #endif
407   }
408 
409   /* Factor the matrix. */
410   PetscStackCall("SuperLU_DIST:PStatInit",PStatInit(&stat));   /* Initialize the statistics variables. */
411   if (lu->MatInputMode == GLOBAL) { /* global mat input */
412 #if defined(PETSC_USE_COMPLEX)
413     PetscStackCall("SuperLU_DIST:pzgssvx_ABglobal",pzgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo));
414 #else
415     PetscStackCall("SuperLU_DIST:pdgssvx_ABglobal",pdgssvx_ABglobal(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, M, 0,&lu->grid, &lu->LUstruct, berr, &stat, &sinfo));
416 #endif
417   } else { /* distributed mat input */
418 #if defined(PETSC_USE_COMPLEX)
419     PetscStackCall("SuperLU_DIST:pzgssvx",pzgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
420     if (sinfo) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pzgssvx fails, info: %d\n",sinfo);
421 #else
422     PetscStackCall("SuperLU_DIST:pdgssvx",pdgssvx(&lu->options, &lu->A_sup, &lu->ScalePermstruct, 0, m, 0, &lu->grid,&lu->LUstruct, &lu->SOLVEstruct, berr, &stat, &sinfo));
423     if (sinfo) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"pdgssvx fails, info: %d\n",sinfo);
424 #endif
425   }
426 
427   if (lu->MatInputMode == GLOBAL && size > 1) {
428     ierr = MatDestroy(&A_seq);CHKERRQ(ierr);
429   }
430 
431   if (lu->options.PrintStat) {
432     PStatPrint(&lu->options, &stat, &lu->grid);  /* Print the statistics. */
433   }
434   PStatFree(&stat);
435   if (size > 1) {
436     F_diag            = ((Mat_MPIAIJ*)(F)->data)->A;
437     F_diag->assembled = PETSC_TRUE;
438   }
439   (F)->assembled    = PETSC_TRUE;
440   (F)->preallocated = PETSC_TRUE;
441   lu->options.Fact  = FACTORED; /* The factored form of A is supplied. Local option used by this func. only */
442   PetscFunctionReturn(0);
443 }
444 
445 /* Note the Petsc r and c permutations are ignored */
446 #undef __FUNCT__
447 #define __FUNCT__ "MatLUFactorSymbolic_SuperLU_DIST"
448 PetscErrorCode MatLUFactorSymbolic_SuperLU_DIST(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
449 {
450   Mat_SuperLU_DIST *lu = (Mat_SuperLU_DIST*)F->spptr;
451   PetscInt         M   = A->rmap->N,N=A->cmap->N;
452 
453   PetscFunctionBegin;
454   /* Initialize the SuperLU process grid. */
455   PetscStackCall("SuperLU_DIST:superlu_gridinit",superlu_gridinit(lu->comm_superlu, lu->nprow, lu->npcol, &lu->grid));
456 
457   /* Initialize ScalePermstruct and LUstruct. */
458   PetscStackCall("SuperLU_DIST:ScalePermstructInit",ScalePermstructInit(M, N, &lu->ScalePermstruct));
459   PetscStackCall("SuperLU_DIST:LUstructInit",LUstructInit(N, &lu->LUstruct));
460   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU_DIST;
461   F->ops->solve           = MatSolve_SuperLU_DIST;
462   F->ops->matsolve        = MatMatSolve_SuperLU_DIST;
463   lu->CleanUpSuperLU_Dist = PETSC_TRUE;
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "MatFactorGetSolverPackage_aij_superlu_dist"
469 PetscErrorCode MatFactorGetSolverPackage_aij_superlu_dist(Mat A,const MatSolverPackage *type)
470 {
471   PetscFunctionBegin;
472   *type = MATSOLVERSUPERLU_DIST;
473   PetscFunctionReturn(0);
474 }
475 
476 #undef __FUNCT__
477 #define __FUNCT__ "MatGetFactor_aij_superlu_dist"
478 PETSC_EXTERN PetscErrorCode MatGetFactor_aij_superlu_dist(Mat A,MatFactorType ftype,Mat *F)
479 {
480   Mat               B;
481   Mat_SuperLU_DIST  *lu;
482   PetscErrorCode    ierr;
483   PetscInt          M=A->rmap->N,N=A->cmap->N,indx;
484   PetscMPIInt       size;
485   superlu_options_t options;
486   PetscBool         flg;
487   const char        *colperm[]     = {"NATURAL","MMD_AT_PLUS_A","MMD_ATA","METIS_AT_PLUS_A","PARMETIS"};
488   const char        *rowperm[]     = {"LargeDiag","NATURAL"};
489   const char        *factPattern[] = {"SamePattern","SamePattern_SameRowPerm"};
490   PetscBool         set;
491 
492   PetscFunctionBegin;
493   /* Create the factorization matrix */
494   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
495   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,M,N);CHKERRQ(ierr);
496   ierr = MatSetType(B,((PetscObject)A)->type_name);CHKERRQ(ierr);
497   ierr = MatSeqAIJSetPreallocation(B,0,NULL);
498   ierr = MatMPIAIJSetPreallocation(B,0,NULL,0,NULL);CHKERRQ(ierr);
499 
500   B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU_DIST;
501   B->ops->view             = MatView_SuperLU_DIST;
502   B->ops->destroy          = MatDestroy_SuperLU_DIST;
503   B->ops->getdiagonal      = MatGetDiagonal_SuperLU_DIST;
504 
505   ierr = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_aij_superlu_dist);CHKERRQ(ierr);
506 
507   B->factortype = MAT_FACTOR_LU;
508 
509   ierr     = PetscNewLog(B,&lu);CHKERRQ(ierr);
510   B->spptr = lu;
511 
512   /* Set the default input options:
513      options.Fact              = DOFACT;
514      options.Equil             = YES;
515      options.ParSymbFact       = NO;
516      options.ColPerm           = METIS_AT_PLUS_A;
517      options.RowPerm           = LargeDiag;
518      options.ReplaceTinyPivot  = YES;
519      options.IterRefine        = DOUBLE;
520      options.Trans             = NOTRANS;
521      options.SolveInitialized  = NO; -hold the communication pattern used MatSolve() and MatMatSolve()
522      options.RefineInitialized = NO;
523      options.PrintStat         = YES;
524   */
525   set_default_options_dist(&options);
526 
527   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)A),&(lu->comm_superlu));CHKERRQ(ierr);
528   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
529   /* Default num of process columns and rows */
530   lu->npcol = (int_t) (0.5 + PetscSqrtReal((PetscReal)size));
531   if (!lu->npcol) lu->npcol = 1;
532   while (lu->npcol > 0) {
533     lu->nprow = (int_t) (size/lu->npcol);
534     if (size == lu->nprow * lu->npcol) break;
535     lu->npcol--;
536   }
537 
538   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU_Dist Options","Mat");CHKERRQ(ierr);
539   ierr = PetscOptionsInt("-mat_superlu_dist_r","Number rows in processor partition","None",lu->nprow,(PetscInt*)&lu->nprow,NULL);CHKERRQ(ierr);
540   ierr = PetscOptionsInt("-mat_superlu_dist_c","Number columns in processor partition","None",lu->npcol,(PetscInt*)&lu->npcol,NULL);CHKERRQ(ierr);
541   if (size != lu->nprow * lu->npcol) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number of processes %d must equal to nprow %d * npcol %d",size,lu->nprow,lu->npcol);
542 
543   lu->MatInputMode = DISTRIBUTED;
544 
545   ierr = PetscOptionsEnum("-mat_superlu_dist_matinput","Matrix input mode (global or distributed)","None",SuperLU_MatInputModes,(PetscEnum)lu->MatInputMode,(PetscEnum*)&lu->MatInputMode,NULL);CHKERRQ(ierr);
546   if (lu->MatInputMode == DISTRIBUTED && size == 1) lu->MatInputMode = GLOBAL;
547 
548   ierr = PetscOptionsBool("-mat_superlu_dist_equil","Equilibrate matrix","None",options.Equil ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
549   if (set && !flg) options.Equil = NO;
550 
551   ierr = PetscOptionsEList("-mat_superlu_dist_rowperm","Row permutation","None",rowperm,2,rowperm[0],&indx,&flg);CHKERRQ(ierr);
552   if (flg) {
553     switch (indx) {
554     case 0:
555       options.RowPerm = LargeDiag;
556       break;
557     case 1:
558       options.RowPerm = NOROWPERM;
559       break;
560     }
561   }
562 
563   ierr = PetscOptionsEList("-mat_superlu_dist_colperm","Column permutation","None",colperm,5,colperm[3],&indx,&flg);CHKERRQ(ierr);
564   if (flg) {
565     switch (indx) {
566     case 0:
567       options.ColPerm = NATURAL;
568       break;
569     case 1:
570       options.ColPerm = MMD_AT_PLUS_A;
571       break;
572     case 2:
573       options.ColPerm = MMD_ATA;
574       break;
575     case 3:
576       options.ColPerm = METIS_AT_PLUS_A;
577       break;
578     case 4:
579       options.ColPerm = PARMETIS;   /* only works for np>1 */
580       break;
581     default:
582       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
583     }
584   }
585 
586   ierr = PetscOptionsBool("-mat_superlu_dist_replacetinypivot","Replace tiny pivots","None",options.ReplaceTinyPivot ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
587   if (set && !flg) options.ReplaceTinyPivot = NO;
588 
589   options.ParSymbFact = NO;
590 
591   ierr = PetscOptionsBool("-mat_superlu_dist_parsymbfact","Parallel symbolic factorization","None",PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
592   if (set && flg && size>1) {
593     if (lu->MatInputMode == GLOBAL) {
594 #if defined(PETSC_USE_INFO)
595       ierr = PetscInfo(A,"Warning: '-mat_superlu_dist_parsymbfact' is ignored because MatInputMode=GLOBAL\n");CHKERRQ(ierr);
596 #endif
597     } else {
598 #if defined(PETSC_HAVE_PARMETIS)
599       options.ParSymbFact = YES;
600       options.ColPerm     = PARMETIS;   /* in v2.2, PARMETIS is forced for ParSymbFact regardless of user ordering setting */
601 #else
602       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"parsymbfact needs PARMETIS");
603 #endif
604     }
605   }
606 
607   lu->FactPattern = SamePattern_SameRowPerm;
608 
609   ierr = PetscOptionsEList("-mat_superlu_dist_fact","Sparsity pattern for repeated matrix factorization","None",factPattern,2,factPattern[1],&indx,&flg);CHKERRQ(ierr);
610   if (flg) {
611     switch (indx) {
612     case 0:
613       lu->FactPattern = SamePattern;
614       break;
615     case 1:
616       lu->FactPattern = SamePattern_SameRowPerm;
617       break;
618     }
619   }
620 
621   options.IterRefine = NOREFINE;
622   ierr               = PetscOptionsBool("-mat_superlu_dist_iterrefine","Use iterative refinement","None",options.IterRefine == NOREFINE ? PETSC_FALSE : PETSC_TRUE ,&flg,&set);CHKERRQ(ierr);
623   if (set) {
624     if (flg) options.IterRefine = SLU_DOUBLE;
625     else options.IterRefine = NOREFINE;
626   }
627 
628   if (PetscLogPrintInfo) options.PrintStat = YES;
629   else options.PrintStat = NO;
630   ierr = PetscOptionsBool("-mat_superlu_dist_statprint","Print factorization information","None",(PetscBool)options.PrintStat,(PetscBool*)&options.PrintStat,NULL);CHKERRQ(ierr);
631   PetscOptionsEnd();
632 
633   lu->options              = options;
634   lu->options.Fact         = DOFACT;
635   lu->matsolve_iscalled    = PETSC_FALSE;
636   lu->matmatsolve_iscalled = PETSC_FALSE;
637 
638   *F = B;
639   PetscFunctionReturn(0);
640 }
641 
642 #undef __FUNCT__
643 #define __FUNCT__ "MatSolverPackageRegister_SuperLU_DIST"
644 PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU_DIST(void)
645 {
646   PetscErrorCode ierr;
647   PetscFunctionBegin;
648   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU_DIST,MATMPIAIJ,  MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr);
649   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU_DIST,MATSEQAIJ,  MAT_FACTOR_LU,MatGetFactor_aij_superlu_dist);CHKERRQ(ierr);
650   PetscFunctionReturn(0);
651 }
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "MatFactorInfo_SuperLU_DIST"
655 PetscErrorCode MatFactorInfo_SuperLU_DIST(Mat A,PetscViewer viewer)
656 {
657   Mat_SuperLU_DIST  *lu=(Mat_SuperLU_DIST*)A->spptr;
658   superlu_options_t options;
659   PetscErrorCode    ierr;
660 
661   PetscFunctionBegin;
662   /* check if matrix is superlu_dist type */
663   if (A->ops->solve != MatSolve_SuperLU_DIST) PetscFunctionReturn(0);
664 
665   options = lu->options;
666   ierr    = PetscViewerASCIIPrintf(viewer,"SuperLU_DIST run parameters:\n");CHKERRQ(ierr);
667   ierr    = PetscViewerASCIIPrintf(viewer,"  Process grid nprow %D x npcol %D \n",lu->nprow,lu->npcol);CHKERRQ(ierr);
668   ierr    = PetscViewerASCIIPrintf(viewer,"  Equilibrate matrix %s \n",PetscBools[options.Equil != NO]);CHKERRQ(ierr);
669   ierr    = PetscViewerASCIIPrintf(viewer,"  Matrix input mode %d \n",lu->MatInputMode);CHKERRQ(ierr);
670   ierr    = PetscViewerASCIIPrintf(viewer,"  Replace tiny pivots %s \n",PetscBools[options.ReplaceTinyPivot != NO]);CHKERRQ(ierr);
671   ierr    = PetscViewerASCIIPrintf(viewer,"  Use iterative refinement %s \n",PetscBools[options.IterRefine == SLU_DOUBLE]);CHKERRQ(ierr);
672   ierr    = PetscViewerASCIIPrintf(viewer,"  Processors in row %d col partition %d \n",lu->nprow,lu->npcol);CHKERRQ(ierr);
673   ierr    = PetscViewerASCIIPrintf(viewer,"  Row permutation %s \n",(options.RowPerm == NOROWPERM) ? "NATURAL" : "LargeDiag");CHKERRQ(ierr);
674 
675   switch (options.ColPerm) {
676   case NATURAL:
677     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation NATURAL\n");CHKERRQ(ierr);
678     break;
679   case MMD_AT_PLUS_A:
680     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_AT_PLUS_A\n");CHKERRQ(ierr);
681     break;
682   case MMD_ATA:
683     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation MMD_ATA\n");CHKERRQ(ierr);
684     break;
685   case METIS_AT_PLUS_A:
686     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation METIS_AT_PLUS_A\n");CHKERRQ(ierr);
687     break;
688   case PARMETIS:
689     ierr = PetscViewerASCIIPrintf(viewer,"  Column permutation PARMETIS\n");CHKERRQ(ierr);
690     break;
691   default:
692     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown column permutation");
693   }
694 
695   ierr = PetscViewerASCIIPrintf(viewer,"  Parallel symbolic factorization %s \n",PetscBools[options.ParSymbFact != NO]);CHKERRQ(ierr);
696 
697   if (lu->FactPattern == SamePattern) {
698     ierr = PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern\n");CHKERRQ(ierr);
699   } else {
700     ierr = PetscViewerASCIIPrintf(viewer,"  Repeated factorization SamePattern_SameRowPerm\n");CHKERRQ(ierr);
701   }
702   PetscFunctionReturn(0);
703 }
704 
705 #undef __FUNCT__
706 #define __FUNCT__ "MatView_SuperLU_DIST"
707 PetscErrorCode MatView_SuperLU_DIST(Mat A,PetscViewer viewer)
708 {
709   PetscErrorCode    ierr;
710   PetscBool         iascii;
711   PetscViewerFormat format;
712 
713   PetscFunctionBegin;
714   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
715   if (iascii) {
716     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
717     if (format == PETSC_VIEWER_ASCII_INFO) {
718       ierr = MatFactorInfo_SuperLU_DIST(A,viewer);CHKERRQ(ierr);
719     }
720   }
721   PetscFunctionReturn(0);
722 }
723 
724 
725 /*MC
726   MATSOLVERSUPERLU_DIST - Parallel direct solver package for LU factorization
727 
728   Use ./configure --download-superlu_dist --download-parmetis --download-metis --download-ptscotch  to have PETSc installed with SuperLU_DIST
729 
730   Use -pc_type lu -pc_factor_mat_solver_package superlu_dist to us this direct solver
731 
732    Works with AIJ matrices
733 
734   Options Database Keys:
735 + -mat_superlu_dist_r <n> - number of rows in processor partition
736 . -mat_superlu_dist_c <n> - number of columns in processor partition
737 . -mat_superlu_dist_matinput <0,1> - matrix input mode; 0=global, 1=distributed
738 . -mat_superlu_dist_equil - equilibrate the matrix
739 . -mat_superlu_dist_rowperm <LargeDiag,NATURAL> - row permutation
740 . -mat_superlu_dist_colperm <MMD_AT_PLUS_A,MMD_ATA,NATURAL> - column permutation
741 . -mat_superlu_dist_replacetinypivot - replace tiny pivots
742 . -mat_superlu_dist_fact <SamePattern> - (choose one of) SamePattern SamePattern_SameRowPerm
743 . -mat_superlu_dist_iterrefine - use iterative refinement
744 - -mat_superlu_dist_statprint - print factorization information
745 
746    Level: beginner
747 
748 .seealso: PCLU
749 
750 .seealso: PCFactorSetMatSolverPackage(), MatSolverPackage
751 
752 M*/
753 
754