Actual source code: slepccontour.h
slepc-3.18.3 2023-03-24
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: #if !defined(SLEPCCONTOUR_H)
12: #define SLEPCCONTOUR_H
14: #include <slepc/private/slepcimpl.h>
15: #include <petscksp.h>
17: /*
18: CISS_BlockHankel - Builds a block Hankel matrix from the contents of Mu.
19: */
20: static inline PetscErrorCode CISS_BlockHankel(PetscScalar *Mu,PetscInt s,PetscInt L,PetscInt M,PetscScalar *H)
21: {
22: PetscInt i,j,k;
24: for (k=0;k<L*M;k++)
25: for (j=0;j<M;j++)
26: for (i=0;i<L;i++)
27: H[j*L+i+k*L*M] = Mu[i+k*L+(j+s)*L*L];
28: return 0;
29: }
31: SLEPC_EXTERN PetscErrorCode SlepcCISS_isGhost(Mat,PetscInt,PetscReal*,PetscReal,PetscBool*);
32: SLEPC_EXTERN PetscErrorCode SlepcCISS_BH_SVD(PetscScalar*,PetscInt,PetscReal,PetscReal*,PetscInt*);
34: /* Data structures and functions for contour integral methods (used in several classes) */
35: struct _n_SlepcContourData {
36: PetscObject parent; /* parent object */
37: PetscSubcomm subcomm; /* subcommunicator for top level parallelization */
38: PetscInt npoints; /* number of integration points assigned to the local subcomm */
39: KSP *ksp; /* ksp array for storing factorizations at integration points */
40: Mat *pA; /* redundant copies of the matrices in the local subcomm */
41: Mat *pP; /* redundant copies of the matrices (preconditioner) */
42: PetscInt nmat; /* number of matrices in pA */
43: Vec xsub; /* aux vector with parallel layout as redundant Mat */
44: Vec xdup; /* aux vector with parallel layout as original Mat (with contiguous order) */
45: VecScatter scatterin; /* to scatter from regular vector to xdup */
46: };
47: typedef struct _n_SlepcContourData* SlepcContourData;
49: SLEPC_EXTERN PetscErrorCode SlepcContourDataCreate(PetscInt,PetscInt,PetscObject,SlepcContourData*);
50: SLEPC_EXTERN PetscErrorCode SlepcContourDataReset(SlepcContourData);
51: SLEPC_EXTERN PetscErrorCode SlepcContourDataDestroy(SlepcContourData*);
53: SLEPC_EXTERN PetscErrorCode SlepcContourRedundantMat(SlepcContourData,PetscInt,Mat*,Mat*);
54: SLEPC_EXTERN PetscErrorCode SlepcContourScatterCreate(SlepcContourData,Vec);
56: #endif