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