Actual source code: spooles.c

  1: #define PETCSMAT_DLL

  3: /* 
  4:    Provides an interface to the Spooles serial sparse solver
  5: */
 6:  #include src/mat/impls/aij/seq/aij.h
 7:  #include src/mat/impls/sbaij/seq/sbaij.h
 8:  #include src/mat/impls/aij/seq/spooles/spooles.h

 10: /* make sun CC happy */
 11: static void (*f)(void);

 16: PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_Spooles_Base(Mat A,MatType type,MatReuse reuse,Mat *newmat) {
 17:   /* This routine is only called to convert an unfactored PETSc-Spooles matrix */
 18:   /* to its base PETSc type, so we will ignore 'MatType type'. */
 20:   Mat            B=*newmat;
 21:   Mat_Spooles    *lu=(Mat_Spooles*)A->spptr;

 24:   if (reuse == MAT_INITIAL_MATRIX) {
 25:     MatDuplicate(A,MAT_COPY_VALUES,&B);
 26:   }
 27:   /* Reset the stashed function pointers set by inherited routines */
 28:   B->ops->duplicate              = lu->MatDuplicate;
 29:   B->ops->choleskyfactorsymbolic = lu->MatCholeskyFactorSymbolic;
 30:   B->ops->lufactorsymbolic       = lu->MatLUFactorSymbolic;
 31:   B->ops->view                   = lu->MatView;
 32:   B->ops->assemblyend            = lu->MatAssemblyEnd;
 33:   B->ops->destroy                = lu->MatDestroy;

 35:   PetscObjectQueryFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C",&f);
 36:   if (f) {
 37:     PetscObjectComposeFunction((PetscObject)B,"MatMPISBAIJSetPreallocation_C","",(PetscVoidFunction)lu->MatPreallocate);
 38:   }
 39:   PetscFree(lu);

 41:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaijspooles_seqaij_C","",PETSC_NULL);
 42:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_seqaijspooles_C","",PETSC_NULL);
 43:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaijspooles_mpiaij_C","",PETSC_NULL);
 44:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpiaij_mpiaijspooles_C","",PETSC_NULL);
 45:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaijspooles_seqsbaij_C","",PETSC_NULL);
 46:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqsbaij_seqsbaijspooles_C","",PETSC_NULL);
 47:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaijspooles_mpisbaij_C","",PETSC_NULL);
 48:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_mpisbaij_mpisbaijspooles_C","",PETSC_NULL);

 50:   PetscObjectChangeTypeName((PetscObject)B,type);
 51:   *newmat = B;
 52:   return(0);
 53: }

 58: PetscErrorCode MatDestroy_SeqAIJSpooles(Mat A)
 59: {
 60:   Mat_Spooles    *lu = (Mat_Spooles*)A->spptr;
 62: 
 64:   if (lu->CleanUpSpooles) {
 65:     FrontMtx_free(lu->frontmtx);
 66:     IV_free(lu->newToOldIV);
 67:     IV_free(lu->oldToNewIV);
 68:     InpMtx_free(lu->mtxA);
 69:     ETree_free(lu->frontETree);
 70:     IVL_free(lu->symbfacIVL);
 71:     SubMtxManager_free(lu->mtxmanager);
 72:     Graph_free(lu->graph);
 73:   }
 74:   MatConvert_Spooles_Base(A,lu->basetype,MAT_REUSE_MATRIX,&A);
 75:   (*A->ops->destroy)(A);
 76:   return(0);
 77: }

 81: PetscErrorCode MatSolve_SeqAIJSpooles(Mat A,Vec b,Vec x)
 82: {
 83:   Mat_Spooles      *lu = (Mat_Spooles*)A->spptr;
 84:   PetscScalar      *array;
 85:   DenseMtx         *mtxY, *mtxX ;
 86:   PetscErrorCode   ierr;
 87:   PetscInt         irow,neqns=A->cmap.n,nrow=A->rmap.n,*iv;
 88: #if defined(PETSC_USE_COMPLEX)
 89:   double           x_real,x_imag;
 90: #else
 91:   double           *entX;
 92: #endif

 95:   mtxY = DenseMtx_new();
 96:   DenseMtx_init(mtxY, lu->options.typeflag, 0, 0, nrow, 1, 1, nrow); /* column major */
 97:   VecGetArray(b,&array);

 99:   if (lu->options.useQR) {   /* copy b to mtxY */
100:     for ( irow = 0 ; irow < nrow; irow++ )
101: #if !defined(PETSC_USE_COMPLEX)
102:       DenseMtx_setRealEntry(mtxY, irow, 0, *array++);
103: #else
104:       DenseMtx_setComplexEntry(mtxY, irow, 0, PetscRealPart(array[irow]), PetscImaginaryPart(array[irow]));
105: #endif
106:   } else {                   /* copy permuted b to mtxY */
107:     iv = IV_entries(lu->oldToNewIV);
108:     for ( irow = 0 ; irow < nrow; irow++ )
109: #if !defined(PETSC_USE_COMPLEX)
110:       DenseMtx_setRealEntry(mtxY, *iv++, 0, *array++);
111: #else
112:       DenseMtx_setComplexEntry(mtxY,*iv++,0,PetscRealPart(array[irow]),PetscImaginaryPart(array[irow]));
113: #endif
114:   }
115:   VecRestoreArray(b,&array);

117:   mtxX = DenseMtx_new();
118:   DenseMtx_init(mtxX, lu->options.typeflag, 0, 0, neqns, 1, 1, neqns);
119:   if (lu->options.useQR) {
120:     FrontMtx_QR_solve(lu->frontmtx, lu->mtxA, mtxX, mtxY, lu->mtxmanager,
121:                   lu->cpus, lu->options.msglvl, lu->options.msgFile);
122:   } else {
123:     FrontMtx_solve(lu->frontmtx, mtxX, mtxY, lu->mtxmanager,
124:                  lu->cpus, lu->options.msglvl, lu->options.msgFile);
125:   }
126:   if ( lu->options.msglvl > 2 ) {
127:     PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n right hand side matrix after permutation");
128:     DenseMtx_writeForHumanEye(mtxY, lu->options.msgFile);
129:     PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n solution matrix in new ordering");
130:     DenseMtx_writeForHumanEye(mtxX, lu->options.msgFile);
131:     fflush(lu->options.msgFile);
132:   }

134:   /* permute solution into original ordering, then copy to x */
135:   DenseMtx_permuteRows(mtxX, lu->newToOldIV);
136:   VecGetArray(x,&array);

138: #if !defined(PETSC_USE_COMPLEX)
139:   entX = DenseMtx_entries(mtxX);
140:   DVcopy(neqns, array, entX);
141: #else
142:   for (irow=0; irow<nrow; irow++){
143:     DenseMtx_complexEntry(mtxX,irow,0,&x_real,&x_imag);
144:     array[irow] = x_real+x_imag*PETSC_i;
145:   }
146: #endif

148:   VecRestoreArray(x,&array);
149: 
150:   /* free memory */
151:   DenseMtx_free(mtxX);
152:   DenseMtx_free(mtxY);
153:   return(0);
154: }

