#include "slepcsys.h" PetscErrorCode SlepcUpdateStrideVectors(PetscInt n_,Vec *V,PetscInt s,PetscInt d,PetscInt e,const PetscScalar *Q,PetscInt ldq_,PetscTruth qtrans)Collective on Vec
n | - number of vectors in V | |
s | - first column of V to be overwritten | |
d | - stride | |
e | - first column of V not to be overwritten | |
Q | - matrix containing the coefficients of the update | |
ldq | - leading dimension of Q | |
qtrans | - flag indicating if Q is to be transposed |
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.
If qtrans=PETSC_TRUE, the operation is V*Q'.
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.
Location: src/sys/slepcutil.c
Index of all sys routines
Table of Contents for all manual pages
Index of all manual pages