slepc-3.20.2 2024-03-15
Report Typos and Errors

BVOrthogonalizeSomeColumn

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

Synopsis

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

Input Parameters

bv  - the basis vectors context
j  - index of column to be orthogonalized
which  - logical array indicating selected columns

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 is similar to BVOrthogonalizeColumn(), but V[j] is orthogonalized only against columns V[i] having which[i]=PETSC_TRUE. The length of array which must be j at least.

The use of this operation is restricted to MGS orthogonalization type.

In the case of an indefinite inner product, the lindep parameter is not computed (set to false).

See Also

BVOrthogonalizeColumn(), BVSetOrthogonalization()

Level

advanced

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