Actual source code: bvec1.c

  1: /*$Id: bvec1.c,v 1.43 2001/09/11 16:31:59 bsmith Exp $*/

  3: /*
  4:    Defines the BLAS based vector operations. Code shared by parallel
  5:   and sequential vectors.
  6: */

 8:  #include vecimpl.h
 9:  #include src/vec/impls/dvecimpl.h
 10:  #include petscblaslapack.h

 14: int VecDot_Seq(Vec xin,Vec yin,PetscScalar *z)
 15: {
 16:   PetscScalar *ya,*xa;
 17:   int         ierr;
 18: #if !defined(PETSC_USE_COMPLEX)
 19:   int         one = 1;
 20: #endif

 23:   VecGetArray(xin,&xa);
 24:   if (xin != yin) {VecGetArray(yin,&ya);}
 25:   else ya = xa;
 26: #if defined(PETSC_USE_COMPLEX)
 27:   /* cannot use BLAS dot for complex because compiler/linker is 
 28:      not happy about returning a double complex */
 29:   {
 30:     int         i;
 31:     PetscScalar sum = 0.0;
 32:     for (i=0; i<xin->n; i++) {
 33:       sum += xa[i]*PetscConj(ya[i]);
 34:     }
 35:     *z = sum;
 36:   }
 37: #else
 38:   *z = BLdot_(&xin->n,xa,&one,ya,&one);
 39: #endif
 40:   VecRestoreArray(xin,&xa);
 41:   if (xin != yin) {VecRestoreArray(yin,&ya);}
 42:   PetscLogFlops(2*xin->n-1);
 43:   return(0);
 44: }

 48: int VecTDot_Seq(Vec xin,Vec yin,PetscScalar *z)
 49: {
 50:   PetscScalar *ya,*xa;
 51:   int         ierr;
 52: #if !defined(PETSC_USE_COMPLEX)
 53:  int     one = 1;
 54: #endif

 57:   VecGetArray(xin,&xa);
 58:   if (xin != yin) {VecGetArray(yin,&ya);}
 59:   else ya = xa;
 60: #if defined(PETSC_USE_COMPLEX)
 61:   /* cannot use BLAS dot for complex because compiler/linker is 
 62:      not happy about returning a double complex */
 63:   int         i;
 64:   PetscScalar sum = 0.0;
 65:   for (i=0; i<xin->n; i++) {
 66:     sum += xa[i]*ya[i];
 67:   }
 68:   *z = sum;
 69: #else
 70:   *z = BLdot_(&xin->n,xa,&one,ya,&one);
 71: #endif
 72:   VecRestoreArray(xin,&xa);
 73:   if (xin != yin) {VecRestoreArray(yin,&ya);}
 74:   PetscLogFlops(2*xin->n-1);
 75:   return(0);
 76: }

 80: int VecScale_Seq(const PetscScalar *alpha,Vec xin)
 81: {
 82:   Vec_Seq *x = (Vec_Seq*)xin->data;
 83:   int     one = 1;

 86:   BLscal_(&xin->n,(PetscScalar *)alpha,x->array,&one);
 87:   PetscLogFlops(xin->n);
 88:   return(0);
 89: }

 93: int VecCopy_Seq(Vec xin,Vec yin)
 94: {
 95:   Vec_Seq     *x = (Vec_Seq *)xin->data;
 96:   PetscScalar *ya;
 97:   int         ierr;

100:   if (xin != yin) {
101:     VecGetArray(yin,&ya);
102:     PetscMemcpy(ya,x->array,xin->n*sizeof(PetscScalar));
103:     VecRestoreArray(yin,&ya);
104:   }
105:   return(0);
106: }

110: int VecSwap_Seq(Vec xin,Vec yin)
111: {
112:   Vec_Seq     *x = (Vec_Seq *)xin->data;
113:   PetscScalar *ya;
114:   int         ierr,one = 1;

117:   if (xin != yin) {
118:     VecGetArray(yin,&ya);
119:     BLswap_(&xin->n,x->array,&one,ya,&one);
120:     VecRestoreArray(yin,&ya);
121:   }
122:   return(0);
123: }

127: int VecAXPY_Seq(const PetscScalar *alpha,Vec xin,Vec yin)
128: {
129:   Vec_Seq      *x = (Vec_Seq *)xin->data;
130:   int          one = 1,ierr;
131:   PetscScalar  *yarray;

134:   VecGetArray(yin,&yarray);
135:   BLaxpy_(&xin->n,(PetscScalar *)alpha,x->array,&one,yarray,&one);
136:   VecRestoreArray(yin,&yarray);
137:   PetscLogFlops(2*xin->n);
138:   return(0);
139: }

143: int VecAXPBY_Seq(const PetscScalar *alpha,const PetscScalar *beta,Vec xin,Vec yin)
144: {
145:   Vec_Seq      *x = (Vec_Seq *)xin->data;
146:   int          n = xin->n,i,ierr;
147:   PetscScalar  *xx = x->array,*yy ,a = *alpha,b = *beta;

150:   VecGetArray(yin,&yy);
151:   for (i=0; i<n; i++) {
152:     yy[i] = a*xx[i] + b*yy[i];
153:   }
154:   VecRestoreArray(yin,&yy);

156:   PetscLogFlops(3*xin->n);
157:   return(0);
158: }