VecDot#

Computes the vector dot product.

Synopsis#

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

Collective

Input Parameters#

  • x - first vector

  • y - second vector

Output Parameter#

  • val - the dot product

Performance Issues#

    per-processor memory bandwidth
    interprocessor latency
    work load imbalance 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#

Vectors and Parallel Data, Vec, VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDotRealPart()

Level#

intermediate

Location#

src/vec/vec/interface/rvector.c

Examples#

src/ksp/ksp/tutorials/ex49.c
src/snes/tutorials/ex15.c
src/snes/tutorials/ex69.c
src/tao/bound/tutorials/jbearing2.c
src/tao/constrained/tutorials/maros.c
src/tao/constrained/tutorials/tomographyADMM.c
src/tao/leastsquares/tutorials/tomography.c
src/tao/pde_constrained/tutorials/elliptic.c
src/tao/pde_constrained/tutorials/hyperbolic.c
src/tao/pde_constrained/tutorials/parabolic.c
src/tao/tutorials/ex3.c

Implementations#

VecDot_MPIKokkos in src/vec/vec/impls/mpi/kokkos/mpikok.kokkos.cxx
VecDot_MPIViennaCL in src/vec/vec/impls/mpi/mpiviennacl/mpiviennacl.cxx
VecDot_MPI in src/vec/vec/impls/mpi/pvec2.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_SeqViennaCL in src/vec/vec/impls/seq/seqviennacl/vecviennacl.cxx


Edit on GitLab

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