Actual source code: bvec1.c

  1: #define PETSCVEC_DLL
  2: /*
  3:    Defines the BLAS based vector operations. Code shared by parallel
  4:   and sequential vectors.
  5: */

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

 13: PetscErrorCode VecDot_Seq(Vec xin,Vec yin,PetscScalar *z)
 14: {
 15:   PetscScalar    *ya,*xa;
 17: #if !defined(PETSC_USE_COMPLEX)
 18:   PetscBLASInt   one = 1,bn = PetscBLASIntCast(xin->map->n);
 19: #endif

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

 49: PetscErrorCode VecTDot_Seq(Vec xin,Vec yin,PetscScalar *z)
 50: {
 51:   PetscScalar    *ya,*xa;
 53: #if !defined(PETSC_USE_COMPLEX)
 54:   PetscBLASInt    one = 1, bn = PetscBLASIntCast(xin->map->n);
 55: #endif

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

 85: PetscErrorCode VecScale_Seq(Vec xin, PetscScalar alpha)
 86: {
 87:   Vec_Seq        *x = (Vec_Seq*)xin->data;
 89:   PetscBLASInt   one = 1,bn = PetscBLASIntCast(xin->map->n);

 92:   if (alpha == 0.0) {
 93:     VecSet_Seq(xin,alpha);
 94:   } else if (alpha != 1.0) {
 95:     PetscScalar a = alpha;
 96:     BLASscal_(&bn,&a,x->array,&one);
 97:     PetscLogFlops(xin->map->n);
 98:   }
 99:   return(0);
100: }

104: PetscErrorCode VecCopy_Seq(Vec xin,Vec yin)
105: {
106:   Vec_Seq        *x = (Vec_Seq *)xin->data;
107:   PetscScalar    *ya;

111:   if (xin != yin) {
112:     VecGetArray(yin,&ya);
113:     PetscMemcpy(ya,x->array,xin->map->n*sizeof(PetscScalar));
114:     VecRestoreArray(yin,&ya);
115:   }
116:   return(0);
117: }

121: PetscErrorCode VecSwap_Seq(Vec xin,Vec yin)
122: {
123:   Vec_Seq        *x = (Vec_Seq *)xin->data;
124:   PetscScalar    *ya;
126:   PetscBLASInt   one = 1,bn = PetscBLASIntCast(xin->map->n);

129:   if (xin != yin) {
130:     VecGetArray(yin,&ya);
131:     BLASswap_(&bn,x->array,&one,ya,&one);
132:     VecRestoreArray(yin,&ya);
133:   }
134:   return(0);
135: }

139: PetscErrorCode VecAXPY_Seq(Vec yin,PetscScalar alpha,Vec xin)
140: {
141:   Vec_Seq        *y = (Vec_Seq *)yin->data;
143:   PetscScalar    *xarray;
144:   PetscBLASInt   one = 1,bn = PetscBLASIntCast(yin->map->n);

147:   /* assume that the BLAS handles alpha == 1.0 efficiently since we have no fast code for it */
148:   if (alpha != 0.0) {
149:     VecGetArray(xin,&xarray);
150:     BLASaxpy_(&bn,&alpha,xarray,&one,y->array,&one);
151:     VecRestoreArray(xin,&xarray);
152:     PetscLogFlops(2.0*yin->map->n);
153:   }
154:   return(0);
155: }

159: PetscErrorCode VecAXPBY_Seq(Vec yin,PetscScalar alpha,PetscScalar beta,Vec xin)
160: {
161:   Vec_Seq           *y = (Vec_Seq *)yin->data;
162:   PetscErrorCode    ierr;
163:   PetscInt          n = yin->map->n,i;
164:   PetscScalar       *yy = y->array,a = alpha,b = beta;
165:   const PetscScalar *xx;

168:   if (a == 0.0) {
169:     VecScale_Seq(yin,beta);
170:   } else if (b == 1.0) {
171:     VecAXPY_Seq(yin,alpha,xin);
172:   } else if (a == 1.0) {
173:     VecAYPX_Seq(yin,beta,xin);
174:   } else if (b == 0.0) {
175:     VecGetArray(xin,(PetscScalar**)&xx);
176:     for (i=0; i<n; i++) {
177:       yy[i] = a*xx[i];
178:     }
179:     VecRestoreArray(xin,(PetscScalar**)&xx);
180:     PetscLogFlops(xin->map->n);
181:   } else {
182:     VecGetArray(xin,(PetscScalar**)&xx);
183:     for (i=0; i<n; i++) {
184:       yy[i] = a*xx[i] + b*yy[i];
185:     }
186:     VecRestoreArray(xin,(PetscScalar**)&xx);
187:     PetscLogFlops(3.0*xin->map->n);
188:   }
189:   return(0);
190: }

194: PetscErrorCode VecAXPBYPCZ_Seq(Vec zin,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec xin,Vec yin)
195: {
196:   PetscErrorCode     ierr;
197:   PetscInt           n = zin->map->n,i;
198:   const PetscScalar  *yy,*xx;
199:   PetscScalar        *zz;

202:   VecGetArray(xin,(PetscScalar**)&xx);
203:   VecGetArray(yin,(PetscScalar**)&yy);
204:   VecGetArray(zin,&zz);
205:   if (alpha == 1.0) {
206:    for (i=0; i<n; i++) {
207:       zz[i] = xx[i] + beta*yy[i] + gamma*zz[i];
208:     }
209:     PetscLogFlops(4.0*n);
210:   } else if (gamma == 1.0) {
211:     for (i=0; i<n; i++) {
212:       zz[i] = alpha*xx[i] + beta*yy[i] + zz[i];
213:     }
214:     PetscLogFlops(4.0*n);
215:   } else {
216:     for (i=0; i<n; i++) {
217:       zz[i] = alpha*xx[i] + beta*yy[i] + gamma*zz[i];
218:     }
219:     PetscLogFlops(5.0*n);
220:   }
221:   VecRestoreArray(xin,(PetscScalar**)&xx);
222:   VecRestoreArray(yin,(PetscScalar**)&yy);
223:   VecRestoreArray(zin,&zz);
224:   return(0);
225: }