Actual source code: mpispooles.c
1: /*$Id: mpispooles.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2: /*
3: Provides an interface to the Spooles parallel sparse solver (MPI SPOOLES)
4: */
6: #include src/mat/impls/aij/seq/aij.h
7: #include src/mat/impls/sbaij/seq/sbaij.h
8: #include src/mat/impls/baij/seq/baij.h
9: #include src/mat/impls/aij/mpi/mpiaij.h
10: #include src/mat/impls/sbaij/mpi/mpisbaij.h
12: #if defined(PETSC_HAVE_SPOOLES) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
13: #include src/mat/impls/aij/seq/spooles.h
15: extern int SetSpoolesOptions(Mat, Spooles_options *);
16: extern int MatDestroy_MPIAIJ(Mat);
18: #undef __FUNCT__
20: int MatDestroy_MPIAIJ_Spooles(Mat A)
21: {
22: Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
23: int ierr;
24:
26:
27: FrontMtx_free(lu->frontmtx) ;
28: IV_free(lu->newToOldIV) ;
29: IV_free(lu->oldToNewIV) ;
30: IV_free(lu->vtxmapIV) ;
31: InpMtx_free(lu->mtxA) ;
32: ETree_free(lu->frontETree) ;
33: IVL_free(lu->symbfacIVL) ;
34: SubMtxManager_free(lu->mtxmanager) ;
35: DenseMtx_free(lu->mtxX) ;
36: DenseMtx_free(lu->mtxY) ;
37:
38: VecDestroy(lu->vec_spooles);
39: ISDestroy(lu->iden);
40: ISDestroy(lu->is_petsc);
41: VecScatterDestroy(lu->scat);
43: PetscFree(lu);
44: MatDestroy_MPIAIJ(A);
46: return(0);
47: }
49: #undef __FUNCT__
51: int MatSolve_MPIAIJ_Spooles(Mat A,Vec b,Vec x)
52: {
53: Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
54: int ierr,size,rank,m=A->m,irow,*rowindY;
55: PetscScalar *array;
56: DenseMtx *newY ;
57: SubMtxManager *solvemanager ;
60: MPI_Comm_size(PETSC_COMM_WORLD,&size);
61: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
62:
63: /* copy b into spooles' rhs mtxY */
64: DenseMtx_init(lu->mtxY, SPOOLES_REAL, 0, 0, m, 1, 1, m) ;
65: VecGetArray(b,&array);
67: /* doesn't work!
68: PetscMalloc(m*sizeof(int),&rowindY);
69: for ( irow = 0 ; irow < m ; irow++ ) rowindY[irow] = irow + lu->rstart;
70: int colind=0;
71: DenseMtx_initWithPointers(lu->mtxY,SPOOLES_REAL,0,0,m,1,1,m,rowindY,&colind,array);
72: */
74: DenseMtx_rowIndices(lu->mtxY, &m, &rowindY) ; /* get m, rowind */
75: for ( irow = 0 ; irow < m ; irow++ ) {
76: rowindY[irow] = irow + lu->rstart; /* global rowind */
77: DenseMtx_setRealEntry(lu->mtxY, irow, 0, *array++) ;
78: }
79: /* DenseMtx_column(lu->mtxY, 0, &array); doesn't work! */
80: VecRestoreArray(b,&array);
81:
82: if ( lu->options.msglvl > 2 ) {
83: fprintf(lu->options.msgFile, "nn 1 matrix in original ordering") ;
84: DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
85: fflush(lu->options.msgFile) ;
86: }
87:
88: /* permute and redistribute Y if necessary */
89: DenseMtx_permuteRows(lu->mtxY, lu->oldToNewIV) ;
90: if ( lu->options.msglvl > 2 ) {
91: fprintf(lu->options.msgFile, "nn rhs matrix in new ordering") ;
92: DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
93: fflush(lu->options.msgFile) ;
94: }
96: MPI_Barrier(MPI_COMM_WORLD) ; /* for initializing firsttag, because the num. of tags used
97: by FrontMtx_MPI_split() is unknown */
98: lu->firsttag = 0;
99: newY = DenseMtx_MPI_splitByRows(lu->mtxY, lu->vtxmapIV, lu->stats, lu->options.msglvl,
100: lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
101: DenseMtx_free(lu->mtxY) ;
102: lu->mtxY = newY ;
103: lu->firsttag += size ;
104: if ( lu->options.msglvl > 2 ) {
105: fprintf(lu->options.msgFile, "nn split DenseMtx Y") ;
106: DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
107: fflush(lu->options.msgFile) ;
108: }
110: if ( FRONTMTX_IS_PIVOTING(lu->frontmtx) ) {
111: /* pivoting has taken place, redistribute the right hand side
112: to match the final rows and columns in the fronts */
113: IV *rowmapIV ;
114: rowmapIV = FrontMtx_MPI_rowmapIV(lu->frontmtx, lu->ownersIV, lu->options.msglvl,
115: lu->options.msgFile, MPI_COMM_WORLD) ;
116: newY = DenseMtx_MPI_splitByRows(lu->mtxY, rowmapIV, lu->stats, lu->options.msglvl,
117: lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
118: DenseMtx_free(lu->mtxY) ;
119: lu->mtxY = newY ;
120: IV_free(rowmapIV) ;
121: lu->firsttag += size;
122: }
123: if ( lu->options.msglvl > 2 ) {
124: fprintf(lu->options.msgFile, "nn rhs matrix after split") ;
125: DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
126: fflush(lu->options.msgFile) ;
127: }
129: if ( lu->nmycol > 0 ) IVcopy(lu->nmycol,lu->rowindX,IV_entries(lu->ownedColumnsIV)); /* must do for each solve */
130:
131: /* solve the linear system */
132: solvemanager = SubMtxManager_new() ;
133: SubMtxManager_init(solvemanager, NO_LOCK, 0) ;
134: FrontMtx_MPI_solve(lu->frontmtx, lu->mtxX, lu->mtxY, solvemanager, lu->solvemap, lu->cpus,
135: lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
136: SubMtxManager_free(solvemanager) ;
137: if ( lu->options.msglvl > 2 ) {
138: fprintf(lu->options.msgFile, "n solution in new ordering") ;
139: DenseMtx_writeForHumanEye(lu->mtxX, lu->options.msgFile) ;
140: }
142: /* permute the solution into the original ordering */
143: DenseMtx_permuteRows(lu->mtxX, lu->newToOldIV) ;
144: if ( lu->options.msglvl > 2 ) {
145: fprintf(lu->options.msgFile, "nn solution in old ordering") ;
146: DenseMtx_writeForHumanEye(lu->mtxX, lu->options.msgFile) ;
147: fflush(lu->options.msgFile) ;
148: }
149:
150: /* scatter local solution mtxX into mpi vector x */
151: if( !lu->scat ){ /* create followings once for each numfactorization */
152: VecCreateSeqWithArray(PETSC_COMM_SELF,lu->nmycol,lu->entX,&lu->vec_spooles); /* vec_spooles <- mtxX */
153: ISCreateStride(PETSC_COMM_SELF,lu->nmycol,0,1,&lu->iden);
154: ISCreateGeneral(PETSC_COMM_SELF,lu->nmycol,lu->rowindX,&lu->is_petsc);
155: VecScatterCreate(lu->vec_spooles,lu->iden,x,lu->is_petsc,&lu->scat);
156: }
158: VecScatterBegin(lu->vec_spooles,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
159: VecScatterEnd(lu->vec_spooles,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
160:
161: return(0);
162: }
164: #undef __FUNCT__
166: int MatFactorNumeric_MPIAIJ_Spooles(Mat A,Mat *F)
167: {
168: Mat_Spooles *lu = (Mat_Spooles*)(*F)->spptr;
169: int rank,size,ierr,lookahead=0;
170: ChvManager *chvmanager ;
171: Chv *rootchv ;
172: Graph *graph ;
173: IVL *adjIVL;
174: DV *cumopsDV ;
175: double droptol=0.0,*opcounts,minops,cutoff,*val;
176: InpMtx *newA ;
177: PetscScalar *av, *bv;
178: int *ai, *aj, *bi,*bj, nz, *ajj, *bjj, *garray,
179: i,j,irow,jcol,countA,countB,jB,*row,*col,colA_start,jj;
180: int M=A->M,m=A->m,root,nedges,tagbound,lasttag;
181:
183: MPI_Comm_size(PETSC_COMM_WORLD,&size);
184: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
185:
186: /* copy A to Spooles' InpMtx object */
187: if ( lu->options.symflag == SPOOLES_NONSYMMETRIC ) {
188: Mat_MPIAIJ *mat = (Mat_MPIAIJ*)A->data;
189: Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(mat->A)->data;
190: Mat_SeqAIJ *bb=(Mat_SeqAIJ*)(mat->B)->data;
191: ai=aa->i; aj=aa->j; av=aa->a;
192: bi=bb->i; bj=bb->j; bv=bb->a;
193: lu->rstart = mat->rstart;
194: nz = aa->nz + bb->nz;
195: garray = mat->garray;
196: } else { /* SPOOLES_SYMMETRIC */
197: Mat_MPISBAIJ *mat = (Mat_MPISBAIJ*)A->data;
198: Mat_SeqSBAIJ *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
199: Mat_SeqBAIJ *bb=(Mat_SeqBAIJ*)(mat->B)->data;
200: ai=aa->i; aj=aa->j; av=aa->a;
201: bi=bb->i; bj=bb->j; bv=bb->a;
202: lu->rstart = mat->rstart;
203: nz = aa->s_nz + bb->nz;
204: garray = mat->garray;
205: }
206:
207: if(lu->flg == DIFFERENT_NONZERO_PATTERN) { lu->mtxA = InpMtx_new() ; }
208: InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, SPOOLES_REAL, nz, 0) ;
209: row = InpMtx_ivec1(lu->mtxA);
210: col = InpMtx_ivec2(lu->mtxA);
211: val = InpMtx_dvec(lu->mtxA);
212:
213: jj = 0; jB = 0; irow = lu->rstart;
214: for ( i=0; i<m; i++ ) {
215: ajj = aj + ai[i]; /* ptr to the beginning of this row */
216: countA = ai[i+1] - ai[i];
217: countB = bi[i+1] - bi[i];
218: bjj = bj + bi[i];
219:
220: if (lu->options.symflag == SPOOLES_NONSYMMETRIC ){
221: /* B part, smaller col index */
222: colA_start = lu->rstart + ajj[0]; /* the smallest col index for A */
223: for (j=0; j<countB; j++){
224: jcol = garray[bjj[j]];
225: if (jcol > colA_start) {
226: jB = j;
227: break;
228: }
229: row[jj] = irow; col[jj] = jcol; val[jj++] = *bv++;
230: if (j==countB-1) jB = countB;
231: }
232: }
233: /* A part */
234: for (j=0; j<countA; j++){
235: row[jj] = irow; col[jj] = lu->rstart + ajj[j]; val[jj++] = *av++;
236: }
237: /* B part, larger col index */
238: for (j=jB; j<countB; j++){
239: row[jj] = irow; col[jj] = garray[bjj[j]]; val[jj++] = *bv++;
240: }
241: irow++;
242: }
244: InpMtx_inputRealTriples(lu->mtxA, nz, row, col, val);
245: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
246: /* fprintf() terminates T3E. Error msg from totalview: no parameters.
247: if ( lu->options.msglvl > 2 ) {
248: fprintf(lu->options.msgFile, "nn input matrix") ;
249: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
250: fflush(lu->options.msgFile) ;
251: }
252: */
253: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
255: (*F)->ops->solve = MatSolve_MPIAIJ_Spooles;
256: (*F)->ops->destroy = MatDestroy_MPIAIJ_Spooles;
257: (*F)->assembled = PETSC_TRUE;
259: /* to be used by MatSolve() */
260: lu->mtxY = DenseMtx_new() ;
261: lu->mtxX = DenseMtx_new() ;
262: lu->scat = PETSC_NULL;
264: IVzero(20, lu->stats) ;
265: DVzero(20, lu->cpus) ;
266:
267: /* get input parameters */
268: SetSpoolesOptions(A, &lu->options);
270: /*
271: find a low-fill ordering
272: (1) create the Graph object
273: (2) order the graph using multiple minimum degree
274: (3) find out who has the best ordering w.r.t. op count,
275: and broadcast that front tree object
276: */
277: graph = Graph_new() ;
278: adjIVL = InpMtx_MPI_fullAdjacency(lu->mtxA, lu->stats,
279: lu->options.msglvl, lu->options.msgFile, MPI_COMM_WORLD) ;
280: nedges = IVL_tsize(adjIVL) ;
281: Graph_init2(graph, 0, M, 0, nedges, M, nedges, adjIVL, NULL, NULL) ;
282: if ( lu->options.msglvl > 2 ) {
283: fprintf(lu->options.msgFile, "nn graph of the input matrix") ;
284: Graph_writeForHumanEye(graph, lu->options.msgFile) ;
285: fflush(lu->options.msgFile) ;
286: }
288: switch (lu->options.ordering) {
289: case 0:
290: lu->frontETree = orderViaBestOfNDandMS(graph,
291: lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize,
292: lu->options.seed + rank, lu->options.msglvl, lu->options.msgFile); break;
293: case 1:
294: lu->frontETree = orderViaMMD(graph,lu->options.seed + rank,lu->options.msglvl,lu->options.msgFile); break;
295: case 2:
296: lu->frontETree = orderViaMS(graph, lu->options.maxdomainsize,
297: lu->options.seed + rank,lu->options.msglvl,lu->options.msgFile); break;
298: case 3:
299: lu->frontETree = orderViaND(graph, lu->options.maxdomainsize,
300: lu->options.seed + rank,lu->options.msglvl,lu->options.msgFile); break;
301: default:
302: SETERRQ(1,"Unknown Spooles's ordering");
303: }
305: Graph_free(graph) ;
306: if ( lu->options.msglvl > 2 ) {
307: fprintf(lu->options.msgFile, "nn front tree from ordering") ;
308: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
309: fflush(lu->options.msgFile) ;
310: }
312: opcounts = DVinit(size, 0.0) ;
313: opcounts[rank] = ETree_nFactorOps(lu->frontETree, SPOOLES_REAL, lu->options.symflag) ;
314: MPI_Allgather((void *) &opcounts[rank], 1, MPI_DOUBLE,
315: (void *) opcounts, 1, MPI_DOUBLE, MPI_COMM_WORLD) ;
316: minops = DVmin(size, opcounts, &root) ;
317: DVfree(opcounts) ;
318:
319: lu->frontETree = ETree_MPI_Bcast(lu->frontETree, root,
320: lu->options.msglvl, lu->options.msgFile, MPI_COMM_WORLD) ;
321: if ( lu->options.msglvl > 2 ) {
322: fprintf(lu->options.msgFile, "nn best front tree") ;
323: ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
324: fflush(lu->options.msgFile) ;
325: }
326:
327: /* get the permutations, permute the front tree, permute the matrix */
328: lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree) ;
329: lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree) ;
331: ETree_permuteVertices(lu->frontETree, lu->oldToNewIV) ;
333: InpMtx_permute(lu->mtxA, IV_entries(lu->oldToNewIV), IV_entries(lu->oldToNewIV)) ;
334:
335: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;
337: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
338: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
340: /* generate the owners map IV object and the map from vertices to owners */
341: cutoff = 1./(2*size) ;
342: cumopsDV = DV_new() ;
343: DV_init(cumopsDV, size, NULL) ;
344: lu->ownersIV = ETree_ddMap(lu->frontETree,
345: SPOOLES_REAL, lu->options.symflag, cumopsDV, cutoff) ;
346: DV_free(cumopsDV) ;
347: lu->vtxmapIV = IV_new() ;
348: IV_init(lu->vtxmapIV, M, NULL) ;
349: IVgather(M, IV_entries(lu->vtxmapIV),
350: IV_entries(lu->ownersIV), ETree_vtxToFront(lu->frontETree)) ;
351: if ( lu->options.msglvl > 2 ) {
352: fprintf(lu->options.msgFile, "nn map from fronts to owning processes") ;
353: IV_writeForHumanEye(lu->ownersIV, lu->options.msgFile) ;
354: fprintf(lu->options.msgFile, "nn map from vertices to owning processes") ;
355: IV_writeForHumanEye(lu->vtxmapIV, lu->options.msgFile) ;
356: fflush(lu->options.msgFile) ;
357: }
359: /* redistribute the matrix */
360: lu->firsttag = 0 ;
361: newA = InpMtx_MPI_split(lu->mtxA, lu->vtxmapIV, lu->stats,
362: lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
363: lu->firsttag += size ;
365: InpMtx_free(lu->mtxA) ;
366: lu->mtxA = newA ;
367: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
368: if ( lu->options.msglvl > 2 ) {
369: fprintf(lu->options.msgFile, "nn split InpMtx") ;
370: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
371: fflush(lu->options.msgFile) ;
372: }
373:
374: /* compute the symbolic factorization */
375: lu->symbfacIVL = SymbFac_MPI_initFromInpMtx(lu->frontETree, lu->ownersIV, lu->mtxA,
376: lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
377: lu->firsttag += lu->frontETree->nfront ;
378: if ( lu->options.msglvl > 2 ) {
379: fprintf(lu->options.msgFile, "nn local symbolic factorization") ;
380: IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile) ;
381: fflush(lu->options.msgFile) ;
382: }
384: lu->mtxmanager = SubMtxManager_new() ;
385: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;
386: lu->frontmtx = FrontMtx_new() ;
388: } else { /* new num factorization using previously computed symbolic factor */
389: if (lu->options.pivotingflag) { /* different FrontMtx is required */
390: FrontMtx_free(lu->frontmtx) ;
391: lu->frontmtx = FrontMtx_new() ;
392: }
394: SubMtxManager_free(lu->mtxmanager) ;
395: lu->mtxmanager = SubMtxManager_new() ;
396: SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;
398: /* permute mtxA */
399: InpMtx_permute(lu->mtxA, IV_entries(lu->oldToNewIV), IV_entries(lu->oldToNewIV)) ;
400: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;
401:
402: InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
403: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
405: /* redistribute the matrix */
406: MPI_Barrier(MPI_COMM_WORLD) ;
407: lu->firsttag = 0;
408: newA = InpMtx_MPI_split(lu->mtxA, lu->vtxmapIV, lu->stats,
409: lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
410: lu->firsttag += size ;
412: InpMtx_free(lu->mtxA) ;
413: lu->mtxA = newA ;
414: InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
415: if ( lu->options.msglvl > 2 ) {
416: fprintf(lu->options.msgFile, "nn split InpMtx") ;
417: InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
418: fflush(lu->options.msgFile) ;
419: }
420: } /* end of if ( lu->flg == DIFFERENT_NONZERO_PATTERN) */
422: FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, SPOOLES_REAL, lu->options.symflag,
423: FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, rank,
424: lu->ownersIV, lu->mtxmanager, lu->options.msglvl, lu->options.msgFile) ;
426: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
427: if ( lu->options.patchAndGoFlag == 1 ) {
428: lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
429: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
430: lu->options.storeids, lu->options.storevalues) ;
431: } else if ( lu->options.patchAndGoFlag == 2 ) {
432: lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
433: PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
434: lu->options.storeids, lu->options.storevalues) ;
435: }
436: }
438: /* numerical factorization */
439: chvmanager = ChvManager_new() ;
440: ChvManager_init(chvmanager, NO_LOCK, 0) ;
442: tagbound = maxTagMPI(MPI_COMM_WORLD) ;
443: lasttag = lu->firsttag + 3*lu->frontETree->nfront + 2;
444: /* if(!rank) PetscPrintf(PETSC_COMM_SELF,"n firsttag: %d, nfront: %dn",lu->firsttag, lu->frontETree->nfront);*/
445: if ( lasttag > tagbound ) {
446: SETERRQ3(1,"fatal error in FrontMtx_MPI_factorInpMtx(), tag range is [%d,%d], tag_bound = %d",447: lu->firsttag, lasttag, tagbound) ;
448: }
449: rootchv = FrontMtx_MPI_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, droptol,
450: chvmanager, lu->ownersIV, lookahead, &ierr, lu->cpus,
451: lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
452: ChvManager_free(chvmanager) ;
453: lu->firsttag = lasttag;
454: if ( lu->options.msglvl > 2 ) {
455: fprintf(lu->options.msgFile, "nn numeric factorization") ;
456: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
457: fflush(lu->options.msgFile) ;
458: }
460: if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
461: if ( lu->options.patchAndGoFlag == 1 ) {
462: if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
463: if (lu->options.msglvl > 0 ){
464: fprintf(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: fprintf(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: fprintf(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: }
483: if ( ierr >= 0 ) SETERRQ2(1,"n proc %d : factorization error at front %d", rank, ierr) ;
484:
485: /* post-process the factorization and split
486: the factor matrices into submatrices */
487: lasttag = lu->firsttag + 5*size;
488: if ( lasttag > tagbound ) {
489: SETERRQ3(1,"fatal error in FrontMtx_MPI_postProcess(), tag range is [%d,%d], tag_bound = %d",490: lu->firsttag, lasttag, tagbound) ;
491: }
492: FrontMtx_MPI_postProcess(lu->frontmtx, lu->ownersIV, lu->stats, lu->options.msglvl,
493: lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
494: lu->firsttag += 5*size ;
495: if ( lu->options.msglvl > 2 ) {
496: fprintf(lu->options.msgFile, "nn numeric factorization after post-processing");
497: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
498: fflush(lu->options.msgFile) ;
499: }
500:
501: /* create the solve map object */
502: lu->solvemap = SolveMap_new() ;
503: SolveMap_ddMap(lu->solvemap, lu->frontmtx->symmetryflag,
504: FrontMtx_upperBlockIVL(lu->frontmtx),
505: FrontMtx_lowerBlockIVL(lu->frontmtx),
506: size, lu->ownersIV, FrontMtx_frontTree(lu->frontmtx),
507: lu->options.seed, lu->options.msglvl, lu->options.msgFile);
508: if ( lu->options.msglvl > 2 ) {
509: SolveMap_writeForHumanEye(lu->solvemap, lu->options.msgFile) ;
510: fflush(lu->options.msgFile) ;
511: }
513: /* redistribute the submatrices of the factors */
514: FrontMtx_MPI_split(lu->frontmtx, lu->solvemap,
515: lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
516: if ( lu->options.msglvl > 2 ) {
517: fprintf(lu->options.msgFile, "nn numeric factorization after split") ;
518: FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
519: fflush(lu->options.msgFile) ;
520: }
522: /* create a solution DenseMtx object */
523: lu->ownedColumnsIV = FrontMtx_ownedColumnsIV(lu->frontmtx, rank, lu->ownersIV,
524: lu->options.msglvl, lu->options.msgFile) ;
525: lu->nmycol = IV_size(lu->ownedColumnsIV) ;
526: if ( lu->nmycol > 0) {
527: DenseMtx_init(lu->mtxX, SPOOLES_REAL, 0, 0, lu->nmycol, 1, 1, lu->nmycol) ;
528: /* get pointers rowindX and entX */
529: DenseMtx_rowIndices(lu->mtxX, &lu->nmycol, &lu->rowindX);
530: lu->entX = DenseMtx_entries(lu->mtxX) ;
531: } else { /* lu->nmycol == 0 */
532: lu->entX = 0;
533: lu->rowindX = 0;
534: }
536: if ( lu->scat ){
537: VecDestroy(lu->vec_spooles);
538: ISDestroy(lu->iden);
539: ISDestroy(lu->is_petsc);
540: VecScatterDestroy(lu->scat);
541: }
542: lu->scat = PETSC_NULL;
543:
544: lu->flg = SAME_NONZERO_PATTERN;
545: return(0);
546: }
548: #endif