slepc-3.23.1 2025-05-01
   
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 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()
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