VecDot

Computes the vector dot product.

Synopsis

#include "petscvec.h"   
PetscErrorCode  VecDot(Vec x,Vec y,PetscScalar *val)
Collective on Vec

Input Parameters

x, y - the vectors

Output Parameter

val - the dot product

Performance Issues

   per-processor memory bandwidth
   interprocessor latency
   work load inbalance that causes certain processes to arrive much earlier than others

Notes for Users of Complex Numbers

For complex vectors, VecDot() computes
    val = (x,y) = y^H x,
where y^H denotes the conjugate transpose of y. Note that this corresponds to the usual "mathematicians" complex inner product where the SECOND argument gets the complex conjugate. Since the BLASdot() complex conjugates the first first argument we call the BLASdot() with the arguments reversed.

Use VecTDot() for the indefinite form

    val = (x,y) = y^T x,
where y^T denotes the transpose of y.

See Also

VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDotRealPart()

Level

intermediate

Location

src/vec/vec/interface/rvector.c

Examples

src/vec/vec/tutorials/ex1.c.html
src/vec/vec/tutorials/ex18.c.html
src/vec/vec/tutorials/performance.c.html
src/vec/vec/tutorials/ex1f.F90.html
src/vec/vec/tutorials/ex1f90.F90.html
src/vec/vec/tutorials/ex20f90.F90.html
src/ksp/ksp/tutorials/ex49.c.html
src/snes/tutorials/ex15.c.html
src/tao/unconstrained/tutorials/spectraladjointassimilation.c.html
src/tao/constrained/tutorials/maros.c.html
src/tao/constrained/tutorials/tomographyADMM.c.html

Implementations

VecDot_MPIKokkos in src/vec/vec/impls/mpi/kokkos/mpikok.kokkos.cxx
VecDot_MPICUDA in src/vec/vec/impls/mpi/mpicuda/mpicuda.cu
VecDot_MPIViennaCL in src/vec/vec/impls/mpi/mpiviennacl/mpiviennacl.cxx
VecDot_MPI in src/vec/vec/impls/mpi/pbvec.c
VecDot_Nest in src/vec/vec/impls/nest/vecnest.c
VecDot_Seq in src/vec/vec/impls/seq/bvec1.c
VecDot_SeqKokkos in src/vec/vec/impls/seq/kokkos/veckok.kokkos.cxx
VecDot_SeqCUDA in src/vec/vec/impls/seq/seqcuda/veccuda2.cu
VecDot_SeqViennaCL in src/vec/vec/impls/seq/seqviennacl/vecviennacl.cxx

Index of all Vec routines
Table of Contents for all manual pages
Index of all manual pages