Actual source code: spooles.c
1: /*$Id: spooles.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2: /*
3: Provides an interface to the Spooles serial sparse solver
4: */
5: #include src/mat/impls/aij/seq/aij.h
6: #include src/mat/impls/sbaij/seq/sbaij.h
7: #include src/mat/impls/aij/seq/spooles/spooles.h
9: EXTERN_C_BEGIN
12: int MatConvert_Spooles_Base(Mat A,const MatType type,Mat *newmat) {
13: /* This routine is only called to convert an unfactored PETSc-Spooles matrix */
14: /* to its base PETSc type, so we will ignore 'MatType type'. */
15: int ierr;
16: Mat B=*newmat;
17: Mat_Spooles *lu=(Mat_Spooles*)A->spptr;
19: if (B != A) {
20: MatDuplicate(A,MAT_COPY_VALUES,&B);
21: }
22: /* Reset the stashed function pointers set by inherited routines */
23: B->ops->duplicate = lu->MatDuplicate;
24: B->ops->choleskyfactorsymbolic = lu->MatCholeskyFactorSymbolic;
25: B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
26: B->ops->view = lu->MatView;
27: B->ops->assemblyend = lu->MatAssemblyEnd;
28: B->ops->destroy = lu->MatDestroy;
30: PetscObjectChangeTypeName((PetscObject)B,lu->basetype);
31: PetscFree(lu);
32: *newmat = B;
33: return(0);
34: }
35: EXTERN_C_END
39: int MatDestroy_SeqAIJSpooles(Mat A)
40: {
41: Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
42: int ierr;
43:
45:
46: if (lu->CleanUpSpooles) {
47: FrontMtx_free(lu->frontmtx) ;
48: IV_free(lu->newToOldIV) ;
49: IV_free(lu->oldToNewIV) ;
50: InpMtx_free(lu->mtxA) ;
51: ETree_free(lu->frontETree) ;
52: IVL_free(lu->symbfacIVL) ;
53: SubMtxManager_free(lu->mtxmanager) ;
54: Graph_free(lu->graph);
55: }
56: MatConvert_Spooles_Base(A,lu->basetype,&A);
57: (*A->ops->destroy)(A);
58: return(0);
59: }
63: int MatSolve_SeqAIJSpooles(Mat A,Vec b,Vec x)
64: {
65: Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
66: PetscScalar *array;
67: DenseMtx *mtxY, *mtxX ;
68: int ierr,irow,neqns=A->n,nrow=A->m,*iv;
69: #if defined(PETSC_USE_COMPLEX)
70: double x_real,x_imag;
71: #else
72: double *entX;
73: #endif
77: mtxY = DenseMtx_new() ;
78: DenseMtx_init(mtxY, lu->options.typeflag, 0, 0, nrow, 1, 1, nrow) ; /* column major */
79: VecGetArray(b,&array);
81: if (lu->options.useQR) { /* copy b to mtxY */
82: for ( irow = 0 ; irow < nrow; irow++ )
83: #if !defined(PETSC_USE_COMPLEX)
84: DenseMtx_setRealEntry(mtxY, irow, 0, *array++) ;
85: #else
86: DenseMtx_setComplexEntry(mtxY, irow, 0, PetscRealPart(array[irow]), PetscImaginaryPart(array[irow]));
87: #endif
88: } else { /* copy permuted b to mtxY */
89: iv = IV_entries(lu->oldToNewIV);
90: for ( irow = 0 ; irow < nrow; irow++ )
91: #if !defined(PETSC_USE_COMPLEX)
92: DenseMtx_setRealEntry(mtxY, *iv++, 0, *array++) ;
93: #else
94: DenseMtx_setComplexEntry(mtxY,*iv++,0,PetscRealPart(array[irow]),PetscImaginaryPart(array[irow]));
95: #endif
96: }
97: VecRestoreArray(b,&array);
99: mtxX = DenseMtx_new() ;
100: DenseMtx_init(mtxX, lu->options.typeflag, 0, 0, neqns, 1, 1, neqns) ;
101: if (lu->options.useQR) {
102: FrontMtx_QR_solve(lu->frontmtx, lu->mtxA, mtxX, mtxY, lu->mtxmanager,
103: lu->cpus, lu->options.msglvl, lu->options.msgFile) ;
104: } else {
105: FrontMtx_solve(lu->frontmtx, mtxX, mtxY, lu->mtxmanager,
106: lu->cpus, lu->options.msglvl, lu->options.msgFile) ;
107: }
108: if ( lu->options.msglvl > 2 ) {
109: fprintf(lu->options.msgFile, "\n\n right hand side matrix after permutation") ;
110: DenseMtx_writeForHumanEye(mtxY, lu->options.msgFile) ;
111: fprintf(lu->options.msgFile, "\n\n solution matrix in new ordering") ;
112: DenseMtx_writeForHumanEye(mtxX, lu->options.msgFile) ;
113: fflush(lu->options.msgFile) ;
114: }
116: /* permute solution into original ordering, then copy to x */
117: DenseMtx_permuteRows(mtxX, lu->newToOldIV);
118: VecGetArray(x,&array);
120: #if !defined(PETSC_USE_COMPLEX)
121: entX = DenseMtx_entries(mtxX);
122: DVcopy(neqns, array, entX);
123: #else
124: for (irow=0; irow<nrow; irow++){
125: DenseMtx_complexEntry(mtxX,irow,0,&x_real,&x_imag);
126: array[irow] = x_real+x_imag*PETSC_i;
127: }
128: #endif
130: VecRestoreArray(x,&array);
131:
132: /* free memory */
133: DenseMtx_free(mtxX) ;
134: DenseMtx_free(mtxY) ;
136: return(0);
137: }
141: int MatFactorNumeric_SeqAIJSpooles(Mat A,Mat *F)
142: {
143: Mat_Spooles *lu = (Mat_Spooles*)(*F)->spptr;
144: ChvManager *chvmanager ;
145: Chv *rootchv ;
146: IVL *adjIVL;
147: int ierr,nz,nrow=A->m,irow,nedges,neqns=A->n,*ai,*aj,i,*diag=0;
148: PetscScalar *av;
149: double cputotal,facops;
150: #if defined(PETSC_USE_COMPLEX)
151: int nz_row,*aj_tmp;
152: PetscScalar *av_tmp;
153: #else
154: int *ivec1,*ivec2,j;
155: double *dvec;
156: #endif
157: PetscTruth isAIJ;
158:
160: if (lu->flg == DIFFERENT_NONZERO_PATTERN) { /* first numeric factorization */
161: (*F)->ops->solve = MatSolve_SeqAIJSpooles;
162: (*F)->ops->destroy = MatDestroy_SeqAIJSpooles;
163: (*F)->assembled = PETSC_TRUE;
164:
165: /* set Spooles options */
166: SetSpoolesOptions(A, &lu->options);
168: lu->mtxA = InpMtx_new() ;
169: }
171: /* copy A to Spooles' InpMtx object */
172: PetscTypeCompare((PetscObject)A,MATSEQAIJSPOOLES,&isAIJ);
173: if (isAIJ){
174: Mat_SeqAIJ *mat = (Mat_SeqAIJ*)A->data;
175: ai=mat->i; aj=mat->j; av=mat->a;
176: if (lu->options.symflag == SPOOLES_NONSYMMETRIC) {
177: nz=mat->nz;
178: } else { /* SPOOLES_SYMMETRIC || SPOOLES_HERMITIAN */
179: nz=(mat->nz + A->m)/2;
180: if (!mat->diag){
181: MatMarkDiagonal_SeqAIJ(A);
182: }
183: diag=mat->diag;
184: }
185: } else { /* A is SBAIJ */
186: Mat_SeqSBAIJ *mat = (Mat_SeqSBAIJ*)A->data;
187: ai=mat->i; aj=mat->j; av=mat->a;
188: nz=mat->nz;
189: }
190: InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, lu->options.typeflag, nz, 0) ;
191:
192: #if defined(PETSC_USE_COMPLEX)
193: for (irow=0; irow<nrow; irow++) {
194: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
195: nz_row = ai[irow+1] - ai[irow];
196: aj_tmp = aj + ai[irow];
197: av_tmp = av + ai[irow];
198: } else {
199: nz_row = ai[irow+1] - diag[irow];
200: aj_tmp = aj + diag[irow];
201: av_tmp = av + diag[irow];
202: }
203: for (i=0; i<nz_row; i++){
204: InpMtx_inputComplexEntry(lu->mtxA, irow, *aj_tmp++,PetscRealPart(*av_tmp),PetscImaginaryPart(*av_tmp));
205: av_tmp++;
206: }
207: }
208: #else
209: ivec1 = InpMtx_ivec1(lu->mtxA);
210: ivec2 = InpMtx_ivec2(lu->mtxA);
211: dvec = InpMtx_dvec(lu->mtxA);
212: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
213: for (irow = 0; irow < nrow; irow++){
214: for (i = ai[irow]; i<ai[irow+1]; i++) ivec1[i] = irow;
215: }
216: IVcopy(nz, ivec2, aj);
217: DVcopy(nz, dvec, av);
218: } else {
219: nz = 0;
220: for (irow = 0; irow < nrow; irow++){
221: for (j = diag[irow]; j<ai[irow+1]; j++) {
222: ivec1[nz] = irow;
223: ivec2[nz] = aj[j];
224: dvec[nz] = av[j];
225: nz++;
226: }
227: }
228: }
229: InpMtx_inputRealTriples(lu->mtxA, nz, ivec1, ivec2, dvec);
230: #endif
232: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
233: if ( lu->options.msglvl > 0 ) {
234: printf("\n\n input matrix") ;
235: fprintf(lu->options.msgFile, "\n\n input matrix") ;
236: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
237: fflush(lu->options.msgFile) ;
238: }
240: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
241: /*---------------------------------------------------
242: find a low-fill ordering
243: (1) create the Graph object
244: (2) order the graph
245: -------------------------------------------------------*/
246: if (lu->options.useQR){
247: adjIVL = InpMtx_adjForATA(lu->mtxA) ;
248: } else {
249: adjIVL = InpMtx_fullAdjacency(lu->mtxA) ;
250: }
251: nedges = IVL_tsize(adjIVL) ;
253: lu->graph = Graph_new() ;
254: Graph_init2(lu->graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL) ;
255: if ( lu->options.msglvl > 2 ) {
256: if (lu->options.useQR){
257: fprintf(lu->options.msgFile, "\n\n graph of A^T A") ;
258: } else {
259: fprintf(lu->options.msgFile, "\n\n graph of the input matrix") ;
260: }
261: Graph_writeForHumanEye(lu->graph, lu->options.msgFile) ;
262: fflush(lu->options.msgFile) ;
263: }
265: switch (lu->options.ordering) {
266: case 0:
267: lu->frontETree = orderViaBestOfNDandMS(lu->graph,
268: lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize,
269: lu->options.seed, lu->options.msglvl, lu->options.msgFile); break;
270: case 1:
271: lu->frontETree = orderViaMMD(lu->graph,lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
272: case 2:
273: lu->frontETree = orderViaMS(lu->graph, lu->options.maxdomainsize,
274: lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
275: case 3:
276: lu->frontETree = orderViaND(lu->graph, lu->options.maxdomainsize,
277: lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
278: default:
279: SETERRQ(1,"Unknown Spooles's ordering");
280: }
282: if ( lu->options.msglvl > 0 ) {
283: fprintf(lu->options.msgFile, "\n\n front tree from ordering") ;
284: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
285: fflush(lu->options.msgFile) ;
286: }
287:
288: /* get the permutation, permute the front tree */
289: lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree) ;
290: lu->oldToNew = IV_entries(lu->oldToNewIV) ;
291: lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree) ;
292: if (!lu->options.useQR) ETree_permuteVertices(lu->frontETree, lu->oldToNewIV) ;
294: /* permute the matrix */
295: if (lu->options.useQR){
296: InpMtx_permute(lu->mtxA, NULL, lu->oldToNew) ;
297: } else {
298: InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew) ;
299: if ( lu->options.symflag == SPOOLES_SYMMETRIC) {
300: InpMtx_mapToUpperTriangle(lu->mtxA) ;
301: }
302: #if defined(PETSC_USE_COMPLEX)
303: if ( lu->options.symflag == SPOOLES_HERMITIAN ) {
304: InpMtx_mapToUpperTriangleH(lu->mtxA) ;
305: }
306: #endif
307: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
308: }
309: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
311: /* get symbolic factorization */
312: if (lu->options.useQR){
313: lu->symbfacIVL = SymbFac_initFromGraph(lu->frontETree, lu->graph) ;
314: IVL_overwrite(lu->symbfacIVL, lu->oldToNewIV) ;
315: IVL_sortUp(lu->symbfacIVL) ;
316: ETree_permuteVertices(lu->frontETree, lu->oldToNewIV) ;
317: } else {
318: lu->symbfacIVL = SymbFac_initFromInpMtx(lu->frontETree, lu->mtxA) ;
319: }
320: if ( lu->options.msglvl > 2 ) {
321: fprintf(lu->options.msgFile, "\n\n old-to-new permutation vector") ;
322: IV_writeForHumanEye(lu->oldToNewIV, lu->options.msgFile) ;
323: fprintf(lu->options.msgFile, "\n\n new-to-old permutation vector") ;
324: IV_writeForHumanEye(lu->newToOldIV, lu->options.msgFile) ;
325: fprintf(lu->options.msgFile, "\n\n front tree after permutation") ;
326: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
327: fprintf(lu->options.msgFile, "\n\n input matrix after permutation") ;
328: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
329: fprintf(lu->options.msgFile, "\n\n symbolic factorization") ;
330: IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile) ;
331: fflush(lu->options.msgFile) ;
332: }
334: lu->frontmtx = FrontMtx_new() ;
335: lu->mtxmanager = SubMtxManager_new() ;
336: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;
338: } else { /* new num factorization using previously computed symbolic factor */
340: if (lu->options.pivotingflag) { /* different FrontMtx is required */
341: FrontMtx_free(lu->frontmtx) ;
342: lu->frontmtx = FrontMtx_new() ;
343: } else {
344: FrontMtx_clearData (lu->frontmtx);
345: }
347: SubMtxManager_free(lu->mtxmanager) ;
348: lu->mtxmanager = SubMtxManager_new() ;
349: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;
351: /* permute mtxA */
352: if (lu->options.useQR){
353: InpMtx_permute(lu->mtxA, NULL, lu->oldToNew) ;
354: } else {
355: InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew) ;
356: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
357: InpMtx_mapToUpperTriangle(lu->mtxA) ;
358: }
359: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
360: }
361: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
362: if ( lu->options.msglvl > 2 ) {
363: fprintf(lu->options.msgFile, "\n\n input matrix after permutation") ;
364: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
365: }
366: } /* end of if( lu->flg == DIFFERENT_NONZERO_PATTERN) */
367:
368: if (lu->options.useQR){
369: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag,
370: SPOOLES_SYMMETRIC, FRONTMTX_DENSE_FRONTS,
371: SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL,
372: lu->mtxmanager, lu->options.msglvl, lu->options.msgFile) ;
373: } else {
374: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag, lu->options.symflag,
375: FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, 0, NULL,
376: lu->mtxmanager, lu->options.msglvl, lu->options.msgFile) ;
377: }
379: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
380: if ( lu->options.patchAndGoFlag == 1 ) {
381: lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
382: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
383: lu->options.storeids, lu->options.storevalues) ;
384: } else if ( lu->options.patchAndGoFlag == 2 ) {
385: lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
386: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
387: lu->options.storeids, lu->options.storevalues) ;
388: }
389: }
391: /* numerical factorization */
392: chvmanager = ChvManager_new() ;
393: ChvManager_init(chvmanager, NO_LOCK, 1) ;
394: DVfill(10, lu->cpus, 0.0) ;
395: if (lu->options.useQR){
396: facops = 0.0 ;
397: FrontMtx_QR_factor(lu->frontmtx, lu->mtxA, chvmanager,
398: lu->cpus, &facops, lu->options.msglvl, lu->options.msgFile) ;
399: if ( lu->options.msglvl > 1 ) {
400: fprintf(lu->options.msgFile, "\n\n factor matrix") ;
401: fprintf(lu->options.msgFile, "\n facops = %9.2f", facops) ;
402: }
403: } else {
404: IVfill(20, lu->stats, 0) ;
405: rootchv = FrontMtx_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, 0.0,
406: chvmanager, &ierr, lu->cpus,lu->stats,lu->options.msglvl,lu->options.msgFile) ;
407: if ( rootchv != NULL ) SETERRQ(1,"\n matrix found to be singular");
408: if ( ierr >= 0 ) SETERRQ1(1,"\n error encountered at front %d", ierr);
409:
410: if(lu->options.FrontMtxInfo){
411: PetscPrintf(PETSC_COMM_SELF,"\n %8d pivots, %8d pivot tests, %8d delayed rows and columns\n",\
412: lu->stats[0], lu->stats[1], lu->stats[2]);
413: cputotal = lu->cpus[8] ;
414: if ( cputotal > 0.0 ) {
415: PetscPrintf(PETSC_COMM_SELF,
416: "\n cpus cpus/totaltime"
417: "\n initialize fronts %8.3f %6.2f"
418: "\n load original entries %8.3f %6.2f"
419: "\n update fronts %8.3f %6.2f"
420: "\n assemble postponed data %8.3f %6.2f"
421: "\n factor fronts %8.3f %6.2f"
422: "\n extract postponed data %8.3f %6.2f"
423: "\n store factor entries %8.3f %6.2f"
424: "\n miscellaneous %8.3f %6.2f"
425: "\n total time %8.3f \n",
426: lu->cpus[0], 100.*lu->cpus[0]/cputotal,
427: lu->cpus[1], 100.*lu->cpus[1]/cputotal,
428: lu->cpus[2], 100.*lu->cpus[2]/cputotal,
429: lu->cpus[3], 100.*lu->cpus[3]/cputotal,
430: lu->cpus[4], 100.*lu->cpus[4]/cputotal,
431: lu->cpus[5], 100.*lu->cpus[5]/cputotal,
432: lu->cpus[6], 100.*lu->cpus[6]/cputotal,
433: lu->cpus[7], 100.*lu->cpus[7]/cputotal, cputotal) ;
434: }
435: }
436: }
437: ChvManager_free(chvmanager) ;
439: if ( lu->options.msglvl > 0 ) {
440: fprintf(lu->options.msgFile, "\n\n factor matrix") ;
441: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
442: fflush(lu->options.msgFile) ;
443: }
445: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
446: if ( lu->options.patchAndGoFlag == 1 ) {
447: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
448: if (lu->options.msglvl > 0 ){
449: fprintf(lu->options.msgFile, "\n small pivots found at these locations") ;
450: IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
451: }
452: }
453: PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
454: } else if ( lu->options.patchAndGoFlag == 2 ) {
455: if (lu->options.msglvl > 0 ){
456: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
457: fprintf(lu->options.msgFile, "\n small pivots found at these locations") ;
458: IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
459: }
460: if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
461: fprintf(lu->options.msgFile, "\n perturbations") ;
462: DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile) ;
463: }
464: }
465: PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
466: }
467: }
469: /* post-process the factorization */
470: FrontMtx_postProcess(lu->frontmtx, lu->options.msglvl, lu->options.msgFile) ;
471: if ( lu->options.msglvl > 2 ) {
472: fprintf(lu->options.msgFile, "\n\n factor matrix after post-processing") ;
473: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
474: fflush(lu->options.msgFile) ;
475: }
477: lu->flg = SAME_NONZERO_PATTERN;
478: lu->CleanUpSpooles = PETSC_TRUE;
479: return(0);
480: }
482: EXTERN_C_BEGIN
485: int MatConvert_SeqAIJ_SeqAIJSpooles(Mat A,const MatType type,Mat *newmat) {
486: /* This routine is only called to convert a MATSEQAIJ matrix */
487: /* to a MATSEQAIJSPOOLES matrix, so we will ignore 'MatType type'. */
488: int ierr;
489: Mat B=*newmat;
490: Mat_Spooles *lu;
493: if (B != A) {
494: /* This routine is inherited, so we know the type is correct. */
495: MatDuplicate(A,MAT_COPY_VALUES,&B);
496: }
497: PetscNew(Mat_Spooles,&lu);
498: B->spptr = (void*)lu;
500: lu->basetype = MATSEQAIJ;
501: lu->useQR = PETSC_FALSE;
502: lu->CleanUpSpooles = PETSC_FALSE;
503: lu->MatDuplicate = A->ops->duplicate;
504: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
505: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
506: lu->MatView = A->ops->view;
507: lu->MatAssemblyEnd = A->ops->assemblyend;
508: lu->MatDestroy = A->ops->destroy;
509: B->ops->duplicate = MatDuplicate_SeqAIJSpooles;
510: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJSpooles;
511: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJSpooles;
512: B->ops->view = MatView_SeqAIJSpooles;
513: B->ops->assemblyend = MatAssemblyEnd_SeqAIJSpooles;
514: B->ops->destroy = MatDestroy_SeqAIJSpooles;
516: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaijspooles_seqaij_C",
517: "MatConvert_Spooles_Base",MatConvert_Spooles_Base);
518: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqaijspooles_C",
519: "MatConvert_SeqAIJ_SeqAIJSpooles",MatConvert_SeqAIJ_SeqAIJSpooles);
520: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJSPOOLES);
521: *newmat = B;
522: return(0);
523: }
524: EXTERN_C_END
528: int MatDuplicate_SeqAIJSpooles(Mat A, MatDuplicateOption op, Mat *M) {
530: Mat_Spooles *lu=(Mat_Spooles *)A->spptr;
533: (*lu->MatDuplicate)(A,op,M);
534: PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Spooles));
535: return(0);
536: }
538: /*MC
539: MATSEQAIJSPOOLES - MATSEQAIJSPOOLES = "seqaijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential matrices
540: via the external package SEQAIJSPOOLES.
542: If SEQAIJSPOOLES is installed (see the manual for
543: instructions on how to declare the existence of external packages),
544: a matrix type can be constructed which invokes SPOOLES solvers.
545: After calling MatCreate(...,A), simply call MatSetType(A,MATSEQAIJSPOOLES).
546: This matrix type is only supported for double precision real.
548: This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is
549: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
550: the MATSEQAIJ type without data copy.
552: Options Database Keys:
553: + -mat_type seqaijspooles - sets the matrix type to "seqaijspooles" during a call to MatSetFromOptions()
554: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
555: . -mat_spooles_seed <seed> - random number seed used for ordering
556: . -mat_spooles_msglvl <msglvl> - message output level
557: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
558: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
559: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
560: . -mat_spooles_maxsize <n> - maximum size of a supernode
561: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
562: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
563: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
564: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
565: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
566: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
567: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
569: Level: beginner
571: .seealso: PCLU
572: M*/
574: EXTERN_C_BEGIN
577: int MatCreate_SeqAIJSpooles(Mat A) {
581: /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SeqAIJSpooles types */
582: PetscObjectChangeTypeName((PetscObject)A,MATSEQAIJSPOOLES);
583: MatSetType(A,MATSEQAIJ);
584: MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,&A);
585: return(0);
586: }
587: EXTERN_C_END