xref: /petsc/src/mat/impls/aij/seq/essl/essl.c (revision 460c49030a5280f0e471e5c8a7d2f060f5488c87)
1e8aa55a4SKris Buschelman /*$Id: essl.c,v 1.49 2001/08/07 03:02:47 balay Exp $*/
2e8aa55a4SKris Buschelman 
3e8aa55a4SKris Buschelman /*
4e8aa55a4SKris Buschelman         Provides an interface to the IBM RS6000 Essl sparse solver
5e8aa55a4SKris Buschelman 
6e8aa55a4SKris Buschelman */
7e8aa55a4SKris Buschelman #include "src/mat/impls/aij/seq/aij.h"
8e8aa55a4SKris Buschelman /* #include <essl.h> This doesn't work!  */
9e8aa55a4SKris Buschelman 
10e8aa55a4SKris Buschelman typedef struct {
11e8aa55a4SKris Buschelman   int         n,nz;
12e8aa55a4SKris Buschelman   PetscScalar *a;
13e8aa55a4SKris Buschelman   int         *ia;
14e8aa55a4SKris Buschelman   int         *ja;
15e8aa55a4SKris Buschelman   int         lna;
16e8aa55a4SKris Buschelman   int         iparm[5];
17e8aa55a4SKris Buschelman   PetscReal   rparm[5];
18e8aa55a4SKris Buschelman   PetscReal   oparm[5];
19e8aa55a4SKris Buschelman   PetscScalar *aux;
20e8aa55a4SKris Buschelman   int         naux;
21e8aa55a4SKris Buschelman 
22*460c4903SKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
23*460c4903SKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
24e8aa55a4SKris Buschelman   int (*MatDestroy)(Mat);
25e8aa55a4SKris Buschelman } Mat_SeqAIJ_Essl;
26e8aa55a4SKris Buschelman 
27*460c4903SKris Buschelman EXTERN int MatLUFactorSymbolic_SeqAIJ_Essl(Mat,IS,IS,MatFactorInfo*,Mat*);
28*460c4903SKris Buschelman 
29*460c4903SKris Buschelman EXTERN_C_BEGIN
30*460c4903SKris Buschelman #undef __FUNCT__
31*460c4903SKris Buschelman #define __FUNCT__ "MatConvert_Essl_SeqAIJ"
32*460c4903SKris Buschelman int MatConvert_Essl_SeqAIJ(Mat A,MatType type,Mat *newmat) {
33*460c4903SKris Buschelman   int             ierr;
34*460c4903SKris Buschelman   Mat             B=*newmat;
35*460c4903SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
36*460c4903SKris Buschelman 
37*460c4903SKris Buschelman   PetscFunctionBegin;
38*460c4903SKris Buschelman   if (B != A) {
39*460c4903SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);
40*460c4903SKris Buschelman   } else {
41*460c4903SKris Buschelman     /* free the Essl datastructures */
42*460c4903SKris Buschelman     ierr = PetscFree(essl->a);CHKERRQ(ierr);
43*460c4903SKris Buschelman     ierr = PetscFree(essl);CHKERRQ(ierr);
44*460c4903SKris Buschelman     ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
45*460c4903SKris Buschelman   }
46*460c4903SKris Buschelman   *newmat = B;
47*460c4903SKris Buschelman   PetscFunctionReturn(0);
48*460c4903SKris Buschelman }
49*460c4903SKris Buschelman EXTERN_C_END
50*460c4903SKris Buschelman 
51e8aa55a4SKris Buschelman #undef __FUNCT__
52e8aa55a4SKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_Essl"
53e8aa55a4SKris Buschelman int MatDestroy_SeqAIJ_Essl(Mat A)
54e8aa55a4SKris Buschelman {
55*460c4903SKris Buschelman   int             ierr;
56e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
57*460c4903SKris Buschelman   int             (*destroy)(Mat)=essl->MatDestroy;
58e8aa55a4SKris Buschelman 
59e8aa55a4SKris Buschelman   PetscFunctionBegin;
60*460c4903SKris Buschelman   ierr = MatConvert_Essl_SeqAIJ(A,MATSEQAIJ,&A);
61e8aa55a4SKris Buschelman   ierr = (*destroy)(A);CHKERRQ(ierr);
62e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
63e8aa55a4SKris Buschelman }
64e8aa55a4SKris Buschelman 
65e8aa55a4SKris Buschelman #undef __FUNCT__
66*460c4903SKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_Essl"
67*460c4903SKris Buschelman int MatAssemblyEnd_SeqAIJ_Essl(Mat A,MatAssemblyType mode) {
68*460c4903SKris Buschelman   int             ierr;
69*460c4903SKris Buschelman   Mat_SeqAIJ_Essl *essl=(Mat_SeqAIJ_Essl*)(A->spptr);
70*460c4903SKris Buschelman 
71*460c4903SKris Buschelman   PetscFunctionBegin;
72*460c4903SKris Buschelman   ierr = (*essl->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
73*460c4903SKris Buschelman 
74*460c4903SKris Buschelman   essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
75*460c4903SKris Buschelman   A->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ_Essl;
76*460c4903SKris Buschelman   PetscLogInfo(0,"Using ESSL for SeqAIJ LU factorization and solves");
77*460c4903SKris Buschelman   PetscFunctionReturn(0);
78*460c4903SKris Buschelman }
79*460c4903SKris Buschelman 
80*460c4903SKris Buschelman #undef __FUNCT__
81e8aa55a4SKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_Essl"
82e8aa55a4SKris Buschelman int MatSolve_SeqAIJ_Essl(Mat A,Vec b,Vec x)
83e8aa55a4SKris Buschelman {
84e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
85e8aa55a4SKris Buschelman   PetscScalar     *xx;
86e8aa55a4SKris Buschelman   int             ierr,m,zero = 0;
87e8aa55a4SKris Buschelman 
88e8aa55a4SKris Buschelman   PetscFunctionBegin;
89e8aa55a4SKris Buschelman   ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr);
90e8aa55a4SKris Buschelman   ierr = VecCopy(b,x);CHKERRQ(ierr);
91e8aa55a4SKris Buschelman   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
92e8aa55a4SKris Buschelman   dgss(&zero,&A->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
93e8aa55a4SKris Buschelman   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
94e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
95e8aa55a4SKris Buschelman }
96e8aa55a4SKris Buschelman 
97e8aa55a4SKris Buschelman #undef __FUNCT__
98e8aa55a4SKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_Essl"
99e8aa55a4SKris Buschelman int MatLUFactorNumeric_SeqAIJ_Essl(Mat A,Mat *F)
100e8aa55a4SKris Buschelman {
101e8aa55a4SKris Buschelman   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)(*F)->data;
102e8aa55a4SKris Buschelman   Mat_SeqAIJ      *aa = (Mat_SeqAIJ*)(A)->data;
103e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)(*F)->spptr;
104e8aa55a4SKris Buschelman   int             i,ierr,one = 1;
105e8aa55a4SKris Buschelman 
106e8aa55a4SKris Buschelman   PetscFunctionBegin;
107e8aa55a4SKris Buschelman   /* copy matrix data into silly ESSL data structure (1-based Frotran style) */
108e8aa55a4SKris Buschelman   for (i=0; i<A->m+1; i++) essl->ia[i] = aa->i[i] + 1;
109e8aa55a4SKris Buschelman   for (i=0; i<aa->nz; i++) essl->ja[i]  = aa->j[i] + 1;
110e8aa55a4SKris Buschelman 
111e8aa55a4SKris Buschelman   ierr = PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));CHKERRQ(ierr);
112e8aa55a4SKris Buschelman 
113e8aa55a4SKris Buschelman   /* set Essl options */
114e8aa55a4SKris Buschelman   essl->iparm[0] = 1;
115e8aa55a4SKris Buschelman   essl->iparm[1] = 5;
116e8aa55a4SKris Buschelman   essl->iparm[2] = 1;
117e8aa55a4SKris Buschelman   essl->iparm[3] = 0;
118e8aa55a4SKris Buschelman   essl->rparm[0] = 1.e-12;
119e8aa55a4SKris Buschelman   essl->rparm[1] = A->lupivotthreshold;
120e8aa55a4SKris Buschelman 
121e8aa55a4SKris Buschelman   dgsf(&one,&A->m,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,
122e8aa55a4SKris Buschelman                essl->rparm,essl->oparm,essl->aux,&essl->naux);
123e8aa55a4SKris Buschelman 
124e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
125e8aa55a4SKris Buschelman }
126e8aa55a4SKris Buschelman 
127e8aa55a4SKris Buschelman #undef __FUNCT__
128e8aa55a4SKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_Essl"
129e8aa55a4SKris Buschelman int MatLUFactorSymbolic_SeqAIJ_Essl(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
130e8aa55a4SKris Buschelman {
131e8aa55a4SKris Buschelman   Mat             B;
132e8aa55a4SKris Buschelman   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data,*b;
133e8aa55a4SKris Buschelman   int             ierr,*ridx,*cidx,i,len;
134e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl;
135e8aa55a4SKris Buschelman   PetscReal       f = 1.0;
136e8aa55a4SKris Buschelman 
137e8aa55a4SKris Buschelman   PetscFunctionBegin;
138e8aa55a4SKris Buschelman   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
139e8aa55a4SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,A->m,A->n,&B);CHKERRQ(ierr);
140e8aa55a4SKris Buschelman   ierr = MatSetType(B,MATESSL);CHKERRQ(ierr);
141e8aa55a4SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
142e8aa55a4SKris Buschelman 
143e8aa55a4SKris Buschelman   B->ops->solve           = MatSolve_SeqAIJ_Essl;
144e8aa55a4SKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Essl;
145e8aa55a4SKris Buschelman   B->factor               = FACTOR_LU;
146e8aa55a4SKris Buschelman 
147e8aa55a4SKris Buschelman   essl = (Mat_SeqAIJ_Essl *)(B->spptr);
148e8aa55a4SKris Buschelman 
149e8aa55a4SKris Buschelman   /* allocate the work arrays required by ESSL */
150e8aa55a4SKris Buschelman   f = info->fill;
151e8aa55a4SKris Buschelman   essl->nz   = a->nz;
152e8aa55a4SKris Buschelman   essl->lna  = (int)a->nz*f;
153e8aa55a4SKris Buschelman   essl->naux = 100 + 10*A->m;
154e8aa55a4SKris Buschelman 
155e8aa55a4SKris Buschelman   /* since malloc is slow on IBM we try a single malloc */
156e8aa55a4SKris Buschelman   len        = essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar);
157e8aa55a4SKris Buschelman   ierr       = PetscMalloc(len,&essl->a);CHKERRQ(ierr);
158e8aa55a4SKris Buschelman   essl->aux  = essl->a + essl->lna;
159e8aa55a4SKris Buschelman   essl->ia   = (int*)(essl->aux + essl->naux);
160e8aa55a4SKris Buschelman   essl->ja   = essl->ia + essl->lna;
161e8aa55a4SKris Buschelman 
162e8aa55a4SKris Buschelman   PetscLogObjectMemory(B,len+sizeof(Mat_SeqAIJ_Essl));
163e8aa55a4SKris Buschelman   *F = B;
164e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
165e8aa55a4SKris Buschelman }
166e8aa55a4SKris Buschelman 
167*460c4903SKris Buschelman EXTERN_C_BEGIN
168e8aa55a4SKris Buschelman #undef __FUNCT__
169*460c4903SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_Essl"
170*460c4903SKris Buschelman int MatConvert_SeqAIJ_Essl(Mat A,MatType type,Mat *newmat) {
171*460c4903SKris Buschelman   Mat             B=*newmat;
172*460c4903SKris Buschelman   int             ierr;
173*460c4903SKris Buschelman   Mat_SeqAIJ_Essl *essl;
174*460c4903SKris Buschelman 
175e8aa55a4SKris Buschelman   PetscFunctionBegin;
176*460c4903SKris Buschelman 
177*460c4903SKris Buschelman   if (B != A) {
178*460c4903SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
179*460c4903SKris Buschelman   }
180*460c4903SKris Buschelman 
181*460c4903SKris Buschelman   ierr                      = PetscNew(Mat_SeqAIJ_Essl,&essl);CHKERRQ(ierr);
182*460c4903SKris Buschelman   essl->MatAssemblyEnd      = A->ops->assemblyend;
183*460c4903SKris Buschelman   essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
184*460c4903SKris Buschelman   essl->MatDestroy          = A->ops->destroy;
185*460c4903SKris Buschelman   B->spptr                  = (void *)essl;
186*460c4903SKris Buschelman 
187*460c4903SKris Buschelman   B->ops->assemblyend       = MatAssemblyEnd_SeqAIJ_Essl;
188*460c4903SKris Buschelman   B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ_Essl;
189*460c4903SKris Buschelman   B->ops->destroy           = MatDestroy_SeqAIJ_Essl;
190*460c4903SKris Buschelman 
191*460c4903SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_essl_C",
192*460c4903SKris Buschelman                                            "MatConvert_SeqAIJ_Essl",MatConvert_SeqAIJ_Essl);CHKERRQ(ierr);
193*460c4903SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_essl_seqaij_C",
194*460c4903SKris Buschelman                                            "MatConvert_Essl_SeqAIJ",MatConvert_Essl_SeqAIJ);CHKERRQ(ierr);
195*460c4903SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
196*460c4903SKris Buschelman   *newmat = B;
197e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
198e8aa55a4SKris Buschelman }
199*460c4903SKris Buschelman EXTERN_C_END
200e8aa55a4SKris Buschelman 
201e8aa55a4SKris Buschelman EXTERN_C_BEGIN
202e8aa55a4SKris Buschelman #undef __FUNCT__
203e8aa55a4SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_Essl"
204e8aa55a4SKris Buschelman int MatCreate_SeqAIJ_Essl(Mat A) {
205e8aa55a4SKris Buschelman   int             ierr;
206e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl;
207e8aa55a4SKris Buschelman 
208e8aa55a4SKris Buschelman   PetscFunctionBegin;
209e8aa55a4SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);
210*460c4903SKris Buschelman   ierr = MatConvert_SeqAIJ_Essl(A,MATESSL,&A);CHKERRQ(ierr);
211e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
212e8aa55a4SKris Buschelman }
2134eb8e494SKris Buschelman EXTERN_C_END
214