:orphan:
# VecAXPBY
Computes `y = alpha x + beta y`. 
## Synopsis
```
#include "petscvec.h"   
PetscErrorCode VecAXPBY(Vec y, PetscScalar alpha, PetscScalar beta, Vec x)
```
Logically Collective


## Input Parameters

- ***alpha -*** first scalar
- ***beta -*** second scalar
- ***x -*** the first scaled vector
- ***y -*** the second scaled vector



## Output Parameter

- ***y -*** output vector





## Developer Note
The implementation is optimized for `alpha` and/or `beta` values of 0.0 and 1.0


## See Also
 [](ch_vectors), `Vec`, `VecAYPX()`, `VecMAXPY()`, `VecWAXPY()`, `VecAXPY()`, `VecAXPBYPCZ()`

## Level
intermediate

## Location
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/interface/rvector.c.html#VecAXPBY">src/vec/vec/interface/rvector.c</A>

## Implementations

<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/impls/nest/vecnest.c.html#VecAXPBY_Nest">VecAXPBY_Nest in src/vec/vec/impls/nest/vecnest.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/impls/seq/bvec1.c.html#VecAXPBY_Seq">VecAXPBY_Seq in src/vec/vec/impls/seq/bvec1.c</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/impls/seq/kokkos/veckok.kokkos.cxx.html#VecAXPBY_SeqKokkos">VecAXPBY_SeqKokkos in src/vec/vec/impls/seq/kokkos/veckok.kokkos.cxx</A><BR>
<A HREF="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/impls/seq/seqviennacl/vecviennacl.cxx.html#VecAXPBY_SeqViennaCL">VecAXPBY_SeqViennaCL in src/vec/vec/impls/seq/seqviennacl/vecviennacl.cxx</A><BR>


---
[Edit on GitLab](https://gitlab.com/petsc/petsc/-/edit/release/src/vec/vec/interface/rvector.c)


[Index of all Vec routines](index.md)  
[Table of Contents for all manual pages](/manualpages/index.md)  
[Index of all manual pages](/manualpages/singleindex.md)  
