Actual source code: essl.c
1: #define PETSCMAT_DLL
3: /*
4: Provides an interface to the IBM RS6000 Essl sparse solver
6: */
7: #include ../src/mat/impls/aij/seq/aij.h
9: /* #include <essl.h> This doesn't work! */
12: void dgss(int*,int*,double*,int*,int*,int*,double*,double*,int*);
13: void dgsf(int*,int*,int*,double*,int*,int*,int*,int*,double*,double*,double*,int*);
16: typedef struct {
17: int n,nz;
18: PetscScalar *a;
19: int *ia;
20: int *ja;
21: int lna;
22: int iparm[5];
23: PetscReal rparm[5];
24: PetscReal oparm[5];
25: PetscScalar *aux;
26: int naux;
28: PetscTruth CleanUpESSL;
29: } Mat_Essl;
33: PetscErrorCode MatDestroy_Essl(Mat A)
34: {
36: Mat_Essl *essl=(Mat_Essl*)A->spptr;
39: if (essl->CleanUpESSL) {
40: PetscFree(essl->a);
41: }
42: PetscFree(essl);
43: MatDestroy_SeqAIJ(A);
44: return(0);
45: }
49: PetscErrorCode MatSolve_Essl(Mat A,Vec b,Vec x)
50: {
51: Mat_Essl *essl = (Mat_Essl*)A->spptr;
52: PetscScalar *xx;
54: int m,zero = 0;
57: VecGetLocalSize(b,&m);
58: VecCopy(b,x);
59: VecGetArray(x,&xx);
60: dgss(&zero,&A->cmap->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
61: VecRestoreArray(x,&xx);
62: return(0);
63: }
67: PetscErrorCode MatLUFactorNumeric_Essl(Mat F,Mat A,const MatFactorInfo *info)
68: {
69: Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(A)->data;
70: Mat_Essl *essl=(Mat_Essl*)(F)->spptr;
72: int i,one = 1;
75: /* copy matrix data into silly ESSL data structure (1-based Frotran style) */
76: for (i=0; i<A->rmap->n+1; i++) essl->ia[i] = aa->i[i] + 1;
77: for (i=0; i<aa->nz; i++) essl->ja[i] = aa->j[i] + 1;
78:
79: PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));
80:
81: /* set Essl options */
82: essl->iparm[0] = 1;
83: essl->iparm[1] = 5;
84: essl->iparm[2] = 1;
85: essl->iparm[3] = 0;
86: essl->rparm[0] = 1.e-12;
87: essl->rparm[1] = 1.0;
88: PetscOptionsGetReal(((PetscObject)A)->prefix,"-matessl_lu_threshold",&essl->rparm[1],PETSC_NULL);
90: dgsf(&one,&A->rmap->n,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,essl->rparm,essl->oparm,essl->aux,&essl->naux);
92: F->ops->solve = MatSolve_Essl;
93: (F)->assembled = PETSC_TRUE;
94: (F)->preallocated = PETSC_TRUE;
95: return(0);
96: }
103: PetscErrorCode MatLUFactorSymbolic_Essl(Mat B,Mat A,IS r,IS c,const MatFactorInfo *info)
104: {
105: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
107: Mat_Essl *essl;
108: PetscReal f = 1.0;
111: essl = (Mat_Essl *)(B->spptr);
113: /* allocate the work arrays required by ESSL */
114: f = info->fill;
115: essl->nz = a->nz;
116: essl->lna = (int)a->nz*f;
117: essl->naux = 100 + 10*A->rmap->n;
119: /* since malloc is slow on IBM we try a single malloc */
120: PetscMalloc(essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar),&essl->a);
121: essl->aux = essl->a + essl->lna;
122: essl->ia = (int*)(essl->aux + essl->naux);
123: essl->ja = essl->ia + essl->lna;
124: essl->CleanUpESSL = PETSC_TRUE;
126: PetscLogObjectMemory(B,essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar));
127: B->ops->lufactornumeric = MatLUFactorNumeric_Essl;
128: return(0);
129: }
134: PetscErrorCode MatFactorGetSolverPackage_essl(Mat A,const MatSolverPackage *type)
135: {
137: *type = MAT_SOLVER_ESSL;
138: return(0);
139: }
141: /*MC
142: MAT_SOLVER_ESSL - "essl" - Provides direct solvers (LU) for sequential matrices
143: via the external package ESSL.
145: If ESSL is installed (see the manual for
146: instructions on how to declare the existence of external packages),
148: Works with MATSEQAIJ matrices
150: Level: beginner
152: .seealso: PCLU, PCFactorSetMatSolverPackage(), MatSolverPackage
153: M*/
157: PetscErrorCode MatGetFactor_seqaij_essl(Mat A,MatFactorType ftype,Mat *F)
158: {
159: Mat B;
161: Mat_Essl *essl;
164: if (A->cmap->N != A->rmap->N) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
165: MatCreate(((PetscObject)A)->comm,&B);
166: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,A->rmap->n,A->cmap->n);
167: MatSetType(B,((PetscObject)A)->type_name);
168: MatSeqAIJSetPreallocation(B,0,PETSC_NULL);
170: PetscNewLog(B,Mat_Essl,&essl);
171: B->spptr = essl;
172: B->ops->lufactorsymbolic = MatLUFactorSymbolic_Essl;
173: PetscObjectComposeFunctionDynamic((PetscObject)B,"MatFactorGetSolverPackage_C","MatFactorGetSolverPackage_essl",MatFactorGetSolverPackage_essl);
174: B->factor = MAT_FACTOR_LU;
175: *F = B;
176: return(0);
177: }