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: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n right hand side matrix after permutation");
128: DenseMtx_writeForHumanEye(mtxY, lu->options.msgFile);
129: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n solution matrix in new ordering");
130: DenseMtx_writeForHumanEye(mtxX, lu->options.msgFile);
131: fflush(lu->options.msgFile);
132: }
134: /* permute solution into original ordering, then copy to x */
135: DenseMtx_permuteRows(mtxX, lu->newToOldIV);
136: VecGetArray(x,&array);
138: #if !defined(PETSC_USE_COMPLEX)
139: entX = DenseMtx_entries(mtxX);
140: DVcopy(neqns, array, entX);
141: #else
142: for (irow=0; irow<nrow; irow++){
143: DenseMtx_complexEntry(mtxX,irow,0,&x_real,&x_imag);
144: array[irow] = x_real+x_imag*PETSC_i;
145: }
146: #endif
148: VecRestoreArray(x,&array);
149:
150: /* free memory */
151: DenseMtx_free(mtxX);
152: DenseMtx_free(mtxY);
153: return(0);
154: }
158: PetscErrorCode MatFactorNumeric_SeqSpooles(Mat A,MatFactorInfo *info,Mat *F)
159: {
160: Mat_Spooles *lu = (Mat_Spooles*)(*F)->spptr;
161: ChvManager *chvmanager ;
162: Chv *rootchv ;
163: IVL *adjIVL;
164: PetscErrorCode ierr;
165: PetscInt nz,nrow=A->rmap.n,irow,nedges,neqns=A->cmap.n,*ai,*aj,i,*diag=0,fierr;
166: PetscScalar *av;
167: double cputotal,facops;
168: #if defined(PETSC_USE_COMPLEX)
169: PetscInt nz_row,*aj_tmp;
170: PetscScalar *av_tmp;
171: #else
172: PetscInt *ivec1,*ivec2,j;
173: double *dvec;
174: #endif
175: PetscTruth isAIJ,isSeqAIJ;
176:
178: if (lu->flg == DIFFERENT_NONZERO_PATTERN) { /* first numeric factorization */
179: (*F)->ops->solve = MatSolve_SeqSpooles;
180: (*F)->ops->destroy = MatDestroy_SeqAIJSpooles;
181: (*F)->assembled = PETSC_TRUE;
182:
183: /* set Spooles options */
184: SetSpoolesOptions(A, &lu->options);
186: lu->mtxA = InpMtx_new();
187: }
189: /* copy A to Spooles' InpMtx object */
190: PetscTypeCompare((PetscObject)A,MATSEQAIJSPOOLES,&isSeqAIJ);
191: PetscTypeCompare((PetscObject)A,MATAIJSPOOLES,&isAIJ);
192: if (isSeqAIJ || isAIJ){
193: Mat_SeqAIJ *mat = (Mat_SeqAIJ*)A->data;
194: ai=mat->i; aj=mat->j; av=mat->a;
195: if (lu->options.symflag == SPOOLES_NONSYMMETRIC) {
196: nz=mat->nz;
197: } else { /* SPOOLES_SYMMETRIC || SPOOLES_HERMITIAN */
198: nz=(mat->nz + A->rmap.n)/2;
199: diag=mat->diag;
200: }
201: } else { /* A is SBAIJ */
202: Mat_SeqSBAIJ *mat = (Mat_SeqSBAIJ*)A->data;
203: ai=mat->i; aj=mat->j; av=mat->a;
204: nz=mat->nz;
205: }
206: InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, lu->options.typeflag, nz, 0);
207:
208: #if defined(PETSC_USE_COMPLEX)
209: for (irow=0; irow<nrow; irow++) {
210: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
211: nz_row = ai[irow+1] - ai[irow];
212: aj_tmp = aj + ai[irow];
213: av_tmp = av + ai[irow];
214: } else {
215: nz_row = ai[irow+1] - diag[irow];
216: aj_tmp = aj + diag[irow];
217: av_tmp = av + diag[irow];
218: }
219: for (i=0; i<nz_row; i++){
220: InpMtx_inputComplexEntry(lu->mtxA, irow, *aj_tmp++,PetscRealPart(*av_tmp),PetscImaginaryPart(*av_tmp));
221: av_tmp++;
222: }
223: }
224: #else
225: ivec1 = InpMtx_ivec1(lu->mtxA);
226: ivec2 = InpMtx_ivec2(lu->mtxA);
227: dvec = InpMtx_dvec(lu->mtxA);
228: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
229: for (irow = 0; irow < nrow; irow++){
230: for (i = ai[irow]; i<ai[irow+1]; i++) ivec1[i] = irow;
231: }
232: IVcopy(nz, ivec2, aj);
233: DVcopy(nz, dvec, av);
234: } else {
235: nz = 0;
236: for (irow = 0; irow < nrow; irow++){
237: for (j = diag[irow]; j<ai[irow+1]; j++) {
238: ivec1[nz] = irow;
239: ivec2[nz] = aj[j];
240: dvec[nz] = av[j];
241: nz++;
242: }
243: }
244: }
245: InpMtx_inputRealTriples(lu->mtxA, nz, ivec1, ivec2, dvec);
246: #endif
248: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
249: if ( lu->options.msglvl > 0 ) {
250: printf("\n\n input matrix");
251: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix");
252: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
253: fflush(lu->options.msgFile);
254: }
256: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
257: /*---------------------------------------------------
258: find a low-fill ordering
259: (1) create the Graph object
260: (2) order the graph
261: -------------------------------------------------------*/
262: if (lu->options.useQR){
263: adjIVL = InpMtx_adjForATA(lu->mtxA);
264: } else {
265: adjIVL = InpMtx_fullAdjacency(lu->mtxA);
266: }
267: nedges = IVL_tsize(adjIVL);
269: lu->graph = Graph_new();
270: Graph_init2(lu->graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL);
271: if ( lu->options.msglvl > 2 ) {
272: if (lu->options.useQR){
273: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of A^T A");
274: } else {
275: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of the input matrix");
276: }
277: Graph_writeForHumanEye(lu->graph, lu->options.msgFile);
278: fflush(lu->options.msgFile);
279: }
281: switch (lu->options.ordering) {
282: case 0:
283: lu->frontETree = orderViaBestOfNDandMS(lu->graph,
284: lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize,
285: lu->options.seed, lu->options.msglvl, lu->options.msgFile); break;
286: case 1:
287: lu->frontETree = orderViaMMD(lu->graph,lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
288: case 2:
289: lu->frontETree = orderViaMS(lu->graph, lu->options.maxdomainsize,
290: lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
291: case 3:
292: lu->frontETree = orderViaND(lu->graph, lu->options.maxdomainsize,
293: lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
294: default:
295: SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown Spooles's ordering");
296: }
298: if ( lu->options.msglvl > 0 ) {
299: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree from ordering");
300: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile);
301: fflush(lu->options.msgFile);
302: }
303:
304: /* get the permutation, permute the front tree */
305: lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree);
306: lu->oldToNew = IV_entries(lu->oldToNewIV);
307: lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree);
308: if (!lu->options.useQR) ETree_permuteVertices(lu->frontETree, lu->oldToNewIV);
310: /* permute the matrix */
311: if (lu->options.useQR){
312: InpMtx_permute(lu->mtxA, NULL, lu->oldToNew);
313: } else {
314: InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew);
315: if ( lu->options.symflag == SPOOLES_SYMMETRIC) {
316: InpMtx_mapToUpperTriangle(lu->mtxA);
317: }
318: #if defined(PETSC_USE_COMPLEX)
319: if ( lu->options.symflag == SPOOLES_HERMITIAN ) {
320: InpMtx_mapToUpperTriangleH(lu->mtxA);
321: }
322: #endif
323: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS);
324: }
325: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
327: /* get symbolic factorization */
328: if (lu->options.useQR){
329: lu->symbfacIVL = SymbFac_initFromGraph(lu->frontETree, lu->graph);
330: IVL_overwrite(lu->symbfacIVL, lu->oldToNewIV);
331: IVL_sortUp(lu->symbfacIVL);
332: ETree_permuteVertices(lu->frontETree, lu->oldToNewIV);
333: } else {
334: lu->symbfacIVL = SymbFac_initFromInpMtx(lu->frontETree, lu->mtxA);
335: }
336: if ( lu->options.msglvl > 2 ) {
337: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n old-to-new permutation vector");
338: IV_writeForHumanEye(lu->oldToNewIV, lu->options.msgFile);
339: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n new-to-old permutation vector");
340: IV_writeForHumanEye(lu->newToOldIV, lu->options.msgFile);
341: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree after permutation");
342: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile);
343: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");
344: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
345: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n symbolic factorization");
346: IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile);
347: fflush(lu->options.msgFile);
348: }
350: lu->frontmtx = FrontMtx_new();
351: lu->mtxmanager = SubMtxManager_new();
352: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0);
354: } else { /* new num factorization using previously computed symbolic factor */
356: if (lu->options.pivotingflag) { /* different FrontMtx is required */
357: FrontMtx_free(lu->frontmtx);
358: lu->frontmtx = FrontMtx_new();
359: } else {
360: FrontMtx_clearData (lu->frontmtx);
361: }
363: SubMtxManager_free(lu->mtxmanager);
364: lu->mtxmanager = SubMtxManager_new();
365: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0);
367: /* permute mtxA */
368: if (lu->options.useQR){
369: InpMtx_permute(lu->mtxA, NULL, lu->oldToNew);
370: } else {
371: InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew);
372: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
373: InpMtx_mapToUpperTriangle(lu->mtxA);
374: }
375: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS);
376: }
377: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
378: if ( lu->options.msglvl > 2 ) {
379: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");
380: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
381: }
382: } /* end of if( lu->flg == DIFFERENT_NONZERO_PATTERN) */
383:
384: if (lu->options.useQR){
385: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag,
386: SPOOLES_SYMMETRIC, FRONTMTX_DENSE_FRONTS,
387: SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL,
388: lu->mtxmanager, lu->options.msglvl, lu->options.msgFile);
389: } else {
390: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag, lu->options.symflag,
391: FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, 0, NULL,
392: lu->mtxmanager, lu->options.msglvl, lu->options.msgFile);
393: }
395: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
396: if ( lu->options.patchAndGoFlag == 1 ) {
397: lu->frontmtx->patchinfo = PatchAndGoInfo_new();
398: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
399: lu->options.storeids, lu->options.storevalues);
400: } else if ( lu->options.patchAndGoFlag == 2 ) {
401: lu->frontmtx->patchinfo = PatchAndGoInfo_new();
402: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
403: lu->options.storeids, lu->options.storevalues);
404: }
405: }
407: /* numerical factorization */
408: chvmanager = ChvManager_new();
409: ChvManager_init(chvmanager, NO_LOCK, 1);
410: DVfill(10, lu->cpus, 0.0);
411: if (lu->options.useQR){
412: facops = 0.0 ;
413: FrontMtx_QR_factor(lu->frontmtx, lu->mtxA, chvmanager,
414: lu->cpus, &facops, lu->options.msglvl, lu->options.msgFile);
415: if ( lu->options.msglvl > 1 ) {
416: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");
417: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n facops = %9.2f", facops);
418: }
419: } else {
420: IVfill(20, lu->stats, 0);
421: rootchv = FrontMtx_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, 0.0,
422: chvmanager, &fierr, lu->cpus,lu->stats,lu->options.msglvl,lu->options.msgFile);
423: if (rootchv) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"\n matrix found to be singular");
424: if (fierr >= 0) SETERRQ1(PETSC_ERR_LIB,"\n error encountered at front %D", fierr);
425:
426: if(lu->options.FrontMtxInfo){
427: 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]);
428: cputotal = lu->cpus[8] ;
429: if ( cputotal > 0.0 ) {
430: PetscPrintf(PETSC_COMM_SELF,
431: "\n cpus cpus/totaltime"
432: "\n initialize fronts %8.3f %6.2f"
433: "\n load original entries %8.3f %6.2f"
434: "\n update fronts %8.3f %6.2f"
435: "\n assemble postponed data %8.3f %6.2f"
436: "\n factor fronts %8.3f %6.2f"
437: "\n extract postponed data %8.3f %6.2f"
438: "\n store factor entries %8.3f %6.2f"
439: "\n miscellaneous %8.3f %6.2f"
440: "\n total time %8.3f \n",
441: lu->cpus[0], 100.*lu->cpus[0]/cputotal,
442: lu->cpus[1], 100.*lu->cpus[1]/cputotal,
443: lu->cpus[2], 100.*lu->cpus[2]/cputotal,
444: lu->cpus[3], 100.*lu->cpus[3]/cputotal,
445: lu->cpus[4], 100.*lu->cpus[4]/cputotal,
446: lu->cpus[5], 100.*lu->cpus[5]/cputotal,
447: lu->cpus[6], 100.*lu->cpus[6]/cputotal,
448: lu->cpus[7], 100.*lu->cpus[7]/cputotal, cputotal);
449: }
450: }
451: }
452: ChvManager_free(chvmanager);
454: if ( lu->options.msglvl > 0 ) {
455: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");
456: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile);
457: fflush(lu->options.msgFile);
458: }
460: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
461: if ( lu->options.patchAndGoFlag == 1 ) {
462: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
463: if (lu->options.msglvl > 0 ){
464: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");
465: IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile);
466: }
467: }
468: PatchAndGoInfo_free(lu->frontmtx->patchinfo);
469: } else if ( lu->options.patchAndGoFlag == 2 ) {
470: if (lu->options.msglvl > 0 ){
471: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
472: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");
473: IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile);
474: }
475: if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
476: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n perturbations");
477: DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile);
478: }
479: }
480: PatchAndGoInfo_free(lu->frontmtx->patchinfo);
481: }
482: }
484: /* post-process the factorization */
485: FrontMtx_postProcess(lu->frontmtx, lu->options.msglvl, lu->options.msgFile);
486: if ( lu->options.msglvl > 2 ) {
487: PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix after post-processing");
488: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile);
489: fflush(lu->options.msgFile);
490: }
492: lu->flg = SAME_NONZERO_PATTERN;
493: lu->CleanUpSpooles = PETSC_TRUE;
494: return(0);
495: }
500: PetscErrorCode MatConvert_SeqAIJ_SeqAIJSpooles(Mat A,MatType type,MatReuse reuse,Mat *newmat)
501: {
503: Mat B=*newmat;
504: Mat_Spooles *lu;
507: PetscNew(Mat_Spooles,&lu);
508: if (reuse == MAT_INITIAL_MATRIX) {
509: /* This routine is inherited, so we know the type is correct. */
510: MatDuplicate(A,MAT_COPY_VALUES,&B);
511: lu->MatDuplicate = B->ops->duplicate;
512: lu->MatCholeskyFactorSymbolic = B->ops->choleskyfactorsymbolic;
513: lu->MatLUFactorSymbolic = B->ops->lufactorsymbolic;
514: lu->MatView = B->ops->view;
515: lu->MatAssemblyEnd = B->ops->assemblyend;
516: lu->MatDestroy = B->ops->destroy;
517: } else {
518: lu->MatDuplicate = A->ops->duplicate;
519: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
520: lu->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
521: lu->MatView = A->ops->view;
522: lu->MatAssemblyEnd = A->ops->assemblyend;
523: lu->MatDestroy = A->ops->destroy;
524: }
525: B->spptr = (void*)lu;
526: lu->basetype = MATSEQAIJ;
527: lu->useQR = PETSC_FALSE;
528: lu->CleanUpSpooles = PETSC_FALSE;
530: B->ops->duplicate = MatDuplicate_Spooles;
531: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJSpooles;
532: B->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJSpooles;
533: B->ops->view = MatView_Spooles;
534: B->ops->assemblyend = MatAssemblyEnd_SeqAIJSpooles;
535: B->ops->destroy = MatDestroy_SeqAIJSpooles;
537: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaijspooles_seqaij_C",
538: "MatConvert_Spooles_Base",MatConvert_Spooles_Base);
539: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqaijspooles_C",
540: "MatConvert_SeqAIJ_SeqAIJSpooles",MatConvert_SeqAIJ_SeqAIJSpooles);
541: /* PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJSPOOLES); */
542: PetscObjectChangeTypeName((PetscObject)B,type);
543: *newmat = B;
544: return(0);
545: }
550: PetscErrorCode MatDuplicate_Spooles(Mat A, MatDuplicateOption op, Mat *M) {
552: Mat_Spooles *lu=(Mat_Spooles *)A->spptr;
555: (*lu->MatDuplicate)(A,op,M);
556: PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Spooles));
557: return(0);
558: }
560: /*MC
561: MATSEQAIJSPOOLES - MATSEQAIJSPOOLES = "seqaijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential matrices
562: via the external package SPOOLES.
564: If SPOOLES is installed (see the manual for
565: instructions on how to declare the existence of external packages),
566: a matrix type can be constructed which invokes SPOOLES solvers.
567: After calling MatCreate(...,A), simply call MatSetType(A,MATSEQAIJSPOOLES).
569: This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is
570: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
571: the MATSEQAIJ type without data copy.
573: Options Database Keys:
574: + -mat_type seqaijspooles - sets the matrix type to "seqaijspooles" during a call to MatSetFromOptions()
575: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
576: . -mat_spooles_seed <seed> - random number seed used for ordering
577: . -mat_spooles_msglvl <msglvl> - message output level
578: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
579: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
580: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
581: . -mat_spooles_maxsize <n> - maximum size of a supernode
582: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
583: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
584: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
585: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
586: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
587: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
588: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
590: Level: beginner
592: .seealso: PCLU
593: M*/
598: PetscErrorCode MatCreate_SeqAIJSpooles(Mat A)
599: {
603: MatSetType(A,MATSEQAIJ);
604: MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,MAT_REUSE_MATRIX,&A);
605: return(0);
606: }
609: /*MC
610: MATAIJSPOOLES - MATAIJSPOOLES = "aijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential and parellel matrices
611: via the external package SPOOLES.
613: If SPOOLES is installed (see the manual for
614: instructions on how to declare the existence of external packages),
615: a matrix type can be constructed which invokes SPOOLES solvers.
616: After calling MatCreate(...,A), simply call MatSetType(A,MATAIJSPOOLES).
617: This matrix type is supported for double precision real and complex.
619: This matrix inherits from MATAIJ. As a result, MatSeqAIJSetPreallocation and MatMPIAIJSetPreallocation are
620: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
621: the MATAIJ type without data copy.
623: Options Database Keys:
624: + -mat_type aijspooles - sets the matrix type to "aijspooles" during a call to MatSetFromOptions()
625: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
626: . -mat_spooles_seed <seed> - random number seed used for ordering
627: . -mat_spooles_msglvl <msglvl> - message output level
628: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
629: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
630: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
631: . -mat_spooles_maxsize <n> - maximum size of a supernode
632: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
633: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
634: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
635: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
636: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
637: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
638: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object
640: Level: beginner
642: .seealso: PCLU
643: M*/
647: PetscErrorCode MatCreate_AIJSpooles(Mat A)
648: {
650: PetscMPIInt size;
653: MPI_Comm_size(A->comm,&size);
654: if (size == 1) {
655: MatSetType(A,MATSEQAIJ);
656: MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,MAT_REUSE_MATRIX,&A);
657: } else {
658: MatSetType(A,MATMPIAIJ);
659: MatConvert_MPIAIJ_MPIAIJSpooles(A,MATMPIAIJSPOOLES,MAT_REUSE_MATRIX,&A);
660: }
661: return(0);
662: }