Actual source code: essl.c

  2: /* 
  3:         Provides an interface to the IBM RS6000 Essl sparse solver

  5: */
 6:  #include src/mat/impls/aij/seq/aij.h
  7: /* #include <essl.h> This doesn't work!  */

  9: typedef struct {
 10:   int         n,nz;
 11:   PetscScalar *a;
 12:   int         *ia;
 13:   int         *ja;
 14:   int         lna;
 15:   int         iparm[5];
 16:   PetscReal   rparm[5];
 17:   PetscReal   oparm[5];
 18:   PetscScalar *aux;
 19:   int         naux;

 21:   PetscErrorCode (*MatDuplicate)(Mat,MatDuplicateOption,Mat*);
 22:   PetscErrorCode (*MatAssemblyEnd)(Mat,MatAssemblyType);
 23:   PetscErrorCode (*MatLUFactorSymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
 24:   PetscErrorCode (*MatDestroy)(Mat);
 25:   PetscTruth CleanUpESSL;
 26: } Mat_Essl;

 28: EXTERN PetscErrorCode MatDuplicate_Essl(Mat,MatDuplicateOption,Mat*);

 33: PetscErrorCode MatConvert_Essl_SeqAIJ(Mat A,const MatType type,Mat *newmat) {
 35:   Mat      B=*newmat;
 36:   Mat_Essl *essl=(Mat_Essl*)A->spptr;
 37: 
 39:   if (B != A) {
 40:     MatDuplicate(A,MAT_COPY_VALUES,&B);
 41:   }
 42:   B->ops->duplicate        = essl->MatDuplicate;
 43:   B->ops->assemblyend      = essl->MatAssemblyEnd;
 44:   B->ops->lufactorsymbolic = essl->MatLUFactorSymbolic;
 45:   B->ops->destroy          = essl->MatDestroy;

 47:   /* free the Essl datastructures */
 48:   PetscFree(essl);

 50:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqaij_essl_C","",PETSC_NULL);
 51:   PetscObjectComposeFunction((PetscObject)B,"MatConvert_essl_seqaij_C","",PETSC_NULL);

 53:   PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);
 54:   *newmat = B;
 55:   return(0);
 56: }

 61: PetscErrorCode MatDestroy_Essl(Mat A)
 62: {
 64:   Mat_Essl *essl=(Mat_Essl*)A->spptr;

 67:   if (essl->CleanUpESSL) {
 68:     PetscFree(essl->a);
 69:   }
 70:   MatConvert_Essl_SeqAIJ(A,MATSEQAIJ,&A);
 71:   (*A->ops->destroy)(A);
 72:   return(0);
 73: }

 77: PetscErrorCode MatSolve_Essl(Mat A,Vec b,Vec x) {
 78:   Mat_Essl    *essl = (Mat_Essl*)A->spptr;
 79:   PetscScalar *xx;
 81:   int         m,zero = 0;

 84:   VecGetLocalSize(b,&m);
 85:   VecCopy(b,x);
 86:   VecGetArray(x,&xx);
 87:   dgss(&zero,&A->n,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
 88:   VecRestoreArray(x,&xx);
 89:   return(0);
 90: }

 94: PetscErrorCode MatLUFactorNumeric_Essl(Mat A,Mat *F) {
 95:   Mat_SeqAIJ *aa=(Mat_SeqAIJ*)(A)->data;
 96:   Mat_Essl   *essl=(Mat_Essl*)(*F)->spptr;
 98:   int        i,one = 1;

101:   /* copy matrix data into silly ESSL data structure (1-based Frotran style) */
102:   for (i=0; i<A->m+1; i++) essl->ia[i] = aa->i[i] + 1;
103:   for (i=0; i<aa->nz; i++) essl->ja[i]  = aa->j[i] + 1;
104: 
105:   PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));
106: 
107:   /* set Essl options */
108:   essl->iparm[0] = 1;
109:   essl->iparm[1] = 5;
110:   essl->iparm[2] = 1;
111:   essl->iparm[3] = 0;
112:   essl->rparm[0] = 1.e-12;
113:   essl->rparm[1] = A->lupivotthreshold;

115:   dgsf(&one,&A->m,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,
116:                essl->rparm,essl->oparm,essl->aux,&essl->naux);

118:   return(0);
119: }

