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: int MatDestroy_MPIAIJ_Spooles(Mat A)
 19: {
 20:   Mat_Spooles   *lu = (Mat_Spooles*)A->spptr;
 21:   int           ierr;
 22: 
 24: 
 25:   FrontMtx_free(lu->frontmtx) ;
 26:   IV_free(lu->newToOldIV) ;
 27:   IV_free(lu->oldToNewIV) ;
 28:   IV_free(lu->vtxmapIV) ;
 29:   InpMtx_free(lu->mtxA) ;
 30:   ETree_free(lu->frontETree) ;
 31:   IVL_free(lu->symbfacIVL) ;
 32:   SubMtxManager_free(lu->mtxmanager) ;
 33:   DenseMtx_free(lu->mtxX) ;
 34:   DenseMtx_free(lu->mtxY) ;
 35: 
 36:   VecDestroy(lu->vec_spooles);
 37:   ISDestroy(lu->iden);
 38:   ISDestroy(lu->is_petsc);
 39:   VecScatterDestroy(lu->scat);

 41:   PetscFree(lu);
 42:   MatDestroy_MPIAIJ(A);

 44:   return(0);
 45: }

 47: int MatSolve_MPIAIJ_Spooles(Mat A,Vec b,Vec x)
 48: {
 49:   Mat_Spooles      *lu = (Mat_Spooles*)A->spptr;
 50:   int              ierr,size,rank,m=A->m,irow,*rowindY;
 51:   PetscScalar      *array;
 52:   DenseMtx         *newY ;
 53:   SubMtxManager    *solvemanager ;

 56:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 57:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 58: 
 59:   /* copy b into spooles' rhs mtxY */
 60:   DenseMtx_init(lu->mtxY, SPOOLES_REAL, 0, 0, m, 1, 1, m) ;
 61:   VecGetArray(b,&array);

 63:   /* doesn't work! 
 64:   PetscMalloc(m*sizeof(int),&rowindY);
 65:   for ( irow = 0 ; irow < m ; irow++ ) rowindY[irow] = irow + lu->rstart; 
 66:   int colind=0;
 67:   DenseMtx_initWithPointers(lu->mtxY,SPOOLES_REAL,0,0,m,1,1,m,rowindY,&colind,array); 
 68:   */

 70:   DenseMtx_rowIndices(lu->mtxY, &m, &rowindY) ;  /* get m, rowind */
 71:   for ( irow = 0 ; irow < m ; irow++ ) {
 72:     rowindY[irow] = irow + lu->rstart;           /* global rowind */
 73:     DenseMtx_setRealEntry(lu->mtxY, irow, 0, *array++) ;
 74:   }
 75:   /* DenseMtx_column(lu->mtxY, 0, &array);  doesn't work! */
 76:   VecRestoreArray(b,&array);
 77: 
 78:   if ( lu->options.msglvl > 2 ) {
 79:     fprintf(lu->options.msgFile, "nn 1 matrix in original ordering") ;
 80:     DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
 81:     fflush(lu->options.msgFile) ;
 82:   }
 83: 
 84:   /* permute and redistribute Y if necessary */
 85:   DenseMtx_permuteRows(lu->mtxY, lu->oldToNewIV) ;
 86:   if ( lu->options.msglvl > 2 ) {
 87:     fprintf(lu->options.msgFile, "nn rhs matrix in new ordering") ;
 88:     DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
 89:    fflush(lu->options.msgFile) ;
 90:   }
 91: 
 92:   newY = DenseMtx_MPI_splitByRows(lu->mtxY, lu->vtxmapIV, lu->stats, lu->options.msglvl,
 93:                                 lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
 94:   DenseMtx_free(lu->mtxY) ;
 95:   lu->mtxY = newY ;
 96:   lu->firsttag += size ;
 97:   if ( lu->options.msglvl > 2 ) {
 98:     fprintf(lu->options.msgFile, "nn split DenseMtx Y") ;
 99:     DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
100:     fflush(lu->options.msgFile) ;
101:   }

103:   if ( FRONTMTX_IS_PIVOTING(lu->frontmtx) ) {
104:     /*   pivoting has taken place, redistribute the right hand side
105:          to match the final rows and columns in the fronts             */
106:     IV *rowmapIV ;
107:     rowmapIV = FrontMtx_MPI_rowmapIV(lu->frontmtx, lu->ownersIV, lu->options.msglvl,
108:                                     lu->options.msgFile, MPI_COMM_WORLD) ;
109:     newY = DenseMtx_MPI_splitByRows(lu->mtxY, rowmapIV, lu->stats, lu->options.msglvl,
110:                                    lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
111:     DenseMtx_free(lu->mtxY) ;
112:     lu->mtxY = newY ;
113:     IV_free(rowmapIV) ;
114:   }
115:   if ( lu->options.msglvl > 2 ) {
116:     fprintf(lu->options.msgFile, "nn rhs matrix after split") ;
117:     DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
118:     fflush(lu->options.msgFile) ;
119:   }

121:   if ( lu->nmycol > 0 ) IVcopy(lu->nmycol,lu->rowindX,IV_entries(lu->ownedColumnsIV)); /* must do for each solve */
122: 
123:   /* solve the linear system */
124:   solvemanager = SubMtxManager_new() ;
125:   SubMtxManager_init(solvemanager, NO_LOCK, 0) ;
126:   FrontMtx_MPI_solve(lu->frontmtx, lu->mtxX, lu->mtxY, solvemanager, lu->solvemap, lu->cpus,
127:                    lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
128:   SubMtxManager_free(solvemanager) ;
129:   if ( lu->options.msglvl > 2 ) {
130:     fprintf(lu->options.msgFile, "n solution in new ordering") ;
131:     DenseMtx_writeForHumanEye(lu->mtxX, lu->options.msgFile) ;
132:   }

134:   /* permute the solution into the original ordering */
135:   DenseMtx_permuteRows(lu->mtxX, lu->newToOldIV) ;
136:   if ( lu->options.msglvl > 2 ) {
137:     fprintf(lu->options.msgFile, "nn solution in old ordering") ;
138:     DenseMtx_writeForHumanEye(lu->mtxX, lu->options.msgFile) ;
139:     fflush(lu->options.msgFile) ;
140:   }
141: 
142:   /* scatter local solution mtxX into mpi vector x */
143:   if( !lu->scat ){ /* create followings once for each numfactorization */
144:     VecCreateSeqWithArray(PETSC_COMM_SELF,lu->nmycol,lu->entX,&lu->vec_spooles); /* vec_spooles <- mtxX */
145:     ISCreateStride(PETSC_COMM_SELF,lu->nmycol,0,1,&lu->iden);
146:     ISCreateGeneral(PETSC_COMM_SELF,lu->nmycol,lu->rowindX,&lu->is_petsc);
147:     VecScatterCreate(lu->vec_spooles,lu->iden,x,lu->is_petsc,&lu->scat);
148:   }

150:   VecScatterBegin(lu->vec_spooles,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
151:   VecScatterEnd(lu->vec_spooles,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
152: 
153:   return(0);
154: }

156: int MatFactorNumeric_MPIAIJ_Spooles(Mat A,Mat *F)
157: {
158:   Mat_Spooles     *lu = (Mat_Spooles*)(*F)->spptr;
159:   int             rank,size,ierr,lookahead=0;
160:   ChvManager      *chvmanager ;
161:   Chv             *rootchv ;
162:   Graph           *graph ;
163:   IVL             *adjIVL;
164:   DV              *cumopsDV ;
165:   double          droptol=0.0,*opcounts,minops,cutoff,*val;
166:   InpMtx          *newA ;
167:   PetscScalar     *av, *bv;
168:   int             *ai, *aj, *bi,*bj, nz, *ajj, *bjj, *garray,
169:                   i,j,irow,jcol,countA,countB,jB,*row,*col,colA_start,jj;
170:   int             M=A->M,m=A->m,root,nedges;
171: 
173:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
174:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
175: 
176:   /* copy A to Spooles' InpMtx object */
177:   if ( lu->options.symflag == SPOOLES_NONSYMMETRIC ) {
178:     Mat_MPIAIJ  *mat =  (Mat_MPIAIJ*)A->data;
179:     Mat_SeqAIJ  *aa=(Mat_SeqAIJ*)(mat->A)->data;
180:     Mat_SeqAIJ  *bb=(Mat_SeqAIJ*)(mat->B)->data;
181:     ai=aa->i; aj=aa->j; av=aa->a;
182:     bi=bb->i; bj=bb->j; bv=bb->a;
183:     lu->rstart = mat->rstart;
184:     nz         = aa->nz + bb->nz;
185:     garray     = mat->garray;
186:   } else {         /* SPOOLES_SYMMETRIC  */
187:     Mat_MPISBAIJ  *mat = (Mat_MPISBAIJ*)A->data;
188:     Mat_SeqSBAIJ  *aa=(Mat_SeqSBAIJ*)(mat->A)->data;
189:     Mat_SeqBAIJ    *bb=(Mat_SeqBAIJ*)(mat->B)->data;
190:     ai=aa->i; aj=aa->j; av=aa->a;
191:     bi=bb->i; bj=bb->j; bv=bb->a;
192:     lu->rstart = mat->rstart;
193:     nz         = aa->s_nz + bb->nz;
194:     garray     = mat->garray;
195:   }
196: 
197:   if(lu->flg == DIFFERENT_NONZERO_PATTERN) { lu->mtxA = InpMtx_new() ; }
198:   InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, SPOOLES_REAL, nz, 0) ;
199:   row   = InpMtx_ivec1(lu->mtxA);
200:   col   = InpMtx_ivec2(lu->mtxA);
201:   val   = InpMtx_dvec(lu->mtxA);
202: 
203:   jj = 0; jB = 0; irow = lu->rstart;
204:   for ( i=0; i<m; i++ ) {
205:     ajj = aj + ai[i];                 /* ptr to the beginning of this row */
206:     countA = ai[i+1] - ai[i];
207:     countB = bi[i+1] - bi[i];
208:     bjj = bj + bi[i];
209: 
210:     if (lu->options.symflag == SPOOLES_NONSYMMETRIC ){
211:       /* B part, smaller col index */
212:       colA_start = lu->rstart + ajj[0]; /* the smallest col index for A */
213:       for (j=0; j<countB; j++){
214:         jcol = garray[bjj[j]];
215:         if (jcol > colA_start) {
216:           jB = j;
217:           break;
218:         }
219:         row[jj] = irow; col[jj] = jcol; val[jj++] = *bv++;
220:         if (j==countB-1) jB = countB;
221:       }
222:     }
223:     /* A part */
224:     for (j=0; j<countA; j++){
225:       row[jj] = irow; col[jj] = lu->rstart + ajj[j]; val[jj++] = *av++;
226:     }
227:     /* B part, larger col index */
228:     for (j=jB; j<countB; j++){
229:       row[jj] = irow; col[jj] = garray[bjj[j]]; val[jj++] = *bv++;
230:     }
231:     irow++;
232:   }

234:   InpMtx_inputRealTriples(lu->mtxA, nz, row, col, val);
235:   InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
236:   if ( lu->options.msglvl > 2 ) {
237:     fprintf(lu->options.msgFile, "nn input matrix") ;
238:     InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
239:     fflush(lu->options.msgFile) ;
240:   }

242:   if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */

244:     (*F)->ops->solve   = MatSolve_MPIAIJ_Spooles;
245:     (*F)->ops->destroy = MatDestroy_MPIAIJ_Spooles;
246:     (*F)->assembled    = PETSC_TRUE;

248:     /* to be used by MatSolve() */
249:     lu->mtxY = DenseMtx_new() ;
250:     lu->mtxX = DenseMtx_new() ;
251:     lu->scat = PETSC_NULL;

253:     IVzero(20, lu->stats) ;
254:     DVzero(20, lu->cpus) ;
255: 
256:     /* get input parameters */
257:     SetSpoolesOptions(A, &lu->options);

259:     /*
260:       find a low-fill ordering
261:       (1) create the Graph object
262:       (2) order the graph using multiple minimum degree
263:       (3) find out who has the best ordering w.r.t. op count,
264:           and broadcast that front tree object
265:     */
266:     graph = Graph_new() ;
267:     adjIVL = InpMtx_MPI_fullAdjacency(lu->mtxA, lu->stats,
268:               lu->options.msglvl, lu->options.msgFile, MPI_COMM_WORLD) ;
269:     nedges = IVL_tsize(adjIVL) ;
270:     Graph_init2(graph, 0, M, 0, nedges, M, nedges, adjIVL, NULL, NULL) ;
271:     if ( lu->options.msglvl > 2 ) {
272:       fprintf(lu->options.msgFile, "nn graph of the input matrix") ;
273:       Graph_writeForHumanEye(graph, lu->options.msgFile) ;
274:       fflush(lu->options.msgFile) ;
275:     }

277:     switch (lu->options.ordering) {
278:     case 0:
279:       lu->frontETree = orderViaBestOfNDandMS(graph,
280:                      lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize,
281:                      lu->options.seed + rank, lu->options.msglvl, lu->options.msgFile); break;
282:     case 1:
283:       lu->frontETree = orderViaMMD(graph,lu->options.seed + rank,lu->options.msglvl,lu->options.msgFile); break;
284:     case 2:
285:       lu->frontETree = orderViaMS(graph, lu->options.maxdomainsize,
286:                      lu->options.seed + rank,lu->options.msglvl,lu->options.msgFile); break;
287:     case 3:
288:       lu->frontETree = orderViaND(graph, lu->options.maxdomainsize,
289:                      lu->options.seed + rank,lu->options.msglvl,lu->options.msgFile); break;
290:     default:
291:       SETERRQ(1,"Unknown Spooles's ordering");
292:     }

294:     Graph_free(graph) ;
295:     if ( lu->options.msglvl > 2 ) {
296:       fprintf(lu->options.msgFile, "nn front tree from ordering") ;
297:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
298:       fflush(lu->options.msgFile) ;
299:     }

301:     opcounts = DVinit(size, 0.0) ;
302:     opcounts[rank] = ETree_nFactorOps(lu->frontETree, SPOOLES_REAL, lu->options.symflag) ;
303:     MPI_Allgather((void *) &opcounts[rank], 1, MPI_DOUBLE,
304:               (void *) opcounts, 1, MPI_DOUBLE, MPI_COMM_WORLD) ;
305:     minops = DVmin(size, opcounts, &root) ;
306:     DVfree(opcounts) ;
307: 
308:     lu->frontETree = ETree_MPI_Bcast(lu->frontETree, root,
309:                              lu->options.msglvl, lu->options.msgFile, MPI_COMM_WORLD) ;
310:     if ( lu->options.msglvl > 2 ) {
311:       fprintf(lu->options.msgFile, "nn best front tree") ;
312:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
313:       fflush(lu->options.msgFile) ;
314:     }
315: 
316:     /* get the permutations, permute the front tree, permute the matrix */
317:     lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree) ;
318:     lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree) ;

320:     ETree_permuteVertices(lu->frontETree, lu->oldToNewIV) ;

322:     InpMtx_permute(lu->mtxA, IV_entries(lu->oldToNewIV), IV_entries(lu->oldToNewIV)) ;
323: 
324:     if (  lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;

326:     InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
327:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;

329:     /* generate the owners map IV object and the map from vertices to owners */
330:     cutoff   = 1./(2*size) ;
331:     cumopsDV = DV_new() ;
332:     DV_init(cumopsDV, size, NULL) ;
333:     lu->ownersIV = ETree_ddMap(lu->frontETree,
334:                        SPOOLES_REAL, lu->options.symflag, cumopsDV, cutoff) ;
335:     DV_free(cumopsDV) ;
336:     lu->vtxmapIV = IV_new() ;
337:     IV_init(lu->vtxmapIV, M, NULL) ;
338:     IVgather(M, IV_entries(lu->vtxmapIV),
339:              IV_entries(lu->ownersIV), ETree_vtxToFront(lu->frontETree)) ;
340:     if ( lu->options.msglvl > 2 ) {
341:       fprintf(lu->options.msgFile, "nn map from fronts to owning processes") ;
342:       IV_writeForHumanEye(lu->ownersIV, lu->options.msgFile) ;
343:       fprintf(lu->options.msgFile, "nn map from vertices to owning processes") ;
344:       IV_writeForHumanEye(lu->vtxmapIV, lu->options.msgFile) ;
345:       fflush(lu->options.msgFile) ;
346:     }

348:     /* redistribute the matrix */
349:     lu->firsttag = 0 ;
350:     newA = InpMtx_MPI_split(lu->mtxA, lu->vtxmapIV, lu->stats,
351:                         lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
352:     lu->firsttag++ ;

354:     InpMtx_free(lu->mtxA) ;
355:     lu->mtxA = newA ;
356:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
357:     if ( lu->options.msglvl > 2 ) {
358:       fprintf(lu->options.msgFile, "nn split InpMtx") ;
359:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
360:       fflush(lu->options.msgFile) ;
361:     }
362: 
363:     /* compute the symbolic factorization */
364:     lu->symbfacIVL = SymbFac_MPI_initFromInpMtx(lu->frontETree, lu->ownersIV, lu->mtxA,
365:                      lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
366:     lu->firsttag += lu->frontETree->nfront ;
367:     if ( lu->options.msglvl > 2 ) {
368:       fprintf(lu->options.msgFile, "nn local symbolic factorization") ;
369:       IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile) ;
370:       fflush(lu->options.msgFile) ;
371:     }

373:     lu->mtxmanager = SubMtxManager_new() ;
374:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;
375:     lu->frontmtx = FrontMtx_new() ;

377:   } else { /* new num factorization using previously computed symbolic factor */
378:     if (lu->options.pivotingflag) {                  /* different FrontMtx is required */
379:       FrontMtx_free(lu->frontmtx) ;
380:       lu->frontmtx   = FrontMtx_new() ;
381:     }

383:     SubMtxManager_free(lu->mtxmanager) ;
384:     lu->mtxmanager = SubMtxManager_new() ;
385:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;

387:     /* permute mtxA */
388:     InpMtx_permute(lu->mtxA, IV_entries(lu->oldToNewIV), IV_entries(lu->oldToNewIV)) ;
389:     if ( lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;
390: 
391:     InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
392:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;

394:     /* redistribute the matrix */
395:     newA = InpMtx_MPI_split(lu->mtxA, lu->vtxmapIV, lu->stats,
396:                         lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
397:     lu->firsttag++ ;

399:     InpMtx_free(lu->mtxA) ;
400:     lu->mtxA = newA ;
401:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
402:     if ( lu->options.msglvl > 2 ) {
403:       fprintf(lu->options.msgFile, "nn split InpMtx") ;
404:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
405:       fflush(lu->options.msgFile) ;
406:     }
407:   } /* end of if ( lu->flg == DIFFERENT_NONZERO_PATTERN) */

409:   FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, SPOOLES_REAL, lu->options.symflag,
410:               FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, rank,
411:               lu->ownersIV, lu->mtxmanager, lu->options.msglvl, lu->options.msgFile) ;

413:     if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
414:     if ( lu->options.patchAndGoFlag == 1 ) {
415:       lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
416:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
417:                        lu->options.storeids, lu->options.storevalues) ;
418:     } else if ( lu->options.patchAndGoFlag == 2 ) {
419:       lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
420:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
421:                        lu->options.storeids, lu->options.storevalues) ;
422:     }
423:   }

425:   /* numerical factorization */
426:   chvmanager = ChvManager_new() ;
427:   ChvManager_init(chvmanager, NO_LOCK, 0) ;
428:   rootchv = FrontMtx_MPI_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, droptol,
429:                      chvmanager, lu->ownersIV, lookahead, &ierr, lu->cpus,
430:                      lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
431:   ChvManager_free(chvmanager) ;
432:   lu->firsttag += 3*lu->frontETree->nfront + 2 ;
433:   if ( lu->options.msglvl > 2 ) {
434:     fprintf(lu->options.msgFile, "nn numeric factorization") ;
435:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
436:     fflush(lu->options.msgFile) ;
437:   }

439:   if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
440:     if ( lu->options.patchAndGoFlag == 1 ) {
441:       if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
442:         if (lu->options.msglvl > 0 ){
443:           fprintf(lu->options.msgFile, "n small pivots found at these locations") ;
444:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
445:         }
446:       }
447:       PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
448:     } else if ( lu->options.patchAndGoFlag == 2 ) {
449:       if (lu->options.msglvl > 0 ){
450:         if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
451:           fprintf(lu->options.msgFile, "n small pivots found at these locations") ;
452:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
453:         }
454:         if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
455:           fprintf(lu->options.msgFile, "n perturbations") ;
456:           DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile) ;
457:         }
458:       }
459:       PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
460:     }
461:   }
462:   if ( ierr >= 0 ) SETERRQ2(1,"n proc %d : factorization error at front %d", rank, ierr) ;
463: 
464:   /*  post-process the factorization and split 
465:       the factor matrices into submatrices */
466:   FrontMtx_MPI_postProcess(lu->frontmtx, lu->ownersIV, lu->stats, lu->options.msglvl,
467:                          lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
468:   lu->firsttag += 5*size ;
469:   if ( lu->options.msglvl > 2 ) {
470:     fprintf(lu->options.msgFile, "nn numeric factorization after post-processing");
471:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
472:     fflush(lu->options.msgFile) ;
473:   }
474: 
475:   /* create the solve map object */
476:   lu->solvemap = SolveMap_new() ;
477:   SolveMap_ddMap(lu->solvemap, lu->frontmtx->symmetryflag,
478:                FrontMtx_upperBlockIVL(lu->frontmtx),
479:                FrontMtx_lowerBlockIVL(lu->frontmtx),
480:                size, lu->ownersIV, FrontMtx_frontTree(lu->frontmtx),
481:                lu->options.seed, lu->options.msglvl, lu->options.msgFile);
482:   if ( lu->options.msglvl > 2 ) {
483:     SolveMap_writeForHumanEye(lu->solvemap, lu->options.msgFile) ;
484:     fflush(lu->options.msgFile) ;
485:   }

487:   /* redistribute the submatrices of the factors */
488:   FrontMtx_MPI_split(lu->frontmtx, lu->solvemap,
489:                    lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
490:   if ( lu->options.msglvl > 2 ) {
491:     fprintf(lu->options.msgFile, "nn numeric factorization after split") ;
492:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
493:     fflush(lu->options.msgFile) ;
494:   }

496:   /* create a solution DenseMtx object */
497:   lu->ownedColumnsIV = FrontMtx_ownedColumnsIV(lu->frontmtx, rank, lu->ownersIV,
498:                                          lu->options.msglvl, lu->options.msgFile) ;
499:   lu->nmycol = IV_size(lu->ownedColumnsIV) ;
500:   if ( lu->nmycol > 0) {
501:     DenseMtx_init(lu->mtxX, SPOOLES_REAL, 0, 0, lu->nmycol, 1, 1, lu->nmycol) ;
502:     /* get pointers rowindX and entX */
503:     DenseMtx_rowIndices(lu->mtxX, &lu->nmycol, &lu->rowindX);
504:     lu->entX = DenseMtx_entries(lu->mtxX) ;
505:   } else { /* lu->nmycol == 0 */
506:     lu->entX    = 0;
507:     lu->rowindX = 0;
508:   }

510:   if ( lu->scat ){
511:     VecDestroy(lu->vec_spooles);
512:     ISDestroy(lu->iden);
513:     ISDestroy(lu->is_petsc);
514:     VecScatterDestroy(lu->scat);
515:   }
516:   lu->scat = PETSC_NULL;
517: 
518:   lu->flg = SAME_NONZERO_PATTERN;
519:   return(0);
520: }

522: #endif