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: int MatDestroy_SeqAIJ_Essl(Mat A)
29: {
30: Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
31: int ierr;
34: /* free the Essl datastructures */
35: PetscFree(essl->a);
36: MatDestroy_SeqAIJ(A);
37: return(0);
38: }
40: int MatSolve_SeqAIJ_Essl(Mat A,Vec b,Vec x)
41: {
42: Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
43: PetscScalar *xx;
44: int ierr,m,zero = 0;
47: VecGetLocalSize(b,&m);
48: VecCopy(b,x);
49: VecGetArray(x,&xx);
50: dgss(&zero,&A->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
51: VecRestoreArray(x,&xx);
52: return(0);
53: }
55: int MatLUFactorNumeric_SeqAIJ_Essl(Mat A,Mat *F)
56: {
57: Mat_SeqAIJ *a = (Mat_SeqAIJ*)(*F)->data;
58: Mat_SeqAIJ *aa = (Mat_SeqAIJ*)(A)->data;
59: Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)(*F)->spptr;
60: int i,ierr,one = 1;
63: /* copy matrix data into silly ESSL data structure */
64: if (!a->indexshift) {
65: for (i=0; i<A->m+1; i++) essl->ia[i] = aa->i[i] + 1;
66: for (i=0; i<aa->nz; i++) essl->ja[i] = aa->j[i] + 1;
67: } else {
68: PetscMemcpy(essl->ia,aa->i,(A->m+1)*sizeof(int));
69: PetscMemcpy(essl->ja,aa->j,(aa->nz)*sizeof(int));
70: }
71: PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));
72:
73: /* set Essl options */
74: essl->iparm[0] = 1;
75: essl->iparm[1] = 5;
76: essl->iparm[2] = 1;
77: essl->iparm[3] = 0;
78: essl->rparm[0] = 1.e-12;
79: essl->rparm[1] = A->lupivotthreshold;
81: dgsf(&one,&A->m,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,
82: essl->rparm,essl->oparm,essl->aux,&essl->naux);
84: return(0);
85: }
87: int MatLUFactorSymbolic_SeqAIJ_Essl(Mat A,IS r,IS c,MatLUInfo *info,Mat *F)
88: {
89: Mat B;
90: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b;
91: int ierr,*ridx,*cidx,i,len;
92: Mat_SeqAIJ_Essl *essl;
93: PetscReal f = 1.0;
96: if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
97: ierr = MatCreateSeqAIJ(A->comm,A->m,A->n,0,PETSC_NULL,F);
98: B = *F;
99: B->ops->solve = MatSolve_SeqAIJ_Essl;
100: B->ops->destroy = MatDestroy_SeqAIJ_Essl;
101: B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Essl;
102: B->factor = FACTOR_LU;
103:
104: ierr = PetscNew(Mat_SeqAIJ_Essl,&essl);
105: B->spptr = (void*)essl;
107: /* allocate the work arrays required by ESSL */
108: if (info) f = info->fill;
109: essl->nz = a->nz;
110: essl->lna = (int)a->nz*f;
111: essl->naux = 100 + 10*A->m;
113: /* since malloc is slow on IBM we try a single malloc */
114: len = essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar);
115: ierr = PetscMalloc(len,&essl->a);
116: essl->aux = essl->a + essl->lna;
117: essl->ia = (int*)(essl->aux + essl->naux);
118: essl->ja = essl->ia + essl->lna;
120: PetscLogObjectMemory(B,len+sizeof(Mat_SeqAIJ_Essl));
121: return(0);
122: }
124: int MatUseEssl_SeqAIJ(Mat A)
125: {
127: A->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ_Essl;
128: PetscLogInfo(0,"Using ESSL for SeqAIJ LU factorization and solves");
129: return(0);
130: }
132: #else
134: int MatUseEssl_SeqAIJ(Mat A)
135: {
137: return(0);
138: }
140: #endif