BVBiorthonormalizeColumn#

Bi-orthonormalize a column of two BV objects.

Synopsis#

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

Collective

Input Parameters#

  • V - first basis vectors context

  • W - second basis vectors context

  • j - index of column to be bi-orthonormalized

Output Parameter#

  • 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#

BV: Basis Vectors, BVOrthonormalizeColumn(), BVBiorthogonalizeColumn()

Level#

advanced

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