158: PetscErrorCode MatFactorNumeric_SeqAIJSpooles(Mat A,MatFactorInfo *info,Mat *F)
159: {
160:   Mat_Spooles        *lu = (Mat_Spooles*)(*F)->spptr;
161:   ChvManager         *chvmanager ;
162:   Chv                *rootchv ;
163:   IVL                *adjIVL;
164:   PetscErrorCode     ierr;
165:   PetscInt           nz,nrow=A->rmap.n,irow,nedges,neqns=A->cmap.n,*ai,*aj,i,*diag=0,fierr;
166:   PetscScalar        *av;
167:   double             cputotal,facops;
168: #if defined(PETSC_USE_COMPLEX)
169:   PetscInt           nz_row,*aj_tmp;
170:   PetscScalar        *av_tmp;
171: #else
172:   PetscInt           *ivec1,*ivec2,j;
173:   double             *dvec;
174: #endif
175:   PetscTruth         isAIJ,isSeqAIJ;
176: 
178:   if (lu->flg == DIFFERENT_NONZERO_PATTERN) { /* first numeric factorization */
179:     (*F)->ops->solve   = MatSolve_SeqAIJSpooles;
180:     (*F)->ops->destroy = MatDestroy_SeqAIJSpooles;
181:     (*F)->assembled    = PETSC_TRUE;
182: 
183:     /* set Spooles options */
184:     SetSpoolesOptions(A, &lu->options);

186:     lu->mtxA = InpMtx_new();
187:   }

189:   /* copy A to Spooles' InpMtx object */
190:   PetscTypeCompare((PetscObject)A,MATSEQAIJSPOOLES,&isSeqAIJ);
191:   PetscTypeCompare((PetscObject)A,MATAIJSPOOLES,&isAIJ);
192:   if (isSeqAIJ || isAIJ){
193:     Mat_SeqAIJ   *mat = (Mat_SeqAIJ*)A->data;
194:     ai=mat->i; aj=mat->j; av=mat->a;
195:     if (lu->options.symflag == SPOOLES_NONSYMMETRIC) {
196:       nz=mat->nz;
197:     } else { /* SPOOLES_SYMMETRIC || SPOOLES_HERMITIAN */
198:       nz=(mat->nz + A->rmap.n)/2;
199:       if (!mat->diag){
200:         MatMarkDiagonal_SeqAIJ(A);
201:       }
202:       diag=mat->diag;
203:     }
204:   } else { /* A is SBAIJ */
205:       Mat_SeqSBAIJ *mat = (Mat_SeqSBAIJ*)A->data;
206:       ai=mat->i; aj=mat->j; av=mat->a;
207:       nz=mat->nz;
208:   }
209:   InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, lu->options.typeflag, nz, 0);
210: 
211: #if defined(PETSC_USE_COMPLEX)
212:     for (irow=0; irow<nrow; irow++) {
213:       if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
214:         nz_row = ai[irow+1] - ai[irow];
215:         aj_tmp = aj + ai[irow];
216:         av_tmp = av + ai[irow];
217:       } else {
218:         nz_row = ai[irow+1] - diag[irow];
219:         aj_tmp = aj + diag[irow];
220:         av_tmp = av + diag[irow];
221:       }
222:       for (i=0; i<nz_row; i++){
223:         InpMtx_inputComplexEntry(lu->mtxA, irow, *aj_tmp++,PetscRealPart(*av_tmp),PetscImaginaryPart(*av_tmp));
224:         av_tmp++;
225:       }
226:     }
227: #else
228:     ivec1 = InpMtx_ivec1(lu->mtxA);
229:     ivec2 = InpMtx_ivec2(lu->mtxA);
230:     dvec  = InpMtx_dvec(lu->mtxA);
231:     if ( lu->options.symflag == SPOOLES_NONSYMMETRIC || !isAIJ){
232:       for (irow = 0; irow < nrow; irow++){
233:         for (i = ai[irow]; i<ai[irow+1]; i++) ivec1[i] = irow;
234:       }
235:       IVcopy(nz, ivec2, aj);
236:       DVcopy(nz, dvec, av);
237:     } else {
238:       nz = 0;
239:       for (irow = 0; irow < nrow; irow++){
240:         for (j = diag[irow]; j<ai[irow+1]; j++) {
241:           ivec1[nz] = irow;
242:           ivec2[nz] = aj[j];
243:           dvec[nz]  = av[j];
244:           nz++;
245:         }
246:       }
247:     }
248:     InpMtx_inputRealTriples(lu->mtxA, nz, ivec1, ivec2, dvec);
249: #endif

