File: BVOrthonormalizeColumn.html

package info (click to toggle)
slepc 3.23.1%2Bdfsg1-1exp2
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 48,276 kB
  • sloc: ansic: 103,363; python: 6,078; f90: 3,286; cpp: 1,528; makefile: 772; sh: 311
file content (65 lines) | stat: -rw-r--r-- 4,435 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> <link rel="canonical" href="https://slepc.upv.es/documentation/current//Users/jroman/tmp/slepc-3.23.1/docs/manualpages/BV/BVOrthonormalizeColumn.html" />
<META NAME="GENERATOR" CONTENT="DOCTEXT">
<link rel="stylesheet" href="/slepc.css" type="text/css">
<TITLE>BVOrthonormalizeColumn</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/BV/BVOrthonormalizeColumn.html "><small>Report Typos and Errors</small></a></div>
<H1>BVOrthonormalizeColumn</H1>
Orthonormalize one of the column vectors with respect to the previous ones. 
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcbv.h"   
<A HREF="https://petsc.org/release/manualpages/Sys/PetscErrorCode.html#PetscErrorCode">PetscErrorCode</A> <A HREF="../BV/BVOrthonormalizeColumn.html#BVOrthonormalizeColumn">BVOrthonormalizeColumn</A>(<A HREF="../BV/BV.html#BV">BV</A> bv,<A HREF="https://petsc.org/release/manualpages/Sys/PetscInt.html#PetscInt">PetscInt</A> j,<A HREF="https://petsc.org/release/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</A> replace,<A HREF="https://petsc.org/release/manualpages/Sys/PetscReal.html#PetscReal">PetscReal</A> *norm,<A HREF="https://petsc.org/release/manualpages/Sys/PetscBool.html#PetscBool">PetscBool</A> *lindep)
</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>bv      </B></TD><TD>&nbsp;- the basis vectors context
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>j       </B></TD><TD>&nbsp;- index of column to be orthonormalized
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>replace </B></TD><TD>&nbsp;- whether it is allowed to set the vector randomly
</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>norm    </B></TD><TD>&nbsp;- (optional) norm of the vector after orthogonalization and before normalization
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>lindep  </B></TD><TD>&nbsp;- (optional) flag indicating that linear dependence was determined during
orthogonalization
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
This is equivalent to a call to <A HREF="../BV/BVOrthogonalizeColumn.html#BVOrthogonalizeColumn">BVOrthogonalizeColumn</A>() followed by a
call to <A HREF="../BV/BVScaleColumn.html#BVScaleColumn">BVScaleColumn</A>() with the reciprocal of the norm.
<P>
This function first orthogonalizes vector V[j] with respect to V[0..j-1],
where V[.] are the vectors of <A HREF="../BV/BV.html#BV">BV</A>. A byproduct of this computation is norm,
the norm of the vector after orthogonalization. Secondly, it scales the
vector with 1/norm, so that the resulting vector has unit norm.
<P>
If after orthogonalization the vector V[j] is exactly zero, it cannot be normalized
because norm=0. In that case, it could be left as zero or replaced by a random
vector that is then orthonormalized. The latter is achieved by setting the
argument replace to TRUE. The vector will be replaced by a random vector also
if lindep was set to TRUE, even if the norm is not exactly zero.
<P>
If the vector has been replaced by a random vector, the output arguments norm and
lindep will be set according to the orthogonalization of this new vector.
<P>

<P>
<H3><FONT COLOR="#883300">See Also</FONT></H3>
 <A HREF="../BV/BVOrthogonalizeColumn.html#BVOrthogonalizeColumn">BVOrthogonalizeColumn</A>(), <A HREF="../BV/BVScaleColumn.html#BVScaleColumn">BVScaleColumn</A>()
<BR><P><B></B><H3><FONT COLOR="#883300">Level</FONT></H3>advanced<BR>
<H3><FONT COLOR="#883300">Location</FONT></H3>
</B><A HREF="../../../src/sys/classes/bv/interface/bvorthog.c.html#BVOrthonormalizeColumn">src/sys/classes/bv/interface/bvorthog.c</A>
<BR><BR><A HREF="./index.html">Index of all BV 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>