File: IPOrthogonalize.html

package info (click to toggle)
slepc 3.2-p5-1
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 11,356 kB
  • sloc: ansic: 34,162; makefile: 2,041; python: 1,411; fortran: 486; f90: 184; sh: 9
file content (65 lines) | stat: -rw-r--r-- 3,038 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
60
61
62
63
64
65
<!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>IPOrthogonalize</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">

<H1>IPOrthogonalize</H1>
Orthogonalize a vector with respect to two set of vectors. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcip.h" 
PetscErrorCode IPOrthogonalize(IP ip,PetscInt nds,Vec *DS,PetscInt n,PetscBool *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
</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 (<A HREF="../IP/IP.html#IP">IP</A>) context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>nds    </B></TD><TD> - number of columns of DS
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>DS     </B></TD><TD> - first set of vectors
</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>which  </B></TD><TD> - logical array indicating columns of V to be used
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>V      </B></TD><TD> - second set of vectors
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Input/Output Parameter</FONT></H3>
<DT><B>v      </B> - (input) vector to be orthogonalized and (output) result of 
orthogonalization
<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 with V
</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>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>lindep </B></TD><TD> - flag indicating that refinement did not improve the quality
of orthogonalization
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
This function applies an orthogonal projector to project vector v onto the
orthogonal complement of the span of the columns of DS and V.
<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/IPBiOrthogonalize.html#IPBiOrthogonalize">IPBiOrthogonalize</A>()
<BR><P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/ip/iporthog.c.html#IPOrthogonalize">src/ip/iporthog.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>