:orphan:
# VecNorm
Computes the vector norm.
## Synopsis
```
#include "petscvec.h"
PetscErrorCode VecNorm(Vec x, NormType type, PetscReal *val)
```
Collective
## Input Parameters
- ***x -*** the vector
- ***type -*** the type of the norm requested
## Output Parameter
- ***val -*** the norm
## Values of NormType
- ***`NORM_1` -*** sum_i |x[i]|
- ***`NORM_2` -*** sqrt(sum_i |x[i]|^2)
- ***`NORM_INFINITY` -*** max_i |x[i]|
- ***`NORM_1_AND_2` -*** computes efficiently both `NORM_1` and `NORM_2` and stores them each in an output array
## Notes
For complex numbers `NORM_1` will return the traditional 1 norm of the 2 norm of the complex numbers; that is the 1
norm of the absolute values of the complex entries. In PETSc 3.6 and earlier releases it returned the 1 norm of
the 1 norm of the complex entries (what is returned by the BLAS routine asum()). Both are valid norms but most
people expect the former.
This routine stashes the computed norm value, repeated calls before the vector entries are changed are then rapid since the
precomputed value is immediately available. Certain vector operations such as `VecSet()` store the norms so the value is
immediately available and does not need to be explicitly computed. `VecScale()` updates any stashed norm values, thus calls after `VecScale()`
do not need to explicitly recompute the norm.
## Performance Issues
- ***per-*** processor memory bandwidth - limits the speed of the computation of local portion of the norm
- ***interprocessor latency -*** limits the accumulation of the result across ranks, .i.e. MPI_Allreduce() time
- ***number of ranks -*** the time for the result will grow with the log base 2 of the number of ranks sharing the vector
- ***work load imbalance -*** the rank with the largest number of vector entries will limit the speed up
## See Also
[](ch_vectors), `Vec`, `NormType`, `VecDot()`, `VecTDot()`, `VecDotBegin()`, `VecDotEnd()`, `VecNormAvailable()`,
`VecNormBegin()`, `VecNormEnd()`, `NormType()`
## Level
intermediate
## Location
src/vec/vec/interface/rvector.c
## Examples
src/dm/impls/stag/tutorials/ex1.c
src/dm/impls/stag/tutorials/ex2.c
src/dm/impls/stag/tutorials/ex3.c
src/dm/impls/stag/tutorials/ex4.c
src/dm/impls/swarm/tutorials/ex1.c
src/dm/impls/swarm/tutorials/ex1f90.F90
src/dm/tutorials/ex10.c
src/ksp/ksp/tutorials/bench_kspsolve.c
src/ksp/ksp/tutorials/ex1.c
src/ksp/ksp/tutorials/ex100.c
src/ksp/ksp/tutorials/ex100f.F90
## Implementations
VecNorm_MPIKokkos in src/vec/vec/impls/mpi/kokkos/mpikok.kokkos.cxx
VecNorm_MPIViennaCL in src/vec/vec/impls/mpi/mpiviennacl/mpiviennacl.cxx
VecNorm_MPI in src/vec/vec/impls/mpi/pvec2.c
VecNorm_Nest in src/vec/vec/impls/nest/vecnest.c
VecNorm_Seq in src/vec/vec/impls/seq/bvec2.c
VecNorm_SeqKokkos in src/vec/vec/impls/seq/kokkos/veckok.kokkos.cxx
VecNorm_SeqViennaCL in src/vec/vec/impls/seq/seqviennacl/vecviennacl.cxx
---
[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)