Actual source code: essl.c
1: /*$Id: essl.c,v 1.49 2001/08/07 03:02:47 balay Exp $*/
3: /*
4: Provides an interface to the IBM RS6000 Essl sparse solver
6: */
7: #include src/mat/impls/aij/seq/aij.h
9: #if defined(PETSC_HAVE_ESSL) && !defined(__cplusplus)
10: /* #include <essl.h> This doesn't work! */
12: typedef struct {
13: int n,nz;
14: PetscScalar *a;
15: int *ia;
16: int *ja;
17: int lna;
18: int iparm[5];
19: PetscReal rparm[5];
20: PetscReal oparm[5];
21: PetscScalar *aux;
22: int naux;
23: } Mat_SeqAIJ_Essl;
26: EXTERN int MatDestroy_SeqAIJ(Mat);
28: #undef __FUNCT__
30: int MatDestroy_SeqAIJ_Essl(Mat A)
31: {
32: Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
33: int ierr;
36: /* free the Essl datastructures */
37: PetscFree(essl->a);
38: MatDestroy_SeqAIJ(A);
39: return(0);
40: }
42: #undef __FUNCT__
44: int MatSolve_SeqAIJ_Essl(Mat A,Vec b,Vec x)
45: {
46: Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
47: PetscScalar *xx;
48: int ierr,m,zero = 0;
51: VecGetLocalSize(b,&m);
52: VecCopy(b,x);
53: VecGetArray(x,&xx);
54: dgss(&zero,&A->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
55: VecRestoreArray(x,&xx);
56: return(0);
57: }
59: #undef __FUNCT__
61: int MatLUFactorNumeric_SeqAIJ_Essl(Mat A,Mat *F)
62: {
63: Mat_SeqAIJ *a = (Mat_SeqAIJ*)(*F)->data;
64: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data;
65: Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)(*F)->spptr;
66: int i,ierr,one = 1;
69: /* copy matrix data into silly ESSL data structure */
70: if (!a->indexshift) {
71: for (i=0; i<A->m+1; i++) essl->ia[i] = aa->i[i] + 1;
72: for (i=0; i<aa->nz; i++) essl->ja[i] = aa->j[i] + 1;
73: } else {
74: PetscMemcpy(essl->ia,aa->i,(A->m+1)*sizeof(int));
75: PetscMemcpy(essl->ja,aa->j,(aa->nz)*sizeof(int));
76: }
77: PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));
78:
79: /* set Essl options */
80: essl->iparm[0] = 1;
81: essl->iparm[1] = 5;
82: essl->iparm[2] = 1;
83: essl->iparm[3] = 0;
84: essl->rparm[0] = 1.e-12;
85: essl->rparm[1] = A->lupivotthreshold;
87: dgsf(&one,&A->m,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,
88: essl->rparm,essl->oparm,essl->aux,&essl->naux);
90: return(0);
91: }
93: #undef __FUNCT__
95: int MatLUFactorSymbolic_SeqAIJ_Essl(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
96: {
97: Mat B;
98: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
99: int ierr,*ridx,*cidx,i,len;
100: Mat_SeqAIJ_Essl *essl;
101: PetscReal f = 1.0;
104: if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
105: ierr = MatCreateSeqAIJ(A->comm,A->m,A->n,0,PETSC_NULL,F);
106: B = *F;
107: B->ops->solve = MatSolve_SeqAIJ_Essl;
108: B->ops->destroy = MatDestroy_SeqAIJ_Essl;
109: B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Essl;
110: B->factor = FACTOR_LU;
111:
112: ierr = PetscNew(Mat_SeqAIJ_Essl,&essl);
113: B->spptr = (void*)essl;
115: /* allocate the work arrays required by ESSL */
116: f = info->fill;
117: essl->nz = a->nz;
118: essl->lna = (int)a->nz*f;
119: essl->naux = 100 + 10*A->m;
121: /* since malloc is slow on IBM we try a single malloc */
122: len = essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar);
123: ierr = PetscMalloc(len,&essl->a);
124: essl->aux = essl->a + essl->lna;
125: essl->ia = (int*)(essl->aux + essl->naux);
126: essl->ja = essl->ia + essl->lna;
128: PetscLogObjectMemory(B,len+sizeof(Mat_SeqAIJ_Essl));
129: return(0);
130: }
132: #undef __FUNCT__
134: int MatUseEssl_SeqAIJ(Mat A)
135: {
137: A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_Essl;
138: PetscLogInfo(0,"Using ESSL for SeqAIJ LU factorization and solves");
139: return(0);
140: }
142: #else
144: #undef __FUNCT__
146: int MatUseEssl_SeqAIJ(Mat A)
147: {
149: return(0);
150: }
152: #endif