Actual source code: dsc.c
1: /*$Id: dsc.c,v 1.7 2001/08/06 21:15:14 bsmith Exp $*/
2: /*
3: Provides an interface to the DSCPACK-S
4: */
6: #include src/mat/impls/aij/seq/aij.h
7: #if defined(PETSC_HAVE_DSCPACKS) && !defined(__cplusplus) && !defined(PETSC_USE_SINGLE)
8: EXTERN_C_BEGIN
9: #include "dscmain.h"
10: extern int Initialize_A_Nonz(int,int*,int*,real_number_type*,int,int*,int*,real_number_type**);
11: EXTERN_C_END
13: extern int MatDestroy_SeqAIJ(Mat);
15: /* golbal data for DSCPACK communcation between reordering and factorization */
16: int dsc_s_nz = 0; /* num of nonzeros in lower/upper half of the matrix */
17: int dsc_pass = 0; /* num of numeric factorizations for a single symbolic factorization */
19: #undef __FUNCT__
21: int MatDestroy_SeqAIJ_DSC_Fac(Mat A)
22: {
26: MatDestroy_SeqAIJ(A);
27: DSC_Final_Free_All(); /* free Cholesky factor and other relevant data structures */
28: /* DSC_Do_Stats(); */
29: DSC_Clean_Up(); /* terminate DSC solver */
31: return(0);
32: }
34: #undef __FUNCT__
36: int MatSolve_SeqAIJ_DSC(Mat A,Vec b,Vec x)
37: {
38: double *rhs_vec, *solution_vec;
39: int ierr;
40:
42: VecGetArray(x, &solution_vec);
43: VecGetArray(b, &rhs_vec);
44:
45: DSC_Input_Rhs(rhs_vec, A->m);
46: DSC_N_Solve();
47: if (DSC_STATUS.cont_or_stop == DSC_STOP_TYPE) goto ERROR_HANDLE;
49: DSC_Get_Solution(solution_vec);
50: if (DSC_STATUS.cont_or_stop == DSC_STOP_TYPE) goto ERROR_HANDLE;
51:
52: VecRestoreArray(x, &solution_vec);
53: VecRestoreArray(b, &rhs_vec);
55: ERROR_HANDLE:
56: if (DSC_STATUS.error_code != DSC_NO_ERROR) {
57: DSC_Error_Display();
58: SETERRQ(PETSC_ERR_ARG_SIZ, "DSC ERROR");
59: }
61: return(0);
62: }
64: #undef __FUNCT__
66: int MatCholeskyFactorNumeric_SeqAIJ_DSC(Mat A, Mat *F)
67: {
68: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data, *fac=(Mat_SeqAIJ*)(*F)->data;
69: IS iscol = fac->col,isicol = fac->icol;
70: PetscTruth flg;
71: int m,ierr;
72: int *ai = a->i, *aj = a->j;
73: int *perm, *iperm;
74: real_number_type *a_nonz = a->a, *s_a_nonz;
77: m = A->m;
78: if (dsc_pass == 0){ /* check the arguments */
79: if (m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ, "matrix must be square");
80: PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
81: if (!flg) SETERRQ(PETSC_ERR_ARG_SIZ, "matrix must be Seq_AIJ");
82: if (m != (*F)->m) SETERRQ(PETSC_ERR_ARG_SIZ, "factorization struct inconsistent");
84: } else { /* frees up Cholesky factor used by previous numeric factorization */
85: DSC_N_Fact_Free();
86: DSC_Re_Init();
87: }
89: ISGetIndices(iscol,&perm);
90: ISGetIndices(isicol,&iperm);
91: Initialize_A_Nonz(m,ai,aj,a_nonz,dsc_s_nz,perm,iperm, &s_a_nonz);
92: if (ierr <0) SETERRQ(PETSC_ERR_ARG_SIZ, "Error setting up permuted nonzero vector");
93:
94: DSC_N_Fact(s_a_nonz);
96: free ((char *) s_a_nonz);
97: dsc_pass++;
98: return(0);
99: }
101: #undef __FUNCT__
103: int MatCholeskyFactorSymbolic_SeqAIJ_DSC(Mat A,IS perm,PetscReal f,Mat *F)
104: {
105: /************************************************************************/
106: /* Input */
107: /* A - matrix to factor */
108: /* perm - row/col permutation (ignored) */
109: /* f - fill (ignored) */
110: /* */
111: /* Output */
112: /* F - matrix storing partial information for factorization */
113: /************************************************************************/
115: int ierr,m;
116: int max_mem_estimate, max_single_malloc_blk,MAX_MEM_ALLOWED=800;
117: PetscTruth flg;
118: IS iperm;
119: Mat_SeqAIJ *b;
120:
122: m = A->m;
123: if (m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ, "matrix must be square");
124: PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
125: if (!flg) SETERRQ(PETSC_ERR_ARG_SIZ, "matrix must be Seq_AIJ");
127: /* Create the factorization */
128: MatCreateSeqAIJ(A->comm, m, m, 0, PETSC_NULL, F);
129:
130: (*F)->ops->destroy = MatDestroy_SeqAIJ_DSC_Fac;
131: (*F)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_DSC;
132: (*F)->ops->solve = MatSolve_SeqAIJ_DSC;
133: (*F)->factor = FACTOR_CHOLESKY;
135: b = (Mat_SeqAIJ*)(*F)->data;
136: ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
137: (b)->col = perm;
138: (b)->icol = iperm;
140: /* Symbolic factorization */
141: DSC_S_Fact (&max_mem_estimate, &max_single_malloc_blk, MAX_MEM_ALLOWED);
142: if (DSC_STATUS.cont_or_stop == DSC_STOP_TYPE) goto ERROR_HANDLE;
144: ERROR_HANDLE:
145: if (DSC_STATUS.error_code != DSC_NO_ERROR) {
146: DSC_Error_Display();
147: SETERRQ(PETSC_ERR_ARG_SIZ, "DSC_ERROR");
148: }
150: return(0);
151: }
153: #undef __FUNCT__
155: int MatSeqAIJUseDSC(Mat A)
156: {
157: int ierr;
158: PetscTruth flg;
159:
162: if (A->m != A->n) SETERRQ(PETSC_ERR_ARG_SIZ, "matrix must be square");
164: PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);
165: if (!flg) SETERRQ(PETSC_ERR_ARG_SIZ, "matrix must be SeqAIJ");
166:
167: A->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ_DSC;
168: PetscLogInfo(0,"Using DSC for SeqAIJ Cholesky factorization and solve.");
169: return(0);
170: }
172: #else
174: #undef __FUNCT__
176: int MatSeqAIJUseDSC(Mat A)
177: {
180: PetscLogInfo(0,"DSCPACK not istalled. Not using DSC.");
181: return(0);
182: }
184: #endif