Actual source code: dscpack.c

  1: #define PETSCMAT_DLL

  3: /* 
  4:         Provides an interface to the DSCPACK (Domain-Separator Codes) sparse direct solver
  5: */

 7:  #include ../src/mat/impls/baij/seq/baij.h
 8:  #include ../src/mat/impls/baij/mpi/mpibaij.h

 11: #include "dscmain.h"

 14: typedef struct {
 15:   DSC_Solver        My_DSC_Solver;
 16:   PetscInt      num_local_strucs, *local_struc_old_num,
 17:                 num_local_cols, num_local_nonz,
 18:                 *global_struc_new_col_num,
 19:                 *global_struc_new_num, *global_struc_owner,
 20:                 dsc_id,bs,*local_cols_old_num,*replication;
 21:   PetscInt      order_code,scheme_code,factor_type, stat,
 22:                 LBLASLevel,DBLASLevel,max_mem_allowed;
 23:   MatStructure  flg;
 24:   IS            my_cols,iden,iden_dsc;
 25:   Vec           vec_dsc;
 26:   VecScatter    scat;
 27:   MPI_Comm      comm_dsc;

 29:   /* A few inheritance details */
 30:   PetscMPIInt    size;

 32:   PetscTruth CleanUpDSCPACK;
 33: } Mat_DSCPACK;


 36: /* DSC function */
 39: void isort2(PetscInt size, PetscInt *list, PetscInt *idx_dsc) {
 40:   /* in increasing order */
 41:   /* idx_dsc will contain indices such that */
 42:   /* list can be accessed in sorted order */
 43:   PetscInt i, j, x, y;
 44: 
 45:   for (i=0; i<size; i++) idx_dsc[i] =i;

 47:   for (i=1; i<size; i++){
 48:     y= idx_dsc[i];
 49:     x=list[idx_dsc[i]];
 50:     for (j=i-1; ((j>=0) && (x<list[idx_dsc[j]])); j--)
 51:       idx_dsc[j+1]=idx_dsc[j];
 52:     idx_dsc[j+1]=y;
 53:   }
 54: }/*end isort2*/

 58: PetscErrorCode  BAIJtoMyANonz( PetscInt *AIndex, PetscInt *AStruct, PetscInt bs,
 59:                     RealNumberType *ANonz, PetscInt NumLocalStructs, 
 60:                     PetscInt NumLocalNonz,  PetscInt *GlobalStructNewColNum,                
 61:                     PetscInt *LocalStructOldNum,
 62:                     PetscInt *LocalStructLocalNum,
 63:                     RealNumberType **adr_MyANonz)
 64: /* 
 65:    Extract non-zero values of lower triangular part
 66:    of the permuted matrix that belong to this processor.

 68:    Only output parameter is adr_MyANonz -- is malloced and changed.
 69:    Rest are input parameters left unchanged.

 71:    When LocalStructLocalNum == PETSC_NULL,
 72:         AIndex, AStruct, and ANonz contain entire original matrix A 
 73:         in PETSc SeqBAIJ format,
 74:         otherwise,
 75:         AIndex, AStruct, and ANonz are indeces for the submatrix
 76:         of A whose colomns (in increasing order) belong to this processor.

 78:    Other variables supply information on ownership of columns
 79:    and the new numbering in a fill-reducing permutation

 81:    This information is used to setup lower half of A nonzeroes
 82:    for columns owned by this processor
 83:  */
 84: {
 86:   PetscInt            i, j, k, iold,inew, jj, kk, bs2=bs*bs,
 87:                  *idx, *NewColNum,
 88:                  MyANonz_last, max_struct=0, struct_size;
 89:   RealNumberType *MyANonz;


 93:   /* loop: to find maximum number of subscripts over columns
 94:      assigned to this processor */
 95:   for (i=0; i <NumLocalStructs; i++) {
 96:     /* for each struct i (local) assigned to this processor */
 97:     if (LocalStructLocalNum){
 98:       iold = LocalStructLocalNum[i];
 99:     } else {
100:       iold = LocalStructOldNum[i];
101:     }
102: 
103:     struct_size = AIndex[iold+1] - AIndex[iold];
104:     if ( max_struct <= struct_size) max_struct = struct_size;
105:   }

107:   /* allocate tmp arrays large enough to hold densest struct */
108:   PetscMalloc((2*max_struct+1)*sizeof(PetscInt),&NewColNum);
109:   idx = NewColNum + max_struct;
110: 
111:   PetscMalloc(NumLocalNonz*sizeof(RealNumberType),&MyANonz);
112:   *adr_MyANonz = MyANonz;

114:   /* loop to set up nonzeroes in MyANonz */
115:   MyANonz_last = 0 ; /* points to first empty space in MyANonz */
116:   for (i=0; i <NumLocalStructs; i++) {

118:     /* for each struct i (local) assigned to this processor */
119:     if (LocalStructLocalNum){
120:       iold = LocalStructLocalNum[i];
121:     } else {
122:       iold = LocalStructOldNum[i];
123:     }

125:     struct_size = AIndex[iold+1] - AIndex[iold];
126:     for (k=0, j=AIndex[iold]; j<AIndex[iold+1]; j++){
127:       NewColNum[k] = GlobalStructNewColNum[AStruct[j]];
128:       k++;
129:     }
130:     isort2(struct_size, NewColNum, idx);
131: 
132:     kk = AIndex[iold]*bs2; /* points to 1st element of iold block col in ANonz */
133:     inew = GlobalStructNewColNum[LocalStructOldNum[i]];

135:     for (jj = 0; jj < bs; jj++) {
136:       for (j=0; j<struct_size; j++){
137:         for ( k = 0; k<bs; k++){
138:           if (NewColNum[idx[j]] + k >= inew)
139:             MyANonz[MyANonz_last++] = ANonz[kk + idx[j]*bs2 + k*bs + jj];
140:         }
141:       }
142:       inew++;
143:     }
144:   } /* end outer loop for i */

146:   PetscFree(NewColNum);
147:   if (MyANonz_last != NumLocalNonz) SETERRQ2(PETSC_ERR_PLIB,"MyANonz_last %d != NumLocalNonz %d\n",MyANonz_last, NumLocalNonz);
148:   return(0);
149: }

