slepc-3.14.1 2020-12-08
Report Typos and Errors

BVBiorthonormalizeColumn

Bi-orthonormalize a column of two BV objects.

Synopsis

#include "slepcbv.h"   
PetscErrorCode BVBiorthonormalizeColumn(BV V,BV W,PetscInt j,PetscReal *delta)
Collective on V

Input Parameters

V,W  - two basis vectors contexts
j  - index of column to be bi-orthonormalized

Output Parameters

delta  - (optional) value used for normalization

Notes

This function first bi-orthogonalizes vectors V[j],W[j] against W[0..j-1], and V[0..j-1], respectively. Then, it scales the vectors with 1/delta, so that the resulting vectors satisfy W[j]'*V[j] = 1.

See Also

BVOrthonormalizeColumn(), BVBiorthogonalizeColumn()

Location: src/sys/classes/bv/interface/bvbiorthog.c
Index of all BV routines
Table of Contents for all manual pages
Index of all manual pages