xref: /petsc/src/mat/impls/aij/seq/essl/essl.c (revision 2f71e70437ff1f4f201a034ed19298a7e2c6d424)
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 
22460c4903SKris Buschelman   int (*MatAssemblyEnd)(Mat,MatAssemblyType);
23460c4903SKris Buschelman   int (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
24e8aa55a4SKris Buschelman   int (*MatDestroy)(Mat);
25*2f71e704SKris Buschelman   PetscTruth CleanUpESSL;
26e8aa55a4SKris Buschelman } Mat_SeqAIJ_Essl;
27e8aa55a4SKris Buschelman 
28460c4903SKris Buschelman EXTERN_C_BEGIN
29460c4903SKris Buschelman #undef __FUNCT__
30460c4903SKris Buschelman #define __FUNCT__ "MatConvert_Essl_SeqAIJ"
31460c4903SKris Buschelman int MatConvert_Essl_SeqAIJ(Mat A,MatType type,Mat *newmat) {
32460c4903SKris Buschelman   int             ierr;
33460c4903SKris Buschelman   Mat             B=*newmat;
34460c4903SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
35460c4903SKris Buschelman 
36460c4903SKris Buschelman   PetscFunctionBegin;
37460c4903SKris Buschelman   if (B != A) {
38460c4903SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);
39460c4903SKris Buschelman   } else {
40*2f71e704SKris Buschelman     B->ops->matassemblyend   = essl->MatAssemblyEnd;
41*2f71e704SKris Buschelman     B->ops->lufactorsymbolic = essl->MatLUFactorSymbolic;
42*2f71e704SKris Buschelman     B->ops->destroy          = essl->MatDestroy;
43*2f71e704SKris Buschelman 
44460c4903SKris Buschelman     /* free the Essl datastructures */
45460c4903SKris Buschelman     ierr = PetscFree(essl);CHKERRQ(ierr);
46460c4903SKris Buschelman     ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr);
47460c4903SKris Buschelman   }
48460c4903SKris Buschelman   *newmat = B;
49460c4903SKris Buschelman   PetscFunctionReturn(0);
50460c4903SKris Buschelman }
51460c4903SKris Buschelman EXTERN_C_END
52460c4903SKris Buschelman 
53e8aa55a4SKris Buschelman #undef __FUNCT__
54e8aa55a4SKris Buschelman #define __FUNCT__ "MatDestroy_SeqAIJ_Essl"
55e8aa55a4SKris Buschelman int MatDestroy_SeqAIJ_Essl(Mat A)
56e8aa55a4SKris Buschelman {
57460c4903SKris Buschelman   int             ierr;
58e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
59e8aa55a4SKris Buschelman 
60e8aa55a4SKris Buschelman   PetscFunctionBegin;
61*2f71e704SKris Buschelman   if (essl->CleanUpESSL) {
62*2f71e704SKris Buschelman     ierr = PetscFree(essl->a);CHKERRQ(ierr);
63e8aa55a4SKris Buschelman   }
64*2f71e704SKris Buschelman   ierr = MatConvert_Essl_SeqAIJ(A,MATSEQAIJ,&A);
65*2f71e704SKris Buschelman   ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
66460c4903SKris Buschelman   PetscFunctionReturn(0);
67460c4903SKris Buschelman }
68460c4903SKris Buschelman 
69460c4903SKris Buschelman #undef __FUNCT__
70e8aa55a4SKris Buschelman #define __FUNCT__ "MatSolve_SeqAIJ_Essl"
71e8aa55a4SKris Buschelman int MatSolve_SeqAIJ_Essl(Mat A,Vec b,Vec x)
72e8aa55a4SKris Buschelman {
73e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)A->spptr;
74e8aa55a4SKris Buschelman   PetscScalar     *xx;
75e8aa55a4SKris Buschelman   int             ierr,m,zero = 0;
76e8aa55a4SKris Buschelman 
77e8aa55a4SKris Buschelman   PetscFunctionBegin;
78e8aa55a4SKris Buschelman   ierr = VecGetLocalSize(b,&m);CHKERRQ(ierr);
79e8aa55a4SKris Buschelman   ierr = VecCopy(b,x);CHKERRQ(ierr);
80e8aa55a4SKris Buschelman   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
81e8aa55a4SKris Buschelman   dgss(&zero,&A->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
82e8aa55a4SKris Buschelman   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
83e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
84e8aa55a4SKris Buschelman }
85e8aa55a4SKris Buschelman 
86e8aa55a4SKris Buschelman #undef __FUNCT__
87e8aa55a4SKris Buschelman #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_Essl"
88e8aa55a4SKris Buschelman int MatLUFactorNumeric_SeqAIJ_Essl(Mat A,Mat *F)
89e8aa55a4SKris Buschelman {
90e8aa55a4SKris Buschelman   Mat_SeqAIJ      *aa = (Mat_SeqAIJ*)(A)->data;
91e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl = (Mat_SeqAIJ_Essl*)(*F)->spptr;
92e8aa55a4SKris Buschelman   int             i,ierr,one = 1;
93e8aa55a4SKris Buschelman 
94e8aa55a4SKris Buschelman   PetscFunctionBegin;
95e8aa55a4SKris Buschelman   /* copy matrix data into silly ESSL data structure (1-based Frotran style) */
96e8aa55a4SKris Buschelman   for (i=0; i<A->m+1; i++) essl->ia[i] = aa->i[i] + 1;
97e8aa55a4SKris Buschelman   for (i=0; i<aa->nz; i++) essl->ja[i]  = aa->j[i] + 1;
98e8aa55a4SKris Buschelman 
99e8aa55a4SKris Buschelman   ierr = PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));CHKERRQ(ierr);
100e8aa55a4SKris Buschelman 
101e8aa55a4SKris Buschelman   /* set Essl options */
102e8aa55a4SKris Buschelman   essl->iparm[0] = 1;
103e8aa55a4SKris Buschelman   essl->iparm[1] = 5;
104e8aa55a4SKris Buschelman   essl->iparm[2] = 1;
105e8aa55a4SKris Buschelman   essl->iparm[3] = 0;
106e8aa55a4SKris Buschelman   essl->rparm[0] = 1.e-12;
107e8aa55a4SKris Buschelman   essl->rparm[1] = A->lupivotthreshold;
108e8aa55a4SKris Buschelman 
109e8aa55a4SKris Buschelman   dgsf(&one,&A->m,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,
110e8aa55a4SKris Buschelman                essl->rparm,essl->oparm,essl->aux,&essl->naux);
111e8aa55a4SKris Buschelman 
112e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
113e8aa55a4SKris Buschelman }
114e8aa55a4SKris Buschelman 
115e8aa55a4SKris Buschelman #undef __FUNCT__
116e8aa55a4SKris Buschelman #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ_Essl"
117e8aa55a4SKris Buschelman int MatLUFactorSymbolic_SeqAIJ_Essl(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F)
118e8aa55a4SKris Buschelman {
119e8aa55a4SKris Buschelman   Mat             B;
120eef49c96SKris Buschelman   Mat_SeqAIJ      *a = (Mat_SeqAIJ*)A->data;
121eef49c96SKris Buschelman   int             ierr,len;
122e8aa55a4SKris Buschelman   Mat_SeqAIJ_Essl *essl;
123e8aa55a4SKris Buschelman   PetscReal       f = 1.0;
124e8aa55a4SKris Buschelman 
125e8aa55a4SKris Buschelman   PetscFunctionBegin;
126e8aa55a4SKris Buschelman   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
127e8aa55a4SKris Buschelman   ierr = MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,A->m,A->n,&B);CHKERRQ(ierr);
128e8aa55a4SKris Buschelman   ierr = MatSetType(B,MATESSL);CHKERRQ(ierr);
129e8aa55a4SKris Buschelman   ierr = MatSeqAIJSetPreallocation(B,0,PETSC_NULL);CHKERRQ(ierr);
130e8aa55a4SKris Buschelman 
131e8aa55a4SKris Buschelman   B->ops->solve           = MatSolve_SeqAIJ_Essl;
132e8aa55a4SKris Buschelman   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Essl;
133e8aa55a4SKris Buschelman   B->factor               = FACTOR_LU;
134e8aa55a4SKris Buschelman 
135e8aa55a4SKris Buschelman   essl = (Mat_SeqAIJ_Essl *)(B->spptr);
136e8aa55a4SKris Buschelman 
137e8aa55a4SKris Buschelman   /* allocate the work arrays required by ESSL */
138e8aa55a4SKris Buschelman   f = info->fill;
139e8aa55a4SKris Buschelman   essl->nz   = a->nz;
140e8aa55a4SKris Buschelman   essl->lna  = (int)a->nz*f;
141e8aa55a4SKris Buschelman   essl->naux = 100 + 10*A->m;
142e8aa55a4SKris Buschelman 
143e8aa55a4SKris Buschelman   /* since malloc is slow on IBM we try a single malloc */
144e8aa55a4SKris Buschelman   len               = essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar);
145e8aa55a4SKris Buschelman   ierr              = PetscMalloc(len,&essl->a);CHKERRQ(ierr);
146e8aa55a4SKris Buschelman   essl->aux         = essl->a + essl->lna;
147e8aa55a4SKris Buschelman   essl->ia          = (int*)(essl->aux + essl->naux);
148e8aa55a4SKris Buschelman   essl->ja          = essl->ia + essl->lna;
149*2f71e704SKris Buschelman   essl->CleanUpESSL = PETSC_TRUE;
150e8aa55a4SKris Buschelman 
151e8aa55a4SKris Buschelman   PetscLogObjectMemory(B,len+sizeof(Mat_SeqAIJ_Essl));
152e8aa55a4SKris Buschelman   *F = B;
153e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
154e8aa55a4SKris Buschelman }
155e8aa55a4SKris Buschelman 
156*2f71e704SKris Buschelman #undef __FUNCT__
157*2f71e704SKris Buschelman #define __FUNCT__ "MatAssemblyEnd_SeqAIJ_Essl"
158*2f71e704SKris Buschelman int MatAssemblyEnd_SeqAIJ_Essl(Mat A,MatAssemblyType mode) {
159*2f71e704SKris Buschelman   int             ierr;
160*2f71e704SKris Buschelman   Mat_SeqAIJ_Essl *essl=(Mat_SeqAIJ_Essl*)(A->spptr);
161*2f71e704SKris Buschelman 
162*2f71e704SKris Buschelman   PetscFunctionBegin;
163*2f71e704SKris Buschelman   ierr = (*essl->MatAssemblyEnd)(A,mode);CHKERRQ(ierr);
164*2f71e704SKris Buschelman 
165*2f71e704SKris Buschelman   essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
166*2f71e704SKris Buschelman   A->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ_Essl;
167*2f71e704SKris Buschelman   PetscLogInfo(0,"Using ESSL for SeqAIJ LU factorization and solves");
168*2f71e704SKris Buschelman   PetscFunctionReturn(0);
169*2f71e704SKris Buschelman }
170*2f71e704SKris Buschelman 
171460c4903SKris Buschelman EXTERN_C_BEGIN
172e8aa55a4SKris Buschelman #undef __FUNCT__
173460c4903SKris Buschelman #define __FUNCT__ "MatConvert_SeqAIJ_Essl"
174460c4903SKris Buschelman int MatConvert_SeqAIJ_Essl(Mat A,MatType type,Mat *newmat) {
175460c4903SKris Buschelman   Mat             B=*newmat;
176460c4903SKris Buschelman   int             ierr;
177460c4903SKris Buschelman   Mat_SeqAIJ_Essl *essl;
178460c4903SKris Buschelman 
179e8aa55a4SKris Buschelman   PetscFunctionBegin;
180460c4903SKris Buschelman 
181460c4903SKris Buschelman   if (B != A) {
182460c4903SKris Buschelman     ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
183460c4903SKris Buschelman   }
184460c4903SKris Buschelman 
185460c4903SKris Buschelman   ierr                      = PetscNew(Mat_SeqAIJ_Essl,&essl);CHKERRQ(ierr);
186460c4903SKris Buschelman   essl->MatAssemblyEnd      = A->ops->assemblyend;
187460c4903SKris Buschelman   essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
188460c4903SKris Buschelman   essl->MatDestroy          = A->ops->destroy;
189*2f71e704SKris Buschelman   essl->CleanUpESSL         = PETSC_FALSE;
190460c4903SKris Buschelman 
191*2f71e704SKris Buschelman   B->spptr                  = (void *)essl;
192460c4903SKris Buschelman   B->ops->assemblyend       = MatAssemblyEnd_SeqAIJ_Essl;
193460c4903SKris Buschelman   B->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ_Essl;
194460c4903SKris Buschelman   B->ops->destroy           = MatDestroy_SeqAIJ_Essl;
195460c4903SKris Buschelman 
196460c4903SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_essl_C",
197460c4903SKris Buschelman                                            "MatConvert_SeqAIJ_Essl",MatConvert_SeqAIJ_Essl);CHKERRQ(ierr);
198460c4903SKris Buschelman   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_essl_seqaij_C",
199460c4903SKris Buschelman                                            "MatConvert_Essl_SeqAIJ",MatConvert_Essl_SeqAIJ);CHKERRQ(ierr);
200460c4903SKris Buschelman   ierr = PetscObjectChangeTypeName((PetscObject)B,type);CHKERRQ(ierr);
201460c4903SKris Buschelman   *newmat = B;
202e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
203e8aa55a4SKris Buschelman }
204460c4903SKris Buschelman EXTERN_C_END
205e8aa55a4SKris Buschelman 
206*2f71e704SKris Buschelman /*MC
207*2f71e704SKris Buschelman   MATESSL - a matrix type providing direct solvers (LU) for sequential matrices
208*2f71e704SKris Buschelman   via the external package ESSL.
209*2f71e704SKris Buschelman 
210*2f71e704SKris Buschelman   If ESSL is installed (see the manual for
211*2f71e704SKris Buschelman   instructions on how to declare the existence of external packages),
212*2f71e704SKris Buschelman   a matrix type can be constructed which invokes ESSL solvers.
213*2f71e704SKris Buschelman   After calling MatCreate(...,A), simply call MatSetType(A,MATESSL).
214*2f71e704SKris Buschelman   This matrix type is only supported for double precision real.
215*2f71e704SKris Buschelman 
216*2f71e704SKris Buschelman   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is
217*2f71e704SKris Buschelman   supported for this matrix type.
218*2f71e704SKris Buschelman 
219*2f71e704SKris Buschelman   Options Database Keys:
220*2f71e704SKris Buschelman . -mat_type essl - sets the matrix type to essl during a call to MatSetFromOptions()
221*2f71e704SKris Buschelman 
222*2f71e704SKris Buschelman    Level: beginner
223*2f71e704SKris Buschelman 
224*2f71e704SKris Buschelman .seealso: PCLU
225*2f71e704SKris Buschelman M*/
226*2f71e704SKris Buschelman 
227e8aa55a4SKris Buschelman EXTERN_C_BEGIN
228e8aa55a4SKris Buschelman #undef __FUNCT__
229e8aa55a4SKris Buschelman #define __FUNCT__ "MatCreate_SeqAIJ_Essl"
230e8aa55a4SKris Buschelman int MatCreate_SeqAIJ_Essl(Mat A) {
231e8aa55a4SKris Buschelman   int             ierr;
232e8aa55a4SKris Buschelman 
233e8aa55a4SKris Buschelman   PetscFunctionBegin;
234e8aa55a4SKris Buschelman   ierr = MatSetType(A,MATSEQAIJ);
235460c4903SKris Buschelman   ierr = MatConvert_SeqAIJ_Essl(A,MATESSL,&A);CHKERRQ(ierr);
236e8aa55a4SKris Buschelman   PetscFunctionReturn(0);
237e8aa55a4SKris Buschelman }
2384eb8e494SKris Buschelman EXTERN_C_END
239