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,*idx, *NewColNum, MyANonz_last, max_struct=0, struct_size;
87: RealNumberType *MyANonz;
91: /* loop: to find maximum number of subscripts over columns
92: assigned to this processor */
93: for (i=0; i <NumLocalStructs; i++) {
94: /* for each struct i (local) assigned to this processor */
95: if (LocalStructLocalNum){
96: iold = LocalStructLocalNum[i];
97: } else {
98: iold = LocalStructOldNum[i];
99: }
100:
101: struct_size = AIndex[iold+1] - AIndex[iold];
102: if ( max_struct <= struct_size) max_struct = struct_size;
103: }
105: /* allocate tmp arrays large enough to hold densest struct */
106: PetscMalloc2(max_struct,PetscInt,&NewColNum,max_struct,PetscInt,&idx);
107:
108: PetscMalloc(NumLocalNonz*sizeof(RealNumberType),&MyANonz);
109: *adr_MyANonz = MyANonz;
111: /* loop to set up nonzeroes in MyANonz */
112: MyANonz_last = 0 ; /* points to first empty space in MyANonz */
113: for (i=0; i <NumLocalStructs; i++) {
115: /* for each struct i (local) assigned to this processor */
116: if (LocalStructLocalNum){
117: iold = LocalStructLocalNum[i];
118: } else {
119: iold = LocalStructOldNum[i];
120: }
122: struct_size = AIndex[iold+1] - AIndex[iold];
123: for (k=0, j=AIndex[iold]; j<AIndex[iold+1]; j++){
124: NewColNum[k] = GlobalStructNewColNum[AStruct[j]];
125: k++;
126: }
127: isort2(struct_size, NewColNum, idx);
128:
129: kk = AIndex[iold]*bs2; /* points to 1st element of iold block col in ANonz */
130: inew = GlobalStructNewColNum[LocalStructOldNum[i]];
132: for (jj = 0; jj < bs; jj++) {
133: for (j=0; j<struct_size; j++){
134: for ( k = 0; k<bs; k++){
135: if (NewColNum[idx[j]] + k >= inew)
136: MyANonz[MyANonz_last++] = ANonz[kk + idx[j]*bs2 + k*bs + jj];
137: }
138: }
139: inew++;
140: }
141: } /* end outer loop for i */
143: PetscFree2(NewColNum);
144: if (MyANonz_last != NumLocalNonz) SETERRQ2(PETSC_ERR_PLIB,"MyANonz_last %d != NumLocalNonz %d\n",MyANonz_last, NumLocalNonz);
145: return(0);
146: }
150: PetscErrorCode MatDestroy_DSCPACK(Mat A)
151: {
152: Mat_DSCPACK *lu=(Mat_DSCPACK*)A->spptr;
154:
156: if (lu->CleanUpDSCPACK) {
157: if (lu->dsc_id != -1) {
158: if(lu->stat) DSC_DoStats(lu->My_DSC_Solver);
159: DSC_FreeAll(lu->My_DSC_Solver);
160: DSC_Close0(lu->My_DSC_Solver);
161:
162: PetscFree(lu->local_cols_old_num);
163: }
164: DSC_End(lu->My_DSC_Solver);
165:
166: MPI_Comm_free(&lu->comm_dsc);
167: ISDestroy(lu->my_cols);
168: PetscFree(lu->replication);
169: VecDestroy(lu->vec_dsc);
170: ISDestroy(lu->iden_dsc);
171: VecScatterDestroy(lu->scat);
172: if (lu->size >1 && lu->iden) {ISDestroy(lu->iden);}
173: }
174: if (lu->size == 1) {
175: MatDestroy_SeqBAIJ(A);
176: } else {
177: MatDestroy_MPIBAIJ(A);
178: }
179: return(0);
180: }
184: PetscErrorCode MatSolve_DSCPACK(Mat A,Vec b,Vec x)
185: {
186: Mat_DSCPACK *lu= (Mat_DSCPACK*)A->spptr;
188: RealNumberType *solution_vec,*rhs_vec;
191: /* scatter b into seq vec_dsc */
192: if ( !lu->scat ) {
193: VecScatterCreate(b,lu->my_cols,lu->vec_dsc,lu->iden_dsc,&lu->scat);
194: }
195: VecScatterBegin(lu->scat,b,lu->vec_dsc,INSERT_VALUES,SCATTER_FORWARD);
196: VecScatterEnd(lu->scat,b,lu->vec_dsc,INSERT_VALUES,SCATTER_FORWARD);
198: if (lu->dsc_id != -1){
199: VecGetArray(lu->vec_dsc,&rhs_vec);
200: DSC_InputRhsLocalVec(lu->My_DSC_Solver, rhs_vec, lu->num_local_cols);
201: VecRestoreArray(lu->vec_dsc,&rhs_vec);
202:
203: DSC_Solve(lu->My_DSC_Solver);
204: if (ierr != DSC_NO_ERROR) {
205: DSC_ErrorDisplay(lu->My_DSC_Solver);
206: SETERRQ(PETSC_ERR_LIB,"Error in calling DSC_Solve");
207: }
209: /* get the permuted local solution */
210: VecGetArray(lu->vec_dsc,&solution_vec);
211: DSC_GetLocalSolution(lu->My_DSC_Solver,solution_vec, lu->num_local_cols);
212: VecRestoreArray(lu->vec_dsc,&solution_vec);
214: } /* end of if (lu->dsc_id != -1) */
216: /* put permuted local solution solution_vec into x in the original order */
217: VecScatterBegin(lu->scat,lu->vec_dsc,x,INSERT_VALUES,SCATTER_REVERSE);
218: VecScatterEnd(lu->scat,lu->vec_dsc,x,INSERT_VALUES,SCATTER_REVERSE);
220: return(0);
221: }
225: PetscErrorCode MatCholeskyFactorNumeric_DSCPACK(Mat F,Mat A,const MatFactorInfo *info)
226: {
227: Mat_SeqBAIJ *a_seq;
228: Mat_DSCPACK *lu=(Mat_DSCPACK*)(F)->spptr;
229: Mat *tseq,A_seq=PETSC_NULL;
230: RealNumberType *my_a_nonz;
232: PetscMPIInt size;
233: PetscInt M=A->rmap->N,Mbs=M/lu->bs,max_mem_estimate,max_single_malloc_blk,
234: number_of_procs,i,j,next,iold,*idx,*iidx=0,*itmp;
235: IS my_cols_sorted;
236: Mat F_diag;
237:
239: MPI_Comm_size(((PetscObject)A)->comm,&size);
240: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
241: /* convert A to A_seq */
242: if (size > 1) {
243: if (!lu->iden){
244: ISCreateStride(PETSC_COMM_SELF,M,0,1,&lu->iden);
245: }
246: MatGetSubMatrices(A,1,&lu->iden,&lu->iden,MAT_INITIAL_MATRIX,&tseq);
247: A_seq = tseq[0];
248: a_seq = (Mat_SeqBAIJ*)A_seq->data;
249: } else {
250: a_seq = (Mat_SeqBAIJ*)A->data;
251: }
252:
253: PetscMalloc(Mbs*sizeof(PetscInt),&lu->replication);
254: for (i=0; i<Mbs; i++) lu->replication[i] = lu->bs;
256: number_of_procs = DSC_Analyze(Mbs, a_seq->i, a_seq->j, lu->replication);
257:
258: i = size;
259: if ( number_of_procs < i ) i = number_of_procs;
260: number_of_procs = 1;
261: while ( i > 1 ){
262: number_of_procs *= 2; i /= 2;
263: }
265: /* DSC_Solver starts */
266: DSC_Open0( lu->My_DSC_Solver, number_of_procs, &lu->dsc_id, lu->comm_dsc );
268: if (lu->dsc_id != -1) {
269: DSC_Order(lu->My_DSC_Solver,lu->order_code,Mbs,a_seq->i,a_seq->j,lu->replication,
270: &M,&lu->num_local_strucs,
271: &lu->num_local_cols, &lu->num_local_nonz, &lu->global_struc_new_col_num,
272: &lu->global_struc_new_num, &lu->global_struc_owner,
273: &lu->local_struc_old_num);
274: if (ierr != DSC_NO_ERROR) {
275: DSC_ErrorDisplay(lu->My_DSC_Solver);
276: SETERRQ(PETSC_ERR_LIB,"Error when use DSC_Order()");
277: }
279: DSC_SFactor(lu->My_DSC_Solver,&max_mem_estimate,&max_single_malloc_blk,
280: lu->max_mem_allowed, lu->LBLASLevel, lu->DBLASLevel);
281: if (ierr != DSC_NO_ERROR) {
282: DSC_ErrorDisplay(lu->My_DSC_Solver);
283: SETERRQ(PETSC_ERR_LIB,"Error when use DSC_Order");
284: }
286: BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
287: lu->num_local_strucs, lu->num_local_nonz,
288: lu->global_struc_new_col_num,
289: lu->local_struc_old_num,
290: PETSC_NULL,
291: &my_a_nonz);
292: if (ierr <0) {
293: DSC_ErrorDisplay(lu->My_DSC_Solver);
294: SETERRQ1(PETSC_ERR_LIB,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
295: }
297: /* get local_cols_old_num and IS my_cols to be used later */
298: PetscMalloc(lu->num_local_cols*sizeof(PetscInt),&lu->local_cols_old_num);
299: for (next = 0, i=0; i<lu->num_local_strucs; i++){
300: iold = lu->bs*lu->local_struc_old_num[i];
301: for (j=0; j<lu->bs; j++)
302: lu->local_cols_old_num[next++] = iold++;
303: }
304: ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,lu->local_cols_old_num,&lu->my_cols);
305:
306: } else { /* lu->dsc_id == -1 */
307: lu->num_local_cols = 0;
308: lu->local_cols_old_num = 0;
309: ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,lu->local_cols_old_num,&lu->my_cols);
310: }
311: /* generate vec_dsc and iden_dsc to be used later */
312: VecCreateSeq(PETSC_COMM_SELF,lu->num_local_cols,&lu->vec_dsc);
313: ISCreateStride(PETSC_COMM_SELF,lu->num_local_cols,0,1,&lu->iden_dsc);
314: lu->scat = PETSC_NULL;
316: if ( size>1 ) {
317: MatDestroyMatrices(1,&tseq);
318: }
319: } else { /* use previously computed symbolic factor */
320: /* convert A to my A_seq */
321: if (size > 1) {
322: if (lu->dsc_id == -1) {
323: itmp = 0;
324: } else {
325: PetscMalloc2(lu->num_local_strucs,PetscInt,&idx,lu->num_local_strucs,PetscInt,&iidx);
326: PetscMalloc(lu->num_local_cols*sizeof(PetscInt),&itmp);
327:
328: isort2(lu->num_local_strucs, lu->local_struc_old_num, idx);
329: for (next=0, i=0; i< lu->num_local_strucs; i++) {
330: iold = lu->bs*lu->local_struc_old_num[idx[i]];
331: for (j=0; j<lu->bs; j++){
332: itmp[next++] = iold++; /* sorted local_cols_old_num */
333: }
334: }
335: for (i=0; i< lu->num_local_strucs; i++) {
336: iidx[idx[i]] = i; /* inverse of idx */
337: }
338: } /* end of (lu->dsc_id == -1) */
339: ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,itmp,&my_cols_sorted);
340: MatGetSubMatrices(A,1,&my_cols_sorted,&lu->iden,MAT_INITIAL_MATRIX,&tseq);
341: ISDestroy(my_cols_sorted);
342: A_seq = tseq[0];
343:
344: if (lu->dsc_id != -1) {
345: DSC_ReFactorInitialize(lu->My_DSC_Solver);
347: a_seq = (Mat_SeqBAIJ*)A_seq->data;
348: BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
349: lu->num_local_strucs, lu->num_local_nonz,
350: lu->global_struc_new_col_num,
351: lu->local_struc_old_num,
352: iidx,
353: &my_a_nonz);
354: if (ierr <0) {
355: DSC_ErrorDisplay(lu->My_DSC_Solver);
356: SETERRQ1(PETSC_ERR_LIB,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
357: }
358: PetscFree2(idx,iidex);
359: PetscFree(itmp);
360: } /* end of if(lu->dsc_id != -1) */
361: } else { /* size == 1 */
362: a_seq = (Mat_SeqBAIJ*)A->data;
363:
364: BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
365: lu->num_local_strucs, lu->num_local_nonz,
366: lu->global_struc_new_col_num,
367: lu->local_struc_old_num,
368: PETSC_NULL,
369: &my_a_nonz);
370: if (ierr <0) {
371: DSC_ErrorDisplay(lu->My_DSC_Solver);
372: SETERRQ1(PETSC_ERR_LIB,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
373: }
374: }
375: if ( size>1 ) {MatDestroyMatrices(1,&tseq); }
376: }
377:
378: if (lu->dsc_id != -1) {
379: DSC_NFactor(lu->My_DSC_Solver, lu->scheme_code, my_a_nonz, lu->factor_type, lu->LBLASLevel, lu->DBLASLevel);
380: PetscFree(my_a_nonz);
381: }
382:
383: if (size > 1) {
384: F_diag = ((Mat_MPIBAIJ *)(F)->data)->A;
385: F_diag->assembled = PETSC_TRUE;
386: }
387: F->assembled = PETSC_TRUE;
388: lu->flg = SAME_NONZERO_PATTERN;
389: F->ops->solve = MatSolve_DSCPACK;
391: return(0);
392: }
394: /* Note the Petsc permutation r is ignored */
397: PetscErrorCode MatCholeskyFactorSymbolic_DSCPACK(Mat F,Mat A,IS r,const MatFactorInfo *info)
398: {
399: Mat_DSCPACK *lu = (Mat_DSCPACK*)(F)->spptr;
403: lu->My_DSC_Solver = DSC_Begin();
404: lu->CleanUpDSCPACK = PETSC_TRUE;
405: (F)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_DSCPACK;
406: return(0);
407: }
413: PetscErrorCode MatFactorGetSolverPackage_seqaij_dscpack(Mat A,const MatSolverPackage *type)
414: {
416: *type = MAT_SOLVER_DSCPACK;
417: return(0);
418: }
420:
423: PetscErrorCode MatGetFactor_seqbaij_dscpack(Mat A,MatFactorType ftype,Mat *F)
424: {
425: Mat B;
426: Mat_DSCPACK *lu;
428: PetscInt bs,indx;
429: PetscTruth flg;
430: const char *ftype[]={"LDLT","LLT"},*ltype[]={"LBLAS1","LBLAS2","LBLAS3"},*dtype[]={"DBLAS1","DBLAS2"};
434: /* Create the factorization matrix F */
435: MatGetBlockSize(A,&bs);
436: MatCreate(((PetscObject)A)->comm,&B);
437: MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
438: MatSetType(B,((PetscObject)A)->type_name);
439: MatSeqBAIJSetPreallocation(B,bs,0,PETSC_NULL);
440: MatMPIBAIJSetPreallocation(B,bs,0,PETSC_NULL,0,PETSC_NULL);
441: PetscNewLog(B,Mat_DSCPACKPACK,&lu);
443: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_DSCPACK;
444: B->ops->destroy = MatDestroy_DSCPACK;
445: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_seqaij_dscpack",MatFactorGetSolverPackage_seqaij_dscpack);
446: B->factor = MAT_FACTOR_CHOLESKY;
448: /* Set the default input options */
449: lu->order_code = 2;
450: lu->scheme_code = 1;
451: lu->factor_type = 2;
452: lu->stat = 0; /* do not display stats */
453: lu->LBLASLevel = DSC_LBLAS3;
454: lu->DBLASLevel = DSC_DBLAS2;
455: lu->max_mem_allowed = 256;
456: MPI_Comm_dup(((PetscObject)A)->comm,&lu->comm_dsc);
457: /* Get the runtime input options */
458: PetscOptionsBegin(((PetscObject)A)->comm,((PetscObject)A)->prefix,"DSCPACK Options","Mat");
460: PetscOptionsInt("-mat_dscpack_order","order_code: \n\
461: 1 = ND, 2 = Hybrid with Minimum Degree, 3 = Hybrid with Minimum Deficiency", \
462: "None",
463: lu->order_code,&lu->order_code,PETSC_NULL);
465: PetscOptionsInt("-mat_dscpack_scheme","scheme_code: \n\
466: 1 = standard factorization, 2 = factorization + selective inversion", \
467: "None",
468: lu->scheme_code,&lu->scheme_code,PETSC_NULL);
469:
470: PetscOptionsEList("-mat_dscpack_factor","factor_type","None",ftype,2,ftype[0],&indx,&flg);
471: if (flg) {
472: switch (indx) {
473: case 0:
474: lu->factor_type = DSC_LDLT;
475: break;
476: case 1:
477: lu->factor_type = DSC_LLT;
478: break;
479: }
480: }
481: PetscOptionsInt("-mat_dscpack_MaxMemAllowed","in Mbytes","None",
482: lu->max_mem_allowed,&lu->max_mem_allowed,PETSC_NULL);
484: PetscOptionsInt("-mat_dscpack_stats","display stats: 0 = no display, 1 = display",
485: "None", lu->stat,&lu->stat,PETSC_NULL);
486:
487: PetscOptionsEList("-mat_dscpack_LBLAS","BLAS level used in the local phase","None",ltype,3,ltype[2],&indx,&flg);
488: if (flg) {
489: switch (indx) {
490: case 0:
491: lu->LBLASLevel = DSC_LBLAS1;
492: break;
493: case 1:
494: lu->LBLASLevel = DSC_LBLAS2;
495: break;
496: case 2:
497: lu->LBLASLevel = DSC_LBLAS3;
498: break;
499: }
500: }
502: PetscOptionsEList("-mat_dscpack_DBLAS","BLAS level used in the distributed phase","None",dtype,2,dtype[1],&indx,&flg);
503: if (flg) {
504: switch (indx) {
505: case 0:
506: lu->DBLASLevel = DSC_DBLAS1;
507: break;
508: case 1:
509: lu->DBLASLevel = DSC_DBLAS2;
510: break;
511: }
512: }
513: PetscOptionsEnd();
514: lu->flg = DIFFERENT_NONZERO_PATTERN;
515: *F = B;
516: return(0);
517: }
521: PetscErrorCode MatFactorInfo_DSCPACK(Mat A,PetscViewer viewer)
522: {
523: Mat_DSCPACK *lu=(Mat_DSCPACK*)A->spptr;
525: const char *s=0;
526:
528: PetscViewerASCIIPrintf(viewer,"DSCPACK run parameters:\n");
530: switch (lu->order_code) {
531: case 1: s = "ND"; break;
532: case 2: s = "Hybrid with Minimum Degree"; break;
533: case 3: s = "Hybrid with Minimum Deficiency"; break;
534: }
535: PetscViewerASCIIPrintf(viewer," order_code: %s \n",s);
537: switch (lu->scheme_code) {
538: case 1: s = "standard factorization"; break;
539: case 2: s = "factorization + selective inversion"; break;
540: }
541: PetscViewerASCIIPrintf(viewer," scheme_code: %s \n",s);
543: switch (lu->stat) {
544: case 0: s = "NO"; break;
545: case 1: s = "YES"; break;
546: }
547: PetscViewerASCIIPrintf(viewer," display stats: %s \n",s);
548:
549: if ( lu->factor_type == DSC_LLT) {
550: s = "LLT";
551: } else if ( lu->factor_type == DSC_LDLT){
552: s = "LDLT";
553: } else if (lu->factor_type == 0) {
554: s = "None";
555: } else {
556: SETERRQ(PETSC_ERR_PLIB,"Unknown factor type");
557: }
558: PetscViewerASCIIPrintf(viewer," factor type: %s \n",s);
560: if ( lu->LBLASLevel == DSC_LBLAS1) {
561: s = "BLAS1";
562: } else if ( lu->LBLASLevel == DSC_LBLAS2){
563: s = "BLAS2";
564: } else if ( lu->LBLASLevel == DSC_LBLAS3){
565: s = "BLAS3";
566: } else if (lu->LBLASLevel == 0) {
567: s = "None";
568: } else {
569: SETERRQ(PETSC_ERR_PLIB,"Unknown local phase BLAS level");
570: }
571: PetscViewerASCIIPrintf(viewer," local phase BLAS level: %s \n",s);
572:
573: if ( lu->DBLASLevel == DSC_DBLAS1) {
574: s = "BLAS1";
575: } else if ( lu->DBLASLevel == DSC_DBLAS2){
576: s = "BLAS2";
577: } else if (lu->DBLASLevel == 0) {
578: s = "None";
579: } else {
580: SETERRQ(PETSC_ERR_PLIB,"Unknown distributed phase BLAS level");
581: }
582: PetscViewerASCIIPrintf(viewer," distributed phase BLAS level: %s \n",s);
583: return(0);
584: }
586: EXTERN PetscErrorCode MatView_SeqBAIJ(Mat,PetscViewer);
587: EXTERN PetscErrorCode MatView_MPIBAIJ(Mat,PetscViewer);
592: PetscErrorCode MatView_DSCPACK(Mat A,PetscViewer viewer)
593: {
594: PetscErrorCode ierr;
595: PetscMPIInt size;
596: PetscTruth iascii;
597: PetscViewerFormat format;
598: Mat_DSCPACK *lu=(Mat_DSCPACK*)A->spptr;
602: /* This convertion ugliness is because MatView for BAIJ types calls MatConvert to AIJ */
603: size = lu->size;
604: if (size==1) {
605: MatView_SeqBAIJ(A,viewer);
606: } else {
607: MatView_MPIBAIJ(A,viewer);
608: }
610: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
611: if (iascii) {
612: PetscViewerGetFormat(viewer,&format);
613: if (format == PETSC_VIEWER_ASCII_INFO) {
614: MatFactorInfo_DSCPACK(A,viewer);
615: }
616: }
617: return(0);
618: }
621: /*MC
622: MAT_SOLVER_DSCPACK - "dscpack" - Provides direct solvers (Cholesky) for sequential
623: or distributed matrices via the external package DSCPACK.
626: Options Database Keys:
627: + -mat_dscpack_order <1,2,3> - DSCPACK ordering, 1:ND, 2:Hybrid with Minimum Degree, 3:Hybrid with Minimum Deficiency
628: . -mat_dscpack_scheme <1,2> - factorization scheme, 1:standard factorization, 2: factorization with selective inversion
629: . -mat_dscpack_factor <LLT,LDLT> - the type of factorization to be performed.
630: . -mat_dscpack_MaxMemAllowed <n> - the maximum memory to be used during factorization
631: . -mat_dscpack_stats <0,1> - display stats of the factorization and solves during MatDestroy(), 0: no display, 1: display
632: . -mat_dscpack_LBLAS <LBLAS1,LBLAS2,LBLAS3> - BLAS level used in the local phase
633: - -mat_dscpack_DBLAS <DBLAS1,DBLAS2> - BLAS level used in the distributed phase
635: Level: beginner
637: .seealso: PCCHOLESKY, PCFactorSetMatSolverPackage(), MatSolverPackage
639: M*/