petsc-3.14.6 2021-03-30
VecMin
Determines the vector component with minimum real part and its location.
Synopsis
#include "petscvec.h"
PetscErrorCode VecMin(Vec x,PetscInt *p,PetscReal *val)
Collective on Vec
Input Parameters
Output Parameter
| p | - the location of val (pass NULL if you don't want this location)
|
| val | - the minimum component
|
Notes
Returns the value PETSC_MAX_REAL and p = -1 if the vector is of length 0.
This returns the smallest index with the minumum value
See Also
VecMax()
Level
intermediate
Location
src/vec/vec/interface/rvector.c
Examples
src/vec/vec/tutorials/ex1.c.html
src/dm/tutorials/ex15.c.html
src/ksp/ksp/tutorials/ex72.c.html
src/ts/tutorials/ex9.c.html
src/ts/tutorials/ex10.c.html
src/ts/tutorials/ex17.c.html
Implementations
VecMin_MPIKokkos in src/vec/vec/impls/mpi/kokkos/mpikok.kokkos.cxx
VecMin_MPI in src/vec/vec/impls/mpi/pvec2.c
VecMin_Nest_Recursive in src/vec/vec/impls/nest/vecnest.c
VecMin_Nest in src/vec/vec/impls/nest/vecnest.c
VecMin_Node in src/vec/vec/impls/node/vecnode.c
VecMin_Seq in src/vec/vec/impls/seq/dvec2.c
VecMin_SeqKokkos in src/vec/vec/impls/seq/kokkos/veckok.kokkos.cxx
Index of all Vec routines
Table of Contents for all manual pages
Index of all manual pages