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

DSTranslateHarmonic

Computes a translation of the dense system.

Synopsis

#include "slepcds.h" 
PetscErrorCode DSTranslateHarmonic(DS ds,PetscScalar tau,PetscReal beta,PetscBool recover,PetscScalar *g,PetscReal *gamma)
Logically Collective

Input Parameters

ds  - the direct solver context
tau  - the translation amount
beta  - last component of vector b
recover  - boolean flag to indicate whether to recover or not

Output Parameters

g  - the computed vector (optional)
gamma  - scale factor (optional)

Notes

This function is intended for use in the context of Krylov methods only. It computes a translation of a Krylov decomposition in order to extract eigenpair approximations by harmonic Rayleigh-Ritz. The matrix is updated as A + g*b' where g = (A-tau*eye(n))'\b and vector b is assumed to be beta*e_n^T.

The gamma factor is defined as sqrt(1+g'*g) and can be interpreted as the factor by which the residual of the Krylov decomposition is scaled.

If the recover flag is activated, the computed translation undoes the translation done previously. In that case, parameter tau is ignored.

See Also

DSTranslateRKS()

Level

developer

Location

src/sys/classes/ds/interface/dsops.c

Index of all DS routines
Table of Contents for all manual pages
Index of all manual pages