File: IPBiOrthogonalize.html

package info (click to toggle)
slepc 3.4.2.dfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 15,584 kB
  • ctags: 66,037
  • sloc: ansic: 51,275; makefile: 2,762; python: 1,577; fortran: 785; f90: 186; sh: 9
file content (58 lines) | stat: -rw-r--r-- 2,469 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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML3.2 EN">
<HTML>
<HEAD>
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc/slepc.css" type="text/css">
<TITLE>IPBiOrthogonalize</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">

<H1>IPBiOrthogonalize</H1>
Bi-orthogonalize a vector with respect to a set of vectors. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcip.h" 
PetscErrorCode IPBiOrthogonalize(IP ip,PetscInt n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
</PRE>
Collective on <A HREF="../IP/IP.html#IP">IP</A> and Vec
<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>ip </B></TD><TD> - the inner product context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>n </B></TD><TD> - number of columns of V
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>V </B></TD><TD> - set of vectors
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>W </B></TD><TD> - set of vectors
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Input/Output Parameter</FONT></H3>
<DT><B>v </B> - vector to be orthogonalized
<br>
<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>H  </B></TD><TD> - coefficients computed during orthogonalization
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>norm </B></TD><TD> - norm of the vector after being orthogonalized
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
This function applies an oblique projector to project vector v onto the
span of the columns of V along the orthogonal complement of the column
space of W.
<P>
On exit, v0 = [V v]*H, where v0 is the original vector v.
<P>
This routine does not normalize the resulting vector.
<P>

<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
 <A HREF="../IP/IPSetOrthogonalization.html#IPSetOrthogonalization">IPSetOrthogonalization</A>(), <A HREF="../IP/IPOrthogonalize.html#IPOrthogonalize">IPOrthogonalize</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/ip/ipbiorthog.c.html#IPBiOrthogonalize">src/ip/ipbiorthog.c</A>
<BR><A HREF="./index.html">Index of all IP routines</A>
<BR><A HREF="../../index.html">Table of Contents for all manual pages</A>
<BR><A HREF="../singleindex.html">Index of all manual pages</A>
</BODY></HTML>