251:   InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
252:   if ( lu->options.msglvl > 0 ) {
253:     printf("\n\n input matrix");
254:     PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix");
255:     InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
256:     fflush(lu->options.msgFile);
257:   }

259:   if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
260:     /*---------------------------------------------------
261:     find a low-fill ordering
262:          (1) create the Graph object
263:          (2) order the graph 
264:     -------------------------------------------------------*/
265:     if (lu->options.useQR){
266:       adjIVL = InpMtx_adjForATA(lu->mtxA);
267:     } else {
268:       adjIVL = InpMtx_fullAdjacency(lu->mtxA);
269:     }
270:     nedges = IVL_tsize(adjIVL);

272:     lu->graph = Graph_new();
273:     Graph_init2(lu->graph, 0, neqns, 0, nedges, neqns, nedges, adjIVL, NULL, NULL);
274:     if ( lu->options.msglvl > 2 ) {
275:       if (lu->options.useQR){
276:         PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of A^T A");
277:       } else {
278:         PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n graph of the input matrix");
279:       }
280:       Graph_writeForHumanEye(lu->graph, lu->options.msgFile);
281:       fflush(lu->options.msgFile);
282:     }

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

301:     if ( lu->options.msglvl > 0 ) {
302:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree from ordering");
303:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile);
304:       fflush(lu->options.msgFile);
305:     }
306: 
307:     /* get the permutation, permute the front tree */
308:     lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree);
309:     lu->oldToNew   = IV_entries(lu->oldToNewIV);
310:     lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree);
311:     if (!lu->options.useQR) ETree_permuteVertices(lu->frontETree, lu->oldToNewIV);

