Actual source code: essl.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: /*
  3:         Provides an interface to the IBM RS6000 Essl sparse solver

  5: */
  6: #include <../src/mat/impls/aij/seq/aij.h>

  8: /* #include <essl.h> This doesn't work!  */

 10: PETSC_EXTERN void dgss(int*,int*,double*,int*,int*,int*,double*,double*,int*);
 11: PETSC_EXTERN void dgsf(int*,int*,int*,double*,int*,int*,int*,int*,double*,double*,double*,int*);

 13: typedef struct {
 14:   int         n,nz;
 15:   PetscScalar *a;
 16:   int         *ia;
 17:   int         *ja;
 18:   int         lna;
 19:   int         iparm[5];
 20:   PetscReal   rparm[5];
 21:   PetscReal   oparm[5];
 22:   PetscScalar *aux;
 23:   int         naux;

 25:   PetscBool CleanUpESSL;
 26: } Mat_Essl;

 30: PetscErrorCode MatDestroy_Essl(Mat A)
 31: {
 33:   Mat_Essl       *essl=(Mat_Essl*)A->spptr;

 36:   if (essl && essl->CleanUpESSL) {
 37:     PetscFree4(essl->a,essl->aux,essl->ia,essl->ja);
 38:   }
 39:   PetscFree(A->spptr);
 40:   MatDestroy_SeqAIJ(A);
 41:   return(0);
 42: }

 46: PetscErrorCode MatSolve_Essl(Mat A,Vec b,Vec x)
 47: {
 48:   Mat_Essl       *essl = (Mat_Essl*)A->spptr;
 49:   PetscScalar    *xx;
 51:   int            nessl,zero = 0;

 54:   PetscBLASIntCast(A->cmap->n,&nessl);
 55:   VecCopy(b,x);
 56:   VecGetArray(x,&xx);
 57:   dgss(&zero,&nessl,essl->a,essl->ia,essl->ja,&essl->lna,xx,essl->aux,&essl->naux);
 58:   VecRestoreArray(x,&xx);
 59:   return(0);
 60: }

 64: PetscErrorCode MatLUFactorNumeric_Essl(Mat F,Mat A,const MatFactorInfo *info)
 65: {
 66:   Mat_SeqAIJ     *aa  =(Mat_SeqAIJ*)(A)->data;
 67:   Mat_Essl       *essl=(Mat_Essl*)(F)->spptr;
 69:   int            nessl,i,one = 1;

 72:   PetscBLASIntCast(A->rmap->n,&nessl);
 73:   /* copy matrix data into silly ESSL data structure (1-based Frotran style) */
 74:   for (i=0; i<A->rmap->n+1; i++) essl->ia[i] = aa->i[i] + 1;
 75:   for (i=0; i<aa->nz; i++) essl->ja[i] = aa->j[i] + 1;

 77:   PetscMemcpy(essl->a,aa->a,(aa->nz)*sizeof(PetscScalar));

 79:   /* set Essl options */
 80:   essl->iparm[0] = 1;
 81:   essl->iparm[1] = 5;
 82:   essl->iparm[2] = 1;
 83:   essl->iparm[3] = 0;
 84:   essl->rparm[0] = 1.e-12;
 85:   essl->rparm[1] = 1.0;

 87:   PetscOptionsGetReal(NULL,((PetscObject)A)->prefix,"-matessl_lu_threshold",&essl->rparm[1],NULL);

 89:   dgsf(&one,&nessl,&essl->nz,essl->a,essl->ia,essl->ja,&essl->lna,essl->iparm,essl->rparm,essl->oparm,essl->aux,&essl->naux);

 91:   F->ops->solve     = MatSolve_Essl;
 92:   (F)->assembled    = PETSC_TRUE;
 93:   (F)->preallocated = PETSC_TRUE;
 94:   return(0);
 95: }




102: PetscErrorCode MatLUFactorSymbolic_Essl(Mat B,Mat A,IS r,IS c,const MatFactorInfo *info)
103: {
104:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
106:   Mat_Essl       *essl;
107:   PetscReal      f = 1.0;

110:   essl = (Mat_Essl*)(B->spptr);

112:   /* allocate the work arrays required by ESSL */
113:   f    = info->fill;
114:   PetscBLASIntCast(a->nz,&essl->nz);
115:   PetscBLASIntCast((PetscInt)(a->nz*f),&essl->lna);
116:   PetscBLASIntCast(100 + 10*A->rmap->n,&essl->naux);

118:   /* since malloc is slow on IBM we try a single malloc */
119:   PetscMalloc4(essl->lna,&essl->a,essl->naux,&essl->aux,essl->lna,&essl->ia,essl->lna,&essl->ja);

121:   essl->CleanUpESSL = PETSC_TRUE;

123:   PetscLogObjectMemory((PetscObject)B,essl->lna*(2*sizeof(int)+sizeof(PetscScalar)) + essl->naux*sizeof(PetscScalar));

125:   B->ops->lufactornumeric = MatLUFactorNumeric_Essl;
126:   return(0);
127: }

131: PetscErrorCode MatFactorGetSolverPackage_essl(Mat A,const MatSolverPackage *type)
132: {
134:   *type = MATSOLVERESSL;
135:   return(0);
136: }

138: /*MC
139:   MATSOLVERESSL - "essl" - Provides direct solvers (LU) for sequential matrices
140:                               via the external package ESSL.

142:   If ESSL is installed (see the manual for
143:   instructions on how to declare the existence of external packages),

145:   Works with MATSEQAIJ matrices

147:    Level: beginner

149: .seealso: PCLU, PCFactorSetMatSolverPackage(), MatSolverPackage
150: M*/

154: PETSC_EXTERN PetscErrorCode MatGetFactor_seqaij_essl(Mat A,MatFactorType ftype,Mat *F)
155: {
156:   Mat            B;
158:   Mat_Essl       *essl;

161:   if (A->cmap->N != A->rmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"matrix must be square");
162:   MatCreate(PetscObjectComm((PetscObject)A),&B);
163:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,A->rmap->n,A->cmap->n);
164:   MatSetType(B,((PetscObject)A)->type_name);
165:   MatSeqAIJSetPreallocation(B,0,NULL);

167:   PetscNewLog(B,&essl);

169:   B->spptr                 = essl;
170:   B->ops->lufactorsymbolic = MatLUFactorSymbolic_Essl;

172:   PetscObjectComposeFunction((PetscObject)B,"MatFactorGetSolverPackage_C",MatFactorGetSolverPackage_essl);

174:   B->factortype = MAT_FACTOR_LU;
175:   PetscFree(B->solvertype);
176:   PetscStrallocpy(MATSOLVERESSL,&B->solvertype);

178:   *F            = B;
179:   return(0);
180: }

184: PETSC_EXTERN PetscErrorCode MatSolverPackageRegister_Essl(void)
185: {
188:   MatSolverPackageRegister(MATSOLVERESSL,MATSEQAIJ,          MAT_FACTOR_LU,MatGetFactor_seqaij_essl);
189:   return(0);
190: }