petsc-3.7.7 2017-09-25
Report Typos and Errors

VecStrideMin

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  VecStrideMin(Vec v,PetscInt start,PetscInt *idex,PetscReal *nrm)
Collective on Vec

Input Parameter

v - the vector
start - starting point of the subvector (defined by a stride)

Output Parameter

idex - the location where the minimum occurred. (pass NULL if not required)
nrm - the minimum value in the subvector

Notes

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

If xa is the array representing the vector x, then this computes the min of the array (xa[start],xa[start+stride],xa[start+2*stride], ....)

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