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
|
<!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>SlepcUpdateStrideVectors</TITLE>
</HEAD>
<BODY BGCOLOR="FFFFFF">
<H1>SlepcUpdateStrideVectors</H1>
Update a set of vectors V as V(:,s:d:e-1) = V*Q(:,s:e-1).
<H3><FONT COLOR="#883300">Synopsis</FONT></H3>
<PRE>
#include "slepcvec.h"
PetscErrorCode SlepcUpdateStrideVectors(PetscInt n_,Vec *V,PetscInt s,PetscInt d,PetscInt e,const PetscScalar *Q,PetscInt ldq_,PetscBool qtrans)
</PRE>
Not 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>n </B></TD><TD> - number of vectors in V
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>s </B></TD><TD> - first column of V to be overwritten
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>d </B></TD><TD> - stride
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>e </B></TD><TD> - first column of V not to be overwritten
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>Q </B></TD><TD> - matrix containing the coefficients of the update
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>ldq </B></TD><TD> - leading dimension of Q
</TD></TR>
<TR><TD WIDTH=40></TD><TD ALIGN=LEFT VALIGN=TOP><B>qtrans </B></TD><TD> - flag indicating if Q is to be transposed
</TD></TR></TABLE>
<P>
<H3><FONT COLOR="#883300">Input/Output parameter</FONT></H3>
<DT><B>V </B> - set of vectors
<br>
<P>
<H3><FONT COLOR="#883300">Notes</FONT></H3>
This function computes V(:,s:d:e-1) = V*Q(:,s:e-1), that is, given a set
of vectors V, columns from s to e-1 are overwritten with columns from s to
e-1 of the matrix-matrix product V*Q.
<P>
Matrix V is represented as an array of Vec, whereas Q is represented as
a column-major dense array of leading dimension ldq. Only columns s to e-1
of Q are referenced.
<P>
If qtrans=PETSC_TRUE, the operation is V*Q'.
<P>
This routine is implemented with a call to BLAS, therefore V is an array
of Vec which have the data stored contiguously in memory as a Fortran matrix.
PETSc does not create such arrays by default.
<P>
<P>
<P><B><FONT COLOR="#883300">Location: </FONT></B><A HREF="../../../src/vec/contiguous.c.html#SlepcUpdateStrideVectors">src/vec/contiguous.c</A>
<BR><A HREF="./index.html">Index of all sys 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>
|