313:     /* permute the matrix */
314:     if (lu->options.useQR){
315:       InpMtx_permute(lu->mtxA, NULL, lu->oldToNew);
316:     } else {
317:       InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew);
318:       if ( lu->options.symflag == SPOOLES_SYMMETRIC) {
319:         InpMtx_mapToUpperTriangle(lu->mtxA);
320:       }
321: #if defined(PETSC_USE_COMPLEX)
322:       if ( lu->options.symflag == SPOOLES_HERMITIAN ) {
323:         InpMtx_mapToUpperTriangleH(lu->mtxA);
324:       }
325: #endif
326:       InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS);
327:     }
328:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);

330:     /* get symbolic factorization */
331:     if (lu->options.useQR){
332:       lu->symbfacIVL = SymbFac_initFromGraph(lu->frontETree, lu->graph);
333:       IVL_overwrite(lu->symbfacIVL, lu->oldToNewIV);
334:       IVL_sortUp(lu->symbfacIVL);
335:       ETree_permuteVertices(lu->frontETree, lu->oldToNewIV);
336:     } else {
337:       lu->symbfacIVL = SymbFac_initFromInpMtx(lu->frontETree, lu->mtxA);
338:     }
339:     if ( lu->options.msglvl > 2 ) {
340:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n old-to-new permutation vector");
341:       IV_writeForHumanEye(lu->oldToNewIV, lu->options.msgFile);
342:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n new-to-old permutation vector");
343:       IV_writeForHumanEye(lu->newToOldIV, lu->options.msgFile);
344:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n front tree after permutation");
345:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile);
346:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");
347:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
348:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n symbolic factorization");
349:       IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile);
350:       fflush(lu->options.msgFile);
351:     }

353:     lu->frontmtx   = FrontMtx_new();
354:     lu->mtxmanager = SubMtxManager_new();
355:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0);

357:   } else { /* new num factorization using previously computed symbolic factor */

359:     if (lu->options.pivotingflag) { /* different FrontMtx is required */
360:       FrontMtx_free(lu->frontmtx);
361:       lu->frontmtx   = FrontMtx_new();
362:     } else {
363:       FrontMtx_clearData (lu->frontmtx);
364:     }

366:     SubMtxManager_free(lu->mtxmanager);
367:     lu->mtxmanager = SubMtxManager_new();
368:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0);

370:     /* permute mtxA */
371:     if (lu->options.useQR){
372:       InpMtx_permute(lu->mtxA, NULL, lu->oldToNew);
373:     } else {
374:       InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew);
375:       if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
376:         InpMtx_mapToUpperTriangle(lu->mtxA);
377:       }
378:       InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS);
379:     }
380:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS);
381:     if ( lu->options.msglvl > 2 ) {
382:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n input matrix after permutation");
383:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile);
384:     }
385:   } /* end of if( lu->flg == DIFFERENT_NONZERO_PATTERN) */
386: 
387:   if (lu->options.useQR){
388:     FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag,
389:                  SPOOLES_SYMMETRIC, FRONTMTX_DENSE_FRONTS,
390:                  SPOOLES_NO_PIVOTING, NO_LOCK, 0, NULL,
391:                  lu->mtxmanager, lu->options.msglvl, lu->options.msgFile);
392:   } else {
393:     FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, lu->options.typeflag, lu->options.symflag,
394:                 FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, 0, NULL,
395:                 lu->mtxmanager, lu->options.msglvl, lu->options.msgFile);
396:   }