123: PetscErrorCode MatLUFactorSymbolic_Essl(Mat A,IS r,IS c,MatFactorInfo *info,Mat *F) {
124:   Mat        B;
125:   Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
127:   int        len;
128:   Mat_Essl   *essl;
129:   PetscReal  f = 1.0;

132:   if (A->N != A->M) SETERRQ(PETSC_ERR_ARG_SIZ,"matrix must be square");
133:   MatCreate(A->comm,PETSC_DECIDE,PETSC_DECIDE,A->m,A->n,&B);
134:   MatSetType(B,A->type_name);
135:   MatSeqAIJSetPreallocation(B,0,PETSC_NULL);

137:   B->ops->solve           = MatSolve_Essl;
138:   B->ops->lufactornumeric = MatLUFactorNumeric_Essl;
139:   B->factor               = FACTOR_LU;
140: 
141:   essl = (Mat_Essl *)(B->spptr);

143:   /* allocate the work arrays required by ESSL */
144:   f = info->fill;
145:   essl->nz   = a->nz;
146:   essl->lna  = (int)a->nz*f;
147:   essl->naux = 100 + 10*A->m;

149:   /* since malloc is slow on IBM we try a single malloc */
150:   len               = essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar);
151:   PetscMalloc(len,&essl->a);
152:   essl->aux         = essl->a + essl->lna;
153:   essl->ia          = (int*)(essl->aux + essl->naux);
154:   essl->ja          = essl->ia + essl->lna;
155:   essl->CleanUpESSL = PETSC_TRUE;

157:   PetscLogObjectMemory(B,len+sizeof(Mat_Essl));
158:   *F = B;
159:   return(0);
160: }

164: PetscErrorCode MatAssemblyEnd_Essl(Mat A,MatAssemblyType mode) {
166:   Mat_Essl *essl=(Mat_Essl*)(A->spptr);

169:   (*essl->MatAssemblyEnd)(A,mode);

171:   essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
172:   A->ops->lufactorsymbolic  = MatLUFactorSymbolic_Essl;
173:   PetscLogInfo(0,"Using ESSL for LU factorization and solves");
174:   return(0);
175: }

180: PetscErrorCode MatConvert_SeqAIJ_Essl(Mat A,const MatType type,Mat *newmat)
181: {
182:   Mat      B=*newmat;
184:   Mat_Essl *essl;

187:   if (B != A) {
188:     MatDuplicate(A,MAT_COPY_VALUES,&B);
189:   }

191:   PetscNew(Mat_Essl,&essl);
192:   essl->MatDuplicate        = A->ops->duplicate;
193:   essl->MatAssemblyEnd      = A->ops->assemblyend;
194:   essl->MatLUFactorSymbolic = A->ops->lufactorsymbolic;
195:   essl->MatDestroy          = A->ops->destroy;
196:   essl->CleanUpESSL         = PETSC_FALSE;

198:   B->spptr                  = (void*)essl;
199:   B->ops->duplicate         = MatDuplicate_Essl;
200:   B->ops->assemblyend       = MatAssemblyEnd_Essl;
201:   B->ops->lufactorsymbolic  = MatLUFactorSymbolic_Essl;
202:   B->ops->destroy           = MatDestroy_Essl;

204:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_essl_C",
205:                                            "MatConvert_SeqAIJ_Essl",MatConvert_SeqAIJ_Essl);
206:   PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_essl_seqaij_C",
207:                                            "MatConvert_Essl_SeqAIJ",MatConvert_Essl_SeqAIJ);
208:   PetscObjectChangeTypeName((PetscObject)B,type);
209:   *newmat = B;
210:   return(0);
211: }

216: PetscErrorCode MatDuplicate_Essl(Mat A, MatDuplicateOption op, Mat *M) {
218:   Mat_Essl *lu = (Mat_Essl *)A->spptr;

221:   (*lu->MatDuplicate)(A,op,M);
222:   PetscMemcpy((*M)->spptr,lu,sizeof(Mat_Essl));
223:   return(0);
224: }

226: /*MC
227:   MATESSL - MATESSL = "essl" - A matrix type providing direct solvers (LU) for sequential matrices 
228:   via the external package ESSL.

230:   If ESSL is installed (see the manual for
231:   instructions on how to declare the existence of external packages),
232:   a matrix type can be constructed which invokes ESSL solvers.
233:   After calling MatCreate(...,A), simply call MatSetType(A,MATESSL).
234:   This matrix type is only supported for double precision real.

236:   This matrix inherits from MATSEQAIJ.  As a result, MatSeqAIJSetPreallocation is 
237:   supported for this matrix type.  One can also call MatConvert for an inplace conversion to or from 
238:   the MATSEQAIJ type without data copy.

240:   Options Database Keys:
241: . -mat_type essl - sets the matrix type to "essl" during a call to MatSetFromOptions()

243:    Level: beginner

245: .seealso: PCLU
246: M*/

251: PetscErrorCode MatCreate_Essl(Mat A)
252: {

256:   /* Change type name before calling MatSetType to force proper construction of SeqAIJ and Essl types */
257:   PetscObjectChangeTypeName((PetscObject)A,MATESSL);
258:   MatSetType(A,MATSEQAIJ);
259:   MatConvert_SeqAIJ_Essl(A,MATESSL,&A);
260:   return(0);
261: }