153: PetscErrorCode MatDestroy_DSCPACK(Mat A)
154: {
155:   Mat_DSCPACK    *lu=(Mat_DSCPACK*)A->spptr;
157: 
159:   if (lu->CleanUpDSCPACK) {
160:     if (lu->dsc_id != -1) {
161:       if(lu->stat) DSC_DoStats(lu->My_DSC_Solver);
162:       DSC_FreeAll(lu->My_DSC_Solver);
163:       DSC_Close0(lu->My_DSC_Solver);
164: 
165:       PetscFree(lu->local_cols_old_num);
166:     }
167:     DSC_End(lu->My_DSC_Solver);
168: 
169:     MPI_Comm_free(&lu->comm_dsc);
170:     ISDestroy(lu->my_cols);
171:     PetscFree(lu->replication);
172:     VecDestroy(lu->vec_dsc);
173:     ISDestroy(lu->iden_dsc);
174:     VecScatterDestroy(lu->scat);
175:     if (lu->size >1 && lu->iden) {ISDestroy(lu->iden);}
176:   }
177:   if (lu->size == 1) {
178:     MatDestroy_SeqBAIJ(A);
179:   } else {
180:     MatDestroy_MPIBAIJ(A);
181:   }
182:   return(0);
183: }

187: PetscErrorCode MatSolve_DSCPACK(Mat A,Vec b,Vec x)
188: {
189:   Mat_DSCPACK    *lu= (Mat_DSCPACK*)A->spptr;
191:   RealNumberType *solution_vec,*rhs_vec;

194:   /* scatter b into seq vec_dsc */
195:   if ( !lu->scat ) {
196:     VecScatterCreate(b,lu->my_cols,lu->vec_dsc,lu->iden_dsc,&lu->scat);
197:   }
198:   VecScatterBegin(lu->scat,b,lu->vec_dsc,INSERT_VALUES,SCATTER_FORWARD);
199:   VecScatterEnd(lu->scat,b,lu->vec_dsc,INSERT_VALUES,SCATTER_FORWARD);

201:   if (lu->dsc_id != -1){
202:     VecGetArray(lu->vec_dsc,&rhs_vec);
203:     DSC_InputRhsLocalVec(lu->My_DSC_Solver, rhs_vec, lu->num_local_cols);
204:     VecRestoreArray(lu->vec_dsc,&rhs_vec);
205: 
206:     DSC_Solve(lu->My_DSC_Solver);
207:     if (ierr !=  DSC_NO_ERROR) {
208:       DSC_ErrorDisplay(lu->My_DSC_Solver);
209:       SETERRQ(PETSC_ERR_LIB,"Error in calling DSC_Solve");
210:     }

212:     /* get the permuted local solution */
213:     VecGetArray(lu->vec_dsc,&solution_vec);
214:     DSC_GetLocalSolution(lu->My_DSC_Solver,solution_vec, lu->num_local_cols);
215:     VecRestoreArray(lu->vec_dsc,&solution_vec);

217:   } /* end of if (lu->dsc_id != -1) */

219:   /* put permuted local solution solution_vec into x in the original order */
220:   VecScatterBegin(lu->scat,lu->vec_dsc,x,INSERT_VALUES,SCATTER_REVERSE);
221:   VecScatterEnd(lu->scat,lu->vec_dsc,x,INSERT_VALUES,SCATTER_REVERSE);

223:   return(0);
224: }

