Actual source code: spooles.c
1: #define PETCSMAT_DLL
3: /*
4: Provides an interface to the Spooles serial sparse solver
5: */
6: #include src/mat/impls/aij/seq/aij.h
7: #include src/mat/impls/sbaij/seq/sbaij.h
8: #include src/mat/impls/aij/seq/spooles/spooles.h
10: /* make sun CC happy */
11: static void (*f)(void);
16: PetscErrorCode MatConvert_Spooles_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat)
17: {
19: Mat B=*newmat;
20: Mat_Spooles *lu=(Mat_Spooles*)A->spptr;
23: if (reuse == MAT_INITIAL_MATRIX) {
24: MatDuplicate(A,MAT_COPY_VALUES,&B);
25: }
26: /* Reset the stashed function pointers set by inherited routines */
27: B->ops->duplicate = lu->MatDuplicate;
28: B->ops->choleskyfactorsymbolic = lu->MatCholeskyFactorSymbolic;
29: B->ops->lufactorsymbolic = lu->MatLUFactorSymbolic;
30: B->ops->view = lu->MatView;
31: B->ops->assemblyend = lu->MatAssemblyEnd;
32: B->ops->destroy = lu->MatDestroy;
34: PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);
35: if (f) {
36: PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)lu->MatPreallocate);
37: }
38: PetscFree(lu);
39: A->spptr = PETSC_NULL;
41: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijspooles_seqaij_C","",PETSC_NULL);
42: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijspooles_C","",PETSC_NULL);
43: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaijspooles_mpiaij_C","",PETSC_NULL);
44: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_mpiaijspooles_C","",PETSC_NULL);
45: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaijspooles_seqsbaij_C","",PETSC_NULL);
46: PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbaijspooles_C","",PETSC_NULL);
47: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaijspooles_mpisbaij_C","",PETSC_NULL);
48: PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbaijspooles_C","",PETSC_NULL);
50: PetscObjectChangeTypeName((PetscObject)B,type);
51: *newmat = B;
52: return(0);
53: }
58: PetscErrorCode MatDestroy_SeqAIJSpooles(Mat A)
59: {
60: Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
62:
64: if (lu->CleanUpSpooles) {
65: FrontMtx_free(lu->frontmtx);
66: IV_free(lu->newToOldIV);
67: IV_free(lu->oldToNewIV);
68: InpMtx_free(lu->mtxA);
69: ETree_free(lu->frontETree);
70: IVL_free(lu->symbfacIVL);
71: SubMtxManager_free(lu->mtxmanager);
72: Graph_free(lu->graph);
73: }
74: MatConvert_Spooles_Base(A,lu->basetype,MAT_REUSE_MATRIX,&A);
75: (*A->ops->destroy)(A);
76: return(0);
77: }
81: PetscErrorCode MatSolve_SeqSpooles(Mat A,Vec b,Vec x)
82: {
83: Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
84: PetscScalar *array;
85: DenseMtx *mtxY, *mtxX ;
86: PetscErrorCode ierr;
87: PetscInt irow,neqns=A->cmap.n,nrow=A->rmap.n,*iv;
88: #if defined(PETSC_USE_COMPLEX)
89: double x_real,x_imag;
90: #else
91: double *entX;
92: #endif
95: mtxY = DenseMtx_new();
96: DenseMtx_init(mtxY, lu->options.typeflag, 0, 0, nrow, 1, 1, nrow); /* column major */
97: VecGetArray(b,&array);
99: if (lu->options.useQR) { /* copy b to mtxY */
100: for ( irow = 0 ; irow < nrow; irow++ )
101: #if !defined(PETSC_USE_COMPLEX)
102: DenseMtx_setRealEntry(mtxY, irow, 0, *array++);
103: #else
104: DenseMtx_setComplexEntry(mtxY, irow, 0, PetscRealPart(array[irow]), PetscImaginaryPart(array[irow]));
105: #endif
106: } else { /* copy permuted b to mtxY */
107: iv = IV_entries(lu->oldToNewIV);
108: for ( irow = 0 ; irow < nrow; irow++ )
109: #if !defined(PETSC_USE_COMPLEX)
110: DenseMtx_setRealEntry(mtxY, *iv++, 0, *array++);
111: #else
112: DenseMtx_setComplexEntry(mtxY,*iv++,0,PetscRealPart(array[irow]),PetscImaginaryPart(array[irow]));
113: #endif
114: }
115: VecRestoreArray(b,&array);
117: mtxX = DenseMtx_new();
118: DenseMtx_init(mtxX, lu->options.typeflag, 0, 0, neqns, 1, 1, neqns);
119: if (lu->options.useQR) {
120: FrontMtx_QR_solve(lu->frontmtx, lu->mtxA, mtxX, mtxY, lu->mtxmanager,
121: lu->cpus, lu->options.msglvl, lu->options.msgFile);
122: } else {
123: FrontMtx_solve(lu->frontmtx, mtxX, mtxY, lu->mtxmanager,
124: lu->cpus, lu->options.msglvl, lu->options.msgFile);
125: }
126: if ( lu->options.msglvl > 2 ) {
127: int err;
128: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n right hand side matrix after permutation");
129: DenseMtx_writeForHumanEye(mtxY, lu->options.msgFile);
130: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n solution matrix in new ordering");
131: DenseMtx_writeForHumanEye(mtxX, lu->options.msgFile);
132: err = fflush(lu->options.msgFile);
133: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
134: }
136: /* permute solution into original ordering, then copy to x */
137: DenseMtx_permuteRows(mtxX, lu->newToOldIV);
138: VecGetArray(x,&array);
140: #if !defined(PETSC_USE_COMPLEX)
141: entX = DenseMtx_entries(mtxX);
142: DVcopy(neqns, array, entX);
143: #else
144: for (irow=0; irow<nrow; irow++){
145: DenseMtx_complexEntry(mtxX,irow,0,&x_real,&x_imag);
146: array[irow] = x_real+x_imag*PETSC_i;
147: }
148: #endif
150: VecRestoreArray(x,&array);
151:
152: /* free memory */
153: DenseMtx_free(mtxX);
154: DenseMtx_free(mtxY);
155: return(0);
156: }
160: PetscErrorCode MatFactorNumeric_SeqSpooles(Mat A,MatFactorInfo *info,Mat *F)
161: {
162: Mat_Spooles *lu = (Mat_Spooles*)(*F)->spptr;
163: ChvManager *chvmanager ;
164: Chv *rootchv ;
165: IVL *adjIVL;
166: PetscErrorCode ierr;
167: PetscInt nz,nrow=A->rmap.n,irow,nedges,neqns=A->cmap.n,*ai,*aj,i,*diag=0,fierr;
168: PetscScalar *av;
169: double cputotal,facops;
170: #if defined(PETSC_USE_COMPLEX)
171: PetscInt nz_row,*aj_tmp;
172: PetscScalar *av_tmp;
173: #else
174: PetscInt *ivec1,*ivec2,j;
175: double *dvec;
176: #endif
177: PetscTruth isAIJ,isSeqAIJ;
178:
180: if (lu->flg == DIFFERENT_NONZERO_PATTERN) { /* first numeric factorization */
181: (*F)->ops->solve = MatSolve_SeqSpooles;
182: (*F)->ops->destroy = MatDestroy_SeqAIJSpooles;
183: (*F)->assembled = PETSC_TRUE;
184:
185: /* set Spooles options */
186: SetSpoolesOptions(A, &lu->options);
188: lu->mtxA = InpMtx_new();
189: }
191: /* copy A to Spooles' InpMtx object */
192: PetscTypeCompare((PetscObject)A,MATSEQAIJSPOOLES,&isSeqAIJ);
193: PetscTypeCompare((PetscObject)A,MATAIJSPOOLES,&isAIJ);
194: if (isSeqAIJ || isAIJ){
195: Mat_SeqAIJ *mat = (Mat_SeqAIJ*)A->data;
196: ai=mat->i; aj=mat->j; av=mat->a;
197: if (lu->options.symflag == SPOOLES_NONSYMMETRIC) {
198: nz=mat->nz;
199: } else { /* SPOOLES_SYMMETRIC || SPOOLES_HERMITIAN */
200: nz=(mat->nz + A->rmap.n)/2;
201: diag=mat->diag;
202: }
203: } else { /* A is SBAIJ */
204: Mat_SeqSBAIJ *mat = (Mat_SeqSBAIJ*)A->data;
205: ai=mat->i; aj=mat->j; av=mat->a;
206: nz=mat->nz;
207: }
208: InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, lu->options.typeflag, nz, 0);
209:
210: #if defined(PETSC_USE_COMPLEX)
211: for (irow=0; irow<nrow; irow++) {
212: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
213: nz_row = ai[irow+1] - ai[irow];
214: aj_tmp = aj + ai[irow];
215: av_tmp = av + ai[irow];
216: } else {
217: nz_row = ai[irow+1] - diag[irow];
218: aj_tmp = aj + diag[irow];
219: av_tmp = av + diag[irow];
220: }
221: for (i=0; i<nz_row; i++){
222: InpMtx_inputComplexEntry(lu->mtxA, irow, *aj_tmp++,PetscRealPart(*av_tmp),PetscImaginaryPart(*av_tmp));
223: av_tmp++;
224: }
225: }
226: #else
227: ivec1 = InpMtx_ivec1(lu->mtxA);
228: ivec2 = InpMtx_ivec2(lu->mtxA);
229: dvec = InpMtx_dvec(lu->mtxA);
230: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
231: for (irow = 0; irow < nrow; irow++){
232: for (i = ai[irow]; i<ai[irow+1]; i++) ivec1[i] = irow;
233: }
234: IVcopy(nz, ivec2, aj);
235: DVcopy(nz, dvec, av);
236: } else {
237: nz = 0;
238: for (irow = 0; irow < nrow; irow++){
239: for (j = diag[irow]; j<ai[irow+1]; j++) {
240: ivec1[nz] = irow;
241: ivec2[nz] = aj[j];
242: dvec[nz] = av[j];
243: nz++;
244: }
245: }
246: }
247: InpMtx_inputRealTriples(lu->mtxA, nz, ivec1, ivec2, dvec);
248: #endif
250: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
251: if ( lu->options.msglvl > 0 ) {
252: int err;
253: printf("\n\n input matrix");
254: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix");
255: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
256: err = fflush(lu->options.msgFile);
257: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
258: }
260: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
261: /*---------------------------------------------------
262: find a low-fill ordering
263: (1) create the Graph object
264: (2) order the graph
265: -------------------------------------------------------*/
266: if (lu->options.useQR){
267: adjIVL = InpMtx_adjForATA(lu->mtxA);
268: } else {
269: adjIVL = InpMtx_fullAdjacency(lu->mtxA);
270: }
271: nedges = IVL_tsize(adjIVL);
273: lu->graph = Graph_new();
274: Graph_init2(lu->graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL);
275: if ( lu->options.msglvl > 2 ) {
276: int err;
278: if (lu->options.useQR){
279: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of A^T A");
280: } else {
281: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of the input matrix");
282: }
283: Graph_writeForHumanEye(lu->graph, lu->options.msgFile);
284: err = fflush(lu->options.msgFile);
285: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
286: }
288: switch (lu->options.ordering) {
289: case 0:
290: lu->frontETree = orderViaBestOfNDandMS(lu->graph,
291: lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize,
292: lu->options.seed, lu->options.msglvl, lu->options.msgFile); break;
293: case 1:
294: lu->frontETree = orderViaMMD(lu->graph,lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
295: case 2:
296: lu->frontETree = orderViaMS(lu->graph, lu->options.maxdomainsize,
297: lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
298: case 3:
299: lu->frontETree = orderViaND(lu->graph, lu->options.maxdomainsize,
300: lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
301: default:
302: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown Spooles's ordering");
303: }
305: if ( lu->options.msglvl > 0 ) {
306: int err;
308: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree from ordering");
309: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile);
310: err = fflush(lu->options.msgFile);
311: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
312: }
313:
314: /* get the permutation, permute the front tree */
315: lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree);
316: lu->oldToNew = IV_entries(lu->oldToNewIV);
317: lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree);
318: if (!lu->options.useQR) ETree_permuteVertices(lu->frontETree, lu->oldToNewIV);
320: /* permute the matrix */
321: if (lu->options.useQR){
322: InpMtx_permute(lu->mtxA, NULL, lu->oldToNew);
323: } else {
324: InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew);
325: if ( lu->options.symflag == SPOOLES_SYMMETRIC) {
326: InpMtx_mapToUpperTriangle(lu->mtxA);
327: }
328: #if defined(PETSC_USE_COMPLEX)
329: if ( lu->options.symflag == SPOOLES_HERMITIAN ) {
330: InpMtx_mapToUpperTriangleH(lu->mtxA);
331: }
332: #endif
333: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS);
334: }
335: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
337: /* get symbolic factorization */
338: if (lu->options.useQR){
339: lu->symbfacIVL = SymbFac_initFromGraph(lu->frontETree, lu->graph);
340: IVL_overwrite(lu->symbfacIVL, lu->oldToNewIV);
341: IVL_sortUp(lu->symbfacIVL);
342: ETree_permuteVertices(lu->frontETree, lu->oldToNewIV);
343: } else {
344: lu->symbfacIVL = SymbFac_initFromInpMtx(lu->frontETree, lu->mtxA);
345: }
346: if ( lu->options.msglvl > 2 ) {
347: int err;
349: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n old-to-new permutation vector");
350: IV_writeForHumanEye(lu->oldToNewIV, lu->options.msgFile);
351: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n new-to-old permutation vector");
352: IV_writeForHumanEye(lu->newToOldIV, lu->options.msgFile);
353: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree after permutation");
354: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile);
355: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");
356: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
357: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n symbolic factorization");
358: IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile);
359: err = fflush(lu->options.msgFile);
360: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
361: }
363: lu->frontmtx = FrontMtx_new();
364: lu->mtxmanager = SubMtxManager_new();
365: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0);
367: } else { /* new num factorization using previously computed symbolic factor */
369: if (lu->options.pivotingflag) { /* different FrontMtx is required */
370: FrontMtx_free(lu->frontmtx);
371: lu->frontmtx = FrontMtx_new();
372: } else {
373: FrontMtx_clearData (lu->frontmtx);
374: }
376: SubMtxManager_free(lu->mtxmanager);
377: lu->mtxmanager = SubMtxManager_new();
378: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0);
380: /* permute mtxA */
381: if (lu->options.useQR){
382: InpMtx_permute(lu->mtxA, NULL, lu->oldToNew);
383: } else {
384: InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew);
385: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
386: InpMtx_mapToUpperTriangle(lu->mtxA);
387: }
388: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS);
389: }
390: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
391: if ( lu->options.msglvl > 2 ) {
392: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");
393: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
394: }
395: } /* end of if( lu->flg == DIFFERENT_NONZERO_PATTERN) */
396:
397: if (lu->options.useQR){
398: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag,
399: SPOOLES_SYMMETRIC, FRONTMTX_DENSE_FRONTS,
400: SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL,
401: lu->mtxmanager, lu->options.msglvl, lu->options.msgFile);
402: } else {
403: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag, lu->options.symflag,
404: FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, 0, NULL,
405: lu->mtxmanager, lu->options.msglvl, lu->options.msgFile);
406: }
408: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
409: if ( lu->options.patchAndGoFlag == 1 ) {
410: lu->frontmtx->patchinfo = PatchAndGoInfo_new();
411: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
412: lu->options.storeids, lu->options.storevalues);
413: } else if ( lu->options.patchAndGoFlag == 2 ) {
414: lu->frontmtx->patchinfo = PatchAndGoInfo_new();
415: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
416: lu->options.storeids, lu->options.storevalues);
417: }
418: }
420: /* numerical factorization */
421: chvmanager = ChvManager_new();
422: ChvManager_init(chvmanager, NO_LOCK, 1);
423: DVfill(10, lu->cpus, 0.0);
424: if (lu->options.useQR){
425: facops = 0.0 ;
426: FrontMtx_QR_factor(lu->frontmtx, lu->mtxA, chvmanager,
427: lu->cpus, &facops, lu->options.msglvl, lu->options.msgFile);
428: if ( lu->options.msglvl > 1 ) {
429: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");
430: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n facops = %9.2f", facops);
431: }
432: } else {
433: IVfill(20, lu->stats, 0);
434: rootchv = FrontMtx_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, 0.0,
435: chvmanager, &fierr, lu->cpus,lu->stats,lu->options.msglvl,lu->options.msgFile);
436: if (rootchv) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"\n matrix found to be singular");
437: if (fierr >= 0) SETERRQ1(PETSC_ERR_LIB,"\n error encountered at front %D", fierr);
438:
439: if(lu->options.FrontMtxInfo){
440: PetscPrintf(PETSC_COMM_SELF,"\n %8d pivots, %8d pivot tests, %8d delayed rows and columns\n",lu->stats[0], lu->stats[1], lu->stats[2]);
441: cputotal = lu->cpus[8] ;
442: if ( cputotal > 0.0 ) {
443: PetscPrintf(PETSC_COMM_SELF,
444: "\n cpus cpus/totaltime"
445: "\n initialize fronts %8.3f %6.2f"
446: "\n load original entries %8.3f %6.2f"
447: "\n update fronts %8.3f %6.2f"
448: "\n assemble postponed data %8.3f %6.2f"
449: "\n factor fronts %8.3f %6.2f"
450: "\n extract postponed data %8.3f %6.2f"
451: "\n store factor entries %8.3f %6.2f"
452: "\n miscellaneous %8.3f %6.2f"
453: "\n total time %8.3f \n",
454: lu->cpus[0], 100.*lu->cpus[0]/cputotal,
455: lu->cpus[1], 100.*lu->cpus[1]/cputotal,
456: lu->cpus[2], 100.*lu->cpus[2]/cputotal,
457: lu->cpus[3], 100.*lu->cpus[3]/cputotal,
458: lu->cpus[4], 100.*lu->cpus[4]/cputotal,
459: lu->cpus[5], 100.*lu->cpus[5]/cputotal,
460: lu->cpus[6], 100.*lu->cpus[6]/cputotal,
461: lu->cpus[7], 100.*lu->cpus[7]/cputotal, cputotal);
462: }
463: }
464: }
465: ChvManager_free(chvmanager);
467: if ( lu->options.msglvl > 0 ) {
468: int err;
470: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");
471: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile);
472: err = fflush(lu->options.msgFile);
473: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
474: }
476: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
477: if ( lu->options.patchAndGoFlag == 1 ) {
478: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
479: if (lu->options.msglvl > 0 ){
480: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");
481: IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile);
482: }
483: }
484: PatchAndGoInfo_free(lu->frontmtx->patchinfo);
485: } else if ( lu->options.patchAndGoFlag == 2 ) {
486: if (lu->options.msglvl > 0 ){
487: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
488: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");
489: IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile);
490: }
491: if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
492: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n perturbations");
493: DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile);
494: }
495: }
496: PatchAndGoInfo_free(lu->frontmtx->patchinfo);
497: }
498: }
500: /* post-process the factorization */
501: FrontMtx_postProcess(lu->frontmtx, lu->options.msglvl, lu->options.msgFile);
502: if ( lu->options.msglvl > 2 ) {
503: int err;
505: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix after post-processing");
506: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile);
507: err = fflush(lu->options.msgFile);
508: if (err) SETERRQ(PETSC_ERR_SYS,"fflush() failed on file");
509: }
511: lu->flg = SAME_NONZERO_PATTERN;
512: lu->CleanUpSpooles = PETSC_TRUE;
513: return(0);
514: }
519: PetscErrorCode MatConvert_SeqAIJ_SeqAIJSpooles(Mat A,MatType type,MatReuse reuse,Mat *newmat)
520: {
522: Mat B=*newmat;
523: Mat_Spooles *lu;
526: PetscNewLog(B,Mat_Spooles,&lu);
527: if (reuse == MAT_INITIAL_MATRIX) {
528: /* This routine is inherited, so we know the type is correct. */
529: MatDuplicate(A,MAT_COPY_VALUES,&B);
530: lu->MatDuplicate = B->ops->duplicate;
531: lu->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic;
532: lu->MatLUFactorSymbolic = B->ops->lufactorsymbolic;
533: lu->MatView = B->ops->view;
534: lu->MatAssemblyEnd = B->ops->assemblyend;
535: lu->MatDestroy = B->ops->destroy;
536: } else {
537: lu->MatDuplicate = A->ops->duplicate;
538: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
539: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
540: lu->MatView = A->ops->view;
541: lu->MatAssemblyEnd = A->ops->assemblyend;
542: lu->MatDestroy = A->ops->destroy;
543: }
544: B->spptr = (void*)lu;
545: lu->basetype = MATSEQAIJ;
546: lu->useQR = PETSC_FALSE;
547: lu->CleanUpSpooles = PETSC_FALSE;
549: B->ops->duplicate = MatDuplicate_Spooles;
550: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJSpooles;
551: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJSpooles;
552: B->ops->view = MatView_Spooles;
553: B->ops->assemblyend = MatAssemblyEnd_SeqAIJSpooles;
554: B->ops->destroy = MatDestroy_SeqAIJSpooles;
556: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaijspooles_seqaij_C",
557: "MatConvert_Spooles_Base",MatConvert_Spooles_Base);
558: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqaijspooles_C",
559: "MatConvert_SeqAIJ_SeqAIJSpooles",MatConvert_SeqAIJ_SeqAIJSpooles);
560: /* PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJSPOOLES); */
561: PetscObjectChangeTypeName((PetscObject)B,type);
562: *newmat = B;
563: return(0);
564: }
569: PetscErrorCode MatDuplicate_Spooles(Mat A, MatDuplicateOption op, Mat *M) {
571: Mat_Spooles *lu=(Mat_Spooles *)A->spptr;
574: (*lu->MatDuplicate)(A,op,M);
575: PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Spooles));
576: return(0);
577: }
579: /*MC
580: MATSEQAIJSPOOLES - MATSEQAIJSPOOLES = "seqaijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential matrices
581: via the external package SPOOLES.
583: If SPOOLES is installed (see the manual for
584: instructions on how to declare the existence of external packages),
585: a matrix type can be constructed which invokes SPOOLES solvers.
586: After calling MatCreate(...,A), simply call MatSetType(A,MATSEQAIJSPOOLES), then
587: optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT
588: call MatCreateSeqAIJ() directly or the preallocation information will be LOST!
590: This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation() is
591: supported for this matrix type. One can also call MatConvert() for an inplace conversion to or from
592: the MATSEQAIJ type without data copy AFTER the matrix values have been set.
594: Options Database Keys:
595: + -mat_type seqaijspooles - sets the matrix type to "seqaijspooles" during a call to MatSetFromOptions()
596: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
597: . -mat_spooles_seed <seed> - random number seed used for ordering
598: . -mat_spooles_msglvl <msglvl> - message output level
599: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
600: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
601: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
602: . -mat_spooles_maxsize <n> - maximum size of a supernode
603: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
604: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
605: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
606: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
607: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
608: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
609: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
611: Level: beginner
613: .seealso: PCLU
614: M*/
619: PetscErrorCode MatCreate_SeqAIJSpooles(Mat A)
620: {
624: MatSetType(A,MATSEQAIJ);
625: MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,MAT_REUSE_MATRIX,&A);
626: return(0);
627: }
630: /*MC
631: MATAIJSPOOLES - MATAIJSPOOLES = "aijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential and parellel matrices
632: via the external package SPOOLES.
634: If SPOOLES is installed (see the manual for
635: instructions on how to declare the existence of external packages),
636: a matrix type can be constructed which invokes SPOOLES solvers.
637: After calling MatCreate(...,A), simply call MatSetType(A,MATAIJSPOOLES), then
638: optionally call MatSeqAIJSetPreallocation() or MatMPIAIJSetPreallocation() etc DO NOT
639: call MatCreateSeqAIJ/MPIAIJ() directly or the preallocation information will be LOST!
640: This matrix type is supported for double precision real and complex.
642: This matrix inherits from MATAIJ. As a result, MatSeqAIJSetPreallocation() and MatMPIAIJSetPreallocation() are
643: supported for this matrix type. One can also call MatConvert() for an inplace conversion to or from
644: the MATAIJ type without data copy AFTER the matrix values have been set.
646: Options Database Keys:
647: + -mat_type aijspooles - sets the matrix type to "aijspooles" during a call to MatSetFromOptions()
648: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
649: . -mat_spooles_seed <seed> - random number seed used for ordering
650: . -mat_spooles_msglvl <msglvl> - message output level
651: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
652: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
653: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
654: . -mat_spooles_maxsize <n> - maximum size of a supernode
655: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
656: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
657: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
658: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
659: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
660: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
661: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
663: Level: beginner
665: .seealso: PCLU
666: M*/
670: PetscErrorCode MatCreate_AIJSpooles(Mat A)
671: {
673: PetscMPIInt size;
676: MPI_Comm_size(((PetscObject)A)->comm,&size);
677: if (size == 1) {
678: MatSetType(A,MATSEQAIJ);
679: MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,MAT_REUSE_MATRIX,&A);
680: } else {
681: MatSetType(A,MATMPIAIJ);
682: MatConvert_MPIAIJ_MPIAIJSpooles(A,MATMPIAIJSPOOLES,MAT_REUSE_MATRIX,&A);
683: }
684: return(0);
685: }