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