228: PetscErrorCode MatCholeskyFactorNumeric_DSCPACK(Mat F,Mat A,const MatFactorInfo *info)
229: {
230:   Mat_SeqBAIJ    *a_seq;
231:   Mat_DSCPACK    *lu=(Mat_DSCPACK*)(F)->spptr;
232:   Mat            *tseq,A_seq=PETSC_NULL;
233:   RealNumberType *my_a_nonz;
235:   PetscMPIInt    size;
236:   PetscInt       M=A->rmap->N,Mbs=M/lu->bs,max_mem_estimate,max_single_malloc_blk,
237:                  number_of_procs,i,j,next,iold,*idx,*iidx=0,*itmp;
238:   IS             my_cols_sorted;
239:   Mat            F_diag;
240: 
242:   MPI_Comm_size(((PetscObject)A)->comm,&size);
243:   if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
244:     /* convert A to A_seq */
245:     if (size > 1) {
246:       if (!lu->iden){
247:         ISCreateStride(PETSC_COMM_SELF,M,0,1,&lu->iden);
248:       }
249:       MatGetSubMatrices(A,1,&lu->iden,&lu->iden,MAT_INITIAL_MATRIX,&tseq);
250:       A_seq = tseq[0];
251:       a_seq = (Mat_SeqBAIJ*)A_seq->data;
252:     } else {
253:       a_seq = (Mat_SeqBAIJ*)A->data;
254:     }
255: 
256:     PetscMalloc(Mbs*sizeof(PetscInt),&lu->replication);
257:     for (i=0; i<Mbs; i++) lu->replication[i] = lu->bs;

259:     number_of_procs = DSC_Analyze(Mbs, a_seq->i, a_seq->j, lu->replication);
260: 
261:     i = size;
262:     if ( number_of_procs < i ) i = number_of_procs;
263:     number_of_procs = 1;
264:     while ( i > 1 ){
265:       number_of_procs  *= 2; i /= 2;
266:     }

268:     /* DSC_Solver starts */
269:     DSC_Open0( lu->My_DSC_Solver, number_of_procs, &lu->dsc_id, lu->comm_dsc );

271:     if (lu->dsc_id != -1) {
272:       DSC_Order(lu->My_DSC_Solver,lu->order_code,Mbs,a_seq->i,a_seq->j,lu->replication,
273:                    &M,&lu->num_local_strucs,
274:                    &lu->num_local_cols, &lu->num_local_nonz,  &lu->global_struc_new_col_num,
275:                    &lu->global_struc_new_num, &lu->global_struc_owner,
276:                    &lu->local_struc_old_num);
277:       if (ierr !=  DSC_NO_ERROR) {
278:         DSC_ErrorDisplay(lu->My_DSC_Solver);
279:         SETERRQ(PETSC_ERR_LIB,"Error when use DSC_Order()");
280:       }

282:       DSC_SFactor(lu->My_DSC_Solver,&max_mem_estimate,&max_single_malloc_blk,
283:                      lu->max_mem_allowed, lu->LBLASLevel, lu->DBLASLevel);
284:       if (ierr !=  DSC_NO_ERROR) {
285:         DSC_ErrorDisplay(lu->My_DSC_Solver);
286:         SETERRQ(PETSC_ERR_LIB,"Error when use DSC_Order");
287:       }

289:       BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
290:                        lu->num_local_strucs, lu->num_local_nonz,
291:                        lu->global_struc_new_col_num,
292:                        lu->local_struc_old_num,
293:                        PETSC_NULL,
294:                        &my_a_nonz);
295:       if (ierr <0) {
296:           DSC_ErrorDisplay(lu->My_DSC_Solver);
297:           SETERRQ1(PETSC_ERR_LIB,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
298:       }

300:       /* get local_cols_old_num and IS my_cols to be used later */
301:       PetscMalloc(lu->num_local_cols*sizeof(PetscInt),&lu->local_cols_old_num);
302:       for (next = 0, i=0; i<lu->num_local_strucs; i++){
303:         iold = lu->bs*lu->local_struc_old_num[i];
304:         for (j=0; j<lu->bs; j++)
305:           lu->local_cols_old_num[next++] = iold++;
306:       }
307:       ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,lu->local_cols_old_num,&lu->my_cols);
308: 
309:     } else {    /* lu->dsc_id == -1 */
310:       lu->num_local_cols = 0;
311:       lu->local_cols_old_num = 0;
312:       ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,lu->local_cols_old_num,&lu->my_cols);
313:     }
314:     /* generate vec_dsc and iden_dsc to be used later */
315:     VecCreateSeq(PETSC_COMM_SELF,lu->num_local_cols,&lu->vec_dsc);
316:     ISCreateStride(PETSC_COMM_SELF,lu->num_local_cols,0,1,&lu->iden_dsc);
317:     lu->scat = PETSC_NULL;