398:   if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {  /* || SPOOLES_HERMITIAN ? */
399:     if ( lu->options.patchAndGoFlag == 1 ) {
400:       lu->frontmtx->patchinfo = PatchAndGoInfo_new();
401:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
402:                        lu->options.storeids, lu->options.storevalues);
403:     } else if ( lu->options.patchAndGoFlag == 2 ) {
404:       lu->frontmtx->patchinfo = PatchAndGoInfo_new();
405:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
406:                        lu->options.storeids, lu->options.storevalues);
407:     }
408:   }

410:   /* numerical factorization */
411:   chvmanager = ChvManager_new();
412:   ChvManager_init(chvmanager, NO_LOCK, 1);
413:   DVfill(10, lu->cpus, 0.0);
414:   if (lu->options.useQR){
415:     facops = 0.0 ;
416:     FrontMtx_QR_factor(lu->frontmtx, lu->mtxA, chvmanager,
417:                    lu->cpus, &facops, lu->options.msglvl, lu->options.msgFile);
418:     if ( lu->options.msglvl > 1 ) {
419:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");
420:       PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n facops = %9.2f", facops);
421:     }
422:   } else {
423:     IVfill(20, lu->stats, 0);
424:     rootchv = FrontMtx_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, 0.0,
425:             chvmanager, &fierr, lu->cpus,lu->stats,lu->options.msglvl,lu->options.msgFile);
426:     if (rootchv) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"\n matrix found to be singular");
427:     if (fierr >= 0) SETERRQ1(PETSC_ERR_LIB,"\n error encountered at front %D", fierr);
428: 
429:     if(lu->options.FrontMtxInfo){
430:       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]);
431:       cputotal = lu->cpus[8] ;
432:       if ( cputotal > 0.0 ) {
433:         PetscPrintf(PETSC_COMM_SELF,
434:            "\n                               cpus   cpus/totaltime"
435:            "\n    initialize fronts       %8.3f %6.2f"
436:            "\n    load original entries   %8.3f %6.2f"
437:            "\n    update fronts           %8.3f %6.2f"
438:            "\n    assemble postponed data %8.3f %6.2f"
439:            "\n    factor fronts           %8.3f %6.2f"
440:            "\n    extract postponed data  %8.3f %6.2f"
441:            "\n    store factor entries    %8.3f %6.2f"
442:            "\n    miscellaneous           %8.3f %6.2f"
443:            "\n    total time              %8.3f \n",
444:            lu->cpus[0], 100.*lu->cpus[0]/cputotal,
445:            lu->cpus[1], 100.*lu->cpus[1]/cputotal,
446:            lu->cpus[2], 100.*lu->cpus[2]/cputotal,
447:            lu->cpus[3], 100.*lu->cpus[3]/cputotal,
448:            lu->cpus[4], 100.*lu->cpus[4]/cputotal,
449:            lu->cpus[5], 100.*lu->cpus[5]/cputotal,
450:            lu->cpus[6], 100.*lu->cpus[6]/cputotal,
451:            lu->cpus[7], 100.*lu->cpus[7]/cputotal, cputotal);
452:       }
453:     }
454:   }
455:   ChvManager_free(chvmanager);

457:   if ( lu->options.msglvl > 0 ) {
458:     PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix");
459:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile);
460:     fflush(lu->options.msgFile);
461:   }

463:   if ( lu->options.symflag == SPOOLES_SYMMETRIC ) { /* || SPOOLES_HERMITIAN ? */
464:     if ( lu->options.patchAndGoFlag == 1 ) {
465:       if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
466:         if (lu->options.msglvl > 0 ){
467:           PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");
468:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile);
469:         }
470:       }
471:       PatchAndGoInfo_free(lu->frontmtx->patchinfo);
472:     } else if ( lu->options.patchAndGoFlag == 2 ) {
473:       if (lu->options.msglvl > 0 ){
474:         if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
475:           PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n small pivots found at these locations");
476:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile);
477:         }
478:         if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
479:           PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n perturbations");
480:           DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile);
481:         }
482:       }
483:       PatchAndGoInfo_free(lu->frontmtx->patchinfo);
484:     }
485:   }

