petsc-3.4.5 2014-06-29

VecStrideMinAll

Computes the minimum of subvector of a vector defined by a starting point and a stride and optionally its location.

Synopsis

#include "petscvec.h" 
PetscErrorCode  VecStrideMinAll(Vec v,PetscInt idex[],PetscReal nrm[])
Collective on Vec

Input Parameter

v -the vector

Output Parameter

idex - the location where the minimum occurred (not supported, pass NULL, if you need this, send mail to [email protected] to request it)
nrm - the minimums

Notes

One must call VecSetBlockSize() before this routine to set the stride information, or use a vector created from a multicomponent DMDA.

This is useful for computing, say the minimum of the pressure variable when the pressure is stored (interlaced) with other variables, e.g., density, etc. This will only work if the desire subvector is a stride subvector.

See Also

VecMin(), VecStrideNorm(), VecStrideGather(), VecStrideScatter(), VecStrideMax()

Level:advanced
Location:
src/vec/vec/utils/vinv.c
Index of all Vec routines
Table of Contents for all manual pages
Index of all manual pages