Actual source code: nepresolv.c
slepc-3.18.0 2022-10-01
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: */
10: /*
11: NEP routines related to resolvent T^{-1}(z) = sum_i (z-lambda_i)^{-1} x_i y_i'
12: */
14: #include <slepc/private/nepimpl.h>
16: typedef struct {
17: NEP nep;
18: RG rg;
19: PetscScalar omega;
20: PetscScalar *nfactor; /* normalization factors y_i'*T'(lambda_i)*x_i */
21: PetscBool *nfactor_avail;
22: PetscScalar *dots; /* inner products y_i'*v */
23: PetscBool *dots_avail;
24: PetscObjectId vid;
25: PetscObjectState vstate;
26: } NEP_RESOLVENT_MATSHELL;
28: static PetscErrorCode MatMult_Resolvent(Mat M,Vec v,Vec r)
29: {
30: NEP_RESOLVENT_MATSHELL *ctx;
31: NEP nep;
32: PetscInt i,inside=1;
33: PetscScalar alpha;
34: Vec x,y,z,w;
36: MatShellGetContext(M,&ctx);
37: nep = ctx->nep;
38: w = nep->work[0];
39: z = nep->work[1];
40: if (((PetscObject)v)->id != ctx->vid || ((PetscObject)v)->state != ctx->vstate) {
41: PetscArrayzero(ctx->dots_avail,ctx->nep->nconv);
42: PetscObjectGetId((PetscObject)v,&ctx->vid);
43: PetscObjectStateGet((PetscObject)v,&ctx->vstate);
44: }
45: VecSet(r,0.0);
46: for (i=0;i<nep->nconv;i++) {
47: if (ctx->rg) RGCheckInside(ctx->rg,1,&nep->eigr[i],&nep->eigi[i],&inside);
48: if (inside>=0) {
49: BVGetColumn(nep->V,i,&x);
50: BVGetColumn(nep->W,i,&y);
51: NEPApplyJacobian(nep,nep->eigr[i],x,z,w,NULL);
52: if (!ctx->dots_avail[i]) {
53: VecDot(v,y,&ctx->dots[i]);
54: ctx->dots_avail[i] = PETSC_TRUE;
55: }
56: if (!ctx->nfactor_avail[i]) {
57: VecDot(w,y,&ctx->nfactor[i]);
58: ctx->nfactor_avail[i] = PETSC_TRUE;
59: }
60: alpha = ctx->dots[i]/(ctx->nfactor[i]*(ctx->omega-nep->eigr[i]));
61: VecAXPY(r,alpha,x);
62: BVRestoreColumn(nep->V,i,&x);
63: BVRestoreColumn(nep->W,i,&y);
64: }
65: }
66: return 0;
67: }
69: static PetscErrorCode MatDestroy_Resolvent(Mat M)
70: {
71: NEP_RESOLVENT_MATSHELL *ctx;
73: if (M) {
74: MatShellGetContext(M,&ctx);
75: PetscFree4(ctx->nfactor,ctx->nfactor_avail,ctx->dots,ctx->dots_avail);
76: PetscFree(ctx);
77: }
78: return 0;
79: }
81: /*@
82: NEPApplyResolvent - Applies the resolvent T^{-1}(z) to a given vector.
84: Collective on nep
86: Input Parameters:
87: + nep - eigensolver context obtained from NEPCreate()
88: . rg - optional region
89: . omega - value where the resolvent must be evaluated
90: - v - input vector
92: Output Parameter:
93: . r - result vector
95: Notes:
96: The resolvent T^{-1}(z) = sum_i (z-lambda_i)^{-1}*x_i*y_i' is evaluated at
97: z=omega and the matrix-vector multiplication r = T^{-1}(omega)*v is computed.
98: Vectors x_i and y_i are right and left eigenvectors, respectively, normalized
99: so that y_i'*T'(lambda_i)*x_i=1. The sum contains only eigenvectors that have
100: been previously computed with NEPSolve(), and if a region rg is given then only
101: those corresponding to eigenvalues inside the region are considered.
103: Level: intermediate
105: .seealso: NEPGetLeftEigenvector(), NEPSolve()
106: @*/
107: PetscErrorCode NEPApplyResolvent(NEP nep,RG rg,PetscScalar omega,Vec v,Vec r)
108: {
109: NEP_RESOLVENT_MATSHELL *ctx;
115: NEPCheckSolved(nep,1);
117: PetscLogEventBegin(NEP_Resolvent,nep,0,0,0);
118: if (!nep->resolvent) {
119: PetscNew(&ctx);
120: ctx->nep = nep;
121: PetscCalloc4(nep->nconv,&ctx->nfactor,nep->nconv,&ctx->nfactor_avail,nep->nconv,&ctx->dots,nep->nconv,&ctx->dots_avail);
122: MatCreateShell(PetscObjectComm((PetscObject)nep),nep->nloc,nep->nloc,nep->n,nep->n,ctx,&nep->resolvent);
123: MatShellSetOperation(nep->resolvent,MATOP_MULT,(void(*)(void))MatMult_Resolvent);
124: MatShellSetOperation(nep->resolvent,MATOP_DESTROY,(void(*)(void))MatDestroy_Resolvent);
125: } else MatShellGetContext(nep->resolvent,&ctx);
126: NEPComputeVectors(nep);
127: NEPSetWorkVecs(nep,2);
128: ctx->rg = rg;
129: ctx->omega = omega;
130: MatMult(nep->resolvent,v,r);
131: PetscLogEventEnd(NEP_Resolvent,nep,0,0,0);
132: return 0;
133: }