Actual source code: spooles.c

  1: /*$Id: spooles.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
  2: /* 
  3:    Provides an interface to the Spooles serial sparse solver
  4: */

 6:  #include src/mat/impls/aij/seq/aij.h
 7:  #include src/mat/impls/sbaij/seq/sbaij.h

  9: #if defined(PETSC_HAVE_SPOOLES) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_COMPLEX)
 10:  #include src/mat/impls/aij/seq/spooles.h

 12: extern int MatDestroy_SeqAIJ(Mat);
 13: int MatDestroy_SeqAIJ_Spooles(Mat A)
 14: {
 15:   Mat_Spooles *lu = (Mat_Spooles*)A->spptr;
 16:   int         ierr;
 17: 
 19: 
 20:   FrontMtx_free(lu->frontmtx) ;
 21:   IV_free(lu->newToOldIV) ;
 22:   IV_free(lu->oldToNewIV) ;
 23:   InpMtx_free(lu->mtxA) ;
 24:   ETree_free(lu->frontETree) ;
 25:   IVL_free(lu->symbfacIVL) ;
 26:   SubMtxManager_free(lu->mtxmanager) ;
 27: 
 28:   PetscFree(lu);
 29:   MatDestroy_SeqAIJ(A);

 31:   return(0);
 32: }

 34: int MatSolve_SeqAIJ_Spooles(Mat A,Vec b,Vec x)
 35: {
 36:   Mat_Spooles      *lu = (Mat_Spooles*)A->spptr;
 37:   PetscScalar      *array;
 38:   DenseMtx         *mtxY, *mtxX ;
 39:   double           *entX;
 40:   int              ierr,irow,neqns=A->m,*iv;

 43:   /* copy permuted b to mtxY */
 44:   mtxY = DenseMtx_new() ;
 45:   DenseMtx_init(mtxY, SPOOLES_REAL, 0, 0, neqns, 1, 1, neqns) ; /* column major */
 46:   iv = IV_entries(lu->oldToNewIV);
 47:   VecGetArray(b,&array);
 48:   for ( irow = 0 ; irow < neqns ; irow++ ) DenseMtx_setRealEntry(mtxY, *iv++, 0, *array++) ;
 49:   VecRestoreArray(b,&array);

 51:   mtxX = DenseMtx_new() ;
 52:   DenseMtx_init(mtxX, SPOOLES_REAL, 0, 0, neqns, 1, 1, neqns) ;
 53:   DenseMtx_zero(mtxX) ;
 54:   FrontMtx_solve(lu->frontmtx, mtxX, mtxY, lu->mtxmanager,
 55:                  lu->cpus, lu->options.msglvl, lu->options.msgFile) ;
 56:   if ( lu->options.msglvl > 2 ) {
 57:     fprintf(lu->options.msgFile, "nn right hand side matrix after permutation") ;
 58:     DenseMtx_writeForHumanEye(mtxY, lu->options.msgFile) ;
 59:     fprintf(lu->options.msgFile, "nn solution matrix in new ordering") ;
 60:     DenseMtx_writeForHumanEye(mtxX, lu->options.msgFile) ;
 61:     fflush(lu->options.msgFile) ;
 62:   }

 64:   /* permute solution into original ordering, then copy to x */
 65:   DenseMtx_permuteRows(mtxX, lu->newToOldIV);
 66:   VecGetArray(x,&array);
 67:   entX = DenseMtx_entries(mtxX);
 68:   DVcopy(neqns, array, entX);
 69:   VecRestoreArray(x,&array);
 70: 
 71:   /* free memory */
 72:   DenseMtx_free(mtxX) ;
 73:   DenseMtx_free(mtxY) ;
 74: 
 75:   return(0);
 76: }

 78: int MatFactorNumeric_SeqAIJ_Spooles(Mat A,Mat *F)
 79: {
 80:   Mat_Spooles        *lu = (Mat_Spooles*)(*F)->spptr;
 81:   ChvManager         *chvmanager ;
 82:   Chv                *rootchv ;
 83:   Graph              *graph ;
 84:   IVL                *adjIVL;
 85:   int                ierr,nz,m=A->m,irow,nedges,
 86:                      *ai,*aj,*ivec1, *ivec2, i;
 87:   PetscScalar        *av;
 88:   double             *dvec,cputotal;
 89: 
 91:   /* copy A to Spooles' InpMtx object */
 92:   if ( lu->options.symflag == SPOOLES_NONSYMMETRIC ) {
 93:     Mat_SeqAIJ   *mat = (Mat_SeqAIJ*)A->data;
 94:     ai=mat->i; aj=mat->j; av=mat->a;
 95:     nz=mat->nz;
 96:   } else {
 97:     Mat_SeqSBAIJ *mat = (Mat_SeqSBAIJ*)A->data;
 98:     ai=mat->i; aj=mat->j; av=mat->a;
 99:     nz=mat->s_nz;
100:   }
101:   if (lu->flg == DIFFERENT_NONZERO_PATTERN) lu->mtxA = InpMtx_new() ;
102:   InpMtx_init(lu->mtxA, INPMTX_BY_ROWS, SPOOLES_REAL, nz, m) ;
103:   ivec1 = InpMtx_ivec1(lu->mtxA);
104:   ivec2 = InpMtx_ivec2(lu->mtxA);
105:   dvec  = InpMtx_dvec(lu->mtxA);
106:   for (irow = 0; irow < m; irow++){
107:     for (i = ai[irow]; i<ai[irow+1]; i++) ivec1[i] = irow;
108:   }
109:   IVcopy(nz, ivec2, aj);
110:   DVcopy(nz, dvec, av);
111:   InpMtx_inputRealTriples(lu->mtxA, nz, ivec1, ivec2, dvec);
112:   InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;

114:   if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
115: 
116:     (*F)->ops->solve   = MatSolve_SeqAIJ_Spooles;
117:     (*F)->ops->destroy = MatDestroy_SeqAIJ_Spooles;
118:     (*F)->assembled    = PETSC_TRUE;
119: 
120:     SetSpoolesOptions(A, &lu->options);

122:     /*---------------------------------------------------
123:     find a low-fill ordering
124:          (1) create the Graph object
125:          (2) order the graph using multiple minimum degree
126:     -------------------------------------------------------*/
127:     graph = Graph_new() ;
128:     adjIVL = InpMtx_fullAdjacency(lu->mtxA) ;
129:     nedges = IVL_tsize(adjIVL) ;
130:     Graph_init2(graph, 0, m, 0, nedges, m, nedges, adjIVL,NULL, NULL) ;
131:     if ( lu->options.msglvl > 2 ) {
132:       fprintf(lu->options.msgFile, "nn graph of the input matrix") ;
133:       Graph_writeForHumanEye(graph, lu->options.msgFile) ;
134:       fflush(lu->options.msgFile) ;
135:     }

137:     switch (lu->options.ordering) {
138:     case 0:
139:       lu->frontETree = orderViaBestOfNDandMS(graph,
140:                      lu->options.maxdomainsize, lu->options.maxzeros, lu->options.maxsize,
141:                      lu->options.seed, lu->options.msglvl, lu->options.msgFile); break;
142:     case 1:
143:       lu->frontETree = orderViaMMD(graph,lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
144:     case 2:
145:       lu->frontETree = orderViaMS(graph, lu->options.maxdomainsize,
146:                      lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
147:     case 3:
148:       lu->frontETree = orderViaND(graph, lu->options.maxdomainsize,
149:                      lu->options.seed,lu->options.msglvl,lu->options.msgFile); break;
150:     default:
151:       SETERRQ(1,"Unknown Spooles's ordering");
152:     }
153:     Graph_free(graph) ;

155:     if ( lu->options.msglvl > 0 ) {
156:       fprintf(lu->options.msgFile, "nn front tree from ordering") ;
157:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
158:       fflush(lu->options.msgFile) ;
159:     }
160: 
161:     /* get the permutation, permute the front tree, permute the matrix */
162:     lu->oldToNewIV = ETree_oldToNewVtxPerm(lu->frontETree) ;
163:     lu->oldToNew   = IV_entries(lu->oldToNewIV) ;
164:     lu->newToOldIV = ETree_newToOldVtxPerm(lu->frontETree) ;
165:     ETree_permuteVertices(lu->frontETree, lu->oldToNewIV) ;

167:     InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew) ;
168:     if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
169:       InpMtx_mapToUpperTriangle(lu->mtxA) ;
170:     }
171:     InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
172:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;

174:     /* get symbolic factorization */
175:     lu->symbfacIVL = SymbFac_initFromInpMtx(lu->frontETree, lu->mtxA) ;

177:     if ( lu->options.msglvl > 2 ) {
178:       fprintf(lu->options.msgFile, "nn old-to-new permutation vector") ;
179:       IV_writeForHumanEye(lu->oldToNewIV, lu->options.msgFile) ;
180:       fprintf(lu->options.msgFile, "nn new-to-old permutation vector") ;
181:       IV_writeForHumanEye(lu->newToOldIV, lu->options.msgFile) ;
182:       fprintf(lu->options.msgFile, "nn front tree after permutation") ;
183:       ETree_writeForHumanEye(lu->frontETree, lu->options.msgFile) ;
184:       fprintf(lu->options.msgFile, "nn input matrix after permutation") ;
185:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
186:       fprintf(lu->options.msgFile, "nn symbolic factorization") ;
187:       IVL_writeForHumanEye(lu->symbfacIVL, lu->options.msgFile) ;
188:       fflush(lu->options.msgFile) ;
189:     }

191:     lu->frontmtx   = FrontMtx_new() ;
192:     lu->mtxmanager = SubMtxManager_new() ;
193:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;

195:   } else { /* new num factorization using previously computed symbolic factor */
196:     if (lu->options.pivotingflag) {              /* different FrontMtx is required */
197:       FrontMtx_free(lu->frontmtx) ;
198:       lu->frontmtx   = FrontMtx_new() ;
199:     }

201:     SubMtxManager_free(lu->mtxmanager) ;
202:     lu->mtxmanager = SubMtxManager_new() ;
203:     SubMtxManager_init(lu->mtxmanager, NO_LOCK, 0) ;

205:     /* permute mtxA */
206:     InpMtx_permute(lu->mtxA, lu->oldToNew, lu->oldToNew) ;
207:     if ( lu->options.symflag == SPOOLES_SYMMETRIC ) InpMtx_mapToUpperTriangle(lu->mtxA) ;
208:     InpMtx_changeCoordType(lu->mtxA, INPMTX_BY_CHEVRONS) ;
209:     InpMtx_changeStorageMode(lu->mtxA, INPMTX_BY_VECTORS) ;
210:     if ( lu->options.msglvl > 2 ) {
211:       fprintf(lu->options.msgFile, "nn input matrix after permutation") ;
212:       InpMtx_writeForHumanEye(lu->mtxA, lu->options.msgFile) ;
213:     }
214:   } /* end of if( lu->flg == DIFFERENT_NONZERO_PATTERN) */