319:     if ( size>1 ) {
320:       MatDestroyMatrices(1,&tseq);
321:     }
322:   } else { /* use previously computed symbolic factor */
323:     /* convert A to my A_seq */
324:     if (size > 1) {
325:       if (lu->dsc_id == -1) {
326:         itmp = 0;
327:       } else {
328:         PetscMalloc(2*lu->num_local_strucs*sizeof(PetscInt),&idx);
329:         iidx = idx + lu->num_local_strucs;
330:         PetscMalloc(lu->num_local_cols*sizeof(PetscInt),&itmp);
331: 
332:         isort2(lu->num_local_strucs, lu->local_struc_old_num, idx);
333:         for (next=0, i=0; i< lu->num_local_strucs; i++) {
334:           iold = lu->bs*lu->local_struc_old_num[idx[i]];
335:           for (j=0; j<lu->bs; j++){
336:             itmp[next++] = iold++; /* sorted local_cols_old_num */
337:           }
338:         }
339:         for (i=0; i< lu->num_local_strucs; i++) {
340:           iidx[idx[i]] = i;       /* inverse of idx */
341:         }
342:       } /* end of (lu->dsc_id == -1) */
343:       ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,itmp,&my_cols_sorted);
344:       MatGetSubMatrices(A,1,&my_cols_sorted,&lu->iden,MAT_INITIAL_MATRIX,&tseq);
345:       ISDestroy(my_cols_sorted);
346:       A_seq = tseq[0];
347: 
348:       if (lu->dsc_id != -1) {
349:         DSC_ReFactorInitialize(lu->My_DSC_Solver);

351:         a_seq = (Mat_SeqBAIJ*)A_seq->data;
352:         BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
353:                        lu->num_local_strucs, lu->num_local_nonz,
354:                        lu->global_struc_new_col_num,
355:                        lu->local_struc_old_num,
356:                        iidx,
357:                        &my_a_nonz);
358:         if (ierr <0) {
359:           DSC_ErrorDisplay(lu->My_DSC_Solver);
360:           SETERRQ1(PETSC_ERR_LIB,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
361:         }
362:         PetscFree(idx);
363:         PetscFree(itmp);
364:       } /* end of if(lu->dsc_id != -1)  */
365:     } else { /* size == 1 */
366:       a_seq = (Mat_SeqBAIJ*)A->data;
367: 
368:       BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
369:                        lu->num_local_strucs, lu->num_local_nonz,
370:                        lu->global_struc_new_col_num,
371:                        lu->local_struc_old_num,
372:                        PETSC_NULL,
373:                        &my_a_nonz);
374:       if (ierr <0) {
375:         DSC_ErrorDisplay(lu->My_DSC_Solver);
376:         SETERRQ1(PETSC_ERR_LIB,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
377:       }
378:     }
379:     if ( size>1 ) {MatDestroyMatrices(1,&tseq); }
380:   }
381: 
382:   if (lu->dsc_id != -1) {
383:     DSC_NFactor(lu->My_DSC_Solver, lu->scheme_code, my_a_nonz, lu->factor_type, lu->LBLASLevel, lu->DBLASLevel);
384:     PetscFree(my_a_nonz);
385:   }
386: 
387:   if (size > 1) {
388:     F_diag = ((Mat_MPIBAIJ *)(F)->data)->A;
389:     F_diag->assembled = PETSC_TRUE;
390:   }
391:   F->assembled   = PETSC_TRUE;
392:   lu->flg           = SAME_NONZERO_PATTERN;
393:   F->ops->solve                  = MatSolve_DSCPACK;

