Actual source code: nepresolv.c

slepc-3.21.0 2024-03-30
Report Typos and Errors
  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:   PetscFunctionBegin;
 37:   PetscCall(MatShellGetContext(M,&ctx));
 38:   nep = ctx->nep;
 39:   w = nep->work[0];
 40:   z = nep->work[1];
 41:   if (((PetscObject)v)->id != ctx->vid || ((PetscObject)v)->state != ctx->vstate) {
 42:     PetscCall(PetscArrayzero(ctx->dots_avail,ctx->nep->nconv));
 43:     PetscCall(PetscObjectGetId((PetscObject)v,&ctx->vid));
 44:     PetscCall(PetscObjectStateGet((PetscObject)v,&ctx->vstate));
 45:   }
 46:   PetscCall(VecSet(r,0.0));
 47:   for (i=0;i<nep->nconv;i++) {
 48:     if (ctx->rg) PetscCall(RGCheckInside(ctx->rg,1,&nep->eigr[i],&nep->eigi[i],&inside));
 49:     if (inside>=0) {
 50:       PetscCall(BVGetColumn(nep->V,i,&x));
 51:       PetscCall(BVGetColumn(nep->W,i,&y));
 52:       PetscCall(NEPApplyJacobian(nep,nep->eigr[i],x,z,w,NULL));
 53:       if (!ctx->dots_avail[i]) {
 54:         PetscCall(VecDot(v,y,&ctx->dots[i]));
 55:         ctx->dots_avail[i] = PETSC_TRUE;
 56:       }
 57:       if (!ctx->nfactor_avail[i]) {
 58:         PetscCall(VecDot(w,y,&ctx->nfactor[i]));
 59:         ctx->nfactor_avail[i] = PETSC_TRUE;
 60:       }
 61:       alpha = ctx->dots[i]/(ctx->nfactor[i]*(ctx->omega-nep->eigr[i]));
 62:       PetscCall(VecAXPY(r,alpha,x));
 63:       PetscCall(BVRestoreColumn(nep->V,i,&x));
 64:       PetscCall(BVRestoreColumn(nep->W,i,&y));
 65:     }
 66:   }
 67:   PetscFunctionReturn(PETSC_SUCCESS);
 68: }

 70: static PetscErrorCode MatDestroy_Resolvent(Mat M)
 71: {
 72:   NEP_RESOLVENT_MATSHELL *ctx;

 74:   PetscFunctionBegin;
 75:   if (M) {
 76:     PetscCall(MatShellGetContext(M,&ctx));
 77:     PetscCall(PetscFree4(ctx->nfactor,ctx->nfactor_avail,ctx->dots,ctx->dots_avail));
 78:     PetscCall(PetscFree(ctx));
 79:   }
 80:   PetscFunctionReturn(PETSC_SUCCESS);
 81: }

 83: /*@
 84:    NEPApplyResolvent - Applies the resolvent T^{-1}(z) to a given vector.

 86:    Collective

 88:    Input Parameters:
 89: +  nep   - eigensolver context obtained from NEPCreate()
 90: .  rg    - optional region
 91: .  omega - value where the resolvent must be evaluated
 92: -  v     - input vector

 94:    Output Parameter:
 95: .  r     - result vector

 97:    Notes:
 98:    The resolvent T^{-1}(z) = sum_i (z-lambda_i)^{-1}*x_i*y_i' is evaluated at
 99:    z=omega and the matrix-vector multiplication r = T^{-1}(omega)*v is computed.
100:    Vectors x_i and y_i are right and left eigenvectors, respectively, normalized
101:    so that y_i'*T'(lambda_i)*x_i=1. The sum contains only eigenvectors that have
102:    been previously computed with NEPSolve(), and if a region rg is given then only
103:    those corresponding to eigenvalues inside the region are considered.

105:    Level: intermediate

107: .seealso: NEPGetLeftEigenvector(), NEPSolve()
108: @*/
109: PetscErrorCode NEPApplyResolvent(NEP nep,RG rg,PetscScalar omega,Vec v,Vec r)
110: {
111:   NEP_RESOLVENT_MATSHELL *ctx;

113:   PetscFunctionBegin;
118:   NEPCheckSolved(nep,1);

120:   PetscCall(PetscLogEventBegin(NEP_Resolvent,nep,0,0,0));
121:   if (!nep->resolvent) {
122:     PetscCall(PetscNew(&ctx));
123:     ctx->nep = nep;
124:     PetscCall(PetscCalloc4(nep->nconv,&ctx->nfactor,nep->nconv,&ctx->nfactor_avail,nep->nconv,&ctx->dots,nep->nconv,&ctx->dots_avail));
125:     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)nep),nep->nloc,nep->nloc,nep->n,nep->n,ctx,&nep->resolvent));
126:     PetscCall(MatShellSetOperation(nep->resolvent,MATOP_MULT,(void(*)(void))MatMult_Resolvent));
127:     PetscCall(MatShellSetOperation(nep->resolvent,MATOP_DESTROY,(void(*)(void))MatDestroy_Resolvent));
128:   } else PetscCall(MatShellGetContext(nep->resolvent,&ctx));
129:   PetscCall(NEPComputeVectors(nep));
130:   PetscCall(NEPSetWorkVecs(nep,2));
131:   ctx->rg    = rg;
132:   ctx->omega = omega;
133:   PetscCall(MatMult(nep->resolvent,v,r));
134:   PetscCall(PetscLogEventEnd(NEP_Resolvent,nep,0,0,0));
135:   PetscFunctionReturn(PETSC_SUCCESS);
136: }