216:   FrontMtx_init(lu->frontmtx, lu->frontETree, lu->symbfacIVL, SPOOLES_REAL, lu->options.symflag,
217:                 FRONTMTX_DENSE_FRONTS, lu->options.pivotingflag, NO_LOCK, 0, NULL,
218:                 lu->mtxmanager, lu->options.msglvl, lu->options.msgFile) ;
219: 

221:   if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
222:     if ( lu->options.patchAndGoFlag == 1 ) {
223:       lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
224:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 1, lu->options.toosmall, lu->options.fudge,
225:                        lu->options.storeids, lu->options.storevalues) ;
226:     } else if ( lu->options.patchAndGoFlag == 2 ) {
227:       lu->frontmtx->patchinfo = PatchAndGoInfo_new() ;
228:       PatchAndGoInfo_init(lu->frontmtx->patchinfo, 2, lu->options.toosmall, lu->options.fudge,
229:                        lu->options.storeids, lu->options.storevalues) ;
230:     }
231:   }

233:   /* numerical factorization */
234:   chvmanager = ChvManager_new() ;
235:   ChvManager_init(chvmanager, NO_LOCK, 1) ;
236:   DVfill(10, lu->cpus, 0.0) ;
237:   IVfill(20, lu->stats, 0) ;
238:   rootchv = FrontMtx_factorInpMtx(lu->frontmtx, lu->mtxA, lu->options.tau, 0.0,
239:             chvmanager, &ierr, lu->cpus, lu->stats, lu->options.msglvl, lu->options.msgFile) ;
240:   ChvManager_free(chvmanager) ;
241:   if(lu->options.FrontMtxInfo){
242:     PetscPrintf(PETSC_COMM_SELF,"n %8d pivots, %8d pivot tests, %8d delayed rows and columnsn",243:                lu->stats[0], lu->stats[1], lu->stats[2]);
244:     cputotal = lu->cpus[8] ;
245:     if ( cputotal > 0.0 ) {
246:       PetscPrintf(PETSC_COMM_SELF,
247:            "n                               cpus   cpus/totaltime"
248:            "n    initialize fronts       %8.3f %6.2f"
249:            "n    load original entries   %8.3f %6.2f"
250:            "n    update fronts           %8.3f %6.2f"
251:            "n    assemble postponed data %8.3f %6.2f"
252:            "n    factor fronts           %8.3f %6.2f"
253:            "n    extract postponed data  %8.3f %6.2f"
254:            "n    store factor entries    %8.3f %6.2f"
255:            "n    miscellaneous           %8.3f %6.2f"
256:            "n    total time              %8.3f n",
257:            lu->cpus[0], 100.*lu->cpus[0]/cputotal,
258:            lu->cpus[1], 100.*lu->cpus[1]/cputotal,
259:            lu->cpus[2], 100.*lu->cpus[2]/cputotal,
260:            lu->cpus[3], 100.*lu->cpus[3]/cputotal,
261:            lu->cpus[4], 100.*lu->cpus[4]/cputotal,
262:            lu->cpus[5], 100.*lu->cpus[5]/cputotal,
263:            lu->cpus[6], 100.*lu->cpus[6]/cputotal,
264:            lu->cpus[7], 100.*lu->cpus[7]/cputotal, cputotal) ;
265:     }
266:   }
267:   if ( lu->options.msglvl > 0 ) {
268:     fprintf(lu->options.msgFile, "nn factor matrix") ;
269:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
270:     fflush(lu->options.msgFile) ;
271:   }