395:   return(0);
396: }

398: /* Note the Petsc permutation r is ignored */
401: PetscErrorCode MatCholeskyFactorSymbolic_DSCPACK(Mat F,Mat A,IS r,const MatFactorInfo *info)
402: {
403:   Mat_DSCPACK    *lu = (Mat_DSCPACK*)(F)->spptr;

407:   lu->My_DSC_Solver = DSC_Begin();
408:   lu->CleanUpDSCPACK = PETSC_TRUE;
409:   (F)->ops->choleskyfactornumeric  = MatCholeskyFactorNumeric_DSCPACK;
410:   return(0);
411: }


417: PetscErrorCode MatFactorGetSolverPackage_seqaij_dscpack(Mat A,const MatSolverPackage *type)
418: {
420:   *type = MAT_SOLVER_DSCPACK;
421:   return(0);
422: }
424: 
427: PetscErrorCode MatGetFactor_seqbaij_dscpack(Mat A,MatFactorType ftype,Mat *F)
428: {
429:   Mat            B;
430:   Mat_DSCPACK    *lu;
432:   PetscInt       bs,indx;
433:   PetscTruth     flg;
434:   const char     *ftype[]={"LDLT","LLT"},*ltype[]={"LBLAS1","LBLAS2","LBLAS3"},*dtype[]={"DBLAS1","DBLAS2"};


438:   /* Create the factorization matrix F */
439:   MatGetBlockSize(A,&bs);
440:   MatCreate(((PetscObject)A)->comm,&B);
441:   MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
442:   MatSetType(B,((PetscObject)A)->type_name);
443:   MatSeqBAIJSetPreallocation(B,bs,0,PETSC_NULL);
444:   MatMPIBAIJSetPreallocation(B,bs,0,PETSC_NULL,0,PETSC_NULL);
445:   PetscNewLog(B,Mat_DSCPACKPACK,&lu);

447:   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_DSCPACK;
448:   B->ops->destroy                = MatDestroy_DSCPACK;
449:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_dscpack",MatFactorGetSolverPackage_seqaij_dscpack);
450:   B->factor                      = MAT_FACTOR_CHOLESKY;

452:   /* Set the default input options */
453:   lu->order_code  = 2;
454:   lu->scheme_code = 1;
455:   lu->factor_type = 2;
456:   lu->stat        = 0; /* do not display stats */
457:   lu->LBLASLevel  = DSC_LBLAS3;
458:   lu->DBLASLevel  = DSC_DBLAS2;
459:   lu->max_mem_allowed = 256;
460:   MPI_Comm_dup(((PetscObject)A)->comm,&lu->comm_dsc);
461:   /* Get the runtime input options */
462:   PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"DSCPACK Options","Mat");

