File: SVDComputeError.html

package info (click to toggle)
slepc 3.22.2%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 49,376 kB
  • sloc: ansic: 118,012; python: 4,887; f90: 3,620; cpp: 1,526; makefile: 811; sh: 311
file content (59 lines) | stat: -rw-r--r-- 3,543 bytes parent folder | download
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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD> <link rel="canonical" href="https://slepc.upv.es/documentation/current/docs/manualpages/SVD/SVDComputeError.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>SVDComputeError</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
   <div id="version" align=right><b>slepc-3.22.2 2024-12-02</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.22.2 v3.22.2 docs/manualpages/SVD/SVDComputeError.html "><small>Report Typos and Errors</small></a></div>
<H1>SVDComputeError</H1>
Computes the error (based on the residual norm) associated with the i-th singular triplet. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcsvd.h" 
<A HREF="https://petsc.org/release/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../SVD/SVDComputeError.html#SVDComputeError">SVDComputeError</A>(<A HREF="../SVD/SVD.html#SVD">SVD</A> svd,<A HREF="https://petsc.org/release/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> i,<A HREF="../SVD/SVDErrorType.html#SVDErrorType">SVDErrorType</A> type,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> *error)
</PRE>
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>svd  </B></TD><TD>&nbsp;- the singular value solver context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>i    </B></TD><TD>&nbsp;- the solution index
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>type </B></TD><TD>&nbsp;- the type of error to compute
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Output Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>error </B></TD><TD>&nbsp;- the error
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
The error can be computed in various ways, all of them based on the residual
norm obtained as sqrt(n1^2+n2^2) with n1 = ||A*v-sigma*u||_2 and
n2 = ||A^T*u-sigma*v||_2, where sigma is the singular value, u is the left
singular vector and v is the right singular vector.
<P>
In the case of the GSVD, the two components of the residual norm are
n1 = ||s^2*A'*u-c*B'*B*x||_2 and n2 = ||c^2*B'*v-s*A'*A*x||_2, where [u;v]
are the left singular vectors and x is the right singular vector, with
sigma=c/s.
<P>

<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
 <A HREF="../SVD/SVDErrorType.html#SVDErrorType">SVDErrorType</A>, <A HREF="../SVD/SVDSolve.html#SVDSolve">SVDSolve</A>()
<BR><P><B></B><H3><FONT COLOR="#883300">Level</FONT></H3>beginner<BR>
<H3><FONT COLOR="#883300">Location</FONT></H3>
</B><A HREF="../../../src/svd/interface/svdsolve.c.html#SVDComputeError">src/svd/interface/svdsolve.c</A>
<P><H3><FONT COLOR="#883300">Examples</FONT></H3>
<A HREF="../../../src/svd/tutorials/ex15.c.html">src/svd/tutorials/ex15.c</A><BR>
<A HREF="../../../src/svd/tutorials/ex15f.F90.html">src/svd/tutorials/ex15f.F90</A><BR>
<A HREF="../../../src/svd/tutorials/ex45.c.html">src/svd/tutorials/ex45.c</A><BR>
<BR><BR><A HREF="./index.html">Index of all SVD 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>