Actual source code: dscpack.c
1: /*$Id: dscpack.c,v 1.10 2001/08/15 15:56:50 bsmith Exp $*/
2: /*
3: Provides an interface to the DSCPACK (Domain-Separator Codes) sparse direct solver
4: */
6: #include src/mat/impls/baij/seq/baij.h
7: #include src/mat/impls/baij/mpi/mpibaij.h
9: EXTERN_C_BEGIN
10: #include "dscmain.h"
11: EXTERN_C_END
13: typedef struct {
14: DSC_Solver My_DSC_Solver;
15: int num_local_strucs, *local_struc_old_num,
16: num_local_cols, num_local_nonz,
17: *global_struc_new_col_num,
18: *global_struc_new_num, *global_struc_owner,
19: dsc_id,bs,*local_cols_old_num,*replication;
20: int order_code,scheme_code,factor_type, stat,
21: LBLASLevel,DBLASLevel,max_mem_allowed;
22: MatStructure flg;
23: IS my_cols,iden,iden_dsc;
24: Vec vec_dsc;
25: VecScatter scat;
26: MPI_Comm comm_dsc;
28: /* A few inheritance details */
29: int size;
30: int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
31: int (*MatView)(Mat,PetscViewer);
32: int (*MatAssemblyEnd)(Mat,MatAssemblyType);
33: int (*MatCholeskyFactorSymbolic)(Mat,IS,MatFactorInfo*,Mat*);
34: int (*MatDestroy)(Mat);
35: /* Clean up flag for destructor */
36: PetscTruth CleanUpDSCPACK;
37: } Mat_DSC;
39: EXTERN int MatDuplicate_DSCPACK(Mat,MatDuplicateOption,Mat*);
40: /* DSC function */
43: void isort2(int size, int *list, int *index) {
44: /* in increasing order */
45: /* index will contain indices such that */
46: /* list can be accessed in sorted order */
47: int i, j, x, y;
48:
49: for (i=0; i<size; i++) index[i] =i;
51: for (i=1; i<size; i++){
52: y= index[i];
53: x=list[index[i]];
54: for (j=i-1; ((j>=0) && (x<list[index[j]])); j--)
55: index[j+1]=index[j];
56: index[j+1]=y;
57: }
58: }/*end isort2*/
62: int BAIJtoMyANonz( int *AIndex, int *AStruct, int bs,
63: RealNumberType *ANonz, int NumLocalStructs,
64: int NumLocalNonz, int *GlobalStructNewColNum,
65: int *LocalStructOldNum,
66: int *LocalStructLocalNum,
67: RealNumberType **adr_MyANonz)
68: /*
69: Extract non-zero values of lower triangular part
70: of the permuted matrix that belong to this processor.
72: Only output parameter is adr_MyANonz -- is malloced and changed.
73: Rest are input parameters left unchanged.
75: When LocalStructLocalNum == PETSC_NULL,
76: AIndex, AStruct, and ANonz contain entire original matrix A
77: in PETSc SeqBAIJ format,
78: otherwise,
79: AIndex, AStruct, and ANonz are indeces for the submatrix
80: of A whose colomns (in increasing order) belong to this processor.
82: Other variables supply information on ownership of columns
83: and the new numbering in a fill-reducing permutation
85: This information is used to setup lower half of A nonzeroes
86: for columns owned by this processor
87: */
88: {
89: int i, j, k, iold,inew, jj, kk,ierr, bs2=bs*bs,
90: *idx, *NewColNum,
91: MyANonz_last, max_struct=0, struct_size;
92: RealNumberType *MyANonz;
96: /* loop: to find maximum number of subscripts over columns
97: assigned to this processor */
98: for (i=0; i <NumLocalStructs; i++) {
99: /* for each struct i (local) assigned to this processor */
100: if (LocalStructLocalNum){
101: iold = LocalStructLocalNum[i];
102: } else {
103: iold = LocalStructOldNum[i];
104: }
105:
106: struct_size = AIndex[iold+1] - AIndex[iold];
107: if ( max_struct <= struct_size) max_struct = struct_size;
108: }
110: /* allocate tmp arrays large enough to hold densest struct */
111: PetscMalloc((2*max_struct+1)*sizeof(int),&NewColNum);
112: idx = NewColNum + max_struct;
113:
114: PetscMalloc(NumLocalNonz*sizeof(RealNumberType),&MyANonz);
115: *adr_MyANonz = MyANonz;
117: /* loop to set up nonzeroes in MyANonz */
118: MyANonz_last = 0 ; /* points to first empty space in MyANonz */
119: for (i=0; i <NumLocalStructs; i++) {
121: /* for each struct i (local) assigned to this processor */
122: if (LocalStructLocalNum){
123: iold = LocalStructLocalNum[i];
124: } else {
125: iold = LocalStructOldNum[i];
126: }
128: struct_size = AIndex[iold+1] - AIndex[iold];
129: for (k=0, j=AIndex[iold]; j<AIndex[iold+1]; j++){
130: NewColNum[k] = GlobalStructNewColNum[AStruct[j]];
131: k++;
132: }
133: isort2(struct_size, NewColNum, idx);
134:
135: kk = AIndex[iold]*bs2; /* points to 1st element of iold block col in ANonz */
136: inew = GlobalStructNewColNum[LocalStructOldNum[i]];
138: for (jj = 0; jj < bs; jj++) {
139: for (j=0; j<struct_size; j++){
140: for ( k = 0; k<bs; k++){
141: if (NewColNum[idx[j]] + k >= inew)
142: MyANonz[MyANonz_last++] = ANonz[kk + idx[j]*bs2 + k*bs + jj];
143: }
144: }
145: inew++;
146: }
147: } /* end outer loop for i */
149: PetscFree(NewColNum);
150: if (MyANonz_last != NumLocalNonz)
151: SETERRQ2(1,"MyANonz_last %d != NumLocalNonz %d\n",MyANonz_last, NumLocalNonz);
152: return(0);
153: }
155: EXTERN_C_BEGIN
158: int MatConvert_DSCPACK_Base(Mat A,MatType type,Mat *newmat) {
159: int ierr;
160: Mat B=*newmat;
161: Mat_DSC *lu=(Mat_DSC*)A->spptr;
162:
164: if (B != A) {
165: MatDuplicate(A,MAT_COPY_VALUES,&B);
166: }
167: /* Reset the original function pointers */
168: B->ops->duplicate = lu->MatDuplicate;
169: B->ops->view = lu->MatView;
170: B->ops->assemblyend = lu->MatAssemblyEnd;
171: B->ops->choleskyfactorsymbolic = lu->MatCholeskyFactorSymbolic;
172: B->ops->destroy = lu->MatDestroy;
174: PetscObjectChangeTypeName((PetscObject)B,type);
175: PetscFree(lu);
176: *newmat = B;
178: return(0);
179: }
180: EXTERN_C_END
184: int MatDestroy_DSCPACK(Mat A) {
185: Mat_DSC *lu=(Mat_DSC*)A->spptr;
186: int ierr;
187:
189: if (lu->CleanUpDSCPACK) {
190: if (lu->dsc_id != -1) {
191: if(lu->stat) DSC_DoStats(lu->My_DSC_Solver);
192: DSC_FreeAll(lu->My_DSC_Solver);
193: DSC_Close0(lu->My_DSC_Solver);
194:
195: PetscFree(lu->local_cols_old_num);
196: }
197: DSC_End(lu->My_DSC_Solver);
198:
199: MPI_Comm_free(&(lu->comm_dsc));
200: ISDestroy(lu->my_cols);
201: PetscFree(lu->replication);
202: VecDestroy(lu->vec_dsc);
203: ISDestroy(lu->iden_dsc);
204: VecScatterDestroy(lu->scat);
205:
206: if (lu->size >1) ISDestroy(lu->iden);
207: }
208: if (lu->size == 1) {
209: MatConvert_DSCPACK_Base(A,MATSEQBAIJ,&A);
210: } else {
211: MatConvert_DSCPACK_Base(A,MATMPIBAIJ,&A);
212: }
213: (*A->ops->destroy)(A);
214: return(0);
215: }
219: int MatSolve_DSCPACK(Mat A,Vec b,Vec x) {
220: Mat_DSC *lu= (Mat_DSC*)A->spptr;
221: int ierr;
222: RealNumberType *solution_vec,*rhs_vec;
225: /* scatter b into seq vec_dsc */
226: if ( !lu->scat ) {
227: VecScatterCreate(b,lu->my_cols,lu->vec_dsc,lu->iden_dsc,&lu->scat);
228: }
229: VecScatterBegin(b,lu->vec_dsc,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
230: VecScatterEnd(b,lu->vec_dsc,INSERT_VALUES,SCATTER_FORWARD,lu->scat);
232: if (lu->dsc_id != -1){
233: VecGetArray(lu->vec_dsc,&rhs_vec);
234: DSC_InputRhsLocalVec(lu->My_DSC_Solver, rhs_vec, lu->num_local_cols);
235: VecRestoreArray(lu->vec_dsc,&rhs_vec);
236:
237: DSC_Solve(lu->My_DSC_Solver);
238: if (ierr != DSC_NO_ERROR) {
239: DSC_ErrorDisplay(lu->My_DSC_Solver);
240: SETERRQ(1,"Error in calling DSC_Solve");
241: }
243: /* get the permuted local solution */
244: VecGetArray(lu->vec_dsc,&solution_vec);
245: DSC_GetLocalSolution(lu->My_DSC_Solver,solution_vec, lu->num_local_cols);
246: VecRestoreArray(lu->vec_dsc,&solution_vec);
248: } /* end of if (lu->dsc_id != -1) */
250: /* put permuted local solution solution_vec into x in the original order */
251: VecScatterBegin(lu->vec_dsc,x,INSERT_VALUES,SCATTER_REVERSE,lu->scat);
252: VecScatterEnd(lu->vec_dsc,x,INSERT_VALUES,SCATTER_REVERSE,lu->scat);
254: return(0);
255: }
259: int MatCholeskyFactorNumeric_DSCPACK(Mat A,Mat *F) {
260: Mat_SeqBAIJ *a_seq;
261: Mat_DSC *lu=(Mat_DSC*)(*F)->spptr;
262: Mat *tseq,A_seq;
263: RealNumberType *my_a_nonz;
264: int ierr, M=A->M, Mbs=M/lu->bs, size,
265: max_mem_estimate, max_single_malloc_blk,
266: number_of_procs,i,j,next,iold,
267: *idx,*iidx,*itmp;
268: IS my_cols_sorted;
269:
271: MPI_Comm_size(A->comm,&size);
272:
273: if ( lu->flg == DIFFERENT_NONZERO_PATTERN){ /* first numeric factorization */
275: /* convert A to A_seq */
276: if (size > 1) {
277: ISCreateStride(PETSC_COMM_SELF,M,0,1,&lu->iden);
278: MatGetSubMatrices(A,1,&lu->iden,&lu->iden,MAT_INITIAL_MATRIX,&tseq);
279:
280: A_seq = *tseq;
281: PetscFree(tseq);
282: a_seq = (Mat_SeqBAIJ*)A_seq->data;
283: } else {
284: a_seq = (Mat_SeqBAIJ*)A->data;
285: }
286:
287: PetscMalloc(Mbs*sizeof(int),&lu->replication);
288: for (i=0; i<Mbs; i++) lu->replication[i] = lu->bs;
290: number_of_procs = DSC_Analyze(Mbs, a_seq->i, a_seq->j, lu->replication);
291:
292: i = size;
293: if ( number_of_procs < i ) i = number_of_procs;
294: number_of_procs = 1;
295: while ( i > 1 ){
296: number_of_procs *= 2; i /= 2;
297: }
299: /* DSC_Solver starts */
300: DSC_Open0( lu->My_DSC_Solver, number_of_procs, &lu->dsc_id, lu->comm_dsc );
302: if (lu->dsc_id != -1) {
303: DSC_Order(lu->My_DSC_Solver,lu->order_code,Mbs,a_seq->i,a_seq->j,lu->replication,
304: &M,&lu->num_local_strucs,
305: &lu->num_local_cols, &lu->num_local_nonz, &lu->global_struc_new_col_num,
306: &lu->global_struc_new_num, &lu->global_struc_owner,
307: &lu->local_struc_old_num);
308: if (ierr != DSC_NO_ERROR) {
309: DSC_ErrorDisplay(lu->My_DSC_Solver);
310: SETERRQ(1,"Error when use DSC_Order()");
311: }
313: DSC_SFactor(lu->My_DSC_Solver,&max_mem_estimate,&max_single_malloc_blk,
314: lu->max_mem_allowed, lu->LBLASLevel, lu->DBLASLevel);
315: if (ierr != DSC_NO_ERROR) {
316: DSC_ErrorDisplay(lu->My_DSC_Solver);
317: SETERRQ(1,"Error when use DSC_Order");
318: }
320: BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
321: lu->num_local_strucs, lu->num_local_nonz,
322: lu->global_struc_new_col_num,
323: lu->local_struc_old_num,
324: PETSC_NULL,
325: &my_a_nonz);
326: if (ierr <0) {
327: DSC_ErrorDisplay(lu->My_DSC_Solver);
328: SETERRQ1(1,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
329: }
331: /* get local_cols_old_num and IS my_cols to be used later */
332: PetscMalloc(lu->num_local_cols*sizeof(int),&lu->local_cols_old_num);
333: for (next = 0, i=0; i<lu->num_local_strucs; i++){
334: iold = lu->bs*lu->local_struc_old_num[i];
335: for (j=0; j<lu->bs; j++)
336: lu->local_cols_old_num[next++] = iold++;
337: }
338: ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,lu->local_cols_old_num,&lu->my_cols);
339:
340: } else { /* lu->dsc_id == -1 */
341: lu->num_local_cols = 0;
342: lu->local_cols_old_num = 0;
343: ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,lu->local_cols_old_num,&lu->my_cols);
344: }
345: /* generate vec_dsc and iden_dsc to be used later */
346: VecCreateSeq(PETSC_COMM_SELF,lu->num_local_cols,&lu->vec_dsc);
347: ISCreateStride(PETSC_COMM_SELF,lu->num_local_cols,0,1,&lu->iden_dsc);
348: lu->scat = PETSC_NULL;
350: if ( size>1 ) {MatDestroy(A_seq); }
352: } else { /* use previously computed symbolic factor */
353: /* convert A to my A_seq */
354: if (size > 1) {
355: if (lu->dsc_id == -1) {
356: itmp = 0;
357: } else {
358: PetscMalloc(2*lu->num_local_strucs*sizeof(int),&idx);
359: iidx = idx + lu->num_local_strucs;
360: PetscMalloc(lu->num_local_cols*sizeof(int),&itmp);
361:
362: isort2(lu->num_local_strucs, lu->local_struc_old_num, idx);
363: for (next=0, i=0; i< lu->num_local_strucs; i++) {
364: iold = lu->bs*lu->local_struc_old_num[idx[i]];
365: for (j=0; j<lu->bs; j++){
366: itmp[next++] = iold++; /* sorted local_cols_old_num */
367: }
368: }
369: for (i=0; i< lu->num_local_strucs; i++) {
370: iidx[idx[i]] = i; /* inverse of idx */
371: }
372: } /* end of (lu->dsc_id == -1) */
373: ISCreateGeneral(PETSC_COMM_SELF,lu->num_local_cols,itmp,&my_cols_sorted);
374: MatGetSubMatrices(A,1,&my_cols_sorted,&lu->iden,MAT_INITIAL_MATRIX,&tseq);
375: ISDestroy(my_cols_sorted);
376:
377: A_seq = *tseq;
378: PetscFree(tseq);
379:
380: if (lu->dsc_id != -1) {
381: DSC_ReFactorInitialize(lu->My_DSC_Solver);
383: a_seq = (Mat_SeqBAIJ*)A_seq->data;
384: BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
385: lu->num_local_strucs, lu->num_local_nonz,
386: lu->global_struc_new_col_num,
387: lu->local_struc_old_num,
388: iidx,
389: &my_a_nonz);
390: if (ierr <0) {
391: DSC_ErrorDisplay(lu->My_DSC_Solver);
392: SETERRQ1(1,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
393: }
394:
395: PetscFree(idx);
396: PetscFree(itmp);
397: } /* end of if(lu->dsc_id != -1) */
398: } else { /* size == 1 */
399: a_seq = (Mat_SeqBAIJ*)A->data;
400:
401: BAIJtoMyANonz(a_seq->i, a_seq->j, lu->bs, a_seq->a,
402: lu->num_local_strucs, lu->num_local_nonz,
403: lu->global_struc_new_col_num,
404: lu->local_struc_old_num,
405: PETSC_NULL,
406: &my_a_nonz);
407: if (ierr <0) {
408: DSC_ErrorDisplay(lu->My_DSC_Solver);
409: SETERRQ1(1,"Error setting local nonzeroes at processor %d \n", lu->dsc_id);
410: }
411: }
412: if ( size>1 ) {MatDestroy(A_seq); }
413: }
414:
415: if (lu->dsc_id != -1) {
416: DSC_NFactor(lu->My_DSC_Solver, lu->scheme_code, my_a_nonz, lu->factor_type, lu->LBLASLevel, lu->DBLASLevel);
417: PetscFree(my_a_nonz);
418: }
419:
420: (*F)->assembled = PETSC_TRUE;
421: lu->flg = SAME_NONZERO_PATTERN;
423: return(0);
424: }
426: /* Note the Petsc permutation r is ignored */
429: int MatCholeskyFactorSymbolic_DSCPACK(Mat A,IS r,MatFactorInfo *info,Mat *F) {
430: Mat B;
431: Mat_DSC *lu;
432: int ierr,bs;
433: PetscTruth flg;
434: char buff[32],*ftype[]={"LDLT","LLT"},*ltype[]={"LBLAS1","LBLAS2","LBLAS3"},*dtype[]={"DBLAS1","DBLAS2"};
438: /* Create the factorization matrix F */
439: MatGetBlockSize(A,&bs);
440: MatCreate(A->comm,A->m,A->n,A->M,A->N,&B);
441: MatSetType(B,MATDSCPACK);
442: MatSeqBAIJSetPreallocation(B,bs,0,PETSC_NULL);
443: MatMPIBAIJSetPreallocation(B,bs,0,PETSC_NULL,0,PETSC_NULL);
444:
445: lu = (Mat_DSC*)B->spptr;
446: lu->bs = bs;
448: B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_DSCPACK;
449: B->ops->solve = MatSolve_DSCPACK;
450: B->factor = 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(A->comm,&(lu->comm_dsc));
461: /* Get the runtime input options */
462: PetscOptionsBegin(A->comm,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",
475: ftype,2,ftype[0],buff,32,&flg);
476: while (flg) {
477: PetscStrcmp(buff,"LLT",&flg);
478: if (flg) {
479: lu->factor_type = DSC_LLT;
480: break;
481: }
482: PetscStrcmp(buff,"LDLT",&flg);
483: if (flg) {
484: lu->factor_type = DSC_LDLT;
485: break;
486: }
487: SETERRQ1(1,"Unknown factor type %s",buff);
488: }
489: PetscOptionsInt("-mat_dscpack_MaxMemAllowed","", \
490: "None",
491: lu->max_mem_allowed,&lu->max_mem_allowed,PETSC_NULL);
493: PetscOptionsInt("-mat_dscpack_stats","display stats: 0 = no display, 1 = display",
494: "None", lu->stat,&lu->stat,PETSC_NULL);
495:
496: PetscOptionsEList("-mat_dscpack_LBLAS","BLAS level used in the local phase","None",
497: ltype,3,ltype[2],buff,32,&flg);
498: while (flg) {
499: PetscStrcmp(buff,"LBLAS1",&flg);
500: if (flg) {
501: lu->LBLASLevel = DSC_LBLAS1;
502: break;
503: }
504: PetscStrcmp(buff,"LBLAS2",&flg);
505: if (flg) {
506: lu->LBLASLevel = DSC_LBLAS2;
507: break;
508: }
509: PetscStrcmp(buff,"LBLAS3",&flg);
510: if (flg) {
511: lu->LBLASLevel = DSC_LBLAS3;
512: break;
513: }
514: SETERRQ1(1,"Unknown local phase BLAS level %s",buff);
515: }
517: PetscOptionsEList("-mat_dscpack_DBLAS","BLAS level used in the distributed phase","None",
518: dtype,2,dtype[1],buff,32,&flg);
519: while (flg) {
520: PetscStrcmp(buff,"DBLAS1",&flg);
521: if (flg) {
522: lu->DBLASLevel = DSC_DBLAS1;
523: break;
524: }
525: PetscStrcmp(buff,"DBLAS2",&flg);
526: if (flg) {
527: lu->DBLASLevel = DSC_DBLAS2;
528: break;
529: }
530: SETERRQ1(1,"Unknown distributed phase BLAS level %s",buff);
531: }
533: PetscOptionsEnd();
534:
535: lu->flg = DIFFERENT_NONZERO_PATTERN;
537: lu->My_DSC_Solver = DSC_Begin();
538: lu->CleanUpDSCPACK = PETSC_TRUE;
539: *F = B;
540: return(0);
541: }
545: int MatAssemblyEnd_DSCPACK(Mat A,MatAssemblyType mode) {
546: int ierr;
547: Mat_DSC *lu=(Mat_DSC*)A->spptr;
550: (*lu->MatAssemblyEnd)(A,mode);
551: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
552: A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_DSCPACK;
553: return(0);
554: }
558: int MatFactorInfo_DSCPACK(Mat A,PetscViewer viewer)
559: {
560: Mat_DSC *lu=(Mat_DSC*)A->spptr;
561: int ierr;
562: char *s;
563:
565: PetscViewerASCIIPrintf(viewer,"DSCPACK run parameters:\n");
567: switch (lu->order_code) {
568: case 1: s = "ND"; break;
569: case 2: s = "Hybrid with Minimum Degree"; break;
570: case 3: s = "Hybrid with Minimum Deficiency"; break;
571: }
572: PetscViewerASCIIPrintf(viewer," order_code: %s \n",s);
574: switch (lu->scheme_code) {
575: case 1: s = "standard factorization"; break;
576: case 2: s = "factorization + selective inversion"; break;
577: }
578: PetscViewerASCIIPrintf(viewer," scheme_code: %s \n",s);
580: switch (lu->stat) {
581: case 0: s = "NO"; break;
582: case 1: s = "YES"; break;
583: }
584: PetscViewerASCIIPrintf(viewer," display stats: %s \n",s);
585:
586: if ( lu->factor_type == DSC_LLT) {
587: s = "LLT";
588: } else if ( lu->factor_type == DSC_LDLT){
589: s = "LDLT";
590: } else {
591: SETERRQ(1,"Unknown factor type");
592: }
593: PetscViewerASCIIPrintf(viewer," factor type: %s \n",s);
595: if ( lu->LBLASLevel == DSC_LBLAS1) {
596: s = "BLAS1";
597: } else if ( lu->LBLASLevel == DSC_LBLAS2){
598: s = "BLAS2";
599: } else if ( lu->LBLASLevel == DSC_LBLAS3){
600: s = "BLAS3";
601: } else {
602: SETERRQ(1,"Unknown local phase BLAS level");
603: }
604: PetscViewerASCIIPrintf(viewer," local phase BLAS level: %s \n",s);
605:
606: if ( lu->DBLASLevel == DSC_DBLAS1) {
607: s = "BLAS1";
608: } else if ( lu->DBLASLevel == DSC_DBLAS2){
609: s = "BLAS2";
610: } else {
611: SETERRQ(1,"Unknown distributed phase BLAS level");
612: }
613: PetscViewerASCIIPrintf(viewer," distributed phase BLAS level: %s \n",s);
614: return(0);
615: }
619: int MatView_DSCPACK(Mat A,PetscViewer viewer) {
620: int ierr;
621: PetscTruth isascii;
622: PetscViewerFormat format;
623: Mat_DSC *lu=(Mat_DSC*)A->spptr;
626: (*lu->MatView)(A,viewer);
628: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
629: if (isascii) {
630: PetscViewerGetFormat(viewer,&format);
631: if (format == PETSC_VIEWER_ASCII_FACTOR_INFO) {
632: MatFactorInfo_DSCPACK(A,viewer);
633: }
634: }
635: return(0);
636: }
638: EXTERN_C_BEGIN
641: int MatConvert_Base_DSCPACK(Mat A,MatType type,Mat *newmat) {
642: /* This routine is only called to convert to MATDSCPACK */
643: /* from MATSEQBAIJ if A has a single process communicator */
644: /* or MATMPIBAIJ otherwise, so we will ignore 'MatType type'. */
645: int ierr;
646: MPI_Comm comm;
647: Mat B=*newmat;
648: Mat_DSC *lu;
651: if (B != A) {
652: MatDuplicate(A,MAT_COPY_VALUES,&B);
653: }
655: PetscObjectGetComm((PetscObject)A,&comm);
656: PetscNew(Mat_DSC,&lu);
658: lu->MatView = A->ops->view;
659: lu->MatAssemblyEnd = A->ops->assemblyend;
660: lu->MatCholeskyFactorSymbolic = A->ops->choleskyfactorsymbolic;
661: lu->MatDestroy = A->ops->destroy;
662: lu->CleanUpDSCPACK = PETSC_FALSE;
664: B->spptr = (void*)lu;
665: B->ops->duplicate = MatDuplicate_DSCPACK;
666: B->ops->view = MatView_DSCPACK;
667: B->ops->assemblyend = MatAssemblyEnd_DSCPACK;
668: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_DSCPACK;
669: B->ops->destroy = MatDestroy_DSCPACK;
671: MPI_Comm_size(comm,&(lu->size));
672: if (lu->size == 1) {
673: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqbaij_dscpack_C",
674: "MatConvert_Base_DSCPACK",MatConvert_Base_DSCPACK);
675: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_dscpack_seqbaij_C",
676: "MatConvert_DSCPACK_Base",MatConvert_DSCPACK_Base);
677: } else {
678: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_mpibaij_dscpack_C",
679: "MatConvert_Base_DSCPACK",MatConvert_Base_DSCPACK);
680: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_dscpack_mpibaij_C",
681: "MatConvert_DSCPACK_Base",MatConvert_DSCPACK_Base);
682: }
683: PetscObjectChangeTypeName((PetscObject)B,MATDSCPACK);
684: *newmat = B;
685: return(0);
686: }
687: EXTERN_C_END
691: int MatDuplicate_DSCPACK(Mat A, MatDuplicateOption op, Mat *M) {
694: (*A->ops->duplicate)(A,op,M);
695: MatConvert_Base_DSCPACK(*M,MATDSCPACK,M);
696: return(0);
697: }
699: /*MC
700: MATDSCPACK - MATDSCPACK = "dscpack" - A matrix type providing direct solvers (Cholesky) for sequential
701: or distributed matrices via the external package DSCPACK.
703: If DSCPACK is installed (see the manual for
704: instructions on how to declare the existence of external packages),
705: a matrix type can be constructed which invokes DSCPACK solvers.
706: After calling MatCreate(...,A), simply call MatSetType(A,MATDSCPACK).
707: This matrix type is only supported for double precision real.
709: This matrix inherits from MATSEQBAIJ if constructed with a single process communicator,
710: and from MATMPIBAIJ otherwise. As a result, for sequential matrices, MatSeqBAIJSetPreallocation is
711: supported, and similarly MatMPIBAIJSetPreallocation is supported for distributed matrices. It is
712: recommended that you call both of the above preallocation routines for simplicity. Also,
713: MatConvert can be called to perform inplace conversion to and from MATSEQBAIJ or MATMPIBAIJ
714: for sequential or distributed matrices respectively.
716: Options Database Keys:
717: + -mat_type dscpack - sets the matrix type to dscpack during a call to MatSetFromOptions()
718: . -mat_dscpack_order <1,2,3> - DSCPACK ordering, 1:ND, 2:Hybrid with Minimum Degree, 3:Hybrid with Minimum Deficiency
719: . -mat_dscpack_scheme <1,2> - factorization scheme, 1:standard factorization, 2: factorization with selective inversion
720: . -mat_dscpack_factor <LLT,LDLT> - the type of factorization to be performed.
721: . -mat_dscpack_MaxMemAllowed <n> - the maximum memory to be used during factorization
722: . -mat_dscpack_stats <0,1> - display stats of the factorization and solves during MatDestroy(), 0: no display, 1: display
723: . -mat_dscpack_LBLAS <LBLAS1,LBLAS2,LBLAS3> - BLAS level used in the local phase
724: - -mat_dscpack_DBLAS <DBLAS1,DBLAS2> - BLAS level used in the distributed phase
726: Level: beginner
728: .seealso: PCCHOLESKY
729: M*/
731: EXTERN_C_BEGIN
734: int MatCreate_DSCPACK(Mat A) {
735: int ierr,size;
738: /* Change type name before calling MatSetType to force proper construction of SeqBAIJ or MPIBAIJ */
739: /* and DSCPACK types */
740: PetscObjectChangeTypeName((PetscObject)A,MATDSCPACK);
741: MPI_Comm_size(A->comm,&size);
742: if (size == 1) {
743: MatSetType(A,MATSEQBAIJ);
744: } else {
745: MatSetType(A,MATMPIBAIJ);
746: }
747: MatConvert_Base_DSCPACK(A,MATDSCPACK,&A);
748: return(0);
749: }
750: EXTERN_C_END