slepc-main 2024-11-09
Report Typos and Errors

BVOrthogonalize

Orthogonalize all columns (starting from the leading ones), that is, compute the QR decomposition.

Synopsis

#include "slepcbv.h"   
PetscErrorCode BVOrthogonalize(BV V,Mat R)
Collective

Input Parameters

V  - basis vectors to be orthogonalized (or B-orthogonalized), modified on output
R  - a sequential dense matrix (or NULL), on output the triangular factor of the QR decomposition

Notes

On input, matrix R must be a square sequential dense Mat, with at least as many rows and columns as the number of active columns of V. The output satisfies V0 = V*R (where V0 represent the input V) and V'*V = I (or V'*B*V = I if an inner product matrix B has been specified with BVSetMatrix()).

If V has leading columns, then they are not modified (are assumed to be already orthonormal) and the leading columns of R are not referenced. Let the decomposition be

   [ V01 V02 ] = [ V1 V2 ] [ R11 R12 ]
                           [  0  R22 ]
then V1 is left unchanged (equal to V01) as well as R11 (it should satisfy V01 = V1*R11).

Can pass NULL if R is not required.

The method to be used for block orthogonalization can be set with BVSetOrthogonalization(). If set to GS, the computation is done column by column with successive calls to BVOrthogonalizeColumn(). Note that in the SVQB method the R factor is not upper triangular.

If V is rank-deficient or very ill-conditioned, that is, one or more columns are (almost) linearly dependent with respect to the rest, then the algorithm may break down or result in larger numerical error. Linearly dependent columns are essentially replaced by random directions, and the corresponding diagonal entry in R is set to (nearly) zero.

See Also

BVOrthogonalizeColumn(), BVOrthogonalizeVec(), BVSetMatrix(), BVSetActiveColumns(), BVSetOrthogonalization(), BVOrthogBlockType

Level

intermediate

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