Actual source code: bvec1.c
2: /*
3: Defines the BLAS based vector operations. Code shared by parallel
4: and sequential vectors.
5: */
7: #include vecimpl.h
8: #include src/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 bn = (PetscBLASInt)xin->n, one = 1;
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: int i;
30: PetscScalar sum = 0.0;
31: for (i=0; i<xin->n; i++) {
32: sum += xa[i]*PetscConj(ya[i]);
33: }
34: *z = sum;
35: }
36: #else
37: *z = BLdot_(&bn,xa,&one,ya,&one);
38: #endif
39: VecRestoreArray(xin,&xa);
40: if (xin != yin) {VecRestoreArray(yin,&ya);}
41: PetscLogFlops(2*xin->n-1);
42: return(0);
43: }
47: PetscErrorCode VecTDot_Seq(Vec xin,Vec yin,PetscScalar *z)
48: {
49: PetscScalar *ya,*xa;
51: #if !defined(PETSC_USE_COMPLEX)
52: PetscBLASInt bn = (PetscBLASInt)xin->n, one = 1;
53: #endif
56: VecGetArray(xin,&xa);
57: if (xin != yin) {VecGetArray(yin,&ya);}
58: else ya = xa;
59: #if defined(PETSC_USE_COMPLEX)
60: /* cannot use BLAS dot for complex because compiler/linker is
61: not happy about returning a double complex */
62: int i;
63: PetscScalar sum = 0.0;
64: for (i=0; i<xin->n; i++) {
65: sum += xa[i]*ya[i];
66: }
67: *z = sum;
68: #else
69: *z = BLdot_(&bn,xa,&one,ya,&one);
70: #endif
71: VecRestoreArray(xin,&xa);
72: if (xin != yin) {VecRestoreArray(yin,&ya);}
73: PetscLogFlops(2*xin->n-1);
74: return(0);
75: }
79: PetscErrorCode VecScale_Seq(const PetscScalar *alpha,Vec xin)
80: {
81: Vec_Seq *x = (Vec_Seq*)xin->data;
82: PetscBLASInt bn = (PetscBLASInt)xin->n, one = 1;
85: BLscal_(&bn,(PetscScalar *)alpha,x->array,&one);
86: PetscLogFlops(xin->n);
87: return(0);
88: }
92: PetscErrorCode VecCopy_Seq(Vec xin,Vec yin)
93: {
94: Vec_Seq *x = (Vec_Seq *)xin->data;
95: PetscScalar *ya;
99: if (xin != yin) {
100: VecGetArray(yin,&ya);
101: PetscMemcpy(ya,x->array,xin->n*sizeof(PetscScalar));
102: VecRestoreArray(yin,&ya);
103: }
104: return(0);
105: }
109: PetscErrorCode VecSwap_Seq(Vec xin,Vec yin)
110: {
111: Vec_Seq *x = (Vec_Seq *)xin->data;
112: PetscScalar *ya;
114: PetscBLASInt bn = (PetscBLASInt)xin->n, one = 1;
117: if (xin != yin) {
118: VecGetArray(yin,&ya);
119: BLswap_(&bn,x->array,&one,ya,&one);
120: VecRestoreArray(yin,&ya);
121: }
122: return(0);
123: }
127: PetscErrorCode VecAXPY_Seq(const PetscScalar *alpha,Vec xin,Vec yin)
128: {
129: Vec_Seq *x = (Vec_Seq *)xin->data;
131: PetscBLASInt bn = (PetscBLASInt)xin->n, one = 1;
132: PetscScalar *yarray;
135: VecGetArray(yin,&yarray);
136: BLaxpy_(&bn,(PetscScalar *)alpha,x->array,&one,yarray,&one);
137: VecRestoreArray(yin,&yarray);
138: PetscLogFlops(2*xin->n);
139: return(0);
140: }
144: PetscErrorCode VecAXPBY_Seq(const PetscScalar *alpha,const PetscScalar *beta,Vec xin,Vec yin)
145: {
146: Vec_Seq *x = (Vec_Seq *)xin->data;
148: int n = xin->n,i;
149: PetscScalar *xx = x->array,*yy ,a = *alpha,b = *beta;
152: VecGetArray(yin,&yy);
153: for (i=0; i<n; i++) {
154: yy[i] = a*xx[i] + b*yy[i];
155: }
156: VecRestoreArray(yin,&yy);
158: PetscLogFlops(3*xin->n);
159: return(0);
160: }