SlepcCheckOrthogonality

Checks (or prints) the level of orthogonality of a set of vectors.

Synopsis

#include "slepcsys.h" 
PetscErrorCode SlepcCheckOrthogonality(Vec *V,PetscInt nv,Vec *W,PetscInt nw,Mat B,PetscViewer viewer,PetscReal *lev)
Collective on Vec

Input parameters

V - a set of vectors
nv - number of V vectors
W - an alternative set of vectors (optional)
nw - number of W vectors
B - matrix defining the inner product (optional)
viewer - optional visualization context

Output parameter

lev - level of orthogonality (optional)

Notes

This function computes W'*V and prints the result. It is intended to check the level of bi-orthogonality of the vectors in the two sets. If W is equal to NULL then V is used, thus checking the orthogonality of the V vectors.

If matrix B is provided then the check uses the B-inner product, W'*B*V.

If lev is not NULL, it will contain the maximum entry of matrix W'*V - I (in absolute value). Otherwise, the matrix W'*V is printed.

Location: src/sys/slepcutil.c
Index of all sys routines
Table of Contents for all manual pages
Index of all manual pages

Examples

src/eps/examples/tutorials/ex12.c.html