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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

586:    Level: beginner

588: .seealso: PCLU
589: M*/

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

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

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

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

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

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

638:    Level: beginner

640: .seealso: PCLU
641: M*/
645: PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_AIJSpooles(Mat A)
646: {
648:   PetscMPIInt    size;

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