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
8: /* #include <essl.h> This doesn't work! */
10: typedef struct {
11: int n,nz;
12: PetscScalar *a;
13: int *ia;
14: int *ja;
15: int lna;
16: int iparm[5];
17: PetscReal rparm[5];
18: PetscReal oparm[5];
19: PetscScalar *aux;
20: int naux;
22: int (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
23: int (*MatAssemblyEnd)(Mat,MatAssemblyType);
24: int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
25: int (*MatDestroy)(Mat);
26: PetscTruth CleanUpESSL;
27: } Mat_Essl;
29: EXTERN int MatDuplicate_Essl(Mat,MatDuplicateOption,Mat*);
31: EXTERN_C_BEGIN
34: int MatConvert_Essl_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
35: int ierr;
36: Mat B=*newmat;
37: Mat_Essl *essl=(Mat_Essl*)A->spptr;
38:
40: if (B != A) {
41: MatDuplicate(A,MAT_COPY_VALUES,&B);
42: }
43: B->ops->duplicate = essl->MatDuplicate;
44: B->ops->assemblyend = essl->MatAssemblyEnd;
45: B->ops->lufactorsymbolic = essl->MatLUFactorSymbolic;
46: B->ops->destroy = essl->MatDestroy;
48: /* free the Essl datastructures */
49: PetscFree(essl);
50: PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
51: *newmat = B;
52: return(0);
53: }
54: EXTERN_C_END
58: int MatDestroy_Essl(Mat A) {
59: int ierr;
60: Mat_Essl *essl=(Mat_Essl*)A->spptr;
63: if (essl->CleanUpESSL) {
64: PetscFree(essl->a);
65: }
66: MatConvert_Essl_SeqAIJ(A,MATSEQAIJ,&A);
67: (*A->ops->destroy)(A);
68: return(0);
69: }
73: int MatSolve_Essl(Mat A,Vec b,Vec x) {
74: Mat_Essl *essl = (Mat_Essl*)A->spptr;
75: PetscScalar *xx;
76: int ierr,m,zero = 0;
79: VecGetLocalSize(b,&m);
80: VecCopy(b,x);
81: VecGetArray(x,&xx);
82: dgss(&zero,&A->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
83: VecRestoreArray(x,&xx);
84: return(0);
85: }
89: int MatLUFactorNumeric_Essl(Mat A,Mat *F) {
90: Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(A)->data;
91: Mat_Essl *essl=(Mat_Essl*)(*F)->spptr;
92: int i,ierr,one = 1;
95: /* copy matrix data into silly ESSL data structure (1-based Frotran style) */
96: for (i=0; i<A->m+1; i++) essl->ia[i] = aa->i[i] + 1;
97: for (i=0; i<aa->nz; i++) essl->ja[i] = aa->j[i] + 1;
98:
99: PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));
100:
101: /* set Essl options */
102: essl->iparm[0] = 1;
103: essl->iparm[1] = 5;
104: essl->iparm[2] = 1;
105: essl->iparm[3] = 0;
106: essl->rparm[0] = 1.e-12;
107: essl->rparm[1] = A->lupivotthreshold;
109: dgsf(&one,&A->m,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,
110: essl->rparm,essl->oparm,essl->aux,&essl->naux);
112: return(0);
113: }
117: int MatLUFactorSymbolic_Essl(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
118: Mat B;
119: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
120: int ierr,len;
121: Mat_Essl *essl;
122: PetscReal f = 1.0;
125: if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
126: MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,A->m,A->n,&B);
127: MatSetType(B,A->type_name);
128: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
130: B->ops->solve = MatSolve_Essl;
131: B->ops->lufactornumeric = MatLUFactorNumeric_Essl;
132: B->factor = FACTOR_LU;
133:
134: essl = (Mat_Essl *)(B->spptr);
136: /* allocate the work arrays required by ESSL */
137: f = info->fill;
138: essl->nz = a->nz;
139: essl->lna = (int)a->nz*f;
140: essl->naux = 100 + 10*A->m;
142: /* since malloc is slow on IBM we try a single malloc */
143: len = essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar);
144: PetscMalloc(len,&essl->a);
145: essl->aux = essl->a + essl->lna;
146: essl->ia = (int*)(essl->aux + essl->naux);
147: essl->ja = essl->ia + essl->lna;
148: essl->CleanUpESSL = PETSC_TRUE;
150: PetscLogObjectMemory(B,len+sizeof(Mat_Essl));
151: *F = B;
152: return(0);
153: }
157: int MatAssemblyEnd_Essl(Mat A,MatAssemblyType mode) {
158: int ierr;
159: Mat_Essl *essl=(Mat_Essl*)(A->spptr);
162: (*essl->MatAssemblyEnd)(A,mode);
164: essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
165: A->ops->lufactorsymbolic = MatLUFactorSymbolic_Essl;
166: PetscLogInfo(0,"Using ESSL for LU factorization and solves");
167: return(0);
168: }
170: EXTERN_C_BEGIN
173: int MatConvert_SeqAIJ_Essl(Mat A,const MatType type,Mat *newmat) {
174: Mat B=*newmat;
175: int ierr;
176: Mat_Essl *essl;
180: if (B != A) {
181: MatDuplicate(A,MAT_COPY_VALUES,&B);
182: }
184: PetscNew(Mat_Essl,&essl);
185: essl->MatDuplicate = A->ops->duplicate;
186: essl->MatAssemblyEnd = A->ops->assemblyend;
187: essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
188: essl->MatDestroy = A->ops->destroy;
189: essl->CleanUpESSL = PETSC_FALSE;
191: B->spptr = (void *)essl;
192: B->ops->duplicate = MatDuplicate_Essl;
193: B->ops->assemblyend = MatAssemblyEnd_Essl;
194: B->ops->lufactorsymbolic = MatLUFactorSymbolic_Essl;
195: B->ops->destroy = MatDestroy_Essl;
197: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_essl_C",
198: "MatConvert_SeqAIJ_Essl",MatConvert_SeqAIJ_Essl);
199: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_essl_seqaij_C",
200: "MatConvert_Essl_SeqAIJ",MatConvert_Essl_SeqAIJ);
201: PetscObjectChangeTypeName((PetscObject)B,type);
202: *newmat = B;
203: return(0);
204: }
205: EXTERN_C_END
209: int MatDuplicate_Essl(Mat A, MatDuplicateOption op, Mat *M) {
210: int ierr;
211: Mat_Essl *lu = (Mat_Essl *)A->spptr;
214: (*lu->MatDuplicate)(A,op,M);
215: PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Essl));
216: return(0);
217: }
219: /*MC
220: MATESSL - MATESSL = "essl" - A matrix type providing direct solvers (LU) for sequential matrices
221: via the external package ESSL.
223: If ESSL is installed (see the manual for
224: instructions on how to declare the existence of external packages),
225: a matrix type can be constructed which invokes ESSL solvers.
226: After calling MatCreate(...,A), simply call MatSetType(A,MATESSL).
227: This matrix type is only supported for double precision real.
229: This matrix inherits from MATSEQAIJ. As a result, MatSeqAIJSetPreallocation is
230: supported for this matrix type. One can also call MatConvert for an inplace conversion to or from
231: the MATSEQAIJ type without data copy.
233: Options Database Keys:
234: . -mat_type essl - sets the matrix type to "essl" during a call to MatSetFromOptions()
236: Level: beginner
238: .seealso: PCLU
239: M*/
241: EXTERN_C_BEGIN
244: int MatCreate_Essl(Mat A) {
248: /* Change type name before calling MatSetType to force proper construction of SeqAIJ and Essl types */
249: PetscObjectChangeTypeName((PetscObject)A,MATESSL);
250: MatSetType(A,MATSEQAIJ);
251: MatConvert_SeqAIJ_Essl(A,MATESSL,&A);
252: return(0);
253: }
254: EXTERN_C_END