487:   /* post-process the factorization */
488:   FrontMtx_postProcess(lu->frontmtx, lu->options.msglvl, lu->options.msgFile);
489:   if ( lu->options.msglvl > 2 ) {
490:     PetscFPrintf(PETSC_COMM_SELF,lu->options.msgFile, "\n\n factor matrix after post-processing");
491:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile);
492:     fflush(lu->options.msgFile);
493:   }

495:   lu->flg = SAME_NONZERO_PATTERN;
496:   lu->CleanUpSpooles = PETSC_TRUE;
497:   return(0);
498: }

503: PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqAIJSpooles(Mat A,MatType type,MatReuse reuse,Mat *newmat) {
504:   /* This routine is only called to convert a MATSEQAIJ matrix */
505:   /* to a MATSEQAIJSPOOLES matrix, so we will ignore 'MatType type'. */
507:   Mat            B=*newmat;
508:   Mat_Spooles    *lu;

511:   if (reuse == MAT_INITIAL_MATRIX) {
512:     /* This routine is inherited, so we know the type is correct. */
513:     MatDuplicate(A,MAT_COPY_VALUES,&B);
514:   }
515:   PetscNew(Mat_Spooles,&lu);
516:   B->spptr = (void*)lu;

518:   lu->basetype                   = MATSEQAIJ;
519:   lu->useQR                      = PETSC_FALSE;
520:   lu->CleanUpSpooles             = PETSC_FALSE;
521:   lu->MatDuplicate               = A->ops->duplicate;
522:   lu->MatCholeskyFactorSymbolic  = A->ops->choleskyfactorsymbolic;
523:   lu->MatLUFactorSymbolic        = A->ops->lufactorsymbolic;
524:   lu->MatView                    = A->ops->view;
525:   lu->MatAssemblyEnd             = A->ops->assemblyend;
526:   lu->MatDestroy                 = A->ops->destroy;
527:   B->ops->duplicate              = MatDuplicate_Spooles;
528:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJSpooles;
529:   B->ops->lufactorsymbolic       = MatLUFactorSymbolic_SeqAIJSpooles;
530:   B->ops->view                   = MatView_SeqAIJSpooles;
531:   B->ops->assemblyend            = MatAssemblyEnd_SeqAIJSpooles;
532:   B->ops->destroy                = MatDestroy_SeqAIJSpooles;

534:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaijspooles_seqaij_C",
535:                                            "MatConvert_Spooles_Base",MatConvert_Spooles_Base);
536:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqaijspooles_C",
537:                                            "MatConvert_SeqAIJ_SeqAIJSpooles",MatConvert_SeqAIJ_SeqAIJSpooles);
538:   /* PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJSPOOLES); */
539:   PetscObjectChangeTypeName((PetscObject)B,type);
540:   *newmat = B;
541:   return(0);
542: }

547: PetscErrorCode MatDuplicate_Spooles(Mat A, MatDuplicateOption op, Mat *M) {
549:   Mat_Spooles    *lu=(Mat_Spooles *)A->spptr;

552:   (*lu->MatDuplicate)(A,op,M);
553:   PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Spooles));
554:   return(0);
555: }

557: /*MC
558:   MATSEQAIJSPOOLES - MATSEQAIJSPOOLES = "seqaijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential matrices 
559:   via the external package SPOOLES.

561:   If SPOOLES is installed (see the manual for
562:   instructions on how to declare the existence of external packages),
563:   a matrix type can be constructed which invokes SPOOLES solvers.
564:   After calling MatCreate(...,A), simply call MatSetType(A,MATSEQAIJSPOOLES).

566:   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is 
567:   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from 
568:   the MATSEQAIJ type without data copy.

570:   Options Database Keys:
571: + -mat_type seqaijspooles - sets the matrix type to "seqaijspooles" during a call to MatSetFromOptions()
572: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
573: . -mat_spooles_seed <seed> - random number seed used for ordering
574: . -mat_spooles_msglvl <msglvl> - message output level
575: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
576: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
577: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
578: . -mat_spooles_maxsize <n> - maximum size of a supernode
579: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
580: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
581: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
582: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
583: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
584: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
585: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object

587:    Level: beginner

589: .seealso: PCLU
590: M*/

