Actual source code: shift.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:    Shift spectral transformation, applies (A + sigma I) as operator, or
 12:    inv(B)(A + sigma B) for generalized problems
 13: */

 15: #include <slepc/private/stimpl.h>

 17: PetscErrorCode STApply_Shift(ST st,Vec x,Vec y)
 18: {

 22:   if (st->nmat>1) {
 23:     /* generalized eigenproblem: y = B^-1 (A - sB) x */
 24:     MatMult(st->T[0],x,st->work[0]);
 25:     STMatSolve(st,st->work[0],y);
 26:   } else {
 27:     /* standard eigenproblem: y = (A - sI) x */
 28:     MatMult(st->T[0],x,y);
 29:   }
 30:   return(0);
 31: }

 33: PetscErrorCode STApplyTranspose_Shift(ST st,Vec x,Vec y)
 34: {

 38:   if (st->nmat>1) {
 39:     /* generalized eigenproblem: y = (A - sB)^T B^-T  x */
 40:     STMatSolveTranspose(st,x,st->work[0]);
 41:     MatMultTranspose(st->T[0],st->work[0],y);
 42:   } else {
 43:     /* standard eigenproblem: y = (A^T - sI) x */
 44:     MatMultTranspose(st->T[0],x,y);
 45:   }
 46:   return(0);
 47: }

 49: PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 50: {
 51:   PetscInt j;

 54:   for (j=0;j<n;j++) {
 55:     eigr[j] += st->sigma;
 56:   }
 57:   return(0);
 58: }

 60: PetscErrorCode STPostSolve_Shift(ST st)
 61: {

 65:   if (st->shift_matrix == ST_MATMODE_INPLACE) {
 66:     if (st->nmat>1) {
 67:       MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
 68:     } else {
 69:       MatShift(st->A[0],st->sigma);
 70:     }
 71:     st->Astate[0] = ((PetscObject)st->A[0])->state;
 72:     st->state = ST_STATE_INITIAL;
 73:   }
 74:   return(0);
 75: }

 77: PetscErrorCode STSetUp_Shift(ST st)
 78: {
 80:   PetscInt       k,nc,nmat=PetscMax(st->nmat,2);
 81:   PetscScalar    *coeffs=NULL;

 84:   if (st->nmat>1) {
 85:     STSetWorkVecs(st,1);
 86:   }
 87:   if (nmat<3 || st->transform) {
 88:     if (nmat>2) {
 89:       nc = (nmat*(nmat+1))/2;
 90:       PetscMalloc1(nc,&coeffs);
 91:       /* Compute coeffs */
 92:       STCoeffs_Monomial(st,coeffs);
 93:     }
 94:     /* T[n] = A_n */
 95:     k = nmat-1;
 96:     PetscObjectReference((PetscObject)st->A[k]);
 97:     MatDestroy(&st->T[k]);
 98:     st->T[k] = st->A[k];
 99:     for (k=0;k<nmat-1;k++) {
100:       STMatMAXPY_Private(st,nmat>2?st->sigma:-st->sigma,0.0,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[k]);
101:     }
102:      if (nmat>2) { PetscFree(coeffs); }
103:   } else {
104:     for (k=0;k<nmat;k++) {
105:       PetscObjectReference((PetscObject)st->A[k]);
106:       MatDestroy(&st->T[k]);
107:       st->T[k] = st->A[k];
108:     }
109:   }
110:   if (nmat>=2 && st->transform) {
111:     PetscObjectReference((PetscObject)st->T[nmat-1]);
112:     MatDestroy(&st->P);
113:     st->P = st->T[nmat-1];
114:   }
115:   if (st->P) {
116:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
117:     STCheckFactorPackage(st);
118:     KSPSetOperators(st->ksp,st->P,st->P);
119:     KSPSetUp(st->ksp);
120:   }
121:   return(0);
122: }

124: PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
125: {
127:   PetscInt       k,nc,nmat=PetscMax(st->nmat,2);
128:   PetscScalar    *coeffs=NULL;

131:   if (st->transform) {
132:     if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
133:       nc = (nmat*(nmat+1))/2;
134:       PetscMalloc1(nc,&coeffs);
135:       /* Compute coeffs */
136:       STCoeffs_Monomial(st,coeffs);
137:     }
138:     for (k=0;k<nmat-1;k++) {
139:       STMatMAXPY_Private(st,nmat>2?newshift:-newshift,nmat>2?st->sigma:-st->sigma,k,coeffs?coeffs+((nmat-k)*(nmat-k-1))/2:NULL,PETSC_FALSE,&st->T[k]);
140:     }
141:     if (st->shift_matrix == ST_MATMODE_COPY && nmat>2) {
142:         PetscFree(coeffs);
143:     }
144:   }
145:   return(0);
146: }

148: SLEPC_EXTERN PetscErrorCode STCreate_Shift(ST st)
149: {
151:   st->ops->apply           = STApply_Shift;
152:   st->ops->getbilinearform = STGetBilinearForm_Default;
153:   st->ops->applytrans      = STApplyTranspose_Shift;
154:   st->ops->postsolve       = STPostSolve_Shift;
155:   st->ops->backtransform   = STBackTransform_Shift;
156:   st->ops->setup           = STSetUp_Shift;
157:   st->ops->setshift        = STSetShift_Shift;
158:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
159:   return(0);
160: }