273:   if ( lu->options.symflag == SPOOLES_SYMMETRIC ) {
274:     if ( lu->options.patchAndGoFlag == 1 ) {
275:       if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
276:         if (lu->options.msglvl > 0 ){
277:           fprintf(lu->options.msgFile, "n small pivots found at these locations") ;
278:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
279:         }
280:       }
281:       PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
282:     } else if ( lu->options.patchAndGoFlag == 2 ) {
283:       if (lu->options.msglvl > 0 ){
284:         if ( lu->frontmtx->patchinfo->fudgeIV != NULL ) {
285:           fprintf(lu->options.msgFile, "n small pivots found at these locations") ;
286:           IV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeIV, lu->options.msgFile) ;
287:         }
288:         if ( lu->frontmtx->patchinfo->fudgeDV != NULL ) {
289:           fprintf(lu->options.msgFile, "n perturbations") ;
290:           DV_writeForHumanEye(lu->frontmtx->patchinfo->fudgeDV, lu->options.msgFile) ;
291:         }
292:       }
293:       PatchAndGoInfo_free(lu->frontmtx->patchinfo) ;
294:     }
295:   }

297:   if ( rootchv != NULL ) SETERRQ(1,"n matrix found to be singular");
298:   if ( ierr >= 0 ) SETERRQ1(1,"n error encountered at front %d", ierr);

300:   /* post-process the factorization */
301:   FrontMtx_postProcess(lu->frontmtx, lu->options.msglvl, lu->options.msgFile) ;
302:   if ( lu->options.msglvl > 2 ) {
303:     fprintf(lu->options.msgFile, "nn factor matrix after post-processing") ;
304:     FrontMtx_writeForHumanEye(lu->frontmtx, lu->options.msgFile) ;
305:     fflush(lu->options.msgFile) ;
306:   }

308:   lu->flg         = SAME_NONZERO_PATTERN;
309: 
310:   return(0);
311: }

313: #endif