464:   PetscOptionsInt("-mat_dscpack_order","order_code: \n\
465:          1 = ND, 2 = Hybrid with Minimum Degree, 3 = Hybrid with Minimum Deficiency", \
466:          "None",
467:          lu->order_code,&lu->order_code,PETSC_NULL);

469:   PetscOptionsInt("-mat_dscpack_scheme","scheme_code: \n\
470:          1 = standard factorization,  2 = factorization + selective inversion", \
471:          "None",
472:          lu->scheme_code,&lu->scheme_code,PETSC_NULL);
473: 
474:   PetscOptionsEList("-mat_dscpack_factor","factor_type","None",ftype,2,ftype[0],&indx,&flg);
475:   if (flg) {
476:     switch (indx) {
477:     case 0:
478:       lu->factor_type = DSC_LDLT;
479:       break;
480:     case 1:
481:       lu->factor_type = DSC_LLT;
482:       break;
483:     }
484:   }
485:   PetscOptionsInt("-mat_dscpack_MaxMemAllowed","in Mbytes","None",
486:          lu->max_mem_allowed,&lu->max_mem_allowed,PETSC_NULL);

488:   PetscOptionsInt("-mat_dscpack_stats","display stats: 0 = no display,  1 = display",
489:          "None", lu->stat,&lu->stat,PETSC_NULL);
490: 
491:   PetscOptionsEList("-mat_dscpack_LBLAS","BLAS level used in the local phase","None",ltype,3,ltype[2],&indx,&flg);
492:   if (flg) {
493:     switch (indx) {
494:     case 0:
495:       lu->LBLASLevel = DSC_LBLAS1;
496:       break;
497:     case 1:
498:       lu->LBLASLevel = DSC_LBLAS2;
499:       break;
500:     case 2:
501:       lu->LBLASLevel = DSC_LBLAS3;
502:       break;
503:     }
504:   }

506:   PetscOptionsEList("-mat_dscpack_DBLAS","BLAS level used in the distributed phase","None",dtype,2,dtype[1],&indx,&flg);
507:   if (flg) {
508:     switch (indx) {
509:     case 0:
510:       lu->DBLASLevel = DSC_DBLAS1;
511:       break;
512:     case 1:
513:       lu->DBLASLevel = DSC_DBLAS2;
514:       break;
515:     }
516:   }
517:   PetscOptionsEnd();
518:   lu->flg = DIFFERENT_NONZERO_PATTERN;
519:   *F = B;
520:   return(0);
521: }

