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 src/vec/vecimpl.h
9: #include src/vec/impls/dvecimpl.h
10: #include petscblaslapack.h
12: int VecDot_Seq(Vec xin,Vec yin,PetscScalar *z)
13: {
14: PetscScalar *ya,*xa;
15: int ierr;
16: #if !defined(PETSC_USE_COMPLEX)
17: int one = 1;
18: #endif
21: VecGetArrayFast(xin,&xa);
22: if (xin != yin) {VecGetArrayFast(yin,&ya);}
23: else ya = xa;
24: #if defined(PETSC_USE_COMPLEX)
25: /* cannot use BLAS dot for complex because compiler/linker is
26: not happy about returning a double complex */
27: {
28: int i;
29: PetscScalar sum = 0.0;
30: for (i=0; i<xin->n; i++) {
31: sum += xa[i]*PetscConj(ya[i]);
32: }
33: *z = sum;
34: }
35: #else
36: *z = BLdot_(&xin->n,xa,&one,ya,&one);
37: #endif
38: VecRestoreArrayFast(xin,&xa);
39: if (xin != yin) {VecRestoreArrayFast(yin,&ya);}
40: PetscLogFlops(2*xin->n-1);
41: return(0);
42: }
44: int VecTDot_Seq(Vec xin,Vec yin,PetscScalar *z)
45: {
46: PetscScalar *ya,*xa;
47: int ierr;
48: #if !defined(PETSC_USE_COMPLEX)
49: int one = 1;
50: #endif
53: VecGetArrayFast(xin,&xa);
54: if (xin != yin) {VecGetArrayFast(yin,&ya);}
55: else ya = xa;
56: #if defined(PETSC_USE_COMPLEX)
57: /* cannot use BLAS dot for complex because compiler/linker is
58: not happy about returning a double complex */
59: int i;
60: PetscScalar sum = 0.0;
61: for (i=0; i<xin->n; i++) {
62: sum += xa[i]*ya[i];
63: }
64: *z = sum;
65: #else
66: *z = BLdot_(&xin->n,xa,&one,ya,&one);
67: #endif
68: VecRestoreArrayFast(xin,&xa);
69: if (xin != yin) {VecRestoreArrayFast(yin,&ya);}
70: PetscLogFlops(2*xin->n-1);
71: return(0);
72: }
74: int VecScale_Seq(const PetscScalar *alpha,Vec xin)
75: {
76: Vec_Seq *x = (Vec_Seq*)xin->data;
77: int one = 1;
80: BLscal_(&xin->n,(PetscScalar *)alpha,x->array,&one);
81: PetscLogFlops(xin->n);
82: return(0);
83: }
85: int VecCopy_Seq(Vec xin,Vec yin)
86: {
87: Vec_Seq *x = (Vec_Seq *)xin->data;
88: PetscScalar *ya;
89: int ierr;
92: if (xin != yin) {
93: VecGetArrayFast(yin,&ya);
94: PetscMemcpy(ya,x->array,xin->n*sizeof(PetscScalar));
95: VecRestoreArrayFast(yin,&ya);
96: }
97: return(0);
98: }
100: int VecSwap_Seq(Vec xin,Vec yin)
101: {
102: Vec_Seq *x = (Vec_Seq *)xin->data;
103: PetscScalar *ya;
104: int ierr,one = 1;
107: if (xin != yin) {
108: VecGetArrayFast(yin,&ya);
109: BLswap_(&xin->n,x->array,&one,ya,&one);
110: VecRestoreArrayFast(yin,&ya);
111: }
112: return(0);
113: }
115: int VecAXPY_Seq(const PetscScalar *alpha,Vec xin,Vec yin)
116: {
117: Vec_Seq *x = (Vec_Seq *)xin->data;
118: int one = 1,ierr;
119: PetscScalar *yarray;
122: VecGetArrayFast(yin,&yarray);
123: BLaxpy_(&xin->n,(PetscScalar *)alpha,x->array,&one,yarray,&one);
124: VecRestoreArrayFast(yin,&yarray);
125: PetscLogFlops(2*xin->n);
126: return(0);
127: }
129: int VecAXPBY_Seq(const PetscScalar *alpha,const PetscScalar *beta,Vec xin,Vec yin)
130: {
131: Vec_Seq *x = (Vec_Seq *)xin->data;
132: int n = xin->n,i,ierr;
133: PetscScalar *xx = x->array,*yy ,a = *alpha,b = *beta;
136: VecGetArrayFast(yin,&yy);
137: for (i=0; i<n; i++) {
138: yy[i] = a*xx[i] + b*yy[i];
139: }
140: VecRestoreArrayFast(yin,&yy);
142: PetscLogFlops(3*xin->n);
143: return(0);
144: }