Actual source code: nepresolv.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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>       /*I "slepcnep.h" I*/

 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: } ResolventCtx;

 28: static PetscErrorCode MatMult_Resolvent(Mat M,Vec v,Vec r)
 29: {
 31:   ResolventCtx   *ctx;
 32:   NEP            nep;
 33:   PetscInt       i,inside=1;
 34:   PetscScalar    alpha;
 35:   Vec            x,y,z,w;

 38:   MatShellGetContext(M,(void**)&ctx);
 39:   nep = ctx->nep;
 40:   w = nep->work[0];
 41:   z = nep->work[1];
 42:   if (((PetscObject)v)->id != ctx->vid || ((PetscObject)v)->state != ctx->vstate) {
 43:     PetscMemzero(ctx->dots_avail,ctx->nep->nconv*sizeof(PetscBool));
 44:     PetscObjectGetId((PetscObject)v,&ctx->vid);
 45:     PetscObjectStateGet((PetscObject)v,&ctx->vstate);
 46:   }
 47:   VecSet(r,0.0);
 48:   for (i=0;i<nep->nconv;i++) {
 49:     if (ctx->rg) {
 50:       RGCheckInside(ctx->rg,1,&nep->eigr[i],&nep->eigi[i],&inside);
 51:     }
 52:     if (inside>=0) {
 53:       BVGetColumn(nep->V,i,&x);
 54:       BVGetColumn(nep->W,i,&y);
 55:       NEPApplyJacobian(nep,nep->eigr[i],x,z,w,NULL);
 56:       if (!ctx->dots_avail[i]) {
 57:         VecDot(v,y,&ctx->dots[i]);
 58:         ctx->dots_avail[i] = PETSC_TRUE;
 59:       }
 60:       if (!ctx->nfactor_avail[i]) {
 61:         VecDot(w,y,&ctx->nfactor[i]);
 62:         ctx->nfactor_avail[i] = PETSC_TRUE;
 63:       }
 64:       alpha = ctx->dots[i]/(ctx->nfactor[i]*(ctx->omega-nep->eigr[i]));
 65:       VecAXPY(r,alpha,x);
 66:       BVRestoreColumn(nep->V,i,&x);
 67:       BVRestoreColumn(nep->W,i,&y);
 68:     }
 69:   }
 70:   return(0);
 71: }

 73: static PetscErrorCode MatDestroy_Resolvent(Mat M)
 74: {
 76:   ResolventCtx   *ctx;

 79:   if (M) {
 80:     MatShellGetContext(M,(void**)&ctx);
 81:     PetscFree4(ctx->nfactor,ctx->nfactor_avail,ctx->dots,ctx->dots_avail);
 82:     PetscFree(ctx);
 83:   }
 84:   return(0);
 85: }

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

 90:    Collective on NEP

 92:    Input Parameters:
 93: +  nep   - eigensolver context obtained from NEPCreate()
 94: .  rg    - optional region
 95: .  omega - value where the resolvent must be evaluated
 96: -  v     - input vector

 98:    Output Parameter:
 99: .  r     - result vector

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

109:    Level: intermediate

111: .seealso: NEPGetLeftEigenvector(), NEPSolve()
112: @*/
113: PetscErrorCode NEPApplyResolvent(NEP nep,RG rg,PetscScalar omega,Vec v,Vec r)
114: {
116:   ResolventCtx   *ctx;

123:   NEPCheckSolved(nep,1);

125:   PetscLogEventBegin(NEP_Resolvent,nep,0,0,0);
126:   if (!nep->resolvent) {
127:     PetscNew(&ctx);
128:     ctx->nep = nep;
129:     PetscCalloc4(nep->nconv,&ctx->nfactor,nep->nconv,&ctx->nfactor_avail,nep->nconv,&ctx->dots,nep->nconv,&ctx->dots_avail);
130:     MatCreateShell(PetscObjectComm((PetscObject)nep),nep->nloc,nep->nloc,nep->n,nep->n,ctx,&nep->resolvent);
131:     MatShellSetOperation(nep->resolvent,MATOP_MULT,(void(*)(void))MatMult_Resolvent);
132:     MatShellSetOperation(nep->resolvent,MATOP_DESTROY,(void(*)(void))MatDestroy_Resolvent);
133:   } else {
134:     MatShellGetContext(nep->resolvent,(void**)&ctx);
135:   }
136:   NEPComputeVectors(nep);
137:   NEPSetWorkVecs(nep,2);
138:   ctx->rg    = rg;
139:   ctx->omega = omega;
140:   MatMult(nep->resolvent,v,r);
141:   PetscLogEventEnd(NEP_Resolvent,nep,0,0,0);
142:   return(0);
143: }