595: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqAIJSpooles(Mat A)
596: {

600:   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and SeqAIJSpooles types */
601:   PetscObjectChangeTypeName((PetscObject)A,MATSEQAIJSPOOLES);
602:   MatSetType(A,MATSEQAIJ);
603:   MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,MAT_REUSE_MATRIX,&A);
604:   return(0);
605: }

608: /*MC
609:   MATAIJSPOOLES - MATAIJSPOOLES = "aijspooles" - A matrix type providing direct solvers (LU or Cholesky) for sequential and parellel matrices 
610:   via the external package SPOOLES.

612:   If SPOOLES is installed (see the manual for
613:   instructions on how to declare the existence of external packages),
614:   a matrix type can be constructed which invokes SPOOLES solvers.
615:   After calling MatCreate(...,A), simply call MatSetType(A,MATAIJSPOOLES).
616:   This matrix type is supported for double precision real and complex.

618:   This matrix inherits from MATAIJ.  As a result, MatSeqAIJSetPreallocation and MatMPIAIJSetPreallocation are
619:   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from 
620:   the MATAIJ type without data copy.

622:   Options Database Keys:
623: + -mat_type aijspooles - sets the matrix type to "aijspooles" during a call to MatSetFromOptions()
624: . -mat_spooles_tau <tau> - upper bound on the magnitude of the largest element in L or U
625: . -mat_spooles_seed <seed> - random number seed used for ordering
626: . -mat_spooles_msglvl <msglvl> - message output level
627: . -mat_spooles_ordering <BestOfNDandMS,MMD,MS,ND> - ordering used
628: . -mat_spooles_maxdomainsize <n> - maximum subgraph size used by Spooles orderings
629: . -mat_spooles_maxzeros <n> - maximum number of zeros inside a supernode
630: . -mat_spooles_maxsize <n> - maximum size of a supernode
631: . -mat_spooles_FrontMtxInfo <true,fase> - print Spooles information about the computed factorization
632: . -mat_spooles_symmetryflag <0,1,2> - 0: SPOOLES_SYMMETRIC, 1: SPOOLES_HERMITIAN, 2: SPOOLES_NONSYMMETRIC
633: . -mat_spooles_patchAndGoFlag <0,1,2> - 0: no patch, 1: use PatchAndGo strategy 1, 2: use PatchAndGo strategy 2
634: . -mat_spooles_toosmall <dt> - drop tolerance for PatchAndGo strategy 1
635: . -mat_spooles_storeids <bool integer> - if nonzero, stores row and col numbers where patches were applied in an IV object
636: . -mat_spooles_fudge <delta> - fudge factor for rescaling diagonals with PatchAndGo strategy 2
637: - -mat_spooles_storevalues <bool integer> - if nonzero and PatchAndGo strategy 2 is used, store change in diagonal value in a DV object

639:    Level: beginner

641: .seealso: PCLU
642: M*/
646: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJSpooles(Mat A)
647: {
649:   PetscMPIInt    size;

652:   /* Change type name before calling MatSetType to force proper construction of SeqAIJSpooles or MPIAIJSpooles */
653:   PetscObjectChangeTypeName((PetscObject)A,MATAIJSPOOLES);
654:   MPI_Comm_size(A->comm,&size);
655:   if (size == 1) {
656:     MatSetType(A,MATSEQAIJ);
657:     MatConvert_SeqAIJ_SeqAIJSpooles(A,MATSEQAIJSPOOLES,MAT_REUSE_MATRIX,&A);
658:   } else {
659:     MatSetType(A,MATMPIAIJ);
660:     MatConvert_MPIAIJ_MPIAIJSpooles(A,MATMPIAIJSPOOLES,MAT_REUSE_MATRIX,&A);
661:   }
662:   return(0);
663: }