| 12
 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
 66
 67
 68
 69
 70
 71
 72
 73
 
 | <!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>IPBOrthogonalize</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<H1>IPBOrthogonalize</H1>
B-Orthogonalize a vector with respect to two set of vectors. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcip.h" 
PetscErrorCode IPBOrthogonalize(IP ip,PetscInt nds,Vec *DS, Vec *BDS,PetscInt n,PetscBool *which,Vec *V,Vec *BV,Vec v,Vec Bv,PetscScalar *H,PetscReal *norm,PetscBool *lindep)
</PRE>
Collective on <A HREF="../IP/IP.html#IP">IP</A>
<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>BDS    </B></TD><TD> - B * DS
</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>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>BV     </B></TD><TD> - B * V
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Input/Output Parameter</FONT></H3>
<TABLE border="0" cellpadding="0" cellspacing="0">
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>v      </B></TD><TD> - (input) vector to be orthogonalized and (output) result of 
orthogonalization
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Bv     </B></TD><TD> - (input/output) B * v
</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>H      </B></TD><TD> - coefficients computed during orthogonalization with V, of size nds+n
if norm == PETSC_NULL, and nds+n+1 otherwise.
</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/ipborthog.c.html#IPBOrthogonalize">src/ip/ipborthog.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>
 |