slepc-3.7.4 2017-05-17
Report Typos and Errors

BVOrthogonalizeColumn

Orthogonalize one of the column vectors with respect to the previous ones.

Synopsis

#include "slepcbv.h"   
PetscErrorCode BVOrthogonalizeColumn(BV bv,PetscInt j,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
Collective on BV

Input Parameters

bv  - the basis vectors context
j  - index of column to be orthogonalized

Output Parameters

H  - (optional) coefficients computed during orthogonalization
norm  - (optional) norm of the vector after being orthogonalized
lindep  - (optional) flag indicating that refinement did not improve the quality of orthogonalization

Notes

This function applies an orthogonal projector to project vector V[j] onto the orthogonal complement of the span of the columns of V[0..j-1], where V[.] are the vectors of BV. The columns V[0..j-1] are assumed to be mutually orthonormal.

Leading columns V[0..l-1] also participate in the orthogonalization.

If a non-standard inner product has been specified with BVSetMatrix(), then the vector is B-orthogonalized, using the non-standard inner product defined by matrix B. The output vector satisfies V[j]'*B*V[0..j-1] = 0.

This routine does not normalize the resulting vector.

See Also

BVSetOrthogonalization(), BVSetMatrix(), BVSetActiveColumns(), BVOrthogonalize(), BVOrthogonalizeVec()

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