slepc-3.20.2 2024-03-15
Report Typos and Errors

BVDotQuadrature

Computes the projection terms required in the quadrature rule to approximate the contour integral.

Synopsis

#include "slepcbv.h" 
PetscErrorCode BVDotQuadrature(BV Y,BV V,PetscScalar *Mu,PetscInt M,PetscInt L,PetscInt L_max,PetscScalar *w,PetscScalar *zn,PetscSubcomm subcomm,PetscInt npoints,PetscBool useconj)
Collective

Input Parameters

Y  - first basis vectors
V  - second basis vectors
M  - number of moments
L  - block size
L_max  - maximum block size
w  - quadrature weights
zn  - normalized quadrature points
subcomm  - subcommunicator layout
npoints  - number of points to process by the subcommunicator
useconj  - whether conjugate points can be used or not

Output Parameter

Mu  - computed result

Notes

This is a generalization of BVDot(). The resulting matrix Mu consists of M blocks of size LxL (placed horizontally), each of them computed as Mu_k = sum_j w_j*zn_j^k*V'*Y_j, where Y_j is the j-th panel of Y containing the result of solving T(z_j)^{-1}*X for each integration point j. L_max is the width of the panels in Y.

When using subcommunicators, Y is stored in the subcommunicators for a subset of integration points. In that case, the computation is done in the subcomm and then the final result is combined via reduction. The value npoints is the number of points to be processed in this subcomm and the flag useconj indicates whether symmetric points can be reused.

See Also

BVDot(), BVScatter(), BVSumQuadrature(), RGComputeQuadrature(), RGCanUseConjugates()

Level

developer

Location

src/sys/classes/bv/interface/bvcontour.c

Index of all BV routines
Table of Contents for all manual pages
Index of all manual pages