slepc-3.21.0 2024-03-30
Report Typos and Errors

BVSumQuadrature

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

Synopsis

#include "slepcbv.h" 
PetscErrorCode BVSumQuadrature(BV S,BV Y,PetscInt M,PetscInt L,PetscInt L_max,PetscScalar *w,PetscScalar *zn,VecScatter scat,PetscSubcomm subcomm,PetscInt npoints,PetscBool useconj)
Collective

Input Parameters

Y  - input basis vectors
M  - number of moments
L  - block size
L_max  - maximum block size
w  - quadrature weights
zn  - normalized quadrature points
scat  - (optional) VecScatter object to communicate between subcommunicators
subcomm  - subcommunicator layout
npoints  - number of points to process by the subcommunicator
useconj  - whether conjugate points can be used or not

Output Parameter

S  - output basis vectors

Notes

This is a generalization of BVMult(). The resulting matrix S consists of M panels of L columns, and the following formula is computed for each panel S_k = sum_j w_j*zn_j^k*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 scattered to the whole communicator in S using the VecScatter scat. 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

BVMult(), BVScatter(), BVDotQuadrature(), 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