slepc-main 2024-12-17
BVMultInPlace
Update a set of vectors as V(:,s:e-1) = V*Q(:,s:e-1).
Synopsis
#include "slepcbv.h"
PetscErrorCode BVMultInPlace(BV V,Mat Q,PetscInt s,PetscInt e)
Logically Collective
Input Parameters
| Q | - a sequential dense matrix
|
| s | - first column of V to be overwritten
|
| e | - first column of V not to be overwritten
|
Input/Output Parameter
Notes
The matrix Q must be a sequential dense Mat, with all entries equal on
all processes (otherwise each process will compute a different update).
This function computes V(:,s: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. Only columns s to e-1 of Q are
referenced.
See Also
BVMult(), BVMultVec(), BVMultInPlaceHermitianTranspose(), BVSetActiveColumns()
Level
intermediate
Location
src/sys/classes/bv/interface/bvops.c
Index of all BV routines
Table of Contents for all manual pages
Index of all manual pages