1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61
|
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="https://slepc.upv.es/documentation/current//Users/jroman/tmp/slepc-3.23.1/docs/manualpages/DS/DSTranslateHarmonic.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>DSTranslateHarmonic</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<div id="version" align=right><b>slepc-3.23.1 2025-05-01</b></div>
<div id="bugreport" align=right><a href="mailto:slepc-maint@upv.es?subject=Typo or Error in Documentation &body=Please describe the typo or error in the documentation: slepc-3.23.1 v3.23.1 /Users/jroman/tmp/slepc-3.23.1/docs/manualpages/DS/DSTranslateHarmonic.html "><small>Report Typos and Errors</small></a></div>
<H1>DSTranslateHarmonic</H1>
Computes a translation of the dense system.
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcds.h"
<A HREF="https://petsc.org/release/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../DS/DSTranslateHarmonic.html#DSTranslateHarmonic">DSTranslateHarmonic</A>(<A HREF="../DS/DS.html#DS">DS</A> ds,<A HREF="https://petsc.org/release/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> tau,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> beta,<A HREF="https://petsc.org/release/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</A> recover,<A HREF="https://petsc.org/release/manualpages/Sys/PetscScalar.html#PetscScalar">PetscScalar</A> *g,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> *gamma)
</PRE>
Logically Collective
<P>
<H3><FONT COLOR="#883300">Input Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ds </B></TD><TD> - the direct solver context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>tau </B></TD><TD> - the translation amount
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>beta </B></TD><TD> - last component of vector b
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>recover </B></TD><TD> - boolean flag to indicate whether to recover or not
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Output Parameters</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>g </B></TD><TD> - the computed vector (optional)
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>gamma </B></TD><TD> - scale factor (optional)
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
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.
<P>
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.
<P>
If the recover flag is activated, the computed translation undoes the
translation done previously. In that case, parameter tau is ignored.
<P>
<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
<A HREF="../DS/DSTranslateRKS.html#DSTranslateRKS">DSTranslateRKS</A>()
<BR><P><B></B><H3><FONT COLOR="#883300">Level</FONT></H3>developer<BR>
<H3><FONT COLOR="#883300">Location</FONT></H3>
</B><A HREF="../../../src/sys/classes/ds/interface/dsops.c.html#DSTranslateHarmonic">src/sys/classes/ds/interface/dsops.c</A>
<BR><BR><A HREF="./index.html">Index of all DS routines</A>
<BR><A HREF="../../../docs/manual.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>
|