525: PetscErrorCode MatFactorInfo_DSCPACK(Mat A,PetscViewer viewer)
526: {
527:   Mat_DSCPACK    *lu=(Mat_DSCPACK*)A->spptr;
529:   const char     *s=0;
530: 
532:   PetscViewerASCIIPrintf(viewer,"DSCPACK run parameters:\n");

534:   switch (lu->order_code) {
535:   case 1: s = "ND"; break;
536:   case 2: s = "Hybrid with Minimum Degree"; break;
537:   case 3: s = "Hybrid with Minimum Deficiency"; break;
538:   }
539:   PetscViewerASCIIPrintf(viewer,"  order_code: %s \n",s);

541:   switch (lu->scheme_code) {
542:   case 1: s = "standard factorization"; break;
543:   case 2: s = "factorization + selective inversion"; break;
544:   }
545:   PetscViewerASCIIPrintf(viewer,"  scheme_code: %s \n",s);

547:   switch (lu->stat) {
548:   case 0: s = "NO"; break;
549:   case 1: s = "YES"; break;
550:   }
551:   PetscViewerASCIIPrintf(viewer,"  display stats: %s \n",s);
552: 
553:   if ( lu->factor_type == DSC_LLT) {
554:     s = "LLT";
555:   } else if ( lu->factor_type == DSC_LDLT){
556:     s = "LDLT";
557:   } else if (lu->factor_type == 0) {
558:     s = "None";
559:   } else {
560:     SETERRQ(PETSC_ERR_PLIB,"Unknown factor type");
561:   }
562:   PetscViewerASCIIPrintf(viewer,"  factor type: %s \n",s);

564:   if ( lu->LBLASLevel == DSC_LBLAS1) {
565:     s = "BLAS1";
566:   } else if ( lu->LBLASLevel == DSC_LBLAS2){
567:     s = "BLAS2";
568:   } else if ( lu->LBLASLevel == DSC_LBLAS3){
569:     s = "BLAS3";
570:   } else if (lu->LBLASLevel == 0) {
571:     s = "None";
572:   } else {
573:     SETERRQ(PETSC_ERR_PLIB,"Unknown local phase BLAS level");
574:   }
575:   PetscViewerASCIIPrintf(viewer,"  local phase BLAS level: %s \n",s);
576: 
577:   if ( lu->DBLASLevel == DSC_DBLAS1) {
578:     s = "BLAS1";
579:   } else if ( lu->DBLASLevel == DSC_DBLAS2){
580:     s = "BLAS2";
581:   } else if (lu->DBLASLevel == 0) {
582:     s = "None";
583:   } else {
584:     SETERRQ(PETSC_ERR_PLIB,"Unknown distributed phase BLAS level");
585:   }
586:   PetscViewerASCIIPrintf(viewer,"  distributed phase BLAS level: %s \n",s);
587:   return(0);
588: }

590: EXTERN PetscErrorCode MatView_SeqBAIJ(Mat,PetscViewer);
591: EXTERN PetscErrorCode MatView_MPIBAIJ(Mat,PetscViewer);


596: PetscErrorCode MatView_DSCPACK(Mat A,PetscViewer viewer)
597: {
598:   PetscErrorCode    ierr;
599:   PetscMPIInt       size;
600:   PetscTruth        iascii;
601:   PetscViewerFormat format;
602:   Mat_DSCPACK       *lu=(Mat_DSCPACK*)A->spptr;


606:   /* This convertion ugliness is because MatView for BAIJ types calls MatConvert to AIJ */
607:   size = lu->size;
608:   if (size==1) {
609:     MatView_SeqBAIJ(A,viewer);
610:   } else {
611:     MatView_MPIBAIJ(A,viewer);
612:   }

614:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
615:   if (iascii) {
616:     PetscViewerGetFormat(viewer,&format);
617:     if (format == PETSC_VIEWER_ASCII_INFO) {
618:       MatFactorInfo_DSCPACK(A,viewer);
619:     }
620:   }
621:   return(0);
622: }


625: /*MC
626:   MAT_SOLVER_DSCPACK -  "dscpack" - Provides direct solvers (Cholesky) for sequential 
627:   or distributed matrices via the external package DSCPACK.


630:   Options Database Keys:
631: + -mat_dscpack_order <1,2,3> - DSCPACK ordering, 1:ND, 2:Hybrid with Minimum Degree, 3:Hybrid with Minimum Deficiency
632: . -mat_dscpack_scheme <1,2> - factorization scheme, 1:standard factorization,  2: factorization with selective inversion
633: . -mat_dscpack_factor <LLT,LDLT> - the type of factorization to be performed.
634: . -mat_dscpack_MaxMemAllowed <n> - the maximum memory to be used during factorization
635: . -mat_dscpack_stats <0,1> - display stats of the factorization and solves during MatDestroy(), 0: no display,  1: display
636: . -mat_dscpack_LBLAS <LBLAS1,LBLAS2,LBLAS3> - BLAS level used in the local phase
637: - -mat_dscpack_DBLAS <DBLAS1,DBLAS2> - BLAS level used in the distributed phase

639:    Level: beginner

641: .seealso: PCCHOLESKY, PCFactorSetMatSolverPackage(), MatSolverPackage

643: M*/