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:   }

 92:   MPI_Barrier(MPI_COMM_WORLD) ; /* for initializing firsttag, because the num. of tags used
 93:                                    by FrontMtx_MPI_split() is unknown */
 94:   lu->firsttag = 0;
 95:   newY = DenseMtx_MPI_splitByRows(lu->mtxY, lu->vtxmapIV, lu->stats, lu->options.msglvl,
 96:                                 lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
 97:   DenseMtx_free(lu->mtxY) ;
 98:   lu->mtxY = newY ;
 99:   lu->firsttag += size ;
100:   if ( lu->options.msglvl > 2 ) {
101:     fprintf(lu->options.msgFile, "nn split DenseMtx Y") ;
102:     DenseMtx_writeForHumanEye(lu->mtxY, lu->options.msgFile) ;
103:     fflush(lu->options.msgFile) ;
104:   }

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

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

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

154:   VecScatterBegin(lu->vec_spooles,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
155:   VecScatterEnd(lu->vec_spooles,x,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
156: 
157:   return(0);
158: }

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

238:   InpMtx_inputRealTriples(lu->mtxA, nz, row, col, val);
239:   InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
240:   /* fprintf() terminates T3E. Error msg from totalview: no parameters.
241:   if ( lu->options.msglvl > 2 ) {
242:     fprintf(lu->options.msgFile, "nn input matrix") ;
243:     InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
244:     fflush(lu->options.msgFile) ;
245:   }
246:   */
247:   if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */

249:     (*F)->ops->solve   = MatSolve_MPIAIJ_Spooles;
250:     (*F)->ops->destroy = MatDestroy_MPIAIJ_Spooles;
251:     (*F)->assembled    = PETSC_TRUE;

253:     /* to be used by MatSolve() */
254:     lu->mtxY = DenseMtx_new() ;
255:     lu->mtxX = DenseMtx_new() ;
256:     lu->scat = PETSC_NULL;

258:     IVzero(20, lu->stats) ;
259:     DVzero(20, lu->cpus) ;
260: 
261:     /* get input parameters */
262:     SetSpoolesOptions(A, &lu->options);

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

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

299:     Graph_free(graph) ;
300:     if ( lu->options.msglvl > 2 ) {
301:       fprintf(lu->options.msgFile, "nn front tree from ordering") ;
302:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
303:       fflush(lu->options.msgFile) ;
304:     }

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

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

327:     InpMtx_permute(lu->mtxA, IV_entries(lu->oldToNewIV), IV_entries(lu->oldToNewIV)) ;
328: 
329:     if (  lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;

331:     InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
332:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;

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

353:     /* redistribute the matrix */
354:     lu->firsttag = 0 ;
355:     newA = InpMtx_MPI_split(lu->mtxA, lu->vtxmapIV, lu->stats,
356:                         lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
357:     lu->firsttag += size ;

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

378:     lu->mtxmanager = SubMtxManager_new() ;
379:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;
380:     lu->frontmtx = FrontMtx_new() ;

382:   } else { /* new num factorization using previously computed symbolic factor */
383:     if (lu->options.pivotingflag) {                  /* different FrontMtx is required */
384:       FrontMtx_free(lu->frontmtx) ;
385:       lu->frontmtx   = FrontMtx_new() ;
386:     }

388:     SubMtxManager_free(lu->mtxmanager) ;
389:     lu->mtxmanager = SubMtxManager_new() ;
390:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;

392:     /* permute mtxA */
393:     InpMtx_permute(lu->mtxA, IV_entries(lu->oldToNewIV), IV_entries(lu->oldToNewIV)) ;
394:     if ( lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;
395: 
396:     InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
397:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;

399:     /* redistribute the matrix */
400:     MPI_Barrier(MPI_COMM_WORLD) ;
401:     lu->firsttag = 0;
402:     newA = InpMtx_MPI_split(lu->mtxA, lu->vtxmapIV, lu->stats,
403:                         lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
404:     lu->firsttag += size ;

406:     InpMtx_free(lu->mtxA) ;
407:     lu->mtxA = newA ;
408:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
409:     if ( lu->options.msglvl > 2 ) {
410:       fprintf(lu->options.msgFile, "nn split InpMtx") ;
411:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
412:       fflush(lu->options.msgFile) ;
413:     }
414:   } /* end of if ( lu->flg == DIFFERENT_NONZERO_PATTERN) */

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

420:     if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
421:     if ( lu->options.patchAndGoFlag == 1 ) {
422:       lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
423:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
424:                        lu->options.storeids, lu->options.storevalues) ;
425:     } else if ( lu->options.patchAndGoFlag == 2 ) {
426:       lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
427:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
428:                        lu->options.storeids, lu->options.storevalues) ;
429:     }
430:   }

432:   /* numerical factorization */
433:   chvmanager = ChvManager_new() ;
434:   ChvManager_init(chvmanager, NO_LOCK, 0) ;

436:   tagbound = maxTagMPI(MPI_COMM_WORLD) ;
437:   lasttag  = lu->firsttag + 3*lu->frontETree->nfront + 2;
438:   /* if(!rank) PetscPrintf(PETSC_COMM_SELF,"n firsttag: %d, nfront: %dn",lu->firsttag, lu->frontETree->nfront);*/
439:   if ( lasttag > tagbound ) {
440:       SETERRQ3(1,"fatal error in FrontMtx_MPI_factorInpMtx(), tag range is [%d,%d], tag_bound = %d",441:                lu->firsttag, lasttag, tagbound) ;
442:   }
443:   rootchv = FrontMtx_MPI_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, droptol,
444:                      chvmanager, lu->ownersIV, lookahead, &ierr, lu->cpus,
445:                      lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
446:   ChvManager_free(chvmanager) ;
447:   lu->firsttag = lasttag;
448:   if ( lu->options.msglvl > 2 ) {
449:     fprintf(lu->options.msgFile, "nn numeric factorization") ;
450:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
451:     fflush(lu->options.msgFile) ;
452:   }

454:   if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
455:     if ( lu->options.patchAndGoFlag == 1 ) {
456:       if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
457:         if (lu->options.msglvl > 0 ){
458:           fprintf(lu->options.msgFile, "n small pivots found at these locations") ;
459:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
460:         }
461:       }
462:       PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
463:     } else if ( lu->options.patchAndGoFlag == 2 ) {
464:       if (lu->options.msglvl > 0 ){
465:         if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
466:           fprintf(lu->options.msgFile, "n small pivots found at these locations") ;
467:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
468:         }
469:         if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
470:           fprintf(lu->options.msgFile, "n perturbations") ;
471:           DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile) ;
472:         }
473:       }
474:       PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
475:     }
476:   }
477:   if ( ierr >= 0 ) SETERRQ2(1,"n proc %d : factorization error at front %d", rank, ierr) ;
478: 
479:   /*  post-process the factorization and split 
480:       the factor matrices into submatrices */
481:   lasttag  = lu->firsttag + 5*size;
482:   if ( lasttag > tagbound ) {
483:       SETERRQ3(1,"fatal error in FrontMtx_MPI_postProcess(), tag range is [%d,%d], tag_bound = %d",484:                lu->firsttag, lasttag, tagbound) ;
485:   }
486:   FrontMtx_MPI_postProcess(lu->frontmtx, lu->ownersIV, lu->stats, lu->options.msglvl,
487:                          lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
488:   lu->firsttag += 5*size ;
489:   if ( lu->options.msglvl > 2 ) {
490:     fprintf(lu->options.msgFile, "nn numeric factorization after post-processing");
491:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
492:     fflush(lu->options.msgFile) ;
493:   }
494: 
495:   /* create the solve map object */
496:   lu->solvemap = SolveMap_new() ;
497:   SolveMap_ddMap(lu->solvemap, lu->frontmtx->symmetryflag,
498:                FrontMtx_upperBlockIVL(lu->frontmtx),
499:                FrontMtx_lowerBlockIVL(lu->frontmtx),
500:                size, lu->ownersIV, FrontMtx_frontTree(lu->frontmtx),
501:                lu->options.seed, lu->options.msglvl, lu->options.msgFile);
502:   if ( lu->options.msglvl > 2 ) {
503:     SolveMap_writeForHumanEye(lu->solvemap, lu->options.msgFile) ;
504:     fflush(lu->options.msgFile) ;
505:   }

507:   /* redistribute the submatrices of the factors */
508:   FrontMtx_MPI_split(lu->frontmtx, lu->solvemap,
509:                    lu->stats, lu->options.msglvl, lu->options.msgFile, lu->firsttag, MPI_COMM_WORLD) ;
510:   if ( lu->options.msglvl > 2 ) {
511:     fprintf(lu->options.msgFile, "nn numeric factorization after split") ;
512:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
513:     fflush(lu->options.msgFile) ;
514:   }

516:   /* create a solution DenseMtx object */
517:   lu->ownedColumnsIV = FrontMtx_ownedColumnsIV(lu->frontmtx, rank, lu->ownersIV,
518:                                          lu->options.msglvl, lu->options.msgFile) ;
519:   lu->nmycol = IV_size(lu->ownedColumnsIV) ;
520:   if ( lu->nmycol > 0) {
521:     DenseMtx_init(lu->mtxX, SPOOLES_REAL, 0, 0, lu->nmycol, 1, 1, lu->nmycol) ;
522:     /* get pointers rowindX and entX */
523:     DenseMtx_rowIndices(lu->mtxX, &lu->nmycol, &lu->rowindX);
524:     lu->entX = DenseMtx_entries(lu->mtxX) ;
525:   } else { /* lu->nmycol == 0 */
526:     lu->entX    = 0;
527:     lu->rowindX = 0;
528:   }

530:   if ( lu->scat ){
531:     VecDestroy(lu->vec_spooles);
532:     ISDestroy(lu->iden);
533:     ISDestroy(lu->is_petsc);
534:     VecScatterDestroy(lu->scat);
535:   }
536:   lu->scat = PETSC_NULL;
537: 
538:   lu->flg = SAME_NONZERO_PATTERN;
539:   return(0);
540: }

542: #endif