Actual source code: superlu.c
1: /*$Id: superlu.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
3: /*
4: Provides an interface to the SuperLU sparse solver
5: Modified for SuperLU 2.0 by Matthew Knepley
6: */
8: #include src/mat/impls/aij/seq/aij.h
10: EXTERN_C_BEGIN
11: #if defined(PETSC_USE_COMPLEX)
12: #include "zsp_defs.h"
13: #else
14: #include "dsp_defs.h"
15: #endif
16: #include "util.h"
17: EXTERN_C_END
19: typedef struct {
20: SuperMatrix A,B,AC,L,U;
21: int *perm_r,*perm_c,ispec,relax,panel_size;
22: double pivot_threshold;
23: NCformat *store;
24: MatStructure flg;
25: PetscTruth SuperluMatOdering;
27: /* A few function pointers for inheritance */
28: int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
29: int (*MatView)(Mat,PetscViewer);
30: int (*MatAssemblyEnd)(Mat,MatAssemblyType);
31: int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
32: int (*MatDestroy)(Mat);
34: /* Flag to clean up (non-global) SuperLU objects during Destroy */
35: PetscTruth CleanUpSuperLU;
36: } Mat_SuperLU;
39: EXTERN int MatFactorInfo_SuperLU(Mat,PetscViewer);
40: EXTERN int MatLUFactorSymbolic_SuperLU(Mat,IS,IS,MatFactorInfo*,Mat*);
42: EXTERN_C_BEGIN
43: EXTERN int MatConvert_SuperLU_SeqAIJ(Mat,MatType,Mat*);
44: EXTERN int MatConvert_SeqAIJ_SuperLU(Mat,MatType,Mat*);
45: EXTERN_C_END
49: int MatDestroy_SuperLU(Mat A)
50: {
51: int ierr;
52: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
55: if (lu->CleanUpSuperLU) {
56: /* We have to free the global data or SuperLU crashes (sucky design)*/
57: /* Since we don't know if more solves on other matrices may be done
58: we cannot free the yucky SuperLU global data
59: StatFree();
60: */
61:
62: /* Free the SuperLU datastructures */
63: Destroy_CompCol_Permuted(&lu->AC);
64: Destroy_SuperNode_Matrix(&lu->L);
65: Destroy_CompCol_Matrix(&lu->U);
66: PetscFree(lu->B.Store);
67: PetscFree(lu->A.Store);
68: PetscFree(lu->perm_r);
69: PetscFree(lu->perm_c);
70: }
71: MatConvert_SuperLU_SeqAIJ(A,MATSEQAIJ,&A);
72: (*A->ops->destroy)(A);
73: return(0);
74: }
78: int MatView_SuperLU(Mat A,PetscViewer viewer)
79: {
80: int ierr;
81: PetscTruth isascii;
82: PetscViewerFormat format;
83: Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr);
86: (*lu->MatView)(A,viewer);
88: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
89: if (isascii) {
90: PetscViewerGetFormat(viewer,&format);
91: if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
92: MatFactorInfo_SuperLU(A,viewer);
93: }
94: }
95: return(0);
96: }
100: int MatAssemblyEnd_SuperLU(Mat A,MatAssemblyType mode) {
101: int ierr;
102: Mat_SuperLU *lu=(Mat_SuperLU*)(A->spptr);
105: (*lu->MatAssemblyEnd)(A,mode);
107: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
108: A->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
109: return(0);
110: }
112: #include src/mat/impls/dense/seq/dense.h
115: int MatCreateNull_SuperLU(Mat A,Mat *nullMat)
116: {
117: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
118: int numRows = A->m,numCols = A->n;
119: SCformat *Lstore;
120: int numNullCols,size;
121: #if defined(PETSC_USE_COMPLEX)
122: doublecomplex *nullVals,*workVals;
123: #else
124: PetscScalar *nullVals,*workVals;
125: #endif
126: int row,newRow,col,newCol,block,b,ierr;
131: if (!A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Unfactored matrix");
132: numNullCols = numCols - numRows;
133: if (numNullCols < 0) SETERRQ(PETSC_ERR_ARG_WRONG,"Function only applies to underdetermined problems");
134: /* Create the null matrix */
135: MatCreateSeqDense(A->comm,numRows,numNullCols,PETSC_NULL,nullMat);
136: if (numNullCols == 0) {
137: MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);
138: MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);
139: return(0);
140: }
141: #if defined(PETSC_USE_COMPLEX)
142: nullVals = (doublecomplex*)((Mat_SeqDense*)(*nullMat)->data)->v;
143: #else
144: nullVals = ((Mat_SeqDense*)(*nullMat)->data)->v;
145: #endif
146: /* Copy in the columns */
147: Lstore = (SCformat*)lu->L.Store;
148: for(block = 0; block <= Lstore->nsuper; block++) {
149: newRow = Lstore->sup_to_col[block];
150: size = Lstore->sup_to_col[block+1] - Lstore->sup_to_col[block];
151: for(col = Lstore->rowind_colptr[newRow]; col < Lstore->rowind_colptr[newRow+1]; col++) {
152: newCol = Lstore->rowind[col];
153: if (newCol >= numRows) {
154: for(b = 0; b < size; b++)
155: #if defined(PETSC_USE_COMPLEX)
156: nullVals[(newCol-numRows)*numRows+newRow+b] = ((doublecomplex*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
157: #else
158: nullVals[(newCol-numRows)*numRows+newRow+b] = ((double*)Lstore->nzval)[Lstore->nzval_colptr[newRow+b]+col];
159: #endif
160: }
161: }
162: }
163: /* Permute rhs to form P^T_c B */
164: PetscMalloc(numRows*sizeof(double),&workVals);
165: for(b = 0; b < numNullCols; b++) {
166: for(row = 0; row < numRows; row++) workVals[lu->perm_c[row]] = nullVals[b*numRows+row];
167: for(row = 0; row < numRows; row++) nullVals[b*numRows+row] = workVals[row];
168: }
169: /* Backward solve the upper triangle A x = b */
170: for(b = 0; b < numNullCols; b++) {
171: #if defined(PETSC_USE_COMPLEX)
172: sp_ztrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr);
173: #else
174: sp_dtrsv("L","T","U",&lu->L,&lu->U,&nullVals[b*numRows],&ierr);
175: #endif
176: if (ierr < 0)
177: SETERRQ1(PETSC_ERR_ARG_WRONG,"The argument %d was invalid",-ierr);
178: }
179: PetscFree(workVals);
181: MatAssemblyBegin(*nullMat,MAT_FINAL_ASSEMBLY);
182: MatAssemblyEnd(*nullMat,MAT_FINAL_ASSEMBLY);
183: return(0);
184: }
188: int MatSolve_SuperLU(Mat A,Vec b,Vec x)
189: {
190: Mat_SuperLU *lu = (Mat_SuperLU*)A->spptr;
191: PetscScalar *array;
192: int m,ierr;
195: VecGetLocalSize(b,&m);
196: VecCopy(b,x);
197: VecGetArray(x,&array);
198: /* Create the Rhs */
199: lu->B.Stype = SLU_DN;
200: lu->B.Mtype = SLU_GE;
201: lu->B.nrow = m;
202: lu->B.ncol = 1;
203: ((DNformat*)lu->B.Store)->lda = m;
204: ((DNformat*)lu->B.Store)->nzval = array;
205: #if defined(PETSC_USE_COMPLEX)
206: lu->B.Dtype = SLU_Z;
207: zgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr);
208: #else
209: lu->B.Dtype = SLU_D;
210: dgstrs("T",&lu->L,&lu->U,lu->perm_r,lu->perm_c,&lu->B,&ierr);
211: #endif
212: if (ierr < 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr);
213: VecRestoreArray(x,&array);
214: return(0);
215: }
217: static int StatInitCalled = 0;
221: int MatLUFactorNumeric_SuperLU(Mat A,Mat *F)
222: {
223: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data;
224: Mat_SuperLU *lu = (Mat_SuperLU*)(*F)->spptr;
225: int *etree,ierr;
226: PetscTruth flag;
229: /* Create the SuperMatrix for A^T:
230: Since SuperLU only likes column-oriented matrices,we pass it the transpose,
231: and then solve A^T X = B in MatSolve().
232: */
233:
234: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numerical factorization */
235: lu->A.Stype = SLU_NC;
236: #if defined(PETSC_USE_COMPLEX)
237: lu->A.Dtype = SLU_Z;
238: #else
239: lu->A.Dtype = SLU_D;
240: #endif
241: lu->A.Mtype = SLU_GE;
242: lu->A.nrow = A->n;
243: lu->A.ncol = A->m;
244:
245: PetscMalloc(sizeof(NCformat),&lu->store);
246: PetscMalloc(sizeof(DNformat),&lu->B.Store);
247: }
248: lu->store->nnz = aa->nz;
249: lu->store->colptr = aa->i;
250: lu->store->rowind = aa->j;
251: lu->store->nzval = aa->a;
252: lu->A.Store = lu->store;
253:
254: /* Set SuperLU options */
255: lu->relax = sp_ienv(2);
256: lu->panel_size = sp_ienv(1);
257: /* We have to initialize global data or SuperLU crashes (sucky design) */
258: if (!StatInitCalled) {
259: StatInit(lu->panel_size,lu->relax);
260: }
261: StatInitCalled++;
263: PetscOptionsBegin(A->comm,A->prefix,"SuperLU Options","Mat");
264: /* use SuperLU mat ordeing */
265: PetscOptionsInt("-mat_superlu_ordering","SuperLU ordering type (one of 0, 1, 2, 3)\n 0: natural ordering;\n 1: MMD applied to A'*A;\n 2: MMD applied to A'+A;\n 3: COLAMD, approximate minimum degree column ordering","None",lu->ispec,&lu->ispec,&flag);
266: if (flag) {
267: get_perm_c(lu->ispec, &lu->A, lu->perm_c);
268: lu->SuperluMatOdering = PETSC_TRUE;
269: }
270: PetscOptionsEnd();
272: /* Create the elimination tree */
273: PetscMalloc(A->n*sizeof(int),&etree);
274: sp_preorder("N",&lu->A,lu->perm_c,etree,&lu->AC);
275: /* Factor the matrix */
276: #if defined(PETSC_USE_COMPLEX)
277: zgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr);
278: #else
279: dgstrf("N",&lu->AC,lu->pivot_threshold,0.0,lu->relax,lu->panel_size,etree,PETSC_NULL,0,lu->perm_r,lu->perm_c,&lu->L,&lu->U,&ierr);
280: #endif
281: if (ierr < 0) {
282: SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element of row %d was invalid",-ierr);
283: } else if (ierr > 0) {
284: if (ierr <= A->m) {
285: SETERRQ1(PETSC_ERR_ARG_WRONG,"The diagonal element %d of U is exactly zero",ierr);
286: } else {
287: SETERRQ1(PETSC_ERR_ARG_WRONG,"Memory allocation failure after %d bytes were allocated",ierr-A->m);
288: }
289: }
291: /* Cleanup */
292: PetscFree(etree);
294: lu->flg = SAME_NONZERO_PATTERN;
295: return(0);
296: }
298: /*
299: Note the r permutation is ignored
300: */
303: int MatLUFactorSymbolic_SuperLU(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
304: {
305: Mat B;
306: Mat_SuperLU *lu;
307: int ierr,*ca;
310:
311: MatCreate(A->comm,A->m,A->n,PETSC_DETERMINE,PETSC_DETERMINE,&B);
312: MatSetType(B,MATSUPERLU);
313: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
315: B->ops->lufactornumeric = MatLUFactorNumeric_SuperLU;
316: B->ops->solve = MatSolve_SuperLU;
317: B->factor = FACTOR_LU;
318: B->assembled = PETSC_TRUE; /* required by -sles_view */
319:
320: lu = (Mat_SuperLU*)(B->spptr);
321: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatCreateNull","MatCreateNull_SuperLU",
322: (void(*)(void))MatCreateNull_SuperLU);
324: /* Allocate the work arrays required by SuperLU (notice sizes are for the transpose) */
325: PetscMalloc(A->n*sizeof(int),&lu->perm_r);
326: PetscMalloc(A->m*sizeof(int),&lu->perm_c);
328: /* use PETSc mat ordering */
329: ISGetIndices(c,&ca);
330: PetscMemcpy(lu->perm_c,ca,A->m*sizeof(int));
331: ISRestoreIndices(c,&ca);
332: lu->SuperluMatOdering = PETSC_FALSE;
334: lu->pivot_threshold = info->dtcol;
335: PetscLogObjectMemory(B,(A->m+A->n)*sizeof(int)+sizeof(Mat_SuperLU));
337: lu->flg = DIFFERENT_NONZERO_PATTERN;
338: lu->CleanUpSuperLU = PETSC_TRUE;
340: *F = B;
341: return(0);
342: }
344: /* used by -sles_view */
347: int MatFactorInfo_SuperLU(Mat A,PetscViewer viewer)
348: {
349: Mat_SuperLU *lu= (Mat_SuperLU*)A->spptr;
350: int ierr;
353: PetscViewerASCIIPrintf(viewer,"SuperLU run parameters:\n");
354: if(lu->SuperluMatOdering) {
355: PetscViewerASCIIPrintf(viewer," SuperLU mat ordering: %d\n",lu->ispec);
356: }
357: return(0);
358: }
362: int MatDuplicate_SuperLU(Mat A, MatDuplicateOption op, Mat *M) {
365: (*A->ops->duplicate)(A,op,M);
366: MatConvert_SeqAIJ_SuperLU(*M,MATSUPERLU,M);
367: return(0);
368: }
370: EXTERN_C_BEGIN
373: int MatConvert_SuperLU_SeqAIJ(Mat A,MatType type,Mat *newmat) {
374: /* This routine is only called to convert an unfactored PETSc-SuperLU matrix */
375: /* to its base PETSc type, so we will ignore 'MatType type'. */
376: int ierr;
377: Mat B=*newmat;
378: Mat_SuperLU *lu=(Mat_SuperLU *)A->spptr;
381: if (B != A) {
382: MatDuplicate(A,MAT_COPY_VALUES,&B);
383: }
384: /* Reset the original function pointers */
385: B->ops->duplicate = lu->MatDuplicate;
386: B->ops->view = lu->MatView;
387: B->ops->assemblyend = lu->MatAssemblyEnd;
388: B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
389: B->ops->destroy = lu->MatDestroy;
390: /* lu is only a function pointer stash unless we've factored the matrix, which we haven't! */
391: PetscFree(lu);
392: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
393: *newmat = B;
394: return(0);
395: }
396: EXTERN_C_END
398: EXTERN_C_BEGIN
401: int MatConvert_SeqAIJ_SuperLU(Mat A,MatType type,Mat *newmat) {
402: /* This routine is only called to convert to MATSUPERLU */
403: /* from MATSEQAIJ, so we will ignore 'MatType type'. */
404: int ierr;
405: Mat B=*newmat;
406: Mat_SuperLU *lu;
409: if (B != A) {
410: MatDuplicate(A,MAT_COPY_VALUES,&B);
411: }
413: PetscNew(Mat_SuperLU,&lu);
414: lu->MatDuplicate = A->ops->duplicate;
415: lu->MatView = A->ops->view;
416: lu->MatAssemblyEnd = A->ops->assemblyend;
417: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
418: lu->MatDestroy = A->ops->destroy;
419: lu->CleanUpSuperLU = PETSC_FALSE;
421: B->spptr = (void*)lu;
422: B->ops->duplicate = MatDuplicate_SuperLU;
423: B->ops->view = MatView_SuperLU;
424: B->ops->assemblyend = MatAssemblyEnd_SuperLU;
425: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SuperLU;
426: B->ops->destroy = MatDestroy_SuperLU;
428: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_superlu_C",
429: "MatConvert_SeqAIJ_SuperLU",MatConvert_SeqAIJ_SuperLU);
430: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_superlu_seqaij_C",
431: "MatConvert_SuperLU_SeqAIJ",MatConvert_SuperLU_SeqAIJ);
432: PetscLogInfo(0,"Using SuperLU for SeqAIJ LU factorization and solves.");
433: PetscObjectChangeTypeName((PetscObject)B,MATSUPERLU);
434: *newmat = B;
435: return(0);
436: }
437: EXTERN_C_END
439: /*MC
440: MATSUPERLU - MATSUPERLU = "superlu" - A matrix type providing direct solvers (LU) for sequential matrices
441: via the external package SuperLU.
443: If SuperLU is installed (see the manual for
444: instructions on how to declare the existence of external packages),
445: a matrix type can be constructed which invokes SuperLU solvers.
446: After calling MatCreate(...,A), simply call MatSetType(A,MATSUPERLU).
447: This matrix type is only supported for double precision real.
449: This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is
450: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
451: the MATSEQAIJ type without data copy.
453: Options Database Keys:
454: + -mat_type superlu - sets the matrix type to "superlu" during a call to MatSetFromOptions()
455: - -mat_superlu_ordering <0,1,2,3> - 0: natural ordering,
456: 1: MMD applied to A'*A,
457: 2: MMD applied to A'+A,
458: 3: COLAMD, approximate minimum degree column ordering
460: Level: beginner
462: .seealso: PCLU
463: M*/
465: EXTERN_C_BEGIN
468: int MatCreate_SuperLU(Mat A) {
472: /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SUPERLU types */
473: PetscObjectChangeTypeName((PetscObject)A,MATSUPERLU);
474: MatSetType(A,MATSEQAIJ);
475: MatConvert_SeqAIJ_SuperLU(A,MATSUPERLU,&A);
476: return(0);
477: }
478: EXTERN_C_END