xref: /petsc/src/mat/impls/aij/seq/superlu/superlu.c (revision 603e8f969c662e39ea604366202a4e5ffbbd2908)
114b396bbSKris Buschelman 
25a46d3faSBarry Smith /*  --------------------------------------------------------------------
35a46d3faSBarry Smith 
45a46d3faSBarry Smith      This file implements a subclass of the SeqAIJ matrix class that uses
5c2b89b5dSBarry Smith      the SuperLU sparse solver.
614b396bbSKris Buschelman */
714b396bbSKris Buschelman 
85a46d3faSBarry Smith /*
95a46d3faSBarry Smith      Defines the data structure for the base matrix type (SeqAIJ)
105a46d3faSBarry Smith */
11c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>    /*I "petscmat.h" I*/
1214b396bbSKris Buschelman 
135a46d3faSBarry Smith /*
145a46d3faSBarry Smith      SuperLU include files
155a46d3faSBarry Smith */
1614b396bbSKris Buschelman EXTERN_C_BEGIN
1714b396bbSKris Buschelman #if defined(PETSC_USE_COMPLEX)
183cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
193cb270beSHong Zhang #include <slu_cdefs.h>
203cb270beSHong Zhang #else
21c6db04a5SJed Brown #include <slu_zdefs.h>
223cb270beSHong Zhang #endif
233cb270beSHong Zhang #else
243cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
253cb270beSHong Zhang #include <slu_sdefs.h>
2614b396bbSKris Buschelman #else
27c6db04a5SJed Brown #include <slu_ddefs.h>
2814b396bbSKris Buschelman #endif
293cb270beSHong Zhang #endif
30c6db04a5SJed Brown #include <slu_util.h>
3114b396bbSKris Buschelman EXTERN_C_END
3214b396bbSKris Buschelman 
335a46d3faSBarry Smith /*
34245fece9SBarry Smith      This is the data that defines the SuperLU factored matrix type
355a46d3faSBarry Smith */
3614b396bbSKris Buschelman typedef struct {
3775af56d4SHong Zhang   SuperMatrix       A,L,U,B,X;
3875af56d4SHong Zhang   superlu_options_t options;
39da7a1d00SHong Zhang   PetscInt          *perm_c; /* column permutation vector */
40da7a1d00SHong Zhang   PetscInt          *perm_r; /* row permutations from partial pivoting */
41da7a1d00SHong Zhang   PetscInt          *etree;
42da7a1d00SHong Zhang   PetscReal         *R, *C;
4375af56d4SHong Zhang   char              equed[1];
44da7a1d00SHong Zhang   PetscInt          lwork;
4575af56d4SHong Zhang   void              *work;
46da7a1d00SHong Zhang   PetscReal         rpg, rcond;
4775af56d4SHong Zhang   mem_usage_t       mem_usage;
4875af56d4SHong Zhang   MatStructure      flg;
495d8b2efaSHong Zhang   SuperLUStat_t     stat;
50cae5a23dSHong Zhang   Mat               A_dup;
51cae5a23dSHong Zhang   PetscScalar       *rhs_dup;
52c147c726SHong Zhang   GlobalLU_t        Glu;
5314b396bbSKris Buschelman 
5414b396bbSKris Buschelman   /* Flag to clean up (non-global) SuperLU objects during Destroy */
55ace3abfcSBarry Smith   PetscBool CleanUpSuperLU;
56f0c56d0fSKris Buschelman } Mat_SuperLU;
5714b396bbSKris Buschelman 
585a46d3faSBarry Smith /*
595a46d3faSBarry Smith     Utility function
605a46d3faSBarry Smith */
615a46d3faSBarry Smith #undef __FUNCT__
625a46d3faSBarry Smith #define __FUNCT__ "MatFactorInfo_SuperLU"
63245fece9SBarry Smith static PetscErrorCode MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
645a46d3faSBarry Smith {
65245fece9SBarry Smith   Mat_SuperLU       *lu= (Mat_SuperLU*)A->data;
665a46d3faSBarry Smith   PetscErrorCode    ierr;
675a46d3faSBarry Smith   superlu_options_t options;
685a46d3faSBarry Smith 
695a46d3faSBarry Smith   PetscFunctionBegin;
705a46d3faSBarry Smith   options = lu->options;
712205254eSKarl Rupp 
725a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");CHKERRQ(ierr);
735a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  Equil: %s\n",(options.Equil != NO) ? "YES" : "NO");CHKERRQ(ierr);
745a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ColPerm: %D\n",options.ColPerm);CHKERRQ(ierr);
755a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  IterRefine: %D\n",options.IterRefine);CHKERRQ(ierr);
765a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  SymmetricMode: %s\n",(options.SymmetricMode != NO) ? "YES" : "NO");CHKERRQ(ierr);
775a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  DiagPivotThresh: %g\n",options.DiagPivotThresh);CHKERRQ(ierr);
785a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PivotGrowth: %s\n",(options.PivotGrowth != NO) ? "YES" : "NO");CHKERRQ(ierr);
795a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ConditionNumber: %s\n",(options.ConditionNumber != NO) ? "YES" : "NO");CHKERRQ(ierr);
805a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  RowPerm: %D\n",options.RowPerm);CHKERRQ(ierr);
815a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  ReplaceTinyPivot: %s\n",(options.ReplaceTinyPivot != NO) ? "YES" : "NO");CHKERRQ(ierr);
825a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  PrintStat: %s\n",(options.PrintStat != NO) ? "YES" : "NO");CHKERRQ(ierr);
835a46d3faSBarry Smith   ierr = PetscViewerASCIIPrintf(viewer,"  lwork: %D\n",lu->lwork);CHKERRQ(ierr);
84d5f3da31SBarry Smith   if (A->factortype == MAT_FACTOR_ILU) {
85cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropTol: %g\n",options.ILU_DropTol);CHKERRQ(ierr);
86cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillTol: %g\n",options.ILU_FillTol);CHKERRQ(ierr);
87cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_FillFactor: %g\n",options.ILU_FillFactor);CHKERRQ(ierr);
88cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_DropRule: %D\n",options.ILU_DropRule);CHKERRQ(ierr);
89cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_Norm: %D\n",options.ILU_Norm);CHKERRQ(ierr);
90cffbb591SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  ILU_MILU: %D\n",options.ILU_MILU);CHKERRQ(ierr);
91cffbb591SHong Zhang   }
925a46d3faSBarry Smith   PetscFunctionReturn(0);
935a46d3faSBarry Smith }
945a46d3faSBarry Smith 
95245fece9SBarry Smith #undef __FUNCT__
96245fece9SBarry Smith #define __FUNCT__ "MatSolve_SuperLU_Private"
97245fece9SBarry Smith PetscErrorCode MatSolve_SuperLU_Private(Mat A,Vec b,Vec x)
98245fece9SBarry Smith {
99245fece9SBarry Smith   Mat_SuperLU       *lu = (Mat_SuperLU*)A->data;
100245fece9SBarry Smith   const PetscScalar *barray;
101245fece9SBarry Smith   PetscScalar       *xarray;
102245fece9SBarry Smith   PetscErrorCode    ierr;
103245fece9SBarry Smith   PetscInt          info,i,n;
104245fece9SBarry Smith   PetscReal         ferr,berr;
105245fece9SBarry Smith   static PetscBool  cite = PETSC_FALSE;
106245fece9SBarry Smith 
107245fece9SBarry Smith   PetscFunctionBegin;
108245fece9SBarry Smith   if (lu->lwork == -1) PetscFunctionReturn(0);
109245fece9SBarry Smith   ierr = PetscCitationsRegister("@article{superlu99,\n  author  = {James W. Demmel and Stanley C. Eisenstat and\n             John R. Gilbert and Xiaoye S. Li and Joseph W. H. Liu},\n  title = {A supernodal approach to sparse partial pivoting},\n  journal = {SIAM J. Matrix Analysis and Applications},\n  year = {1999},\n  volume  = {20},\n  number = {3},\n  pages = {720-755}\n}\n",&cite);CHKERRQ(ierr);
110245fece9SBarry Smith 
111245fece9SBarry Smith   ierr = VecGetLocalSize(x,&n);CHKERRQ(ierr);
112245fece9SBarry Smith   lu->B.ncol = 1;   /* Set the number of right-hand side */
113245fece9SBarry Smith   if (lu->options.Equil && !lu->rhs_dup) {
114245fece9SBarry Smith     /* superlu overwrites b when Equil is used, thus create rhs_dup to keep user's b unchanged */
115245fece9SBarry Smith     ierr = PetscMalloc1(n,&lu->rhs_dup);CHKERRQ(ierr);
116245fece9SBarry Smith   }
117245fece9SBarry Smith   if (lu->options.Equil) {
118245fece9SBarry Smith     /* Copy b into rsh_dup */
119245fece9SBarry Smith     ierr   = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
120245fece9SBarry Smith     ierr   = PetscMemcpy(lu->rhs_dup,barray,n*sizeof(PetscScalar));CHKERRQ(ierr);
121245fece9SBarry Smith     ierr   = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
122245fece9SBarry Smith     barray = lu->rhs_dup;
123245fece9SBarry Smith   } else {
124245fece9SBarry Smith     ierr = VecGetArrayRead(b,&barray);CHKERRQ(ierr);
125245fece9SBarry Smith   }
126245fece9SBarry Smith   ierr = VecGetArray(x,&xarray);CHKERRQ(ierr);
127245fece9SBarry Smith 
128245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX)
129245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
130245fece9SBarry Smith   ((DNformat*)lu->B.Store)->nzval = (singlecomplex*)barray;
131245fece9SBarry Smith   ((DNformat*)lu->X.Store)->nzval = (singlecomplex*)xarray;
132245fece9SBarry Smith #else
133245fece9SBarry Smith   ((DNformat*)lu->B.Store)->nzval = (doublecomplex*)barray;
134245fece9SBarry Smith   ((DNformat*)lu->X.Store)->nzval = (doublecomplex*)xarray;
135245fece9SBarry Smith #endif
136245fece9SBarry Smith #else
137245fece9SBarry Smith   ((DNformat*)lu->B.Store)->nzval = (void*)barray;
138245fece9SBarry Smith   ((DNformat*)lu->X.Store)->nzval = xarray;
139245fece9SBarry Smith #endif
140245fece9SBarry Smith 
141245fece9SBarry Smith   lu->options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
142245fece9SBarry Smith   if (A->factortype == MAT_FACTOR_LU) {
143245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX)
144245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
145245fece9SBarry Smith     PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
146245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
147245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
148245fece9SBarry Smith #else
149245fece9SBarry Smith     PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
150245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
151245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
152245fece9SBarry Smith #endif
153245fece9SBarry Smith #else
154245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
155245fece9SBarry Smith     PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
156245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
157245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
158245fece9SBarry Smith #else
159245fece9SBarry Smith     PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
160245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
161245fece9SBarry Smith                                      &lu->Glu,&lu->mem_usage, &lu->stat, &info));
162245fece9SBarry Smith #endif
163245fece9SBarry Smith #endif
164245fece9SBarry Smith   } else if (A->factortype == MAT_FACTOR_ILU) {
165245fece9SBarry Smith #if defined(PETSC_USE_COMPLEX)
166245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
167245fece9SBarry Smith     PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
168245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
169245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
170245fece9SBarry Smith #else
171245fece9SBarry Smith     PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
172245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
173245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
174245fece9SBarry Smith #endif
175245fece9SBarry Smith #else
176245fece9SBarry Smith #if defined(PETSC_USE_REAL_SINGLE)
177245fece9SBarry Smith     PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
178245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
179245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
180245fece9SBarry Smith #else
181245fece9SBarry Smith     PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
182245fece9SBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
183245fece9SBarry Smith                                      &lu->Glu, &lu->mem_usage, &lu->stat, &info));
184245fece9SBarry Smith #endif
185245fece9SBarry Smith #endif
186245fece9SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
187245fece9SBarry Smith   if (!lu->options.Equil) {
188245fece9SBarry Smith     ierr = VecRestoreArrayRead(b,&barray);CHKERRQ(ierr);
189245fece9SBarry Smith   }
190245fece9SBarry Smith   ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr);
191245fece9SBarry Smith 
192245fece9SBarry Smith   if (!info || info == lu->A.ncol+1) {
193245fece9SBarry Smith     if (lu->options.IterRefine) {
194245fece9SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"Iterative Refinement:\n");
195245fece9SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  %8s%8s%16s%16s\n", "rhs", "Steps", "FERR", "BERR");
196245fece9SBarry Smith       for (i = 0; i < 1; ++i) {
197245fece9SBarry Smith         ierr = PetscPrintf(PETSC_COMM_SELF,"  %8d%8d%16e%16e\n", i+1, lu->stat.RefineSteps, ferr, berr);
198245fece9SBarry Smith       }
199245fece9SBarry Smith     }
200245fece9SBarry Smith   } else if (info > 0) {
201245fece9SBarry Smith     if (lu->lwork == -1) {
202245fece9SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  ** Estimated memory: %D bytes\n", info - lu->A.ncol);
203245fece9SBarry Smith     } else {
204245fece9SBarry Smith       ierr = PetscPrintf(PETSC_COMM_SELF,"  Warning: gssvx() returns info %D\n",info);
205245fece9SBarry Smith     }
206245fece9SBarry Smith   } else if (info < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", info,-info);
207245fece9SBarry Smith 
208245fece9SBarry Smith   if (lu->options.PrintStat) {
209245fece9SBarry Smith     ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve__SuperLU():\n");
210245fece9SBarry Smith     PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
211245fece9SBarry Smith   }
212245fece9SBarry Smith   PetscFunctionReturn(0);
213245fece9SBarry Smith }
214245fece9SBarry Smith 
215245fece9SBarry Smith #undef __FUNCT__
216245fece9SBarry Smith #define __FUNCT__ "MatSolve_SuperLU"
217245fece9SBarry Smith PetscErrorCode MatSolve_SuperLU(Mat A,Vec b,Vec x)
218245fece9SBarry Smith {
219245fece9SBarry Smith   Mat_SuperLU    *lu = (Mat_SuperLU*)A->data;
220245fece9SBarry Smith   PetscErrorCode ierr;
221245fece9SBarry Smith 
222245fece9SBarry Smith   PetscFunctionBegin;
223*603e8f96SBarry Smith   if (A->factorerrortype) {
224245fece9SBarry Smith     ierr = PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n");CHKERRQ(ierr);
225245fece9SBarry Smith     ierr = VecSetInf(x);CHKERRQ(ierr);
226245fece9SBarry Smith     PetscFunctionReturn(0);
227245fece9SBarry Smith   }
228245fece9SBarry Smith 
229245fece9SBarry Smith   lu->options.Trans = TRANS;
230245fece9SBarry Smith   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
231245fece9SBarry Smith   PetscFunctionReturn(0);
232245fece9SBarry Smith }
233245fece9SBarry Smith 
234245fece9SBarry Smith #undef __FUNCT__
235245fece9SBarry Smith #define __FUNCT__ "MatSolveTranspose_SuperLU"
236245fece9SBarry Smith PetscErrorCode MatSolveTranspose_SuperLU(Mat A,Vec b,Vec x)
237245fece9SBarry Smith {
238245fece9SBarry Smith   Mat_SuperLU    *lu = (Mat_SuperLU*)A->data;
239245fece9SBarry Smith   PetscErrorCode ierr;
240245fece9SBarry Smith 
241245fece9SBarry Smith   PetscFunctionBegin;
242*603e8f96SBarry Smith   if (A->factorerrortype) {
243245fece9SBarry Smith     ierr = PetscInfo(A,"MatSolve is called with singular matrix factor, skip\n");CHKERRQ(ierr);
244245fece9SBarry Smith     ierr = VecSetInf(x);CHKERRQ(ierr);
245245fece9SBarry Smith     PetscFunctionReturn(0);
246245fece9SBarry Smith   }
247245fece9SBarry Smith 
248245fece9SBarry Smith   lu->options.Trans = NOTRANS;
249245fece9SBarry Smith   ierr = MatSolve_SuperLU_Private(A,b,x);CHKERRQ(ierr);
250245fece9SBarry Smith   PetscFunctionReturn(0);
251245fece9SBarry Smith }
252245fece9SBarry Smith 
2535a46d3faSBarry Smith #undef __FUNCT__
2545a46d3faSBarry Smith #define __FUNCT__ "MatLUFactorNumeric_SuperLU"
255245fece9SBarry Smith static PetscErrorCode MatLUFactorNumeric_SuperLU(Mat F,Mat A,const MatFactorInfo *info)
2565a46d3faSBarry Smith {
257245fece9SBarry Smith   Mat_SuperLU    *lu = (Mat_SuperLU*)F->data;
258cae5a23dSHong Zhang   Mat_SeqAIJ     *aa;
2595a46d3faSBarry Smith   PetscErrorCode ierr;
2605a46d3faSBarry Smith   PetscInt       sinfo;
2615a46d3faSBarry Smith   PetscReal      ferr, berr;
2625a46d3faSBarry Smith   NCformat       *Ustore;
2635a46d3faSBarry Smith   SCformat       *Lstore;
2645a46d3faSBarry Smith 
2655a46d3faSBarry Smith   PetscFunctionBegin;
2665a46d3faSBarry Smith   if (lu->flg == SAME_NONZERO_PATTERN) { /* successing numerical factorization */
2675a46d3faSBarry Smith     lu->options.Fact = SamePattern;
2685a46d3faSBarry Smith     /* Ref: ~SuperLU_3.0/EXAMPLE/dlinsolx2.c */
2695a46d3faSBarry Smith     Destroy_SuperMatrix_Store(&lu->A);
270cae5a23dSHong Zhang     if (lu->options.Equil) {
271cae5a23dSHong Zhang       ierr = MatCopy_SeqAIJ(A,lu->A_dup,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
272cae5a23dSHong Zhang     }
2735a46d3faSBarry Smith     if (lu->lwork >= 0) {
274d387c056SBarry Smith       PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
275d387c056SBarry Smith       PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
2765a46d3faSBarry Smith       lu->options.Fact = SamePattern;
2775a46d3faSBarry Smith     }
2785a46d3faSBarry Smith   }
2795a46d3faSBarry Smith 
2805a46d3faSBarry Smith   /* Create the SuperMatrix for lu->A=A^T:
2815a46d3faSBarry Smith        Since SuperLU likes column-oriented matrices,we pass it the transpose,
2825a46d3faSBarry Smith        and then solve A^T X = B in MatSolve(). */
283cae5a23dSHong Zhang   if (lu->options.Equil) {
284cae5a23dSHong Zhang     aa = (Mat_SeqAIJ*)(lu->A_dup)->data;
285cae5a23dSHong Zhang   } else {
286cae5a23dSHong Zhang     aa = (Mat_SeqAIJ*)(A)->data;
287cae5a23dSHong Zhang   }
2885a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
2893cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
290d387c056SBarry Smith   PetscStackCall("SuperLU:cCreate_CompCol_Matrix",cCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(singlecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_C,SLU_GE));
2913cb270beSHong Zhang #else
292d387c056SBarry Smith   PetscStackCall("SuperLU:zCreate_CompCol_Matrix",zCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,(doublecomplex*)aa->a,aa->j,aa->i,SLU_NC,SLU_Z,SLU_GE));
2933cb270beSHong Zhang #endif
2943cb270beSHong Zhang #else
2953cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
296d387c056SBarry Smith   PetscStackCall("SuperLU:sCreate_CompCol_Matrix",sCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_S,SLU_GE));
2975a46d3faSBarry Smith #else
298d387c056SBarry Smith   PetscStackCall("SuperLU:dCreate_CompCol_Matrix",dCreate_CompCol_Matrix(&lu->A,A->cmap->n,A->rmap->n,aa->nz,aa->a,aa->j,aa->i,SLU_NC,SLU_D,SLU_GE));
2995a46d3faSBarry Smith #endif
3003cb270beSHong Zhang #endif
3015a46d3faSBarry Smith 
3025a46d3faSBarry Smith   /* Numerical factorization */
3035a46d3faSBarry Smith   lu->B.ncol = 0;  /* Indicate not to solve the system */
304d5f3da31SBarry Smith   if (F->factortype == MAT_FACTOR_LU) {
3055a46d3faSBarry Smith #if defined(PETSC_USE_COMPLEX)
3063cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
307d387c056SBarry Smith     PetscStackCall("SuperLU:cgssvx",cgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3083cb270beSHong Zhang                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3095db81dd2SSatish Balay                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
3103cb270beSHong Zhang #else
311d387c056SBarry Smith     PetscStackCall("SuperLU:zgssvx",zgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3125a46d3faSBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3135db81dd2SSatish Balay                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
3143cb270beSHong Zhang #endif
3153cb270beSHong Zhang #else
3163cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
317d387c056SBarry Smith     PetscStackCall("SuperLU:sgssvx",sgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3183cb270beSHong Zhang                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
3195db81dd2SSatish Balay                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
3205a46d3faSBarry Smith #else
321d387c056SBarry Smith     PetscStackCall("SuperLU:dgssvx",dgssvx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3225a46d3faSBarry Smith                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond, &ferr, &berr,
323c147c726SHong Zhang                                      &lu->Glu,&lu->mem_usage, &lu->stat, &sinfo));
3245a46d3faSBarry Smith #endif
3253cb270beSHong Zhang #endif
326d5f3da31SBarry Smith   } else if (F->factortype == MAT_FACTOR_ILU) {
327cffbb591SHong Zhang     /* Compute the incomplete factorization, condition number and pivot growth */
328cffbb591SHong Zhang #if defined(PETSC_USE_COMPLEX)
3293cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
330d387c056SBarry Smith     PetscStackCall("SuperLU:cgsisx",cgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
3313cb270beSHong Zhang                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3325db81dd2SSatish Balay                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
3333cb270beSHong Zhang #else
334d387c056SBarry Smith     PetscStackCall("SuperLU:zgsisx",zgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r,lu->etree, lu->equed, lu->R, lu->C,
335cffbb591SHong Zhang                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3365db81dd2SSatish Balay                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
3373cb270beSHong Zhang #endif
3383cb270beSHong Zhang #else
3393cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
340d387c056SBarry Smith     PetscStackCall("SuperLU:sgsisx",sgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
3413cb270beSHong Zhang                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
3425db81dd2SSatish Balay                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
343cffbb591SHong Zhang #else
344d387c056SBarry Smith     PetscStackCall("SuperLU:dgsisx",dgsisx(&lu->options, &lu->A, lu->perm_c, lu->perm_r, lu->etree, lu->equed, lu->R, lu->C,
345cffbb591SHong Zhang                                      &lu->L, &lu->U, lu->work, lu->lwork, &lu->B, &lu->X, &lu->rpg, &lu->rcond,
346c147c726SHong Zhang                                      &lu->Glu, &lu->mem_usage, &lu->stat, &sinfo));
347cffbb591SHong Zhang #endif
3483cb270beSHong Zhang #endif
349f23aa3ddSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
3505a46d3faSBarry Smith   if (!sinfo || sinfo == lu->A.ncol+1) {
3512205254eSKarl Rupp     if (lu->options.PivotGrowth) {
352675d1226SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. pivot growth = %e\n", lu->rpg);CHKERRQ(ierr);
3532205254eSKarl Rupp     }
3542205254eSKarl Rupp     if (lu->options.ConditionNumber) {
355675d1226SHong Zhang       ierr = PetscPrintf(PETSC_COMM_SELF,"  Recip. condition number = %e\n", lu->rcond);CHKERRQ(ierr);
3562205254eSKarl Rupp     }
3575a46d3faSBarry Smith   } else if (sinfo > 0) {
358675d1226SHong Zhang     if (A->erroriffailure) {
359675d1226SHong Zhang       SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot in row %D",sinfo);
360675d1226SHong Zhang     } else {
361675d1226SHong Zhang       if (sinfo <= lu->A.ncol) {
362675d1226SHong Zhang         if (lu->options.ILU_FillTol == 0.0) {
363*603e8f96SBarry Smith           F->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
364675d1226SHong Zhang         }
365675d1226SHong Zhang         ierr = PetscInfo2(F,"Number of zero pivots %D, ILU_FillTol %g\n",sinfo,lu->options.ILU_FillTol);CHKERRQ(ierr);
366675d1226SHong Zhang       } else if (sinfo == lu->A.ncol + 1) {
367675d1226SHong Zhang         /*
368675d1226SHong Zhang          U is nonsingular, but RCOND is less than machine
369675d1226SHong Zhang  		      precision, meaning that the matrix is singular to
370675d1226SHong Zhang  		      working precision. Nevertheless, the solution and
371675d1226SHong Zhang  		      error bounds are computed because there are a number
372675d1226SHong Zhang  		      of situations where the computed solution can be more
373675d1226SHong Zhang  		      accurate than the value of RCOND would suggest.
374675d1226SHong Zhang          */
375675d1226SHong Zhang         ierr = PetscInfo1(F,"Matrix factor U is nonsingular, but is singular to working precision. The solution is computed. info %D",sinfo);CHKERRQ(ierr);
376675d1226SHong Zhang       } else { /* sinfo > lu->A.ncol + 1 */
377*603e8f96SBarry Smith         F->factorerrortype = MAT_FACTOR_OUTMEMORY;
378675d1226SHong Zhang         ierr = PetscInfo1(F,"Number of bytes allocated when memory allocation fails %D\n",sinfo);CHKERRQ(ierr);
379675d1226SHong Zhang       }
380675d1226SHong Zhang     }
381f23aa3ddSBarry Smith   } else SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB, "info = %D, the %D-th argument in gssvx() had an illegal value", sinfo,-sinfo);
3825a46d3faSBarry Smith 
3835a46d3faSBarry Smith   if (lu->options.PrintStat) {
384675d1226SHong Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"MatLUFactorNumeric_SuperLU():\n");CHKERRQ(ierr);
385d387c056SBarry Smith     PetscStackCall("SuperLU:StatPrint",StatPrint(&lu->stat));
3865a46d3faSBarry Smith     Lstore = (SCformat*) lu->L.Store;
3875a46d3faSBarry Smith     Ustore = (NCformat*) lu->U.Store;
3885a46d3faSBarry Smith     ierr   = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor L = %D\n", Lstore->nnz);
3895a46d3faSBarry Smith     ierr   = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in factor U = %D\n", Ustore->nnz);
3905a46d3faSBarry Smith     ierr   = PetscPrintf(PETSC_COMM_SELF,"  No of nonzeros in L+U = %D\n", Lstore->nnz + Ustore->nnz - lu->A.ncol);
3916da386baSHong Zhang     ierr   = PetscPrintf(PETSC_COMM_SELF,"  L\\U MB %.3f\ttotal MB needed %.3f\n",
3926da386baSHong Zhang                          lu->mem_usage.for_lu/1e6, lu->mem_usage.total_needed/1e6);
3935a46d3faSBarry Smith   }
3945a46d3faSBarry Smith 
3955a46d3faSBarry Smith   lu->flg                = SAME_NONZERO_PATTERN;
3961d5ca7f3SHong Zhang   F->ops->solve          = MatSolve_SuperLU;
3971d5ca7f3SHong Zhang   F->ops->solvetranspose = MatSolveTranspose_SuperLU;
3981aef8b4cSStefano Zampini   F->ops->matsolve       = NULL;
3995a46d3faSBarry Smith   PetscFunctionReturn(0);
4005a46d3faSBarry Smith }
4015a46d3faSBarry Smith 
40214b396bbSKris Buschelman #undef __FUNCT__
403f0c56d0fSKris Buschelman #define __FUNCT__ "MatDestroy_SuperLU"
404245fece9SBarry Smith static PetscErrorCode MatDestroy_SuperLU(Mat A)
40514b396bbSKris Buschelman {
406dfbe8321SBarry Smith   PetscErrorCode ierr;
407245fece9SBarry Smith   Mat_SuperLU    *lu=(Mat_SuperLU*)A->data;
40814b396bbSKris Buschelman 
40914b396bbSKris Buschelman   PetscFunctionBegin;
410245fece9SBarry Smith   if (lu->CleanUpSuperLU) { /* Free the SuperLU datastructures */
411d387c056SBarry Smith     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->A));
412d387c056SBarry Smith     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->B));
413d387c056SBarry Smith     PetscStackCall("SuperLU:Destroy_SuperMatrix_Store",Destroy_SuperMatrix_Store(&lu->X));
414d387c056SBarry Smith     PetscStackCall("SuperLU:StatFree",StatFree(&lu->stat));
4150e742b69SHong Zhang     if (lu->lwork >= 0) {
416d387c056SBarry Smith       PetscStackCall("SuperLU:Destroy_SuperNode_Matrix",Destroy_SuperNode_Matrix(&lu->L));
417d387c056SBarry Smith       PetscStackCall("SuperLU:Destroy_CompCol_Matrix",Destroy_CompCol_Matrix(&lu->U));
4180e742b69SHong Zhang     }
4190e742b69SHong Zhang   }
4209ce50919SHong Zhang   ierr = PetscFree(lu->etree);CHKERRQ(ierr);
4219ce50919SHong Zhang   ierr = PetscFree(lu->perm_r);CHKERRQ(ierr);
4229ce50919SHong Zhang   ierr = PetscFree(lu->perm_c);CHKERRQ(ierr);
4239ce50919SHong Zhang   ierr = PetscFree(lu->R);CHKERRQ(ierr);
4249ce50919SHong Zhang   ierr = PetscFree(lu->C);CHKERRQ(ierr);
425bf0cc555SLisandro Dalcin   ierr = PetscFree(lu->rhs_dup);CHKERRQ(ierr);
426bf0cc555SLisandro Dalcin   ierr = MatDestroy(&lu->A_dup);CHKERRQ(ierr);
427245fece9SBarry Smith   ierr = PetscFree(A->data);CHKERRQ(ierr);
4280e742b69SHong Zhang 
429d954db57SHong Zhang   /* clear composed functions */
430bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatFactorGetSolverPackage_C",NULL);CHKERRQ(ierr);
431bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatSuperluSetILUDropTol_C",NULL);CHKERRQ(ierr);
43214b396bbSKris Buschelman   PetscFunctionReturn(0);
43314b396bbSKris Buschelman }
43414b396bbSKris Buschelman 
43514b396bbSKris Buschelman #undef __FUNCT__
436f0c56d0fSKris Buschelman #define __FUNCT__ "MatView_SuperLU"
437245fece9SBarry Smith static PetscErrorCode MatView_SuperLU(Mat A,PetscViewer viewer)
43814b396bbSKris Buschelman {
439dfbe8321SBarry Smith   PetscErrorCode    ierr;
440ace3abfcSBarry Smith   PetscBool         iascii;
44114b396bbSKris Buschelman   PetscViewerFormat format;
44214b396bbSKris Buschelman 
44314b396bbSKris Buschelman   PetscFunctionBegin;
444251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
44532077d6dSBarry Smith   if (iascii) {
44614b396bbSKris Buschelman     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
4472f59403fSHong Zhang     if (format == PETSC_VIEWER_ASCII_INFO) {
448f0c56d0fSKris Buschelman       ierr = MatFactorInfo_SuperLU(A,viewer);CHKERRQ(ierr);
44914b396bbSKris Buschelman     }
45014b396bbSKris Buschelman   }
45114b396bbSKris Buschelman   PetscFunctionReturn(0);
45214b396bbSKris Buschelman }
45314b396bbSKris Buschelman 
454e0b74bf9SHong Zhang #undef __FUNCT__
455e0b74bf9SHong Zhang #define __FUNCT__ "MatMatSolve_SuperLU"
456e0b74bf9SHong Zhang PetscErrorCode MatMatSolve_SuperLU(Mat A,Mat B,Mat X)
457e0b74bf9SHong Zhang {
458245fece9SBarry Smith   Mat_SuperLU    *lu = (Mat_SuperLU*)A->data;
459cd723cd1SBarry Smith   PetscBool      flg;
460cd723cd1SBarry Smith   PetscErrorCode ierr;
461e0b74bf9SHong Zhang 
462e0b74bf9SHong Zhang   PetscFunctionBegin;
4630298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
464ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix");
4650298fd71SBarry Smith   ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr);
466ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix");
4672205254eSKarl Rupp   lu->options.Trans = TRANS;
468e0b74bf9SHong Zhang   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatSolve_SuperLU() is not implemented yet");
469e0b74bf9SHong Zhang   PetscFunctionReturn(0);
470e0b74bf9SHong Zhang }
471e0b74bf9SHong Zhang 
47214b396bbSKris Buschelman /*
47314b396bbSKris Buschelman    Note the r permutation is ignored
47414b396bbSKris Buschelman */
47514b396bbSKris Buschelman #undef __FUNCT__
476f0c56d0fSKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SuperLU"
477245fece9SBarry Smith static PetscErrorCode MatLUFactorSymbolic_SuperLU(Mat F,Mat A,IS r,IS c,const MatFactorInfo *info)
47814b396bbSKris Buschelman {
479245fece9SBarry Smith   Mat_SuperLU *lu = (Mat_SuperLU*)(F->data);
480b24902e0SBarry Smith 
481b24902e0SBarry Smith   PetscFunctionBegin;
482b24902e0SBarry Smith   lu->flg                 = DIFFERENT_NONZERO_PATTERN;
483b24902e0SBarry Smith   lu->CleanUpSuperLU      = PETSC_TRUE;
4841d5ca7f3SHong Zhang   F->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
485b24902e0SBarry Smith   PetscFunctionReturn(0);
486b24902e0SBarry Smith }
487b24902e0SBarry Smith 
48835bd34faSBarry Smith #undef __FUNCT__
4895ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol_SuperLU"
490b2573a8aSBarry Smith static PetscErrorCode MatSuperluSetILUDropTol_SuperLU(Mat F,PetscReal dtol)
4915ccb76cbSHong Zhang {
492245fece9SBarry Smith   Mat_SuperLU *lu= (Mat_SuperLU*)F->data;
4935ccb76cbSHong Zhang 
4945ccb76cbSHong Zhang   PetscFunctionBegin;
4955ccb76cbSHong Zhang   lu->options.ILU_DropTol = dtol;
4965ccb76cbSHong Zhang   PetscFunctionReturn(0);
4975ccb76cbSHong Zhang }
4985ccb76cbSHong Zhang 
4995ccb76cbSHong Zhang #undef __FUNCT__
5005ccb76cbSHong Zhang #define __FUNCT__ "MatSuperluSetILUDropTol"
5015ccb76cbSHong Zhang /*@
5025ccb76cbSHong Zhang   MatSuperluSetILUDropTol - Set SuperLU ILU drop tolerance
5035ccb76cbSHong Zhang    Logically Collective on Mat
5045ccb76cbSHong Zhang 
5055ccb76cbSHong Zhang    Input Parameters:
5065ccb76cbSHong Zhang +  F - the factored matrix obtained by calling MatGetFactor() from PETSc-SuperLU interface
5075ccb76cbSHong Zhang -  dtol - drop tolerance
5085ccb76cbSHong Zhang 
5095ccb76cbSHong Zhang   Options Database:
5105ccb76cbSHong Zhang .   -mat_superlu_ilu_droptol <dtol>
5115ccb76cbSHong Zhang 
5125ccb76cbSHong Zhang    Level: beginner
5135ccb76cbSHong Zhang 
51496a0c994SBarry Smith    References:
51596a0c994SBarry Smith .      SuperLU Users' Guide
5165ccb76cbSHong Zhang 
5175ccb76cbSHong Zhang .seealso: MatGetFactor()
5185ccb76cbSHong Zhang @*/
5195ccb76cbSHong Zhang PetscErrorCode MatSuperluSetILUDropTol(Mat F,PetscReal dtol)
5205ccb76cbSHong Zhang {
5215ccb76cbSHong Zhang   PetscErrorCode ierr;
5225ccb76cbSHong Zhang 
5235ccb76cbSHong Zhang   PetscFunctionBegin;
5245ccb76cbSHong Zhang   PetscValidHeaderSpecific(F,MAT_CLASSID,1);
52569fbec6eSBarry Smith   PetscValidLogicalCollectiveReal(F,dtol,2);
5265ccb76cbSHong Zhang   ierr = PetscTryMethod(F,"MatSuperluSetILUDropTol_C",(Mat,PetscReal),(F,dtol));CHKERRQ(ierr);
5275ccb76cbSHong Zhang   PetscFunctionReturn(0);
5285ccb76cbSHong Zhang }
5295ccb76cbSHong Zhang 
5305ccb76cbSHong Zhang #undef __FUNCT__
53135bd34faSBarry Smith #define __FUNCT__ "MatFactorGetSolverPackage_seqaij_superlu"
53235bd34faSBarry Smith PetscErrorCode MatFactorGetSolverPackage_seqaij_superlu(Mat A,const MatSolverPackage *type)
53335bd34faSBarry Smith {
53435bd34faSBarry Smith   PetscFunctionBegin;
5352692d6eeSBarry Smith   *type = MATSOLVERSUPERLU;
53635bd34faSBarry Smith   PetscFunctionReturn(0);
53735bd34faSBarry Smith }
53835bd34faSBarry Smith 
539b24902e0SBarry Smith /*MC
540ba992d64SSatish Balay   MATSOLVERSUPERLU = "superlu" - A solver package providing solvers LU and ILU for sequential matrices
541b24902e0SBarry Smith   via the external package SuperLU.
542b24902e0SBarry Smith 
543e2e64c6bSBarry Smith   Use ./configure --download-superlu to have PETSc installed with SuperLU
544b24902e0SBarry Smith 
545c2b89b5dSBarry Smith   Use -pc_type lu -pc_factor_mat_solver_package superlu to us this direct solver
546c2b89b5dSBarry Smith 
547b24902e0SBarry Smith   Options Database Keys:
548e08999f5SMatthew G Knepley + -mat_superlu_equil <FALSE>            - Equil (None)
549e08999f5SMatthew G Knepley . -mat_superlu_colperm <COLAMD>         - (choose one of) NATURAL MMD_ATA MMD_AT_PLUS_A COLAMD
550e08999f5SMatthew G Knepley . -mat_superlu_iterrefine <NOREFINE>    - (choose one of) NOREFINE SINGLE DOUBLE EXTRA
551e08999f5SMatthew G Knepley . -mat_superlu_symmetricmode: <FALSE>   - SymmetricMode (None)
552e08999f5SMatthew G Knepley . -mat_superlu_diagpivotthresh <1>      - DiagPivotThresh (None)
553e08999f5SMatthew G Knepley . -mat_superlu_pivotgrowth <FALSE>      - PivotGrowth (None)
554e08999f5SMatthew G Knepley . -mat_superlu_conditionnumber <FALSE>  - ConditionNumber (None)
555e08999f5SMatthew G Knepley . -mat_superlu_rowperm <NOROWPERM>      - (choose one of) NOROWPERM LargeDiag
556e08999f5SMatthew G Knepley . -mat_superlu_replacetinypivot <FALSE> - ReplaceTinyPivot (None)
557e08999f5SMatthew G Knepley . -mat_superlu_printstat <FALSE>        - PrintStat (None)
558e08999f5SMatthew G Knepley . -mat_superlu_lwork <0>                - size of work array in bytes used by factorization (None)
559e08999f5SMatthew G Knepley . -mat_superlu_ilu_droptol <0>          - ILU_DropTol (None)
560e08999f5SMatthew G Knepley . -mat_superlu_ilu_filltol <0>          - ILU_FillTol (None)
561e08999f5SMatthew G Knepley . -mat_superlu_ilu_fillfactor <0>       - ILU_FillFactor (None)
562e08999f5SMatthew G Knepley . -mat_superlu_ilu_droprull <0>         - ILU_DropRule (None)
563e08999f5SMatthew G Knepley . -mat_superlu_ilu_norm <0>             - ILU_Norm (None)
564e08999f5SMatthew G Knepley - -mat_superlu_ilu_milu <0>             - ILU_MILU (None)
565b24902e0SBarry Smith 
5662692d6eeSBarry Smith    Notes: Do not confuse this with MATSOLVERSUPERLU_DIST which is for parallel sparse solves
5675c9eb25fSBarry Smith 
568b24902e0SBarry Smith    Level: beginner
569b24902e0SBarry Smith 
570d45987f3SHong Zhang .seealso: PCLU, PCILU, MATSOLVERSUPERLU_DIST, MATSOLVERMUMPS, PCFactorSetMatSolverPackage(), MatSolverPackage
571b24902e0SBarry Smith M*/
572b24902e0SBarry Smith 
573b24902e0SBarry Smith #undef __FUNCT__
574b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqaij_superlu"
575db87b0f2SBarry Smith static PetscErrorCode MatGetFactor_seqaij_superlu(Mat A,MatFactorType ftype,Mat *F)
576b24902e0SBarry Smith {
57714b396bbSKris Buschelman   Mat            B;
578f0c56d0fSKris Buschelman   Mat_SuperLU    *lu;
5796849ba73SBarry Smith   PetscErrorCode ierr;
5805d8b2efaSHong Zhang   PetscInt       indx,m=A->rmap->n,n=A->cmap->n;
5818afaa268SBarry Smith   PetscBool      flg,set;
5823cb270beSHong Zhang   PetscReal      real_input;
5835ca28756SHong Zhang   const char     *colperm[]   ={"NATURAL","MMD_ATA","MMD_AT_PLUS_A","COLAMD"}; /* MY_PERMC - not supported by the petsc interface yet */
5845ca28756SHong Zhang   const char     *iterrefine[]={"NOREFINE", "SINGLE", "DOUBLE", "EXTRA"};
5855ca28756SHong Zhang   const char     *rowperm[]   ={"NOROWPERM", "LargeDiag"}; /* MY_PERMC - not supported by the petsc interface yet */
58614b396bbSKris Buschelman 
58714b396bbSKris Buschelman   PetscFunctionBegin;
588ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr);
589d0f46423SBarry Smith   ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
590245fece9SBarry Smith   ierr = PetscStrallocpy("superlu",&((PetscObject)B)->type_name);CHKERRQ(ierr);
591245fece9SBarry Smith   ierr = MatSetUp(B);CHKERRQ(ierr);
592cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) {
593b24902e0SBarry Smith     B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SuperLU;
594cffbb591SHong Zhang     B->ops->ilufactorsymbolic = MatLUFactorSymbolic_SuperLU;
595b3a44c85SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
596cffbb591SHong Zhang 
59700c67f3bSHong Zhang   ierr = PetscFree(B->solvertype);CHKERRQ(ierr);
59800c67f3bSHong Zhang   ierr = PetscStrallocpy(MATSOLVERSUPERLU,&B->solvertype);CHKERRQ(ierr);
59900c67f3bSHong Zhang 
600b9c12af5SBarry Smith   B->ops->getinfo     = MatGetInfo_External;
601b24902e0SBarry Smith   B->ops->destroy     = MatDestroy_SuperLU;
6023519fcdcSHong Zhang   B->ops->view        = MatView_SuperLU;
603d5f3da31SBarry Smith   B->factortype       = ftype;
60494b7f48cSBarry Smith   B->assembled        = PETSC_TRUE;           /* required by -ksp_view */
6055c9eb25fSBarry Smith   B->preallocated     = PETSC_TRUE;
60614b396bbSKris Buschelman 
607b00a9115SJed Brown   ierr = PetscNewLog(B,&lu);CHKERRQ(ierr);
608cae5a23dSHong Zhang 
609cffbb591SHong Zhang   if (ftype == MAT_FACTOR_LU) {
6109ce50919SHong Zhang     set_default_options(&lu->options);
6113d6581fbSHong Zhang     /* Comments from SuperLU_4.0/SRC/dgssvx.c:
6123d6581fbSHong Zhang       "Whether or not the system will be equilibrated depends on the
6133d6581fbSHong Zhang        scaling of the matrix A, but if equilibration is used, A is
6143d6581fbSHong Zhang        overwritten by diag(R)*A*diag(C) and B by diag(R)*B
6153d6581fbSHong Zhang        (if options->Trans=NOTRANS) or diag(C)*B (if options->Trans = TRANS or CONJ)."
6163d6581fbSHong Zhang      We set 'options.Equil = NO' as default because additional space is needed for it.
6173d6581fbSHong Zhang     */
6183d6581fbSHong Zhang     lu->options.Equil = NO;
619cffbb591SHong Zhang   } else if (ftype == MAT_FACTOR_ILU) {
6200c74a584SJed Brown     /* Set the default input options of ilu: */
621d387c056SBarry Smith     PetscStackCall("SuperLU:ilu_set_default_options",ilu_set_default_options(&lu->options));
622cffbb591SHong Zhang   }
6239ce50919SHong Zhang   lu->options.PrintStat = NO;
6241a2e2f44SHong Zhang 
6255d8b2efaSHong Zhang   /* Initialize the statistics variables. */
626d387c056SBarry Smith   PetscStackCall("SuperLU:StatInit",StatInit(&lu->stat));
6278914a3f7SHong Zhang   lu->lwork = 0;   /* allocate space internally by system malloc */
6289ce50919SHong Zhang 
629ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"SuperLU Options","Mat");CHKERRQ(ierr);
6308afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_superlu_equil","Equil","None",(PetscBool)lu->options.Equil,(PetscBool*)&lu->options.Equil,NULL);CHKERRQ(ierr);
6318914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_colperm","ColPerm","None",colperm,4,colperm[3],&indx,&flg);CHKERRQ(ierr);
6322205254eSKarl Rupp   if (flg) lu->options.ColPerm = (colperm_t)indx;
6338914a3f7SHong Zhang   ierr = PetscOptionsEList("-mat_superlu_iterrefine","IterRefine","None",iterrefine,4,iterrefine[0],&indx,&flg);CHKERRQ(ierr);
6342205254eSKarl Rupp   if (flg) lu->options.IterRefine = (IterRefine_t)indx;
6358afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_superlu_symmetricmode","SymmetricMode","None",(PetscBool)lu->options.SymmetricMode,&flg,&set);CHKERRQ(ierr);
6368afaa268SBarry Smith   if (set && flg) lu->options.SymmetricMode = YES;
6373cb270beSHong Zhang   ierr = PetscOptionsReal("-mat_superlu_diagpivotthresh","DiagPivotThresh","None",lu->options.DiagPivotThresh,&real_input,&flg);CHKERRQ(ierr);
6383cb270beSHong Zhang   if (flg) lu->options.DiagPivotThresh = (double) real_input;
6398afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_superlu_pivotgrowth","PivotGrowth","None",(PetscBool)lu->options.PivotGrowth,&flg,&set);CHKERRQ(ierr);
6408afaa268SBarry Smith   if (set && flg) lu->options.PivotGrowth = YES;
6418afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_superlu_conditionnumber","ConditionNumber","None",(PetscBool)lu->options.ConditionNumber,&flg,&set);CHKERRQ(ierr);
6428afaa268SBarry Smith   if (set && flg) lu->options.ConditionNumber = YES;
643d7ebd59bSHong Zhang   ierr = PetscOptionsEList("-mat_superlu_rowperm","rowperm","None",rowperm,2,rowperm[lu->options.RowPerm],&indx,&flg);CHKERRQ(ierr);
6442205254eSKarl Rupp   if (flg) lu->options.RowPerm = (rowperm_t)indx;
6458afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_superlu_replacetinypivot","ReplaceTinyPivot","None",(PetscBool)lu->options.ReplaceTinyPivot,&flg,&set);CHKERRQ(ierr);
6468afaa268SBarry Smith   if (set && flg) lu->options.ReplaceTinyPivot = YES;
6478afaa268SBarry Smith   ierr = PetscOptionsBool("-mat_superlu_printstat","PrintStat","None",(PetscBool)lu->options.PrintStat,&flg,&set);CHKERRQ(ierr);
6488afaa268SBarry Smith   if (set && flg) lu->options.PrintStat = YES;
6490298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_superlu_lwork","size of work array in bytes used by factorization","None",lu->lwork,&lu->lwork,NULL);CHKERRQ(ierr);
6505fe6150dSHong Zhang   if (lu->lwork > 0) {
651d87de817SBarry Smith     /* lwork is in bytes, hence PetscMalloc() is used here, not PetscMalloc1()*/
6525fe6150dSHong Zhang     ierr = PetscMalloc(lu->lwork,&lu->work);CHKERRQ(ierr);
6535fe6150dSHong Zhang   } else if (lu->lwork != 0 && lu->lwork != -1) {
65477431f27SBarry Smith     ierr      = PetscPrintf(PETSC_COMM_SELF,"   Warning: lwork %D is not supported by SUPERLU. The default lwork=0 is used.\n",lu->lwork);
6558914a3f7SHong Zhang     lu->lwork = 0;
6568914a3f7SHong Zhang   }
657cffbb591SHong Zhang   /* ilu options */
6583cb270beSHong Zhang   ierr = PetscOptionsReal("-mat_superlu_ilu_droptol","ILU_DropTol","None",lu->options.ILU_DropTol,&real_input,&flg);CHKERRQ(ierr);
6593cb270beSHong Zhang   if (flg) lu->options.ILU_DropTol = (double) real_input;
6603cb270beSHong Zhang   ierr = PetscOptionsReal("-mat_superlu_ilu_filltol","ILU_FillTol","None",lu->options.ILU_FillTol,&real_input,&flg);CHKERRQ(ierr);
6613cb270beSHong Zhang   if (flg) lu->options.ILU_FillTol = (double) real_input;
6623cb270beSHong Zhang   ierr = PetscOptionsReal("-mat_superlu_ilu_fillfactor","ILU_FillFactor","None",lu->options.ILU_FillFactor,&real_input,&flg);CHKERRQ(ierr);
6633cb270beSHong Zhang   if (flg) lu->options.ILU_FillFactor = (double) real_input;
6640298fd71SBarry Smith   ierr = PetscOptionsInt("-mat_superlu_ilu_droprull","ILU_DropRule","None",lu->options.ILU_DropRule,&lu->options.ILU_DropRule,NULL);CHKERRQ(ierr);
665cffbb591SHong Zhang   ierr = PetscOptionsInt("-mat_superlu_ilu_norm","ILU_Norm","None",lu->options.ILU_Norm,&indx,&flg);CHKERRQ(ierr);
6662205254eSKarl Rupp   if (flg) lu->options.ILU_Norm = (norm_t)indx;
667cffbb591SHong Zhang   ierr = PetscOptionsInt("-mat_superlu_ilu_milu","ILU_MILU","None",lu->options.ILU_MILU,&indx,&flg);CHKERRQ(ierr);
6682205254eSKarl Rupp   if (flg) lu->options.ILU_MILU = (milu_t)indx;
669245fece9SBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
67038abfdaeSHong Zhang   if (lu->options.Equil == YES) {
67138abfdaeSHong Zhang     /* superlu overwrites input matrix and rhs when Equil is used, thus create A_dup to keep user's A unchanged */
67238abfdaeSHong Zhang     ierr = MatDuplicate_SeqAIJ(A,MAT_COPY_VALUES,&lu->A_dup);CHKERRQ(ierr);
67338abfdaeSHong Zhang   }
6749ce50919SHong Zhang 
6755d8b2efaSHong Zhang   /* Allocate spaces (notice sizes are for the transpose) */
676785e854fSJed Brown   ierr = PetscMalloc1(m,&lu->etree);CHKERRQ(ierr);
677785e854fSJed Brown   ierr = PetscMalloc1(n,&lu->perm_r);CHKERRQ(ierr);
678785e854fSJed Brown   ierr = PetscMalloc1(m,&lu->perm_c);CHKERRQ(ierr);
679785e854fSJed Brown   ierr = PetscMalloc1(n,&lu->R);CHKERRQ(ierr);
680785e854fSJed Brown   ierr = PetscMalloc1(m,&lu->C);CHKERRQ(ierr);
6815d8b2efaSHong Zhang 
6825d8b2efaSHong Zhang   /* create rhs and solution x without allocate space for .Store */
6835d8b2efaSHong Zhang #if defined(PETSC_USE_COMPLEX)
6843cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
685d387c056SBarry Smith   PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
686d387c056SBarry Smith   PetscStackCall("SuperLU:cCreate_Dense_Matrix(",cCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_C, SLU_GE));
6873cb270beSHong Zhang #else
688d387c056SBarry Smith   PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
689d387c056SBarry Smith   PetscStackCall("SuperLU:zCreate_Dense_Matrix",zCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_Z, SLU_GE));
6903cb270beSHong Zhang #endif
6913cb270beSHong Zhang #else
6923cb270beSHong Zhang #if defined(PETSC_USE_REAL_SINGLE)
693d387c056SBarry Smith   PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
694d387c056SBarry Smith   PetscStackCall("SuperLU:sCreate_Dense_Matrix",sCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_S, SLU_GE));
6955d8b2efaSHong Zhang #else
696d387c056SBarry Smith   PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->B, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
697d387c056SBarry Smith   PetscStackCall("SuperLU:dCreate_Dense_Matrix",dCreate_Dense_Matrix(&lu->X, m, 1, NULL, m, SLU_DN, SLU_D, SLU_GE));
6985d8b2efaSHong Zhang #endif
6993cb270beSHong Zhang #endif
7005d8b2efaSHong Zhang 
701bdf89e91SBarry Smith   ierr     = PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_seqaij_superlu);CHKERRQ(ierr);
702bdf89e91SBarry Smith   ierr     = PetscObjectComposeFunction((PetscObject)B,"MatSuperluSetILUDropTol_C",MatSuperluSetILUDropTol_SuperLU);CHKERRQ(ierr);
703245fece9SBarry Smith   B->data  = lu;
70420be8e61SHong Zhang 
705899d7b4fSKris Buschelman   *F       = B;
70614b396bbSKris Buschelman   PetscFunctionReturn(0);
70714b396bbSKris Buschelman }
708d954db57SHong Zhang 
70942c9c57cSBarry Smith #undef __FUNCT__
71042c9c57cSBarry Smith #define __FUNCT__ "MatSolverPackageRegister_SuperLU"
71129b38603SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_SuperLU(void)
71242c9c57cSBarry Smith {
71342c9c57cSBarry Smith   PetscErrorCode ierr;
71442c9c57cSBarry Smith 
71542c9c57cSBarry Smith   PetscFunctionBegin;
71642c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ,       MAT_FACTOR_LU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr);
71742c9c57cSBarry Smith   ierr = MatSolverPackageRegister(MATSOLVERSUPERLU,MATSEQAIJ,       MAT_FACTOR_ILU,MatGetFactor_seqaij_superlu);CHKERRQ(ierr);
71842c9c57cSBarry Smith   PetscFunctionReturn(0);
71942c9c57cSBarry Smith }
720